Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
NonLinElasAssFE_def.hpp
1#ifndef NonLinElasAssFE_def_hpp
2#define NonLinElasAssFE_def_hpp
3#include "NonLinElasAssFE_decl.hpp"
12
13
14namespace FEDD {
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):
17NonLinearProblem<SC,LO,GO,NO>( parameterList, domain->getComm() ),
18u_rep_()
19{
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();
24
25 C_ = this->parameterList_->sublist("Parameter").get("C",1.);
26
27 double density = this->parameterList_->sublist("Parameter").get("Density",1.);
28
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);
32
33 timeSteppingTool_ = Teuchos::rcp(new TimeSteppingTools(sublist(this->parameterList_,"Timestepping Parameter") , this->comm_));
34
35}
36
37template<class SC,class LO,class GO,class NO>
38NonLinElasAssFE<SC,LO,GO,NO>::~NonLinElasAssFE(){
39
40}
41
42template<class SC,class LO,class GO,class NO>
43void NonLinElasAssFE<SC,LO,GO,NO>::info(){
44 this->infoProblem();
45 this->infoNonlinProblem();
46}
47
48template<class SC,class LO,class GO,class NO>
49void NonLinElasAssFE<SC,LO,GO,NO>::assemble(std::string type) const{
50
51 if (type == ""){
52 if (this->verbose_)
53 std::cout << "-- Assembly nonlinear elasticity ... " << std::flush;
54
55 MatrixPtr_Type A(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
56
57 this->feFactory_->assemblyEmptyMatrix(A);
58
59 this->system_.reset(new BlockMatrix_Type(1));
60 this->system_->addBlock( A, 0, 0 );
61
62 double density = this->parameterList_->sublist("Parameter").get("Density",1000.);
63 std::string sourceType = this->parameterList_->sublist("Parameter").get("Source Type","volume");
64
65 this->assembleSourceTerm( 0. );
66 if(sourceType == "volume")
67 this->sourceTerm_->scale(density);
68
69 this->addToRhs( this->sourceTerm_ );
70
71 this->setBoundariesRHS();
72
73 this->solution_->putScalar(0.);
74
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);
78
79 if (this->verbose_)
80 std::cout << "done -- " << std::endl;
81
82 this->reAssemble("Newton-Residual");
83 }
84 else if(type == "UpdateTime")
85 {
86 if(this->verbose_)
87 std::cout << "-- Reassembly (UpdateTime)" << '\n';
88
89 updateTime();
90 return;
91 }
92 else
93 this->reAssemble(type);
94}
95
96// Damit die richtige timeSteppingTool_->currentTime() genommen wird.
97template<class SC,class LO,class GO,class NO>
98void NonLinElasAssFE<SC,LO,GO,NO>::updateTime() const
99{
100 timeSteppingTool_->t_ = timeSteppingTool_->t_ + timeSteppingTool_->dt_prev_;
101}
102
103template<class SC,class LO,class GO,class NO>
104void NonLinElasAssFE<SC,LO,GO,NO>::reAssemble(std::string type) const {
105
106 std::string material_model = this->parameterList_->sublist("Parameter").get("Structure Model","Neo-Hooke");
107
108 if (this->verbose_)
109 std::cout << "-- Reassembly nonlinear elasticity with structure material model " << material_model <<" ("<<type <<") ... " << std::flush;
110
111 if (type=="Newton-Residual") {
112 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
113 u_rep_->importFromVector(u, true);
114
115 MatrixPtr_Type W = Teuchos::rcp(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
116
117
118 this->system_->addBlock( W, 0, 0 );
119
120 MultiVectorPtr_Type fRep = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated(), 1 ) );
121 fRep->putScalar(0.);
122 this->residualVec_->addBlock( fRep, 0 );
123
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/*call fillComplete*/);
126 else
127 this->feFactory_->assemblyNonLinearElasticity(this->dim_, this->getDomain(0)->getFEType(),2, this->dim_, u_rep_, this->system_, this->residualVec_, this->parameterList_,true);
128 //fRep->scale(-1.);
129
130 MultiVectorPtr_Type fUnique = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldUnique(), 1 ) );
131 fUnique->putScalar(0.);
132 fUnique->exportFromVector( fRep, true, "Add" );
133
134 this->residualVec_->addBlock( fUnique, 0 );
135
136 if(loadStepping_)
137 assembleSourceTermLoadstepping();
138
139
140 }
141 else if(type=="Newton"){ //we already assemble the new tangent when we calculate the stresses above
142
143 }
144 if (this->verbose_)
145 std::cout << "done -- " << std::endl;
146}
147
148template<class SC,class LO,class GO,class NO>
149void NonLinElasAssFE<SC,LO,GO,NO>::reAssemble( MatrixPtr_Type& massmatrix, std::string type ) const
150{
151
152}
153
154template<class SC,class LO,class GO,class NO>
155void NonLinElasAssFE<SC,LO,GO,NO>::reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions){
156
157
158 TEUCHOS_TEST_FOR_EXCEPTION( true, std::logic_error, "Only Newton/NOX implemented for nonlinear material models!");
159
160}
161
162template<class SC,class LO,class GO,class NO>
163void NonLinElasAssFE<SC,LO,GO,NO>::assembleSourceTermLoadstepping(double time) const
164{
165 double dt = timeSteppingTool_->get_dt();
166 double beta = timeSteppingTool_->get_beta();
167 double gamma = timeSteppingTool_->get_gamma();
168
169
170 if( loadStepping_ == true){
171 // The if condition resets the rhs. If we skip it when we attemp loadstepping, the rhs will be updated continously and wrongly increase with each timestep
172 this->getRhs()->scale(0.0);
173 }
174
175
176 if (this->hasSourceTerm())
177 {
178 if(externalForce_){
179
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_;
183 // how can we use different parameters for different blocks here?
184 funcParameter[1] =this->parameterList_->sublist("Parameter").get("Volume force",0.00211);
185
186 funcParameter[3] =this->parameterList_->sublist("Parameter").get("Final time force",1.0);
187 funcParameter[4] =this->parameterList_->sublist("Parameter").get("dt",0.1);
188
189
190 if(nonlinearExternalForce_){
191 MatrixPtr_Type A( new Matrix_Type (this->system_->getBlock(0,0)));
192 //A->print();
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_);
197
198 A->addMatrix(1.,AKext,0.);
199 Kext->addMatrix(1.,AKext,1.);
200
201 AKext->fillComplete(this->getDomain(0)->getMapVecFieldUnique(),this->getDomain(0)->getMapVecFieldUnique());
202
203 this->system_->addBlock(AKext,0,0);
204
205 }
206 else
207 this->feFactory_->assemblySurfaceIntegralExternal(this->dim_, this->getDomain(0)->getFEType(),FERhs, u_rep_, funcParameter, this->rhsFuncVec_[0],this->parameterList_);
208
209
210 this->sourceTerm_->getBlockNonConst(0)->exportFromVector( FERhs, false, "Add" );
211 }
212 else
213 this->assembleSourceTerm( timeSteppingTool_->t_ );
214 //this->problemTimeStructure_->getSourceTerm()->scale(density);
215 // Fuege die rechte Seite der DGL (f bzw. f_{n+1}) der rechten Seite hinzu (skaliert mit coeffSourceTerm)
216 // Die Skalierung mit der Dichte erfolgt schon in der Assemblierungsfunktion!
217
218 // addSourceTermToRHS() aus DAESolverInTime
219 double coeffSourceTermStructure = 1.0;
220 BlockMultiVectorPtr_Type tmpPtr= this->sourceTerm_;
221
222 this->getRhs()->update(coeffSourceTermStructure, *tmpPtr, 1.);
223 this->rhs_->addBlock( this->getRhs()->getBlockNonConst(0), 0 );
224
225 }
226}
227
228template<class SC,class LO,class GO,class NO>
229void NonLinElasAssFE<SC,LO,GO,NO>::calculateNonLinResidualVec(std::string type, double time) const{
230
231
232 this->reAssemble("Newton-Residual");
233
234 if(this->parameterList_->sublist("Parameter").get("SCI",false) == true || this->parameterList_->sublist("Parameter").get("FSCI",false) == true ){
235 //this->feFactory_->assemblyNonLinearElasticity(this->dim_, this->getDomain(0)->getFEType(),2, this->dim_, u_rep_, this->system_, this->residualVec_, this->parameterList_,this->getDomain(0), this->eModVec_,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/*call fillComplete*/);
237 }
238 else{
239 if (!type.compare("standard")){
240 this->residualVec_->update(-1.,*this->rhs_,1.);
241 //if ( !this->sourceTerm_.is_null() )
242 // this->residualVec_->update(-1.,*this->sourceTerm_,1.);
243 }
244 else if(!type.compare("reverse")){
245 this->residualVec_->update(1.,*this->rhs_,-1.); // this = -1*this + 1*rhs
246 //if ( !this->sourceTerm_.is_null() )
247 // this->residualVec_->update(1.,*this->sourceTerm_,1.);
248 }
249 else{
250 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Unknown type for residual computation.");
251 }
252
253 // this might be set again by the TimeProblem after adding of M*u
254 this->bcFactory_->setBCMinusVector( this->residualVec_, this->solution_, time );
255 }
256
257}
258}
259#endif
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