#include <tutorials/adv-0/src/PlasticOpsMonitor.hpp>
|
| 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) |
|
MoFEMErrorCode | preProcess () |
|
MoFEMErrorCode | operator() () |
|
MoFEMErrorCode | postProcess () |
|
- Examples
- plastic.cpp.
Definition at line 9 of file PlasticOpsMonitor.hpp.
◆ Monitor()
PlasticOps::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 |
|
) |
| |
|
inline |
◆ operator()()
MoFEMErrorCode PlasticOps::Monitor::operator() |
( |
| ) |
|
|
inline |
◆ postProcess()
MoFEMErrorCode PlasticOps::Monitor::postProcess |
( |
| ) |
|
|
inline |
- Examples
- PlasticOpsMonitor.hpp.
Definition at line 28 of file PlasticOpsMonitor.hpp.
34 auto make_vtk = [&]() {
40 boost::lexical_cast<std::string>(ts_step) +
47 "out_skin_plastic_" + boost::lexical_cast<std::string>(ts_step) +
53 auto calculate_reaction = [&]() {
61 auto post_proc_residual = [&](
auto dm,
auto f_res,
auto out_name) {
64 boost::make_shared<PostProcBrokenMeshInMoab<DomainEle>>(
66 using OpPPMap = OpPostProcMapInMoab<SPACE_DIM, SPACE_DIM>;
67 auto u_vec = boost::make_shared<MatrixDouble>();
68 post_proc_fe->getOpPtrVector().push_back(
69 new OpCalculateVectorFieldValues<SPACE_DIM>(
"U", u_vec, f_res));
70 post_proc_fe->getOpPtrVector().push_back(
74 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
85 post_proc_fe->writeFile(
"res.h5m");
89 CHKERR post_proc_residual(
dM,
r,
"reaction");
95 auto print_max_min = [&](
auto &tuple,
const std::string msg) {
97 CHKERR VecScatterBegin(std::get<1>(tuple), ts_u, std::get<0>(tuple),
98 INSERT_VALUES, SCATTER_FORWARD);
99 CHKERR VecScatterEnd(std::get<1>(tuple), ts_u, std::get<0>(tuple),
100 INSERT_VALUES, SCATTER_FORWARD);
102 CHKERR VecMax(std::get<0>(tuple), PETSC_NULL, &max);
103 CHKERR VecMin(std::get<0>(tuple), PETSC_NULL, &min);
104 MOFEM_LOG_C(
"PLASTICITY", Sev::inform,
"%s time %3.4e min %3.4e max %3.4e",
105 msg.c_str(), ts_t, min, max);
111 if (!(ts_step % se)) {
115 CHKERR calculate_reaction();
◆ preProcess()
MoFEMErrorCode PlasticOps::Monitor::preProcess |
( |
| ) |
|
|
inline |
◆ dM
SmartPetscObj<DM> PlasticOps::Monitor::dM |
|
private |
◆ postProcFe
boost::shared_ptr<PostProcEle> PlasticOps::Monitor::postProcFe |
|
private |
◆ reactionFe
boost::shared_ptr<DomainEle> PlasticOps::Monitor::reactionFe |
|
private |
◆ skinPostProcFe
◆ uXScatter
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > PlasticOps::Monitor::uXScatter |
|
private |
◆ uYScatter
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > PlasticOps::Monitor::uYScatter |
|
private |
◆ uZScatter
std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > PlasticOps::Monitor::uZScatter |
|
private |
The documentation for this struct was generated from the following file:
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Post post-proc data at points from hash maps.