v0.13.2
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Private Attributes | Friends | List of all members
MoFEM::TsCtx Struct Reference

Interface for Time Stepping (TS) solver. More...

#include <src/petsc/TsCtx.hpp>

Collaboration diagram for MoFEM::TsCtx:
[legend]

Public Types

typedef MoFEM::PairNameFEMethodPtr PairNameFEMethodPtr
 
typedef MoFEM::FEMethodsSequence FEMethodsSequence
 
typedef MoFEM::BasicMethodsSequence BasicMethodsSequence
 

Public Member Functions

 TsCtx (MoFEM::Interface &m_field, const std::string &problem_name)
 
virtual ~TsCtx ()=default
 
FEMethodsSequencegetLoopsIFunction ()
 Get the loops to do IFunction object. More...
 
FEMethodsSequencegetLoopsRHSFunction ()
 Get the loops to do RHSFunction object. More...
 
FEMethodsSequencegetLoopsIJacobian ()
 Get the loops to do IJacobian object. More...
 
FEMethodsSequencegetLoopsRHSJacobian ()
 Get the loops to do RHSJacobian object. More...
 
FEMethodsSequencegetLoopsMonitor ()
 Get the loops to do Monitor object. More...
 
BasicMethodsSequencegetPreProcessIFunction ()
 Get the preProcess to do IFunction object. More...
 
BasicMethodsSequencegetPostProcessIFunction ()
 Get the postProcess to do IFunction object. More...
 
BasicMethodsSequencegetPreProcessIJacobian ()
 Get the preProcess to do IJacobian object. More...
 
BasicMethodsSequencegetPostProcessIJacobian ()
 Get the postProcess to do IJacobian object. More...
 
BasicMethodsSequencegetPreProcessMonitor ()
 Get the preProcess to do Monitor object. More...
 
BasicMethodsSequencegetPostProcessMonitor ()
 Get the postProcess to do Monitor object. More...
 
BasicMethodsSequencegetPreProcessRHSJacobian ()
 Get the preProcess to do RHSJacobian object. More...
 
BasicMethodsSequencegetPostProcessRHSJacobian ()
 Get the postProcess to do RHSJacobian object. More...
 
BasicMethodsSequencegetPreProcessRHSFunction ()
 Get the preProcess to do RHSFunction object. More...
 
BasicMethodsSequencegetPostProcessRHSFunction ()
 Get the postProcess to do RHSFunction object. More...
 
MoFEMErrorCode clearLoops ()
 Clear loops. More...
 

Public Attributes

MoFEM::InterfacemField
 
moab::Interface & moab
 
std::string problemName
 
MoFEMTypes bH
 If set to MF_EXIST check if element exist. More...
 
FEMethodsSequence loopsIJacobian
 
FEMethodsSequence loopsIFunction
 
FEMethodsSequence loopsMonitor
 
FEMethodsSequence loopsRHSJacobian
 
FEMethodsSequence loopsRHSFunction
 
BasicMethodsSequence preProcessIJacobian
 
BasicMethodsSequence postProcessIJacobian
 
BasicMethodsSequence preProcessIFunction
 
BasicMethodsSequence postProcessIFunction
 
BasicMethodsSequence preProcessMonitor
 
BasicMethodsSequence postProcessMonitor
 
BasicMethodsSequence preProcessRHSJacobian
 
BasicMethodsSequence preProcessRHSFunction
 
BasicMethodsSequence postProcessRHSJacobian
 
BasicMethodsSequence postProcessRHSFunction
 
bool zeroMatrix
 

Private Attributes

PetscLogEvent MOFEM_EVENT_TsCtxRHSFunction
 
PetscLogEvent MOFEM_EVENT_TsCtxRHSJacobian
 
PetscLogEvent MOFEM_EVENT_TsCtxIFunction
 
PetscLogEvent MOFEM_EVENT_TsCtxIJacobian
 
PetscLogEvent MOFEM_EVENT_TsCtxMonitor
 
PetscLogEvent MOFEM_EVENT_TsCtxI2Function
 
PetscLogEvent MOFEM_EVENT_TsCtxI2Jacobian
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 

Friends

