v0.14.0
Public Member Functions | Public Attributes | List of all members
FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags Struct Reference

#include <users_modules/fracture_mechanics/src/GriffithForceElement.hpp>

Collaboration diagram for FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags:
[legend]

Public Member Functions

 SaveGriffithForceOnTags (MoFEM::Interface &m_field, Tag th, Range &tris, Tag *th_coords=nullptr)
 
MoFEMErrorCode operator() ()
 

Public Attributes

MoFEM::InterfacemField
 
Tag tH
 
RangetRis
 
Tag * thCoordsPtr
 
AuxFunctions< doubleauxFun
 

Detailed Description

Definition at line 532 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ SaveGriffithForceOnTags()

FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::SaveGriffithForceOnTags ( MoFEM::Interface m_field,
Tag  th,
Range tris,
Tag *  th_coords = nullptr 
)
inline

Definition at line 540 of file GriffithForceElement.hpp.

542  : mField(m_field), tH(th), tRis(tris), thCoordsPtr(th_coords) {}

Member Function Documentation

◆ operator()()

MoFEMErrorCode FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::operator() ( )
inline

Definition at line 546 of file GriffithForceElement.hpp.

546  {
548 
549  Range verts;
550  CHKERR mField.get_moab().get_connectivity(tRis, verts, true);
551  {
552  std::vector<double> data(verts.size() * 3, 0);
553  CHKERR mField.get_moab().tag_set_data(tH, verts, &*data.begin());
554  }
555  MatrixDouble diff_n(3, 2);
556  CHKERR ShapeDiffMBTRI(&*diff_n.data().begin());
557 
558  for (Range::iterator tit = tRis.begin(); tit != tRis.end(); ++tit) {
559  int num_nodes;
560  const EntityHandle *conn;
561  CHKERR mField.get_moab().get_connectivity(*tit, conn, num_nodes, true);
562  VectorDouble9 coords(9);
563 
564  int nb_dofs = num_nodes * 3;
565  auxFun.referenceCoords.resize(9, false);
566  auxFun.currentCoords.resize(9, false);
567 
568  if (thCoordsPtr) {
569  CHKERR mField.get_moab().tag_get_data(*thCoordsPtr, conn, num_nodes,
570  &*coords.begin());
571  } else {
572  CHKERR mField.get_moab().get_coords(conn, num_nodes,
573  &*coords.begin());
574  }
575  for (int dd = 0; dd != nb_dofs; dd++) {
576  auxFun.currentCoords[dd] = coords[dd];
577  }
578  if (thCoordsPtr)
579  CHKERR mField.get_moab().get_coords(conn, num_nodes,
580  &*coords.begin());
581  for (int dd = 0; dd != nb_dofs; dd++) {
582  auxFun.referenceCoords[dd] = coords[dd];
583  }
584 
585  CHKERR auxFun.calculateGriffithForce(1, 0.5, diff_n);
586 
587  double *tag_vals[3];
588  CHKERR mField.get_moab().tag_get_by_ptr(tH, conn, 3,
589  (const void **)tag_vals);
590  for (int nn = 0; nn != 3; ++nn) {
591  for (int dd = 0; dd != 3; ++dd) {
592  tag_vals[nn][dd] += auxFun.griffithForce[3 * nn + dd];
593  }
594  }
595  }
596 
598  }

Member Data Documentation

◆ auxFun

AuxFunctions<double> FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::auxFun

Definition at line 544 of file GriffithForceElement.hpp.

◆ mField

MoFEM::Interface& FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::mField

Definition at line 534 of file GriffithForceElement.hpp.

◆ tH

Tag FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::tH

Definition at line 536 of file GriffithForceElement.hpp.

◆ thCoordsPtr

Tag* FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::thCoordsPtr

Definition at line 538 of file GriffithForceElement.hpp.

◆ tRis

Range& FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::tRis

Definition at line 537 of file GriffithForceElement.hpp.


The documentation for this struct was generated from the following file:
FractureMechanics::GriffithForceElement::AuxFunctions::calculateGriffithForce
MoFEMErrorCode calculateGriffithForce(const double gc, const double beta, const MatrixDouble &diff_n)
Definition: GriffithForceElement.hpp:201
EntityHandle
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::tRis
Range & tRis
Definition: GriffithForceElement.hpp:537
FractureMechanics::GriffithForceElement::AuxFunctions::currentCoords
ublas::vector< TYPE > currentCoords
Definition: GriffithForceElement.hpp:126
FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::auxFun
AuxFunctions< double > auxFun
Definition: GriffithForceElement.hpp:544
ShapeDiffMBTRI
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:194
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
FractureMechanics::GriffithForceElement::AuxFunctions::referenceCoords
ublas::vector< TYPE > referenceCoords
Definition: GriffithForceElement.hpp:125
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::mField
MoFEM::Interface & mField
Definition: GriffithForceElement.hpp:534
Range
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
FractureMechanics::GriffithForceElement::AuxFunctions::griffithForce
ublas::vector< TYPE > griffithForce
Definition: GriffithForceElement.hpp:127
MoFEM::Types::VectorDouble9
VectorBoundedArray< double, 9 > VectorDouble9
Definition: Types.hpp:96
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
FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::thCoordsPtr
Tag * thCoordsPtr
Definition: GriffithForceElement.hpp:538
FractureMechanics::GriffithForceElement::SaveGriffithForceOnTags::tH
Tag tH
Definition: GriffithForceElement.hpp:536