v0.15.4
Loading...
Searching...
No Matches
Static Public Member Functions | Static Private Attributes | List of all members
EshelbianPlasticity::TSElasticPostStep Struct Reference

#include "users_modules/eshelbian_plasticity/src/TSElasticPostStep.hpp"

Collaboration diagram for EshelbianPlasticity::TSElasticPostStep:
[legend]

Static Public Member Functions

static MoFEMErrorCode postStepInitialise (EshelbianCore *ep_ptr)
 
static MoFEMErrorCode postStepDestroy ()
 
static MoFEMErrorCode preStepFun (TS ts)
 
static MoFEMErrorCode postStepFun (TS ts)
 

Static Private Attributes

static SmartPetscObj< KSP > prjKsp
 
static SmartPetscObj< Vec > prjD
 
static SmartPetscObj< Vec > prjF
 
static SmartPetscObj< DM > prjDM
 
static EshelbianCoreepPtr
 
static boost::shared_ptr< FEMethodpreProcRhs
 

Detailed Description

Definition at line 17 of file TSElasticPostStep.hpp.

Member Function Documentation

◆ postStepDestroy()

MoFEMErrorCode EshelbianPlasticity::TSElasticPostStep::postStepDestroy ( )
static

Definition at line 96 of file TSElasticPostStep.cpp.

96 {
98 prjKsp.reset();
99 prjD.reset();
100 prjF.reset();
101 prjDM.reset();
102 preProcRhs.reset();
104};
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
static boost::shared_ptr< FEMethod > preProcRhs

◆ postStepFun()

MoFEMErrorCode EshelbianPlasticity::TSElasticPostStep::postStepFun ( TS  ts)
static

Definition at line 233 of file TSElasticPostStep.cpp.

233 {
235
236 double time;
237 CHKERR TSGetTime(ts, &time);
238
239 MOFEM_LOG("EP", Sev::inform) << "Solve H1 post-step";
240 CHKERR VecZeroEntries(prjF);
241 CHKERR KSPSolve(prjKsp, prjF, prjD);
242 CHKERR VecGhostUpdateBegin(prjD, INSERT_VALUES, SCATTER_FORWARD);
243 CHKERR VecGhostUpdateEnd(prjD, INSERT_VALUES, SCATTER_FORWARD);
244 CHKERR DMoFEMMeshToLocalVector(prjDM, prjD, INSERT_VALUES, SCATTER_REVERSE);
245
248 case GRIFFITH_FORCE:
250 MOFEM_LOG("EP", Sev::inform) << "Calculate Griffith force";
252 break;
253 default:
254 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
255 "Energy release selector not handled");
256 };
258 MOFEM_LOG("EP", Sev::inform) << "Calculate faces orientation";
259
261 }
262 }
263
264 CHKERR VecZeroEntries(epPtr->solTSStep);
265
267};
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define CHKERR
Inline error check.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
#define MOFEM_LOG(channel, severity)
Log.
static PetscBool crackingOn
MoFEMErrorCode calculateFaceMaterialForce(const int tag, TS ts)
static double crackingStartTime
MoFEMErrorCode calculateOrientation(const int tag, bool set_orientation)
static enum EnergyReleaseSelector energyReleaseSelector
SmartPetscObj< Vec > solTSStep

◆ postStepInitialise()

MoFEMErrorCode EshelbianPlasticity::TSElasticPostStep::postStepInitialise ( EshelbianCore ep_ptr)
static

Definition at line 22 of file TSElasticPostStep.cpp.

