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

arc-length element More...

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

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

Classes

struct  OpAssemble_db
 assemble internal residual derivative More...
 
struct  OpIntLambda
 assemble internal residual More...
 

Public Member Functions

 FrontArcLengthControl (int tag, GriffithForceElement::BlockData &block_data, GriffithForceElement::CommonData &common_data, const std::string &lambda_field_name, boost::shared_ptr< ArcLengthCtx > &arc_ptr)
 
virtual ~FrontArcLengthControl ()
 
MoFEMErrorCode preProcess ()
 
double calculateLambdaInt ()
 Calculate internal lambda. More...
 
MoFEMErrorCode operator() ()
 
MoFEMErrorCode postProcess ()
 
MoFEMErrorCode calculateDb ()
 Calculate db. More...
 
MoFEMErrorCode calculateDxAndDlambda (Vec x)
 

Public Attributes

int tAg
 
boost::shared_ptr< ArcLengthCtxarcPtr
 
BlockDatablockData
 
CommonDatacommonData
 
std::string lambdaFieldName
 
Vec ghostIntLambda
 
double intLambdaArray [2]
 
doubleintLambdaDelta
 
doubleintLambdaLambda
 
double intLambda
 
boost::shared_ptr< MyTriangleFEfeIntLambda
 

Detailed Description

arc-length element

Calualte residual for arc length method amd vector db, which is derivative of residual.

Definition at line 1288 of file GriffithForceElement.hpp.

Constructor & Destructor Documentation

◆ FrontArcLengthControl()

FractureMechanics::GriffithForceElement::FrontArcLengthControl::FrontArcLengthControl ( int  tag,
GriffithForceElement::BlockData block_data,
GriffithForceElement::CommonData common_data,
const std::string &  lambda_field_name,
boost::shared_ptr< ArcLengthCtx > &  arc_ptr 
)
inline

Definition at line 1497 of file GriffithForceElement.hpp.

1501  : tAg(tag), arcPtr(arc_ptr), blockData(block_data),
1502  commonData(common_data), lambdaFieldName(lambda_field_name),
1505  int two[] = {0, 1};
1506  int nb_local = arcPtr->mField.get_comm_rank() == 0 ? 2 : 0;
1507  int nb_ghost = nb_local ? 0 : 2;
1508  ierr = VecCreateGhostWithArray(PETSC_COMM_WORLD, nb_local, 2, nb_ghost,
1510  CHKERRABORT(PETSC_COMM_WORLD, ierr);
1511  feIntLambda = boost::shared_ptr<MyTriangleFE>(
1512  new GriffithForceElement::MyTriangleFE(arcPtr->mField));
1513  }

◆ ~FrontArcLengthControl()

virtual FractureMechanics::GriffithForceElement::FrontArcLengthControl::~FrontArcLengthControl ( )
inlinevirtual

Definition at line 1514 of file GriffithForceElement.hpp.

1514  {
1515  ierr = VecDestroy(&ghostIntLambda);
1516  CHKERRABORT(PETSC_COMM_WORLD, ierr);
1517  }

Member Function Documentation

◆ calculateDb()

MoFEMErrorCode FractureMechanics::GriffithForceElement::FrontArcLengthControl::calculateDb ( )
inline

Calculate db.

Definition at line 1588 of file GriffithForceElement.hpp.

