v0.10.0
Classes | Typedefs | Enumerations | Functions | Variables
FractureMechanics Namespace Reference

Classes

struct  AnalyticalDisp
 
struct  AnalyticalForces
 
struct  AnalyticalOptions
 
struct  ConstantArea
 Constant area constrains. More...
 
struct  CPMeshCut
 
struct  CPSolvers
 
struct  CrackFrontSingularBase
 
struct  CrackPropagation
 
struct  FirendVolumeOnSide
 
struct  GetSmoothingElementsSkin
 
struct  GriffithForceElement
 Implementation of Griffith element. More...
 
struct  MWLSApprox
 
struct  OpAleLhsWithDensitySingularElement_dX_dX
 
struct  OpAleLhsWithDensitySingularElement_dx_dX
 
struct  OpAnalyticalMaterialTraction
 
struct  OpAnalyticalSpatialTraction
 
struct  OpGetCrackFrontCommonDataAtGaussPts
 
struct  OpGetCrackFrontDataAtGaussPts
 
struct  OpGetCrackFrontDataGradientAtGaussPts
 
struct  OpGetDensityFieldForTesting
 Op to generate artificial density field. More...
 
struct  OpLhsBoneExplicitDerivariveWithHooke_dX
 Calculate explicit derivative of energy. More...
 
struct  OpLhsBoneExplicitDerivariveWithHooke_dx
 
struct  OpPostProcDisplacements
 
struct  OpPrint
 
struct  OpRhsBoneExplicitDerivariveWithHooke
 Calculate explicit derivative of energy. More...
 
struct  OpSetTagRangeOnSkin
 Mark crack surfaces on skin. More...
 
struct  OpTransfromSingularBaseFunctions
 
class  TwoType
 

Typedefs

typedef CrackFrontSingularBase< NonlinearElasticElement::MyVolumeFE, VolumeElementForcesAndSourcesCoreCrackFrontElement
 

Enumerations

enum  Materials {
  HOOKE, KIRCHHOFF, NEOHOOKEAN, BONEHOOKE,
  LASTOP
}
 
enum  CrackInterfaces { CRACK_PROPAGATION_INTERFACE = 1 << 0, CP_SOLVERS_INTERFACE = 1 << 1, CP_CUT_MESH_INTERFACE = 1 << 2 }
 
enum  AdolcTags {
  ELASTIC_TAG = 1, MATERIAL_TAG, GRIFFITH_FORCE_TAG, GRIFFITH_CONSTRAINS_TAG,
  CONSTANT_AREA_TAG, ARC_LENGTH_TAG, FRONT_TANGENT, EXTERIOR_DERIVATIVE_TAG,
  SMOOTHING_TAG, SURFACE_SLIDING_TAG, EDGE_SLIDING_TAG
}
 Tapes numbers used by ADOL-C. More...
 
enum  tangent_tests { MWLS_STRESS_TAN, MWLS_DENSITY_TAN, MWLS_GRIFFITH_TAN, LASTTAN }
 Names of tangent matrices tests. More...
 

Functions

static double calMax (double a, double b, double r)
 
static double diffCalMax_a (double a, double b, double r)
 
static MoFEMErrorCode elastic_snes_rhs (SNES snes, Vec x, Vec f, void *ctx)
 
static MoFEMErrorCode elastic_snes_mat (SNES snes, Vec x, Mat A, Mat B, void *ctx)
 
static MoFEMErrorCode propagation_snes_rhs (SNES snes, Vec x, Vec f, void *ctx)
 
static MoFEMErrorCode propagation_snes_mat (SNES snes, Vec x, Mat A, Mat B, void *ctx)
 

Variables

static const MOFEMuuid IDD_MOFEMCPMeshCutInterface
 
static const MOFEMuuid IDD_MOFEMCPSolversInterface
 
static const MOFEMuuid IDD_MOFEMCrackPropagationInterface
 
static std::map< long int, MatrixDoublemapRefCoords
 
const char * materials_list []
 

Typedef Documentation

◆ CrackFrontElement

Definition at line 378 of file CrackFrontElement.hpp.

Enumeration Type Documentation

◆ AdolcTags

◆ CrackInterfaces

Three interfaces used to solve fracture problem

Enumerator
CRACK_PROPAGATION_INTERFACE 
CP_SOLVERS_INTERFACE 
CP_CUT_MESH_INTERFACE 

Definition at line 27 of file CrackPropagation.hpp.

◆ Materials

◆ tangent_tests

Names of tangent matrices tests.

Enumerator
MWLS_STRESS_TAN 
MWLS_DENSITY_TAN 
MWLS_GRIFFITH_TAN 
LASTTAN 

Definition at line 62 of file CrackPropagation.hpp.

Function Documentation

◆ calMax()

