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

assemble internal residual More...

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

Inheritance diagram for FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda:
[legend]
Collaboration diagram for FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda:
[legend]

Public Member Functions

 OpIntLambda (MoFEM::Interface &m_field, GriffithForceElement::BlockData &block_data, GriffithForceElement::CommonData &common_data, const std::string &lambda_field_name, FrontArcLengthControl &arc_front)
 
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
 
std::string lambdaFieldName
 
FrontArcLengthControlarcFront
 
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

assemble internal residual

Parameters
block_datagriffith block data
common_datagriffith common data
lambda_field_namelagrange multipliers controling crack surface area
arc_frontreference to FrontArcLengthControl object

Definition at line 1312 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ OpIntLambda()

FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::OpIntLambda ( MoFEM::Interface m_field,
GriffithForceElement::BlockData block_data,
GriffithForceElement::CommonData common_data,
const std::string &  lambda_field_name,
FrontArcLengthControl arc_front 
)
inline

Definition at line 1319 of file GriffithForceElement.hpp.

1325  "MESH_NODE_POSITIONS", lambda_field_name,
1327  false // not symmetric operator
1328  ),
1329  AuxOp(0, block_data, common_data), mField(m_field),
1330  lambdaFieldName(lambda_field_name), arcFront(arc_front) {}

Member Function Documentation

◆ doWork()

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

Definition at line 1334 of file GriffithForceElement.hpp.

1335  {
1337  if (type != MBVERTEX) {
1339  }
1340  // get indices of Lagrange multiplier DOFs
1342  // get dofs increment
1343  const double *dx;
1344  CHKERR VecGetArrayRead(arcFront.arcPtr->dx, &dx);
1345  // get dofs at begining of load step
1346  const double *x0;
1347  CHKERR VecGetArrayRead(arcFront.arcPtr->x0, &x0);
1348  auxFun.referenceCoords.resize(9, false);
1349  auxFun.currentCoords.resize(9, false);
1350  for (int dd = 0; dd != 9; dd++) {
1351  int loc_idx = data.getLocalIndices()[dd];
1352  if (loc_idx != -1) {
1353  auxFun.referenceCoords[dd] = x0[loc_idx];
1354  auxFun.currentCoords[dd] = auxFun.referenceCoords[dd] + dx[loc_idx];
1355  } else {
1356  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Note expected");
1357  }
1358  }
1359  CHKERR VecRestoreArrayRead(arcFront.arcPtr->dx, &dx);
1360  CHKERR VecRestoreArrayRead(arcFront.arcPtr->x0, &x0);
1361 
1362  // calculate crack area increment
1363  CHKERR auxFun.calculateGriffithForce(1, getGaussPts()(2, 0) * 0.5,
1364  data.getDiffN(0));
1365 
1366  // calculate area change
1367  double delta = 0;
1368  for (int nn = 0; nn != 3; nn++) {
1369  if (rowLambdaIndices[nn] == -1)
1370  continue;
1371  auto griffith_force =
1372  getVectorAdaptor(&(auxFun.griffithForce[3 * nn]), 3);
1373  // griffith_force /= norm_2(griffith_force);
1374  for (int dd = 0; dd != 3; dd++) {
1375  const int idx = 3 * nn + dd;
1376  delta += griffith_force[dd] *
1378  }
1379  }
1380  if (delta != delta) {
1381  SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "nan value");
1382  }
1385  };

Member Data Documentation

◆ arcFront

FrontArcLengthControl& FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::arcFront

Definition at line 1317 of file GriffithForceElement.hpp.

◆ auxFun

AuxFunctions<double> FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::auxFun

Definition at line 1332 of file GriffithForceElement.hpp.

◆ lambdaFieldName

std::string FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::lambdaFieldName

Definition at line 1316 of file GriffithForceElement.hpp.

◆ mField

MoFEM::Interface& FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::mField

Definition at line 1315 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::FrontArcLengthControl::arcPtr
boost::shared_ptr< ArcLengthCtx > arcPtr
Definition: GriffithForceElement.hpp:1291
FractureMechanics::GriffithForceElement::AuxOp::AuxOp
AuxOp(int tag, BlockData &block_data, CommonData &common_data)
Definition: GriffithForceElement.hpp:221
FractureMechanics::GriffithForceElement::AuxFunctions::currentCoords
ublas::vector< TYPE > currentCoords
Definition: GriffithForceElement.hpp:126
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
FractureMechanics::GriffithForceElement::AuxFunctions::referenceCoords
ublas::vector< TYPE > referenceCoords
Definition: GriffithForceElement.hpp:125
delta
static constexpr double delta
Definition: prism_elements_from_surface.cpp:18
MoFEM::getVectorAdaptor
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition: Templates.hpp:31
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::FrontArcLengthControl::OpIntLambda::auxFun
AuxFunctions< double > auxFun
Definition: GriffithForceElement.hpp:1332
FractureMechanics::GriffithForceElement::AuxOp::rowLambdaIndices
VectorInt3 rowLambdaIndices
Definition: GriffithForceElement.hpp:227
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::mField
MoFEM::Interface & mField
Definition: GriffithForceElement.hpp:1315
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
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::arcFront
FrontArcLengthControl & arcFront
Definition: GriffithForceElement.hpp:1317
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaDelta
double & intLambdaDelta
Definition: GriffithForceElement.hpp:1298
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
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::lambdaFieldName
std::string lambdaFieldName
Definition: GriffithForceElement.hpp:1316
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567