v0.9.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, DataForcesAndSourcesCore::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...
 
doublegetVolume ()
 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...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHoGaussPtsJac ()
 
MatrixDoublegetHoGaussPtsInvJac ()
 
VectorDoublegetHoGaussPtsDetJac ()
 
auto getFTenosr0HoMeasure ()
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 Get coordinates at integration points taking geometry from field. More...
 
auto getFTensor2HoGaussPtsJac ()
 
auto getFTensor2HoGaussPtsInvJac ()
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
MoFEMErrorCode getDivergenceOfHDivBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
 Get divergence of base functions at integration point. More...
 
MoFEMErrorCode getCurlOfHCurlBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
 Get curl of base functions at integration point. More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, 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...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (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 ()
 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...
 
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
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
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 

Additional Inherited Members

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

Detailed Description

Operator performs automatic differentiation.

Examples
cell_forces.cpp.

Definition at line 386 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 246 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 NonlinearElasticElement::OpJacobianEshelbyStress, and Smoother::OpJacobianSmoother.

Definition at line 258 of file NonLinearElasticElement.cpp.

259  {
261 
262  ierr = dAta.materialAdoublePtr->calculateP_PiolaKirchhoffI(
264  CHKERRG(ierr);
265  if (aLe) {
267  dAta.materialAdoublePtr->detH *
268  prod(dAta.materialAdoublePtr->P, trans(dAta.materialAdoublePtr->invH));
269  }
270  commonData.sTress[gg].resize(3, 3, false);
271  for (int dd1 = 0; dd1 < 3; dd1++) {
272  for (int dd2 = 0; dd2 < 3; dd2++) {
273  dAta.materialAdoublePtr->P(dd1, dd2) >>=
274  (commonData.sTress[gg])(dd1, dd2);
275  }
276  }
277 
279 }
boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble > > materialAdoublePtr
std::vector< MatrixDouble3by3 > sTress
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:544
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
bool aLe
true if arbitrary Lagrangian-Eulerian formulation
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ doWork()

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

Calculate stress or jacobian at gauss points.

Parameters
row_side
row_type
row_data
Returns
error code

Definition at line 379 of file NonLinearElasticElement.cpp.

381  {
383 
384  // do it only once, no need to repeat this for edges,faces or tets
385  if (row_type != MBVERTEX)
387 
388  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
389  dAta.tEts.end()) {
391  }
392 
393  int nb_dofs = row_data.getFieldData().size();
394  if (nb_dofs == 0)
396  dAta.materialAdoublePtr->commonDataPtr = &commonData;
397  dAta.materialAdoublePtr->opPtr = this;
398 
399  int nb_gauss_pts = row_data.getN().size1();
400  commonData.sTress.resize(nb_gauss_pts);
401  commonData.jacStress.resize(nb_gauss_pts);
402 
404  if (aLe) {
406  }
407 
408  for (int gg = 0; gg != nb_gauss_pts; gg++) {
409 
410  dAta.materialAdoublePtr->gG = gg;
411 
412  // Record tag and calculate stress
414  CHKERR recordTag(gg);
415  }
416 
417  // Set active variables vector
418  if (jAcobian || (!recordTagForIntegrationPoint(gg))) {
419  activeVariables.resize(nbActiveVariables, false);
420  if (!aLe) {
421  for (int dd1 = 0; dd1 < 3; dd1++) {
422  for (int dd2 = 0; dd2 < 3; dd2++) {
423  activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
424  }
425  }
426  } else {
427  for (int dd1 = 0; dd1 < 3; dd1++) {
428  for (int dd2 = 0; dd2 < 3; dd2++) {
429  activeVariables(dd1 * 3 + dd2) = (*ptrh)[gg](dd1, dd2);
430  }
431  }
432  for (int dd1 = 0; dd1 < 3; dd1++) {
433  for (int dd2 = 0; dd2 < 3; dd2++) {
434  activeVariables(9 + dd1 * 3 + dd2) = (*ptrH)[gg](dd1, dd2);
435  }
436  }
437  }
438  CHKERR dAta.materialAdoublePtr->setUserActiveVariables(activeVariables);
439 
440  // Play tag and calculate stress or tangent
441  if (jAcobian || (!recordTagForIntegrationPoint(gg))) {
442  CHKERR playTag(gg);
443  }
444  }
445  }
446 
448 }
boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble > > materialAdoublePtr
std::vector< MatrixDouble3by3 > sTress
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Range tEts
constrains elements in block set
virtual bool recordTagForIntegrationPoint(const int gg)
Cgeck if tape is recorded for given integration point.
std::vector< MatrixDouble > jacStress
this is simply material tangent operator
std::map< std::string, std::vector< MatrixDouble > > gradAtGaussPts
#define CHKERR
Inline error check.
Definition: definitions.h:596
bool aLe
true if arbitrary Lagrangian-Eulerian formulation
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MoFEMErrorCode playTag(const int gg)
Play ADOL-C tape.
virtual MoFEMErrorCode recordTag(const int gg)
Record ADOL-C tape.
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

◆ playTag()

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

