v0.15.0
Loading...
Searching...
No Matches
TSElasticPostStep Struct Reference
Collaboration diagram for 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 9 of file TSElasticPostStep.cpp.

Member Function Documentation

◆ postStepDestroy()

static MoFEMErrorCode TSElasticPostStep::postStepDestroy ( )
inlinestatic
Examples
EshelbianPlasticity.cpp.

Definition at line 85 of file TSElasticPostStep.cpp.

85 {
87 prjKsp.reset();
88 prjD.reset();
89 prjF.reset();
90 prjDM.reset();
91 preProcRhs.reset();
93 };
#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 SmartPetscObj< Vec > prjD
static SmartPetscObj< Vec > prjF
static SmartPetscObj< KSP > prjKsp
static boost::shared_ptr< FEMethod > preProcRhs
static SmartPetscObj< DM > prjDM

◆ postStepFun()

static MoFEMErrorCode TSElasticPostStep::postStepFun ( TS ts)
inlinestatic
Examples
EshelbianPlasticity.cpp.

Definition at line 223 of file TSElasticPostStep.cpp.

223 {
225
226 double time;
227 CHKERR TSGetTime(ts, &time);
228
229 MOFEM_LOG("EP", Sev::inform) << "Solve H1 post-step";
230 CHKERR VecZeroEntries(prjF);
231 CHKERR KSPSolve(prjKsp, prjF, prjD);
232 CHKERR VecGhostUpdateBegin(prjD, INSERT_VALUES, SCATTER_FORWARD);
233 CHKERR VecGhostUpdateEnd(prjD, INSERT_VALUES, SCATTER_FORWARD);
234 CHKERR DMoFEMMeshToLocalVector(prjDM, prjD, INSERT_VALUES, SCATTER_REVERSE);
235
238 case GRIFFITH_FORCE:
240 MOFEM_LOG("EP", Sev::inform) << "Calculate Griffith force";
242 break;
243 default:
244 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
245 "Energy release selector not handled");
246 };
248 MOFEM_LOG("EP", Sev::inform) << "Calculate faces orientation";
249
251 }
252 }
253
254
255 // double time;
256 // CHKERR TSGetTime(ts, &time);
257 // preProcRhs->ts_t = time;
258 // CHKERR DMoFEMPreProcessFiniteElements(prjDM, preProcRhs.get());
259 CHKERR VecZeroEntries(epPtr->solTSStep);
260
262 };
@ 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
static EshelbianCore * epPtr

◆ postStepInitialise()

static MoFEMErrorCode TSElasticPostStep::postStepInitialise ( EshelbianCore * ep_ptr)
inlinestatic
Examples
EshelbianPlasticity.cpp.

Definition at line 11 of file TSElasticPostStep.cpp.

