v0.14.0
Loading...
Searching...
No Matches
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, VolumeElementForcesAndSourcesCore > CrackFrontElement
 

Enumerations

enum  Materials {
  HOOKE , KIRCHHOFF , NEOHOOKEAN , BONEHOOKE ,
  LASTOP
}
 
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

MoFEMErrorCode clean_pcomms (moab::Interface &moab, boost::shared_ptr< WrapMPIComm > moab_comm_wrap)
 
MoFEMErrorCode broadcast_entities (moab::Interface &moab, moab::Interface &moab_tmp, const int from_proc, Range &entities, const bool adjacencies, const bool tags)
 
static double calMax (double a, double b, double r)
 
static double diffCalMax_a (double a, double b, double r)
 
MoFEMErrorCode clean_pcomms (moab::Interface &moab, boost::shared_ptr< WrapMPIComm > moab_comm_wrap)
 
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 std::map< long int, MatrixDouble > mapRefCoords
 
const char * materials_list []
 

Typedef Documentation

◆ CrackFrontElement

Definition at line 402 of file CrackFrontElement.hpp.

Enumeration Type Documentation

◆ AdolcTags

Tapes numbers used by ADOL-C.

Enumerator
ELASTIC_TAG 
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 

Definition at line 36 of file CrackPropagation.hpp.

◆ Materials

Enumerator
HOOKE 
KIRCHHOFF 
NEOHOOKEAN 
BONEHOOKE 
LASTOP 

Definition at line 23 of file CrackFrontElement.hpp.

◆ tangent_tests

Names of tangent matrices tests.

Enumerator
MWLS_STRESS_TAN 
MWLS_DENSITY_TAN 
MWLS_GRIFFITH_TAN 
LASTTAN 

Definition at line 54 of file CrackPropagation.hpp.

Function Documentation

◆ broadcast_entities()

MoFEMErrorCode FractureMechanics::broadcast_entities ( moab::Interface &  moab,
moab::Interface &  moab_tmp,
const int  from_proc,
Range entities,
const bool  adjacencies,
const bool  tags 
)

Definition at line 147 of file CrackPropagation.cpp.

150 {
151 ErrorCode result = MB_SUCCESS;
152 int success;
153 int buff_size;
154
155 auto add_verts = [&](Range &sent_ents) {
156 // Get the verts adj to these entities, since we'll have to send those
157 // too
158
159 // First check sets
160 std::pair<Range::const_iterator, Range::const_iterator> set_range =
161 sent_ents.equal_range(MBENTITYSET);
162 ErrorCode result = MB_SUCCESS, tmp_result;
163 for (Range::const_iterator rit = set_range.first; rit != set_range.second;
164 ++rit) {
165 tmp_result = moab_tmp.get_entities_by_type(*rit, MBVERTEX, sent_ents);
166 CHK_MOAB_THROW(tmp_result, "Failed to get contained verts");
167 }
168
169 // Now non-sets
170 Range tmp_ents;
171 std::copy(sent_ents.begin(), set_range.first, range_inserter(tmp_ents));
172 result = moab_tmp.get_adjacencies(tmp_ents, 0, false, sent_ents,
173 moab::Interface::UNION);
174 CHK_MOAB_THROW(result, "Failed to get vertices adj to ghosted ents");
175
176 return result;
177 };
178
179 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab_tmp, MYPCOMM_INDEX);
180 auto &procConfig = pcomm->proc_config();
181 const int MAX_BCAST_SIZE = (1 << 28);
182
183 ParallelComm::Buffer buff(ParallelComm::INITIAL_BUFF_SIZE);
184 buff.reset_ptr(sizeof(int));
185 if ((int)procConfig.proc_rank() == from_proc) {
186
187 result = add_verts(entities);
188 CHK_MOAB_THROW(result, "Failed to add adj vertices");
189
190 buff.reset_ptr(sizeof(int));
191 result = pcomm->pack_buffer(entities, adjacencies, tags, false, -1, &buff);
192 CHK_MOAB_THROW(result,
193 "Failed to compute buffer size in broadcast_entities");
194 buff.set_stored_size();
195 buff_size = buff.buff_ptr - buff.mem_ptr;
196 }
197
198 success =
199 MPI_Bcast(&buff_size, 1, MPI_INT, from_proc, procConfig.proc_comm());
200 if (MPI_SUCCESS != success) {
201 THROW_MESSAGE("MPI_Bcast of buffer size failed");
202 }
203
204 if (!buff_size) // No data
205 return MB_SUCCESS;
206
207 if ((int)procConfig.proc_rank() != from_proc)
208 buff.reserve(buff_size);
209
210 size_t offset = 0;
211 while (buff_size) {
212 int sz = std::min(buff_size, MAX_BCAST_SIZE);
213 success = MPI_Bcast(buff.mem_ptr + offset, sz, MPI_UNSIGNED_CHAR, from_proc,
214 procConfig.proc_comm());
215 if (MPI_SUCCESS != success) {
216 THROW_MESSAGE("MPI_Bcast of buffer failed");
217 }
218
219 offset += sz;
220 buff_size -= sz;
221 }
222
223 std::vector<std::vector<EntityHandle>> dum1a, dum1b;
224 std::vector<std::vector<int>> dum1p;
225 std::vector<EntityHandle> dum2, dum4;
226 std::vector<unsigned int> dum3;
227 buff.reset_ptr(sizeof(int));
228 if ((int)procConfig.proc_rank() == from_proc) {
229 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
230 result = pcomm->unpack_buffer(buff.buff_ptr, false, from_proc, -1, dum1a,
231 dum1b, dum1p, dum2, dum2, dum3, dum4);
232 } else {
233 result = pcomm->unpack_buffer(buff.buff_ptr, false, from_proc, -1, dum1a,
234 dum1b, dum1p, dum2, dum2, dum3, dum4);
235 }
236 CHK_MOAB_THROW(result, "Failed to unpack buffer in broadcast_entities");
237 std::copy(dum4.begin(), dum4.end(), range_inserter(entities));
238
239 return MB_SUCCESS;
240};
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
Definition: definitions.h:576
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
Definition: definitions.h:561

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

