v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Public Attributes | Static Private Member Functions | Private Attributes | List of all members
TSPrePostProc Struct Reference

Set of functions called by PETSc solver used to refine and update mesh. More...

Collaboration diagram for TSPrePostProc:
[legend]

Public Member Functions

 TSPrePostProc ()=default
 
virtual ~TSPrePostProc ()=default
 
MoFEMErrorCode tsSetUp (TS ts)
 Used to setup TS solver.
 
 TSPrePostProc ()=default
 
virtual ~TSPrePostProc ()=default
 
MoFEMErrorCode tsSetUp (TS ts)
 Used to setup TS solver.
 
SmartPetscObj< VecScatter > getScatter (Vec x, Vec y, enum FR fr)
 Get scatter context for vector operations.
 
SmartPetscObj< Vec > getSubVector ()
 Create sub-problem vector.
 
 TSPrePostProc ()=default
 
virtual ~TSPrePostProc ()=default
 
MoFEMErrorCode tsSetUp (TS ts)
 Used to setup TS solver.
 

Static Public Member Functions

static MoFEMErrorCode tsPostStage (TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
 [Boundary condition]
 
static MoFEMErrorCode tsPostStep (TS ts)
 
static MoFEMErrorCode tsPreStep (TS ts)
 

Public Attributes

ExamplefsRawPtr
 
SmartPetscObj< DM > solverSubDM
 
SmartPetscObj< Vec > globSol
 
FreeSurfacefsRawPtr
 
ExampleExRawPtr
 

Static Private Member Functions

static MoFEMErrorCode tsPreProc (TS ts)
 Pre process time step.
 
static MoFEMErrorCode tsPostProc (TS ts)
 Post process time step.
 
static MoFEMErrorCode tsPreStage (TS ts)
 Pre-stage processing for time stepping.
 
static MoFEMErrorCode tsSetIFunction (TS ts, PetscReal t, Vec u, Vec u_t, Vec f, void *ctx)
 Set implicit function for time stepping.
 
static MoFEMErrorCode tsSetIJacobian (TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
 Set implicit Jacobian for time stepping.
 
static MoFEMErrorCode tsMonitor (TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
 Monitor solution during time stepping.
 
static MoFEMErrorCode pcSetup (PC pc)
 Setup preconditioner.
 
static MoFEMErrorCode pcApply (PC pc, Vec pc_f, Vec pc_x)
 Apply preconditioner.
 
static MoFEMErrorCode tsPreProc (TS ts)
 Pre process time step.
 
static MoFEMErrorCode tsPostProc (TS ts)
 Post process time step.
 
static MoFEMErrorCode tsPreStage (TS ts)
 

Private Attributes

SmartPetscObj< Vec > globRes
 
SmartPetscObj< Mat > subB
 
SmartPetscObj< KSP > subKSP
 
boost::shared_ptr< SnesCtxsnesCtxPtr
 
boost::shared_ptr< TsCtxtsCtxPtr
 

Detailed Description

Set of functions called by PETSc solver used to refine and update mesh.

Note
Currently theta method is only handled by this code.
Examples
mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp, mofem/tutorials/vec-5_free_surface/free_surface.cpp, and thermoplastic.cpp.

Definition at line 383 of file dynamic_first_order_con_law.cpp.

Constructor & Destructor Documentation

◆ TSPrePostProc() [1/3]

TSPrePostProc::TSPrePostProc ( )
default

◆ ~TSPrePostProc() [1/3]

virtual TSPrePostProc::~TSPrePostProc ( )
virtualdefault

◆ TSPrePostProc() [2/3]

TSPrePostProc::TSPrePostProc ( )
default

◆ ~TSPrePostProc() [2/3]

virtual TSPrePostProc::~TSPrePostProc ( )
virtualdefault

◆ TSPrePostProc() [3/3]

TSPrePostProc::TSPrePostProc ( )
default

◆ ~TSPrePostProc() [3/3]

virtual TSPrePostProc::~TSPrePostProc ( )
virtualdefault

Member Function Documentation

◆ getScatter()

SmartPetscObj< VecScatter > TSPrePostProc::getScatter ( Vec  x,
Vec  y,
enum FR  fr 
)

Get scatter context for vector operations.

Parameters
xLocal sub-vector
yGlobal vector
frDirection flag (F=forward, R=reverse)
Returns
SmartPetscObj<VecScatter> Scatter context for vector operations

Creates scatter context to transfer data between global and sub-problem vectors

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3578 of file free_surface.cpp.

3578 {
3579 if (auto ptr = tsPrePostProc.lock()) {
3580 auto prb_ptr = ptr->fsRawPtr->mField.get_problem("SUB_SOLVER");
3581 if (auto sub_data = prb_ptr->getSubData()) {
3582 auto is = sub_data->getSmartColIs();
3583 VecScatter s;
3584 if (fr == R) {
3585 CHK_THROW_MESSAGE(VecScatterCreate(x, PETSC_NULLPTR, y, is, &s),
3586 "crate scatter");
3587 } else {
3588 CHK_THROW_MESSAGE(VecScatterCreate(x, is, y, PETSC_NULLPTR, &s),
3589 "crate scatter");
3590 }
3591 return SmartPetscObj<VecScatter>(s);
3592 }
3593 }
3596}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ R
intrusive_ptr for managing petsc objects
static boost::weak_ptr< TSPrePostProc > tsPrePostProc

◆ getSubVector()

SmartPetscObj< Vec > TSPrePostProc::getSubVector ( )

Create sub-problem vector.

Returns
SmartPetscObj<Vec> Vector compatible with sub-problem DM

Creates a vector with proper size and ghost structure for sub-problem

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3598 of file free_surface.cpp.

3598 {
3600}
auto createDMVector(DM dm, RowColData rc=RowColData::COL)
Get smart vector from DM.
Definition DMMoFEM.hpp:1237
SmartPetscObj< DM > solverSubDM

◆ pcApply()

MoFEMErrorCode TSPrePostProc::pcApply ( PC  pc,
Vec  pc_f,
Vec  pc_x 
)
staticprivate

Apply preconditioner.

Parameters
pcPETSc preconditioner object
pc_fInput vector (right-hand side)
pc_xOutput vector (preconditioned solution)
Returns
MoFEMErrorCode Success/failure code

Applies preconditioner by solving sub-problem with KSP

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3558 of file free_surface.cpp.

3558 {
3560 if (auto ptr = tsPrePostProc.lock()) {
3561 auto sub_x = ptr->getSubVector();
3562 auto sub_f = vectorDuplicate(sub_x);
3563 auto scatter = ptr->getScatter(sub_x, pc_x, R);
3564 CHKERR VecScatterBegin(scatter, pc_f, sub_f, INSERT_VALUES,
3565 SCATTER_REVERSE);
3566 CHKERR VecScatterEnd(scatter, pc_f, sub_f, INSERT_VALUES, SCATTER_REVERSE);
3567 CHKERR KSPSetOperators(ptr->subKSP, ptr->subB, ptr->subB);
3568 MOFEM_LOG("FS", Sev::verbose) << "PCShell solve";
3569 CHKERR KSPSolve(ptr->subKSP, sub_f, sub_x);
3570 MOFEM_LOG("FS", Sev::verbose) << "PCShell solve <- done";
3571 CHKERR VecScatterBegin(scatter, sub_x, pc_x, INSERT_VALUES,
3572 SCATTER_FORWARD);
3573 CHKERR VecScatterEnd(scatter, sub_x, pc_x, INSERT_VALUES, SCATTER_FORWARD);
3574 }
3576};
#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()
#define CHKERR
Inline error check.
#define MOFEM_LOG(channel, severity)
Log.
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.

◆ pcSetup()

MoFEMErrorCode TSPrePostProc::pcSetup ( PC  pc)
staticprivate

Setup preconditioner.

Parameters
pcPETSc preconditioner object
Returns
MoFEMErrorCode Success/failure code

Initializes KSP solver for shell preconditioner

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3548 of file free_surface.cpp.

3548 {
3550 if (auto ptr = tsPrePostProc.lock()) {
3551 MOFEM_LOG("FS", Sev::verbose) << "SetUP sub PC";
3552 ptr->subKSP = createKSP(ptr->fsRawPtr->mField.get_comm());
3553 CHKERR KSPSetFromOptions(ptr->subKSP);
3554 }
3556};
auto createKSP(MPI_Comm comm)

◆ tsMonitor()

MoFEMErrorCode TSPrePostProc::tsMonitor ( TS  ts,
PetscInt  step,
PetscReal  t,
Vec  u,
void *  ctx 
)
staticprivate

Monitor solution during time stepping.

Parameters
tsPETSc time stepping object
stepCurrent time step number
tCurrent time value
uCurrent solution vector
ctxUser context (pointer to FreeSurface)
Returns
MoFEMErrorCode Success/failure code

Called after each time step to monitor solution and save output Wrapper for TS monitor

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3523 of file free_surface.cpp.

3524 {
3526 if (auto ptr = tsPrePostProc.lock()) {
3527 auto get_norm = [&](auto x) {
3528 double nrm;
3529 CHKERR VecNorm(x, NORM_2, &nrm);
3530 return nrm;
3531 };
3532
3533 auto sub_u = ptr->getSubVector();
3534 auto scatter = ptr->getScatter(sub_u, u, R);
3535 CHKERR VecScatterBegin(scatter, u, sub_u, INSERT_VALUES, SCATTER_REVERSE);
3536 CHKERR VecScatterEnd(scatter, u, sub_u, INSERT_VALUES, SCATTER_REVERSE);
3537 CHKERR VecGhostUpdateBegin(sub_u, INSERT_VALUES, SCATTER_FORWARD);
3538 CHKERR VecGhostUpdateEnd(sub_u, INSERT_VALUES, SCATTER_FORWARD);
3539
3540 MOFEM_LOG("FS", Sev::verbose)
3541 << "u norm " << get_norm(u) << " u sub nom " << get_norm(sub_u);
3542
3543 CHKERR TsMonitorSet(ts, step, t, sub_u, ptr->tsCtxPtr.get());
3544 }
3546}
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition TsCtx.cpp:263
constexpr double t
plate stiffness
Definition plate.cpp:58

◆ tsPostProc() [1/2]

MoFEMErrorCode TSPrePostProc::tsPostProc ( TS  ts)
staticprivate

Post process time step.

Currently this function does not make anything major

Parameters
tsPETSc time stepping object
Returns
MoFEMErrorCode Success/failure code

Called after each time step completion for cleanup operations

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp, and thermoplastic.cpp.

Definition at line 3458 of file free_surface.cpp.

3458 {
3459 if (auto ptr = tsPrePostProc.lock()) {
3460 auto &m_field = ptr->fsRawPtr->mField;
3461 MOFEM_LOG_CHANNEL("SYNC");
3462 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "FS") << "PostProc done";
3463 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::verbose);
3464 }
3465 return 0;
3466}
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.

◆ tsPostProc() [2/2]

static MoFEMErrorCode TSPrePostProc::tsPostProc ( TS  ts)
staticprivate

Post process time step.

Currently that function do not make anything major

Parameters
ts
Returns
MoFEMErrorCode

◆ tsPostStage()

MoFEMErrorCode TSPrePostProc::tsPostStage ( TS  ts,
PetscReal  stagetime,
PetscInt  stageindex,
Vec *  Y 
)
static

[Boundary condition]

Examples
mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp.

Definition at line 580 of file dynamic_first_order_con_law.cpp.

581 {
583 // cerr << "tsPostStage " <<"\n";
584 if (auto ptr = tsPrePostProc.lock()) {
585 auto &m_field = ptr->fsRawPtr->mField;
586
587 auto fb = m_field.getInterface<FieldBlas>();
588 double dt;
589 CHKERR TSGetTimeStep(ts, &dt);
590 double time;
591 CHKERR TSGetTime(ts, &time);
592 PetscInt num_stages;
593 Vec *stage_solutions;
594
595 CHKERR TSGetStages(ts, &num_stages, &stage_solutions);
596 PetscPrintf(PETSC_COMM_WORLD, "Check timestep %d time %e dt %e\n",
597 num_stages, time, dt);
598
599 const double inv_num_step = (double)num_stages;
600 CHKERR fb->fieldCopy(1., "x_1", "x_2");
601 CHKERR fb->fieldAxpy(dt, "V", "x_2");
602 CHKERR fb->fieldCopy(1., "x_2", "x_1");
603
604 CHKERR fb->fieldCopy(-inv_num_step / dt, "F_0", "F_dot");
605 CHKERR fb->fieldAxpy(inv_num_step / dt, "F", "F_dot");
606 CHKERR fb->fieldCopy(1., "F", "F_0");
607 }
609}
double dt
const FTensor::Tensor2< T, Dim, Dim > Vec
Basic algebra on fields.
Definition FieldBlas.hpp:21

◆ tsPostStep()

MoFEMErrorCode TSPrePostProc::tsPostStep ( TS  ts)
static
Examples
mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp.

Definition at line 611 of file dynamic_first_order_con_law.cpp.

611 {
613
614 if (auto ptr = tsPrePostProc.lock()) {
615 double dt;
616 CHKERR TSGetTimeStep(ts, &dt);
617 double time;
618 CHKERR TSGetTime(ts, &time);
619 }
621}

◆ tsPreProc() [1/2]

MoFEMErrorCode TSPrePostProc::tsPreProc ( TS  ts)
staticprivate

Pre process time step.

Refine mesh and update fields before each time step

Parameters
tsPETSc time stepping object
Returns
MoFEMErrorCode Success/failure code

This function is called before each time step to:

  • Apply phase field cutoff constraints
  • Refine mesh based on interface location
  • Project solution data to new mesh
  • Update solver operators

cut-off values at nodes, i.e. abs("H") <= 1

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp, and thermoplastic.cpp.

Definition at line 3343 of file free_surface.cpp.

3343 {
3345
3346 if (auto ptr = tsPrePostProc.lock()) {
3347
3348 /**
3349 * @brief cut-off values at nodes, i.e. abs("H") <= 1
3350 *
3351 */
3352 auto cut_off_dofs = [&]() {
3354
3355 auto &m_field = ptr->fsRawPtr->mField;
3356
3357 Range current_verts;
3358 auto bit_mng = m_field.getInterface<BitRefManager>();
3360 bit(get_current_bit()), BitRefLevel().set(), MBVERTEX, current_verts);
3361
3362 auto cut_off_verts = [&](boost::shared_ptr<FieldEntity> ent_ptr) {
3364 for (auto &h : ent_ptr->getEntFieldData()) {
3365 h = cut_off(h);
3366 }
3368 };
3369
3370 auto field_blas = m_field.getInterface<FieldBlas>();
3371 CHKERR field_blas->fieldLambdaOnEntities(cut_off_verts, "H",
3372 &current_verts);
3374 };
3375
3376 CHKERR cut_off_dofs();
3377 }
3378
3379 if (auto ptr = tsPrePostProc.lock()) {
3380 MOFEM_LOG("FS", Sev::inform) << "Run step pre proc";
3381
3382 auto &m_field = ptr->fsRawPtr->mField;
3383 auto simple = m_field.getInterface<Simple>();
3384
3385 // get vector norm
3386 auto get_norm = [&](auto x) {
3387 double nrm;
3388 CHKERR VecNorm(x, NORM_2, &nrm);
3389 return nrm;
3390 };
3391
3392 // refine problem and project data, including theta data
3393 auto refine_problem = [&]() {
3395 MOFEM_LOG("FS", Sev::inform) << "Refine problem";
3396 CHKERR ptr->fsRawPtr->refineMesh(refine_overlap);
3397 CHKERR ptr->fsRawPtr->projectData();
3399 };
3400
3401 // set new jacobin operator, since problem and thus tangent matrix size has
3402 // changed
3403 auto set_jacobian_operators = [&]() {
3405 ptr->subB = createDMMatrix(ptr->solverSubDM);
3406 CHKERR KSPReset(ptr->subKSP);
3408 };
3409
3410 // set new solution
3411 auto set_solution = [&]() {
3413 MOFEM_LOG("FS", Sev::inform) << "Set solution";
3414
3415 PetscObjectState state;
3416
3417 // Record the state, and set it again. This is to fool PETSc that solution
3418 // vector is not updated. Otherwise PETSc will treat every step as a first
3419 // step.
3420
3421 // globSol is updated as result mesh refinement - this is not really set
3422 // a new solution.
3423
3424 CHKERR PetscObjectStateGet(getPetscObject(ptr->globSol.get()), &state);
3425 CHKERR DMoFEMMeshToLocalVector(simple->getDM(), ptr->globSol,
3426 INSERT_VALUES, SCATTER_FORWARD);
3427 CHKERR PetscObjectStateSet(getPetscObject(ptr->globSol.get()), state);
3428 MOFEM_LOG("FS", Sev::verbose)
3429 << "Set solution, vector norm " << get_norm(ptr->globSol);
3431 };
3432
3433 PetscBool is_theta;
3434 PetscObjectTypeCompare((PetscObject)ts, TSTHETA, &is_theta);
3435 if (is_theta) {
3436
3437 CHKERR refine_problem(); // refine problem
3438 CHKERR set_jacobian_operators(); // set new jacobian
3439 CHKERR set_solution(); // set solution
3440
3441 } else {
3442 SETERRQ(PETSC_COMM_WORLD, MOFEM_NOT_IMPLEMENTED,
3443 "Sorry, only TSTheta handling is implemented");
3444 }
3445
3446 // Need barriers, some functions in TS solver need are called collectively
3447 // and require the same state of variables
3448 PetscBarrier((PetscObject)ts);
3449
3450 MOFEM_LOG_CHANNEL("SYNC");
3451 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "FS") << "PreProc done";
3452 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), Sev::verbose);
3453 }
3454
3456}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
int refine_overlap
auto cut_off
Phase field cutoff function.
auto get_current_bit
dofs bit used to do calculations
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode, RowColData rc=RowColData::COL)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1194
MoFEMErrorCode getEntitiesByTypeAndRefLevel(const BitRefLevel bit, const BitRefLevel mask, const EntityType type, const EntityHandle meshset, int verb=0) const
add all ents from ref level given by bit to meshset
auto bit
set bit
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
PetscObject getPetscObject(T obj)
double h
Managing BitRefLevels.
MoFEMErrorCode fieldLambdaOnEntities(OneFieldFunctionOnEntities lambda, const std::string field_name, Range *ents_ptr=nullptr)
field lambda
Definition FieldBlas.cpp:50
Simple interface for fast problem set-up.
Definition Simple.hpp:27

