1#ifndef NonLinElasAssFE_def_hpp
2#define NonLinElasAssFE_def_hpp
3#include "NonLinElasAssFE_decl.hpp"
15template<
class SC,
class LO,
class GO,
class NO>
16NonLinElasAssFE<SC,LO,GO,NO>::NonLinElasAssFE(
const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList):
20 this->nonLinearTolerance_ = this->parameterList_->sublist(
"Parameter").get(
"relNonLinTol",1.0e-6);
21 this->initNOXParameters();
22 this->addVariable( domain , FEType ,
"u" , domain->getDimension());
23 this->dim_ = this->getDomain(0)->getDimension();
25 C_ = this->parameterList_->sublist(
"Parameter").get(
"C",1.);
27 double density = this->parameterList_->sublist(
"Parameter").get(
"Density",1.);
29 loadStepping_ = !(parameterList->sublist(
"Timestepping Parameter").get(
"Class",
"Singlestep")).compare(
"Loadstepping");
30 externalForce_ = parameterList->sublist(
"Parameter").get(
"External Force",
false);
31 nonlinearExternalForce_ = parameterList->sublist(
"Parameter").get(
"Nonlinear External Force",
false);
33 timeSteppingTool_ = Teuchos::rcp(
new TimeSteppingTools(sublist(this->parameterList_,
"Timestepping Parameter") , this->comm_));
37template<
class SC,
class LO,
class GO,
class NO>
38NonLinElasAssFE<SC,LO,GO,NO>::~NonLinElasAssFE(){
42template<
class SC,
class LO,
class GO,
class NO>
43void NonLinElasAssFE<SC,LO,GO,NO>::info(){
45 this->infoNonlinProblem();
48template<
class SC,
class LO,
class GO,
class NO>
53 std::cout <<
"-- Assembly nonlinear elasticity ... " << std::flush;
55 MatrixPtr_Type A(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
57 this->feFactory_->assemblyEmptyMatrix(A);
59 this->system_.reset(
new BlockMatrix_Type(1));
60 this->system_->addBlock( A, 0, 0 );
62 double density = this->parameterList_->sublist(
"Parameter").get(
"Density",1000.);
63 std::string sourceType = this->parameterList_->sublist(
"Parameter").get(
"Source Type",
"volume");
65 this->assembleSourceTerm( 0. );
66 if(sourceType ==
"volume")
71 this->setBoundariesRHS();
73 this->solution_->putScalar(0.);
75 u_rep_ = Teuchos::rcp(
new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated() ));
76 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
77 u_rep_->importFromVector(u,
true);
80 std::cout <<
"done -- " << std::endl;
82 this->reAssemble(
"Newton-Residual");
84 else if(type ==
"UpdateTime")
87 std::cout <<
"-- Reassembly (UpdateTime)" <<
'\n';
93 this->reAssemble(type);
97template<
class SC,
class LO,
class GO,
class NO>
98void NonLinElasAssFE<SC,LO,GO,NO>::updateTime()
const
100 timeSteppingTool_->t_ = timeSteppingTool_->t_ + timeSteppingTool_->dt_prev_;
103template<
class SC,
class LO,
class GO,
class NO>
104void NonLinElasAssFE<SC,LO,GO,NO>::reAssemble(std::string type)
const {
106 std::string material_model = this->parameterList_->sublist(
"Parameter").get(
"Structure Model",
"Neo-Hooke");
109 std::cout <<
"-- Reassembly nonlinear elasticity with structure material model " << material_model <<
" ("<<type <<
") ... " << std::flush;
111 if (type==
"Newton-Residual") {
112 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
113 u_rep_->importFromVector(u,
true);
115 MatrixPtr_Type W = Teuchos::rcp(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
118 this->system_->addBlock( W, 0, 0 );
120 MultiVectorPtr_Type fRep = Teuchos::rcp(
new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated(), 1 ) );
122 this->residualVec_->addBlock( fRep, 0 );
124 if(this->parameterList_->sublist(
"Parameter").get(
"SCI",
false) ==
true || this->parameterList_->sublist(
"Parameter").get(
"FSCI",
false) ==
true )
125 this->feFactory_->assemblyAceDeformDiffuBlock(this->dim_, this->getDomain(0)->getFEType(), this->getDomain(0)->getFEType(), 2, 1,this->dim_,concentration_,u_rep_,this->system_,0,0,this->residualVec_,0, this->parameterList_,
"Jacobian",
true);
127 this->feFactory_->assemblyNonLinearElasticity(this->dim_, this->getDomain(0)->getFEType(),2, this->dim_, u_rep_, this->system_, this->residualVec_, this->parameterList_,
true);
130 MultiVectorPtr_Type fUnique = Teuchos::rcp(
new MultiVector_Type( this->getDomain(0)->getMapVecFieldUnique(), 1 ) );
131 fUnique->putScalar(0.);
132 fUnique->exportFromVector( fRep,
true,
"Add" );
134 this->residualVec_->addBlock( fUnique, 0 );
137 assembleSourceTermLoadstepping();
141 else if(type==
"Newton"){
145 std::cout <<
"done -- " << std::endl;
148template<
class SC,
class LO,
class GO,
class NO>
149void NonLinElasAssFE<SC,LO,GO,NO>::reAssemble( MatrixPtr_Type& massmatrix, std::string type )
const
154template<
class SC,
class LO,
class GO,
class NO>
155void NonLinElasAssFE<SC,LO,GO,NO>::reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions){
158 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Only Newton/NOX implemented for nonlinear material models!");
162template<
class SC,
class LO,
class GO,
class NO>
163void NonLinElasAssFE<SC,LO,GO,NO>::assembleSourceTermLoadstepping(
double time)
const
165 double dt = timeSteppingTool_->get_dt();
166 double beta = timeSteppingTool_->get_beta();
167 double gamma = timeSteppingTool_->get_gamma();
170 if( loadStepping_ ==
true){
172 this->getRhs()->scale(0.0);
176 if (this->hasSourceTerm())
180 MultiVectorPtr_Type FERhs = Teuchos::rcp(
new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated() ));
181 vec_dbl_Type funcParameter(4,0.);
182 funcParameter[0] = timeSteppingTool_->t_;
184 funcParameter[1] =this->parameterList_->sublist(
"Parameter").get(
"Volume force",0.00211);
186 funcParameter[3] =this->parameterList_->sublist(
"Parameter").get(
"Final time force",1.0);
187 funcParameter[4] =this->parameterList_->sublist(
"Parameter").get(
"dt",0.1);
190 if(nonlinearExternalForce_){
191 MatrixPtr_Type A(
new Matrix_Type (this->system_->getBlock(0,0)));
193 MatrixPtr_Type AKext(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
194 MatrixPtr_Type Kext(
new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow()*2 ) );
195 MultiVectorPtr_Type Kext_vec;
196 this->feFactory_->assemblyNonlinearSurfaceIntegralExternal(this->dim_, this->getDomain(0)->getFEType(),FERhs, u_rep_,Kext, funcParameter, this->rhsFuncVec_[0],this->parameterList_);
198 A->addMatrix(1.,AKext,0.);
199 Kext->addMatrix(1.,AKext,1.);
201 AKext->fillComplete(this->getDomain(0)->getMapVecFieldUnique(),this->getDomain(0)->getMapVecFieldUnique());
203 this->system_->addBlock(AKext,0,0);
207 this->feFactory_->assemblySurfaceIntegralExternal(this->dim_, this->getDomain(0)->getFEType(),FERhs, u_rep_, funcParameter, this->rhsFuncVec_[0],this->parameterList_);
210 this->sourceTerm_->getBlockNonConst(0)->exportFromVector( FERhs,
false,
"Add" );
213 this->assembleSourceTerm( timeSteppingTool_->t_ );
219 double coeffSourceTermStructure = 1.0;
220 BlockMultiVectorPtr_Type tmpPtr= this->sourceTerm_;
222 this->getRhs()->update(coeffSourceTermStructure, *tmpPtr, 1.);
223 this->rhs_->addBlock( this->getRhs()->getBlockNonConst(0), 0 );
228template<
class SC,
class LO,
class GO,
class NO>
232 this->reAssemble(
"Newton-Residual");
234 if(this->parameterList_->sublist(
"Parameter").get(
"SCI",
false) ==
true || this->parameterList_->sublist(
"Parameter").get(
"FSCI",
false) ==
true ){
236 this->feFactory_->assemblyAceDeformDiffuBlock(this->dim_, this->getDomain(0)->getFEType(), this->getDomain(0)->getFEType(), 2, 1,this->dim_,concentration_,u_rep_,this->system_,0,0,this->residualVec_,0, this->parameterList_,
"Rhs",
true);
239 if (!type.compare(
"standard")){
240 this->residualVec_->update(-1.,*this->rhs_,1.);
244 else if(!type.compare(
"reverse")){
245 this->residualVec_->update(1.,*this->rhs_,-1.);
250 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::runtime_error,
"Unknown type for residual computation.");
254 this->bcFactory_->setBCMinusVector( this->residualVec_, this->solution_, time );
virtual void assemble(std::string type="") const
assemble of type exectuted by the derived specific non-linear problem classes
Definition NonLinElasAssFE_def.hpp:49
virtual void calculateNonLinResidualVec(std::string type, double time=0.) const
Virtual function which is implemented in the specific non-linear problem classes to calculate the non...
Definition NonLinElasAssFE_def.hpp:229
Definition NonLinearProblem_decl.hpp:28
BlockMultiVectorPtr_Type sourceTerm_
Definition Problem_decl.hpp:251
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36