static double FractureMechanics::calMax ( double  a,
double  b,
double  r 
)
static

Definition at line 24 of file GriffithForceElement.hpp.

24  {
25  return (a + b + (1 / r) * pow(fabs(a - b), r)) / 2;
26 }
const double r
rate factor

◆ diffCalMax_a()

static double FractureMechanics::diffCalMax_a ( double  a,
double  b,
double  r 
)
static

Definition at line 28 of file GriffithForceElement.hpp.

28  {
29  double sgn = ((a - b) == 0) ? 0 : (((a - b) < 0) ? -1 : 1);
30  return (1 + pow(fabs(a - b), r - 1) * sgn * (+1)) / 2;
31 }
const double r
rate factor

◆ elastic_snes_mat()

static MoFEMErrorCode FractureMechanics::elastic_snes_mat ( SNES  snes,
Vec  x,
Mat  A,
Mat  B,
void *  ctx 
)
static

Definition at line 7556 of file CrackPropagation.cpp.

7557  {
7558  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
7559  (CrackPropagation::ArcLengthSnesCtx *)ctx;
7560  DMMGViaApproxOrdersCtx *dm_ctx = arc_snes_ctx->dmCtx;
7562  CHKERR SnesMat(snes, x, A, B, ctx);
7563  CHKERR MatDiagonalScale(B, arc_snes_ctx->diagM, PETSC_NULL);
7564  if (dm_ctx->coarseningIS.size() > 1) {
7565 #if PETSC_VERSION_GE(3, 8, 0)
7566  CHKERR MatCreateSubMatrix(B, dm_ctx->coarseningIS[0],
7567  dm_ctx->coarseningIS[0], MAT_REUSE_MATRIX,
7568  &(dm_ctx->kspOperators[0]));
7569 #else
7570  CHKERR MatGetSubMatrix(B, dm_ctx->coarseningIS[0], dm_ctx->coarseningIS[0],
7571  MAT_REUSE_MATRIX, &(dm_ctx->kspOperators[0]));
7572 #endif
7573 
7574  // MatView(dm_ctx->kspOperators[0], PETSC_VIEWER_DRAW_WORLD);
7575  // std::string wait;
7576  // std::cin >> wait;
7577  // MatView(B, PETSC_VIEWER_DRAW_WORLD);
7578  // std::cin >> wait;
7579  }
7581 };
PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition: SnesCtx.cpp:126
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
Structure for DM for multi-grid via approximation orders.
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
std::vector< Mat > kspOperators
Get KSP operators.
std::vector< IS > coarseningIS
Coarsening IS.

◆ elastic_snes_rhs()

static MoFEMErrorCode FractureMechanics::elastic_snes_rhs ( SNES  snes,
Vec  x,
Vec  f,
void *  ctx 
)
static

Definition at line 7547 of file CrackPropagation.cpp.

7547  {
7548  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
7549  (CrackPropagation::ArcLengthSnesCtx *)ctx;
7551  CHKERR SnesRhs(snes, x, f, ctx);
7552  CHKERR VecPointwiseMult(f, f, arc_snes_ctx->diagM);
7554 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition: SnesCtx.cpp:17

◆ propagation_snes_mat()

static MoFEMErrorCode FractureMechanics::propagation_snes_mat ( SNES  snes,
Vec  x,
Mat  A,
Mat  B,
void *  ctx 
)
static

Definition at line 7840 of file CrackPropagation.cpp.

7841  {
7842  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
7843  static_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
7844  DMMGViaApproxOrdersCtx *dm_ctx = arc_snes_ctx->dmCtx;
7845  DMMGViaApproxOrdersCtx *dm_ctx_elastic = arc_snes_ctx->dmCtxElastic;
7847  CHKERR SnesMat(snes, x, A, B, ctx);
7848 // CHKERR MatDiagonalScale(B, arc_snes_ctx->diagM, PETSC_NULL);
7849 // CHKERR VecPointwiseMult(arc_snes_ctx->arcPtr->F_lambda,
7850 // arc_snes_ctx->arcPtr->F_lambda, arc_snes_ctx->diagM);
7851 // if (dm_ctx_elastic->coarseningIS.size() > 1) {
7852 // boost::shared_ptr<ComposedProblemsData> cmp_data_ptr =
7853 // dm_ctx->problemPtr->getComposedProblemsData();
7854 // #if PETSC_VERSION_GE(3, 8, 0)
7855 // CHKERR MatCreateSubMatrix(B, cmp_data_ptr->rowIs[0], cmp_data_ptr->rowIs[0],
7856 // MAT_REUSE_MATRIX,
7857 // &(dm_ctx_elastic->kspOperators.back()));
7858 // CHKERR MatCreateSubMatrix(dm_ctx_elastic->kspOperators.back(),
7859 // dm_ctx_elastic->coarseningIS[0],
7860 // dm_ctx_elastic->coarseningIS[0], MAT_REUSE_MATRIX,
7861 // &(dm_ctx->kspOperators[0]));
7862 // #else
7863 // CHKERR MatGetSubMatrix(B, cmp_data_ptr->rowIs[0], cmp_data_ptr->rowIs[0],
7864 // MAT_REUSE_MATRIX,
7865 // &(dm_ctx_elastic->kspOperators.back()));
7866 // CHKERR MatGetSubMatrix(dm_ctx_elastic->kspOperators.back(),
7867 // dm_ctx_elastic->coarseningIS[0],
7868 // dm_ctx_elastic->coarseningIS[0], MAT_REUSE_MATRIX,
7869 // &(dm_ctx->kspOperators[0]));
7870 // #endif
7871 // }
7872  // MatView(B,PETSC_VIEWER_DRAW_WORLD);
7873  // std::string wait;
7874  // std::cin >> wait;
7876 };
PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx)
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
Definition: SnesCtx.cpp:126
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
Structure for DM for multi-grid via approximation orders.
#define CHKERR
Inline error check.
Definition: definitions.h:604
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415

