v0.14.0
PlasticOpsMonitor.hpp
/** \file PlasticOpsMonitor.hpp
* \example PlasticOpsMonitor.hpp
*/
namespace PlasticOps {
struct Monitor : public FEMethod {
Monitor(SmartPetscObj<DM> dm,
std::pair<boost::shared_ptr<PostProcEle>,
boost::shared_ptr<SkinPostProcEle>>
pair_post_proc_fe,
boost::shared_ptr<DomainEle> reaction_fe,
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> ux_scatter,
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uy_scatter,
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uz_scatter)
: dM(dm), reactionFe(reaction_fe), uXScatter(ux_scatter),
uYScatter(uy_scatter), uZScatter(uz_scatter) {
postProcFe = pair_post_proc_fe.first;
skinPostProcFe = pair_post_proc_fe.second;
};
MoFEMErrorCode preProcess() { return 0; }
MoFEMErrorCode operator()() { return 0; }
MoFEM::Interface *m_field_ptr;
auto make_vtk = [&]() {
if (postProcFe) {
getCacheWeakPtr());
CHKERR postProcFe->writeFile("out_plastic_" +
boost::lexical_cast<std::string>(ts_step) +
".h5m");
}
getCacheWeakPtr());
CHKERR skinPostProcFe->writeFile(
"out_skin_plastic_" + boost::lexical_cast<std::string>(ts_step) +
".h5m");
}
};
auto calculate_reaction = [&]() {
auto r = createDMVector(dM);
reactionFe->f = r;
CHKERR VecZeroEntries(r);
CHKERR DMoFEMLoopFiniteElements(dM, "dFE", reactionFe, getCacheWeakPtr());
#ifndef NDEBUG
auto post_proc_residual = [&](auto dm, auto f_res, auto out_name) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(
*m_field_ptr);
using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
auto u_vec = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_vec, f_res));
post_proc_fe->getOpPtrVector().push_back(
new OpPPMap(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{},
{{"RES", u_vec}},
{}, {})
);
CHKERR DMoFEMLoopFiniteElements(dM, "dFE", post_proc_fe);
post_proc_fe->writeFile("res.h5m");
};
CHKERR post_proc_residual(dM, r, "reaction");
#endif // NDEBUG
};
auto print_max_min = [&](auto &tuple, const std::string msg) {
CHKERR VecScatterBegin(std::get<1>(tuple), ts_u, std::get<0>(tuple),
INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecScatterEnd(std::get<1>(tuple), ts_u, std::get<0>(tuple),
INSERT_VALUES, SCATTER_FORWARD);
double max, min;
CHKERR VecMax(std::get<0>(tuple), PETSC_NULL, &max);
CHKERR VecMin(std::get<0>(tuple), PETSC_NULL, &min);
MOFEM_LOG_C("PLASTICITY", Sev::inform, "%s time %3.4e min %3.4e max %3.4e",
msg.c_str(), ts_t, min, max);
};
int se = 1;
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-save_every", &se, PETSC_NULL);
if (!(ts_step % se)) {
CHKERR make_vtk();
}
CHKERR calculate_reaction();
CHKERR print_max_min(uXScatter, "Ux");
CHKERR print_max_min(uYScatter, "Uy");
if constexpr (SPACE_DIM == 3)
CHKERR print_max_min(uZScatter, "Uz");
}
private:
SmartPetscObj<DM> dM;
boost::shared_ptr<PostProcEle> postProcFe;
boost::shared_ptr<SkinPostProcEle> skinPostProcFe;
boost::shared_ptr<DomainEle> reactionFe;
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uXScatter;
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uYScatter;
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter>> uZScatter;
};
}; // namespace PlasticOps
PlasticOps::Monitor::operator()
MoFEMErrorCode operator()()
Definition: PlasticOpsMonitor.hpp:26
PlasticOps::Monitor::uXScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
Definition: PlasticOpsMonitor.hpp:129
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
PlasticOps::Monitor::dM
SmartPetscObj< DM > dM
Definition: PlasticOpsMonitor.hpp:125
PlasticOps::Monitor::uYScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
Definition: PlasticOpsMonitor.hpp:130
sdf.r
int r
Definition: sdf.py:8
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
PlasticOps::Monitor::skinPostProcFe
boost::shared_ptr< SkinPostProcEle > skinPostProcFe
Definition: PlasticOpsMonitor.hpp:127
FEMethod
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::createDMVector
auto createDMVector(DM dm)
Get smart vector from DM.
Definition: DMMoFEM.hpp:1018
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
PlasticOps::Monitor::Monitor
Monitor(SmartPetscObj< DM > dm, std::pair< boost::shared_ptr< PostProcEle >, boost::shared_ptr< SkinPostProcEle >> pair_post_proc_fe, boost::shared_ptr< DomainEle > reaction_fe, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter >> ux_scatter, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter >> uy_scatter, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter >> uz_scatter)
Definition: PlasticOpsMonitor.hpp:11
PlasticOps::Monitor::uZScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
Definition: PlasticOpsMonitor.hpp:131
OpPPMap
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
Definition: photon_diffusion.cpp:29
PlasticOps
Definition: PlasticNaturalBCs.hpp:13
MoFEM::DMoFEMGetInterfacePtr
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition: DMMoFEM.cpp:418
Monitor
[Push operators to pipeline]
Definition: adolc_plasticity.cpp:333
PlasticOps::Monitor::preProcess
MoFEMErrorCode preProcess()
Definition: PlasticOpsMonitor.hpp:25
PlasticOps::Monitor::postProcFe
boost::shared_ptr< PostProcEle > postProcFe
Definition: PlasticOpsMonitor.hpp:126
PlasticOps::Monitor::reactionFe
boost::shared_ptr< DomainEle > reactionFe
Definition: PlasticOpsMonitor.hpp:128
MoFEM::DMoFEMLoopFiniteElements
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:590
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
PlasticOps::Monitor::postProcess
MoFEMErrorCode postProcess()
Definition: PlasticOpsMonitor.hpp:28
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::OpPostProcMapInMoab
Post post-proc data at points from hash maps.
Definition: PostProcBrokenMeshInMoabBase.hpp:698