v0.15.0
Loading...
Searching...
No Matches
MonitorPostProc Struct Reference
Inheritance diagram for MonitorPostProc:
[legend]
Collaboration diagram for MonitorPostProc:
[legend]

Public Member Functions

 MonitorPostProc (MoFEM::Interface &m_field, std::map< int, NonlinearElasticElement::BlockData > &set_of_blocks, NonlinearElasticElement::MyVolumeFE &fe_elastic_energy, ConvectiveMassElement::MyVolumeFE &fe_kinetic_energy)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop
 
MoFEMErrorCode operator() ()
 function is run for every finite element
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop
 
 MonitorPostProc (MoFEM::Interface &m_field)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop
 
MoFEMErrorCode operator() ()
 function is run for every finite element
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop
 
int getLoopSize () const
 get loop size
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities.
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements.
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements.
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method.
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object.
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data.
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data.
 
- Public Member Functions inherited from MoFEM::TaoMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TaoMethod ()
 
virtual ~TaoMethod ()=default
 
MoFEMErrorCode copyTao (const TaoMethod &tao)
 Copy TAO data.
 

Public Attributes

MoFEM::InterfacemField
 
PostProcVolumeOnRefinedMesh postProc
 
std::map< int, NonlinearElasticElement::BlockData > & setOfBlocks
 
NonlinearElasticElement::MyVolumeFEfeElasticEnergy
 calculate elastic energy
 
ConvectiveMassElement::MyVolumeFEfeKineticEnergy
 calculate elastic energy
 
bool iNit
 
int pRT
 
int * step
 
PostProcFaceOnRefinedMesh skinPostProc
 
PetscBool saveSkin
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element.
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element.
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method
 
int loopSize
 local number oe methods to process
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities.
 
int rAnk
 processor rank
 
int sIze
 number of processors in communicator
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities
 
const ProblemproblemPtr
 raw pointer to problem
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing.
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing.
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator.
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context.
 
KSP ksp
 KSP solver.
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec dx
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver
 
Vec & snes_x
 state vector
 
Vec & snes_dx
 solution update
 
Vec & snes_f
 residual
 
Mat & snes_A
 jacobian matrix
 
Mat & snes_B
 preconditioner of jacobian matrix
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver)
 
PetscReal ts_aa
 shift for U_tt shift for U_tt
 
PetscReal ts_t
 time
 
PetscReal ts_dt
 time step size
 
Vec & ts_u
 state vector
 
Vec & ts_u_t
 time derivative of state vector
 
Vec & ts_u_tt
 second time derivative of state vector
 
Vec & ts_F
 residual vector
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A.
 
- Public Attributes inherited from MoFEM::TaoMethod
TAOContext tao_ctx
 
Tao tao
 tao solver
 
Vec & tao_x
 
Vec & tao_f
 state vector
 
Mat & tao_A
 gradient vector
 
Mat & tao_B
 hessian matrix
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_DX = 1 << 4 , CTX_SET_X_T = 1 << 5 , CTX_SET_X_TT = 1 << 6 ,
  CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset<8>
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 
- Public Types inherited from MoFEM::TaoMethod
enum  TAOContext { CTX_TAO_OBJECTIVE , CTX_TAO_GRADIENT , CTX_TAO_HESSIAN , CTX_TAO_NONE }
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetDX = PetscData::Switches(CTX_SET_DX)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

Constructor & Destructor Documentation

◆ MonitorPostProc() [1/2]

MonitorPostProc::MonitorPostProc ( MoFEM::Interface & m_field,
std::map< int, NonlinearElasticElement::BlockData > & set_of_blocks,
NonlinearElasticElement::MyVolumeFE & fe_elastic_energy,
ConvectiveMassElement::MyVolumeFE & fe_kinetic_energy )
inline
Examples
nonlinear_dynamics.cpp.

Definition at line 45 of file nonlinear_dynamics.cpp.