◆ tsPreProc() [2/2]

static MoFEMErrorCode TSPrePostProc::tsPreProc ( TS  ts)
staticprivate

Pre process time step.

Refine mesh and update fields

Parameters
ts
Returns
MoFEMErrorCode

◆ tsPreStage() [1/2]

static MoFEMErrorCode TSPrePostProc::tsPreStage ( TS  ts)
staticprivate

Pre-stage processing for time stepping.

Parameters
tsPETSc time stepping object
Returns
MoFEMErrorCode Success/failure code
Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp, and thermoplastic.cpp.

◆ tsPreStage() [2/2]

static MoFEMErrorCode TSPrePostProc::tsPreStage ( TS  ts)
staticprivate

◆ tsPreStep()

MoFEMErrorCode TSPrePostProc::tsPreStep ( TS  ts)
static
Examples
mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp.

Definition at line 623 of file dynamic_first_order_con_law.cpp.

623 {
625
626 if (auto ptr = tsPrePostProc.lock()) {
627 double dt;
628 CHKERR TSGetTimeStep(ts, &dt);
629 double time;
630 CHKERR TSGetTime(ts, &time);
631 int step_num;
632 CHKERR TSGetStepNumber(ts, &step_num);
633 }
635}

◆ tsSetIFunction()

MoFEMErrorCode TSPrePostProc::tsSetIFunction ( TS  ts,
PetscReal  t,
Vec  u,
Vec  u_t,
Vec  f,
void *  ctx 
)
staticprivate

