v0.14.0
Public Member Functions | Public Attributes | Private Attributes | List of all members
FractureMechanics::CrackPropagation::BothSurfaceConstrains Struct Reference

Constrains material displacement on both sides. More...

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

Inheritance diagram for FractureMechanics::CrackPropagation::BothSurfaceConstrains:
[legend]
Collaboration diagram for FractureMechanics::CrackPropagation::BothSurfaceConstrains:
[legend]

Public Member Functions

 BothSurfaceConstrains (MoFEM::Interface &m_field, const std::string lambda_field_name="LAMBDA_BOTH_SIDES", const std::string spatial_field_name="MESH_NODE_POSITIONS")
 Construct a new Both Surface Constrains object. More...
 
MoFEMErrorCode getOptions ()
 
MoFEMErrorCode preProcess ()
 
MoFEMErrorCode operator() ()
 
MoFEMErrorCode postProcess ()
 

Public Attributes

double aLpha
 
Range masterNodes
 ! Scaling parameter More...
 

Private Attributes

MoFEM::InterfacemField
 ! Node on which lagrange multiplier is set More...
 
const std::string lambdaFieldName
 
const std::string spatialFieldName
 

Detailed Description

Constrains material displacement on both sides.

Both sides of crack surface have the same material displacements.

Definition at line 957 of file CrackPropagation.hpp.

Constructor & Destructor Documentation

◆ BothSurfaceConstrains()

FractureMechanics::CrackPropagation::BothSurfaceConstrains::BothSurfaceConstrains ( MoFEM::Interface m_field,
const std::string  lambda_field_name = "LAMBDA_BOTH_SIDES",
const std::string  spatial_field_name = "MESH_NODE_POSITIONS" 
)
inline

Construct a new Both Surface Constrains object.

Parameters
m_field
lambda_field_namelagrange multiplayers
spatial_field_nameconstrained spatial field

Definition at line 966 of file CrackPropagation.hpp.

970  : mField(m_field), lambdaFieldName(lambda_field_name),
971  spatialFieldName(spatial_field_name), aLpha(1) {
972  ierr = getOptions();
973  CHKERRABORT(PETSC_COMM_WORLD, ierr);
974  }

Member Function Documentation

◆ getOptions()

MoFEMErrorCode FractureMechanics::CrackPropagation::BothSurfaceConstrains::getOptions ( )
inline

Definition at line 976 of file CrackPropagation.hpp.

976  {
978  ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "",
979  "Get Booth surface constrains element scaling",
980  "none");
981  CHKERRQ(ierr);
982  CHKERR PetscOptionsScalar("-booth_side_alpha", "scaling parameter", "",
983  aLpha, &aLpha, PETSC_NULL);
984  ierr = PetscOptionsEnd();
985  CHKERRQ(ierr);
987  }

◆ operator()()

MoFEMErrorCode FractureMechanics::CrackPropagation::BothSurfaceConstrains::operator() ( )

Definition at line 9802 of file CrackPropagation.cpp.

9802  {
9804 
9805  EntityHandle fe_ent = numeredEntFiniteElementPtr->getEnt();
9806  EntityHandle tri3;
9807  CHKERR mField.get_moab().side_element(fe_ent, 2, 3, tri3);
9808  double area = mField.getInterface<Tools>()->getTriArea(tri3);
9809 
9810  auto constrain_nodes = [&]() {
9812  std::vector<int> lambda_dofs(9, -1);
9813  std::vector<double> lambda_vals(9, 0);
9814 
9815  // get dofs of Lagrange multipliers
9816  int sum_dd = 0;
9817 
9818  auto row_dofs = getRowDofsPtr();
9819  auto lo_uid_lambda = DofEntity::getLoFieldEntityUId(
9820  getFieldBitNumber(lambdaFieldName), get_id_for_min_type<MBVERTEX>());
9821  auto hi_uid_lambda = DofEntity::getHiFieldEntityUId(
9822  getFieldBitNumber(lambdaFieldName), get_id_for_max_type<MBVERTEX>());
9823 
9824  for (auto it = row_dofs->lower_bound(lo_uid_lambda);
9825  it != row_dofs->upper_bound(hi_uid_lambda); ++it) {
9826  int side_number = it->get()->getSideNumberPtr()->side_number;
9827  if (side_number > 2)
9828  side_number -= 3;
9829  const int dof_idx = 3 * side_number + it->get()->getEntDofIdx();
9830  lambda_dofs[dof_idx] = it->get()->getPetscGlobalDofIdx();
9831  lambda_vals[dof_idx] = it->get()->getFieldData();
9832  sum_dd++;
9833  }
9834  if (sum_dd > 9)
9835  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
9836  "Only nine Lagrange multipliers can be on element (%d)", sum_dd);
9837 
9838  std::vector<int> positions_dofs(18, -1);
9839  std::vector<double> positions_vals(18, 0);
9840  std::vector<bool> positions_flags(18, true);
9841 
9842  auto lo_uid_pos = DofEntity::getLoFieldEntityUId(
9843  getFieldBitNumber(spatialFieldName), get_id_for_min_type<MBVERTEX>());
9844  auto hi_uid_pos = DofEntity::getHiFieldEntityUId(
9845  getFieldBitNumber(spatialFieldName), get_id_for_max_type<MBVERTEX>());
9846 
9847  // get dofs of material node positions
9848  for (auto it = row_dofs->lower_bound(lo_uid_pos);
9849  it != row_dofs->upper_bound(hi_uid_pos); ++it) {
9850  int side_number = it->get()->getSideNumberPtr()->side_number;
9851  int dof_idx = 3 * side_number + it->get()->getEntDofIdx();
9852  positions_dofs[dof_idx] = it->get()->getPetscGlobalDofIdx();
9853  positions_vals[dof_idx] = it->get()->getFieldData();
9854  Range vert = Range(it->get()->getEnt(), it->get()->getEnt());
9855  if (masterNodes.contains(vert)) {
9856  positions_flags[dof_idx] = false;
9857  }
9858  }
9859 
9860  const double beta = area * aLpha;
9861  switch (snes_ctx) {
9862  case CTX_SNESSETFUNCTION:
9863  CHKERR VecSetOption(snes_f, VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
9864  for (int ii = 0; ii != 9; ++ii) {
9865  if (lambda_dofs[ii] == -1)
9866  continue;
9867  double gap = positions_vals[0 + ii] - positions_vals[9 + ii];
9868  double val1 = beta * gap;
9869  CHKERR VecSetValue(snes_f, lambda_dofs[ii], val1, ADD_VALUES);
9870  double val2 = beta * lambda_vals[ii];
9871  if (positions_flags[0 + ii]) {
9872  CHKERR VecSetValue(snes_f, positions_dofs[0 + ii], +val2, ADD_VALUES);
9873  }
9874  if (positions_flags[9 + ii]) {
9875  CHKERR VecSetValue(snes_f, positions_dofs[9 + ii], -val2, ADD_VALUES);
9876  }
9877  }
9878  break;
9879  case CTX_SNESSETJACOBIAN:
9880  for (int ii = 0; ii != 9; ++ii) {
9881  if (lambda_dofs[ii] == -1)
9882  continue;
9883  CHKERR MatSetValue(snes_B, lambda_dofs[ii], positions_dofs[0 + ii],
9884  +1 * beta, ADD_VALUES);
9885  CHKERR MatSetValue(snes_B, lambda_dofs[ii], positions_dofs[9 + ii],
9886  -1 * beta, ADD_VALUES);
9887  if (positions_flags[0 + ii]) {
9888  CHKERR MatSetValue(snes_B, positions_dofs[0 + ii], lambda_dofs[ii],
9889  +1 * beta, ADD_VALUES);
9890  }
9891  if (positions_flags[9 + ii]) {
9892  CHKERR MatSetValue(snes_B, positions_dofs[9 + ii], lambda_dofs[ii],
9893  -1 * beta, ADD_VALUES);
9894  }
9895  }
9896  break;
9897  default:
9898  break;
9899  }
9901  };
9902 
9903  auto constrain_edges = [&]() {
9905 
9906  map<EntityHandle, EntityHandle> side_map;
9907  Range ents;
9908 
9909  for (auto s : {0, 1, 2}) {
9910  EntityHandle a, b;
9911  CHKERR mField.get_moab().side_element(fe_ent, 1, s, a);
9912  CHKERR mField.get_moab().side_element(fe_ent, 1, s + 6, b);
9913  side_map[a] = b;
9914  side_map[b] = a;
9915  ents.insert(a);
9916  ents.insert(b);
9917  }
9918  {
9919  EntityHandle a, b;
9920  CHKERR mField.get_moab().side_element(fe_ent, 2, 3, a);
9921  CHKERR mField.get_moab().side_element(fe_ent, 2, 4, b);
9922  side_map[a] = b;
9923  side_map[b] = a;
9924  ents.insert(a);
9925  ents.insert(b);
9926  }
9927 
9928  map<EntityHandle, std::vector<int>> side_map_lambda_dofs;
9929  map<EntityHandle, std::vector<double>> side_map_lambda_vals;
9930  map<EntityHandle, std::vector<int>> side_map_positions_dofs;
9931  map<EntityHandle, std::vector<double>> side_map_positions_vals;
9932  map<EntityHandle, double> side_map_positions_sign;
9933 
9934  for (auto e : ents) {
9935 
9936  auto row_dofs = getRowDofsPtr();
9937  auto lo_uid_lambda = row_dofs->lower_bound(DofEntity::getLoFieldEntityUId(
9938  getFieldBitNumber(lambdaFieldName), e));
9939  auto hi_uid_lambda = row_dofs->upper_bound(DofEntity::getHiFieldEntityUId(
9940  getFieldBitNumber(lambdaFieldName), e));
9941  const int nb_lambda_dofs = std::distance(lo_uid_lambda, hi_uid_lambda);
9942 
9943  std::vector<int> lambda_dofs(nb_lambda_dofs, -1);
9944  std::vector<double> lambda_vals(nb_lambda_dofs, 0);
9945 
9946  for (auto it = lo_uid_lambda; it != hi_uid_lambda; ++it) {
9947  const int dof_idx = it->get()->getEntDofIdx();
9948  lambda_dofs[dof_idx] = it->get()->getPetscGlobalDofIdx();
9949  lambda_vals[dof_idx] = it->get()->getFieldData();
9950  }
9951 
9952  if (nb_lambda_dofs) {
9953  side_map_lambda_dofs[e] = lambda_dofs;
9954  side_map_lambda_vals[e] = lambda_vals;
9955  side_map_lambda_dofs[side_map[e]] = lambda_dofs;
9956  side_map_lambda_vals[side_map[e]] = lambda_vals;
9957  }
9958 
9959  auto lo_uid_pos = row_dofs->lower_bound(DofEntity::getLoFieldEntityUId(
9960  getFieldBitNumber(spatialFieldName), e));
9961  auto hi_uid_pos = row_dofs->upper_bound(DofEntity::getHiFieldEntityUId(
9962  getFieldBitNumber(spatialFieldName), e));
9963  const int nb_dofs = std::distance(lo_uid_pos, hi_uid_pos);
9964 
9965  std::vector<int> positions_dofs(nb_dofs, -1);
9966  std::vector<double> positions_vals(nb_dofs, 0);
9967  side_map_positions_sign[e] = (!nb_lambda_dofs) ? 1 : -1;
9968 
9969  // get dofs of material node positions
9970  for (auto it = lo_uid_pos; it != hi_uid_pos; ++it) {
9971  int dof_idx = it->get()->getEntDofIdx();
9972  positions_dofs[dof_idx] = it->get()->getPetscGlobalDofIdx();
9973  positions_vals[dof_idx] = it->get()->getFieldData();
9974  }
9975 
9976  side_map_positions_dofs[e] = positions_dofs;
9977  side_map_positions_vals[e] = positions_vals;
9978  }
9979 
9980  for (auto m : side_map_positions_dofs) {
9981 
9982  auto e = m.first;
9983  if (side_map_lambda_dofs.find(e) != side_map_lambda_dofs.end()) {
9984 
9985  auto &lambda_dofs = side_map_lambda_dofs.at(e);
9986  auto &lambda_vals = side_map_lambda_vals.at(e);
9987  auto &positions_dofs = side_map_positions_dofs.at(e);
9988  auto &positions_vals = side_map_positions_vals.at(e);
9989  auto sign = side_map_positions_sign.at(e);
9990 
9991  const double beta = area * aLpha;
9992 
9993  switch (snes_ctx) {
9994  case CTX_SNESSETFUNCTION:
9995  for (int ii = 0; ii != lambda_dofs.size(); ++ii) {
9996  double val1 = beta * positions_vals[ii] * sign;
9997  CHKERR VecSetValue(snes_f, lambda_dofs[ii], val1, ADD_VALUES);
9998  }
9999  for (int ii = 0; ii != positions_dofs.size(); ++ii) {
10000  double val2 = sign * beta * lambda_vals[ii];
10001  CHKERR VecSetValue(snes_f, positions_dofs[ii], +val2, ADD_VALUES);
10002  }
10003  break;
10004  case CTX_SNESSETJACOBIAN:
10005  for (int ii = 0; ii != lambda_dofs.size(); ++ii) {
10006  CHKERR MatSetValue(snes_B, lambda_dofs[ii], positions_dofs[ii],
10007  beta * sign, ADD_VALUES);
10008  }
10009  for (int ii = 0; ii != positions_dofs.size(); ++ii) {
10010  CHKERR MatSetValue(snes_B, positions_dofs[ii], lambda_dofs[ii],
10011  beta * sign, ADD_VALUES);
10012  }
10013  break;
10014  default:
10015  break;
10016  }
10017  }
10018  }
10019 
10021  };
10022 
10023  CHKERR constrain_nodes();
10024  CHKERR constrain_edges();
10025 
10027 }

