v0.9.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  OpLhsBoneExplicitDerivariveWithHooke_dx
 
struct  OpLhsBoneExplicitDerivariveWithHooke_dX
 Calculate explicit derivative of energy. More...
 
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...
 

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 376 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

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 }

◆ 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 }

◆ elastic_snes_mat()

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

Definition at line 6859 of file CrackPropagation.cpp.

6860  {
6861  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
6862  (CrackPropagation::ArcLengthSnesCtx *)ctx;
6863  DMMGViaApproxOrdersCtx *dm_ctx = arc_snes_ctx->dmCtx;
6865  CHKERR SnesMat(snes, x, A, B, ctx);
6866  CHKERR MatDiagonalScale(B, arc_snes_ctx->diagM, PETSC_NULL);
6867  if (dm_ctx->coarseningIS.size() > 1) {
6868  CHKERR MatGetSubMatrix(B, dm_ctx->coarseningIS[0], dm_ctx->coarseningIS[0],
6869  MAT_REUSE_MATRIX, &(dm_ctx->kspOperators[0]));
6870  // MatView(dm_ctx->kspOperators[0], PETSC_VIEWER_DRAW_WORLD);
6871  // std::string wait;
6872  // std::cin >> wait;
6873  // MatView(B, PETSC_VIEWER_DRAW_WORLD);
6874  // std::cin >> wait;
6875  }
6877 };
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:121
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
Structure for DM for multi-grid via approximation orders.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
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 6850 of file CrackPropagation.cpp.

6850  {
6851  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
6852  (CrackPropagation::ArcLengthSnesCtx *)ctx;
6854  CHKERR SnesRhs(snes, x, f, ctx);
6855  CHKERR VecPointwiseMult(f, f, arc_snes_ctx->diagM);
6857 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
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:23

◆ propagation_snes_mat()

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

Definition at line 7095 of file CrackPropagation.cpp.

7096  {
7097  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
7098  static_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
7099  DMMGViaApproxOrdersCtx *dm_ctx = arc_snes_ctx->dmCtx;
7100  DMMGViaApproxOrdersCtx *dm_ctx_elastic = arc_snes_ctx->dmCtxElastic;
7102  CHKERR SnesMat(snes, x, A, B, ctx);
7103  CHKERR MatDiagonalScale(B, arc_snes_ctx->diagM, PETSC_NULL);
7104  CHKERR VecPointwiseMult(arc_snes_ctx->arcPtr->F_lambda,
7105  arc_snes_ctx->arcPtr->F_lambda, arc_snes_ctx->diagM);
7106  if (dm_ctx_elastic->coarseningIS.size() > 1) {
7107  boost::shared_ptr<ComposedProblemsData> cmp_data_ptr =
7109  CHKERR MatGetSubMatrix(B, cmp_data_ptr->rowIs[0], cmp_data_ptr->rowIs[0],
7110  MAT_REUSE_MATRIX,
7111  &(dm_ctx_elastic->kspOperators.back()));
7112  CHKERR MatGetSubMatrix(dm_ctx_elastic->kspOperators.back(),
7113  dm_ctx_elastic->coarseningIS[0],
7114  dm_ctx_elastic->coarseningIS[0], MAT_REUSE_MATRIX,
7115  &(dm_ctx->kspOperators[0]));
7116  }
7117  // MatView(B,PETSC_VIEWER_DRAW_WORLD);
7118  // std::string wait;
7119  // std::cin >> wait;
7121 };
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:121
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
const Problem * problemPtr
pinter to problem data structure
Definition: DMMoFEM.hpp:884
Structure for DM for multi-grid via approximation orders.
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
boost::shared_ptr< ComposedProblemsData > & getComposedProblemsData() const
Het composed problems data structure.
std::vector< Mat > kspOperators
Get KSP operators.
std::vector< IS > coarseningIS
Coarsening IS.

◆ propagation_snes_rhs()

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

Definition at line 7086 of file CrackPropagation.cpp.

7086  {
7087  CrackPropagation::ArcLengthSnesCtx *arc_snes_ctx =
7088  static_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
7090  CHKERR SnesRhs(snes, x, f, ctx);
7091  CHKERR VecPointwiseMult(f, f, arc_snes_ctx->diagM);
7093 };
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
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:23

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 53 of file CrackFrontElement.cpp.

◆ materials_list

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

Definition at line 332 of file CrackPropagation.cpp.