Set implicit function for time stepping.

Parameters
tsPETSc time stepping object
tCurrent time
uSolution vector at current time
u_tTime derivative of solution
fOutput residual vector F(t,u,u_t)
ctxUser context (unused)
Returns
MoFEMErrorCode Success/failure code

Wrapper that scatters global vectors to sub-problem and evaluates residual

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3468 of file free_surface.cpp.

3469 {
3471 if (auto ptr = tsPrePostProc.lock()) {
3472 auto sub_u = ptr->getSubVector();
3473 auto sub_u_t = vectorDuplicate(sub_u);
3474 auto sub_f = vectorDuplicate(sub_u);
3475 auto scatter = ptr->getScatter(sub_u, u, R);
3476
3477 auto apply_scatter_and_update = [&](auto x, auto sub_x) {
3479 CHKERR VecScatterBegin(scatter, x, sub_x, INSERT_VALUES, SCATTER_REVERSE);
3480 CHKERR VecScatterEnd(scatter, x, sub_x, INSERT_VALUES, SCATTER_REVERSE);
3481 CHKERR VecGhostUpdateBegin(sub_x, INSERT_VALUES, SCATTER_FORWARD);
3482 CHKERR VecGhostUpdateEnd(sub_x, INSERT_VALUES, SCATTER_FORWARD);
3484 };
3485
3486 CHKERR apply_scatter_and_update(u, sub_u);
3487 CHKERR apply_scatter_and_update(u_t, sub_u_t);
3488
3489 CHKERR TsSetIFunction(ts, t, sub_u, sub_u_t, sub_f, ptr->tsCtxPtr.get());
3490 CHKERR VecScatterBegin(scatter, sub_f, f, INSERT_VALUES, SCATTER_FORWARD);
3491 CHKERR VecScatterEnd(scatter, sub_f, f, INSERT_VALUES, SCATTER_FORWARD);
3492 }
3494}
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition TsCtx.cpp:56

◆ tsSetIJacobian()

MoFEMErrorCode TSPrePostProc::tsSetIJacobian ( TS  ts,
PetscReal  t,
Vec  u,
Vec  u_t,
PetscReal  a,
Mat  A,
Mat  B,
void *  ctx 
)
staticprivate

Set implicit Jacobian for time stepping.

Parameters
tsPETSc time stepping object
tCurrent time
uSolution vector at current time
u_tTime derivative of solution
aShift parameter for implicit methods
AJacobian matrix (input)
BJacobian matrix (output, often same as A)
ctxUser context (unused)
Returns
MoFEMErrorCode Success/failure code

Wrapper that assembles Jacobian matrix for sub-problem Wrapper for SNES Lhs

Examples
mofem/tutorials/vec-5_free_surface/free_surface.cpp.

Definition at line 3496 of file free_surface.cpp.

3498 {
3500 if (auto ptr = tsPrePostProc.lock()) {
3501 auto sub_u = ptr->getSubVector();
3502 auto sub_u_t = vectorDuplicate(sub_u);
3503 auto scatter = ptr->getScatter(sub_u, u, R);
3504
3505 auto apply_scatter_and_update = [&](auto x, auto sub_x) {
3507 CHKERR VecScatterBegin(scatter, x, sub_x, INSERT_VALUES, SCATTER_REVERSE);
3508 CHKERR VecScatterEnd(scatter, x, sub_x, INSERT_VALUES, SCATTER_REVERSE);
3509 CHKERR VecGhostUpdateBegin(sub_x, INSERT_VALUES, SCATTER_FORWARD);
3510 CHKERR VecGhostUpdateEnd(sub_x, INSERT_VALUES, SCATTER_FORWARD);
3512 };
3513
3514 CHKERR apply_scatter_and_update(u, sub_u);
3515 CHKERR apply_scatter_and_update(u_t, sub_u_t);
3516
3517 CHKERR TsSetIJacobian(ts, t, sub_u, sub_u_t, a, ptr->subB, ptr->subB,
3518 ptr->tsCtxPtr.get());
3519 }
3521}
constexpr double a
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.
Definition TsCtx.cpp:169