1588  {
1590 
1591  MOFEM_LOG_CHANNEL("WORLD");
1593  MOFEM_LOG_TAG("WORLD", "Arc");
1594 
1595  double alpha = arcPtr->alpha;
1596  double beta = arcPtr->beta;
1597  double d_lambda = arcPtr->dLambda;
1598 
1599  CHKERR VecZeroEntries(ghostIntLambda);
1600  CHKERR VecGhostUpdateBegin(ghostIntLambda, INSERT_VALUES,
1601  SCATTER_FORWARD);
1602  CHKERR VecGhostUpdateEnd(ghostIntLambda, INSERT_VALUES, SCATTER_FORWARD);
1603  CHKERR VecZeroEntries(arcPtr->db);
1604  CHKERR VecGhostUpdateBegin(arcPtr->db, INSERT_VALUES, SCATTER_FORWARD);
1605  CHKERR VecGhostUpdateEnd(arcPtr->db, INSERT_VALUES, SCATTER_FORWARD);
1606 
1607  // calculate internal lambda
1608  feIntLambda->getOpPtrVector().clear();
1609  feIntLambda->getOpPtrVector().push_back(new OpIntLambda(
1610  arcPtr->mField, blockData, commonData, lambdaFieldName, *this));
1611  feIntLambda->getOpPtrVector().push_back(new OpAssemble_db(
1612  tAg, arcPtr->mField, blockData, commonData, lambdaFieldName, *this));
1613  CHKERR arcPtr->mField.loop_finite_elements(
1614  problemPtr->getName(), "CRACKFRONT_AREA_ELEM", (*feIntLambda));
1615  const double *dx;
1616  const auto bit_field_number = getFieldBitNumber(lambdaFieldName);
1617  CHKERR VecGetArrayRead(arcPtr->dx, &dx);
1619  bit_field_number, dit)) {
1620  if (static_cast<int>(dit->get()->getPart()) !=
1621  arcPtr->mField.get_comm_rank())
1622  continue;
1623  int idx = dit->get()->getPetscLocalDofIdx();
1624  double val = dx[idx];
1625  CHKERR VecSetValue(ghostIntLambda, 1, val, ADD_VALUES);
1626  }
1627  CHKERR VecRestoreArrayRead(arcPtr->dx, &dx);
1628 
1629  CHKERR VecAssemblyBegin(ghostIntLambda);
1630  CHKERR VecAssemblyEnd(ghostIntLambda);
1631  CHKERR VecGhostUpdateBegin(ghostIntLambda, ADD_VALUES, SCATTER_REVERSE);
1632  CHKERR VecGhostUpdateEnd(ghostIntLambda, ADD_VALUES, SCATTER_REVERSE);
1633  CHKERR VecGhostUpdateBegin(ghostIntLambda, INSERT_VALUES,
1634  SCATTER_FORWARD);
1635  CHKERR VecGhostUpdateEnd(ghostIntLambda, INSERT_VALUES, SCATTER_FORWARD);
1636 
1637  CHKERR VecAssemblyBegin(arcPtr->db);
1638  CHKERR VecAssemblyEnd(arcPtr->db);
1639  CHKERR VecGhostUpdateBegin(arcPtr->db, ADD_VALUES, SCATTER_REVERSE);
1640  CHKERR VecGhostUpdateEnd(arcPtr->db, ADD_VALUES, SCATTER_REVERSE);
1641  CHKERR VecGhostUpdateBegin(arcPtr->db, INSERT_VALUES, SCATTER_FORWARD);
1642  CHKERR VecGhostUpdateEnd(arcPtr->db, INSERT_VALUES, SCATTER_FORWARD);
1643 
1644  intLambda = alpha * intLambdaDelta +
1645  alpha * beta * beta * d_lambda * d_lambda * arcPtr->F_lambda2;
1646 
1647  MOFEM_LOG_C(
1648  "WORLD", Sev::noisy,
1649  "\tintLambdaLambda = %9.8e intLambdaDelta = %9.8e intLambda = %9.8e",
1651 
1652  double nrm;
1653  CHKERR VecNorm(arcPtr->db, NORM_2, &nrm);
1654  MOFEM_LOG_C("WORLD", Sev::noisy, "\tdb nrm = %9.8e", nrm);
1655 
1657  }

◆ calculateDxAndDlambda()

MoFEMErrorCode FractureMechanics::GriffithForceElement::FrontArcLengthControl::calculateDxAndDlambda ( Vec  x)
inline

Definition at line 1659 of file GriffithForceElement.hpp.