◆ clean_pcomms() [1/2]

MoFEMErrorCode FractureMechanics::clean_pcomms ( moab::Interface &  moab,
boost::shared_ptr< WrapMPIComm >  moab_comm_wrap 
)

◆ clean_pcomms() [2/2]

MoFEMErrorCode FractureMechanics::clean_pcomms ( moab::Interface &  moab,
boost::shared_ptr< WrapMPIComm moab_comm_wrap 
)

Definition at line 126 of file CrackPropagation.cpp.

127 {
129 std::vector<ParallelComm *> list_pcomms;
130 ParallelComm::get_all_pcomm(&moab, list_pcomms);
131 if (list_pcomms[MYPCOMM_INDEX]) {
132 CHKERR moab.tag_delete(list_pcomms[MYPCOMM_INDEX]->pstatus_tag());
133 CHKERR moab.tag_delete(list_pcomms[MYPCOMM_INDEX]->sharedps_tag());
134 CHKERR moab.tag_delete(list_pcomms[MYPCOMM_INDEX]->sharedhs_tag());
135 CHKERR moab.tag_delete(list_pcomms[MYPCOMM_INDEX]->sharedp_tag());
136 CHKERR moab.tag_delete(list_pcomms[MYPCOMM_INDEX]->sharedh_tag());
137 }
138 for (auto p : list_pcomms) {
139 delete p;
140 }
141 int pcomm_id = MYPCOMM_INDEX;
142 list_pcomms[MYPCOMM_INDEX] =
143 new ParallelComm(&moab, moab_comm_wrap->get_comm(), &pcomm_id);
145}
static Index< 'p', 3 > p
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535

◆ 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}
constexpr double a

◆ elastic_snes_mat()

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

Definition at line 8212 of file CrackPropagation.cpp.

8213 {
8214 auto arc_snes_ctx =
8215 reinterpret_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
8217 CHKERR SnesMat(snes, x, A, B, ctx);
8218 CHKERR MatDiagonalScale(B, arc_snes_ctx->getVecDiagM(), PETSC_NULL);
8220};
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:136

◆ elastic_snes_rhs()

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

Definition at line 8203 of file CrackPropagation.cpp.

8203 {
8204 auto arc_snes_ctx =
8205 reinterpret_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
8207 CHKERR SnesRhs(snes, x, f, ctx);
8208 CHKERR VecPointwiseMult(f, f, arc_snes_ctx->getVecDiagM());
8210}
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:27

◆ propagation_snes_mat()

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

Definition at line 8476 of file CrackPropagation.cpp.

8477 {
8478 auto arc_snes_ctx =
8479 reinterpret_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
8481
8482 CHKERR SnesMat(snes, x, A, B, ctx);
8483
8484 CHKERR MatDiagonalScale(B, arc_snes_ctx->getVecDiagM(), PETSC_NULL);
8485 CHKERR VecPointwiseMult(arc_snes_ctx->getArcPtr()->F_lambda,
8486 arc_snes_ctx->getArcPtr()->F_lambda,
8487 arc_snes_ctx->getVecDiagM());
8488
8490};

◆ propagation_snes_rhs()

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

Definition at line 8467 of file CrackPropagation.cpp.

8467 {
8468 auto arc_snes_ctx =
8469 reinterpret_cast<CrackPropagation::ArcLengthSnesCtx *>(ctx);
8471 CHKERR SnesRhs(snes, x, f, ctx);
8472 CHKERR VecPointwiseMult(f, f, arc_snes_ctx->getVecDiagM());
8474};

Variable Documentation

◆ 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 492 of file CrackPropagation.cpp.