v0.9.1
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...
 
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
 
DEPRECATED Vec getSnesF () const
 
DEPRECATED Vec getSnesX () const
 
DEPRECATED Mat getSnesA () const
 
DEPRECATED 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
 
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)
 
virtual ~DataOperator ()=default
 
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)
 
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 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
 
- 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...
 
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, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
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 387 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 226 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 238 of file NonLinearElasticElement.cpp.

239  {
241 
242  ierr = dAta.materialAdoublePtr->calculateP_PiolaKirchhoffI(
244  CHKERRG(ierr);
245  if (aLe) {
247  dAta.materialAdoublePtr->detH *
248  prod(dAta.materialAdoublePtr->P, trans(dAta.materialAdoublePtr->invH));
249  }
250  commonData.sTress[gg].resize(3, 3, false);
251  for (int dd1 = 0; dd1 < 3; dd1++) {
252  for (int dd2 = 0; dd2 < 3; dd2++) {
253  dAta.materialAdoublePtr->P(dd1, dd2) >>=
254  (commonData.sTress[gg])(dd1, dd2);
255  }
256  }
257 
259 }
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:483
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:550
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:413

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

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

324  {
326 
327  int r;
328 
329  if (fUnction) {
330  commonData.sTress[gg].resize(3, 3, false);
331  // play recorder for values
332  r = ::function(tAg, 9, nbActiveVariables, &activeVariables[0],
333  &commonData.sTress[gg](0, 0));
334  if (r < adlocReturnValue) { // function is locally analytic
335  SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
336  "ADOL-C function evaluation with error r = %d", r);
337  }
338  }
339 
340  if (jAcobian) {
341  commonData.jacStress[gg].resize(9, nbActiveVariables, false);
342  double *jac_ptr[] = {
343  &(commonData.jacStress[gg](0, 0)), &(commonData.jacStress[gg](1, 0)),
344  &(commonData.jacStress[gg](2, 0)), &(commonData.jacStress[gg](3, 0)),
345  &(commonData.jacStress[gg](4, 0)), &(commonData.jacStress[gg](5, 0)),
346  &(commonData.jacStress[gg](6, 0)), &(commonData.jacStress[gg](7, 0)),
347  &(commonData.jacStress[gg](8, 0))};
348  // play recorder for jacobians
349  r = jacobian(tAg, 9, nbActiveVariables, &activeVariables[0], jac_ptr);
350  if (r < adlocReturnValue) {
351  SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
352  "ADOL-C function evaluation with error");
353  }
354  }
355 
357 }
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:507
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
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 262 of file NonLinearElasticElement.cpp.

263  {
265 
266  trace_on(tAg, 0);
267 
268  dAta.materialAdoublePtr->F.resize(3, 3, false);
269 
270  if (!aLe) {
271 
272  nbActiveVariables = 0;
273  for (int dd1 = 0; dd1 < 3; dd1++) {
274  for (int dd2 = 0; dd2 < 3; dd2++) {
275  dAta.materialAdoublePtr->F(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
276  if (fieldDisp) {
277  if (dd1 == dd2) {
278  dAta.materialAdoublePtr->F(dd1, dd2) += 1;
279  }
280  }
282  }
283  }
284 
285  } else {
286 
287  nbActiveVariables = 0;
288 
289  dAta.materialAdoublePtr->h.resize(3, 3, false);
290  for (int dd1 = 0; dd1 < 3; dd1++) {
291  for (int dd2 = 0; dd2 < 3; dd2++) {
292  dAta.materialAdoublePtr->h(dd1, dd2) <<= (*ptrh)[gg](dd1, dd2);
294  }
295  }
296 
297  dAta.materialAdoublePtr->H.resize(3, 3, false);
298  for (int dd1 = 0; dd1 < 3; dd1++) {
299  for (int dd2 = 0; dd2 < 3; dd2++) {
300  dAta.materialAdoublePtr->H(dd1, dd2) <<= (*ptrH)[gg](dd1, dd2);
302  }
303  }
304 
306  dAta.materialAdoublePtr->detH);
307  dAta.materialAdoublePtr->invH.resize(3, 3, false);
310  dAta.materialAdoublePtr->invH);
311  noalias(dAta.materialAdoublePtr->F) =
313  }
314 
315  CHKERR dAta.materialAdoublePtr->setUserActiveVariables(nbActiveVariables);
317 
318  trace_off();
319 
321 }
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:483
#define CHKERR
Inline error check.
Definition: definitions.h:602
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:413

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

◆ adlocReturnValue

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::adlocReturnValue

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

Definition at line 396 of file NonLinearElasticElement.hpp.

◆ aLe

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::aLe

true if arbitrary Lagrangian-Eulerian formulation

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

◆ dAta

BlockData& NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::dAta

Structure keeping data about problem, like material parameters

Definition at line 390 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fieldDisp

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

Definition at line 401 of file NonLinearElasticElement.hpp.

◆ fUnction

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::fUnction

if true stress i calculated

Definition at line 399 of file NonLinearElasticElement.hpp.

◆ jAcobian

bool NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::jAcobian

if true Jacobian is calculated

Definition at line 398 of file NonLinearElasticElement.hpp.

◆ nbActiveVariables

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::nbActiveVariables

Definition at line 422 of file NonLinearElasticElement.hpp.

◆ ptrh

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

Definition at line 424 of file NonLinearElasticElement.hpp.

◆ ptrH

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

Definition at line 425 of file NonLinearElasticElement.hpp.

◆ tAg

int NonlinearElasticElement::OpJacobianPiolaKirchhoffStress::tAg

Definition at line 395 of file NonLinearElasticElement.hpp.


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