11 {
13
14 epPtr = ep_ptr;
15
16 auto create_post_step_ksp = [&]() {
17 auto ksp = createKSP(epPtr->mField.get_comm());
18
19 auto set_up = [&]() {
21 using DomainEle = VolumeElementForcesAndSourcesCore;
22 using DomainEleOp = DomainEle::UserDataOperator;
23 using OpDomainMass = FormsIntegrators<DomainEleOp>::Assembly<
25 using OpRhs = FormsIntegrators<DomainEleOp>::Assembly<
27 auto fe_lhs = boost::make_shared<DomainEle>(ep_ptr->mField);
28 auto fe_rhs = boost::make_shared<DomainEle>(ep_ptr->mField);
29
30 fe_lhs->getUserPolynomialBase() =
31 boost::shared_ptr<BaseFunction>(new CGGUserPolynomialBase());
32 fe_rhs->getUserPolynomialBase() =
33 boost::shared_ptr<BaseFunction>(new CGGUserPolynomialBase());
34 CHKERR
36 fe_lhs->getOpPtrVector(), {L2}, ep_ptr->materialH1Positions,
37 ep_ptr->frontAdjEdges);
38 CHKERR
40 fe_rhs->getOpPtrVector(), {L2}, ep_ptr->materialH1Positions,
41 ep_ptr->frontAdjEdges);
42
43 fe_lhs->getOpPtrVector().push_back(
44 new OpDomainMass(ep_ptr->spatialH1Disp, ep_ptr->spatialH1Disp));
45 auto w_ptr = boost::make_shared<MatrixDouble>();
46 fe_rhs->getOpPtrVector().push_back(
47 new OpCalculateVectorFieldValues<3>(ep_ptr->spatialL2Disp, w_ptr));
48 fe_rhs->getOpPtrVector().push_back(
49 new OpRhs(ep_ptr->spatialH1Disp, w_ptr));
50
52 ep_ptr->elementVolumeName, fe_lhs,
53 nullptr, nullptr);
55 ep_ptr->elementVolumeName, fe_rhs,
56 nullptr, nullptr);
57
58 // preProcRhs = boost::make_shared<FEMethod>();
59 // struct MinusOne : public ScalingMethod {
60 // double getScale(const double time) { return -time; }
61 // };
62 // preProcRhs->preProcessHook = EssentialPreProc<DisplacementCubitBcData>(
63 // ep_ptr->mField, preProcRhs, {boost::make_shared<MinusOne>()});
64
65 CHKERR KSPAppendOptionsPrefix(ksp, "prjspatial_");
66 CHKERR KSPSetFromOptions(ksp);
67 CHKERR KSPSetDM(ksp, ep_ptr->dmPrjSpatial);
68 CHKERR KSPSetUp(ksp);
70 };
71
72 CHK_THROW_MESSAGE(set_up(), "set up");
73
74 return ksp;
75 };
76
77 prjKsp = create_post_step_ksp();
80 prjDM = ep_ptr->dmPrjSpatial;
81
83 };
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
[Define dimension]
#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 compute operator for KSP solver via sub-matrix and IS.
Definition DMMoFEM.cpp:627
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1102
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
auto createKSP(MPI_Comm comm)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
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

◆ preStepFun()

static MoFEMErrorCode TSElasticPostStep::preStepFun ( TS ts)
inlinestatic
Examples
EshelbianPlasticity.cpp.

Definition at line 95 of file TSElasticPostStep.cpp.

