v0.15.0
Loading...
Searching...
No Matches
NonlinearElasticElementInterface Struct Reference

Set of functions declaring elements and setting operators for generic element interface. More...

#include "users_modules/basic_finite_elements/nonlinear_elastic_materials/src/NonlinearElasticElementInterface.hpp"

Inheritance diagram for NonlinearElasticElementInterface:
[legend]
Collaboration diagram for NonlinearElasticElementInterface:
[legend]

Public Member Functions

 NonlinearElasticElementInterface (MoFEM::Interface &m_field, string postion_field, string mesh_posi_field_name="MESH_NODE_POSITIONS", bool is_displacement_field=true, PetscBool is_quasi_static=PETSC_TRUE)
 
 ~NonlinearElasticElementInterface ()
 
MoFEMErrorCode getCommandLineParameters ()
 
MoFEMErrorCode addElementFields ()
 
MoFEMErrorCode createElements ()
 
MoFEMErrorCode setOperators ()
 
BitRefLevel getBitRefLevel ()
 
MoFEMErrorCode addElementsToDM (SmartPetscObj< DM > dm)
 
MoFEMErrorCode setupSolverJacobianSNES ()
 
MoFEMErrorCode setupSolverFunctionSNES ()
 
MoFEMErrorCode setupSolverJacobianTS (const TSType type)
 
MoFEMErrorCode setupSolverFunctionTS (const TSType type)
 
MoFEMErrorCode updateElementVariables ()
 
MoFEMErrorCode postProcessElement (int step)
 
- Public Member Functions inherited from GenericElementInterface
 GenericElementInterface ()
 
virtual ~GenericElementInterface ()
 
virtual MoFEMErrorCode setGlobalBoundaryMarker (BcMarkerPtr mark)
 
virtual BcMarkerPtr getGlobalBoundaryMarker ()
 
virtual MoFEMErrorCode setMonitorPtr (boost::shared_ptr< MoFEM::FEMethod > monitor_ptr)
 
virtual BitRefLevel getBitRefLevelMask ()
 
virtual MoFEMErrorCode opFactoryDomainRhs (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
 
virtual MoFEMErrorCode opFactoryDomainLhs (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
 
virtual MoFEMErrorCode setUpdateElementVariablesOperators ()
 
virtual MoFEMErrorCode updateElementVariables (SmartPetscObj< DM > dm, string fe_name)
 
virtual MoFEMErrorCode postProcessElement (int step, SmartPetscObj< DM > dm, string fe_name)
 
 GenericElementInterface ()
 
virtual ~GenericElementInterface ()
 
virtual MoFEMErrorCode setGlobalBoundaryMarker (BcMarkerPtr mark)
 
virtual BcMarkerPtr getGlobalBoundaryMarker ()
 
virtual MoFEMErrorCode setMonitorPtr (boost::shared_ptr< FEMethod > monitor_ptr)
 
virtual BitRefLevel getBitRefLevelMask ()
 
virtual MoFEMErrorCode opFactoryDomainRhs (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
 
virtual MoFEMErrorCode opFactoryDomainLhs (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
 
virtual MoFEMErrorCode setUpdateElementVariablesOperators ()
 
virtual MoFEMErrorCode updateElementVariables (SmartPetscObj< DM > dm, string fe_name)
 
virtual MoFEMErrorCode postProcessElement (int step, SmartPetscObj< DM > dm, string fe_name)
 

Public Attributes

MoFEM::InterfacemField
 
SmartPetscObj< DM > dM
 
PetscBool isQuasiStatic
 
PetscInt oRder
 
bool isDisplacementField
 
BitRefLevel bIt
 
boost::shared_ptr< NonlinearElasticElementelasticElementPtr
 
boost::shared_ptr< ElasticMaterialselasticMaterialsPtr
 
boost::shared_ptr< PostProcVolumeOnRefinedMeshpostProcMeshPtr
 
string positionField
 
string meshNodeField
 
- Public Attributes inherited from GenericElementInterface
BcMarkerPtr mBoundaryMarker
 
int atomTest
 
int restartRunStep
 
boost::shared_ptr< MoFEM::FEMethodmonitorPtr
 
boost::shared_ptr< FEMethodmonitorPtr
 

Additional Inherited Members

- Public Types inherited from GenericElementInterface
enum  TSType {
  EX , IM , IM2 , IMEX ,
  DEFAULT , EX , IM , IM2 ,
  IMEX , DEFAULT
}
 
enum  TSType {
  EX , IM , IM2 , IMEX ,
  DEFAULT , EX , IM , IM2 ,
  IMEX , DEFAULT
}
 
using BcMarkerPtr = boost::shared_ptr<std::vector<char unsigned>>
 

Detailed Description

Set of functions declaring elements and setting operators for generic element interface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 15 of file NonlinearElasticElementInterface.hpp.

Constructor & Destructor Documentation

◆ NonlinearElasticElementInterface()

NonlinearElasticElementInterface::NonlinearElasticElementInterface ( MoFEM::Interface & m_field,
string postion_field,
string mesh_posi_field_name = "MESH_NODE_POSITIONS",
bool is_displacement_field = true,
PetscBool is_quasi_static = PETSC_TRUE )
inline

◆ ~NonlinearElasticElementInterface()

NonlinearElasticElementInterface::~NonlinearElasticElementInterface ( )
inline

Member Function Documentation

◆ addElementFields()

MoFEMErrorCode NonlinearElasticElementInterface::addElementFields ( )
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 56 of file NonlinearElasticElementInterface.hpp.

56 {
58 auto simple = mField.getInterface<Simple>();
61 3);
62 CHKERR simple->addBoundaryField(positionField, H1,
64 CHKERR simple->setFieldOrder(positionField, oRder);
65 }
68 3);
69 CHKERR simple->setFieldOrder(meshNodeField, 2);
70 }
71
73 };
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
Definition definitions.h:60
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
virtual bool check_field(const std::string &name) const =0
check if field is in database
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ addElementsToDM()

MoFEMErrorCode NonlinearElasticElementInterface::addElementsToDM ( SmartPetscObj< DM > dm)
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 112 of file NonlinearElasticElementInterface.hpp.

112 {
114 this->dM = dm;
115 CHKERR DMMoFEMAddElement(dM, "ELASTIC");
116 mField.getInterface<Simple>()->getOtherFiniteElements().push_back(
117 "ELASTIC");
118
120 };
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488

◆ createElements()

MoFEMErrorCode NonlinearElasticElementInterface::createElements ( )
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 75 of file NonlinearElasticElementInterface.hpp.

75 {
77
78 elasticElementPtr = boost::make_shared<NonlinearElasticElement>(mField, 2);
79 elasticMaterialsPtr = boost::make_shared<ElasticMaterials>(mField, "HOOKE");
80 CHKERR elasticMaterialsPtr->setBlocks(elasticElementPtr->setOfBlocks);
81
83 false, false, false);
85 false, false, false);
87 true, false, false, false);
88 CHKERR elasticElementPtr->addElement("ELASTIC", positionField,
89 meshNodeField, false);
90
91
93 };
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
boost::shared_ptr< ElasticMaterials > elasticMaterialsPtr
boost::shared_ptr< NonlinearElasticElement > elasticElementPtr