◆ propagation_snes_rhs()

static MoFEMErrorCode FractureMechanics::propagation_snes_rhs ( SNES  snes,
Vec  x,
Vec  f,
void *  ctx 
)
static

Definition at line 7804 of file CrackPropagation.cpp.

7804  {
7805  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
7806  static_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
7808  CHKERR SnesRhs(snes, x, f, ctx);
7809  // CHKERR VecPointwiseMult(f, f, arc_snes_ctx->diagM);
7810 
7811  SnesCtx *snes_ctx = (SnesCtx *)ctx;
7812  auto &m_field = snes_ctx->mField;
7813  auto prb = m_field.get_problem("CRACK_PROPAGATION");
7814  auto dofs = prb->getNumeredRowDofs();
7815 
7816  double *array;
7817  CHKERR VecGetArray(f, &array);
7818 
7819  constexpr double eps = 1e12;
7820  EntityHandle out_meshset;
7821  CHKERR m_field.get_moab().create_meshset(MESHSET_SET, out_meshset);
7822  for(auto &dof : *dofs) {
7823  auto loc_idx = dof->getPetscLocalDofIdx();
7824  if (array[loc_idx] > eps || array[loc_idx] != array[loc_idx]) {
7825  auto ent = dof->getEnt();
7826  CHKERR m_field.get_moab().add_entities(out_meshset, &ent, 1);
7827  }
7828  }
7829  CHKERR VecRestoreArray(f, &array);
7830 
7831  CHKERR m_field.get_moab().write_file("entites_with_large_residuals.vtk",
7832  "VTK", "", &out_meshset, 1);
7833  CHKERR m_field.get_moab().delete_entities(&out_meshset, 1);
7834 
7835 
7836 
7838 };
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:29
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:27
#define CHKERR
Inline error check.
Definition: definitions.h:604
virtual const Problem * get_problem(const std::string &problem_name) const =0
Get the problem object.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415
DEPRECATED boost::shared_ptr< NumeredDofEntity_multiIndex > & getNumeredRowDofs() const
get access to numeredRowDofsPtr storing DOFs on rows \deprected Use getNumeredRowDofsPtr
static const double eps
PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx)
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
Definition: SnesCtx.cpp:17

Variable Documentation

◆ IDD_MOFEMCPMeshCutInterface

const MOFEMuuid FractureMechanics::IDD_MOFEMCPMeshCutInterface
static
Initial value:

ID of CPMeshCut interface

Definition at line 27 of file CPMeshCut.hpp.

◆ IDD_MOFEMCPSolversInterface

const MOFEMuuid FractureMechanics::IDD_MOFEMCPSolversInterface
static
Initial value:
=
std::bitset< BITINTERFACEUID_SIZE > BitIntefaceId
Definition: Types.hpp:56

ID for CPSolvers interface

Definition at line 27 of file CPSolvers.hpp.

◆ IDD_MOFEMCrackPropagationInterface

const MOFEMuuid FractureMechanics::IDD_MOFEMCrackPropagationInterface
static
Initial value:

Unknown interface ID variable shared by all structs and used to register interface for CrackPropagation at CrackPropagation constructor

Definition at line 37 of file CrackPropagation.hpp.

◆ mapRefCoords

std::map<long int, MatrixDouble> FractureMechanics::mapRefCoords
static

Definition at line 51 of file CrackFrontElement.cpp.

◆ materials_list

const char* FractureMechanics::materials_list[]
Initial value:
= {"HOOKE", "KIRCHHOFF", "NEOHOOKEAN",
"BONEHOOKE"}

Definition at line 336 of file CrackPropagation.cpp.