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

assemble internal residual derivative More...

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

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

Public Member Functions

 OpAssemble_db (int tag, 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::FrontArcLengthControl::OpIntLambda
 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

int tAg
 
AuxFunctions< adoublea_auxFun
 
VectorDouble gradDelta
 
VectorDouble activeVariables
 
- Public Attributes inherited from FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda
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 derivative

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 1396 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ OpAssemble_db()

FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::OpAssemble_db ( int  tag,
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 1398 of file GriffithForceElement.hpp.

1403  : OpIntLambda(m_field, block_data, common_data, lambda_field_name,
1404  arc_front),
1405  tAg(tag) {}

Member Function Documentation

◆ doWork()

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

Definition at line 1411 of file GriffithForceElement.hpp.

1412  {
1414  if (type != MBVERTEX) {
1416  }
1417  CHKERR setIndices(data);
1418  // get indices of Lagrange multiplier DOFs
1420  double d_delta;
1421  adouble a_delta;
1422  // get dofs increment
1423  const double *dx;
1424  CHKERR VecGetArrayRead(arcFront.arcPtr->dx, &dx);
1425  // get dofs at begining of load step
1426  const double *x0;
1427  CHKERR VecGetArrayRead(arcFront.arcPtr->x0, &x0);
1428  a_auxFun.referenceCoords.resize(9, false);
1429  a_auxFun.currentCoords.resize(9, false);
1430 
1431  activeVariables.resize(18, false);
1432  trace_on(tAg);
1433  for (int dd = 0; dd != 9; dd++) {
1434  int loc_idx = data.getLocalIndices()[dd];
1435  if (loc_idx != -1) {
1436  a_auxFun.referenceCoords[dd] <<= x0[loc_idx];
1437  activeVariables[dd] = x0[loc_idx];
1438  } else {
1439  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Note expected");
1440  }
1441  }
1442  for (int dd = 0; dd != 9; dd++) {
1443  int loc_idx = data.getLocalIndices()[dd];
1444  if (loc_idx != -1) {
1445  a_auxFun.currentCoords[dd] <<= x0[dd] + dx[loc_idx];
1446  activeVariables[9 + dd] = x0[loc_idx] + dx[loc_idx];
1447  } else {
1448  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Note expected");
1449  }
1450  }
1451  // calculate crack area increment
1452 
1453  CHKERR a_auxFun.calculateGriffithForce(1, getGaussPts()(2, 0) * 0.5,
1454  data.getDiffN(0));
1455  // calculate area change
1456  a_delta = 0;
1457  for (int nn = 0; nn != 3; nn++) {
1458  if (rowLambdaIndices[nn] == -1)
1459  continue;
1460  auto griffith_force =
1461  getVectorAdaptor(&(a_auxFun.griffithForce[3 * nn]), 3);
1462  // griffith_force /= sqrt(griffith_force[0] * griffith_force[0] +
1463  // griffith_force[1] * griffith_force[1] +
1464  // griffith_force[2] * griffith_force[2]);
1465  for (int dd = 0; dd != 3; dd++) {
1466  const int idx = 3 * nn + dd;
1467  a_delta += griffith_force[dd] * (a_auxFun.currentCoords[idx] -
1468  a_auxFun.referenceCoords[idx]);
1469  }
1470  }
1471  a_delta >>= d_delta;
1472  trace_off();
1473  CHKERR VecRestoreArrayRead(arcFront.arcPtr->dx, &dx);
1474  CHKERR VecRestoreArrayRead(arcFront.arcPtr->x0, &x0);
1475 
1476  gradDelta.resize(18, false);
1477  {
1478  int r;
1479  // play recorder for jacobian
1480  r = ::gradient(tAg, 18, &activeVariables[0], &gradDelta[0]);
1481  if (r < 3) { // function is locally analytic
1482  SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1483  "ADOL-C function evaluation with error r = %d", r);
1484  }
1485  }
1486 
1487  for (int dd = 0; dd != 9; dd++) {
1488  double val = arcFront.arcPtr->alpha * gradDelta[9 + dd];
1489  CHKERR VecSetValue(arcFront.arcPtr->db, data.getIndices()[dd], val,
1490  ADD_VALUES);
1491  }
1492 
1494  }

Member Data Documentation

◆ a_auxFun

AuxFunctions<adouble> FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::a_auxFun

Definition at line 1407 of file GriffithForceElement.hpp.

◆ activeVariables

VectorDouble FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::activeVariables

Definition at line 1409 of file GriffithForceElement.hpp.

◆ gradDelta

VectorDouble FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::gradDelta

Definition at line 1408 of file GriffithForceElement.hpp.

◆ tAg

int FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::tAg

Definition at line 1397 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::FrontArcLengthControl::OpIntLambda::OpIntLambda
OpIntLambda(MoFEM::Interface &m_field, GriffithForceElement::BlockData &block_data, GriffithForceElement::CommonData &common_data, const std::string &lambda_field_name, FrontArcLengthControl &arc_front)
Definition: GriffithForceElement.hpp:1319
FractureMechanics::GriffithForceElement::AuxOp::setIndices
MoFEMErrorCode setIndices(DataForcesAndSourcesCore::EntData &data)
Definition: GriffithForceElement.hpp:232
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::activeVariables
VectorDouble activeVariables
Definition: GriffithForceElement.hpp:1409
sdf.r
int r
Definition: sdf.py:8
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::gradDelta
VectorDouble gradDelta
Definition: GriffithForceElement.hpp:1408
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
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::AuxOp::rowLambdaIndices
VectorInt3 rowLambdaIndices
Definition: GriffithForceElement.hpp:227
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::a_auxFun
AuxFunctions< adouble > a_auxFun
Definition: GriffithForceElement.hpp:1407
adouble
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
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::arcFront
FrontArcLengthControl & arcFront
Definition: GriffithForceElement.hpp:1317
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpAssemble_db::tAg
int tAg
Definition: GriffithForceElement.hpp:1397
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
FractureMechanics::GriffithForceElement::FrontArcLengthControl::OpIntLambda::lambdaFieldName
std::string lambdaFieldName
Definition: GriffithForceElement.hpp:1316