v0.13.0
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. More...
 
virtual MoFEMErrorCode calculateStress (const int gg)
 Calculate Paola-Kirchhoff I stress. More...
 
virtual MoFEMErrorCode recordTag (const int gg)
 Record ADOL-C tape. More...
 
virtual MoFEMErrorCode playTag (const int gg)
 Play ADOL-C tape. More...
 
virtual bool recordTagForIntegrationPoint (const int gg)
 Cgeck if tape is recorded for given integration point. More...
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 Calculate stress or jacobian at gauss points. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
double & getVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
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 getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
double & getMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
- 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. More...
 
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. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Public Attributes

BlockDatadAta
 
CommonDatacommonData
 
int tAg
 
int adlocReturnValue
 
bool jAcobian
 if true Jacobian is calculated More...
 
bool fUnction
 if true stress i calculated More...
 
bool aLe
 true if arbitrary Lagrangian-Eulerian formulation More...
 
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. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 }
 Controls loop over entities on element. More...
 
- 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::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Operator performs automatic differentiation.

Examples
cell_forces.cpp.

Definition at line 383 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 212 of file NonLinearElasticElement.cpp.

Member Function Documentation

◆ calculateStress()

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

Calculate Paola-Kirchhoff I stress.

Returns
error code

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

Definition at line 225 of file NonLinearElasticElement.cpp.

226  {
228 
229  CHKERR dAta.materialAdoublePtr->calculateP_PiolaKirchhoffI(
231 
232  if (aLe) {
233  auto &t_P = dAta.materialAdoublePtr->t_P;
234  auto &t_invH = dAta.materialAdoublePtr->t_invH;
235  t_P(i, j) = t_P(i, k) * t_invH(j, k);
236  t_P(i, j) *= dAta.materialAdoublePtr->detH;
237  }
238 
239  commonData.sTress[gg].resize(3, 3, false);
240  for (int dd1 = 0; dd1 < 3; dd1++) {
241  for (int dd2 = 0; dd2 < 3; dd2++) {
242  dAta.materialAdoublePtr->P(dd1, dd2) >>=
243  (commonData.sTress[gg])(dd1, dd2);
244  }
245  }
246 
248 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble > > materialAdoublePtr
std::vector< MatrixDouble3by3 > sTress

◆ 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 352 of file NonLinearElasticElement.cpp.

354  {
356 
357  // do it only once, no need to repeat this for edges,faces or tets
358  if (row_type != MBVERTEX)
360 
361  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
362  dAta.tEts.end()) {
364  }
365 
366  int nb_dofs = row_data.getFieldData().size();
367  if (nb_dofs == 0)
369  dAta.materialAdoublePtr->commonDataPtr = &commonData;
370  dAta.materialAdoublePtr->opPtr = this;
371 
372  int nb_gauss_pts = row_data.getN().size1();
373  commonData.sTress.resize(nb_gauss_pts);
374  commonData.jacStress.resize(nb_gauss_pts);
375 
377  if (aLe) {
379  }
380 
381  for (int gg = 0; gg != nb_gauss_pts; gg++) {
382 
383  dAta.materialAdoublePtr->gG = gg;
384 
385  // Record tag and calculate stress
387  CHKERR recordTag(gg);
388  }
389 
390  // Set active variables vector
391  if (jAcobian || (!recordTagForIntegrationPoint(gg))) {
392  activeVariables.resize(nbActiveVariables, false);
393  if (!aLe) {
394  for (int dd1 = 0; dd1 < 3; dd1++) {
395  for (int dd2 = 0; dd2 < 3; dd2++) {
396  activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
397  }
398  }
399  } else {
400  for (int dd1 = 0; dd1 < 3; dd1++) {
401  for (int dd2 = 0; dd2 < 3; dd2++) {
402  activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
403  }
404  }
405  for (int dd1 = 0; dd1 < 3; dd1++) {
406  for (int dd2 = 0; dd2 < 3; dd2++) {
407  activeVariables(9 + dd1 * 3 + dd2) = (*ptrH)[gg](dd1, dd2);
408  }
409  }
410  }
411  CHKERR dAta.materialAdoublePtr->setUserActiveVariables(activeVariables);
412 
413  // Play tag and calculate stress or tangent
414  if (jAcobian || (!recordTagForIntegrationPoint(gg))) {
415  CHKERR playTag(gg);
416  }
417  }
418  }
419 
421 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
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 dofs values
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 317 of file NonLinearElasticElement.cpp.

317  {
319 
320  int r;
321 
322  if (fUnction) {
323  commonData.sTress[gg].resize(3, 3, false);
324  // play recorder for values
325  r = ::function(tAg, 9, nbActiveVariables, &activeVariables[0],
326  &commonData.sTress[gg](0, 0));
327  if (r < adlocReturnValue) { // function is locally analytic
328  SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
329  "ADOL-C function evaluation with error r = %d", r);
330  }
331  }
332 
333  if (jAcobian) {
334  commonData.jacStress[gg].resize(9, nbActiveVariables, false);
335  double *jac_ptr[] = {
336  &(commonData.jacStress[gg](0, 0)), &(commonData.jacStress[gg](1, 0)),
337  &(commonData.jacStress[gg](2, 0)), &(commonData.jacStress[gg](3, 0)),
338  &(commonData.jacStress[gg](4, 0)), &(commonData.jacStress[gg](5, 0)),
339  &(commonData.jacStress[gg](6, 0)), &(commonData.jacStress[gg](7, 0)),
340  &(commonData.jacStress[gg](8, 0))};
341  // play recorder for jacobians
342  r = jacobian(tAg, 9, nbActiveVariables, &activeVariables[0], jac_ptr);
343  if (r < adlocReturnValue) {
344  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
345  "ADOL-C function evaluation with error");
346  }
347  }
348 
350 }
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:47
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
const double r
rate factor

