v0.15.5
Loading...
Searching...
No Matches
SnesCtx.hpp
Go to the documentation of this file.
1/**
2 * \file SnesCtx.hpp
3 * \brief Context for PETSc SNES, i.e. nonlinear solver
4 * \author Anonymous contributor under MIT license
5 */
6
7#ifndef __SNESCTX_HPP__
8#define __SNESCTX_HPP__
9
10namespace MoFEM {
11
12/** \brief Interface for nonlinear (SNES) solver
13 * \ingroup mofem_petsc_solvers
14 */
15struct SnesCtx {
16
20 using HookFunction = boost::function<MoFEMErrorCode(
21 SNES snes, Vec x, Vec F, boost::shared_ptr<CacheTuple> cache_ptr,
22 void *ctx)>;
23
24 SnesCtx(Interface &m_field, std::string problem_name);
25
26 virtual ~SnesCtx();
27
28 /**
29 * @return return reference to vector with FEMethod to calculate tangent
30 * matrix
31 */
33
34 /**
35 * @return return vector to vector with FEMethod to calculate residual
36 */
38
39 /**
40 * The sequence of BasicMethod is executed before residual is calculated. It
41 * can be used to setup data structures, e.g. zero global variable which is
42 * integrated in domain, e.g. for calculation of strain energy.
43 *
44 * @return reference to BasicMethod for preprocessing
45 */
47
48 /**
49 * The sequence of BasicMethod is executed after residual is calculated. It
50 * can be used to setup data structures, e.g. aggregate data from processors
51 * or to apply essential boundary conditions.
52 *
53 * @return reference to BasicMethod for postprocessing
54 */
56
57 /**
58 * @return reference to BasicMethod for preprocessing
59 */
61
62 /**
63 * The sequence of BasicMethod is executed after tangent matrix is calculated.
64 * It can be used to setup data structures, e.g. aggregate data from
65 * processors or to apply essential boundary conditions.
66 *
67 * @return reference to BasicMethod for postprocessing
68 */
70
71 /**
72 * @brief Get the finite element pipeline for FunctionFn, that calculate
73 * tangent of function load for arc length control
74 *
75 * @return FEMethodsSequence&
76 */
78
79 /**
80 * @brief Get the BasicMethod sequence for preprocessing of FunctionFn
81 *
82 * @return BasicMethodsSequence&
83 */
85
86 /**
87 * @brief Get the BasicMethod sequence for postprocessing of FunctionFn
88 *
89 * @return BasicMethodsSequence&
90 */
92
93 /**
94 * @brief Get the Load Tangent Hook function
95 *
96 * Using hook you can debug or post-process while load vector tangent is
97 * evaluated
98 *
99 * @return HookFunction&
100 */
102
103 /**
104 * @brief Get the Right Hand Side Hook function
105 *
106 * Using hook you can debug or post-process while right hand side vector is
107 * evaluated
108 *
109 * @return HookFunction&
110 */
112
113 /**
114 * \brief Copy sequences from other SNES context
115 * @param snes_ctx SNES context from which Sequence is copied from
116 * @return error code
117 */
118 MoFEMErrorCode copyLoops(const SnesCtx &snes_ctx);
119
120 /**
121 * @brief Clear loops
122 *
123 * @return MoFEMErrorCode
124 */
126
127 /**
128 * @brief Get mField reference
129 */
131
132 friend PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx);
133 friend PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx);
134
136 MatAssemblyType type);
137 friend MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh);
138
140 MatAssemblyType type);
141
142 friend MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its,
143 PetscReal fgnorm,
144 SnesCtx *snes_ctx);
145
146 friend MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its,
147 PetscReal fgnorm,
148 SnesCtx *snes_ctx);
149
150 friend PetscErrorCode SnesLoadTangent(SNES snes, Vec u, Vec F, void *ctx);
151
152 struct SnesCtxImpl;
153
154protected:
155 boost::movelib::unique_ptr<SnesCtxImpl> snesCtxImpl;
156};
157
158/**
159 * \brief This is MoFEM implementation for the right hand side (residual vector)
160 * evaluation in SNES solver
161 *
162 * For more information pleas look to PETSc manual, i.e. SNESSetFunction
163 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html>
164 *
165 * @param snes SNES solver
166 * @param x Solution vector at current iteration
167 * @param f The right hand side vector
168 * @param ctx Pointer to context i.e. SnesCtx
169 * @return Error code
170 */
171PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx);
172
173/**
174 * \brief This is MoFEM implementation for the left hand side (tangent matrix)
175 * evaluation in SNES solver
176 *
177 * For more information pleas look to PETSc manual, i.e. SNESSetJacobian
178 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian>
179 *
180 * @param snes SNES solver
181 * @param x Solution vector at current iteration
182 * @param A Tangent matrix
183 * @param B Preconditioner tangent matrix
184 * @param ctx Pointer to context i.e. SnesCtx
185 * @return Error code
186 */
187PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx);
188
189/**
190 * \brief This function calls finite element pipeline to compute tangent of of
191 * load vector in arc length control
192 *
193 * For more information pleas look to PETSc manual, i.e. SNESSetFunction
194 * <https://petsc.org/main/manualpages/SNES/SNESSetFunction/>
195 *
196 * @param snes SNES solver
197 * @param u Solution vector at current iteration
198 * @param F The right hand side vector
199 * @param ctx Pointer to context i.e. SnesCtx
200 * @return Error code
201 */
202PetscErrorCode SnesLoadTangent(SNES snes, Vec u, Vec F, void *ctx);
203
204/**
205 * \brief Set assembly type at the end of SnesMat
206 *
207 * \note Note that tangent matrix need have to have final assembly, you would
208 * use flush assembly in special case that you call SnesMat form other function
209 * set to SNESSetJacobian
210 *
211 * @param snes
212 * @param type type of assembly, either MAT_FLUSH_ASSEMBLY or
213 * MAT_FINAL_ASSEMBLY
214 * @return error code
215 */
216MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type);
217
218/**
219 * \brief Set behavior if finite element in sequence does not exist
220 * @param snes
221 * @param bh If set to MF_EXIST check if element exist, default MF_EXIST.
222 * Otherwise set MF_ZERO
223 * @return error code
224 */
226
227/**
228 * @brief Sens monitor printing residual field by field
229 *
230 */
231MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm,
232 SnesCtx *snes_ctx);
233
234/**
235 * @brief Sens monitor printing residual field by field
236 *
237 */
238MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm,
239 SnesCtx *snes_ctx);
240
241} // namespace MoFEM
242
243#endif // __SNESCTX_HPP__
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
@ F
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
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.
Definition SnesCtx.cpp:491
MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:656
MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:600
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.
Definition SnesCtx.cpp:227
MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition SnesCtx.cpp:572
MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition SnesCtx.cpp:586
std::deque< BasicMethodPtr > BasicMethodsSequence
Definition AuxPETSc.hpp:58
PetscErrorCode SnesLoadTangent(SNES snes, Vec x, Vec f, void *ctx)
This function calls finite element pipeline to compute tangent of of load vector in arc length contro...
Definition SnesCtx.cpp:359
std::deque< PairNameFEMethodPtr > FEMethodsSequence
Definition AuxPETSc.hpp:57
constexpr AssemblyType A
Deprecated interface functions.
Interface for nonlinear (SNES) solver.
Definition SnesCtx.hpp:15
MoFEM::BasicMethodsSequence BasicMethodsSequence
Definition SnesCtx.hpp:19
friend 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 contro...
Definition SnesCtx.cpp:359
FEMethodsSequence & getSetOperators()
Definition SnesCtx.cpp:133
friend MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition SnesCtx.cpp:586
friend MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:656
BasicMethodsSequence & getPreProcLoadTangent()
Get the BasicMethod sequence for preprocessing of FunctionFn.
Definition SnesCtx.cpp:161
BasicMethodsSequence & getPreProcSetOperators()
Definition SnesCtx.cpp:149
MoFEM::Interface & getMField()
Get mField reference.
Definition SnesCtx.cpp:175
friend MoFEMErrorCode SNESMoFEMSetAssmblyType(SNES snes, MatAssemblyType type)
MoFEM::FEMethodsSequence FEMethodsSequence
Definition SnesCtx.hpp:18
BasicMethodsSequence & getPostProcLoadTangent()
Get the BasicMethod sequence for postprocessing of FunctionFn.
Definition SnesCtx.cpp:165
MoFEMErrorCode copyLoops(const SnesCtx &snes_ctx)
Copy sequences from other SNES context.
Definition SnesCtx.cpp:169
boost::function< MoFEMErrorCode(SNES snes, Vec x, Vec F, boost::shared_ptr< CacheTuple > cache_ptr, void *ctx)> HookFunction
Definition SnesCtx.hpp:22
friend 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.
Definition SnesCtx.cpp:227
friend 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.
Definition SnesCtx.cpp:491
BasicMethodsSequence & getPostProcComputeRhs()
Definition SnesCtx.cpp:145
HookFunction & getLoadTangentHook()
Get the Load Tangent Hook function.
Definition SnesCtx.cpp:177
BasicMethodsSequence & getPostProcSetOperators()
Definition SnesCtx.cpp:153
boost::movelib::unique_ptr< SnesCtxImpl > snesCtxImpl
Definition SnesCtx.hpp:155
BasicMethodsSequence & getPreProcComputeRhs()
Definition SnesCtx.cpp:141
FEMethodsSequence & getLoadTangent()
Get the finite element pipeline for FunctionFn, that calculate tangent of function load for arc lengt...
Definition SnesCtx.cpp:157
FEMethodsSequence & getComputeRhs()
Definition SnesCtx.cpp:137
HookFunction & getRhsHook()
Get the Right Hand Side Hook function.
Definition SnesCtx.cpp:181
virtual ~SnesCtx()
friend MoFEMErrorCode MoFEMSNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *snes_ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:600
MoFEMErrorCode clearLoops()
Clear loops.
Definition SnesCtx.cpp:173
friend MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition SnesCtx.cpp:572