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

Operator post-procesing stresses for Hook isotropic material. More...

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

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

Public Member Functions

 PostProcHookStress (MoFEM::Interface &m_field, moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, const std::string field_name, PostProcVolumeOnRefinedMesh::CommonData &common_data, const bool is_field_disp=true)
 
MoFEMErrorCode getMatParameters (double *_lambda, double *_mu, int *_block_id)
 get material parameter More...
 
MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Here real work is done. More...
 
- 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

MoFEM::InterfacemField
 
moab::Interface & postProcMesh
 
std::vector< EntityHandle > & mapGaussPts
 
bool isFieldDisp
 
PostProcVolumeOnRefinedMesh::CommonDatacommonData
 
- 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

Operator post-procesing stresses for Hook isotropic material.

Example how to use it

PostProcVolumeOnRefinedMesh post_proc(m_field);
{
CHKERR post_proc.generateReferenceElementMesh();
CHKERR post_proc.addFieldValuesPostProc("DISPLACEMENT");
CHKERR post_proc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
CHKERR post_proc.addFieldValuesGradientPostProc("DISPLACEMENT");
//add postprocessing for stresses
post_proc.getOpPtrVector().push_back(
m_field,
post_proc.postProcMesh,
post_proc.mapGaussPts,
"DISPLACEMENT",
post_proc.commonData,
&elastic.setOfBlocks
)
);
CHKERR DMoFEMLoopFiniteElements(dm,"ELASTIC",&post_proc);
CHKERR post_proc.writeFile("out.h5m");
}
Examples
cell_forces.cpp, and elasticity.cpp.

Definition at line 40 of file PostProcHookStresses.hpp.

Constructor & Destructor Documentation

◆ PostProcHookStress()

PostProcHookStress::PostProcHookStress ( MoFEM::Interface m_field,
moab::Interface &  post_proc_mesh,
std::vector< EntityHandle > &  map_gauss_pts,
const std::string  field_name,
PostProcVolumeOnRefinedMesh::CommonData common_data,
const bool  is_field_disp = true 
)
inline

Constructor

Definition at line 59 of file PostProcHookStresses.hpp.

69  field_name, ForcesAndSourcesCore::UserDataOperator::OPROW),
70  mField(m_field), postProcMesh(post_proc_mesh),
71  mapGaussPts(map_gauss_pts),
72 #ifdef __NONLINEAR_ELASTIC_HPP
73  setOfBlocksMaterialDataPtr(set_of_block_data_ptr),
74 #endif //__NONLINEAR_ELASTIC_HPP
75  commonData(common_data), isFieldDisp(is_field_disp) {
76  }

Member Function Documentation

◆ doWork()

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

Here real work is done.

Definition at line 140 of file PostProcHookStresses.hpp.

141  {
143 
144  if (type != MBVERTEX)
146  if (data.getFieldData().size() == 0)
148 
149  int id;
150  double lambda, mu;
151  CHKERR getMatParameters(&lambda, &mu, &id);
152 
153  MatrixDouble D_lambda, D_mu, D;
154  D_lambda.resize(6, 6);
155  D_lambda.clear();
156  for (int rr = 0; rr < 3; rr++) {
157  for (int cc = 0; cc < 3; cc++) {
158  D_lambda(rr, cc) = 1;
159  }
160  }
161  D_mu.resize(6, 6);
162  D_mu.clear();
163  for (int rr = 0; rr < 6; rr++) {
164  D_mu(rr, rr) = rr < 3 ? 2 : 1;
165  }
166  D = lambda * D_lambda + mu * D_mu;
167 
168  int tag_length = 9;
169  double def_VAL[tag_length];
170  bzero(def_VAL, tag_length * sizeof(double));
171  Tag th_stress;
172  CHKERR postProcMesh.tag_get_handle("STRESS", 9, MB_TYPE_DOUBLE, th_stress,
173  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
174 
175  Tag th_id;
176  int def_block_id = -1;
177  CHKERR postProcMesh.tag_get_handle("BLOCK_ID", 1, MB_TYPE_INTEGER, th_id,
178  MB_TAG_CREAT | MB_TAG_SPARSE,
179  &def_block_id);
180  Range::iterator tit = commonData.tEts.begin();
181  for (; tit != commonData.tEts.end(); tit++)
182  CHKERR postProcMesh.tag_set_data(th_id, &*tit, 1, &id);
183 
184  VectorDouble strain;
185  VectorDouble stress;
186  MatrixDouble Stress;
187 
188  int nb_gauss_pts = data.getN().size1();
189  if (mapGaussPts.size() != (unsigned int)nb_gauss_pts)
190  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
191 
192  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
193 
194  strain.resize(6);
195  strain[0] = (commonData.gradMap[rowFieldName][gg])(0, 0);
196  strain[1] = (commonData.gradMap[rowFieldName][gg])(1, 1);
197  strain[2] = (commonData.gradMap[rowFieldName][gg])(2, 2);
198  strain[3] = (commonData.gradMap[rowFieldName][gg])(0, 1) +
199  (commonData.gradMap[rowFieldName][gg])(1, 0);
200  strain[4] = (commonData.gradMap[rowFieldName][gg])(1, 2) +
201  (commonData.gradMap[rowFieldName][gg])(2, 1);
202  strain[5] = (commonData.gradMap[rowFieldName][gg])(0, 2) +
203  (commonData.gradMap[rowFieldName][gg])(2, 0);
204 
205  if (!isFieldDisp) {
206  strain[0] -= 1.0;
207  strain[1] -= 1.0;
208  strain[2] -= 1.0;
209  }
210 
211  stress.resize(6);
212  noalias(stress) = prod(D, strain);
213 
214  Stress.resize(3, 3);
215  Stress(0, 0) = stress[0];
216  Stress(1, 1) = stress[1];
217  Stress(2, 2) = stress[2];
218  Stress(0, 1) = Stress(1, 0) = stress[3];
219  Stress(1, 2) = Stress(2, 1) = stress[4];
220  Stress(2, 0) = Stress(0, 2) = stress[5];
221 
222  CHKERR postProcMesh.tag_set_data(th_stress, &mapGaussPts[gg], 1,
223  &Stress(0, 0));
224  }
225 
227  }

◆ getMatParameters()

MoFEMErrorCode PostProcHookStress::getMatParameters ( double _lambda,
double _mu,
int *  _block_id 
)
inline

get material parameter

Material parameters are read form BlockSet, however if block data are present, use data how are set for elastic element operators.

Parameters
_lambdaelastic material constant
_muelastic material constant
_block_idblock id
Returns
error code

Definition at line 91 of file PostProcHookStresses.hpp.

92  {
94 
95  *_lambda = 1;
96  *_mu = 1;
97 
100  mField, BLOCKSET | MAT_ELASTICSET, it)) {
101  Mat_Elastic mydata;
102  CHKERR it->getAttributeDataStructure(mydata);
103 
104  Range meshsets;
105  CHKERR mField.get_moab().get_entities_by_type(it->meshset, MBENTITYSET,
106  meshsets, false);
107  meshsets.insert(it->meshset);
108  for (Range::iterator mit = meshsets.begin(); mit != meshsets.end();
109  mit++) {
110  if (mField.get_moab().contains_entities(*mit, &ent, 1)) {
111  *_lambda = LAMBDA(mydata.data.Young, mydata.data.Poisson);
112  *_mu = MU(mydata.data.Young, mydata.data.Poisson);
113  *_block_id = it->getMeshsetId();
114 #ifdef __NONLINEAR_ELASTIC_HPP
115  if (setOfBlocksMaterialDataPtr) {
116  *_lambda =
117  LAMBDA(setOfBlocksMaterialDataPtr->at(*_block_id).E,
118  setOfBlocksMaterialDataPtr->at(*_block_id).PoissonRatio);
119  *_mu = MU(setOfBlocksMaterialDataPtr->at(*_block_id).E,
120  setOfBlocksMaterialDataPtr->at(*_block_id).PoissonRatio);
121  }
122 #endif //__NONLINEAR_ELASTIC_HPP
124  }
125  }
126  }
127 
128  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
129  "Element is not in elastic block, however you run linear elastic "
130  "analysis with that element\n"
131  "top tip: check if you update block sets after mesh refinements or "
132  "interface insertion");
133 
135  }

Member Data Documentation

◆ commonData

PostProcVolumeOnRefinedMesh::CommonData& PostProcHookStress::commonData

Definition at line 54 of file PostProcHookStresses.hpp.

◆ isFieldDisp

bool PostProcHookStress::isFieldDisp

Definition at line 46 of file PostProcHookStresses.hpp.

◆ mapGaussPts

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

Definition at line 45 of file PostProcHookStresses.hpp.

◆ mField

MoFEM::Interface& PostProcHookStress::mField

Definition at line 43 of file PostProcHookStresses.hpp.

◆ postProcMesh

moab::Interface& PostProcHookStress::postProcMesh

Definition at line 44 of file PostProcHookStresses.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:447
EntityHandle
PostProcHookStress::mapGaussPts
std::vector< EntityHandle > & mapGaussPts
Definition: PostProcHookStresses.hpp:45
PostProcHookStress::isFieldDisp
bool isFieldDisp
Definition: PostProcHookStresses.hpp:46
_IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
Definition: MeshsetsManager.hpp:49
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
LAMBDA
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
PostProcHookStress::postProcMesh
moab::Interface & postProcMesh
Definition: PostProcHookStresses.hpp:44
PostProcHookStress::getMatParameters
MoFEMErrorCode getMatParameters(double *_lambda, double *_mu, int *_block_id)
get material parameter
Definition: PostProcHookStresses.hpp:91
MoFEM::ForcesAndSourcesCore::UserDataOperator::rowFieldName
std::string rowFieldName
Definition: ForcesAndSourcesCore.hpp:577
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
PostProcHookStress::PostProcHookStress
PostProcHookStress(MoFEM::Interface &m_field, moab::Interface &post_proc_mesh, std::vector< EntityHandle > &map_gauss_pts, const std::string field_name, PostProcVolumeOnRefinedMesh::CommonData &common_data, const bool is_field_disp=true)
Definition: PostProcHookStresses.hpp:59
convert.type
type
Definition: convert.py:64
MAT_ELASTICSET
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
Definition: definitions.h:159
MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: VolumeElementForcesAndSourcesCore.hpp:105
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Definition: ForcesAndSourcesCore.hpp:999
PostProcHookStress::commonData
PostProcVolumeOnRefinedMesh::CommonData & commonData
Definition: PostProcHookStresses.hpp:54
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
PostProcCommonOnRefMesh::CommonDataForVolume::tEts
Range tEts
Definition: PostProcOnRefMesh.hpp:38
PostProcHookStress::mField
MoFEM::Interface & mField
Definition: PostProcHookStresses.hpp:43
Range
PostProcVolumeOnRefinedMesh
Post processing.
Definition: PostProcOnRefMesh.hpp:955
BLOCKSET
@ BLOCKSET
Definition: definitions.h:148
mu
double mu
Definition: dynamic_first_order_con_law.cpp:98
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::Types::VectorDouble
UBlasVector< double > VectorDouble
Definition: Types.hpp:68
ReactionDiffusionEquation::D
const double D
diffusivity
Definition: reaction_diffusion.cpp:20
PostProcCommonOnRefMesh::CommonData::gradMap
std::map< std::string, std::vector< MatrixDouble > > gradMap
Definition: PostProcOnRefMesh.hpp:34
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:590
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MU
#define MU(E, NU)
Definition: fem_tools.h:23