PetscErrorCode TsSetIFunction (TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
 Set IFunction for TS solver. More...
 
PetscErrorCode TsSetIJacobian (TS ts, PetscReal t, Vec u, Vec U_t, PetscReal a, Mat A, Mat B, void *ctx)
 Set function evaluating jacobina in TS solver. More...
 
PetscErrorCode TsMonitorSet (TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
 Set monitor for TS solver. More...
 
PetscErrorCode TsSetRHSFunction (TS ts, PetscReal t, Vec u, Vec F, void *ctx)
 TS solver function. More...
 
PetscErrorCode TsSetRHSJacobian (TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
 TS solver function. More...
 
PetscErrorCode TsSetI2Function (TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx)
 Calculation the right hand side for second order PDE in time. More...
 
PetscErrorCode TsSetI2Jacobian (TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat P, void *ctx)
 Calculation Jaconian for second order PDE in time. More...
 

Detailed Description

Interface for Time Stepping (TS) solver.

Examples
Remodeling.cpp, level_set.cpp, and nonlinear_dynamics.cpp.

Definition at line 15 of file TsCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

Definition at line 25 of file TsCtx.hpp.

◆ FEMethodsSequence

Definition at line 24 of file TsCtx.hpp.

◆ PairNameFEMethodPtr

Definition at line 23 of file TsCtx.hpp.

Constructor & Destructor Documentation

◆ TsCtx()

MoFEM::TsCtx::TsCtx ( MoFEM::Interface m_field,
const std::string &  problem_name 
)
inline

Definition at line 45 of file TsCtx.hpp.

46 : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
47 bH(MF_EXIST), zeroMatrix(true) {
48 PetscLogEventRegister("LoopTsIFunction", 0, &MOFEM_EVENT_TsCtxIFunction);
49 PetscLogEventRegister("LoopTsIJacobian", 0, &MOFEM_EVENT_TsCtxIJacobian);
50 PetscLogEventRegister("LoopTsRHSFunction", 0,
52 PetscLogEventRegister("LoopTsRHSJacobian", 0,
54 PetscLogEventRegister("LoopTsMonitor", 0, &MOFEM_EVENT_TsCtxMonitor);
55 PetscLogEventRegister("LoopTsI2Function", 0, &MOFEM_EVENT_TsCtxI2Function);
56 PetscLogEventRegister("LoopTsI2Jacobian", 0, &MOFEM_EVENT_TsCtxI2Jacobian);
57 }
@ MF_EXIST
Definition: definitions.h:100
virtual moab::Interface & get_moab()=0
PetscLogEvent MOFEM_EVENT_TsCtxI2Jacobian
Definition: TsCtx.hpp:224
bool zeroMatrix
Definition: TsCtx.hpp:43
MoFEMTypes bH
If set to MF_EXIST check if element exist.
Definition: TsCtx.hpp:21
moab::Interface & moab
Definition: TsCtx.hpp:18
PetscLogEvent MOFEM_EVENT_TsCtxRHSJacobian
Definition: TsCtx.hpp:219
PetscLogEvent MOFEM_EVENT_TsCtxMonitor
Definition: TsCtx.hpp:222
PetscLogEvent MOFEM_EVENT_TsCtxRHSFunction
Definition: TsCtx.hpp:218
PetscLogEvent MOFEM_EVENT_TsCtxIJacobian
Definition: TsCtx.hpp:221
std::string problemName
Definition: TsCtx.hpp:20
MoFEM::Interface & mField
Definition: TsCtx.hpp:17
PetscLogEvent MOFEM_EVENT_TsCtxIFunction
Definition: TsCtx.hpp:220
PetscLogEvent MOFEM_EVENT_TsCtxI2Function
Definition: TsCtx.hpp:223

◆ ~TsCtx()

virtual MoFEM::TsCtx::~TsCtx ( )
virtualdefault

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode MoFEM::TsCtx::clearLoops ( )

Clear loops.

Returns
MoFEMErrorCode

Definition at line 5 of file TsCtx.cpp.

5 {
7 loopsIJacobian.clear();
8 loopsIFunction.clear();
9 loopsMonitor.clear();
10 loopsRHSJacobian.clear();
11 loopsRHSFunction.clear();
12 preProcessIJacobian.clear();
14 preProcessIFunction.clear();
16 preProcessMonitor.clear();
17 postProcessMonitor.clear();
23}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
FEMethodsSequence loopsIJacobian
Definition: TsCtx.hpp:27
FEMethodsSequence loopsRHSFunction
Definition: TsCtx.hpp:31
BasicMethodsSequence postProcessRHSFunction
Definition: TsCtx.hpp:41
BasicMethodsSequence preProcessMonitor
Definition: TsCtx.hpp:36
BasicMethodsSequence preProcessIJacobian
Definition: TsCtx.hpp:32
BasicMethodsSequence postProcessIFunction
Definition: TsCtx.hpp:35
FEMethodsSequence loopsIFunction
Definition: TsCtx.hpp:28
FEMethodsSequence loopsRHSJacobian
Definition: TsCtx.hpp:30
BasicMethodsSequence preProcessIFunction
Definition: TsCtx.hpp:34
FEMethodsSequence loopsMonitor
Definition: TsCtx.hpp:29
BasicMethodsSequence postProcessMonitor
Definition: TsCtx.hpp:37
BasicMethodsSequence preProcessRHSJacobian
Definition: TsCtx.hpp:38
BasicMethodsSequence preProcessRHSFunction
Definition: TsCtx.hpp:39
BasicMethodsSequence postProcessIJacobian
Definition: TsCtx.hpp:33
BasicMethodsSequence postProcessRHSJacobian
Definition: TsCtx.hpp:40

◆ getLoopsIFunction()

FEMethodsSequence & MoFEM::TsCtx::getLoopsIFunction ( )
inline

Get the loops to do IFunction object.

It is sequence of finite elements used to evaluate the right hand side of implicit time solver.

Returns
FEMethodsSequence&
Examples
nonlinear_dynamics.cpp.

Definition at line 69 of file TsCtx.hpp.

69{ return loopsIFunction; }

◆ getLoopsIJacobian()

FEMethodsSequence & MoFEM::TsCtx::getLoopsIJacobian ( )
inline

Get the loops to do IJacobian object.

It is sequence of finite elements used to evalite the left hand sie of implimcit time solver.

Returns
FEMethodsSequence&

Definition at line 89 of file TsCtx.hpp.

89{ return loopsIJacobian; }

◆ getLoopsMonitor()

FEMethodsSequence & MoFEM::TsCtx::getLoopsMonitor ( )
inline

Get the loops to do Monitor object.

It is sequence used to monitor solution of time solver.

Returns
FEMethodsSequence&
Examples
EshelbianPlasticity.cpp, and nonlinear_dynamics.cpp.

Definition at line 108 of file TsCtx.hpp.

108{ return loopsMonitor; }

◆ getLoopsRHSFunction()

FEMethodsSequence & MoFEM::TsCtx::getLoopsRHSFunction ( )
inline

Get the loops to do RHSFunction object.

It is sequence of finite elements used to evaluate the right hand side of implicit time solver.

Returns
FEMethodsSequence&

Definition at line 79 of file TsCtx.hpp.

79{ return loopsRHSFunction; }

◆ getLoopsRHSJacobian()

FEMethodsSequence & MoFEM::TsCtx::getLoopsRHSJacobian ( )
inline

Get the loops to do RHSJacobian object.

It is sequence of finite elements used to evalite the left hand sie of implimcit time solver.

Returns
FEMethodsSequence&

Definition at line 99 of file TsCtx.hpp.

99{ return loopsRHSJacobian; }

◆ getPostProcessIFunction()

BasicMethodsSequence & MoFEM::TsCtx::getPostProcessIFunction ( )
inline

Get the postProcess to do IFunction object.

Returns
BasicMethodsSequence&
Examples
nonlinear_dynamics.cpp.

Definition at line 122 of file TsCtx.hpp.

122 {
124 }

◆ getPostProcessIJacobian()

BasicMethodsSequence & MoFEM::TsCtx::getPostProcessIJacobian ( )
inline

Get the postProcess to do IJacobian object.

Returns
BasicMethodsSequence&
Examples
nonlinear_dynamics.cpp.

Definition at line 138 of file TsCtx.hpp.

138 {
140 }

◆ getPostProcessMonitor()

BasicMethodsSequence & MoFEM::TsCtx::getPostProcessMonitor ( )
inline

Get the postProcess to do Monitor object.

Returns
BasicMethodsSequence&
Examples
Remodeling.cpp, and UnsaturatedFlow.hpp.

Definition at line 154 of file TsCtx.hpp.

154{ return postProcessMonitor; }

◆ getPostProcessRHSFunction()

BasicMethodsSequence & MoFEM::TsCtx::getPostProcessRHSFunction ( )
inline

Get the postProcess to do RHSFunction object.

Returns
BasicMethodsSequence&

Definition at line 188 of file TsCtx.hpp.

188 {
190 }

◆ getPostProcessRHSJacobian()

BasicMethodsSequence & MoFEM::TsCtx::getPostProcessRHSJacobian ( )
inline

Get the postProcess to do RHSJacobian object.

Returns
BasicMethodsSequence&

Definition at line 170 of file TsCtx.hpp.

170 {
172 }

◆ getPreProcessIFunction()

BasicMethodsSequence & MoFEM::TsCtx::getPreProcessIFunction ( )
inline

Get the preProcess to do IFunction object.

Returns
BasicMethodsSequence&
Examples
nonlinear_dynamics.cpp.

Definition at line 115 of file TsCtx.hpp.

115{ return preProcessIFunction; }

◆ getPreProcessIJacobian()

BasicMethodsSequence & MoFEM::TsCtx::getPreProcessIJacobian ( )
inline

Get the preProcess to do IJacobian object.

Returns
BasicMethodsSequence&
Examples
nonlinear_dynamics.cpp.

Definition at line 131 of file TsCtx.hpp.

131{ return preProcessIJacobian; }

◆ getPreProcessMonitor()

BasicMethodsSequence & MoFEM::TsCtx::getPreProcessMonitor ( )
inline

Get the preProcess to do Monitor object.

Returns
BasicMethodsSequence&

Definition at line 147 of file TsCtx.hpp.

147{ return preProcessMonitor; }

◆ getPreProcessRHSFunction()

BasicMethodsSequence & MoFEM::TsCtx::getPreProcessRHSFunction ( )
inline

Get the preProcess to do RHSFunction object.

Returns
BasicMethodsSequence&

Definition at line 179 of file TsCtx.hpp.

179 {
181 }

◆ getPreProcessRHSJacobian()

BasicMethodsSequence & MoFEM::TsCtx::getPreProcessRHSJacobian ( )
inline

Get the preProcess to do RHSJacobian object.

Returns
BasicMethodsSequence&

Definition at line 161 of file TsCtx.hpp.

161 {
163 }

Friends And Related Function Documentation

◆ TsMonitorSet

PetscErrorCode TsMonitorSet ( TS  ts,
PetscInt  step,
PetscReal  t,
Vec  u,
void *  ctx 
)
friend

Set monitor for TS solver.

See PETSc for details

Parameters
ts
step
t
u
ctx
Returns
PetscErrorCode

Definition at line 219 of file TsCtx.cpp.

220 {
222 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
223 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxMonitor, 0, 0, 0, 0);
224 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
225 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
226 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
227 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
228
229 auto set = [&](auto &fe) {
230 fe.ts = ts;
231 fe.ts_u = u;
232 fe.ts_t = t;
233 fe.ts_step = step;
234 fe.ts_F = PETSC_NULL;
236 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
237 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
239 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
240 };
241
242 auto unset = [&](auto &fe) {
243 fe.ts_ctx = TSMethod::CTX_TSNONE;
244 fe.data_ctx = PetscData::CtxSetNone;
245 };
246
247 // preproces
248 for (auto &bit : ts_ctx->preProcessMonitor) {
249 set(*bit);
251 *bit);
252 unset(*bit);
253 }
254
255 auto cache_ptr = boost::make_shared<CacheTuple>();
257
258 for (auto &lit : ts_ctx->loopsMonitor) {
259 set(*lit.second);
261 *(lit.second), nullptr,
262 ts_ctx->bH, cache_ptr);
263 unset(*lit.second);
264 }
265
266 // post process
267 for (auto &bit : ts_ctx->postProcessMonitor) {
268 set(*bit);
270 *bit);
271 unset(*bit);
272 }
273
274 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxMonitor, 0, 0, 0, 0);
276}
@ COL
Definition: definitions.h:123
#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
virtual MoFEMErrorCode problem_basic_method_postProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
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.
auto bit
set bit
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1876
constexpr double t
plate stiffness
Definition: plate.cpp:59
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
static constexpr Switches CtxSetX
Definition: LoopMethods.hpp:39
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:35
static constexpr Switches CtxSetTime
Definition: LoopMethods.hpp:42
TsCtx(MoFEM::Interface &m_field, const std::string &problem_name)
Definition: TsCtx.hpp:45
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ TsSetI2Function