◆ recordTag()

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

Record ADOL-C tape.

Returns
error code

Definition at line 251 of file NonLinearElasticElement.cpp.

252  {
254 
255  trace_on(tAg, 0);
256 
257  dAta.materialAdoublePtr->F.resize(3, 3, false);
258 
259  if (!aLe) {
260 
261  nbActiveVariables = 0;
262  for (int dd1 = 0; dd1 < 3; dd1++) {
263  for (int dd2 = 0; dd2 < 3; dd2++) {
264  dAta.materialAdoublePtr->F(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
265  if (fieldDisp) {
266  if (dd1 == dd2) {
267  dAta.materialAdoublePtr->F(dd1, dd2) += 1;
268  }
269  }
271  }
272  }
273 
274  } else {
275 
276  nbActiveVariables = 0;
277 
278  dAta.materialAdoublePtr->h.resize(3, 3, false);
279  for (int dd1 = 0; dd1 < 3; dd1++) {
280  for (int dd2 = 0; dd2 < 3; dd2++) {
281  dAta.materialAdoublePtr->h(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
283  }
284  }
285 
286  dAta.materialAdoublePtr->H.resize(3, 3, false);
287  for (int dd1 = 0; dd1 < 3; dd1++) {
288  for (int dd2 = 0; dd2 < 3; dd2++) {
289  dAta.materialAdoublePtr->H(dd1, dd2) <<= (*ptrH)[gg](dd1, dd2);
291  }
292  }
293 
295  dAta.materialAdoublePtr->invH.resize(3, 3, false);
297  dAta.materialAdoublePtr->detH,
298  dAta.materialAdoublePtr->invH);
299 
300  auto &t_F = dAta.materialAdoublePtr->t_F;
301  auto &t_h = dAta.materialAdoublePtr->t_h;
302  auto &t_invH = dAta.materialAdoublePtr->t_invH;
303 
304  t_F(i, j) = t_h(i, k) * t_invH(k, j);
305 
306  }
307 
308  CHKERR dAta.materialAdoublePtr->setUserActiveVariables(nbActiveVariables);
310 
311  trace_off();
312 
314 }
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.
Definition: Templates.hpp:1172
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1161
virtual MoFEMErrorCode calculateStress(const int gg)
Calculate Paola-Kirchhoff I stress.

◆ recordTagForIntegrationPoint()

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

Cgeck if tape is recorded for given integration point.

Parameters
ggintegration point
Returns
true if tag is recorded

Definition at line 450 of file NonLinearElasticElement.hpp.

450  {
451  // return true;
452  if (gg == 0)
453  return true;
454  return false;
455  }

Member Data Documentation

◆ activeVariables

VectorDouble NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::activeVariables

Definition at line 417 of file NonLinearElasticElement.hpp.

◆ adlocReturnValue

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::adlocReturnValue

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

Definition at line 392 of file NonLinearElasticElement.hpp.

◆ aLe

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::aLe

true if arbitrary Lagrangian-Eulerian formulation

Definition at line 396 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 388 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::dAta

Structure keeping data about problem, like material parameters

Definition at line 386 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fieldDisp

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

Definition at line 397 of file NonLinearElasticElement.hpp.

◆ fUnction

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fUnction

if true stress i calculated

Definition at line 395 of file NonLinearElasticElement.hpp.

◆ i

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

Definition at line 423 of file NonLinearElasticElement.hpp.

◆ j

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

Definition at line 424 of file NonLinearElasticElement.hpp.

◆ jAcobian

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::jAcobian

if true Jacobian is calculated

Definition at line 394 of file NonLinearElasticElement.hpp.

◆ k

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

Definition at line 425 of file NonLinearElasticElement.hpp.

◆ nbActiveVariables

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::nbActiveVariables

Definition at line 418 of file NonLinearElasticElement.hpp.

◆ ptrh

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

Definition at line 420 of file NonLinearElasticElement.hpp.

◆ ptrH

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

Definition at line 421 of file NonLinearElasticElement.hpp.

◆ tAg

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::tAg

Definition at line 391 of file NonLinearElasticElement.hpp.


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