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

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

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

Public Member Functions

 OpConstrainsRhs (MoFEM::Interface &m_field, int tag, BlockData &block_data, CommonData &common_data, const std::string &lambda_field_name)
 
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

MoFEM::InterfacemField
 
const std::string lambdaFieldName
 
AuxFunctions< doubleauxFun
 
VectorDouble9 nF
 
- 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

Definition at line 993 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ OpConstrainsRhs()

FractureMechanics::GriffithForceElement::OpConstrainsRhs::OpConstrainsRhs ( MoFEM::Interface m_field,
int  tag,
BlockData block_data,
CommonData common_data,
const std::string &  lambda_field_name 
)
inline

Definition at line 1000 of file GriffithForceElement.hpp.

1004  "MESH_NODE_POSITIONS", UserDataOperator::OPROW),
1005  AuxOp(tag, block_data, common_data), mField(m_field),
1006  lambdaFieldName(lambda_field_name) {}

Member Function Documentation

◆ doWork()

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

Definition at line 1012 of file GriffithForceElement.hpp.

1013  {
1015  if (type != MBVERTEX)
1017 
1018  // get indices of material DOFs
1019  CHKERR setIndices(data);
1022 
1023  auxFun.currentCoords.resize(9, false);
1024  noalias(auxFun.currentCoords) = data.getFieldData();
1025 
1026  CHKERR auxFun.calculateGriffithForce(1, getGaussPts()(2, 0) * 0.5,
1027  data.getDiffN(0));
1028 
1029  // assemble local vector
1030  nF.resize(9, false);
1031  nF.clear();
1032  for (int nn = 0; nn != 3; nn++) {
1033  if (rowLambdaIndices[nn] != -1) {
1034  const auto lambda = lambdaAtNodes[nn];
1035  for (int dd = 0; dd != 3; dd++) {
1036  int idx = 3 * nn + dd;
1037  nF[idx] -= lambda * auxFun.griffithForce[idx];
1038  }
1039  }
1040  }
1041 
1042  // assemble the right hand vectors
1043  CHKERR VecSetValues(getFEMethod()->snes_f, 9, &*rowIndices.data().begin(),
1044  &nF[0], ADD_VALUES);
1045 
1047  }

Member Data Documentation

◆ auxFun

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

Definition at line 1008 of file GriffithForceElement.hpp.

◆ lambdaFieldName

const std::string FractureMechanics::GriffithForceElement::OpConstrainsRhs::lambdaFieldName

Definition at line 998 of file GriffithForceElement.hpp.

◆ mField

MoFEM::Interface& FractureMechanics::GriffithForceElement::OpConstrainsRhs::mField

Definition at line 997 of file GriffithForceElement.hpp.

◆ nF

VectorDouble9 FractureMechanics::GriffithForceElement::OpConstrainsRhs::nF

Definition at line 1009 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::OpConstrainsRhs::mField
MoFEM::Interface & mField
Definition: GriffithForceElement.hpp:997
FractureMechanics::GriffithForceElement::AuxOp::setIndices
MoFEMErrorCode setIndices(DataForcesAndSourcesCore::EntData &data)
Definition: GriffithForceElement.hpp:232
FractureMechanics::GriffithForceElement::OpConstrainsRhs::nF
VectorDouble9 nF
Definition: GriffithForceElement.hpp:1009
FractureMechanics::GriffithForceElement::AuxOp::setLambdaNodes
MoFEMErrorCode setLambdaNodes(FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, const std::string &lambda_field_name)
Definition: GriffithForceElement.hpp:276
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1589
FractureMechanics::GriffithForceElement::AuxFunctions::currentCoords
ublas::vector< TYPE > currentCoords
Definition: GriffithForceElement.hpp:126
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
convert.type
type
Definition: convert.py:64
FractureMechanics::GriffithForceElement::AuxOp::setLambdaIndices
MoFEMErrorCode setLambdaIndices(FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, const std::string &lambda_field_name)
Definition: GriffithForceElement.hpp:301
FractureMechanics::GriffithForceElement::AuxOp::rowLambdaIndices
VectorInt3 rowLambdaIndices
Definition: GriffithForceElement.hpp:227
FractureMechanics::GriffithForceElement::OpConstrainsRhs::lambdaFieldName
const std::string lambdaFieldName
Definition: GriffithForceElement.hpp:998
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::AuxOp::lambdaAtNodes
VectorDouble3 lambdaAtNodes
Definition: GriffithForceElement.hpp:229
FractureMechanics::GriffithForceElement::AuxFunctions::griffithForce
ublas::vector< TYPE > griffithForce
Definition: GriffithForceElement.hpp:127
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
FractureMechanics::GriffithForceElement::AuxOp::rowIndices
ublas::vector< int > rowIndices
Definition: GriffithForceElement.hpp:222
FractureMechanics::GriffithForceElement::OpConstrainsRhs::auxFun
AuxFunctions< double > auxFun
Definition: GriffithForceElement.hpp:1008
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
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