PetscErrorCode TsSetI2Function ( TS  ts,
PetscReal  t,
Vec  U,
Vec  U_t,
Vec  U_tt,
Vec  F,
void *  ctx 
)
friend

Calculation the right hand side for second order PDE in time.

PETSc for details

Parameters
ts
t
u
u_t
u_tt
F
ctx
Returns
PetscErrorCode

Definition at line 563 of file TsCtx.cpp.

564 {
566 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
567 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
568 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
569 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
570 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
571 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
572 CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
573 CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
574 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
575 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
576
577 auto zero_ghost_vec = [](Vec g) {
579 Vec l;
580 CHKERR VecGhostGetLocalForm(g, &l);
581 double *a;
582 CHKERR VecGetArray(l, &a);
583 int s;
584 CHKERR VecGetLocalSize(l, &s);
585 for (int i = 0; i != s; ++i)
586 a[i] = 0;
587 CHKERR VecRestoreArray(l, &a);
588 CHKERR VecGhostRestoreLocalForm(g, &l);
590 };
591 CHKERR zero_ghost_vec(F);
592
593 ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
594
595 int step;
596#if PETSC_VERSION_GE(3, 8, 0)
597 CHKERR TSGetStepNumber(ts, &step);
598#else
599 CHKERR TSGetTimeStepNumber(ts, &step);
600#endif
601
602 auto set = [&](auto &fe) {
603 fe.ts_u = u;
604 fe.ts_u_t = u_t;
605 fe.ts_u_tt = u_tt;
606 fe.ts_F = F;
607 fe.ts_t = t;
608 fe.ts_step = step;
611 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
612 fe.data_ctx = PetscData::CtxSetF | PetscData::CtxSetX |
615 fe.ts = ts;
616 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
617 };
618
619 auto unset = [&](auto &fe) {
620 fe.ts_ctx = TSMethod::CTX_TSNONE;
621 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
622 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
623 fe.data_ctx = PetscData::CtxSetNone;
624 };
625
626 // preprocess
627 for (auto &bit : ts_ctx->preProcessIFunction) {
628 bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
629 set(*bit);
631 *bit);
632 unset(*bit);
633 ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
634 }
635
636 auto cache_ptr = boost::make_shared<CacheTuple>();
638
639 // fe loops
640 for (auto &lit : ts_ctx->loopsIFunction) {
641 lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
642 set(*lit.second);
644 *(lit.second), nullptr,
645 ts_ctx->bH, cache_ptr);
646 unset(*lit.second);
647 ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
648 }
649
650 // post process
651 for (auto &bit : ts_ctx->postProcessIFunction) {
652 bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
653 set(*bit);
655 *bit);
656 unset(*bit);
657 ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
658 }
659
661 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
662 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
663 CHKERR VecAssemblyBegin(F);
664 CHKERR VecAssemblyEnd(F);
665 }
666
667 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
669}
constexpr double a
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
const FTensor::Tensor2< T, Dim, Dim > Vec
constexpr double g
static constexpr Switches CtxSetX_TT
Definition: LoopMethods.hpp:41
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:36
static constexpr Switches CtxSetX_T
Definition: LoopMethods.hpp:40
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: TsCtx.hpp:226

