v0.14.0
Loading...
Searching...
No Matches
TSElasticPostStep.cpp
Go to the documentation of this file.
1/** @file TSElasticPostStep.cpp
2 * @brief
3 * @date 2023-05-13
4 *
5 * @license{This project is released under the MIT License.}
6 *
7 */
8
10
11 static MoFEMErrorCode postStepInitialise(EshelbianCore *ep_ptr) {
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
35 fe_lhs->getOpPtrVector().push_back(
36 new OpDomainMass(ep_ptr->spatialH1Disp, ep_ptr->spatialH1Disp));
37 auto w_ptr = boost::make_shared<MatrixDouble>();
38 fe_rhs->getOpPtrVector().push_back(
39 new OpCalculateVectorFieldValues<3>(ep_ptr->spatialL2Disp, w_ptr));
40 fe_rhs->getOpPtrVector().push_back(
41 new OpRhs(ep_ptr->spatialH1Disp, w_ptr));
42
43 CHKERR DMMoFEMKSPSetComputeOperators(ep_ptr->dmPrjSpatial,
44 ep_ptr->elementVolumeName, fe_lhs,
45 nullptr, nullptr);
46 CHKERR DMMoFEMKSPSetComputeRHS(ep_ptr->dmPrjSpatial,
47 ep_ptr->elementVolumeName, fe_rhs,
48 nullptr, nullptr);
49
50 // preProcRhs = boost::make_shared<FEMethod>();
51 // struct MinusOne : public ScalingMethod {
52 // double getScale(const double time) { return -time; }
53 // };
54 // preProcRhs->preProcessHook = EssentialPreProc<DisplacementCubitBcData>(
55 // ep_ptr->mField, preProcRhs, {boost::make_shared<MinusOne>()});
56
57 CHKERR KSPAppendOptionsPrefix(ksp, "prjspatial_");
58 CHKERR KSPSetFromOptions(ksp);
59 CHKERR KSPSetDM(ksp, ep_ptr->dmPrjSpatial);
60 CHKERR KSPSetUp(ksp);
62 };
63
64 CHK_THROW_MESSAGE(set_up(), "set up");
65
66 return ksp;
67 };
68
69 prjKsp = create_post_step_ksp();
70 prjD = createDMVector(ep_ptr->dmPrjSpatial);
71 prjF = vectorDuplicate(prjD);
72 prjDM = ep_ptr->dmPrjSpatial;
73
75 };
76
77 static MoFEMErrorCode postStepDestroy() {
79 prjKsp.reset();
80 prjD.reset();
81 prjF.reset();
82 prjDM.reset();
83 preProcRhs.reset();
85 };
86
87 static MoFEMErrorCode preStepFun(TS ts) {
89 MOFEM_LOG("EP", Sev::inform) << "Pre step";
90
91 Vec T;
92 CHKERR TSGetSolution(ts, &T);
93
94 auto zero_filed = [&](auto name) {
96 auto is_mng = epPtr->getInterface<ISManager>();
97 SmartPetscObj<IS> is;
98 CHKERR is_mng->isCreateProblemFieldAndRankLocal(
99 "ELASTIC_PROBLEM", ROW, epPtr->rotAxis, 0, MAX_DOFS_ON_ENTITY, is);
100 const int *index_ptr;
101 CHKERR ISGetIndices(is, &index_ptr);
102 int size;
103 CHKERR ISGetLocalSize(is, &size);
104 double *a;
105 CHKERR VecGetArray(T, &a);
106 for (auto i = 0; i != size; i++) {
107 a[index_ptr[i]] = 0;
108 }
109 CHKERR VecRestoreArray(T, &a);
110 CHKERR ISRestoreIndices(is, &index_ptr);
112 };
113
114 switch (EshelbianCore::gradApperoximator) {
115 case LARGE_ROT:
116 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (large) fields";
117 CHKERR zero_filed(epPtr->stretchTensor);
118 CHKERR zero_filed(epPtr->rotAxis);
119 break;
120 case MODERATE_ROT:
121 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (moderate) fields";
122 CHKERR zero_filed(epPtr->stretchTensor);
123 CHKERR zero_filed(epPtr->rotAxis);
124 break;
125 case SMALL_ROT:
126 MOFEM_LOG("EP", Sev::verbose) << "Zeroing (small) fields";
127 CHKERR zero_filed(epPtr->rotAxis);
128 break;
129 }
130
132 }
133
134 static MoFEMErrorCode postStepFun(TS ts) {
136
137 MOFEM_LOG("EP", Sev::inform) << "Post step";
138 CHKERR VecZeroEntries(prjF);
139 CHKERR KSPSolve(prjKsp, prjF, prjD);
140 CHKERR VecGhostUpdateBegin(prjD, INSERT_VALUES, SCATTER_FORWARD);
141 CHKERR VecGhostUpdateEnd(prjD, INSERT_VALUES, SCATTER_FORWARD);
142 CHKERR DMoFEMMeshToLocalVector(prjDM, prjD, INSERT_VALUES, SCATTER_REVERSE);
143 double time;
144 CHKERR TSGetTime(ts, &time);
145 // preProcRhs->ts_t = time;
146 // CHKERR DMoFEMPreProcessFiniteElements(prjDM, preProcRhs.get());
148 };
149
150private:
151 static SmartPetscObj<KSP> prjKsp;
152 static SmartPetscObj<Vec> prjD;
153 static SmartPetscObj<Vec> prjF;
154 static SmartPetscObj<DM> prjDM;
155 static EshelbianCore *epPtr;
156
157 static boost::shared_ptr<FEMethod> preProcRhs;
158};
159
160SmartPetscObj<KSP> TSElasticPostStep::prjKsp;
161SmartPetscObj<Vec> TSElasticPostStep::prjD;
162SmartPetscObj<Vec> TSElasticPostStep::prjF;
163SmartPetscObj<DM> TSElasticPostStep::prjDM;
164EshelbianCore *TSElasticPostStep::epPtr;
165boost::shared_ptr<FEMethod> TSElasticPostStep::preProcRhs;
constexpr double a
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
@ ROW
Definition: definitions.h:123
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
Definition: definitions.h:236
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:595
#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
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
FTensor::Index< 'i', SPACE_DIM > i
const double T
static SmartPetscObj< Vec > prjD
static EshelbianCore * epPtr
static MoFEMErrorCode postStepInitialise(EshelbianCore *ep_ptr)
static MoFEMErrorCode postStepDestroy()
static SmartPetscObj< Vec > prjF
static MoFEMErrorCode preStepFun(TS ts)
static SmartPetscObj< KSP > prjKsp
static boost::shared_ptr< FEMethod > preProcRhs
static MoFEMErrorCode postStepFun(TS ts)
static SmartPetscObj< DM > prjDM