22 {
24
25 epPtr = ep_ptr;
26
27 auto create_post_step_ksp = [&]() {
28 auto ksp = createKSP(epPtr->mField.get_comm());
29
30 auto set_up = [&]() {
33 using DomainEleOp = DomainEle::UserDataOperator;
37 GAUSS>::OpBaseTimesVector<1, 3, 1>;
38 auto fe_lhs = boost::make_shared<DomainEle>(ep_ptr->mField);
39 auto fe_rhs = boost::make_shared<DomainEle>(ep_ptr->mField);
40
41 fe_lhs->getUserPolynomialBase() =
42 boost::shared_ptr<BaseFunction>(new CGGUserPolynomialBase());
43 fe_rhs->getUserPolynomialBase() =
44 boost::shared_ptr<BaseFunction>(new CGGUserPolynomialBase());
45 CHKERR
46 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
47 fe_lhs->getOpPtrVector(), {L2}, ep_ptr->materialH1Positions,
48 ep_ptr->frontAdjEdges);
49 CHKERR
50 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
51 fe_rhs->getOpPtrVector(), {L2}, ep_ptr->materialH1Positions,
52 ep_ptr->frontAdjEdges);
53
54 fe_lhs->getOpPtrVector().push_back(
55 new OpDomainMass(ep_ptr->spatialH1Disp, ep_ptr->spatialH1Disp));
56 auto w_ptr = boost::make_shared<MatrixDouble>();
57 fe_rhs->getOpPtrVector().push_back(
59 fe_rhs->getOpPtrVector().push_back(
60 new OpRhs(ep_ptr->spatialH1Disp, w_ptr));
61
63 ep_ptr->elementVolumeName, fe_lhs,
64 nullptr, nullptr);
66 ep_ptr->elementVolumeName, fe_rhs, nullptr,
67 nullptr);
68
69 // preProcRhs = boost::make_shared<FEMethod>();
70 // struct MinusOne : public ScalingMethod {
71 // double getScale(const double time) { return -time; }
72 // };
73 // preProcRhs->preProcessHook = EssentialPreProc<DisplacementCubitBcData>(
74 // ep_ptr->mField, preProcRhs, {boost::make_shared<MinusOne>()});
75
76 CHKERR KSPAppendOptionsPrefix(ksp, "prjspatial_");
77 CHKERR KSPSetFromOptions(ksp);
78 CHKERR KSPSetDM(ksp, ep_ptr->dmPrjSpatial);
79 CHKERR KSPSetUp(ksp);
81 };
82
83 CHK_THROW_MESSAGE(set_up(), "set up");
84
85 return ksp;
86 };
87
88 prjKsp = create_post_step_ksp();
91 prjDM = ep_ptr->dmPrjSpatial;
92
94};
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set KSP right hand side evaluation function
Definition DMMoFEM.cpp:627
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Definition DMMoFEM.cpp:668
@ GAUSS
Gaussian quadrature integration.
@ PETSC
Standard PETSc assembly.
auto createKSP(MPI_Comm comm)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
CGG User Polynomial Base.
boost::shared_ptr< Range > frontAdjEdges
MoFEM::Interface & mField
const std::string spatialL2Disp
const std::string materialH1Positions
const std::string elementVolumeName
const std::string spatialH1Disp
SmartPetscObj< DM > dmPrjSpatial
Projection spatial displacement.
virtual MPI_Comm & get_comm() const =0
Specialization for double precision vector field values calculation.

◆ preStepFun()

MoFEMErrorCode EshelbianPlasticity::TSElasticPostStep::preStepFun ( TS  ts)
static

Definition at line 106 of file TSElasticPostStep.cpp.

