v0.14.0
Public Member Functions | Public Attributes | List of all members
PostProcStress Struct Reference

#include <users_modules/basic_finite_elements/src/PostProcStresses.hpp>

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

Public Member Functions

 PostProcStress (moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, const std::string field_name, NonlinearElasticElement::BlockData &data, PostProcCommonOnRefMesh::CommonDataForVolume &common_data, const bool field_disp=false, const bool replace_nonanumber_by_max_value=false, const double max_val=1e16, const bool print_cauchy_stress=false)
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::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...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () 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...
 
std::string getFEName () const
 Get name of the element. More...
 
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 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...
 
doublegetMeasure ()
 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, 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. More...
 
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. More...
 
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. More...
 
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. 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

moab::Interface & postProcMesh
 
std::vector< EntityHandle > & mapGaussPts
 
NonlinearElasticElement::BlockDatadAta
 
PostProcCommonOnRefMesh::CommonDataForVolumecommonData
 
const bool fieldDisp
 
const bool replaceNonANumberByMaxValue
 
const double maxVal
 
const bool printCauchy
 
NonlinearElasticElement::CommonData nonLinearElementCommonData
 
- 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,
  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

Examples
nonlinear_dynamics.cpp, NonlinearElasticElementInterface.hpp, Remodeling.cpp, and simple_contact.cpp.

Definition at line 17 of file PostProcStresses.hpp.

Constructor & Destructor Documentation

◆ PostProcStress()

PostProcStress::PostProcStress ( moab::Interface &  post_proc_mesh,
std::vector< EntityHandle > &  map_gauss_pts,
const std::string  field_name,
NonlinearElasticElement::BlockData data,
PostProcCommonOnRefMesh::CommonDataForVolume common_data,
const bool  field_disp = false,
const bool  replace_nonanumber_by_max_value = false,
const double  max_val = 1e16,
const bool  print_cauchy_stress = false 
)
inline

Definition at line 30 of file PostProcStresses.hpp.

40  field_name, ForcesAndSourcesCore::UserDataOperator::OPROW),
41  postProcMesh(post_proc_mesh), mapGaussPts(map_gauss_pts), dAta(data),
42  commonData(common_data), fieldDisp(field_disp),
43  replaceNonANumberByMaxValue(replace_nonanumber_by_max_value),
44  maxVal(max_val), printCauchy(print_cauchy_stress) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode PostProcStress::doWork ( int  side,
EntityType  type,
EntitiesFieldData::EntData data 
)
inline

Definition at line 48 of file PostProcStresses.hpp.

49  {
51 
52  if (type != MBVERTEX)
54  if (data.getIndices().size() == 0)
56  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
57  dAta.tEts.end()) {
59  }
60 
61  const auto &dof_ptr = data.getFieldDofs()[0];
62 
63  int id = dAta.iD;
64 
65  Tag th_id;
66  int def_block_id = -1;
67  CHKERR postProcMesh.tag_get_handle("BLOCK_ID", 1, MB_TYPE_INTEGER, th_id,
68  MB_TAG_CREAT | MB_TAG_SPARSE,
69  &def_block_id);
70  Range::iterator tit = commonData.tEts.begin();
71  for (; tit != commonData.tEts.end(); tit++) {
72  CHKERR postProcMesh.tag_set_data(th_id, &*tit, 1, &id);
73  }
74 
75  string tag_name_piola1 = dof_ptr->getName() + "_PIOLA1_STRESS";
76  string tag_name_energy = dof_ptr->getName() + "_ENERGY_DENSITY";
77 
78  int tag_length = 9;
79  double def_VAL[tag_length];
80  bzero(def_VAL, tag_length * sizeof(double));
81  Tag th_piola1, th_energy, th_cauchy;
82  CHKERR postProcMesh.tag_get_handle(tag_name_piola1.c_str(), tag_length,
83  MB_TYPE_DOUBLE, th_piola1,
84  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
85  CHKERR postProcMesh.tag_get_handle(tag_name_energy.c_str(), 1,
86  MB_TYPE_DOUBLE, th_energy,
87  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
88 
89  if (printCauchy) {
90  string tag_name_cauchy = "MED_" + dof_ptr->getName() + "_CAUCHY_STRESS";
91  CHKERR postProcMesh.tag_get_handle(tag_name_cauchy.c_str(), tag_length,
92  MB_TYPE_DOUBLE, th_cauchy,
93  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
94  }
95 
96  int nb_gauss_pts = data.getN().size1();
97  if (mapGaussPts.size() != (unsigned int)nb_gauss_pts) {
98  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
99  "Nb. of integration points is not equal to number points on "
100  "post-processing mesh");
101  }
102  if (commonData.gradMap[rowFieldName].size() != (unsigned int)nb_gauss_pts) {
103  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
104  "Gradient of field not found, filed <%s> not found",
105  rowFieldName.c_str());
106  }
107 
108  MatrixDouble3by3 H, invH;
109  double detH;
110 
112  dAta.materialDoublePtr->opPtr = this;
113  CHKERR dAta.materialDoublePtr->getDataOnPostProcessor(commonData.fieldMap,
115 
118 
119  MatrixDouble3by3 maxP(3, 3);
120  maxP.clear();
121 
122  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
123 
124  dAta.materialDoublePtr->gG = gg;
125  dAta.materialDoublePtr->F.resize(3, 3);
126  noalias(dAta.materialDoublePtr->F) =
128  if (fieldDisp) {
129  for (int dd = 0; dd != 3; dd++) {
130  dAta.materialDoublePtr->F(dd, dd) += 1;
131  }
132  }
133  if (commonData.gradMap["MESH_NODE_POSITIONS"].size() ==
134  (unsigned int)nb_gauss_pts) {
135  H.resize(3, 3);
136  invH.resize(3, 3);
137  noalias(H) = (commonData.gradMap["MESH_NODE_POSITIONS"])[gg];
138  detH = determinantTensor3by3(H);
139  CHKERR invertTensor3by3(H, detH, invH);
140  noalias(dAta.materialDoublePtr->F) =
141  prod(dAta.materialDoublePtr->F, invH);
142  }
143 
144  int nb_active_variables = 9;
145  CHKERR dAta.materialDoublePtr->setUserActiveVariables(
146  nb_active_variables);
147  CHKERR dAta.materialDoublePtr->calculateP_PiolaKirchhoffI(
149  CHKERR dAta.materialDoublePtr->calculateElasticEnergy(
151  CHKERR postProcMesh.tag_set_data(th_piola1, &mapGaussPts[gg], 1,
152  &dAta.materialDoublePtr->P(0, 0));
153  CHKERR postProcMesh.tag_set_data(th_energy, &mapGaussPts[gg], 1,
154  &dAta.materialDoublePtr->eNergy);
155  if (printCauchy) {
156  dAta.materialDoublePtr->sigmaCauchy.resize(3, 3);
157  CHKERR dAta.materialDoublePtr->calculateCauchyStress(
159  CHKERR postProcMesh.tag_set_data(
160  th_cauchy, &mapGaussPts[gg], 1,
161  &dAta.materialDoublePtr->sigmaCauchy(0, 0));
162  }
163  }
164 
166  MatrixDouble3by3 P(3, 3);
167  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
168  double val_energy;
169  CHKERR postProcMesh.tag_get_data(th_energy, &mapGaussPts[gg], 1,
170  &val_energy);
171  if (!std::isnormal(val_energy)) {
172  CHKERR postProcMesh.tag_set_data(th_energy, &mapGaussPts[gg], 1,
173  &maxVal);
174  CHKERR postProcMesh.tag_get_data(th_piola1, &mapGaussPts[gg], 1,
175  &P(0, 0));
176  for (unsigned int r = 0; r != P.size1(); ++r) {
177  for (unsigned int c = 0; c != P.size2(); ++c) {
178  if (!std::isnormal(P(r, c)))
179  P(r, c) = copysign(maxVal, P(r, c));
180  }
181  }
182  CHKERR postProcMesh.tag_set_data(th_piola1, &mapGaussPts[gg], 1,
183  &P(0, 0));
184  }
185  }
186  }
187 
189  }

Member Data Documentation

◆ commonData

PostProcCommonOnRefMesh::CommonDataForVolume& PostProcStress::commonData

Definition at line 24 of file PostProcStresses.hpp.

◆ dAta

NonlinearElasticElement::BlockData& PostProcStress::dAta

Definition at line 23 of file PostProcStresses.hpp.

◆ fieldDisp

const bool PostProcStress::fieldDisp

Definition at line 25 of file PostProcStresses.hpp.

◆ mapGaussPts

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

Definition at line 21 of file PostProcStresses.hpp.

◆ maxVal

const double PostProcStress::maxVal

Definition at line 27 of file PostProcStresses.hpp.

◆ nonLinearElementCommonData

NonlinearElasticElement::CommonData PostProcStress::nonLinearElementCommonData

Definition at line 46 of file PostProcStresses.hpp.

◆ postProcMesh

moab::Interface& PostProcStress::postProcMesh

Definition at line 20 of file PostProcStresses.hpp.

◆ printCauchy

const bool PostProcStress::printCauchy

Definition at line 28 of file PostProcStresses.hpp.

◆ replaceNonANumberByMaxValue

const bool PostProcStress::replaceNonANumberByMaxValue

Definition at line 26 of file PostProcStresses.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
PostProcCommonOnRefMesh::CommonData::fieldMap
std::map< std::string, std::vector< VectorDouble > > fieldMap
Definition: PostProcOnRefMesh.hpp:33
PostProcStress::commonData
PostProcCommonOnRefMesh::CommonDataForVolume & commonData
Definition: PostProcStresses.hpp:24
MoFEM::ForcesAndSourcesCore::UserDataOperator::rowFieldName
std::string rowFieldName
Definition: ForcesAndSourcesCore.hpp:577
sdf.r
int r
Definition: sdf.py:8
PostProcStress::replaceNonANumberByMaxValue
const bool replaceNonANumberByMaxValue
Definition: PostProcStresses.hpp:26
MoFEM::invertTensor3by3
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:1588
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
PostProcStress::mapGaussPts
std::vector< EntityHandle > & mapGaussPts
Definition: PostProcStresses.hpp:21
PostProcStress::dAta
NonlinearElasticElement::BlockData & dAta
Definition: PostProcStresses.hpp:23
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
convert.type
type
Definition: convert.py:64
EshelbianPlasticity::P
@ P
Definition: EshelbianContact.cpp:197
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:105
PostProcStress::nonLinearElementCommonData
NonlinearElasticElement::CommonData nonLinearElementCommonData
Definition: PostProcStresses.hpp:46
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Definition: ForcesAndSourcesCore.hpp:1000
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
PostProcStress::fieldDisp
const bool fieldDisp
Definition: PostProcStresses.hpp:25
PostProcCommonOnRefMesh::CommonDataForVolume::tEts
Range tEts
Definition: PostProcOnRefMesh.hpp:38
MoFEM::determinantTensor3by3
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
Definition: Templates.hpp:1540
PostProcStress::printCauchy
const bool printCauchy
Definition: PostProcStresses.hpp:28
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
NonlinearElasticElement::BlockData::iD
int iD
Definition: HookeElement.hpp:33
NonlinearElasticElement::BlockData::tEts
Range tEts
constrains elements in block set
Definition: HookeElement.hpp:36
PostProcStress::maxVal
const double maxVal
Definition: PostProcStresses.hpp:27
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::Types::MatrixDouble3by3
MatrixBoundedArray< double, 9 > MatrixDouble3by3
Definition: Types.hpp:105
NonlinearElasticElement::CommonData::dataAtGaussPts
std::map< std::string, std::vector< VectorDouble > > dataAtGaussPts
Definition: NonLinearElasticElement.hpp:107
PostProcCommonOnRefMesh::CommonData::gradMap
std::map< std::string, std::vector< MatrixDouble > > gradMap
Definition: PostProcOnRefMesh.hpp:34
convert.int
int
Definition: convert.py:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
NonlinearElasticElement::CommonData::gradAtGaussPts
std::map< std::string, std::vector< MatrixDouble > > gradAtGaussPts
Definition: NonLinearElasticElement.hpp:108
NonlinearElasticElement::BlockData::materialDoublePtr
boost::shared_ptr< FunctionsToCalculatePiolaKirchhoffI< double > > materialDoublePtr
Definition: NonLinearElasticElement.hpp:94
H
double H
Hardening.
Definition: plastic.cpp:124
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
PostProcStress::postProcMesh
moab::Interface & postProcMesh
Definition: PostProcStresses.hpp:20