Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
NonLinElasticity_def.hpp
1#ifndef NonLinElasticity_def_hpp
2#define NonLinElasticity_def_hpp
3#include "NonLinElasticity_decl.hpp"
12
13
14namespace FEDD {
15template<class SC,class LO,class GO,class NO>
16NonLinElasticity<SC,LO,GO,NO>::NonLinElasticity(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// poissonRatio_ = this->parameterList_->sublist("Parameter").get("Poisson Ratio",0.4);
29// mue_ = this->parameterList_->sublist("Parameter").get("Mu",2.0e+6);
30// // Berechne daraus nun E (Youngsches Modul) und die erste Lamé-Konstante \lambda
31// E_ = mue_*2.*(1. + poissonRatio_);
32// lambda_ = (poissonRatio_*E_)/((1 + poissonRatio_)*(1 - 2*poissonRatio_));
33
34}
35
36template<class SC,class LO,class GO,class NO>
37NonLinElasticity<SC,LO,GO,NO>::~NonLinElasticity(){
38
39}
40
41template<class SC,class LO,class GO,class NO>
42void NonLinElasticity<SC,LO,GO,NO>::info(){
43 this->infoProblem();
44 this->infoNonlinProblem();
45}
46
47template<class SC,class LO,class GO,class NO>
48void NonLinElasticity<SC,LO,GO,NO>::assemble(std::string type) const{
49
50 if (type == ""){
51 if (this->verbose_)
52 std::cout << "-- Assembly nonlinear elasticity ... " << std::flush;
53
54 MatrixPtr_Type A(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
55
56 this->feFactory_->assemblyEmptyMatrix(A);
57
58 this->system_.reset(new BlockMatrix_Type(1));
59 this->system_->addBlock( A, 0, 0 );
60
61 double density = this->parameterList_->sublist("Parameter").get("Density",1000.);
62 std::string sourceType = this->parameterList_->sublist("Parameter").get("Source Type","volume");
63
64 this->assembleSourceTerm( 0. );
65 if(sourceType == "volume")
66 this->sourceTerm_->scale(density);
67
68 this->addToRhs( this->sourceTerm_ );
69
70 this->setBoundariesRHS();
71
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
85 this->reAssemble(type);
86}
87
88template<class SC,class LO,class GO,class NO>
89void NonLinElasticity<SC,LO,GO,NO>::reAssemble(std::string type) const {
90 std::string material_model = this->parameterList_->sublist("Parameter").get("Material model","Neo-Hooke");
91
92 if (this->verbose_)
93 std::cout << "-- Reassembly nonlinear elasticity with material model " << material_model <<" ("<<type <<") ... " << std::flush;
94
95 if (type=="Newton-Residual") {
96 MultiVectorConstPtr_Type u = this->solution_->getBlock(0);
97
98 u_rep_->importFromVector(u, true);
99
100 MultiVectorPtr_Type f = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated(), 1 ) );
101 MatrixPtr_Type W = Teuchos::rcp(new Matrix_Type( this->getDomain(0)->getMapVecFieldUnique(), this->getDomain(0)->getDimension() * this->getDomain(0)->getApproxEntriesPerRow() ) );
102 this->feFactory_->assemblyElasticityJacobianAndStressAceFEM(this->dim_, this->getDomain(0)->getFEType(), W, f, u_rep_, this->parameterList_, C_);
103
104 MultiVectorPtr_Type fUnique = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldUnique(), 1 ) );
105 fUnique->putScalar(0.);
106 fUnique->exportFromVector( f, true, "Add" );
107
108 this->residualVec_->addBlock( fUnique, 0 );
109 this->system_->addBlock( W, 0, 0 );
110
111 // MultiVectorPtr_Type f = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldRepeated(), 1 ) );
112// this->feFactory_->assemblyElasticityStressesAceFEM(this->dim_, this->getDomain(0)->getFEType(), f, u_rep_, material_model, E_, nu_, C_);
113// MultiVectorPtr_Type fUnique = Teuchos::rcp( new MultiVector_Type( this->getDomain(0)->getMapVecFieldUnique(), 1 ) );
114// fUnique->putScalar(0.);
115// fUnique->importFromVector( f, true, "Add" );
116// this->residualVec_->addBlock( fUnique, 0 );
117 }
118 else if(type=="Newton"){ //we already assemble the new tangent when we calculate the stresses above
119
120 }
121 if (this->verbose_)
122 std::cout << "done -- " << std::endl;
123}
124
125template<class SC,class LO,class GO,class NO>
126void NonLinElasticity<SC,LO,GO,NO>::reAssemble( MatrixPtr_Type& massmatrix, std::string type ) const
127{
128
129}
130
131template<class SC,class LO,class GO,class NO>
132void NonLinElasticity<SC,LO,GO,NO>::reAssembleExtrapolation(BlockMultiVectorPtrArray_Type previousSolutions){
133
134
135 TEUCHOS_TEST_FOR_EXCEPTION( true, std::logic_error, "Only Newton/NOX implemented for nonlinear material models!");
136
137}
138
139template<class SC,class LO,class GO,class NO>
140void NonLinElasticity<SC,LO,GO,NO>::evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
141 const Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs
142 ) const
143{
144 using Teuchos::RCP;
145 using Teuchos::rcp;
146 using Teuchos::rcp_dynamic_cast;
147 using Teuchos::rcp_const_cast;
148 using Teuchos::ArrayView;
149 using Teuchos::Array;
150
151 TEUCHOS_TEST_FOR_EXCEPTION( this->solution_->getBlock(0)->getMap()->getUnderlyingLib() != "Tpetra", std::runtime_error, "Use of NOX only supports Tpetra. Epetra support must be implemented.");
152 RCP<Teuchos::FancyOStream> fancy = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
153 TEUCHOS_TEST_FOR_EXCEPTION( inArgs.get_x().is_null(), std::logic_error, "inArgs.get_x() is null.");
154
155 RCP< const Thyra::VectorBase< SC > > vecThyra = inArgs.get_x();
156 RCP<Teuchos::FancyOStream> out = Teuchos::VerboseObjectBase::getDefaultOStream();
157
158 RCP< Thyra::VectorBase< SC > > vecThyraNonConst = rcp_const_cast<Thyra::VectorBase< SC > >(vecThyra);
159
160 this->solution_->fromThyraMultiVector(vecThyraNonConst);
161
162 const RCP<Thyra::MultiVectorBase<SC> > f_out = outArgs.get_f();
163 const RCP<Thyra::LinearOpBase<SC> > W_out = outArgs.get_W_op();
164 const RCP<Thyra::PreconditionerBase<SC> > W_prec_out = outArgs.get_W_prec();
165
166 typedef Thyra::TpetraOperatorVectorExtraction<SC,LO,GO,NO> tpetra_extract;
167 typedef Tpetra::CrsMatrix<SC,LO,GO,NO> TpetraMatrix_Type;
168 typedef RCP<TpetraMatrix_Type> TpetraMatrixPtr_Type;
169 typedef RCP<const TpetraMatrix_Type> TpetraMatrixConstPtr_Type;
170
171 const bool fill_f = nonnull(f_out);
172 const bool fill_W = nonnull(W_out);
173 const bool fill_W_prec = nonnull(W_prec_out);
174
175 if ( fill_f || fill_W || fill_W_prec ) {
176
177 // ****************
178 // Get the underlying xpetra objects
179 // ****************
180 if (fill_f) {
181
182 this->calculateNonLinResidualVec("standard");
183
184 RCP<Thyra::MultiVectorBase<SC> > f_thyra = this->getResidualVector()->getThyraMultiVector();
185 f_out->assign(*f_thyra);
186 }
187
188 TpetraMatrixPtr_Type W;
189 if (fill_W) {
190
191 this->reAssemble("Newton");
192 this->setBoundariesSystem();
193
194 Teuchos::RCP<TpetraOp_Type> W_tpetra = tpetra_extract::getTpetraOperator(W_out);
195 Teuchos::RCP<TpetraMatrix_Type> W_tpetraMat = Teuchos::rcp_dynamic_cast<TpetraMatrix_Type>(W_tpetra);
196
197 TpetraMatrixConstPtr_Type W_systemTpetra = this->getSystem()->getBlock( 0, 0 )->getTpetraMatrix();
198 TpetraMatrixPtr_Type W_systemTpetraNonConst = rcp_const_cast<TpetraMatrix_Type>(W_systemTpetra);
199
200 //Tpetra::CrsMatrixWrap<SC,LO,GO,NO>& crsOp = dynamic_cast<Xpetra::CrsMatrixWrap<SC,LO,GO,NO>&>(*W_systemXpetraNonConst);
201 //Xpetra::TpetraCrsMatrix<SC,LO,GO,NO>& xTpetraMat = dynamic_cast<Xpetra::TpetraCrsMatrix<SC,LO,GO,NO>&>(*crsOp.getCrsMatrix());
202
203 Teuchos::RCP<TpetraMatrix_Type> tpetraMatTpetra = W_systemTpetraNonConst; //xTpetraMat.getTpetra_CrsMatrixNonConst();
204
205
206 W_tpetraMat->resumeFill();
207
208 for (auto i=0; i<tpetraMatTpetra->getMap()->getLocalNumElements(); i++) {
209 typename Tpetra::CrsMatrix<SC,LO,GO,NO>::local_inds_host_view_type indices; //ArrayView< const LO > indices
210 typename Tpetra::CrsMatrix<SC,LO,GO,NO>::values_host_view_type values;
211 tpetraMatTpetra->getLocalRowView( i, indices, values);
212 W_tpetraMat->replaceLocalValues( i, indices, values);
213 }
214 W_tpetraMat->fillComplete();
215
216 }
217
218 if (fill_W_prec) {
219 this->setupPreconditioner( "Monolithic" );
220
221 // ch 26.04.19: After each setup of the preconditioner we check if we use a two-level precondtioner with multiplicative combination between the levels.
222 // If this is the case, we need to pre apply the coarse level to the residual(f_out).
223
224 std::string levelCombination = this->parameterList_->sublist("ThyraPreconditioner").sublist("Preconditioner Types").sublist("FROSch").get("Level Combination","Additive");
225 if (!levelCombination.compare("Multiplicative")) {
226 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Multiplicative Level Combination is not supported for NOX. In general we need to pre-apply the coarse problem. This must be implemented here.");
227 }
228
229 }
230 }
231}
232
233template<class SC,class LO,class GO,class NO>
234Teuchos::RCP<Thyra::LinearOpBase<SC> > NonLinElasticity<SC,LO,GO,NO>::create_W_op() const
235{
236
237 Teuchos::RCP<const Thyra::LinearOpBase<SC> > W_opConst = this->system_->getThyraLinOp();
238 Teuchos::RCP<Thyra::LinearOpBase<SC> > W_op = Teuchos::rcp_const_cast<Thyra::LinearOpBase<SC> >(W_opConst);
239 return W_op;
240}
241
242template<class SC,class LO,class GO,class NO>
243Teuchos::RCP<Thyra::PreconditionerBase<SC> > NonLinElasticity<SC,LO,GO,NO>::create_W_prec() const
244{
245 this->initializeSolverBuilder();
246 this->initializePreconditioner();
247
248 Teuchos::RCP<const Thyra::PreconditionerBase<SC> > thyraPrec = this->getPreconditionerConst()->getThyraPrecConst();
249 Teuchos::RCP<Thyra::PreconditionerBase<SC> > thyraPrecNonConst= Teuchos::rcp_const_cast<Thyra::PreconditionerBase<SC> >(thyraPrec);
250
251 return thyraPrecNonConst;
252}
253
254template<class SC,class LO,class GO,class NO>
255void NonLinElasticity<SC,LO,GO,NO>::calculateNonLinResidualVec(std::string type, double time) const{
256
257 this->reAssemble("Newton-Residual");
258 if (!type.compare("standard")){
259 this->residualVec_->update(-1.,*this->rhs_,1.);
260 //if ( !this->sourceTerm_.is_null() )
261 // this->residualVec_->update(-1.,*this->sourceTerm_,1.);
262 }
263 else if(!type.compare("reverse")){
264 this->residualVec_->update(1.,*this->rhs_,-1.); // this = -1*this + 1*rhs
265 //if ( !this->sourceTerm_.is_null() )
266 // this->residualVec_->update(1.,*this->sourceTerm_,1.);
267 }
268 else{
269 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Unknown type for residual computation.");
270 }
271
272 // this might be set again by the TimeProblem after adding of M*u
273 this->bcFactory_->setBCMinusVector( this->residualVec_, this->solution_, time );
274
275}
276}
277#endif
Definition NonLinearProblem_decl.hpp:24
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5