◆ postProcess()

MoFEMErrorCode FractureMechanics::CrackPropagation::BothSurfaceConstrains::postProcess ( )
inline

Definition at line 991 of file CrackPropagation.hpp.

991 { return 0; }

◆ preProcess()

MoFEMErrorCode FractureMechanics::CrackPropagation::BothSurfaceConstrains::preProcess ( )
inline

Definition at line 989 of file CrackPropagation.hpp.

989 { return 0; }

Member Data Documentation

◆ aLpha

double FractureMechanics::CrackPropagation::BothSurfaceConstrains::aLpha

Definition at line 993 of file CrackPropagation.hpp.

◆ lambdaFieldName

const std::string FractureMechanics::CrackPropagation::BothSurfaceConstrains::lambdaFieldName
private

Definition at line 998 of file CrackPropagation.hpp.

◆ masterNodes

Range FractureMechanics::CrackPropagation::BothSurfaceConstrains::masterNodes

! Scaling parameter

Definition at line 994 of file CrackPropagation.hpp.

◆ mField

MoFEM::Interface& FractureMechanics::CrackPropagation::BothSurfaceConstrains::mField
private

! Node on which lagrange multiplier is set

Definition at line 997 of file CrackPropagation.hpp.

◆ spatialFieldName

const std::string FractureMechanics::CrackPropagation::BothSurfaceConstrains::spatialFieldName
private