1659  {
1661  MOFEM_LOG_CHANNEL("WORLD");
1663  MOFEM_LOG_TAG("WORLD", "Arc");
1664 
1665  // Calculate dx
1666  CHKERR VecGhostUpdateBegin(arcPtr->x0, INSERT_VALUES, SCATTER_FORWARD);
1667  CHKERR VecGhostUpdateEnd(arcPtr->x0, INSERT_VALUES, SCATTER_FORWARD);
1668 
1669  Vec l_x, l_x0, l_dx;
1670  CHKERR VecGhostGetLocalForm(x, &l_x);
1671  CHKERR VecGhostGetLocalForm(arcPtr->x0, &l_x0);
1672  CHKERR VecGhostGetLocalForm(arcPtr->dx, &l_dx);
1673  {
1674  double *x_array, *x0_array, *dx_array;
1675  CHKERR VecGetArray(l_x, &x_array);
1676  CHKERR VecGetArray(l_x0, &x0_array);
1677  CHKERR VecGetArray(l_dx, &dx_array);
1678  int size =
1679  problemPtr->getNbLocalDofsRow() + problemPtr->getNbGhostDofsRow();
1680  for (int i = 0; i != size; ++i) {
1681  dx_array[i] = x_array[i] - x0_array[i];
1682  }
1683  CHKERR VecRestoreArray(l_x, &x_array);
1684  CHKERR VecRestoreArray(l_x0, &x0_array);
1685  CHKERR VecRestoreArray(l_dx, &dx_array);
1686  }
1687  CHKERR VecGhostRestoreLocalForm(x, &l_x);
1688  CHKERR VecGhostRestoreLocalForm(arcPtr->x0, &l_x0);
1689  CHKERR VecGhostRestoreLocalForm(arcPtr->dx, &l_dx);
1690 
1691  // Calculate dlambda
1692  if (arcPtr->getPetscLocalDofIdx() != -1) {
1693  double *array;
1694  CHKERR VecGetArray(arcPtr->dx, &array);
1695  arcPtr->dLambda = array[arcPtr->getPetscLocalDofIdx()];
1696  array[arcPtr->getPetscLocalDofIdx()] = 0;
1697  CHKERR VecRestoreArray(arcPtr->dx, &array);
1698  }
1699  CHKERR VecGhostUpdateBegin(arcPtr->ghosTdLambda, INSERT_VALUES,
1700  SCATTER_FORWARD);
1701  CHKERR VecGhostUpdateEnd(arcPtr->ghosTdLambda, INSERT_VALUES,
1702  SCATTER_FORWARD);
1703 
1704  // Calculate dx2
1705  double x_nrm, x0_nrm;
1706  CHKERR VecNorm(x, NORM_2, &x_nrm);
1707  CHKERR VecNorm(arcPtr->x0, NORM_2, &x0_nrm);
1708  CHKERR VecDot(arcPtr->dx, arcPtr->dx, &arcPtr->dx2);
1709  MOFEM_LOG_C("WORLD", Sev::noisy,
1710  "\tx norm = %6.4e x0 norm = %6.4e dx2 = %6.4e", x_nrm, x0_nrm,
1711  arcPtr->dx2);
1713  }

◆ calculateLambdaInt()

double FractureMechanics::GriffithForceElement::FrontArcLengthControl::calculateLambdaInt ( )
inline

Calculate internal lambda.

Definition at line 1536 of file GriffithForceElement.hpp.

1536  {
1538  return intLambda;
1540  }

◆ operator()()

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

Definition at line 1542 of file GriffithForceElement.hpp.

1542  {
1544  MOFEM_LOG_CHANNEL("WORLD");
1546  MOFEM_LOG_TAG("WORLD", "Arc");
1547 
1548  switch (snes_ctx) {
1549  case CTX_SNESSETFUNCTION: {
1550  arcPtr->res_lambda = calculateLambdaInt() - arcPtr->alpha * arcPtr->s;
1551  CHKERR VecSetValue(snes_f, arcPtr->getPetscGlobalDofIdx(),
1552  arcPtr->res_lambda, ADD_VALUES);
1553  } break;
1554  case CTX_SNESSETJACOBIAN: {
1555  arcPtr->dIag = 2 * arcPtr->dLambda * arcPtr->alpha *
1556  pow(arcPtr->beta, 2) * arcPtr->F_lambda2;
1557  MOFEM_LOG_C("WORLD", Sev::noisy, "diag = %9.8e", arcPtr->dIag);
1558  CHKERR MatSetValue(snes_B, arcPtr->getPetscGlobalDofIdx(),
1559  arcPtr->getPetscGlobalDofIdx(), 1, ADD_VALUES);
1560  } break;
1561  default:
1562  break;
1563  }
1565  }

◆ postProcess()

MoFEMErrorCode FractureMechanics::GriffithForceElement::FrontArcLengthControl::postProcess ( )
inline

Definition at line 1567 of file GriffithForceElement.hpp.

1567  {
1569  switch (snes_ctx) {
1570  case CTX_SNESSETFUNCTION: {
1571  CHKERR VecAssemblyBegin(snes_f);
1572  CHKERR VecAssemblyEnd(snes_f);
1573  } break;
1574  case CTX_SNESSETJACOBIAN: {
1575  CHKERR VecGhostUpdateBegin(arcPtr->ghostDiag, INSERT_VALUES,
1576  SCATTER_FORWARD);
1577  CHKERR VecGhostUpdateEnd(arcPtr->ghostDiag, INSERT_VALUES,
1578  SCATTER_FORWARD);
1579  } break;
1580  default:
1581  break;
1582  }
1584  }

◆ preProcess()

MoFEMErrorCode FractureMechanics::GriffithForceElement::FrontArcLengthControl::preProcess ( )
inline

Definition at line 1519 of file GriffithForceElement.hpp.

1519  {
1521  switch (snes_ctx) {
1522  case CTX_SNESSETFUNCTION: {
1523  CHKERR VecAssemblyBegin(snes_f);
1524  CHKERR VecAssemblyEnd(snes_f);
1525  CHKERR calculateDxAndDlambda(snes_x);
1526  CHKERR calculateDb();
1527  } break;
1528  default:
1529  break;
1530  }
1532  }

