![]() |
v0.15.0 |
Interface for nonlinear (SNES) solver. More...
#include "src/petsc/SnesCtx.hpp"
Classes | |
struct | SnesCtxImpl |
Public Types | |
using | PairNameFEMethodPtr = MoFEM::PairNameFEMethodPtr |
using | FEMethodsSequence = MoFEM::FEMethodsSequence |
using | BasicMethodsSequence = MoFEM::BasicMethodsSequence |
using | HookFunction |
Public Member Functions | |
SnesCtx (Interface &m_field, std::string problem_name) | |
virtual | ~SnesCtx () |
FEMethodsSequence & | getSetOperators () |
FEMethodsSequence & | getComputeRhs () |
BasicMethodsSequence & | getPreProcComputeRhs () |
BasicMethodsSequence & | getPostProcComputeRhs () |
BasicMethodsSequence & | getPreProcSetOperators () |
BasicMethodsSequence & | getPostProcSetOperators () |
FEMethodsSequence & | getLoadTangent () |
Get the finite element pipeline for FunctionFn, that calculate tangent of function load for arc length control. | |
BasicMethodsSequence & | getPreProcLoadTangent () |
Get the BasicMethod sequence for preprocessing of FunctionFn. | |
BasicMethodsSequence & | getPostProcLoadTangent () |
Get the BasicMethod sequence for postprocessing of FunctionFn. | |
HookFunction & | getLoadTangentHook () |
Get the Load Tangent Hook function. | |
HookFunction & | getRhsHook () |
Get the Right Hand Side Hook function. | |
MoFEMErrorCode | copyLoops (const SnesCtx &snes_ctx) |
Copy sequences from other SNES context. | |
MoFEMErrorCode | clearLoops () |
Clear loops. | |
Protected Attributes | |
boost::movelib::unique_ptr< SnesCtxImpl > | snesCtxImpl |
Friends | |
PetscErrorCode | SnesRhs (SNES snes, Vec x, Vec f, void *ctx) |
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver. | |
PetscErrorCode | SnesMat (SNES snes, Vec x, Mat A, Mat B, void *ctx) |
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver. | |
MoFEMErrorCode | SNESMoFEMSetAssmblyType (SNES snes, MatAssemblyType type) |
MoFEMErrorCode | SnesMoFEMSetBehavior (SNES snes, MoFEMTypes bh) |
Set behavior if finite element in sequence does not exist. | |
MoFEMErrorCode | SnesMoFEMSetAssemblyType (SNES snes, MatAssemblyType type) |
Set assembly type at the end of SnesMat. | |
MoFEMErrorCode | MoFEMSNESMonitorFields (SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx) |
Sens monitor printing residual field by field. | |
MoFEMErrorCode | MoFEMSNESMonitorEnergy (SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx) |
Sens monitor printing residual field by field. | |
PetscErrorCode | SnesLoadTangent (SNES snes, Vec u, Vec F, void *ctx) |
This function calls finite element pipeline to compute tangent of of load vector in arc length control. | |
Interface for nonlinear (SNES) solver.
Definition at line 15 of file SnesCtx.hpp.
Definition at line 19 of file SnesCtx.hpp.
Definition at line 18 of file SnesCtx.hpp.
Definition at line 20 of file SnesCtx.hpp.
Definition at line 17 of file SnesCtx.hpp.
MoFEM::SnesCtx::SnesCtx | ( | Interface & | m_field, |
std::string | problem_name ) |
Definition at line 125 of file SnesCtx.cpp.
|
virtualdefault |
MoFEMErrorCode MoFEM::SnesCtx::clearLoops | ( | ) |
MoFEMErrorCode MoFEM::SnesCtx::copyLoops | ( | const SnesCtx & | snes_ctx | ) |
Copy sequences from other SNES context.
snes_ctx | SNES context from which Sequence is copied from |
Definition at line 167 of file SnesCtx.cpp.
FEMethodsSequence & MoFEM::SnesCtx::getComputeRhs | ( | ) |
Definition at line 135 of file SnesCtx.cpp.
FEMethodsSequence & MoFEM::SnesCtx::getLoadTangent | ( | ) |
Get the finite element pipeline for FunctionFn, that calculate tangent of function load for arc length control.
Definition at line 155 of file SnesCtx.cpp.
SnesCtx::HookFunction & MoFEM::SnesCtx::getLoadTangentHook | ( | ) |
Get the Load Tangent Hook function.
Using hook you can debug or post-process while load vector tangent is evaluated
Definition at line 173 of file SnesCtx.cpp.
BasicMethodsSequence & MoFEM::SnesCtx::getPostProcComputeRhs | ( | ) |
The sequence of BasicMethod is executed after residual is calculated. It can be used to setup data structures, e.g. aggregate data from processors or to apply essential boundary conditions.
Definition at line 143 of file SnesCtx.cpp.
BasicMethodsSequence & MoFEM::SnesCtx::getPostProcLoadTangent | ( | ) |
Get the BasicMethod sequence for postprocessing of FunctionFn.
Definition at line 163 of file SnesCtx.cpp.
BasicMethodsSequence & MoFEM::SnesCtx::getPostProcSetOperators | ( | ) |
The sequence of BasicMethod is executed after tangent matrix is calculated. It can be used to setup data structures, e.g. aggregate data from processors or to apply essential boundary conditions.
Definition at line 151 of file SnesCtx.cpp.
BasicMethodsSequence & MoFEM::SnesCtx::getPreProcComputeRhs | ( | ) |
The sequence of BasicMethod is executed before residual is calculated. It can be used to setup data structures, e.g. zero global variable which is integrated in domain, e.g. for calculation of strain energy.
Definition at line 139 of file SnesCtx.cpp.
BasicMethodsSequence & MoFEM::SnesCtx::getPreProcLoadTangent | ( | ) |
Get the BasicMethod sequence for preprocessing of FunctionFn.
Definition at line 159 of file SnesCtx.cpp.
BasicMethodsSequence & MoFEM::SnesCtx::getPreProcSetOperators | ( | ) |
Definition at line 147 of file SnesCtx.cpp.
SnesCtx::HookFunction & MoFEM::SnesCtx::getRhsHook | ( | ) |
Get the Right Hand Side Hook function.
Using hook you can debug or post-process while right hand side vector is evaluated
Definition at line 177 of file SnesCtx.cpp.
FEMethodsSequence & MoFEM::SnesCtx::getSetOperators | ( | ) |
Definition at line 131 of file SnesCtx.cpp.
|
friend |
Sens monitor printing residual field by field.
Definition at line 648 of file SnesCtx.cpp.
|
friend |
Sens monitor printing residual field by field.
Definition at line 592 of file SnesCtx.cpp.
|
friend |
This function calls finite element pipeline to compute tangent of of load vector in arc length control.
For more information pleas look to PETSc manual, i.e. SNESSetFunction https://petsc.org/main/manualpages/SNES/SNESSetFunction/
snes | SNES solver |
u | Solution vector at current iteration |
F | The right hand side vector |
ctx | Pointer to context i.e. SnesCtx |
Definition at line 353 of file SnesCtx.cpp.
|
friend |
This is MoFEM implementation for the left hand side (tangent matrix) evaluation in SNES solver.
For more information pleas look to PETSc manual, i.e. SNESSetJacobian http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian
snes | SNES solver |
x | Solution vector at current iteration |
A | Tangent matrix |
B | Preconditioner tangent matrix |
ctx | Pointer to context i.e. SnesCtx |
Definition at line 483 of file SnesCtx.cpp.
|
friend |
Set assembly type at the end of SnesMat.
snes | |
type | type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY |
Definition at line 564 of file SnesCtx.cpp.
|
friend |
|
friend |
Set behavior if finite element in sequence does not exist.
snes | |
bh | If set to MF_EXIST check if element exist, default MF_EXIST. Otherwise set MF_ZERO |
Definition at line 578 of file SnesCtx.cpp.
|
friend |
This is MoFEM implementation for the right hand side (residual vector) evaluation in SNES solver.
For more information pleas look to PETSc manual, i.e. SNESSetFunction http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html
snes | SNES solver |
x | Solution vector at current iteration |
f | The right hand side vector |
ctx | Pointer to context i.e. SnesCtx |
Definition at line 223 of file SnesCtx.cpp.
|
protected |
Definition at line 150 of file SnesCtx.hpp.