95 {
97 MOFEM_LOG("EP", Sev::inform) << "Pre step";
98
99 double time;
100 CHKERR TSGetTime(ts, &time);
101
103
104 auto debug_crack = [&]() {
106
107 PetscBool debug_crack_mesh = PETSC_FALSE;
108 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-debug_crack_mesh",
109 &debug_crack_mesh, PETSC_NULLPTR);
110 if (debug_crack_mesh) {
111
112 auto get_meshsets_from_block =
113 [](MoFEM::Interface &m_field, std::string block_name, int dim,
114 std::vector<EntityHandle> r) {
115
116 auto mesh_mng = m_field.getInterface<MeshsetsManager>();
117 auto bcs = mesh_mng->getCubitMeshsetPtr(
118
119 std::regex((boost::format("%s(.*)") % block_name).str())
120
121 );
122
123 for (auto bc : bcs) {
124 r.push_back(bc->getMeshset());
125 }
126
127 return r;
128 };
129
130 if (epPtr->mField.get_comm_rank() == 0) {
131
132 if (epPtr->maxMovedFaces) {
133 CHKERR epPtr->mField.getInterface<MeshsetsManager>()
134 ->addEntitiesToMeshset(BLOCKSET, epPtr->addCrackMeshsetId,
136 }
137
138 auto meshset_ptr = get_temp_meshset_ptr(epPtr->mField.get_moab());
139 Range tets;
140 CHKERR epPtr->mField.get_moab().get_entities_by_dimension(
141 *meshset_ptr, 3, tets);
142 CHKERR epPtr->mField.get_moab().add_entities(*meshset_ptr, tets);
143
144 std::vector<EntityHandle> meshsets;
145 meshsets.push_back(*meshset_ptr);
146
147 meshsets =
148 get_meshsets_from_block(epPtr->mField, "CRACK", 2, meshsets);
149 meshsets =
150 get_meshsets_from_block(epPtr->mField, "FRONT", 1, meshsets);
151 meshsets =
152 get_meshsets_from_block(epPtr->mField, "EDGE", 1, meshsets);
153
154 int time_step = 0;
155 CHKERR TSGetStepNumber(ts, &time_step);
156 std::string file_name =
157 "crack_meshsets_" + std::to_string(time_step) + ".h5m";
158 CHKERR epPtr->mField.get_moab().write_file(
159 file_name.c_str(), "MOAB", nullptr, meshsets.data());
160 }
161 }
162
164 };
165
167 CHKERR debug_crack();
170 }
171
172 Vec T;
173 CHKERR TSGetSolution(ts, &T);
174
175 auto zero_filled = [&](auto name) {
177 auto is_mng = epPtr->mField.getInterface<ISManager>();
178 SmartPetscObj<IS> is;
179 CHKERR is_mng->isCreateProblemFieldAndRankLocal(
180 "ELASTIC_PROBLEM", ROW, name, 0, MAX_DOFS_ON_ENTITY, is);
181 const int *index_ptr;
182 CHKERR ISGetIndices(is, &index_ptr);
183 int size;
184 CHKERR ISGetLocalSize(is, &size);
185 double *a;
186 CHKERR VecGetArray(T, &a);
187 for (auto i = 0; i != size; i++) {
188 a[index_ptr[i]] = 0;
189 }
190 CHKERR VecRestoreArray(T, &a);
191 CHKERR ISRestoreIndices(is, &index_ptr);
193 };
194
197 break;
198 case LARGE_ROT:
199 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (large) fields";
200 CHKERR zero_filled(epPtr->stretchTensor);
201 CHKERR zero_filled(epPtr->rotAxis);
202 break;
203 case MODERATE_ROT:
204 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (moderate) fields";
205 CHKERR zero_filled(epPtr->stretchTensor);
206 CHKERR zero_filled(epPtr->rotAxis);
207 break;
208 case SMALL_ROT:
209 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (small) fields";
210 CHKERR zero_filled(epPtr->stretchTensor);
211 CHKERR zero_filled(epPtr->rotAxis);
212 break;
213 }
214
215 CHKERR VecCopy(T, epPtr->solTSStep);
216 CHKERR VecGhostUpdateBegin(epPtr->solTSStep, INSERT_VALUES,
217 SCATTER_FORWARD);
218 CHKERR VecGhostUpdateEnd(epPtr->solTSStep, INSERT_VALUES, SCATTER_FORWARD);
219
221 }
constexpr double a
@ ROW
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
@ BLOCKSET
FTensor::Index< 'i', SPACE_DIM > i
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:8
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.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

Member Data Documentation

◆ epPtr

EshelbianCore * TSElasticPostStep::epPtr
staticprivate

Definition at line 269 of file TSElasticPostStep.cpp.

◆ preProcRhs

boost::shared_ptr< FEMethod > TSElasticPostStep::preProcRhs
staticprivate

Definition at line 271 of file TSElasticPostStep.cpp.

◆ prjD

SmartPetscObj< Vec > TSElasticPostStep::prjD
staticprivate

Definition at line 266 of file TSElasticPostStep.cpp.

◆ prjDM

SmartPetscObj< DM > TSElasticPostStep::prjDM
staticprivate

Definition at line 268 of file TSElasticPostStep.cpp.

◆ prjF

SmartPetscObj< Vec > TSElasticPostStep::prjF
staticprivate

Definition at line 267 of file TSElasticPostStep.cpp.

◆ prjKsp

SmartPetscObj< KSP > TSElasticPostStep::prjKsp
staticprivate

Definition at line 265 of file TSElasticPostStep.cpp.


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