v0.15.5
Loading...
Searching...
No Matches
EshelbianCore.hpp
Go to the documentation of this file.
1/**
2 * @file EshelbianCore.hpp
3 * @author your name (you@domain.com)
4 * @brief
5 * @version 0.1
6 * @date 2024-12-31
7 *
8 * @copyright Copyright (c) 2024
9 *
10 */
11
13
21
22 static inline const char *listSolvers[] = {
23 "time_solver", "dynamic_relaxation", "cohesive", "shape_optimization"};
24
32
33 static inline enum SolverType solverType = TimeSolver;
35 static inline enum SymmetrySelector symmetrySelector = NOT_SYMMETRIC;
36 static inline enum RotSelector rotSelector = LARGE_ROT;
37 static inline enum RotSelector gradApproximator = LARGE_ROT;
38 static inline enum StretchSelector stretchSelector = LOG;
39 static inline PetscBool noStretch = PETSC_FALSE; //< no stretch field
40 static inline PetscBool setSingularity = PETSC_FALSE; //< set singularity
41 static inline PetscBool dynamicRelaxation =
42 PETSC_FALSE; //< switch on/off dynamic relaxation
43 static inline PetscBool crackingOn = PETSC_FALSE; //< cracking on
44 static inline double crackingStartTime = 0; //< time when crack starts to grow
45 static inline double crackingAddTime = 0; //< time to add new crack surface
46 static inline int nbJIntegralContours =
47 0; //< number of contours for J integral evaluation
48 static inline double dynamicTime =
49 0; //< true time used when dynamic relaxation is used.
50 static inline int dynamicStep =
51 0; //< true time used when dynamic relaxation is used.
52 static inline PetscBool l2UserBaseScale = PETSC_TRUE; //< scale L2 user base
53 static inline int addCrackMeshsetId = 1000; //< add crack meshset id
54 static inline double griffithEnergy = 1; ///< Griffith energy
55 static inline double crackingRtol = 1e-10; ///< Cracking relative tolerance
56 static inline double crackingAtol = 1e-12; ///< Cracking absolute tolerance
57 static inline enum EnergyReleaseSelector energyReleaseSelector =
58 GRIFFITH_SKELETON; //< energy release selector
59 static inline std::string internalStressTagName =
60 ""; //< internal stress tag name
61 static inline int internalStressInterpOrder =
62 1; //< internal stress interpolation order
63 static inline PetscBool internalStressVoigt =
64 PETSC_FALSE; //< internal stress index notation
65 static inline PetscBool interfaceCrack =
66 PETSC_FALSE; //< interface crack tracking
67 static inline int interfaceRemoveLevel =
68 0; //< number of levels of elements to remove around interface crack
69
70 static double exponentBase;
71 static boost::function<double(const double)> f;
72 static boost::function<double(const double)> d_f;
73 static boost::function<double(const double)> dd_f;
74 static boost::function<double(const double)> inv_f;
75 static boost::function<double(const double)> inv_d_f;
76 static boost::function<double(const double)> inv_dd_f;
77
78 static inline constexpr bool use_quadratic_exp = true;
79 static inline constexpr double v_max = 24;
80 static inline constexpr double v_min = -v_max;
81
82 static double f_log_e_quadratic(const double v) {
83 if (v > v_max) {
84 double e = static_cast<double>(std::exp(v_max));
85 double dv = v - v_max;
86 return 0.5 * e * dv * dv + e * dv + e;
87 } else {
88 return static_cast<double>(std::exp(v));
89 }
90 }
91
92 static double d_f_log_e_quadratic(const double v) {
93 if (v > v_max) {
94 double e = static_cast<double>(std::exp(v_max));
95 double dv = v - v_max;
96 return e * dv + e;
97 } else {
98 return static_cast<double>(std::exp(v));
99 }
100 }
101
102 static double dd_f_log_e_quadratic(const double v) {
103 if (v > v_max) {
104 return static_cast<double>(std::exp(v_max));
105 } else {
106 return static_cast<double>(std::exp(v));
107 }
108 }
109
110 static double f_log_e(const double v) {
111 if constexpr(use_quadratic_exp) {
112 return f_log_e_quadratic(v);
113 } else {
114 if (v > v_max)
115 // y = exp(v_max) * v + exp(v_max) * (1 - v_max);
116 // y/exp(v_max) = v + (1 - v_max);
117 // y/exp(v_max) - (1 - v_max) = v;
118 return std::exp(v_max) * v + std::exp(v_max) * (1 - v_max);
119 else
120 return std::exp(v);
121 }
122 }
123 static double d_f_log_e(const double v) {
124 if constexpr (use_quadratic_exp) {
125 return d_f_log_e_quadratic(v);
126 } else {
127 if (v > v_max)
128 return std::exp(v_max);
129 else
130 return std::exp(v);
131 }
132 }
133 static double dd_f_log_e(const double v) {
134 if constexpr (use_quadratic_exp) {
135 return dd_f_log_e_quadratic(v);
136 } else {
137 if (v > v_max)
138 return 0.;
139 else
140 return std::exp(v);
141 }
142 }
143 static double inv_f_log_e(const double v) {
144 if (v > std::exp(v_min))
145 return std::log(v);
146 else
147 // y = exp(v_min) * v + exp(v_min) * (1 - v_min);
148 // y/exp(v_min) = v + (1 - v_min);
149 // y/exp(v_min) - (1 - v_min) = v;
150 return v / exp(v_min) - (1 - v_min);
151 }
152 static double inv_d_f_log_e(const double v) {
153 if (v > std::exp(v_min))
154 return 1. / v;
155 else
156 return 1. / exp(v_min);
157 }
158 static double inv_dd_f_log_e(const double v) {
159 if (v > std::exp(v_min))
160 return -1. / (v * v);
161 else
162 return 0.;
163 }
164
165 static double f_log(const double v) {
166 return pow(EshelbianCore::exponentBase, v);
167 }
168 static double d_f_log(const double v) {
169 return pow(exponentBase, v) * log(EshelbianCore::exponentBase);
170 }
171 static double dd_f_log(const double v) {
172 return pow(EshelbianCore::exponentBase, v) *
174 }
175
176 static double inv_f_log(const double v) {
177 return log(v) / log(EshelbianCore::exponentBase);
178 }
179 static double inv_d_f_log(const double v) {
180 return (1. / v) / log(EshelbianCore::exponentBase);
181 }
182 static double inv_dd_f_log(const double v) {
183 return -(1. / (v * v)) / log(EshelbianCore::exponentBase);
184 }
185
186 static double f_linear(const double v) { return v + 1; }
187 static double d_f_linear(const double v) { return 1; }
188 static double dd_f_linear(const double v) { return 0; }
189
190 static double inv_f_linear(const double v) { return v - 1; }
191 static double inv_d_f_linear(const double v) { return 0; }
192 static double inv_dd_f_linear(const double v) { return 0; }
193
194 /**
195 * \brief Getting interface of core database
196 * @param uuid unique ID of interface
197 * @param iface returned pointer to interface
198 * @return error code
199 */
200 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
201 UnknownInterface **iface) const;
202
204
205 boost::shared_ptr<DataAtIntegrationPts> dataAtPts;
206 boost::shared_ptr<PhysicalEquations> physicalEquations;
207 boost::shared_ptr<AnalyticalExprPython> AnalyticalExprPythonPtr;
208
209 boost::shared_ptr<VolumeElementForcesAndSourcesCore> elasticFeRhs;
210 boost::shared_ptr<VolumeElementForcesAndSourcesCore> elasticFeLhs;
211 boost::shared_ptr<FaceElementForcesAndSourcesCore> elasticBcLhs;
212 boost::shared_ptr<FaceElementForcesAndSourcesCore> elasticBcRhs;
213 boost::shared_ptr<ForcesAndSourcesCore>
214 contactTreeRhs; ///< Make a contact tree
215
216 SmartPetscObj<DM> dM; ///< Coupled problem all fields
217 SmartPetscObj<DM> dmElastic; ///< Elastic problem
218 SmartPetscObj<DM> dmMaterial; ///< Material problem
219 SmartPetscObj<DM> dmPrjSpatial; ///< Projection spatial displacement
220
221 const std::string piolaStress = "P";
222 // const std::string eshelbyStress = "S";
223 const std::string spatialL2Disp = "wL2";
224 // const std::string materialL2Disp = "WL2";
225 const std::string spatialH1Disp = "wH1";
226 const std::string materialH1Positions = "XH1";
227 const std::string hybridSpatialDisp = "hybridSpatialDisp";
228 // const std::string hybridMaterialDisp = "hybridMaterialDisp";
229 const std::string contactDisp = "contactDisp";
230 const std::string stretchTensor = "u";
231 const std::string rotAxis = "omega";
232 const std::string bubbleField = "bubble";
233
234 const std::string elementVolumeName = "EP";
235 const std::string naturalBcElement = "NATURAL_BC";
236 const std::string skinElement = "SKIN";
237 const std::string skeletonElement = "SKELETON";
238 const std::string contactElement = "CONTACT";
239
241 virtual ~EshelbianCore();
242
243 int spaceOrder = 2;
244 int spaceH1Order = -1;
246 double alphaU = 0;
247 double alphaW = 0;
248 double alphaOmega = 0;
249 double alphaRho = 0;
250 double alphaTau = 0;
251
252 int contactRefinementLevels = 1; //< refinement levels for contact integration
253 int frontLayers = 3; //< front layers with material element
254
255 MoFEMErrorCode getOptions();
256
257 boost::shared_ptr<BcDispVec> bcSpatialDispVecPtr;
258 boost::shared_ptr<BcRotVec> bcSpatialRotationVecPtr;
259 boost::shared_ptr<TractionBcVec> bcSpatialTractionVecPtr;
260 boost::shared_ptr<TractionFreeBc> bcSpatialFreeTractionVecPtr;
261 boost::shared_ptr<NormalDisplacementBcVec> bcSpatialNormalDisplacementVecPtr;
262 boost::shared_ptr<AnalyticalDisplacementBcVec>
264 boost::shared_ptr<AnalyticalTractionBcVec> bcSpatialAnalyticalTractionVecPtr;
265 boost::shared_ptr<PressureBcVec> bcSpatialPressureVecPtr;
266 boost::shared_ptr<ExternalStrainVec> externalStrainVecPtr;
267
268 std::map<std::string, boost::shared_ptr<ScalingMethod>> timeScaleMap;
269
270 template <typename BC>
271 MoFEMErrorCode getBc(boost::shared_ptr<BC> &bc_vec_ptr,
272 const std::string block_name, const int nb_attributes) {
274 for (auto it :
275 mField.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
276
277 (boost::format("%s(.*)") % block_name).str()
278
279 ))
280
281 ) {
282 std::vector<double> block_attributes;
283 CHKERR it->getAttributes(block_attributes);
284 if (block_attributes.size() < nb_attributes) {
285 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
286 "In block %s expected %d attributes, but given %ld",
287 it->getName().c_str(), nb_attributes, block_attributes.size());
288 }
289 Range faces;
290 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 2, faces,
291 true);
292 bc_vec_ptr->emplace_back(it->getName(), block_attributes, faces);
293 }
295 }
296
297 MoFEMErrorCode getSpatialDispBc();
298
299 inline MoFEMErrorCode getSpatialRotationBc() {
301 bcSpatialRotationVecPtr = boost::make_shared<BcRotVec>();
302 CHKERR getBc(bcSpatialRotationVecPtr, "SPATIAL_ROTATION_BC", 4);
303 CHKERR getBc(bcSpatialRotationVecPtr, "SPATIAL_ROTATION_AXIS_BC", 7);
304
305 auto ts_rotation =
306 boost::make_shared<DynamicRelaxationTimeScale>("rotation_history.txt");
307 for (auto &bc : *bcSpatialRotationVecPtr) {
308 timeScaleMap[bc.blockName] =
309 GetBlockScalingMethod<DynamicRelaxationTimeScale>::get(
310 ts_rotation, "rotation_history", ".txt", bc.blockName);
311 }
312
314 }
315
316 MoFEMErrorCode getSpatialTractionBc();
317
318 /**
319 * @brief Remove all, but entities where kinematic constrains are applied.
320 *
321 * @param meshset
322 * @param bc_ptr
323 * @param disp_block_set_name
324 * @param rot_block_set_name
325 * @param contact_set_name
326 * @return MoFEMErrorCode
327 */
328 MoFEMErrorCode getTractionFreeBc(const EntityHandle meshset,
329 boost::shared_ptr<TractionFreeBc> &bc_ptr,
330 const std::string contact_set_name);
331
332 inline MoFEMErrorCode
335 boost::shared_ptr<TractionFreeBc>(new TractionFreeBc());
336 return getTractionFreeBc(meshset, bcSpatialFreeTractionVecPtr, "CONTACT");
337 }
338
339 MoFEMErrorCode getExternalStrain();
340
341 MoFEMErrorCode createExchangeVectors(Sev sev);
342
343 MoFEMErrorCode addFields(const EntityHandle meshset = 0);
344 MoFEMErrorCode projectGeometry(const EntityHandle meshset = 0,
345 double time = 0);
346 MoFEMErrorCode projectInternalStress(const EntityHandle meshset = 0);
347
348 MoFEMErrorCode addVolumeFiniteElement(const EntityHandle meshset = 0);
349 MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset = 0);
350 MoFEMErrorCode addDMs(const BitRefLevel bit = BitRefLevel().set(0),
351 const EntityHandle meshset = 0);
352
353 MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape,
354 const double lambda,
355 const double mu,
356 const double sigma_y);
357
358 MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha,
359 const double beta,
360 const double lambda,
361 const double sigma_y);
362
363 MoFEMErrorCode addMaterial_HMHNeohookean(const int tape, const double c10,
364 const double K);
365
366 MoFEMErrorCode addMaterial_Hencky(double E, double nu);
367
368 MoFEMErrorCode setBaseVolumeElementOps(
369 const int tag, const bool do_rhs, const bool do_lhs,
370 const bool calc_rates,
371 boost::shared_ptr<VolumeElementForcesAndSourcesCore> fe);
372
373 MoFEMErrorCode setVolumeElementOps(
374 const int tag, const bool add_elastic, const bool add_material,
375 boost::shared_ptr<VolumeElementForcesAndSourcesCore> &fe_rhs,
376 boost::shared_ptr<VolumeElementForcesAndSourcesCore> &fe_lhs);
377
378 MoFEMErrorCode
379 setFaceElementOps(const bool add_elastic, const bool add_material,
380 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_rhs,
381 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_lhs);
382
383 MoFEMErrorCode setFaceInterfaceOps(
384 const bool add_elastic, const bool add_material,
385 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_rhs,
386 boost::shared_ptr<FaceElementForcesAndSourcesCore> &fe_lhs);
387
388 MoFEMErrorCode setContactElementRhsOps(
389
390 boost::shared_ptr<ForcesAndSourcesCore> &fe_contact_tree
391
392 );
393
394 MoFEMErrorCode setElasticElementOps(const int tag);
395 MoFEMErrorCode setElasticElementToTs(DM dm);
396
397 /** \brief Add debug to model
398 *
399 * That prints information every SNES step
400 */
401 MoFEMErrorCode addDebugModel(TS ts);
402
403 friend struct solve_elastic_set_up;
404 MoFEMErrorCode solveElastic(TS ts, Vec x);
405
406 /**
407 * @brief Solve problem using dynamic relaxation method
408 *
409 * @param ts solver time stepper
410 * @param x solution vector
411 * @param start_step starting step number
412 * @param start_time starting time
413 * @return MoFEMErrorCode
414 */
415 MoFEMErrorCode solveDynamicRelaxation(TS ts, Vec x, int start_step,
416 double start_time);
417
418 /**
419 * @brief Solve cohesive crack growth problem
420 *
421 * @param ts
422 * @param x
423 * @return * MoFEMErrorCode
424 */
425 MoFEMErrorCode solveCohesiveCrackGrowth(TS ts, Vec x, int start_step,
426 double start_time);
427
428 /**
429 * @brief Solve cohesive crack growth problem
430 *
431 * @param ts
432 * @param x
433 * @return * MoFEMErrorCode
434 */
435 MoFEMErrorCode solveSchapeOptimisation(TS ts, Vec x, int start_step,
436 double start_time);
437
438 MoFEMErrorCode setBlockTagsOnSkin();
439
440 MoFEMErrorCode postProcessResults(const int tag, const std::string file,
441 Vec f_residual = PETSC_NULLPTR,
442 Vec var_vec = PETSC_NULLPTR,
443 std::vector<Tag> tags_to_transfer = {},
444 TS ts = PETSC_NULLPTR);
445 MoFEMErrorCode
446 postProcessSkeletonResults(const int tag, const std::string file,
447 Vec f_residual = PETSC_NULLPTR,
448 std::vector<Tag> tags_to_transfer = {});
449
451 boost::shared_ptr<double> area_ptr); // calculate crack area
452
454
455 struct SetUpSchur {
456 static boost::shared_ptr<SetUpSchur> createSetUpSchur(
457
458 MoFEM::Interface &m_field, EshelbianCore *ep_core_ptr
459
460 );
461 virtual MoFEMErrorCode setUp(TS) = 0;
462
463 protected:
464 SetUpSchur() = default;
465 };
466
468
469 using TimeScale::TimeScale;
470
471 double getScale(const double time) override {
473 return TimeScale::getScale(EshelbianCore::dynamicTime);
474 else
475 return TimeScale::getScale(time);
476 }
477 };
478
479 MoFEMErrorCode calculateFaceMaterialForce(const int tag, TS ts);
480 MoFEMErrorCode calculateOrientation(const int tag, bool set_orientation);
481 MoFEMErrorCode setNewFrontCoordinates();
482 MoFEMErrorCode addCrackSurfaces(const bool debug = false);
483 MoFEMErrorCode saveOrgCoords();
484 MoFEMErrorCode createCrackSurfaceMeshset();
485
486 boost::shared_ptr<Range> contactFaces;
487 boost::shared_ptr<Range> crackFaces;
488 boost::shared_ptr<Range> frontEdges;
489 boost::shared_ptr<Range> frontAdjEdges;
490 boost::shared_ptr<Range> frontVertices;
491 boost::shared_ptr<Range> skeletonFaces;
492 boost::shared_ptr<Range> maxMovedFaces;
493 boost::shared_ptr<Range> interfaceFaces;
494
495 boost::shared_ptr<ParentFiniteElementAdjacencyFunctionSkeleton<2>>
497
498 BitRefLevel bitAdjParent = BitRefLevel().set(); ///< bit ref level for parent
499 BitRefLevel bitAdjParentMask =
500 BitRefLevel().set(); ///< bit ref level for parent parent
501 BitRefLevel bitAdjEnt = BitRefLevel().set(); ///< bit ref level for parent
502 BitRefLevel bitAdjEntMask =
503 BitRefLevel().set(); ///< bit ref level for parent parent
504
505 SmartPetscObj<Vec> solTSStep;
506
507 CommInterface::EntitiesPetscVector volumeExchange;
508 CommInterface::EntitiesPetscVector faceExchange;
509 CommInterface::EntitiesPetscVector edgeExchange;
510 CommInterface::EntitiesPetscVector vertexExchange;
511
512 std::vector<Tag>
513 listTagsToTransfer; ///< list of tags to transfer to postprocessor
514
515 Mat S = PETSC_NULLPTR; //< Schur complement matrix
516 AO aoS = PETSC_NULLPTR; //< AO for Schur complement matrix
517 SmartPetscObj<IS> crackHybridIs; //< IS for crack hybrid field
518 std::vector<std::string> a00FieldList; //< list of fields for Schur complement
519 std::vector<boost::shared_ptr<Range>>
520 a00RangeList; //< list of ranges for Schur complement
521
522 int nbCrackFaces = 0; //< number of crack faces
523};
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static const bool debug
auto bit
set bit
static double lambda
const double v
phase velocity of light in medium (cm/ns)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
double getScale(const double time) override
virtual MoFEMErrorCode setUp(TS)=0
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field, EshelbianCore *ep_core_ptr)
std::vector< boost::shared_ptr< Range > > a00RangeList
MoFEMErrorCode setElasticElementOps(const int tag)
boost::shared_ptr< ExternalStrainVec > externalStrainVecPtr
MoFEMErrorCode addVolumeFiniteElement(const EntityHandle meshset=0)
MoFEMErrorCode addFields(const EntityHandle meshset=0)
static constexpr bool use_quadratic_exp
static enum StretchSelector stretchSelector
boost::shared_ptr< Range > frontAdjEdges
MoFEMErrorCode createCrackSurfaceMeshset()
static int interfaceRemoveLevel
MoFEMErrorCode addBoundaryFiniteElement(const EntityHandle meshset=0)
const std::string skeletonElement
static double inv_f_linear(const double v)
MoFEMErrorCode getSpatialRotationBc()
boost::shared_ptr< TractionBcVec > bcSpatialTractionVecPtr
boost::shared_ptr< Range > contactFaces
static double dd_f_log_e_quadratic(const double v)
static double dynamicTime
static double dd_f_log(const double v)
BitRefLevel bitAdjEnt
bit ref level for parent
MoFEMErrorCode postProcessResults(const int tag, const std::string file, Vec f_residual=PETSC_NULLPTR, Vec var_vec=PETSC_NULLPTR, std::vector< Tag > tags_to_transfer={}, TS ts=PETSC_NULLPTR)
static boost::function< double(const double)> inv_dd_f
MoFEM::Interface & mField
static constexpr double v_max
const std::string spatialL2Disp
static double inv_d_f_log(const double v)
std::map< std::string, boost::shared_ptr< ScalingMethod > > timeScaleMap
MoFEMErrorCode addMaterial_Hencky(double E, double nu)
static enum SolverType solverType
friend struct solve_elastic_set_up
static PetscBool l2UserBaseScale
SmartPetscObj< DM > dM
Coupled problem all fields.
boost::shared_ptr< FaceElementForcesAndSourcesCore > elasticBcRhs
MoFEMErrorCode projectInternalStress(const EntityHandle meshset=0)
MoFEMErrorCode solveSchapeOptimisation(TS ts, Vec x, int start_step, double start_time)
Solve cohesive crack growth problem.
static int internalStressInterpOrder
SmartPetscObj< IS > crackHybridIs
boost::shared_ptr< TractionFreeBc > bcSpatialFreeTractionVecPtr
static const char * listSolvers[]
const std::string materialH1Positions
static int nbJIntegralContours
MoFEMErrorCode setBlockTagsOnSkin()
std::vector< Tag > listTagsToTransfer
list of tags to transfer to postprocessor
boost::shared_ptr< FaceElementForcesAndSourcesCore > elasticBcLhs
static PetscBool crackingOn
MoFEMErrorCode addMaterial_HMHHStVenantKirchhoff(const int tape, const double lambda, const double mu, const double sigma_y)
MoFEMErrorCode addMaterial_HMHMooneyRivlin(const int tape, const double alpha, const double beta, const double lambda, const double sigma_y)
MoFEMErrorCode getTractionFreeBc(const EntityHandle meshset, boost::shared_ptr< TractionFreeBc > &bc_ptr, const std::string contact_set_name)
Remove all, but entities where kinematic constrains are applied.
MoFEMErrorCode setBaseVolumeElementOps(const int tag, const bool do_rhs, const bool do_lhs, const bool calc_rates, boost::shared_ptr< VolumeElementForcesAndSourcesCore > fe)
MoFEMErrorCode calculateFaceMaterialForce(const int tag, TS ts)
static double griffithEnergy
Griffith energy.
MoFEMErrorCode calculateCrackArea(boost::shared_ptr< double > area_ptr)
boost::shared_ptr< VolumeElementForcesAndSourcesCore > elasticFeRhs
static int dynamicStep
const std::string elementVolumeName
static double dd_f_log_e(const double v)
static enum RotSelector rotSelector
MoFEMErrorCode addDebugModel(TS ts)
Add debug to model.
static enum RotSelector gradApproximator
MoFEMErrorCode getBc(boost::shared_ptr< BC > &bc_vec_ptr, const std::string block_name, const int nb_attributes)
CommInterface::EntitiesPetscVector vertexExchange
boost::shared_ptr< BcRotVec > bcSpatialRotationVecPtr
boost::shared_ptr< Range > maxMovedFaces
static PetscBool dynamicRelaxation
const std::string spatialH1Disp
MoFEMErrorCode solveElastic(TS ts, Vec x)
static double d_f_log(const double v)
boost::shared_ptr< NormalDisplacementBcVec > bcSpatialNormalDisplacementVecPtr
static double crackingStartTime
static enum MaterialModel materialModel
MoFEMErrorCode getOptions()
const std::string piolaStress
MoFEMErrorCode setElasticElementToTs(DM dm)
static double inv_d_f_log_e(const double v)
MoFEMErrorCode setFaceInterfaceOps(const bool add_elastic, const bool add_material, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_rhs, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_lhs)
MoFEMErrorCode gettingNorms()
[Getting norms]
boost::shared_ptr< Range > interfaceFaces
std::vector< std::string > a00FieldList
MoFEMErrorCode setVolumeElementOps(const int tag, const bool add_elastic, const bool add_material, boost::shared_ptr< VolumeElementForcesAndSourcesCore > &fe_rhs, boost::shared_ptr< VolumeElementForcesAndSourcesCore > &fe_lhs)
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Getting interface of core database.
const std::string bubbleField
boost::shared_ptr< AnalyticalDisplacementBcVec > bcSpatialAnalyticalDisplacementVecPtr
SmartPetscObj< DM > dmMaterial
Material problem.
MoFEMErrorCode calculateOrientation(const int tag, bool set_orientation)
static double inv_f_log(const double v)
boost::shared_ptr< VolumeElementForcesAndSourcesCore > elasticFeLhs
static PetscBool noStretch
MoFEMErrorCode setNewFrontCoordinates()
boost::shared_ptr< ParentFiniteElementAdjacencyFunctionSkeleton< 2 > > parentAdjSkeletonFunctionDim2
static double crackingAddTime
static double exponentBase
static double dd_f_linear(const double v)
MoFEMErrorCode setFaceElementOps(const bool add_elastic, const bool add_material, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_rhs, boost::shared_ptr< FaceElementForcesAndSourcesCore > &fe_lhs)
MoFEMErrorCode projectGeometry(const EntityHandle meshset=0, double time=0)
MoFEMErrorCode postProcessSkeletonResults(const int tag, const std::string file, Vec f_residual=PETSC_NULLPTR, std::vector< Tag > tags_to_transfer={})
boost::shared_ptr< AnalyticalExprPython > AnalyticalExprPythonPtr
static double crackingAtol
Cracking absolute tolerance.
boost::shared_ptr< Range > skeletonFaces
static double crackingRtol
Cracking relative tolerance.
boost::shared_ptr< PhysicalEquations > physicalEquations
const std::string rotAxis
static double inv_d_f_linear(const double v)
BitRefLevel bitAdjParentMask
bit ref level for parent parent
MoFEMErrorCode solveDynamicRelaxation(TS ts, Vec x, int start_step, double start_time)
Solve problem using dynamic relaxation method.
static double inv_dd_f_log(const double v)
const std::string contactDisp
static std::string internalStressTagName
static enum SymmetrySelector symmetrySelector
CommInterface::EntitiesPetscVector edgeExchange
SmartPetscObj< DM > dmPrjSpatial
Projection spatial displacement.
static boost::function< double(const double)> f
boost::shared_ptr< BcDispVec > bcSpatialDispVecPtr
boost::shared_ptr< ForcesAndSourcesCore > contactTreeRhs
Make a contact tree.
const std::string skinElement
static PetscBool internalStressVoigt
static double inv_dd_f_linear(const double v)
MoFEMErrorCode getSpatialTractionFreeBc(const EntityHandle meshset=0)
static double inv_dd_f_log_e(const double v)
MoFEMErrorCode getExternalStrain()
MoFEMErrorCode getSpatialTractionBc()
static PetscBool setSingularity
virtual ~EshelbianCore()
static double d_f_log_e(const double v)
boost::shared_ptr< AnalyticalTractionBcVec > bcSpatialAnalyticalTractionVecPtr
static double f_log_e_quadratic(const double v)
MoFEMErrorCode addCrackSurfaces(const bool debug=false)
MoFEMErrorCode addDMs(const BitRefLevel bit=BitRefLevel().set(0), const EntityHandle meshset=0)
MoFEMErrorCode solveCohesiveCrackGrowth(TS ts, Vec x, int start_step, double start_time)
Solve cohesive crack growth problem.
MoFEMErrorCode getSpatialDispBc()
[Getting norms]
BitRefLevel bitAdjParent
bit ref level for parent
MoFEMErrorCode setContactElementRhsOps(boost::shared_ptr< ForcesAndSourcesCore > &fe_contact_tree)
static PetscBool interfaceCrack
static double d_f_log_e_quadratic(const double v)
CommInterface::EntitiesPetscVector volumeExchange
MoFEMErrorCode saveOrgCoords()
const std::string naturalBcElement
static boost::function< double(const double)> dd_f
static double f_log_e(const double v)
static constexpr double v_min
static int addCrackMeshsetId
static double inv_f_log_e(const double v)
MoFEMErrorCode createExchangeVectors(Sev sev)
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
boost::shared_ptr< Range > crackFaces
static boost::function< double(const double)> d_f
boost::shared_ptr< Range > frontVertices
static enum EnergyReleaseSelector energyReleaseSelector
MoFEMErrorCode addMaterial_HMHNeohookean(const int tape, const double c10, const double K)
static boost::function< double(const double)> inv_d_f
boost::shared_ptr< PressureBcVec > bcSpatialPressureVecPtr
static double d_f_linear(const double v)
const std::string hybridSpatialDisp
SmartPetscObj< Vec > solTSStep
static double f_log(const double v)
CommInterface::EntitiesPetscVector faceExchange
SmartPetscObj< DM > dmElastic
Elastic problem.
boost::shared_ptr< Range > frontEdges
static boost::function< double(const double)> inv_f
const std::string stretchTensor
BitRefLevel bitAdjEntMask
bit ref level for parent parent
static double f_linear(const double v)
const std::string contactElement
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.