|
| v0.14.0
|
Interface for Time Stepping (TS) solver.
More...
#include <src/petsc/TsCtx.hpp>
|
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 jacobian 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 Jacobian for second order PDE in time. More...
|
|
Interface for Time Stepping (TS) solver.
- Examples
- level_set.cpp, nonlinear_dynamics.cpp, and Remodeling.cpp.
Definition at line 17 of file TsCtx.hpp.
◆ BasicMethodsSequence
◆ FEMethodsSequence
◆ PairNameFEMethodPtr
◆ TsCtx()
MoFEM::TsCtx::TsCtx |
( |
MoFEM::Interface & |
m_field, |
|
|
const std::string & |
problem_name |
|
) |
| |
Definition at line 5 of file TsCtx.cpp.
17 auto core_log = logging::core::get();
◆ ~TsCtx()
virtual MoFEM::TsCtx::~TsCtx |
( |
| ) |
|
|
virtualdefault |
◆ clearLoops()
Clear loops.
- Returns
- MoFEMErrorCode
Definition at line 36 of file TsCtx.cpp.
◆ getLoopsIFunction()
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 59 of file TsCtx.hpp.
◆ getLoopsIJacobian()
Get the loops to do IJacobian object.
It is sequence of finite elements used to evalite the left hand sie of implicit time solver.
- Returns
- FEMethodsSequence&
Definition at line 79 of file TsCtx.hpp.
◆ getLoopsMonitor()
◆ getLoopsRHSFunction()
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 69 of file TsCtx.hpp.
◆ getLoopsRHSJacobian()
Get the loops to do RHSJacobian object.
It is sequence of finite elements used to evalite the left hand sie of implicit time solver.
- Returns
- FEMethodsSequence&
Definition at line 89 of file TsCtx.hpp.
◆ getPostProcessIFunction()
◆ getPostProcessIJacobian()
◆ getPostProcessMonitor()
◆ getPostProcessRHSFunction()
Get the postProcess to do RHSFunction object.
- Returns
- BasicMethodsSequence&
Definition at line 178 of file TsCtx.hpp.
◆ getPostProcessRHSJacobian()
Get the postProcess to do RHSJacobian object.
- Returns
- BasicMethodsSequence&
Definition at line 160 of file TsCtx.hpp.
◆ getPreProcessIFunction()
◆ getPreProcessIJacobian()
◆ getPreProcessMonitor()
Get the preProcess to do Monitor object.
- Returns
- BasicMethodsSequence&
Definition at line 137 of file TsCtx.hpp.
◆ getPreProcessRHSFunction()
Get the preProcess to do RHSFunction object.
- Returns
- BasicMethodsSequence&
Definition at line 169 of file TsCtx.hpp.
◆ getPreProcessRHSJacobian()
Get the preProcess to do RHSJacobian object.
- Returns
- BasicMethodsSequence&
Definition at line 151 of file TsCtx.hpp.
◆ TsMonitorSet
PetscErrorCode TsMonitorSet |
( |
TS |
ts, |
|
|
PetscInt |
step, |
|
|
PetscReal |
t, |
|
|
Vec |
u, |
|
|
void * |
ctx |
|
) |
| |
|
friend |
Set monitor for TS solver.
See PETSc for details
- Parameters
-
- Returns
- PetscErrorCode
Definition at line 259 of file TsCtx.cpp.
264 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
265 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
269 auto cache_ptr = boost::make_shared<CacheTuple>();
272 auto set = [&](
auto &fe) {
277 fe.ts_F = PETSC_NULL;
283 CHKERR TSGetSNES(ts, &fe.snes);
284 CHKERR SNESGetKSP(fe.snes, &fe.ksp);
286 fe.cacheWeakPtr = cache_ptr;
290 auto unset = [&](
auto &fe) {
306 *(lit.second),
nullptr,
◆ 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
-
- Returns
- PetscErrorCode
Definition at line 612 of file TsCtx.cpp.
617 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
618 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
619 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
620 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
621 CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
622 CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
626 auto zero_ghost_vec = [](
Vec g) {
634 for (
int i = 0;
i != s; ++
i)
637 CHKERR VecGhostRestoreLocalForm(
g, &
l);
643 auto cache_ptr = boost::make_shared<CacheTuple>();
647 #if PETSC_VERSION_GE(3, 8, 0)
648 CHKERR TSGetStepNumber(ts, &step);
650 CHKERR TSGetTimeStepNumber(ts, &step);
653 auto set = [&](
auto &fe) {
667 CHKERR TSGetSNES(ts, &fe.snes);
668 CHKERR SNESGetKSP(fe.snes, &fe.ksp);
671 fe.cacheWeakPtr = cache_ptr;
675 auto unset = [&](
auto &fe) {
697 *(lit.second),
nullptr,
714 CHKERR VecGhostUpdateBegin(
F, ADD_VALUES, SCATTER_REVERSE);
715 CHKERR VecGhostUpdateEnd(
F, ADD_VALUES, SCATTER_REVERSE);
◆ 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 Jacobian for second order PDE in time.
See PETSc for details
- Parameters
-
ts | |
t | time at step/stage being solved |
u | state vectora |
u_t | time derivative of state vector |
u_tt | second time derivative of state vector |
a | shift for u_t |
aa | shift for u_tt |
A | Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt |
B | preconditioning matrix for J, may be same as J |
ctx | TsCtx context for matrix evaluation routine |
- Returns
- PetscErrorCode
Definition at line 511 of file TsCtx.cpp.
518 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
519 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
520 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
521 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
522 CHKERR VecGhostUpdateBegin(u_tt, INSERT_VALUES, SCATTER_FORWARD);
523 CHKERR VecGhostUpdateEnd(u_tt, INSERT_VALUES, SCATTER_FORWARD);
530 #if PETSC_VERSION_GE(3, 8, 0)
531 CHKERR TSGetStepNumber(ts, &step);
533 CHKERR TSGetTimeStepNumber(ts, &step);
538 auto cache_ptr = boost::make_shared<CacheTuple>();
541 auto set = [&](
auto &fe) {
559 CHKERR TSGetSNES(ts, &fe.snes);
560 CHKERR SNESGetKSP(fe.snes, &fe.ksp);
563 fe.cacheWeakPtr = cache_ptr;
567 auto unset = [&](
auto &fe) {
588 *(lit.second),
nullptr,
605 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
606 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
◆ 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
-
- Returns
- PetscErrorCode
Definition at line 56 of file TsCtx.cpp.
61 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
62 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
63 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
64 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
68 auto zero_ghost_vec = [](
Vec g) {
76 for (
int i = 0;
i != s; ++
i)
79 CHKERR VecGhostRestoreLocalForm(
g, &
l);
87 #if PETSC_VERSION_GE(3, 8, 0)
88 CHKERR TSGetStepNumber(ts, &step);
90 CHKERR TSGetTimeStepNumber(ts, &step);
93 auto cache_ptr = boost::make_shared<CacheTuple>();
96 auto set = [&](
auto &fe) {
109 CHKERR TSGetSNES(ts, &fe.snes);
110 CHKERR SNESGetKSP(fe.snes, &fe.ksp);
112 fe.cacheWeakPtr = cache_ptr;
116 auto unset = [&](
auto &fe) {
138 *(lit.second),
nullptr,
155 CHKERR VecGhostUpdateBegin(
F, ADD_VALUES, SCATTER_REVERSE);
156 CHKERR VecGhostUpdateEnd(
F, ADD_VALUES, SCATTER_REVERSE);
◆ TsSetIJacobian
PetscErrorCode TsSetIJacobian |
( |
TS |
ts, |
|
|
PetscReal |
t, |
|
|
Vec |
u, |
|
|
Vec |
U_t, |
|
|
PetscReal |
a, |
|
|
Mat |
A, |
|
|
Mat |
B, |
|
|
void * |
ctx |
|
) |
| |
|
friend |
Set function evaluating jacobian in TS solver.
See PETSc for details
- Parameters
-
- Returns
- PetscErrorCode
Definition at line 165 of file TsCtx.cpp.
171 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
172 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
173 CHKERR VecGhostUpdateBegin(u_t, INSERT_VALUES, SCATTER_FORWARD);
174 CHKERR VecGhostUpdateEnd(u_t, INSERT_VALUES, SCATTER_FORWARD);
181 #if PETSC_VERSION_GE(3, 8, 0)
182 CHKERR TSGetStepNumber(ts, &step);
184 CHKERR TSGetTimeStepNumber(ts, &step);
189 auto cache_ptr = boost::make_shared<CacheTuple>();
192 auto set = [&](
auto &fe) {
207 CHKERR TSGetSNES(ts, &fe.snes);
208 CHKERR SNESGetKSP(fe.snes, &fe.ksp);
210 fe.cacheWeakPtr = cache_ptr;
214 auto unset = [&](
auto &fe) {
235 *(lit.second),
nullptr,
252 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
253 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
◆ TsSetRHSFunction
PetscErrorCode TsSetRHSFunction |
( |
TS |
ts, |
|
|
PetscReal |
t, |
|
|
Vec |
u, |
|
|
Vec |
F, |
|
|
void * |
ctx |
|
) |
| |
|
friend |
TS solver function.
See PETSc for details
- Parameters
-
- Returns
- PetscErrorCode
Definition at line 323 of file TsCtx.cpp.
327 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
328 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
332 auto zero_ghost_vec = [](
Vec g) {
340 for (
int i = 0;
i != s; ++
i)
343 CHKERR VecGhostRestoreLocalForm(
g, &
l);
349 auto cache_ptr = boost::make_shared<CacheTuple>();
353 #if PETSC_VERSION_GE(3, 8, 0)
354 CHKERR TSGetStepNumber(ts, &step);
356 CHKERR TSGetTimeStepNumber(ts, &step);
359 auto set = [&](
auto &fe) {
370 fe.cacheWeakPtr = cache_ptr;
374 auto unset = [&](
auto &fe) {
395 *(lit.second),
nullptr,
412 CHKERR VecGhostUpdateBegin(
F, ADD_VALUES, SCATTER_REVERSE);
413 CHKERR VecGhostUpdateEnd(
F, ADD_VALUES, SCATTER_REVERSE);
◆ TsSetRHSJacobian
PetscErrorCode TsSetRHSJacobian |
( |
TS |
ts, |
|
|
PetscReal |
t, |
|
|
Vec |
u, |
|
|
Mat |
A, |
|
|
Mat |
B, |
|
|
void * |
ctx |
|
) |
| |
|
friend |
TS solver function.
See PETSc for details
- Parameters
-
- Returns
- PetscErrorCode
Definition at line 422 of file TsCtx.cpp.
427 CHKERR VecGhostUpdateBegin(u, INSERT_VALUES, SCATTER_FORWARD);
428 CHKERR VecGhostUpdateEnd(u, INSERT_VALUES, SCATTER_FORWARD);
438 auto cache_ptr = boost::make_shared<CacheTuple>();
442 #if PETSC_VERSION_GE(3, 8, 0)
443 CHKERR TSGetStepNumber(ts, &step);
445 CHKERR TSGetTimeStepNumber(ts, &step);
448 auto set = [&](
auto &fe) {
460 fe.cacheWeakPtr = cache_ptr;
464 auto unset = [&](
auto &fe) {
486 *(lit.second),
nullptr,
503 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
504 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
◆ bH
If set to MF_EXIST check if element exist.
Definition at line 23 of file TsCtx.hpp.
◆ loopsIFunction
◆ loopsIJacobian
◆ loopsMonitor
◆ loopsRHSFunction
◆ loopsRHSJacobian
◆ matAssembleSwitch
boost::movelib::unique_ptr<bool> MoFEM::TsCtx::matAssembleSwitch |
|
private |
◆ mField
◆ moab
moab::Interface& MoFEM::TsCtx::moab |
◆ MOFEM_EVENT_TsCtxI2Function
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxI2Function |
|
private |
◆ MOFEM_EVENT_TsCtxI2Jacobian
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxI2Jacobian |
|
private |
◆ MOFEM_EVENT_TsCtxIFunction
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxIFunction |
|
private |
◆ MOFEM_EVENT_TsCtxIJacobian
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxIJacobian |
|
private |
◆ MOFEM_EVENT_TsCtxMonitor
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxMonitor |
|
private |
◆ MOFEM_EVENT_TsCtxRHSFunction
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxRHSFunction |
|
private |
◆ MOFEM_EVENT_TsCtxRHSJacobian
PetscLogEvent MoFEM::TsCtx::MOFEM_EVENT_TsCtxRHSJacobian |
|
private |
◆ postProcessIFunction
◆ postProcessIJacobian
◆ postProcessMonitor
◆ postProcessRHSFunction
◆ postProcessRHSJacobian
◆ preProcessIFunction
◆ preProcessIJacobian
◆ preProcessMonitor
◆ preProcessRHSFunction
◆ preProcessRHSJacobian
◆ problemName
std::string MoFEM::TsCtx::problemName |
◆ vecAssembleSwitch
boost::movelib::unique_ptr<bool> MoFEM::TsCtx::vecAssembleSwitch |
|
private |
◆ zeroMatrix
bool MoFEM::TsCtx::zeroMatrix |
The documentation for this struct was generated from the following files:
BasicMethodsSequence preProcessRHSJacobian
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
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.
static constexpr Switches CtxSetB
static constexpr Switches CtxSetA
PetscLogEvent MOFEM_EVENT_TsCtxIJacobian
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
PetscLogEvent MOFEM_EVENT_TsCtxI2Function
PetscLogEvent MOFEM_EVENT_TsCtxIFunction
PetscLogEvent MOFEM_EVENT_TsCtxMonitor
boost::movelib::unique_ptr< bool > vecAssembleSwitch
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
FEMethodsSequence loopsIJacobian
MoFEMTypes bH
If set to MF_EXIST check if element exist.
FEMethodsSequence loopsRHSJacobian
virtual MoFEMErrorCode problem_basic_method_preProcess(const Problem *problem_ptr, BasicMethod &method, int verb=DEFAULT_VERBOSITY)=0
Set data for BasicMethod.
BasicMethodsSequence postProcessIJacobian
#define CHKERR
Inline error check.
virtual moab::Interface & get_moab()=0
MoFEM::Interface & mField
BasicMethodsSequence postProcessMonitor
FEMethodsSequence loopsMonitor
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
static constexpr Switches CtxSetX_T
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
virtual MoFEMErrorCode cache_problem_entities(const std::string prb_name, CacheTupleWeakPtr cache_ptr)=0
Cache variables.
boost::movelib::unique_ptr< bool > matAssembleSwitch
BasicMethodsSequence preProcessMonitor
BasicMethodsSequence postProcessIFunction
PetscLogEvent MOFEM_EVENT_TsCtxRHSFunction
BasicMethodsSequence preProcessIFunction
FEMethodsSequence loopsIFunction
BasicMethodsSequence preProcessIJacobian
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
static constexpr Switches CtxSetTime
static constexpr Switches CtxSetNone
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
constexpr double t
plate stiffness
FEMethodsSequence loopsRHSFunction
FTensor::Index< 'i', SPACE_DIM > i
PetscLogEvent MOFEM_EVENT_TsCtxRHSJacobian
static constexpr Switches CtxSetX_TT
BasicMethodsSequence postProcessRHSFunction
static constexpr Switches CtxSetX
PetscLogEvent MOFEM_EVENT_TsCtxI2Jacobian
const FTensor::Tensor2< T, Dim, Dim > Vec
BasicMethodsSequence postProcessRHSJacobian
static constexpr Switches CtxSetF
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
TsCtx(MoFEM::Interface &m_field, const std::string &problem_name)
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
BasicMethodsSequence preProcessRHSFunction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'l', 3 > l