◆ getBitRefLevel()

BitRefLevel NonlinearElasticElementInterface::getBitRefLevel ( )
inlinevirtual

◆ getCommandLineParameters()

MoFEMErrorCode NonlinearElasticElementInterface::getCommandLineParameters ( )
inlinevirtual

Reimplemented from GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 45 of file NonlinearElasticElementInterface.hpp.

45 {
47 isQuasiStatic = PETSC_FALSE;
48 oRder = 2;
49 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "-is_quasi_static", &isQuasiStatic,
50 PETSC_NULLPTR);
51 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "-order", &oRder, PETSC_NULLPTR);
52
54 };
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)

◆ postProcessElement()

MoFEMErrorCode NonlinearElasticElementInterface::postProcessElement ( int step)
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 181 of file NonlinearElasticElementInterface.hpp.

181 {
183
184 if (elasticElementPtr->setOfBlocks.empty())
186
187 if (!postProcMeshPtr) {
188 postProcMeshPtr = boost::make_shared<PostProcVolumeOnRefinedMesh>(mField);
189
190 if (mField.check_field("MESH_NODE_POSITIONS"))
191 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", *postProcMeshPtr, true, false,
192 false, false);
193 CHKERR postProcMeshPtr->generateReferenceElementMesh();
194 CHKERR postProcMeshPtr->addFieldValuesPostProc(positionField);
195 CHKERR postProcMeshPtr->addFieldValuesPostProc(meshNodeField);
196 CHKERR postProcMeshPtr->addFieldValuesGradientPostProc(positionField);
197
198 for (auto &sit : elasticElementPtr->setOfBlocks) {
199 postProcMeshPtr->getOpPtrVector().push_back(new PostProcStress(
200 postProcMeshPtr->postProcMesh, postProcMeshPtr->mapGaussPts,
201 positionField, sit.second, postProcMeshPtr->commonData,
203 }
204 }
205
206 elasticElementPtr->getLoopFeEnergy().snes_ctx = SnesMethod::CTX_SNESNONE;
207 elasticElementPtr->getLoopFeEnergy().eNergy = 0;
208 // MOFEM_LOG("WORLD", Sev::inform) << "Loop energy\n";
210 &elasticElementPtr->getLoopFeEnergy());
211
212 auto E = elasticElementPtr->getLoopFeEnergy().eNergy;
213 // Print elastic energy
214 MOFEM_LOG_C("WORLD", Sev::inform, "%d Time %3.2e Elastic energy %3.2e",
215 step, elasticElementPtr->getLoopFeRhs().ts_t, E);
216
218 auto out_name = "out_vol_" + to_string(step) + ".h5m";
219
220 CHKERR postProcMeshPtr->writeFile(out_name);
221
223 };
#define MOFEM_LOG_C(channel, severity, format,...)
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
boost::shared_ptr< PostProcVolumeOnRefinedMesh > postProcMeshPtr

◆ setOperators()