◆ TsSetI2Jacobian

PetscErrorCode TsSetI2Jacobian ( TS  ts,
PetscReal  t,
Vec  U,
Vec  U_t,
Vec  U_tt,
PetscReal  v,
PetscReal  a,
Mat  J,
Mat  P,
void *  ctx 
)
friend

Calculation Jaconian for second order PDE in time.

See PETSc for details

Parameters
ts
ttime at step/stage being solved
ustate vectora
u_ttime derivative of state vector
u_ttsecond time derivative of state vector
ashift for u_t
aashift for u_tt
AJacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt
Bpreconditioning matrix for J, may be same as J
ctxTsCtx context for matrix evaluation routine
Returns
PetscErrorCode

Definition at line 466 of file TsCtx.cpp.

468 {
470
471 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
472 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxI2Function, 0, 0, 0, 0);
473 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
474 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
475 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
476 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
477 CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
478 CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
479 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
480 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
481 if (ts_ctx->zeroMatrix) {
482 CHKERR MatZeroEntries(B);
483 }
484 int step;
485#if PETSC_VERSION_GE(3, 8, 0)
486 CHKERR TSGetStepNumber(ts, &step);
487#else
488 CHKERR TSGetTimeStepNumber(ts, &step);
489#endif
490
492 boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
493
494 auto set = [&](auto &fe) {
495 fe.ts_u = u;
496 fe.ts_u_t = u_t;
497 fe.ts_u_tt = u_tt;
498 fe.ts_A = A;
499 fe.ts_B = B;
500 fe.ts_t = t;
501 fe.ts_a = a;
502 fe.ts_aa = aa;
503 fe.ts_step = step;
504
507 fe.ksp_ctx = KspMethod::CTX_OPERATORS;
511 fe.ts = ts;
512 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
513 };
514
515 auto unset = [&](auto &fe) {
516 fe.ts_ctx = TSMethod::CTX_TSNONE;
517 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
518 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
519 fe.data_ctx = PetscData::CtxSetNone;
520 };
521
522 // preproces
523 for (auto &bit : ts_ctx->preProcessIJacobian) {
524 bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
525 set(*bit);
527 *bit);
528 unset(*bit);
529 ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
530 }
531
532 auto cache_ptr = boost::make_shared<CacheTuple>();
534
535 for (auto &lit : ts_ctx->loopsIJacobian) {
536 lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
537 set(*lit.second);
539 *(lit.second), nullptr,
540 ts_ctx->bH, cache_ptr);
541 unset(*lit.second);
542 ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
543 }
544
545 // post process
546 for (auto &bit : ts_ctx->postProcessIJacobian) {
547 bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
548 set(*bit);
550 *bit);
551 unset(*bit);
552 ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
553 }
554
556 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
557 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
558 }
559 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxI2Function, 0, 0, 0, 0);
561}
constexpr AssemblyType A
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: TsCtx.hpp:227

◆ TsSetIFunction

PetscErrorCode TsSetIFunction ( TS  ts,
PetscReal  t,
Vec  u,
Vec  u_t,
Vec  F,
void *  ctx 
)
friend

Set IFunction for TS solver.

See petsc for details

Parameters
ts
t
u
u_t
F
ctx
Returns
PetscErrorCode

Definition at line 25 of file TsCtx.cpp.

26 {
28 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
29 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
30 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
31 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
32 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
33 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
34 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
35 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
36
37 auto zero_ghost_vec = [](Vec g) {
39 Vec l;
40 CHKERR VecGhostGetLocalForm(g, &l);
41 double *a;
42 CHKERR VecGetArray(l, &a);
43 int s;
44 CHKERR VecGetLocalSize(l, &s);
45 for (int i = 0; i != s; ++i)
46 a[i] = 0;
47 CHKERR VecRestoreArray(l, &a);
48 CHKERR VecGhostRestoreLocalForm(g, &l);
50 };
51 CHKERR zero_ghost_vec(F);
52
53 ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
54
55 int step;
56#if PETSC_VERSION_GE(3, 8, 0)
57 CHKERR TSGetStepNumber(ts, &step);
58#else
59 CHKERR TSGetTimeStepNumber(ts, &step);
60#endif
61
62 auto set = [&](auto &fe) {
63 fe.ts = ts;
64 fe.ts_u = u;
65 fe.ts_u_t = u_t;
66 fe.ts_F = F;
67 fe.ts_t = t;
68 fe.ts_step = step;
71 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
74 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
75 };
76
77 auto unset = [&](auto &fe) {
78 fe.ts_ctx = TSMethod::CTX_TSNONE;
79 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
80 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
81 fe.data_ctx = PetscData::CtxSetNone;
82 };
83
84 // preprocess
85 for (auto &bit : ts_ctx->preProcessIFunction) {
86 bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
87 set(*bit);
89 *bit);
90 unset(*bit);
91 ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
92 }
93
94 auto cache_ptr = boost::make_shared<CacheTuple>();
96
97 // fe loops
98 for (auto &lit : ts_ctx->loopsIFunction) {
99 lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
100 set(*lit.second);
102 *(lit.second), nullptr,
103 ts_ctx->bH, cache_ptr);
104 unset(*lit.second);
105 ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
106 }
107
108 // post process
109 for (auto &bit : ts_ctx->postProcessIFunction) {
110 bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
111 set(*bit);
113 *bit);
114 unset(*bit);
115 ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
116 }
117
119 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
120 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
121 CHKERR VecAssemblyBegin(F);
122 CHKERR VecAssemblyEnd(F);
123 }
124
125 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
127}