50 : FEMethod(), mField(m_field), postProc(m_field),
51 setOfBlocks(set_of_blocks), feElasticEnergy(fe_elastic_energy),
52 feKineticEnergy(fe_kinetic_energy), iNit(false) {
53
54 double def_t_val = 0;
55 const EntityHandle root_meshset = mField.get_moab().get_root_set();
56
57 Tag th_step;
58 rval = m_field.get_moab().tag_get_handle(
59 "_TsStep_", 1, MB_TYPE_INTEGER, th_step,
60 MB_TAG_CREAT | MB_TAG_EXCL | MB_TAG_MESH, &def_t_val);
61 if (rval == MB_ALREADY_ALLOCATED) {
62 MOAB_THROW(m_field.get_moab().tag_get_by_ptr(th_step, &root_meshset, 1,
63 (const void **)&step));
64 } else {
65 MOAB_THROW(m_field.get_moab().tag_set_data(th_step, &root_meshset, 1,
66 &def_t_val));
67 MOAB_THROW(m_field.get_moab().tag_get_by_ptr(th_step, &root_meshset, 1,
68 (const void **)&step));
69 }
70
71 PetscBool flg = PETSC_TRUE;
72 CHK_THROW_MESSAGE(PetscOptionsGetInt(PETSC_NULLPTR, PETSC_NULLPTR,
73 "-my_output_prt", &pRT, &flg),
74 "Can not get option");
75 if (flg != PETSC_TRUE) {
76 pRT = 10;
77 }
78 }
#define MOAB_THROW(err)
Check error code of MoAB function and throw MoFEM exception.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
virtual moab::Interface & get_moab()=0
FEMethod()=default
NonlinearElasticElement::MyVolumeFE & feElasticEnergy
calculate elastic energy
MoFEM::Interface & mField
std::map< int, NonlinearElasticElement::BlockData > & setOfBlocks
ConvectiveMassElement::MyVolumeFE & feKineticEnergy
calculate elastic energy
PostProcVolumeOnRefinedMesh postProc

◆ MonitorPostProc() [2/2]

MonitorPostProc::MonitorPostProc ( MoFEM::Interface & m_field)
inline

Definition at line 60 of file thermal_unsteady.cpp.

61 : FEMethod(), mField(m_field), postProc(m_field), skinPostProc(m_field),
62 iNit(false) {
63
64 PetscBool flg = PETSC_TRUE;
65 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, PETSC_NULLPTR, "-my_output_prt", &pRT,
66 &flg);
67 CHKERRABORT(PETSC_COMM_WORLD, ierr);
68 if (flg != PETSC_TRUE) {
69 pRT = 1;
70 }
71 saveSkin = PETSC_TRUE;
72 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, PETSC_NULLPTR, "-my_save_skin",
73 &saveSkin, PETSC_NULLPTR);
74 }
#define CHKERR
Inline error check.
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
PostProcFaceOnRefinedMesh skinPostProc

Member Function Documentation

◆ operator()() [1/2]

MoFEMErrorCode MonitorPostProc::operator() ( )
inlinevirtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from MoFEM::BasicMethod.

Examples
nonlinear_dynamics.cpp.

Definition at line 128 of file nonlinear_dynamics.cpp.