MoFEMErrorCode NonlinearElasticElementInterface::setOperators ( )
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 95 of file NonlinearElasticElementInterface.hpp.

95 {
97 auto &pipeline_rhs = elasticElementPtr->feRhs.getOpPtrVector();
98 auto &pipeline_lhs = elasticElementPtr->feLhs.getOpPtrVector();
99
100 pipeline_rhs.push_back(new OpSetBc(positionField, true, mBoundaryMarker));
101 pipeline_lhs.push_back(new OpSetBc(positionField, true, mBoundaryMarker));
102
105
106 pipeline_rhs.push_back(new OpUnSetBc(positionField));
107 pipeline_lhs.push_back(new OpUnSetBc(positionField));
109 }

◆ setupSolverFunctionSNES()

MoFEMErrorCode NonlinearElasticElementInterface::setupSolverFunctionSNES ( )
inlinevirtual

Reimplemented from GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 130 of file NonlinearElasticElementInterface.hpp.

130 {
133 &elasticElementPtr->getLoopFeRhs(),
134 PETSC_NULLPTR, PETSC_NULLPTR);
136 };
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition DMMoFEM.cpp:708

◆ setupSolverFunctionTS()

MoFEMErrorCode NonlinearElasticElementInterface::setupSolverFunctionTS ( const TSType type)
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 159 of file NonlinearElasticElementInterface.hpp.

159 {
161 auto &method = elasticElementPtr->getLoopFeRhs();
162 switch (type) {
163 case IM:
164 CHKERR DMMoFEMTSSetIFunction(dM, "ELASTIC", &method, &method, &method);
165 break;
166 case IM2:
167 CHKERR DMMoFEMTSSetI2Function(dM, "ELASTIC", &method, &method, &method);
168 break;
169 case EX:
170 CHKERR DMMoFEMTSSetRHSFunction(dM, "ELASTIC", &method, &method, &method);
171 break;
172 default:
173 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not implemented");
174 break;
175 }
176
178 };
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set TS implicit function evaluation function
Definition DMMoFEM.cpp:790
PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS implicit function evaluation function
Definition DMMoFEM.cpp:965
PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side function
Definition DMMoFEM.cpp:872

◆ setupSolverJacobianSNES()

MoFEMErrorCode NonlinearElasticElementInterface::setupSolverJacobianSNES ( )
inlinevirtual

Reimplemented from GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 122 of file NonlinearElasticElementInterface.hpp.

122 {
124
126 dM, "ELASTIC", &elasticElementPtr->getLoopFeLhs(), NULL, NULL);
127
129 };
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition DMMoFEM.cpp:749

◆ setupSolverJacobianTS()

MoFEMErrorCode NonlinearElasticElementInterface::setupSolverJacobianTS ( const TSType type)
inlinevirtual

Implements GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 138 of file NonlinearElasticElementInterface.hpp.

138 {
140 auto &method = elasticElementPtr->getLoopFeLhs();
141 switch (type) {
142 case IM:
143 CHKERR DMMoFEMTSSetIJacobian(dM, "ELASTIC", &method, &method, &method);
144 break;
145 case IM2:
146 CHKERR DMMoFEMTSSetI2Jacobian(dM, "ELASTIC", &method, &method, &method);
147 break;
148 case EX:
149 CHKERR DMMoFEMTSSetRHSJacobian(dM, "ELASTIC", &method, &method, &method);
150 break;
151 default:
152 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
153 "This TS is not yet implemented");
154 break;
155 }
157 };
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition DMMoFEM.cpp:843
PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition DMMoFEM.cpp:1007
PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side jacobian
Definition DMMoFEM.cpp:912

◆ updateElementVariables()

MoFEMErrorCode NonlinearElasticElementInterface::updateElementVariables ( )
inlinevirtual

Reimplemented from GenericElementInterface.

Examples
NonlinearElasticElementInterface.hpp.

Definition at line 180 of file NonlinearElasticElementInterface.hpp.

180{ return 0; };

Member Data Documentation

◆ bIt

BitRefLevel NonlinearElasticElementInterface::bIt

◆ dM

SmartPetscObj<DM> NonlinearElasticElementInterface::dM

◆ elasticElementPtr

boost::shared_ptr<NonlinearElasticElement> NonlinearElasticElementInterface::elasticElementPtr

◆ elasticMaterialsPtr

boost::shared_ptr<ElasticMaterials> NonlinearElasticElementInterface::elasticMaterialsPtr

◆ isDisplacementField

bool NonlinearElasticElementInterface::isDisplacementField

◆ isQuasiStatic

PetscBool NonlinearElasticElementInterface::isQuasiStatic

◆ meshNodeField

string NonlinearElasticElementInterface::meshNodeField

◆ mField

MoFEM::Interface& NonlinearElasticElementInterface::mField

◆ oRder

PetscInt NonlinearElasticElementInterface::oRder

◆ positionField

string NonlinearElasticElementInterface::positionField

◆ postProcMeshPtr

boost::shared_ptr<PostProcVolumeOnRefinedMesh> NonlinearElasticElementInterface::postProcMeshPtr

The documentation for this struct was generated from the following file: