v0.9.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)
 
MoFEMErrorCode getMatParameters (double *_lambda, double *_mu, int *_block_id)
 get material parameter More...
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Here real work is done. 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

MoFEM::InterfacemField
 
moab::Interface & postProcMesh
 
std::vector< EntityHandle > & mapGaussPts
 
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
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 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 52 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 
)

Constructor

Definition at line 70 of file PostProcHookStresses.hpp.

81  field_name, ForcesAndSourcesCore::UserDataOperator::OPROW),
82  mField(m_field), postProcMesh(post_proc_mesh),
83  mapGaussPts(map_gauss_pts),
84 #ifdef __NONLINEAR_ELASTIC_HPP
85  setOfBlocksMaterialDataPtr(set_of_block_data_ptr),
86 #endif //__NONLINEAR_ELASTIC_HPP
87  commonData(common_data) {
88  }
PostProcVolumeOnRefinedMesh::CommonData & commonData
moab::Interface & postProcMesh
std::vector< EntityHandle > & mapGaussPts
MoFEM::Interface & mField
VolumeElementForcesAndSourcesCoreBase::UserDataOperator UserDataOperator

Member Function Documentation

◆ doWork()

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

Here real work is done.

Definition at line 152 of file PostProcHookStresses.hpp.

153  {
155 
156  if (type != MBVERTEX)
158  if (data.getFieldData().size() == 0)
160 
161  // const MoFEM::FEDofEntity *dof_ptr = data.getFieldDofs()[0];
162 
163  int id;
164  double lambda, mu;
165  CHKERR getMatParameters(&lambda, &mu, &id);
166 
167  MatrixDouble D_lambda, D_mu, D;
168  D_lambda.resize(6, 6);
169  D_lambda.clear();
170  for (int rr = 0; rr < 3; rr++) {
171  for (int cc = 0; cc < 3; cc++) {
172  D_lambda(rr, cc) = 1;
173  }
174  }
175  D_mu.resize(6, 6);
176  D_mu.clear();
177  for (int rr = 0; rr < 6; rr++) {
178  D_mu(rr, rr) = rr < 3 ? 2 : 1;
179  }
180  D = lambda * D_lambda + mu * D_mu;
181 
182  int tag_length = 9;
183  double def_VAL[tag_length];
184  bzero(def_VAL, tag_length * sizeof(double));
185  Tag th_stress;
186  CHKERR postProcMesh.tag_get_handle("STRESS", 9, MB_TYPE_DOUBLE, th_stress,
187  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
188 
189  Tag th_prin_stress_vect1, th_prin_stress_vect2, th_prin_stress_vect3;
190  CHKERR postProcMesh.tag_get_handle("S1", 3, MB_TYPE_DOUBLE,
191  th_prin_stress_vect1,
192  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
193  CHKERR postProcMesh.tag_get_handle("S2", 3, MB_TYPE_DOUBLE,
194  th_prin_stress_vect2,
195  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
196  CHKERR postProcMesh.tag_get_handle("S3", 3, MB_TYPE_DOUBLE,
197  th_prin_stress_vect3,
198  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
199 
200  Tag th_prin_stress_vals;
201  CHKERR postProcMesh.tag_get_handle("PRINCIPAL_STRESS", 3, MB_TYPE_DOUBLE,
202  th_prin_stress_vals,
203  MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
204 
205  Tag th_id;
206  int def_block_id = -1;
207  CHKERR postProcMesh.tag_get_handle("BLOCK_ID", 1, MB_TYPE_INTEGER, th_id,
208  MB_TAG_CREAT | MB_TAG_SPARSE,
209  &def_block_id);
210  Range::iterator tit = commonData.tEts.begin();
211  for (; tit != commonData.tEts.end(); tit++) {
212  CHKERR postProcMesh.tag_set_data(th_id, &*tit, 1, &id);
213  }
214 
215  VectorDouble strain;
216  VectorDouble stress;
217  MatrixDouble Stress;
218 
219  // Combine eigenvalues and vectors to create principal stress vector
220  MatrixDouble prin_stress_vect(3, 3);
221  VectorDouble prin_vals_vect(3);
222 
223  int nb_gauss_pts = data.getN().size1();
224  if (mapGaussPts.size() != (unsigned int)nb_gauss_pts) {
225  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
226  }
227  for (int gg = 0; gg < nb_gauss_pts; gg++) {
228 
229  strain.resize(6);
230  strain[0] = (commonData.gradMap[rowFieldName][gg])(0, 0);
231  strain[1] = (commonData.gradMap[rowFieldName][gg])(1, 1);
232  strain[2] = (commonData.gradMap[rowFieldName][gg])(2, 2);
233  strain[3] = (commonData.gradMap[rowFieldName][gg])(0, 1) +
234  (commonData.gradMap[rowFieldName][gg])(1, 0);
235  strain[4] = (commonData.gradMap[rowFieldName][gg])(1, 2) +
236  (commonData.gradMap[rowFieldName][gg])(2, 1);
237  strain[5] = (commonData.gradMap[rowFieldName][gg])(0, 2) +
238  (commonData.gradMap[rowFieldName][gg])(2, 0);
239 
240  stress.resize(6);
241  noalias(stress) = prod(D, strain);
242 
243  Stress.resize(3, 3);
244  Stress(0, 0) = stress[0];
245  Stress(1, 1) = stress[1];
246  Stress(2, 2) = stress[2];
247  Stress(0, 1) = Stress(1, 0) = stress[3];
248  Stress(1, 2) = Stress(2, 1) = stress[4];
249  Stress(2, 0) = Stress(0, 2) = stress[5];
250 
251  CHKERR postProcMesh.tag_set_data(th_stress, &mapGaussPts[gg], 1,
252  &Stress(0, 0));
253 
254  MatrixDouble eigen_vectors = Stress;
255  VectorDouble eigen_values(3);
256 
257  // LAPACK - eigenvalues and vectors. Applied twice for initial creates
258  // memory space
259  int n = 3, lda = 3, info, lwork = -1;
260  double wkopt;
261  info = lapack_dsyev('V', 'U', n, &(eigen_vectors.data()[0]), lda,
262  &(eigen_values.data()[0]), &wkopt, lwork);
263  if (info != 0)
264  SETERRQ1(PETSC_COMM_SELF, 1,
265  "is something wrong with lapack_dsyev info = %d", info);
266  lwork = (int)wkopt;
267  double work[lwork];
268  info = lapack_dsyev('V', 'U', n, &(eigen_vectors.data()[0]), lda,
269  &(eigen_values.data()[0]), work, lwork);
270  if (info != 0)
271  SETERRQ1(PETSC_COMM_SELF, 1,
272  "is something wrong with lapack_dsyev info = %d", info);
273 
274  map<double, int> eigen_sort;
275  eigen_sort[eigen_values[0]] = 0;
276  eigen_sort[eigen_values[1]] = 1;
277  eigen_sort[eigen_values[2]] = 2;
278 
279  prin_stress_vect.clear();
280  prin_vals_vect.clear();
281 
282  int ii = 0;
283  for (map<double, int>::reverse_iterator mit = eigen_sort.rbegin();
284  mit != eigen_sort.rend(); mit++) {
285  prin_vals_vect[ii] = eigen_values[mit->second];
286  for (int dd = 0; dd != 3; dd++) {
287  prin_stress_vect(ii, dd) = eigen_vectors.data()[3 * mit->second + dd];
288  }
289  ii++;
290  }
291 
292  // Tag principle stress vectors 1, 2, 3
293  CHKERR postProcMesh.tag_set_data(th_prin_stress_vect1, &mapGaussPts[gg],
294  1, &prin_stress_vect(0, 0));
295  CHKERR postProcMesh.tag_set_data(th_prin_stress_vect2, &mapGaussPts[gg],
296  1, &prin_stress_vect(1, 0));
297  CHKERR postProcMesh.tag_set_data(th_prin_stress_vect3, &mapGaussPts[gg],
298  1, &prin_stress_vect(2, 0));
299  CHKERR postProcMesh.tag_set_data(th_prin_stress_vals, &mapGaussPts[gg], 1,
300  &prin_vals_vect[0]);
301  }
302 
304  }
PostProcVolumeOnRefinedMesh::CommonData & commonData
moab::Interface & postProcMesh
std::map< std::string, std::vector< MatrixDouble > > gradMap
Range tEts
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:74
std::vector< EntityHandle > & mapGaussPts
#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
MoFEMErrorCode getMatParameters(double *_lambda, double *_mu, int *_block_id)
get material parameter
static __CLPK_integer lapack_dsyev(char jobz, char uplo, __CLPK_integer n, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_doublereal *w, __CLPK_doublereal *work, __CLPK_integer lwork)
Definition: lapack_wrap.h:273
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
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
ublas::vector< double, DoubleAllocator > VectorDouble
Definition: Types.hpp:72

◆ getMatParameters()

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

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 103 of file PostProcHookStresses.hpp.

104  {
106 
107  *_lambda = 1;
108  *_mu = 1;
109 
112  mField, BLOCKSET | MAT_ELASTICSET, it)) {
113  Mat_Elastic mydata;
114  CHKERR it->getAttributeDataStructure(mydata);
115 
116  Range meshsets;
117  CHKERR mField.get_moab().get_entities_by_type(it->meshset, MBENTITYSET,
118  meshsets, false);
119  meshsets.insert(it->meshset);
120  for (Range::iterator mit = meshsets.begin(); mit != meshsets.end();
121  mit++) {
122  if (mField.get_moab().contains_entities(*mit, &ent, 1)) {
123  *_lambda = LAMBDA(mydata.data.Young, mydata.data.Poisson);
124  *_mu = MU(mydata.data.Young, mydata.data.Poisson);
125  *_block_id = it->getMeshsetId();
126 #ifdef __NONLINEAR_ELASTIC_HPP
127  if (setOfBlocksMaterialDataPtr) {
128  *_lambda =
129  LAMBDA(setOfBlocksMaterialDataPtr->at(*_block_id).E,
130  setOfBlocksMaterialDataPtr->at(*_block_id).PoissonRatio);
131  *_mu = MU(setOfBlocksMaterialDataPtr->at(*_block_id).E,
132  setOfBlocksMaterialDataPtr->at(*_block_id).PoissonRatio);
133  }
134 #endif //__NONLINEAR_ELASTIC_HPP
136  }
137  }
138  }
139 
140  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
141  "Element is not in elastic block, however you run linear elastic "
142  "analysis with that element\n"
143  "top tip: check if you update block sets after mesh refinements or "
144  "interface insertion");
145 
147  }
#define LAMBDA(E, NU)
Definition: fem_tools.h:32
virtual moab::Interface & get_moab()=0
#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.
MoFEM::Interface & mField
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MU(E, NU)
Definition: fem_tools.h:33
block name is "MAT_ELASTIC"
Definition: definitions.h:222
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

Member Data Documentation

◆ commonData

PostProcVolumeOnRefinedMesh::CommonData& PostProcHookStress::commonData

Definition at line 65 of file PostProcHookStresses.hpp.

◆ mapGaussPts

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

Definition at line 57 of file PostProcHookStresses.hpp.

◆ mField

MoFEM::Interface& PostProcHookStress::mField

Definition at line 55 of file PostProcHookStresses.hpp.

◆ postProcMesh

moab::Interface& PostProcHookStress::postProcMesh

Definition at line 56 of file PostProcHookStresses.hpp.


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