Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
NonLinElasticity_def.hpp
1#ifndef NonLinElasticity_def_hpp
2#define NonLinElasticity_def_hpp
3
12
13#include "feddlib/core/FE/Domain.hpp"
14#include "feddlib//core/FE/FE.hpp"
15#include "feddlib/core/General/BCBuilder.hpp"
16
17
18namespace FEDD {
19template<class SC,class LO,class GO,class NO>
20NonLinElasticity<SC,LO,GO,NO>::NonLinElasticity(const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList):
21NonLinearProblem<SC,LO,GO,NO>( parameterList, domain->getComm() ),
22u_rep_()
23{
24 this->nonLinearTolerance_ = this->parameterList_->sublist("Parameter").get("relNonLinTol",1.0e-6);
25 this->initNOXParameters();
26 this->addVariable( domain , FEType , "u" , domain->getDimension());
27 this->dim_ = this->getDomain(0)->getDimension();
28
29 C_ = this->parameterList_->sublist("Parameter").get("C",1.);
30
31 double density = this->parameterList_->sublist("Parameter").get("Density",1.);
32// poissonRatio_ = this->parameterList_->sublist("Parameter").get("Poisson Ratio",0.4);
33// mue_ = this->parameterList_->sublist("Parameter").get("Mu",2.0e+6);
34// // Berechne daraus nun E (Youngsches Modul) und die erste Lamé-Konstante \lambda
35// E_ = mue_*2.*(1. + poissonRatio_);
36// lambda_ = (poissonRatio_*E_)/((1 + poissonRatio_)*(1 - 2*poissonRatio_));
37 loadStepping_ = !(parameterList->sublist("Timestepping Parameter").get("Class","Singlestep")).compare("Loadstepping");
38 externalForce_ = parameterList->sublist("Parameter").get("External Force",false);
39 nonlinearExternalForce_ = parameterList->sublist("Parameter").get("Nonlinear External Force",false);
40
41 timeSteppingTool_ = Teuchos::rcp(new TimeSteppingTools(sublist(this->parameterList_,"Timestepping Parameter") , this->comm_));
42
43}
44
45template<class SC,class LO,class GO,class NO>
46NonLinElasticity<SC,LO,GO,NO>::~NonLinElasticity(){
47
48}
49
50template<class SC,class LO,class GO,class NO>
51void NonLinElasticity<SC,LO,GO,NO>::info(){
52 this->infoProblem();
54}
55
56template<class SC,class LO,class GO,class NO>
57void NonLinElasticity<SC,LO,GO,NO>::assemble(std::string type) const{
58
59 if (type == ""){
60 if (this->verbose_)
61 std::cout << "-- Assembly nonlinear elasticity ... " << std::flush;
62
63 MatrixPtr_Type A(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
64
65 this->feFactory_->assemblyEmptyMatrix(A);
66
67 this->system_.reset(new BlockMatrix_Type(1));
68 this->system_->addBlock( A, 0, 0 );
69
70 double density = this->parameterList_->sublist("Parameter").get("Density",1000.);
71 std::string sourceType = this->parameterList_->sublist("Parameter").get("Source Type","volume");
72
73 this->assembleSourceTerm( 0. );
74 if(sourceType == "volume")
75 this->sourceTerm_->scale(density);
76
77 this->addToRhs( this->sourceTerm_ );
78
79 this->setBoundariesRHS();
80
81 this->solution_->putScalar(0.);
82
83 u_rep_ = Teuchos::rcp(new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated() ));
84 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
85 u_rep_->importFromVector(u, true);
86
87 if (this->verbose_)
88 std::cout << "done -- " << std::endl;
89
90 this->reAssemble("Newton-Residual");
91 }
92 else if(type == "UpdateTime")
93 {
94 if(this->verbose_)
95 std::cout << "-- Reassembly (UpdateTime)" << '\n';
96
97 updateTime();
98 return;
99 }
100 else
101 this->reAssemble(type);
102}
103
104// Damit die richtige timeSteppingTool_->currentTime() genommen wird.
105template<class SC,class LO,class GO,class NO>
106void NonLinElasticity<SC,LO,GO,NO>::updateTime() const
107{
108 timeSteppingTool_->t_ = timeSteppingTool_->t_ + timeSteppingTool_->dt_prev_;
109}
110
111template<class SC,class LO,class GO,class NO>
112void NonLinElasticity<SC,LO,GO,NO>::reAssemble(std::string type) const {
113 std::string material_model = this->parameterList_->sublist("Parameter").get("Material model","Neo-Hooke");
114
115 if (this->verbose_)
116 std::cout << "-- Reassembly nonlinear elasticity with material model " << material_model <<" ("<<type <<") ... " << std::flush;
117
118 if (type=="Newton-Residual") {
119 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
120 u_rep_->importFromVector(u, true);
121
122 MultiVectorPtr_Type f = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated(), 1 ) );
123
124 MatrixPtr_Type W = Teuchos::rcp(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
125
126 #ifdef FEDD_HAVE_ACEGENINTERFACE
127 bool useInterface = this->parameterList_->sublist("Parameter").get("Use AceGen Interface", true);
128 if(this->getFEType(0) =="P2" && useInterface && this->dim_==3 && material_model == "Neo-Hooke"){
129 this->system_->addBlock( W, 0, 0 );
130
131 f->putScalar(0.);
132 this->residualVec_->addBlock( f, 0 );
133
134 if(this->parameterList_->sublist("Parameter").get("SCI",false) == true || this->parameterList_->sublist("Parameter").get("FSCI",false) == true )
135 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*/);
136 else
137 this->feFactory_->assemblyNonLinearElasticity(this->dim_, this->getDomain(0)->getFEType(),2, this->dim_, u_rep_, this->system_, this->residualVec_, this->parameterList_,true);
138
139 }
140 else
141 this->feFactory_->assemblyElasticityJacobianAndStressAceFEM(this->dim_, this->getDomain(0)->getFEType(), W, f, u_rep_, this->parameterList_, C_);
142 #else
143 this->feFactory_->assemblyElasticityJacobianAndStressAceFEM(this->dim_, this->getDomain(0)->getFEType(), W, f, u_rep_, this->parameterList_, C_);
144 #endif
145
146 MultiVectorPtr_Type fUnique = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldUnique(), 1 ) );
147 fUnique->putScalar(0.);
148 fUnique->exportFromVector( f, true, "Add" );
149
150 this->residualVec_->addBlock( fUnique, 0 );
151 this->system_->addBlock( W, 0, 0 );
152
153 if(loadStepping_)
154 assembleSourceTermLoadstepping();
155
156 }
157 else if(type=="Newton"){ //we already assemble the new tangent when we calculate the stresses above
158
159 }
160 if (this->verbose_)
161 std::cout << "done -- " << std::endl;
162}
163
164template<class SC,class LO,class GO,class NO>
165void NonLinElasticity<SC,LO,GO,NO>::reAssemble( MatrixPtr_Type& massmatrix, std::string type ) const
166{
167
168}
169
170template<class SC,class LO,class GO,class NO>
171void NonLinElasticity<SC,LO,GO,NO>::reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions){
172
173
174 TEUCHOS_TEST_FOR_EXCEPTION( true, std::logic_error, "Only Newton/NOX implemented for nonlinear material models!");
175
176}
177
178
179template<class SC,class LO,class GO,class NO>
180void NonLinElasticity<SC,LO,GO,NO>::assembleSourceTermLoadstepping(double time) const
181{
182 double dt = timeSteppingTool_->get_dt();
183 double beta = timeSteppingTool_->get_beta();
184 double gamma = timeSteppingTool_->get_gamma();
185
186
187 if( loadStepping_ == true){
188 // 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
189 this->getRhs()->scale(0.0);
190 }
191
192
193 if (this->hasSourceTerm())
194 {
195 if(externalForce_){
196
197 MultiVectorPtr_Type FERhs = Teuchos::rcp(new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated() ));
198 vec_dbl_Type funcParameter(4,0.);
199 funcParameter[0] = timeSteppingTool_->t_;
200 // how can we use different parameters for different blocks here?
201 funcParameter[1] =this->parameterList_->sublist("Parameter").get("Volume force",0.00211);
202
203 funcParameter[3] =this->parameterList_->sublist("Parameter").get("Final time force",1.0);
204 funcParameter[4] =this->parameterList_->sublist("Parameter").get("dt",0.1);
205
206
207 if(nonlinearExternalForce_){
208 MatrixPtr_Type A( new Matrix_Type (this->system_->getBlock(0,0)));
209 //A->print();
210 MatrixPtr_Type AKext(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
211 MatrixPtr_Type Kext(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow()*2 ) );
212 MultiVectorPtr_Type Kext_vec;
213 this->feFactory_->assemblyNonlinearSurfaceIntegralExternal(this->dim_, this->getDomain(0)->getFEType(),FERhs, u_rep_,Kext, funcParameter, this->rhsFuncVec_[0],this->parameterList_);
214
215 A->addMatrix(1.,AKext,0.);
216 Kext->addMatrix(1.,AKext,1.);
217
218 AKext->fillComplete(this->getDomain(0)->getMapVecFieldUnique(),this->getDomain(0)->getMapVecFieldUnique());
219
220 this->system_->addBlock(AKext,0,0);
221
222 }
223 else
224 this->feFactory_->assemblySurfaceIntegralExternal(this->dim_, this->getDomain(0)->getFEType(),FERhs, u_rep_, funcParameter, this->rhsFuncVec_[0],this->parameterList_);
225
226
227 this->sourceTerm_->getBlockNonConst(0)->exportFromVector( FERhs, false, "Add" );
228 }
229 else
230 this->assembleSourceTerm( timeSteppingTool_->t_ );
231 //this->problemTimeStructure_->getSourceTerm()->scale(density);
232 // Fuege die rechte Seite der DGL (f bzw. f_{n+1}) der rechten Seite hinzu (skaliert mit coeffSourceTerm)
233 // Die Skalierung mit der Dichte erfolgt schon in der Assemblierungsfunktion!
234
235 // addSourceTermToRHS() aus DAESolverInTime
236 double coeffSourceTermStructure = 1.0;
237 BlockMultiVectorPtr_Type tmpPtr= this->sourceTerm_;
238
239 this->getRhs()->update(coeffSourceTermStructure, *tmpPtr, 1.);
240 this->rhs_->addBlock( this->getRhs()->getBlockNonConst(0), 0 );
241
242 }
243}
244
245
246template<class SC,class LO,class GO,class NO>
247void NonLinElasticity<SC,LO,GO,NO>::calculateNonLinResidualVec(std::string type, double time) const{
248
249 this->reAssemble("Newton-Residual");
250
251#ifdef FEDD_HAVE_ACEGENINTERFACE
252 if(this->getFEType(0) =="P2" && this->parameterList_->sublist("Parameter").get("Use AceGen Interface", true)){
253 if(this->parameterList_->sublist("Parameter").get("SCI",false) == true || this->parameterList_->sublist("Parameter").get("FSCI",false) == true ){
254 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*/);
255 }
256 }
257 #endif
258 if(this->parameterList_->sublist("Parameter").get("SCI",false) == false && this->parameterList_->sublist("Parameter").get("FSCI",false) == false ){
259
260 if (!type.compare("standard")){
261 this->residualVec_->update(-1.,*this->rhs_,1.);
262 //if ( !this->sourceTerm_.is_null() )
263 // this->residualVec_->update(-1.,*this->sourceTerm_,1.);
264 }
265 else if(!type.compare("reverse")){
266 this->residualVec_->update(1.,*this->rhs_,-1.); // this = -1*this + 1*rhs
267 //if ( !this->sourceTerm_.is_null() )
268 // this->residualVec_->update(1.,*this->sourceTerm_,1.);
269 }
270 else{
271 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Unknown type for residual computation.");
272 }
273
274 // this might be set again by the TimeProblem after adding of M*u
275 this->bcFactory_->setBCMinusVector( this->residualVec_, this->solution_, time );
276 }
277}
278}
279#endif
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 NonLinElasticity_def.hpp:247
virtual void assemble(std::string type="") const
assemble of type exectuted by the derived specific non-linear problem classes
Definition NonLinElasticity_def.hpp:57
Definition NonLinearProblem_decl.hpp:28
BlockMultiVectorPtr_Type sourceTerm_
Definition Problem_decl.hpp:251
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement_decl.hpp:36