◆ tsSetUp() [1/3]

MoFEMErrorCode TSPrePostProc::tsSetUp ( TS  ts)

Used to setup TS solver.

Parameters
ts
Returns
MoFEMErrorCode
Examples
mofem/tutorials/adv-4_dynamic_first_order_con_law/dynamic_first_order_con_law.cpp, mofem/tutorials/vec-5_free_surface/free_surface.cpp, and thermoplastic.cpp.

Definition at line 3602 of file free_surface.cpp.

3602 {
3603
3604 auto &m_field = fsRawPtr->mField;
3605 auto simple = m_field.getInterface<Simple>();
3606
3608
3609 auto dm = simple->getDM();
3610
3612 CHKERR TSSetIFunction(ts, globRes, tsSetIFunction, nullptr);
3613 CHKERR TSSetIJacobian(ts, PETSC_NULLPTR, PETSC_NULLPTR, tsSetIJacobian, nullptr);
3614 CHKERR TSMonitorSet(ts, tsMonitor, fsRawPtr, PETSC_NULLPTR);
3615
3616 SNES snes;
3617 CHKERR TSGetSNES(ts, &snes);
3618 auto snes_ctx_ptr = getDMSnesCtx(dm);
3619
3620 auto set_section_monitor = [&](auto snes) {
3622 CHKERR SNESMonitorSet(snes,
3623 (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
3624 void *))MoFEMSNESMonitorFields,
3625 (void *)(snes_ctx_ptr.get()), nullptr);
3627 };
3628
3629 CHKERR set_section_monitor(snes);
3630
3631 auto ksp = createKSP(m_field.get_comm());
3632 CHKERR KSPSetType(ksp, KSPPREONLY); // Run KSP internally in ShellPC
3633 auto sub_pc = createPC(fsRawPtr->mField.get_comm());
3634 CHKERR PCSetType(sub_pc, PCSHELL);
3635 CHKERR PCShellSetSetUp(sub_pc, pcSetup);
3636 CHKERR PCShellSetApply(sub_pc, pcApply);
3637 CHKERR KSPSetPC(ksp, sub_pc);
3638 CHKERR SNESSetKSP(snes, ksp);
3639
3642
3643 CHKERR TSSetPreStep(ts, tsPreProc);
3644 CHKERR TSSetPostStep(ts, tsPostProc);
3645
3647}
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition DMMoFEM.hpp:1279
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:600
auto createPC(MPI_Comm comm)
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1265
MoFEM::Interface & mField
Reference to MoFEM interface.
Definition plastic.cpp:226
virtual MPI_Comm & get_comm() const =0
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:799
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
boost::shared_ptr< TsCtx > tsCtxPtr
SmartPetscObj< Vec > globRes
SmartPetscObj< Vec > globSol
static MoFEMErrorCode tsMonitor(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Monitor solution during time stepping.
static MoFEMErrorCode pcApply(PC pc, Vec pc_f, Vec pc_x)
Apply preconditioner.
static MoFEMErrorCode tsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set implicit Jacobian for time stepping.
boost::shared_ptr< SnesCtx > snesCtxPtr
static MoFEMErrorCode tsPreProc(TS ts)
Pre process time step.
static MoFEMErrorCode tsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec f, void *ctx)
Set implicit function for time stepping.
static MoFEMErrorCode pcSetup(PC pc)
Setup preconditioner.
static MoFEMErrorCode tsPostProc(TS ts)
Post process time step.