128 {
131 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

◆ operator()() [2/2]

MoFEMErrorCode MonitorPostProc::operator() ( )
inlinevirtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from MoFEM::BasicMethod.

Definition at line 81 of file thermal_unsteady.cpp.

◆ postProcess() [1/2]

MoFEMErrorCode MonitorPostProc::postProcess ( )
inlinevirtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::BasicMethod.

Examples
nonlinear_dynamics.cpp.

Definition at line 133 of file nonlinear_dynamics.cpp.

◆ postProcess() [2/2]

MoFEMErrorCode MonitorPostProc::postProcess ( )
inlinevirtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::BasicMethod.

Definition at line 86 of file thermal_unsteady.cpp.

86 {
88
89 if (!iNit) {
90 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", postProc, true, false, false,
91 false);
96 CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
97
98 CHKERR addHOOpsFace3D("MESH_NODE_POSITIONS", skinPostProc, false, false);
101
102 iNit = true;
103 }
104 int step;
105 CHKERR TSGetTimeStepNumber(ts, &step);
106
107 if (pRT && (step) % pRT == 0) {
108 // CHKERR mField.loop_finite_elements("DMTHERMAL","THERMAL_FE",postProc);
109 // std::ostringstream sss;
110 // sss << "out_thermal_" << step << ".h5m";
111 // CHKERR postProc.writeFile(sss.str().c_str());
112 if (saveSkin) {
113 CHKERR mField.loop_finite_elements("DMTHERMAL", "POST_PROC_SKIN",
115 std::ostringstream sss;
116 sss << "out_skin_" << step << ".h5m";
117 CHKERR skinPostProc.writeFile(sss.str().c_str());
118 }
119 }
121 }
#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()
MoFEMErrorCode addFieldValuesPostProc(const std::string field_name, Vec v=PETSC_NULLPTR)
Add operator to post-process L2, H1, Hdiv, Hcurl field value.
MoFEMErrorCode addFieldValuesGradientPostProc(const std::string field_name, Vec v=PETSC_NULLPTR)
Add operator to post-process L2 or H1 field gradient.
MoFEMErrorCode writeFile(const std::string file_name, const char *file_type="MOAB", const char *file_options="PARALLEL=WRITE_PART")
wrote results in (MOAB) format, use "file_name.h5m"
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
TS ts
time solver
MoFEMErrorCode generateReferenceElementMesh()
MoFEMErrorCode generateReferenceElementMesh()
Generate reference mesh on single element.

◆ preProcess() [1/2]

MoFEMErrorCode MonitorPostProc::preProcess ( )
inlinevirtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from MoFEM::BasicMethod.

Examples
nonlinear_dynamics.cpp.

Definition at line 80 of file nonlinear_dynamics.cpp.

80 {
82
83 if (!iNit) {
85 if(mField.check_field("MESH_NODE_POSITIONS"))
86 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", postProc, true, false, false,
87 false);
90 CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
92
93 std::map<int, NonlinearElasticElement::BlockData>::iterator sit =
94 setOfBlocks.begin();
95 for (; sit != setOfBlocks.end(); sit++) {
98 sit->second, postProc.commonData, true));
99 }
100
101 iNit = true;
102 }
103
104 if ((*step) % pRT == 0) {
105 CHKERR mField.loop_finite_elements("DYNAMICS", "MASS_ELEMENT", postProc);
106 std::ostringstream sss;
107 sss << "out_values_" << (*step) << ".h5m";
108 CHKERR postProc.writeFile(sss.str().c_str());
109 }
110
113 CHKERR mField.loop_finite_elements("DYNAMICS", "ELASTIC", feElasticEnergy);
114
116 CHKERR mField.loop_finite_elements("DYNAMICS", "MASS_ELEMENT",
118 double E = feElasticEnergy.eNergy;
119 double T = feKineticEnergy.eNergy;
121 "DYNAMIC", Sev::inform,
122 "%d Time %3.2e Elastic energy %3.2e Kinetic Energy %3.2e Total %3.2e\n",
123 ts_step, ts_t, E, T, E + T);
124
126 }
#define MOFEM_LOG_C(channel, severity, format,...)
virtual bool check_field(const std::string &name) const =0
check if field is in database
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
SNESContext snes_ctx
PetscReal ts_t
time
PetscInt ts_step
time step number
std::vector< EntityHandle > mapGaussPts

◆ preProcess() [2/2]

MoFEMErrorCode MonitorPostProc::preProcess ( )
inlinevirtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from MoFEM::BasicMethod.

Definition at line 76 of file thermal_unsteady.cpp.

Member Data Documentation

◆ feElasticEnergy

NonlinearElasticElement::MyVolumeFE& MonitorPostProc::feElasticEnergy

calculate elastic energy

Examples
nonlinear_dynamics.cpp.

Definition at line 36 of file nonlinear_dynamics.cpp.

◆ feKineticEnergy

ConvectiveMassElement::MyVolumeFE& MonitorPostProc::feKineticEnergy

calculate elastic energy

Examples
nonlinear_dynamics.cpp.

Definition at line 38 of file nonlinear_dynamics.cpp.

◆ iNit

bool MonitorPostProc::iNit
Examples
nonlinear_dynamics.cpp.

Definition at line 40 of file nonlinear_dynamics.cpp.

◆ mField

MoFEM::Interface & MonitorPostProc::mField
Examples
nonlinear_dynamics.cpp.

Definition at line 32 of file nonlinear_dynamics.cpp.

◆ postProc

PostProcVolumeOnRefinedMesh MonitorPostProc::postProc
Examples
nonlinear_dynamics.cpp.

Definition at line 33 of file nonlinear_dynamics.cpp.

◆ pRT

int MonitorPostProc::pRT
Examples
nonlinear_dynamics.cpp.

Definition at line 42 of file nonlinear_dynamics.cpp.

◆ saveSkin

PetscBool MonitorPostProc::saveSkin

Definition at line 58 of file thermal_unsteady.cpp.

◆ setOfBlocks

std::map<int, NonlinearElasticElement::BlockData>& MonitorPostProc::setOfBlocks
Examples
nonlinear_dynamics.cpp.

Definition at line 34 of file nonlinear_dynamics.cpp.

◆ skinPostProc

PostProcFaceOnRefinedMesh MonitorPostProc::skinPostProc

Definition at line 54 of file thermal_unsteady.cpp.

◆ step

int* MonitorPostProc::step
Examples
nonlinear_dynamics.cpp.

Definition at line 43 of file nonlinear_dynamics.cpp.


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