Play ADOL-C tape.

Returns
error code

Definition at line 344 of file NonLinearElasticElement.cpp.

344  {
346 
347  int r;
348 
349  if (fUnction) {
350  commonData.sTress[gg].resize(3, 3, false);
351  // play recorder for values
352  r = ::function(tAg, 9, nbActiveVariables, &activeVariables[0],
353  &commonData.sTress[gg](0, 0));
354  if (r < adlocReturnValue) { // function is locally analytic
355  SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
356  "ADOL-C function evaluation with error r = %d", r);
357  }
358  }
359 
360  if (jAcobian) {
361  commonData.jacStress[gg].resize(9, nbActiveVariables, false);
362  double *jac_ptr[] = {
363  &(commonData.jacStress[gg](0, 0)), &(commonData.jacStress[gg](1, 0)),
364  &(commonData.jacStress[gg](2, 0)), &(commonData.jacStress[gg](3, 0)),
365  &(commonData.jacStress[gg](4, 0)), &(commonData.jacStress[gg](5, 0)),
366  &(commonData.jacStress[gg](6, 0)), &(commonData.jacStress[gg](7, 0)),
367  &(commonData.jacStress[gg](8, 0))};
368  // play recorder for jacobians
369  r = jacobian(tAg, 9, nbActiveVariables, &activeVariables[0], jac_ptr);
370  if (r < adlocReturnValue) {
371  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
372  "ADOL-C function evaluation with error");
373  }
374  }
375 
377 }
std::vector< MatrixDouble3by3 > sTress
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
std::vector< MatrixDouble > jacStress
this is simply material tangent operator

◆ recordTag()

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

Record ADOL-C tape.

Returns
error code

Definition at line 282 of file NonLinearElasticElement.cpp.

283  {
285 
286  trace_on(tAg, 0);
287 
288  dAta.materialAdoublePtr->F.resize(3, 3, false);
289 
290  if (!aLe) {
291 
292  nbActiveVariables = 0;
293  for (int dd1 = 0; dd1 < 3; dd1++) {
294  for (int dd2 = 0; dd2 < 3; dd2++) {
295  dAta.materialAdoublePtr->F(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
296  if (fieldDisp) {
297  if (dd1 == dd2) {
298  dAta.materialAdoublePtr->F(dd1, dd2) += 1;
299  }
300  }
302  }
303  }
304 
305  } else {
306 
307  nbActiveVariables = 0;
308 
309  dAta.materialAdoublePtr->h.resize(3, 3, false);
310  for (int dd1 = 0; dd1 < 3; dd1++) {
311  for (int dd2 = 0; dd2 < 3; dd2++) {
312  dAta.materialAdoublePtr->h(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
314  }
315  }
316 
317  dAta.materialAdoublePtr->H.resize(3, 3, false);
318  for (int dd1 = 0; dd1 < 3; dd1++) {
319  for (int dd2 = 0; dd2 < 3; dd2++) {
320  dAta.materialAdoublePtr->H(dd1, dd2) <<= (*ptrH)[gg](dd1, dd2);
322  }
323  }
324 
326  dAta.materialAdoublePtr->detH);
327  dAta.materialAdoublePtr->invH.resize(3, 3, false);
330  dAta.materialAdoublePtr->invH);
331  noalias(dAta.materialAdoublePtr->F) =
333  }
334 
335  CHKERR dAta.materialAdoublePtr->setUserActiveVariables(nbActiveVariables);
337 
338  trace_off();
339 
341 }
boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< adouble > > materialAdoublePtr
virtual MoFEMErrorCode calculateStress(const int gg)
Calculate Paola-Kirchhoff I stress.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERR
Inline error check.
Definition: definitions.h:596
bool aLe
true if arbitrary Lagrangian-Eulerian formulation
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ 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 449 of file NonLinearElasticElement.hpp.

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

Member Data Documentation

◆ activeVariables

VectorDouble NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::activeVariables

Definition at line 420 of file NonLinearElasticElement.hpp.

◆ adlocReturnValue

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::adlocReturnValue

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

Definition at line 395 of file NonLinearElasticElement.hpp.

◆ aLe

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::aLe

true if arbitrary Lagrangian-Eulerian formulation

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

◆ dAta

BlockData& NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::dAta

Structure keeping data about problem, like material parameters

Definition at line 389 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fieldDisp

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

Definition at line 400 of file NonLinearElasticElement.hpp.

◆ fUnction

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fUnction

if true stress i calculated

Definition at line 398 of file NonLinearElasticElement.hpp.

◆ jAcobian

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::jAcobian

if true Jacobian is calculated

Definition at line 397 of file NonLinearElasticElement.hpp.

◆ nbActiveVariables

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::nbActiveVariables

Definition at line 421 of file NonLinearElasticElement.hpp.

◆ ptrh

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

Definition at line 423 of file NonLinearElasticElement.hpp.

◆ ptrH

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

Definition at line 424 of file NonLinearElasticElement.hpp.

◆ tAg

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::tAg

Definition at line 394 of file NonLinearElasticElement.hpp.


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