◆ TsSetIJacobian

PetscErrorCode TsSetIJacobian ( TS  ts,
PetscReal  t,
Vec  u,
Vec  U_t,
PetscReal  a,
Mat  A,
Mat  B,
void *  ctx 
)
friend

Set function evaluating jacobina in TS solver.

See PETSc for details

Parameters
ts
t
u
u_t
a
A
B
ctx
Returns
PetscErrorCode

Definition at line 129 of file TsCtx.cpp.

130 {
132
133 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
134 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
135 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
136 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
137 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
138 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
139 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
140 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
141 if (ts_ctx->zeroMatrix) {
142 CHKERR MatZeroEntries(B);
143 }
144 int step;
145#if PETSC_VERSION_GE(3, 8, 0)
146 CHKERR TSGetStepNumber(ts, &step);
147#else
148 CHKERR TSGetTimeStepNumber(ts, &step);
149#endif
150
152 boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
153
154 auto set = [&](auto &fe) {
155 fe.ts = ts;
156 fe.ts_u = u;
157 fe.ts_u_t = u_t;
158 fe.ts_A = A;
159 fe.ts_B = B;
160 fe.ts_t = t;
161 fe.ts_a = a;
162 fe.ts_step = step;
165 fe.ksp_ctx = KspMethod::CTX_OPERATORS;
168 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
169 };
170
171 auto unset = [&](auto &fe) {
172 fe.ts_ctx = TSMethod::CTX_TSNONE;
173 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
174 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
175 fe.data_ctx = PetscData::CtxSetNone;
176 };
177
178 // preproces
179 for (auto &bit : ts_ctx->preProcessIJacobian) {
180 bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
181 set(*bit);
183 *bit);
184 unset(*bit);
185 ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
186 }
187
188 auto cache_ptr = boost::make_shared<CacheTuple>();
190
191 for (auto &lit : ts_ctx->loopsIJacobian) {
192 lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
193 set(*lit.second);
195 *(lit.second), nullptr,
196 ts_ctx->bH, cache_ptr);
197 unset(*lit.second);
198 ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
199 }
200
201 // post process
202 for (auto &bit : ts_ctx->postProcessIJacobian) {
203 bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
204 set(*bit);
206 *bit);
207 unset(*bit);
208 ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
209 }
210
212 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
213 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
214 }
215 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxIFunction, 0, 0, 0, 0);
217}

◆ TsSetRHSFunction

PetscErrorCode TsSetRHSFunction ( TS  ts,
PetscReal  t,
Vec  u,
Vec  F,
void *  ctx 
)
friend

TS solver function.

See PETSc for details