Member Data Documentation

◆ arcPtr

boost::shared_ptr<ArcLengthCtx> FractureMechanics::GriffithForceElement::FrontArcLengthControl::arcPtr

Definition at line 1291 of file GriffithForceElement.hpp.

◆ blockData

BlockData& FractureMechanics::GriffithForceElement::FrontArcLengthControl::blockData

Definition at line 1292 of file GriffithForceElement.hpp.

◆ commonData

CommonData& FractureMechanics::GriffithForceElement::FrontArcLengthControl::commonData

Definition at line 1293 of file GriffithForceElement.hpp.

◆ feIntLambda

boost::shared_ptr<MyTriangleFE> FractureMechanics::GriffithForceElement::FrontArcLengthControl::feIntLambda

Definition at line 1302 of file GriffithForceElement.hpp.

◆ ghostIntLambda

Vec FractureMechanics::GriffithForceElement::FrontArcLengthControl::ghostIntLambda

Definition at line 1296 of file GriffithForceElement.hpp.

◆ intLambda

double FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambda

Definition at line 1300 of file GriffithForceElement.hpp.

◆ intLambdaArray

double FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaArray[2]

Definition at line 1297 of file GriffithForceElement.hpp.

◆ intLambdaDelta

double& FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaDelta

Definition at line 1298 of file GriffithForceElement.hpp.

◆ intLambdaLambda

double& FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaLambda

Definition at line 1299 of file GriffithForceElement.hpp.

◆ lambdaFieldName

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

Definition at line 1294 of file GriffithForceElement.hpp.

◆ tAg

int FractureMechanics::GriffithForceElement::FrontArcLengthControl::tAg

Definition at line 1290 of file GriffithForceElement.hpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
FractureMechanics::GriffithForceElement::FrontArcLengthControl::arcPtr
boost::shared_ptr< ArcLengthCtx > arcPtr
Definition: GriffithForceElement.hpp:1291
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
FractureMechanics::GriffithForceElement::FrontArcLengthControl::lambdaFieldName
std::string lambdaFieldName
Definition: GriffithForceElement.hpp:1294
FractureMechanics::GriffithForceElement::FrontArcLengthControl::blockData
BlockData & blockData
Definition: GriffithForceElement.hpp:1292
MOFEM_LOG_FUNCTION
#define MOFEM_LOG_FUNCTION()
Set scope.
Definition: LogManager.hpp:325
FractureMechanics::GriffithForceElement::FrontArcLengthControl::calculateDb
MoFEMErrorCode calculateDb()
Calculate db.
Definition: GriffithForceElement.hpp:1588
FractureMechanics::GriffithForceElement::FrontArcLengthControl::tAg
int tAg
Definition: GriffithForceElement.hpp:1290
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
_IT_NUMEREDDOF_ROW_BY_BITNUMBER_FOR_LOOP_
#define _IT_NUMEREDDOF_ROW_BY_BITNUMBER_FOR_LOOP_(PROBLEMPTR, FIELD_BIT_NUMBER, IT)
Definition: ProblemsMultiIndices.hpp:338
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
FractureMechanics::GriffithForceElement::FrontArcLengthControl::calculateDxAndDlambda
MoFEMErrorCode calculateDxAndDlambda(Vec x)
Definition: GriffithForceElement.hpp:1659
FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambda
double intLambda
Definition: GriffithForceElement.hpp:1300
FractureMechanics::GriffithForceElement::FrontArcLengthControl::feIntLambda
boost::shared_ptr< MyTriangleFE > feIntLambda
Definition: GriffithForceElement.hpp:1302
FractureMechanics::GriffithForceElement::FrontArcLengthControl::ghostIntLambda
Vec ghostIntLambda
Definition: GriffithForceElement.hpp:1296
MOFEM_LOG_TAG
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
Definition: LogManager.hpp:339
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
FractureMechanics::GriffithForceElement::FrontArcLengthControl::commonData
CommonData & commonData
Definition: GriffithForceElement.hpp:1293
FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaLambda
double & intLambdaLambda
Definition: GriffithForceElement.hpp:1299
FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaDelta
double & intLambdaDelta
Definition: GriffithForceElement.hpp:1298
FractureMechanics::GriffithForceElement::FrontArcLengthControl::intLambdaArray
double intLambdaArray[2]
Definition: GriffithForceElement.hpp:1297
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
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::FrontArcLengthControl::calculateLambdaInt
double calculateLambdaInt()
Calculate internal lambda.
Definition: GriffithForceElement.hpp:1536