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

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

Inheritance diagram for FractureMechanics::ConstantArea::OpTangentC:
[legend]
Collaboration diagram for FractureMechanics::ConstantArea::OpTangentC:
[legend]

Public Member Functions

 OpTangentC (int tag, CommonData &common_data)
 
PetscErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 
- Public Member Functions inherited from FractureMechanics::ConstantArea::AuxOp
 AuxOp (int tag, CommonData &common_data)
 
PetscErrorCode setVariables (FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 

Public Attributes

MatrixDouble cMat
 
- Public Attributes inherited from FractureMechanics::ConstantArea::AuxOp
int tAg
 
CommonDatacommonData
 
ublas::vector< int > rowIndices
 
VectorDouble activeVariables
 

Detailed Description

Definition at line 700 of file ConstantArea.hpp.

Constructor & Destructor Documentation

◆ OpTangentC()

FractureMechanics::ConstantArea::OpTangentC::OpTangentC ( int  tag,
CommonData common_data 
)
inline

Definition at line 703 of file ConstantArea.hpp.

705  "LAMBDA_CRACKFRONT_AREA_TANGENT", "MESH_NODE_POSITIONS",
707  AuxOp(tag, common_data) {
708  sYmm = false;
709  }

Member Function Documentation

◆ doWork()

PetscErrorCode FractureMechanics::ConstantArea::OpTangentC::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
DataForcesAndSourcesCore::EntData row_data,
DataForcesAndSourcesCore::EntData col_data 
)
inline

Definition at line 713 of file ConstantArea.hpp.

716  {
718 
719  if (row_type != MBVERTEX) {
721  }
722 
723  VectorInt &row_indices = row_data.getIndices();
724  const int nb_rows = row_indices.size();
725  if (!nb_rows)
727  VectorInt &col_indices = col_data.getIndices();
728  const int nb_cols = col_indices.size();
729  if (!nb_cols)
731  cMat.resize(3, 9, false);
732  cMat.clear();
733 
734  CHKERR setVariables(this, col_side, col_type, col_data);
735 
736  int nb_dofs = col_data.getFieldData().size();
737  commonData.crackFrontTangent.resize(nb_dofs, false);
738  int r;
739  // play recorder for values
740  r = ::function(tAg, nb_dofs, 18, &activeVariables[0],
742  if (r < 3) {
743  SETERRQ1(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
744  "ADOL-C function evaluation with error r = %d", r);
745  }
746 
747  // Set face orientation
748  {
749  // Get orientation in respect to adjacent tet
750  const EntityHandle face = getFEEntityHandle();
751  const BitRefLevel bit = getFEMethod()->problemPtr->getBitRefLevel();
752  Range adj_side_elems;
753  CHKERR commonData.mField.getInterface<BitRefManager>()->getAdjacencies(
754  bit, &face, 1, 3, adj_side_elems);
755  adj_side_elems = adj_side_elems.subset_by_type(MBTET);
756  if (adj_side_elems.size() != 1) {
757  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
758  "expect 1 tet but is %u", adj_side_elems.size());
759  }
760  EntityHandle side_elem = *adj_side_elems.begin();
761  int side_number, sense, offset;
762  CHKERR commonData.mField.get_moab().side_number(
763  side_elem, face, side_number, sense, offset);
764  if (sense == -1) {
766  }
767  // Get orientaton of face in respect to face in the prism
768  int side;
769  CHKERR commonData.mField.get_moab().tag_get_data(
770  commonData.thInterfaceSide, &face, 1, &side);
771  if (side == 1) {
773  }
774  }
775 
776  auto check_isnormal = [](double v) {
778  if (v != 0 && !std::isnormal(v)) {
779  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Not a number");
780  }
782  };
783 
784  for (int rr = 0; rr != nb_rows; ++rr) {
785  if (row_indices[rr] != -1) {
786  // cerr << commonData.crackFrontTangent(3*rr+0) << " "
787  // << commonData.crackFrontTangent(3*rr+1) << " "
788  // << commonData.crackFrontTangent(3*rr+2) << endl;
789  CHKERR check_isnormal(commonData.crackFrontTangent[3 * rr + 0]);
790  CHKERR check_isnormal(commonData.crackFrontTangent[3 * rr + 1]);
791  CHKERR check_isnormal(commonData.crackFrontTangent[3 * rr + 2]);
792  cMat(rr, 3 * rr + 0) = commonData.crackFrontTangent[3 * rr + 0];
793  cMat(rr, 3 * rr + 1) = commonData.crackFrontTangent[3 * rr + 1];
794  cMat(rr, 3 * rr + 2) = commonData.crackFrontTangent[3 * rr + 2];
795  }
796  }
797 
798  if (!commonData.setMapV) {
799  CHKERR MatSetValues(getFEMethod()->snes_B, nb_rows, &row_indices[0],
800  nb_cols, &col_indices[0], &cMat(0, 0), ADD_VALUES);
801  } else {
802  for (int rr = 0; rr != nb_rows; ++rr) {
803  if (row_indices[rr] != -1) {
804  if (commonData.mapV.find(row_indices[rr]) ==
805  commonData.mapV.end()) {
806  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
807  "Vector for given DOF in map not found");
808  }
809  for (int cc = 0; cc != nb_cols; ++cc) {
810  if (col_indices[cc] != -1) {
811  CHKERR check_isnormal(cMat(rr, cc));
812  }
813  }
814  CHKERR VecSetValues(commonData.mapV[row_indices[rr]], nb_cols,
815  &col_indices[0], &cMat(rr, 0), ADD_VALUES);
816  }
817  }
818  }
819 
821  }