Parameters
ts
t
u
F
ctx
Returns
PetscErrorCode

Definition at line 278 of file TsCtx.cpp.

278 {
280 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
281 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxRHSFunction, 0, 0, 0, 0);
282 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
283 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
284 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
285 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
286
287 auto zero_ghost_vec = [](Vec g) {
289 Vec l;
290 CHKERR VecGhostGetLocalForm(g, &l);
291 double *a;
292 CHKERR VecGetArray(l, &a);
293 int s;
294 CHKERR VecGetLocalSize(l, &s);
295 for (int i = 0; i != s; ++i)
296 a[i] = 0;
297 CHKERR VecRestoreArray(l, &a);
298 CHKERR VecGhostRestoreLocalForm(g, &l);
300 };
301 CHKERR zero_ghost_vec(F);
302
303 ts_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
304
305 int step;
306#if PETSC_VERSION_GE(3, 8, 0)
307 CHKERR TSGetStepNumber(ts, &step);
308#else
309 CHKERR TSGetTimeStepNumber(ts, &step);
310#endif
311
312 auto set = [&](auto &fe) {
313 fe.ts_u = u;
314 fe.ts_F = F;
315 fe.ts_t = t;
316 fe.ts = ts;
317 fe.ts_step = step;
320 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
321 fe.data_ctx =
323 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
324 };
325
326 auto unset = [&](auto &fe) {
327 fe.ts_ctx = TSMethod::CTX_TSNONE;
328 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
329 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
330 fe.data_ctx = PetscData::CtxSetNone;
331 };
332
333 for (auto &bit : ts_ctx->preProcessRHSJacobian) {
334 bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
335 set(*bit);
337 *bit);
338 unset(*bit);
339 ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
340 }
341
342 auto cache_ptr = boost::make_shared<CacheTuple>();
344
345 // fe loops
346 for (auto &lit : ts_ctx->loopsRHSFunction) {
347 lit.second->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
348 set(*lit.second);
350 *(lit.second), nullptr,
351 ts_ctx->bH, cache_ptr);
352 unset(*lit.second);
353 ts_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
354 }
355
356 // post process
357 for (auto &bit : ts_ctx->postProcessRHSJacobian) {
358 bit->vecAssembleSwitch = boost::move(ts_ctx->vecAssembleSwitch);
359 set(*bit);
361 *bit);
362 unset(*bit);
363 ts_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
364 }
365
367 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
368 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
369 CHKERR VecAssemblyBegin(F);
370 CHKERR VecAssemblyEnd(F);
371 }
372
373 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxRHSFunction, 0, 0, 0, 0);
375}

◆ TsSetRHSJacobian

PetscErrorCode TsSetRHSJacobian ( TS  ts,
PetscReal  t,
Vec  u,
Mat  A,
Mat  B,
void *  ctx 
)
friend

TS solver function.

See PETSc for details

Parameters
ts
t
u
A
B
ctx
Returns
PetscErrorCode

Definition at line 377 of file TsCtx.cpp.

378 {
380 TsCtx *ts_ctx = static_cast<TsCtx *>(ctx);
381 PetscLogEventBegin(ts_ctx->MOFEM_EVENT_TsCtxRHSJacobian, 0, 0, 0, 0);
382 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
383 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
384 CHKERR ts_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
385 ts_ctx->problemName, COL, u, INSERT_VALUES, SCATTER_REVERSE);
386
387 if (ts_ctx->zeroMatrix) {
388 CHKERR MatZeroEntries(B);
389 }
390
392 boost::movelib::make_unique<bool>(ts_ctx->zeroMatrix);
393
394 int step;
395#if PETSC_VERSION_GE(3, 8, 0)
396 CHKERR TSGetStepNumber(ts, &step);
397#else
398 CHKERR TSGetTimeStepNumber(ts, &step);
399#endif
400
401 auto set = [&](auto &fe) {
402 fe.ts_u = u;
403 fe.ts_A = A;
404 fe.ts_B = B;
405 fe.ts_t = t;
406 fe.ts_step = step;
409 fe.ksp_ctx = KspMethod::CTX_OPERATORS;
412 fe.ts = ts;
413 CHK_THROW_MESSAGE(TSGetTimeStep(ts, &fe.ts_dt), "get time step failed");
414 };
415
416 auto unset = [&](auto &fe) {
417 fe.ts_ctx = TSMethod::CTX_TSNONE;
418 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
419 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
420 fe.data_ctx = PetscData::CtxSetNone;
421 };
422
423 // preprocess
424 for (auto &bit : ts_ctx->preProcessRHSJacobian) {
425 bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
426 set(*bit);
428 *bit);
429 unset(*bit);
430 ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
431 }
432
433 auto cache_ptr = boost::make_shared<CacheTuple>();
435
436 // fe loops
437 for (auto &lit : ts_ctx->loopsRHSJacobian) {
438 lit.second->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
439 set(*lit.second);
441 *(lit.second), nullptr,
442 ts_ctx->bH, cache_ptr);
443 unset(*lit.second);
444 ts_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
445 }
446
447 // post process
448 for (auto &bit : ts_ctx->postProcessRHSJacobian) {
449 bit->matAssembleSwitch = boost::move(ts_ctx->matAssembleSwitch);
450 set(*bit);
452 *bit);
453 unset(*bit);
454 ts_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
455 }
456
458 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
459 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
460 }
461
462 PetscLogEventEnd(ts_ctx->MOFEM_EVENT_TsCtxRHSJacobian, 0, 0, 0, 0);
464}

