v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
NonlinearElasticElement::OpJacobianPiolaKirchhoffStress Struct Reference

Operator performs automatic differentiation. More...

#include "users_modules/basic_finite_elements/src/NonLinearElasticElement.hpp"

Inheritance diagram for NonlinearElasticElement::OpJacobianPiolaKirchhoffStress:
[legend]
Collaboration diagram for NonlinearElasticElement::OpJacobianPiolaKirchhoffStress:
[legend]

Public Member Functions

 OpJacobianPiolaKirchhoffStress (const std::string field_name, BlockData &data, CommonData &common_data, int tag, bool jacobian, bool ale, bool field_disp)
 Construct operator to calculate Piola-Kirchhoff stress or its derivatives over gradient deformation.
 
virtual MoFEMErrorCode calculateStress (const int gg)
 Calculate Paola-Kirchhoff I stress.
 
virtual MoFEMErrorCode recordTag (const int gg)
 Record ADOL-C tape.
 
virtual MoFEMErrorCode playTag (const int gg)
 Play ADOL-C tape.
 
virtual bool recordTagForIntegrationPoint (const int gg)
 Cgeck if tape is recorded for given integration point.
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 Calculate stress or jacobian at gauss points.
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes
 
const EntityHandlegetConn ()
 get element connectivity
 
double getVolume () const
 element volume (linear geometry)
 
doublegetVolume ()
 element volume (linear geometry)
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian
 
VectorDoublegetCoords ()
 nodal coordinates
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 Constructor for operators working on finite element spaces.
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 Constructor for operators working on a single field.
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 Constructor for operators working on two fields (bilinear forms)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEM::InterfacegetMField ()
 
moab::Interface & getMoab ()
 
virtual boost::weak_ptr< ForcesAndSourcesCoregetSubPipelinePtr () const
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side.
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side.
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Public Attributes

BlockDatadAta
 
CommonDatacommonData
 
int tAg
 
int adlocReturnValue
 
bool jAcobian
 if true Jacobian is calculated
 
bool fUnction
 if true stress i calculated
 
bool aLe
 true if arbitrary Lagrangian-Eulerian formulation
 
bool fieldDisp
 
VectorDouble activeVariables
 
int nbActiveVariables
 
std::vector< MatrixDouble > * ptrh
 
std::vector< MatrixDouble > * ptrH
 
FTensor::Index< 'i', 3 > i
 
FTensor::Index< 'j', 3 > j
 
FTensor::Index< 'k', 3 > k
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Operator performs automatic differentiation.

Definition at line 370 of file NonLinearElasticElement.hpp.

Constructor & Destructor Documentation

◆ OpJacobianPiolaKirchhoffStress()

NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::OpJacobianPiolaKirchhoffStress ( const std::string  field_name,
BlockData data,
CommonData common_data,
int  tag,
bool  jacobian,
bool  ale,
bool  field_disp 
)

Construct operator to calculate Piola-Kirchhoff stress or its derivatives over gradient deformation.

Parameters
field_nameapproximation field name of spatial positions or displacements
datareference to block data (what is Young modulus, Poisson ratio or what elements are part of the block)
tagadol-c unique tag of the tape
jacobianif true derivative of Piola Stress is calculated otherwise just stress is calculated
field_dispif true approximation field keeps displacements not spatial positions

Definition at line 198 of file NonLinearElasticElement.cpp.

203 : VolumeElementForcesAndSourcesCore::UserDataOperator(
205 dAta(data), commonData(common_data), tAg(tag), adlocReturnValue(0),
206 jAcobian(jacobian), fUnction(!jacobian), aLe(ale), fieldDisp(field_disp) {
207
208}
constexpr auto field_name
@ OPROW
operator doWork function is executed on FE rows
bool aLe
true if arbitrary Lagrangian-Eulerian formulation

Member Function Documentation

◆ calculateStress()

MoFEMErrorCode NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::calculateStress ( const int  gg)
virtual

Calculate Paola-Kirchhoff I stress.

Returns
error code

Reimplemented in NonlinearElasticElement::OpJacobianEshelbyStress, and Smoother::OpJacobianSmoother.

Definition at line 211 of file NonLinearElasticElement.cpp.

212 {
214
215 CHKERR dAta.materialAdoublePtr->calculateP_PiolaKirchhoffI(
217
218 if (aLe) {
219 auto &t_P = dAta.materialAdoublePtr->t_P;
220 auto &t_invH = dAta.materialAdoublePtr->t_invH;
221 t_P(i, j) = t_P(i, k) * t_invH(j, k);
222 t_P(i, j) *= dAta.materialAdoublePtr->detH;
223 }
224
225 commonData.sTress[gg].resize(3, 3, false);
226 for (int dd1 = 0; dd1 < 3; dd1++) {
227 for (int dd2 = 0; dd2 < 3; dd2++) {
228 dAta.materialAdoublePtr->P(dd1, dd2) >>=
229 (commonData.sTress[gg])(dd1, dd2);
230 }
231 }
232
234}
#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()
#define CHKERR
Inline error check.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble > > materialAdoublePtr

◆ doWork()

MoFEMErrorCode NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::doWork ( int  row_side,
EntityType  row_type,
EntitiesFieldData::EntData row_data 
)

Calculate stress or jacobian at gauss points.

Parameters
row_side
row_type
row_data
Returns
error code

Definition at line 338 of file NonLinearElasticElement.cpp.

340 {
342
343 // do it only once, no need to repeat this for edges,faces or tets
344 if (row_type != MBVERTEX)
346
347 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
348 dAta.tEts.end()) {
350 }
351
352 int nb_dofs = row_data.getFieldData().size();
353 if (nb_dofs == 0)
355 dAta.materialAdoublePtr->commonDataPtr = &commonData;
356 dAta.materialAdoublePtr->opPtr = this;
357
358 int nb_gauss_pts = row_data.getN().size1();
359 commonData.sTress.resize(nb_gauss_pts);
360 commonData.jacStress.resize(nb_gauss_pts);
361
363 if (aLe) {
365 }
366
367 for (int gg = 0; gg != nb_gauss_pts; gg++) {
368
369 dAta.materialAdoublePtr->gG = gg;
370
371 // Record tag and calculate stress
373 CHKERR recordTag(gg);
374 }
375
376 // Set active variables vector
378 activeVariables.resize(nbActiveVariables, false);
379 if (!aLe) {
380 for (int dd1 = 0; dd1 < 3; dd1++) {
381 for (int dd2 = 0; dd2 < 3; dd2++) {
382 activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
383 }
384 }
385 } else {
386 for (int dd1 = 0; dd1 < 3; dd1++) {
387 for (int dd2 = 0; dd2 < 3; dd2++) {
388 activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
389 }
390 }
391 for (int dd1 = 0; dd1 < 3; dd1++) {
392 for (int dd2 = 0; dd2 < 3; dd2++) {
393 activeVariables(9 + dd1 * 3 + dd2) = (*ptrH)[gg](dd1, dd2);
394 }
395 }
396 }
397 CHKERR dAta.materialAdoublePtr->setUserActiveVariables(activeVariables);
398
399 // Play tag and calculate stress or tangent
401 CHKERR playTag(gg);
402 }
403 }
404 }
405
407}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
Get DOF values on entity.
Range tEts
constrains elements in block set
std::vector< MatrixDouble > jacStress
this is simply material tangent operator
std::map< std::string, std::vector< MatrixDouble > > gradAtGaussPts
virtual MoFEMErrorCode recordTag(const int gg)
Record ADOL-C tape.
virtual MoFEMErrorCode playTag(const int gg)
Play ADOL-C tape.
virtual bool recordTagForIntegrationPoint(const int gg)
Cgeck if tape is recorded for given integration point.

◆ playTag()

MoFEMErrorCode NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::playTag ( const int  gg)
virtual

Play ADOL-C tape.

Returns
error code

Definition at line 303 of file NonLinearElasticElement.cpp.

303 {
305
306 int r;
307
308 if (fUnction) {
309 commonData.sTress[gg].resize(3, 3, false);
310 // play recorder for values
311 r = ::function(tAg, 9, nbActiveVariables, &activeVariables[0],
312 &commonData.sTress[gg](0, 0));
313 if (r < adlocReturnValue) { // function is locally analytic
314 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
315 "ADOL-C function evaluation with error r = %d", r);
316 }
317 }
318
319 if (jAcobian) {
320 commonData.jacStress[gg].resize(9, nbActiveVariables, false);
321 double *jac_ptr[] = {
322 &(commonData.jacStress[gg](0, 0)), &(commonData.jacStress[gg](1, 0)),
323 &(commonData.jacStress[gg](2, 0)), &(commonData.jacStress[gg](3, 0)),
324 &(commonData.jacStress[gg](4, 0)), &(commonData.jacStress[gg](5, 0)),
325 &(commonData.jacStress[gg](6, 0)), &(commonData.jacStress[gg](7, 0)),
326 &(commonData.jacStress[gg](8, 0))};
327 // play recorder for jacobians
328 r = jacobian(tAg, 9, nbActiveVariables, &activeVariables[0], jac_ptr);
329 if (r < adlocReturnValue) {
330 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
331 "ADOL-C function evaluation with error");
332 }
333 }
334
336}
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
int r
Definition sdf.py:205

◆ recordTag()

MoFEMErrorCode NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::recordTag ( const int  gg)
virtual

Record ADOL-C tape.

Returns
error code

Definition at line 237 of file NonLinearElasticElement.cpp.

238 {
240
241 trace_on(tAg, 0);
242
243 dAta.materialAdoublePtr->F.resize(3, 3, false);
244
245 if (!aLe) {
246
248 for (int dd1 = 0; dd1 < 3; dd1++) {
249 for (int dd2 = 0; dd2 < 3; dd2++) {
250 dAta.materialAdoublePtr->F(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
251 if (fieldDisp) {
252 if (dd1 == dd2) {
253 dAta.materialAdoublePtr->F(dd1, dd2) += 1;
254 }
255 }
257 }
258 }
259
260 } else {
261
263
264 dAta.materialAdoublePtr->h.resize(3, 3, false);
265 for (int dd1 = 0; dd1 < 3; dd1++) {
266 for (int dd2 = 0; dd2 < 3; dd2++) {
267 dAta.materialAdoublePtr->h(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
269 }
270 }
271
272 dAta.materialAdoublePtr->H.resize(3, 3, false);
273 for (int dd1 = 0; dd1 < 3; dd1++) {
274 for (int dd2 = 0; dd2 < 3; dd2++) {
275 dAta.materialAdoublePtr->H(dd1, dd2) <<= (*ptrH)[gg](dd1, dd2);
277 }
278 }
279
281 dAta.materialAdoublePtr->invH.resize(3, 3, false);
285
286 auto &t_F = dAta.materialAdoublePtr->t_F;
287 auto &t_h = dAta.materialAdoublePtr->t_h;
288 auto &t_invH = dAta.materialAdoublePtr->t_invH;
289
290 t_F(i, j) = t_h(i, k) * t_invH(k, j);
291
292 }
293
294 CHKERR dAta.materialAdoublePtr->setUserActiveVariables(nbActiveVariables);
296
297 trace_off();
298
300}
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
virtual MoFEMErrorCode calculateStress(const int gg)
Calculate Paola-Kirchhoff I stress.

◆ recordTagForIntegrationPoint()

virtual bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::recordTagForIntegrationPoint ( const int  gg)
inlinevirtual

Cgeck if tape is recorded for given integration point.

Parameters
ggintegration point
Returns
true if tag is recorded

Definition at line 437 of file NonLinearElasticElement.hpp.

437 {
438 // return true;
439 if (gg == 0)
440 return true;
441 return false;
442 }

Member Data Documentation

◆ activeVariables

VectorDouble NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::activeVariables

Definition at line 404 of file NonLinearElasticElement.hpp.

◆ adlocReturnValue

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::adlocReturnValue

return value from ADOL-C, if non-zero that is error.

Definition at line 379 of file NonLinearElasticElement.hpp.

◆ aLe

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::aLe

true if arbitrary Lagrangian-Eulerian formulation

Definition at line 383 of file NonLinearElasticElement.hpp.

◆ commonData

CommonData& NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::commonData

Structure keeping data abut this particular element, e.g. gradient of deformation at integration points

Definition at line 375 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::dAta

Structure keeping data about problem, like material parameters

Definition at line 373 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fieldDisp

true if field of displacements is given, usually spatial positions are given.

Definition at line 384 of file NonLinearElasticElement.hpp.

◆ fUnction

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fUnction

if true stress i calculated

Definition at line 382 of file NonLinearElasticElement.hpp.

◆ i

FTensor::Index<'i', 3> NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::i

Definition at line 410 of file NonLinearElasticElement.hpp.

◆ j

FTensor::Index<'j', 3> NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::j

Definition at line 411 of file NonLinearElasticElement.hpp.

◆ jAcobian

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::jAcobian

if true Jacobian is calculated

Definition at line 381 of file NonLinearElasticElement.hpp.

◆ k

FTensor::Index<'k', 3> NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::k

Definition at line 412 of file NonLinearElasticElement.hpp.

◆ nbActiveVariables

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::nbActiveVariables

Definition at line 405 of file NonLinearElasticElement.hpp.

◆ ptrh

std::vector<MatrixDouble>* NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::ptrh

Definition at line 407 of file NonLinearElasticElement.hpp.

◆ ptrH

std::vector<MatrixDouble>* NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::ptrH

Definition at line 408 of file NonLinearElasticElement.hpp.

◆ tAg

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::tAg

Definition at line 378 of file NonLinearElasticElement.hpp.


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