Member Data Documentation

◆ cMat

MatrixDouble FractureMechanics::ConstantArea::OpTangentC::cMat

Definition at line 711 of file ConstantArea.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:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
FractureMechanics::ConstantArea::CommonData::setMapV
bool setMapV
Definition: ConstantArea.hpp:46
MoFEM::MatSetValues
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
Definition: EntitiesFieldData.hpp:1644
EntityHandle
FractureMechanics::ConstantArea::CommonData::mField
MoFEM::Interface & mField
Definition: ConstantArea.hpp:34
FractureMechanics::ConstantArea::CommonData::thInterfaceSide
Tag thInterfaceSide
Definition: ConstantArea.hpp:48
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROWCOL
@ OPROWCOL
operator doWork is executed on FE rows &columns
Definition: ForcesAndSourcesCore.hpp:569
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1589
sdf.r
int r
Definition: sdf.py:8
FractureMechanics::ConstantArea::CommonData::mapV
map< int, Vec > mapV
Definition: ConstantArea.hpp:45
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
FractureMechanics::ConstantArea::AuxOp::tAg
int tAg
Definition: ConstantArea.hpp:481
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
FractureMechanics::ConstantArea::AuxOp::AuxOp
AuxOp(int tag, CommonData &common_data)
Definition: ConstantArea.hpp:483
MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition: definitions.h:34
FractureMechanics::ConstantArea::AuxOp::activeVariables
VectorDouble activeVariables
Definition: ConstantArea.hpp:487
FractureMechanics::ConstantArea::AuxOp::commonData
CommonData & commonData
Definition: ConstantArea.hpp:482
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
FractureMechanics::ConstantArea::AuxOp::setVariables
PetscErrorCode setVariables(FaceElementForcesAndSourcesCore::UserDataOperator *fe_ptr, int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
Definition: ConstantArea.hpp:490
FractureMechanics::ConstantArea::OpTangentC::cMat
MatrixDouble cMat
Definition: ConstantArea.hpp:711
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
MoFEM::Types::VectorInt
UBlasVector< int > VectorInt
Definition: Types.hpp:67
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
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::ConstantArea::CommonData::crackFrontTangent
VectorDouble crackFrontTangent
Definition: ConstantArea.hpp:41