Member Data Documentation

◆ bH

MoFEMTypes MoFEM::TsCtx::bH

If set to MF_EXIST check if element exist.

Definition at line 21 of file TsCtx.hpp.

◆ loopsIFunction

FEMethodsSequence MoFEM::TsCtx::loopsIFunction

Definition at line 28 of file TsCtx.hpp.

◆ loopsIJacobian

FEMethodsSequence MoFEM::TsCtx::loopsIJacobian
Examples
EshelbianPlasticity.cpp.

Definition at line 27 of file TsCtx.hpp.

◆ loopsMonitor

FEMethodsSequence MoFEM::TsCtx::loopsMonitor

Definition at line 29 of file TsCtx.hpp.

◆ loopsRHSFunction

FEMethodsSequence MoFEM::TsCtx::loopsRHSFunction

Definition at line 31 of file TsCtx.hpp.

◆ loopsRHSJacobian

FEMethodsSequence MoFEM::TsCtx::loopsRHSJacobian

Definition at line 30 of file TsCtx.hpp.

◆ matAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::TsCtx::matAssembleSwitch
private

Definition at line 227 of file TsCtx.hpp.

◆ mField

MoFEM::Interface& MoFEM::TsCtx::mField

Definition at line 17 of file TsCtx.hpp.

◆ moab

moab::Interface& MoFEM::TsCtx::moab

Definition at line 18 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxI2Function

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxI2Function
private

Definition at line 223 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxI2Jacobian

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxI2Jacobian
private

Definition at line 224 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxIFunction

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxIFunction
private

Definition at line 220 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxIJacobian

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxIJacobian
private

Definition at line 221 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxMonitor

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxMonitor
private

Definition at line 222 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxRHSFunction

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxRHSFunction
private

Definition at line 218 of file TsCtx.hpp.

◆ MOFEM_EVENT_TsCtxRHSJacobian

PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxRHSJacobian
private

Definition at line 219 of file TsCtx.hpp.

◆ postProcessIFunction

BasicMethodsSequence MoFEM::TsCtx::postProcessIFunction

Definition at line 35 of file TsCtx.hpp.

◆ postProcessIJacobian

BasicMethodsSequence MoFEM::TsCtx::postProcessIJacobian
Examples
EshelbianPlasticity.cpp.

Definition at line 33 of file TsCtx.hpp.

◆ postProcessMonitor

BasicMethodsSequence MoFEM::TsCtx::postProcessMonitor

Definition at line 37 of file TsCtx.hpp.

◆ postProcessRHSFunction

BasicMethodsSequence MoFEM::TsCtx::postProcessRHSFunction

Definition at line 41 of file TsCtx.hpp.

◆ postProcessRHSJacobian

BasicMethodsSequence MoFEM::TsCtx::postProcessRHSJacobian

Definition at line 40 of file TsCtx.hpp.

◆ preProcessIFunction

BasicMethodsSequence MoFEM::TsCtx::preProcessIFunction

Definition at line 34 of file TsCtx.hpp.

◆ preProcessIJacobian

BasicMethodsSequence MoFEM::TsCtx::preProcessIJacobian
Examples
EshelbianPlasticity.cpp.

Definition at line 32 of file TsCtx.hpp.

◆ preProcessMonitor

BasicMethodsSequence MoFEM::TsCtx::preProcessMonitor

Definition at line 36 of file TsCtx.hpp.

◆ preProcessRHSFunction

BasicMethodsSequence MoFEM::TsCtx::preProcessRHSFunction

Definition at line 39 of file TsCtx.hpp.

◆ preProcessRHSJacobian

BasicMethodsSequence MoFEM::TsCtx::preProcessRHSJacobian

Definition at line 38 of file TsCtx.hpp.

◆ problemName

std::string MoFEM::TsCtx::problemName

Definition at line 20 of file TsCtx.hpp.

◆ vecAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::TsCtx::vecAssembleSwitch
private

Definition at line 226 of file TsCtx.hpp.

◆ zeroMatrix

bool MoFEM::TsCtx::zeroMatrix

Definition at line 43 of file TsCtx.hpp.


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