v0.9.0
Public Member Functions | Public Attributes | List of all members
OpPostProcStress Struct Reference

#include <users_modules/basic_finite_elements/elasticity_mixed_formulation/src/ElasticityMixedFormulation.hpp>

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

Public Member Functions

 OpPostProcStress (moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, DataAtIntegrationPts &common_data, BlockData &data)
 
PetscErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 
- 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

DataAtIntegrationPtscommonData
 
moab::Interface & postProcMesh
 
std::vector< EntityHandle > & mapGaussPts
 
BlockDatadAta
 
- 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

Examples
elasticity_mixed_formulation.cpp, ElasticityMixedFormulation.hpp, Remodeling.cpp, and Remodeling.hpp.

Definition at line 451 of file ElasticityMixedFormulation.hpp.

Constructor & Destructor Documentation

◆ OpPostProcStress()

OpPostProcStress::OpPostProcStress ( moab::Interface &  post_proc_mesh,
std::vector< EntityHandle > &  map_gauss_pts,
DataAtIntegrationPts common_data,
BlockData data 
)
Examples
ElasticityMixedFormulation.hpp, and Remodeling.cpp.

Definition at line 458 of file ElasticityMixedFormulation.hpp.

462  "U", UserDataOperator::OPROW),
463  commonData(common_data), postProcMesh(post_proc_mesh),
464  mapGaussPts(map_gauss_pts), dAta(data) {
465  doVertices = true;
466  doEdges = false;
467  doQuads = false;
468  doTris = false;
469  doTets = false;
470  doPrisms = false;
471  }
moab::Interface & postProcMesh
bool doVertices
If false skip vertices.
std::vector< EntityHandle > & mapGaussPts
DataAtIntegrationPts & commonData
bool doEdges
If false skip edges.
ForcesAndSourcesCore::UserDataOperator UserDataOperator
BlockData & dAta

Member Function Documentation

◆ doWork()

PetscErrorCode OpPostProcStress::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
Examples
ElasticityMixedFormulation.hpp, and Remodeling.cpp.

Definition at line 473 of file ElasticityMixedFormulation.hpp.

474  {
476  if (type != MBVERTEX)
477  PetscFunctionReturn(9);
478  double def_VAL[9];
479  bzero(def_VAL, 9 * sizeof(double));
480 
481  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
482  dAta.tEts.end()) {
484  }
486 
487  Tag th_stress;
488  CHKERR postProcMesh.tag_get_handle("STRESS", 9, MB_TYPE_DOUBLE, th_stress,
489  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
490  Tag th_strain;
491  CHKERR postProcMesh.tag_get_handle("STRAIN", 9, MB_TYPE_DOUBLE, th_strain,
492  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
493  Tag th_psi;
494  CHKERR postProcMesh.tag_get_handle("ENERGY", 1, MB_TYPE_DOUBLE, th_psi,
495  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
496 
497  auto grad = getFTensor2FromMat<3, 3>(*commonData.gradDispPtr);
499 
500  const int nb_gauss_pts = commonData.gradDispPtr->size2();
501  const double mu = commonData.mU;
502 
507 
508  for (int gg = 0; gg != nb_gauss_pts; gg++) {
509  strain(i, j) = 0.5 * (grad(i, j) + grad(j, i));
510  double psi = 0.5 * p * p + mu * strain(i, j) * strain(i, j);
511 
512  stress(i, j) = 2 * mu * strain(i, j);
513  stress(1, 1) -= p;
514  stress(0, 0) -= p;
515  stress(2, 2) -= p;
516 
517  CHKERR postProcMesh.tag_set_data(th_psi, &mapGaussPts[gg], 1, &psi);
518  CHKERR postProcMesh.tag_set_data(th_strain, &mapGaussPts[gg], 1,
519  &strain(0, 0));
520  CHKERR postProcMesh.tag_set_data(th_stress, &mapGaussPts[gg], 1,
521  &stress(0, 0));
522  ++p;
523  ++grad;
524  }
525 
527  }
static FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:142
moab::Interface & postProcMesh
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< MatrixDouble > gradDispPtr
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
std::vector< EntityHandle > & mapGaussPts
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
boost::shared_ptr< VectorDouble > pPtr
DataAtIntegrationPts & commonData
#define CHKERR
Inline error check.
Definition: definitions.h:596
MoFEMErrorCode getBlockData(BlockData &data)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
BlockData & dAta

Member Data Documentation

◆ commonData

DataAtIntegrationPts& OpPostProcStress::commonData

◆ dAta

BlockData& OpPostProcStress::dAta

◆ mapGaussPts

std::vector<EntityHandle>& OpPostProcStress::mapGaussPts

◆ postProcMesh

moab::Interface& OpPostProcStress::postProcMesh

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