◆ tsSetUp() [2/3]

MoFEMErrorCode TSPrePostProc::tsSetUp ( TS  ts)

Used to setup TS solver.

Parameters
tsPETSc time stepping object to configure
Returns
MoFEMErrorCode Success/failure code

Sets up time stepping solver with custom preconditioner, monitors, and callback functions for mesh refinement and solution projection

◆ tsSetUp() [3/3]

MoFEMErrorCode TSPrePostProc::tsSetUp ( TS  ts)

Used to setup TS solver.

Parameters
ts
Returns
MoFEMErrorCode

Member Data Documentation

◆ ExRawPtr

Example* TSPrePostProc::ExRawPtr
Examples
thermoplastic.cpp.

Definition at line 601 of file thermoplastic.cpp.

◆ fsRawPtr [1/2]

Example* TSPrePostProc::fsRawPtr

◆ fsRawPtr [2/2]

FreeSurface* TSPrePostProc::fsRawPtr

Definition at line 673 of file free_surface.cpp.

◆ globRes

SmartPetscObj<Vec> TSPrePostProc::globRes
private

◆ globSol

SmartPetscObj<Vec> TSPrePostProc::globSol

◆ snesCtxPtr

boost::shared_ptr<SnesCtx> TSPrePostProc::snesCtxPtr
private

◆ solverSubDM

SmartPetscObj<DM> TSPrePostProc::solverSubDM

◆ subB

SmartPetscObj<Mat> TSPrePostProc::subB
private

◆ subKSP

SmartPetscObj<KSP> TSPrePostProc::subKSP
private

◆ tsCtxPtr

boost::shared_ptr<TsCtx> TSPrePostProc::tsCtxPtr
private

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