Finite Element Domain Decomposition Library
FEDDLib
Loading...
Searching...
No Matches
Elasticity_def.hpp
1#ifndef ELASTICITY_def_hpp
2#define ELASTICITY_def_hpp
3#include "Elasticity_decl.hpp"
12
13namespace FEDD {
14template<class SC,class LO,class GO,class NO>
15Elasticity<SC,LO,GO,NO>::Elasticity(const DomainConstPtr_Type &domain, std::string FEType, ParameterListPtr_Type parameterList)
16{
17 if(parameterList->sublist("Parameter").get("Material model","linear")=="linear")
18 linearProblem_ = Teuchos::rcp( new LinElas_Type( domain, FEType, parameterList ) );
19 else
20 nonLinearProblem_ = Teuchos::rcp( new NonLinElas_Type( domain, FEType, parameterList ) );
21}
22template<class SC,class LO,class GO,class NO>
23Elasticity<SC,LO,GO,NO>::~Elasticity(){
24 linearProblem_.reset();
25 nonLinearProblem_.reset();
26}
27
28template<class SC,class LO,class GO,class NO>
29void Elasticity<SC,LO,GO,NO>::info(){
30 this->infoProblem();
31 if(!nonLinearProblem_.is_null())
32 this->infoNonlinProblem();
33}
34
35template<class SC,class LO,class GO,class NO>
36void Elasticity<SC,LO,GO,NO>::assemble(){
37 if(!nonLinearProblem_.is_null())
38 nonLinearProblem_->assemble();
39 else
40 linearProblem_->assemble();
41}
42
43template<class SC,class LO,class GO,class NO>
44void Elasticity<SC,LO,GO,NO>::reAssemble(std::string type) const {
45
46 if(!nonLinearProblem_.is_null())
47 nonLinearProblem_->reAssemble( type );
48 else{
49 if(this->getVerbose())
50 std::cout << "-- Nothing to reassemble for linear elasticity --" << std::endl;
51 }
52
53}
54
55template<class SC,class LO,class GO,class NO>
56void Elasticity<SC,LO,GO,NO>::reAssembleExtrapolation(){
57
58 TEUCHOS_TEST_FOR_EXCEPTION( true, std::logic_error, "Only Newton/NOX implemented for nonlinear material models!");
59
60}
61
62template<class SC,class LO,class GO,class NO>
63void Elasticity<SC,LO,GO,NO>::evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs<SC> &inArgs,
64 const Thyra::ModelEvaluatorBase::OutArgs<SC> &outArgs
65 ) const
66{
67 if(!nonLinearProblem_.is_null())
68 nonLinearProblem_->evalModelImpl( inArgs, outArgs );
69 else{
70 if(this->getVerbose())
71 std::cout << "-- Nothing to evaluate for linear elasticity --" << std::endl;
72 }
73
74}
75
76template<class SC,class LO,class GO,class NO>
77Teuchos::RCP<Thyra::LinearOpBase<SC> > Elasticity<SC,LO,GO,NO>::create_W_op() const
78{
79
80 if(!nonLinearProblem_.is_null())
81 return nonLinearProblem_->create_W_op( );
82 else{
83 if(this->getVerbose())
84 std::cout << "-- No create_W_op for linear elasticity --" << std::endl;
85 Teuchos::RCP<Thyra::LinearOpBase<SC> > W_opDummy;
86 return W_opDummy;
87 }
88}
89
90template<class SC,class LO,class GO,class NO>
91Teuchos::RCP<Thyra::PreconditionerBase<SC> > Elasticity<SC,LO,GO,NO>::create_W_prec() const
92{
93
94 if(!nonLinearProblem_.is_null())
95 return nonLinearProblem_->create_W_prec( );
96 else{
97 if(this->getVerbose())
98 std::cout << "-- No create_W_prec for linear elasticity --" << std::endl;
99 Teuchos::RCP<Thyra::PreconditionerBase<SC> > W_precDummy;
100 return W_precDummy;
101 }
102}
103
104template<class SC,class LO,class GO,class NO>
105void Elasticity<SC,LO,GO,NO>::calculateNonLinResidualVec(std::string type) const{
106
107 if(!nonLinearProblem_.is_null())
108 return nonLinearProblem_->calculateNonLinResidualVec( type );
109 else{
110 if(this->getVerbose())
111 std::cout << "-- No computation of nonlinear residual for linear elasticity --" << std::endl;
112 }
113}
114}
115#endif
Adaptive Mesh Refinement.
Definition AdaptiveMeshRefinement.cpp:5