106 {
108 MOFEM_LOG("EP", Sev::inform) << "Pre step";
109
110 double time;
111 CHKERR TSGetTime(ts, &time);
112
114
115 auto debug_crack = [&]() {
117
118 PetscBool debug_crack_mesh = PETSC_FALSE;
119 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-debug_crack_mesh",
120 &debug_crack_mesh, PETSC_NULLPTR);
121 if (debug_crack_mesh) {
122
123 auto get_meshsets_from_block = [](MoFEM::Interface &m_field,
124 std::string block_name, int dim,
125 std::vector<EntityHandle> r) {
126 auto mesh_mng = m_field.getInterface<MeshsetsManager>();
127 auto bcs = mesh_mng->getCubitMeshsetPtr(
128
129 std::regex((boost::format("%s(.*)") % block_name).str())
130
131 );
132
133 for (auto bc : bcs) {
134 r.push_back(bc->getMeshset());
135 }
136
137 return r;
138 };
139
140 if (epPtr->mField.get_comm_rank() == 0) {
141
142 if (epPtr->maxMovedFaces) {
144 ->addEntitiesToMeshset(BLOCKSET, epPtr->addCrackMeshsetId,
146 }
147
148 auto meshset_ptr = get_temp_meshset_ptr(epPtr->mField.get_moab());
149 Range tets;
150 CHKERR epPtr->mField.get_moab().get_entities_by_dimension(
151 *meshset_ptr, 3, tets);
152 CHKERR epPtr->mField.get_moab().add_entities(*meshset_ptr, tets);
153
154 std::vector<EntityHandle> meshsets;
155 meshsets.push_back(*meshset_ptr);
156
157 meshsets =
158 get_meshsets_from_block(epPtr->mField, "CRACK", 2, meshsets);
159 meshsets =
160 get_meshsets_from_block(epPtr->mField, "FRONT", 1, meshsets);
161 meshsets =
162 get_meshsets_from_block(epPtr->mField, "EDGE", 1, meshsets);
163
164 int time_step = 0;
165 CHKERR TSGetStepNumber(ts, &time_step);
166 std::string file_name =
167 "crack_meshsets_" + std::to_string(time_step) + ".h5m";
168 CHKERR epPtr->mField.get_moab().write_file(file_name.c_str(), "MOAB",
169 nullptr, meshsets.data());
170 }
171 }
172
174 };
175
177 CHKERR debug_crack();
181 }
182
183 Vec T;
184 CHKERR TSGetSolution(ts, &T);
185
186 auto zero_filled = [&](auto name) {
188 auto is_mng = epPtr->mField.getInterface<ISManager>();
190 CHKERR is_mng->isCreateProblemFieldAndRankLocal(
191 "ELASTIC_PROBLEM", ROW, name, 0, MAX_DOFS_ON_ENTITY, is);
192 const int *index_ptr;
193 CHKERR ISGetIndices(is, &index_ptr);
194 int size;
195 CHKERR ISGetLocalSize(is, &size);
196 double *a;
197 CHKERR VecGetArray(T, &a);
198 for (auto i = 0; i != size; i++) {
199 a[index_ptr[i]] = 0;
200 }
201 CHKERR VecRestoreArray(T, &a);
202 CHKERR ISRestoreIndices(is, &index_ptr);
204 };
205
208 break;
209 case LARGE_ROT:
210 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (large) fields";
211 CHKERR zero_filled(epPtr->stretchTensor);
212 CHKERR zero_filled(epPtr->rotAxis);
213 break;
214 case MODERATE_ROT:
215 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (moderate) fields";
216 CHKERR zero_filled(epPtr->stretchTensor);
217 CHKERR zero_filled(epPtr->rotAxis);
218 break;
219 case SMALL_ROT:
220 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (small) fields";
221 CHKERR zero_filled(epPtr->stretchTensor);
222 CHKERR zero_filled(epPtr->rotAxis);
223 break;
224 }
225
226 CHKERR VecCopy(T, epPtr->solTSStep);
227 CHKERR VecGhostUpdateBegin(epPtr->solTSStep, INSERT_VALUES, SCATTER_FORWARD);
228 CHKERR VecGhostUpdateEnd(epPtr->solTSStep, INSERT_VALUES, SCATTER_FORWARD);
229
231}
constexpr double a
@ ROW
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
@ BLOCKSET
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
const FTensor::Tensor2< T, Dim, Dim > Vec
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
int r
Definition sdf.py:205
MoFEMErrorCode projectInternalStress(const EntityHandle meshset=0)
MoFEMErrorCode projectGeometry(const EntityHandle meshset=0)
static enum RotSelector gradApproximator
boost::shared_ptr< Range > maxMovedFaces
MoFEMErrorCode setNewFrontCoordinates()
const std::string rotAxis
MoFEMErrorCode addCrackSurfaces(const bool debug=false)
static int addCrackMeshsetId
const std::string stretchTensor
virtual moab::Interface & get_moab()=0
virtual int get_comm_rank() const =0
Deprecated interface functions.
Section manager is used to create indexes and sections.
Definition ISManager.hpp:23
Interface for managing meshsets containing materials and boundary conditions.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

Member Data Documentation

◆ epPtr

EshelbianCore* EshelbianPlasticity::TSElasticPostStep::epPtr
inlinestaticprivate

Definition at line 29 of file TSElasticPostStep.hpp.

◆ preProcRhs

boost::shared_ptr<FEMethod> EshelbianPlasticity::TSElasticPostStep::preProcRhs
inlinestaticprivate

Definition at line 30 of file TSElasticPostStep.hpp.

◆ prjD

SmartPetscObj<Vec> EshelbianPlasticity::TSElasticPostStep::prjD
inlinestaticprivate

Definition at line 26 of file TSElasticPostStep.hpp.

◆ prjDM

SmartPetscObj<DM> EshelbianPlasticity::TSElasticPostStep::prjDM
inlinestaticprivate

Definition at line 28 of file TSElasticPostStep.hpp.

◆ prjF

SmartPetscObj<Vec> EshelbianPlasticity::TSElasticPostStep::prjF
inlinestaticprivate

Definition at line 27 of file TSElasticPostStep.hpp.

◆ prjKsp

SmartPetscObj<KSP> EshelbianPlasticity::TSElasticPostStep::prjKsp
inlinestaticprivate

Definition at line 25 of file TSElasticPostStep.hpp.


The documentation for this struct was generated from the following files: