v0.14.0
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 
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;
23  using OpDomainMass = FormsIntegrators<DomainEleOp>::Assembly<
25  using OpRhs = FormsIntegrators<DomainEleOp>::Assembly<
26  PETSC>::LinearForm<GAUSS>::OpBaseTimesVector<1, 3, 1>;
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
35  EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
36  fe_lhs->getOpPtrVector(), {L2}, ep_ptr->materialH1Positions,
37  ep_ptr->frontAdjEdges);
38  CHKERR
39  EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
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();
78  prjD = createDMVector(ep_ptr->dmPrjSpatial);
80  prjDM = ep_ptr->dmPrjSpatial;
81 
83  };
84 
87  prjKsp.reset();
88  prjD.reset();
89  prjF.reset();
90  prjDM.reset();
91  preProcRhs.reset();
93  };
94 
95  static MoFEMErrorCode preStepFun(TS ts) {
97  MOFEM_LOG("EP", Sev::inform) << "Pre step";
98 
103  }
104 
105  Vec T;
106  CHKERR TSGetSolution(ts, &T);
107 
108  auto zero_filled = [&](auto name) {
110  auto is_mng = epPtr->mField.getInterface<ISManager>();
111  SmartPetscObj<IS> is;
112  CHKERR is_mng->isCreateProblemFieldAndRankLocal(
113  "ELASTIC_PROBLEM", ROW, name, 0, MAX_DOFS_ON_ENTITY, is);
114  const int *index_ptr;
115  CHKERR ISGetIndices(is, &index_ptr);
116  int size;
117  CHKERR ISGetLocalSize(is, &size);
118  double *a;
119  CHKERR VecGetArray(T, &a);
120  for (auto i = 0; i != size; i++) {
121  a[index_ptr[i]] = 0;
122  }
123  CHKERR VecRestoreArray(T, &a);
124  CHKERR ISRestoreIndices(is, &index_ptr);
126  };
127 
129  case NO_H1_CONFIGURATION:
130  break;
131  case LARGE_ROT:
132  MOFEM_LOG("EP", Sev::verbose) << "Zeroing (large) fields";
133  CHKERR zero_filled(epPtr->stretchTensor);
134  CHKERR zero_filled(epPtr->rotAxis);
135  break;
136  case MODERATE_ROT:
137  MOFEM_LOG("EP", Sev::verbose) << "Zeroing (moderate) fields";
138  CHKERR zero_filled(epPtr->stretchTensor);
139  CHKERR zero_filled(epPtr->rotAxis);
140  break;
141  case SMALL_ROT:
142  MOFEM_LOG("EP", Sev::verbose) << "Zeroing (small) fields";
143  CHKERR zero_filled(epPtr->stretchTensor);
144  CHKERR zero_filled(epPtr->rotAxis);
145  break;
146  }
147 
148  CHKERR VecCopy(T, epPtr->solTSStep);
149  CHKERR VecGhostUpdateBegin(epPtr->solTSStep, INSERT_VALUES,
150  SCATTER_FORWARD);
151  CHKERR VecGhostUpdateEnd(epPtr->solTSStep, INSERT_VALUES, SCATTER_FORWARD);
152 
154  }
155 
156  static MoFEMErrorCode postStepFun(TS ts) {
159  MOFEM_LOG("EP", Sev::inform) << "Calculate faces energy";
161  MOFEM_LOG("EP", Sev::inform) << "Calculate faces orientation";
163  }
164 
165  MOFEM_LOG("EP", Sev::inform) << "Solve H1 post-step";
166  CHKERR VecZeroEntries(prjF);
167  CHKERR KSPSolve(prjKsp, prjF, prjD);
168  CHKERR VecGhostUpdateBegin(prjD, INSERT_VALUES, SCATTER_FORWARD);
169  CHKERR VecGhostUpdateEnd(prjD, INSERT_VALUES, SCATTER_FORWARD);
170  CHKERR DMoFEMMeshToLocalVector(prjDM, prjD, INSERT_VALUES, SCATTER_REVERSE);
171  // double time;
172  // CHKERR TSGetTime(ts, &time);
173  // preProcRhs->ts_t = time;
174  // CHKERR DMoFEMPreProcessFiniteElements(prjDM, preProcRhs.get());
175  CHKERR VecZeroEntries(epPtr->solTSStep);
176 
178  };
179 
180 private:
181  static SmartPetscObj<KSP> prjKsp;
182  static SmartPetscObj<Vec> prjD;
183  static SmartPetscObj<Vec> prjF;
184  static SmartPetscObj<DM> prjDM;
186 
187  static boost::shared_ptr<FEMethod> preProcRhs;
188 };
189 
190 SmartPetscObj<KSP> TSElasticPostStep::prjKsp;
191 SmartPetscObj<Vec> TSElasticPostStep::prjD;
192 SmartPetscObj<Vec> TSElasticPostStep::prjF;
193 SmartPetscObj<DM> TSElasticPostStep::prjDM;
195 boost::shared_ptr<FEMethod> TSElasticPostStep::preProcRhs;
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
EshelbianCore::setNewFrontCoordinates
MoFEMErrorCode setNewFrontCoordinates()
Definition: EshelbianPlasticity.cpp:4990
MoFEM::DMMoFEMKSPSetComputeRHS
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:637
EshelbianCore::spatialL2Disp
const std::string spatialL2Disp
Definition: EshelbianCore.hpp:102
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
TSElasticPostStep::postStepFun
static MoFEMErrorCode postStepFun(TS ts)
Definition: TSElasticPostStep.cpp:156
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
EshelbianCore::calculateOrientation
MoFEMErrorCode calculateOrientation(const int tag, bool set_orientation)
Definition: EshelbianPlasticity.cpp:4068
EshelbianCore::addCrackSurfaces
MoFEMErrorCode addCrackSurfaces()
Definition: EshelbianPlasticity.cpp:5028
MoFEM::PETSC
@ PETSC
Definition: FormsIntegrators.hpp:105
MoFEM::DMoFEMMeshToLocalVector
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:523
MoFEM::createKSP
auto createKSP(MPI_Comm comm)
Definition: PetscSmartObj.hpp:261
EshelbianPlasticity::NO_H1_CONFIGURATION
@ NO_H1_CONFIGURATION
Definition: EshelbianPlasticity.hpp:45
MoFEM::DMMoFEMKSPSetComputeOperators
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:678
OpDomainMass
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, FIELD_DIM > OpDomainMass
Definition: child_and_parent.cpp:53
TSElasticPostStep
Definition: TSElasticPostStep.cpp:9
ROW
@ ROW
Definition: definitions.h:136
EshelbianCore::gradApproximator
static enum RotSelector gradApproximator
Definition: EshelbianCore.hpp:16
EshelbianCore::calculateEnergyRelease
MoFEMErrorCode calculateEnergyRelease(const int tag, TS ts)
Definition: EshelbianPlasticity.cpp:3638
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1099
TSElasticPostStep::prjKsp
static SmartPetscObj< KSP > prjKsp
Definition: TSElasticPostStep.cpp:178
TSElasticPostStep::preStepFun
static MoFEMErrorCode preStepFun(TS ts)
Definition: TSElasticPostStep.cpp:95
a
constexpr double a
Definition: approx_sphere.cpp:30
EshelbianCore::rotAxis
const std::string rotAxis
Definition: EshelbianCore.hpp:110
EshelbianCore::materialH1Positions
const std::string materialH1Positions
Definition: EshelbianCore.hpp:105
EshelbianCore::elementVolumeName
const std::string elementVolumeName
Definition: EshelbianCore.hpp:113
EshelbianCore::spatialH1Disp
const std::string spatialH1Disp
Definition: EshelbianCore.hpp:104
EshelbianCore::mField
MoFEM::Interface & mField
Definition: EshelbianCore.hpp:84
TSElasticPostStep::epPtr
static EshelbianCore * epPtr
Definition: TSElasticPostStep.cpp:185
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
EshelbianCore::stretchTensor
const std::string stretchTensor
Definition: EshelbianCore.hpp:109
BiLinearForm
TSElasticPostStep::postStepDestroy
static MoFEMErrorCode postStepDestroy()
Definition: TSElasticPostStep.cpp:85
OpRhs
Definition: approx_sphere.cpp:68
TSElasticPostStep::preProcRhs
static boost::shared_ptr< FEMethod > preProcRhs
Definition: TSElasticPostStep.cpp:187
TSElasticPostStep::prjF
static SmartPetscObj< Vec > prjF
Definition: TSElasticPostStep.cpp:183
MAX_DOFS_ON_ENTITY
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
Definition: definitions.h:249
DomainEleOp
EshelbianPlasticity::SMALL_ROT
@ SMALL_ROT
Definition: EshelbianPlasticity.hpp:45
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::vectorDuplicate
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
Definition: PetscSmartObj.hpp:221
EshelbianPlasticity::LARGE_ROT
@ LARGE_ROT
Definition: EshelbianPlasticity.hpp:45
EshelbianCore::frontAdjEdges
boost::shared_ptr< Range > frontAdjEdges
Definition: EshelbianCore.hpp:286
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
TSElasticPostStep::prjD
static SmartPetscObj< Vec > prjD
Definition: TSElasticPostStep.cpp:182
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:64
TSElasticPostStep::postStepInitialise
static MoFEMErrorCode postStepInitialise(EshelbianCore *ep_ptr)
Definition: TSElasticPostStep.cpp:11
EshelbianCore::projectGeometry
MoFEMErrorCode projectGeometry(const EntityHandle meshset=0)
Definition: EshelbianPlasticity.cpp:1251
DomainEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
Definition: child_and_parent.cpp:34
EshelbianCore::crackingOn
static PetscBool crackingOn
Definition: EshelbianCore.hpp:22
EshelbianPlasticity::MODERATE_ROT
@ MODERATE_ROT
Definition: EshelbianPlasticity.hpp:45
EshelbianCore::solTSStep
SmartPetscObj< Vec > solTSStep
Definition: EshelbianCore.hpp:302
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
EshelbianCore::dmPrjSpatial
SmartPetscObj< DM > dmPrjSpatial
Projection spatial displacement.
Definition: EshelbianCore.hpp:98
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
TSElasticPostStep::prjDM
static SmartPetscObj< DM > prjDM
Definition: TSElasticPostStep.cpp:184
EshelbianCore
Definition: EshelbianCore.hpp:12