Definition at line 999 of file CrackPropagation.hpp.


The documentation for this struct was generated from the following files:
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
ContactOps::sign
double sign(double x)
Definition: ContactOps.hpp:580
FractureMechanics::CrackPropagation::BothSurfaceConstrains::spatialFieldName
const std::string spatialFieldName
Definition: CrackPropagation.hpp:999
EntityHandle
MoFEM::Tools
Auxiliary tools.
Definition: Tools.hpp:19
FractureMechanics::CrackPropagation::BothSurfaceConstrains::lambdaFieldName
const std::string lambdaFieldName
Definition: CrackPropagation.hpp:998
FractureMechanics::CrackPropagation::BothSurfaceConstrains::aLpha
double aLpha
Definition: CrackPropagation.hpp:993
FractureMechanics::CrackPropagation::BothSurfaceConstrains::getOptions
MoFEMErrorCode getOptions()
Definition: CrackPropagation.hpp:976
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
a
constexpr double a
Definition: approx_sphere.cpp:30
FractureMechanics::CrackPropagation::BothSurfaceConstrains::mField
MoFEM::Interface & mField
! Node on which lagrange multiplier is set
Definition: CrackPropagation.hpp:997
Range
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
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::CrackPropagation::BothSurfaceConstrains::masterNodes
Range masterNodes
! Scaling parameter
Definition: CrackPropagation.hpp:994