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

calculate griffith force vector More...

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

Inheritance diagram for FractureMechanics::GriffithForceElement::OpCalculateGriffithForce:
[legend]
Collaboration diagram for FractureMechanics::GriffithForceElement::OpCalculateGriffithForce:
[legend]

Public Member Functions

 OpCalculateGriffithForce (int tag, BlockData &block_data, CommonData &common_data)
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 
- Public Member Functions inherited from FractureMechanics::GriffithForceElement::AuxOp
 AuxOp (int tag, BlockData &block_data, CommonData &common_data)
 
MoFEMErrorCode setIndices (DataForcesAndSourcesCore::EntData &data)
 
MoFEMErrorCode setVariables (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, DataForcesAndSourcesCore::EntData &data)
 
MoFEMErrorCode setLambdaNodes (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, const std::string &lambda_field_name)
 
MoFEMErrorCode setLambdaIndices (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, const std::string &lambda_field_name)
 

Public Attributes

AuxFunctions< doubleauxFun
 
- Public Attributes inherited from FractureMechanics::GriffithForceElement::AuxOp
int tAg
 
BlockDatablockData
 
CommonDatacommonData
 
ublas::vector< int > rowIndices
 
ublas::vector< int > rowIndicesLocal
 
VectorInt3 rowLambdaIndices
 
VectorInt3 rowLambdaIndicesLocal
 
VectorDouble3 lambdaAtNodes
 
VectorDouble activeVariables
 

Detailed Description

calculate griffith force vector

\[ \mathbf{f}_{\mathrm{m}}^{\mathrm{h}}=\frac{1}{2}\left(\tilde{\mathbf{A}}_{\Gamma}^{\mathrm{h}}\right)^{\mathrm{T}} \mathbf{g}_{c}(\rho(\mathbf{X})) ]\f The matrix $\mathbf{A}_{\Gamma}^h$ comprising of direction vectors along the crack front that are normal to the crack front and tangent to the crack surface is defined as follows: \f[ A^h_{\Gamma} = \| \mathbf{N}(\tilde{\mathbf{X}}) \| = \left\| \epsilon_{ijk} \frac{\partial \Phi^\alpha_p}{\partial \xi_i} \frac{\partial \Phi^\beta_r}{\partial \xi_j} \tilde{X}^\alpha_p \tilde{X}^\beta_r \right\| \]

Definition at line 378 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ OpCalculateGriffithForce()

FractureMechanics::GriffithForceElement::OpCalculateGriffithForce::OpCalculateGriffithForce ( int  tag,
BlockData block_data,
CommonData common_data 
)
inline

Definition at line 381 of file GriffithForceElement.hpp.

384  "MESH_NODE_POSITIONS", UserDataOperator::OPROW),
385  AuxOp(tag, block_data, common_data) {}

Member Function Documentation

◆ doWork()

MoFEMErrorCode FractureMechanics::GriffithForceElement::OpCalculateGriffithForce::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
inline

Definition at line 387 of file GriffithForceElement.hpp.

388  {
390 
391  if (type != MBVERTEX)
393  CHKERR setIndices(data);
394  CHKERR setVariables(this, data);
395  int nb_dofs = data.getIndices().size();
396  int nb_gauss_pts = data.getN().size1();
397 
398  auxFun.referenceCoords.resize(nb_dofs, false);
399  auxFun.currentCoords.resize(nb_dofs, false);
400 
401  for (int dd = 0; dd != nb_dofs; dd++) {
402  auxFun.referenceCoords[dd] = getCoords()[dd];
403  auxFun.currentCoords[dd] = data.getFieldData()[dd];
404  }
405 
406  for (int gg = 0; gg != nb_gauss_pts; gg++) {
407  double val = getGaussPts()(2, gg) * 0.5;
409  data.getDiffN(gg));
410  }
411  commonData.griffithForce.resize(nb_dofs, false);
412  std::fill(commonData.griffithForce.begin(),
413  commonData.griffithForce.end(), 0);
414  std::copy(auxFun.griffithForce.begin(), auxFun.griffithForce.end(),
415  commonData.griffithForce.begin());
416 
418  }

Member Data Documentation

◆ auxFun

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

Definition at line 386 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
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
FractureMechanics::GriffithForceElement::AuxOp::AuxOp
AuxOp(int tag, BlockData &block_data, CommonData &common_data)
Definition: GriffithForceElement.hpp:221
FractureMechanics::GriffithForceElement::AuxOp::setIndices
MoFEMErrorCode setIndices(DataForcesAndSourcesCore::EntData &data)
Definition: GriffithForceElement.hpp:232
FractureMechanics::GriffithForceElement::BlockData::gc
double gc
Griffith energy.
Definition: GriffithForceElement.hpp:75
FractureMechanics::GriffithForceElement::AuxFunctions::currentCoords
ublas::vector< TYPE > currentCoords
Definition: GriffithForceElement.hpp:126
FractureMechanics::GriffithForceElement::AuxOp::setVariables
MoFEMErrorCode setVariables(FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, DataForcesAndSourcesCore::EntData &data)
Definition: GriffithForceElement.hpp:261
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
FractureMechanics::GriffithForceElement::AuxFunctions::referenceCoords
ublas::vector< TYPE > referenceCoords
Definition: GriffithForceElement.hpp:125
convert.type
type
Definition: convert.py:64
FractureMechanics::GriffithForceElement::OpCalculateGriffithForce::auxFun
AuxFunctions< double > auxFun
Definition: GriffithForceElement.hpp:386
FractureMechanics::GriffithForceElement::AuxOp::blockData
BlockData & blockData
Definition: GriffithForceElement.hpp:218
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
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
FractureMechanics::GriffithForceElement::CommonData::griffithForce
VectorDouble griffithForce
Griffith force at nodes.
Definition: GriffithForceElement.hpp:65
FractureMechanics::GriffithForceElement::AuxOp::commonData
CommonData & commonData
Definition: GriffithForceElement.hpp:219
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567