v0.10.0
SnesCtx.hpp
Go to the documentation of this file.
1 /** \file SnesCtx.hpp
2  * \brief Context for PETSc SNES, i.e. nonlinear solver
3  */
4 
5 /* This file is part of MoFEM.
6  * MoFEM is free software: you can redistribute it and/or modify it under
7  * the terms of the GNU Lesser General Public License as published by the
8  * Free Software Foundation, either version 3 of the License, or (at your
9  * option) any later version.
10  *
11  * MoFEM is distributed in the hope that it will be useful, but WITHOUT
12  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
14  * License for more details.
15  *
16  * You should have received a copy of the GNU Lesser General Public
17  * License along with MoFEM. If not, see <http://www.gnu.org/licenses/>. */
18 
19 #ifndef __SNESCTX_HPP__
20 #define __SNESCTX_HPP__
21 
22 namespace MoFEM {
23 
24 /** \brief Interface for nonlinear (SNES) solver
25  * \ingroup mofem_petsc_solvers
26  */
27 struct SnesCtx {
28 
29  MoFEM::Interface &mField; ///< database Interface
30  moab::Interface &moab; ///< moab Interface
31 
32  std::string problemName; ///< problem name
34  bH; ///< If set to MF_EXIST check if element exist, default MF_EXIST
35  bool zeroPreCondMatrixB; ///< If true zero matrix, otherwise user need to do
36  ///< it, default true
37  MatAssemblyType typeOfAssembly; ///< type of assembly at the end
38  bool vErify; ///< If true verify vector
39 
43 
44  FEMethodsSequence loops_to_do_Mat; ///< Sequence of finite elements instances
45  ///< assembling tangent matrix
46  FEMethodsSequence loops_to_do_Rhs; ///< Sequence of finite elements instances
47  ///< assembling residual vector
48  BasicMethodsSequence preProcess_Mat; ///< Sequence of methods run before
49  ///< tangent matrix is assembled
50  BasicMethodsSequence postProcess_Mat; ///< Sequence of methods run after
51  ///< tangent matrix is assembled
53  preProcess_Rhs; ///< Sequence of methods run before residual is assembled
55  postProcess_Rhs; ///< Sequence of methods run after residual is assembled
56 
57  /**
58  * \brief Copy sequences from other SNES contex
59  * @param snes_ctx SNES contex from which Sequence is copied from
60  * @return error code
61  */
62  MoFEMErrorCode copyLoops(const SnesCtx &snes_ctx) {
66  preProcess_Mat = snes_ctx.preProcess_Mat;
68  preProcess_Rhs = snes_ctx.preProcess_Rhs;
71  }
72 
75  loops_to_do_Mat.clear();
76  loops_to_do_Rhs.clear();
77  preProcess_Mat.clear();
78  postProcess_Mat.clear();
79  preProcess_Rhs.clear();
80  postProcess_Rhs.clear();
82  }
83 
84  SnesCtx(Interface &m_field, const std::string &problem_name)
85  : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
87  typeOfAssembly(MAT_FINAL_ASSEMBLY), vErify(false) {
88  PetscLogEventRegister("LoopSNESRhs", 0, &MOFEM_EVENT_SnesRhs);
89  PetscLogEventRegister("LoopSNESMat", 0, &MOFEM_EVENT_SnesMat);
90  }
91 
92  virtual ~SnesCtx() = default;
93 
94  /**
95  * @return return reference to vector with FEMethod to calculate tangent
96  * matrix
97  */
99 
100  /**
101  * @return return vector to vector with FEMethod to calculate residual
102  */
104 
105  /**
106  * The sequence of BasicMethod is executed before residual is calculated. It
107  * can be used to setup data structures, e.g. zero global variable which is
108  * integrated in domain, e.g. for calculation of strain energy.
109  *
110  * @return reference to BasicMethod for preprocessing
111  */
113 
114  /**
115  * The sequence of BasicMethod is executed after residual is calculated. It
116  * can be used to setup data structures, e.g. aggregate data from processors
117  * or to apply essential boundary conditions.
118  *
119  * @return reference to BasicMethod for postprocessing
120  */
122 
123  /**
124  * @return reference to BasicMethod for preprocessing
125  */
127 
128  /**
129  * The sequence of BasicMethod is executed after tangent matrix is calculated.
130  * It can be used to setup data structures, e.g. aggregate data from
131  * processors or to apply essential boundary conditions.
132  *
133  * @return reference to BasicMethod for postprocessing
134  */
136 
137  friend PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx);
138  friend PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx);
139 
140  friend MoFEMErrorCode SNESMoFEMSetAssmblyType(SNES snes,
141  MatAssemblyType type);
142  friend MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh);
143 
144 private:
145 
146  boost::movelib::unique_ptr<bool> vecAssembleSwitch;
147  boost::movelib::unique_ptr<bool> matAssembleSwitch;
148  PetscLogEvent MOFEM_EVENT_SnesRhs; ///< Log events to assemble residual
149  PetscLogEvent MOFEM_EVENT_SnesMat; ///< Log events to assemble tangent matrix
150 
151 };
152 
153 /**
154  * \brief This is MoFEM implementation for the right hand side (residual vector)
155  * evaluation in SNES solver
156  *
157  * For more information pleas look to PETSc manual, i.e. SNESSetFunction
158  * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetFunction.html>
159  *
160  * @param snes SNES solver
161  * @param x Solution vector at current iteration
162  * @param f The right hand side vector
163  * @param ctx Pointer to context i.e. SnesCtx
164  * @return Error code
165  */
166 PetscErrorCode SnesRhs(SNES snes, Vec x, Vec f, void *ctx);
167 
168 /**
169  * \brief This is MoFEM implementation for the left hand side (tangent matrix)
170  * evaluation in SNES solver
171  *
172  * For more information pleas look to PETSc manual, i.e. SNESSetJacobian
173  * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/SNES/SNESSetJacobian.html#SNESSetJacobian>
174  *
175  * @param snes SNES solver
176  * @param x Solution vector at current iteration
177  * @param A Tangent matrix
178  * @param B Preconditioner tangent matrix
179  * @param ctx Pointer to context i.e. SnesCtx
180  * @return Error code
181  */
182 PetscErrorCode SnesMat(SNES snes, Vec x, Mat A, Mat B, void *ctx);
183 
184 /**
185  * \brief Set assembly type at the end of SnesMat
186  *
187  * \note Note that tangent matrix need have to have final assembly, you would
188  * use flush assembly in special case that you call SnesMat form other function
189  * set to SNESSetJacobian
190  *
191  * @param snes
192  * @param type type of assembly, either MAT_FLUSH_ASSEMBLY or
193  * MAT_FINAL_ASSEMBLY
194  * @return error code
195  */
196 MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type);
197 
198 /**
199  * \brief Set behavior if finite element in sequence does not exist
200  * @param snes
201  * @param bh If set to MF_EXIST check if element exist, default MF_EXIST.
202  * Otherwise set MF_ZERO
203  * @return error code
204  */
206 
207 } // namespace MoFEM
208 
209 #endif // __SNESCTX_HPP__
moab::Interface & moab
moab Interface
Definition: SnesCtx.hpp:30
MatAssemblyType typeOfAssembly
type of assembly at the end
Definition: SnesCtx.hpp:37
Deprecated interface functions.
bool zeroPreCondMatrixB
Definition: SnesCtx.hpp:35
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition: SnesCtx.hpp:149
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:126
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition: SnesCtx.hpp:148
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509
friend MoFEMErrorCode SNESMoFEMSetAssmblyType(SNES snes, MatAssemblyType type)
virtual ~SnesCtx()=default
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: SnesCtx.hpp:147
MoFEM::Interface & mField
database Interface
Definition: SnesCtx.hpp:29
std::string problemName
problem name
Definition: SnesCtx.hpp:32
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:126
MoFEM::BasicMethodsSequence BasicMethodsSequence
Definition: SnesCtx.hpp:42
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
BasicMethodsSequence & get_preProcess_to_do_Mat()
Definition: SnesCtx.hpp:126
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:21
MoFEMErrorCode SnesMoFEMSetAssemblyType(SNES snes, MatAssemblyType type)
Set assembly type at the end of SnesMat.
Definition: SnesCtx.cpp:196
BasicMethodsSequence postProcess_Rhs
Sequence of methods run after residual is assembled.
Definition: SnesCtx.hpp:55
FEMethodsSequence & get_loops_to_do_Mat()
Definition: SnesCtx.hpp:98
MoFEMErrorCode copyLoops(const SnesCtx &snes_ctx)
Copy sequences from other SNES contex.
Definition: SnesCtx.hpp:62
bool vErify
If true verify vector.
Definition: SnesCtx.hpp:38
std::vector< BasicMethodPtr > BasicMethodsSequence
Definition: AuxPETSc.hpp:112
friend MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition: SnesCtx.cpp:205
MoFEMErrorCode clearLoops()
Definition: SnesCtx.hpp:73
FEMethodsSequence loops_to_do_Rhs
Definition: SnesCtx.hpp:46
BasicMethodsSequence & get_preProcess_to_do_Rhs()
Definition: SnesCtx.hpp:112
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:67
FEMethodsSequence loops_to_do_Mat
Definition: SnesCtx.hpp:44
Interface for nonlinear (SNES) solver.
Definition: SnesCtx.hpp:27
BasicMethodsSequence & get_postProcess_to_do_Rhs()
Definition: SnesCtx.hpp:121
BasicMethodsSequence postProcess_Mat
Definition: SnesCtx.hpp:50
MoFEM::FEMethodsSequence FEMethodsSequence
Definition: SnesCtx.hpp:41
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:17
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: SnesCtx.hpp:146
FEMethodsSequence & get_loops_to_do_Rhs()
Definition: SnesCtx.hpp:103
std::vector< PairNameFEMethodPtr > FEMethodsSequence
Definition: AuxPETSc.hpp:111
BasicMethodsSequence & get_postProcess_to_do_Mat()
Definition: SnesCtx.hpp:135
BasicMethodsSequence preProcess_Rhs
Sequence of methods run before residual is assembled.
Definition: SnesCtx.hpp:53
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1943
MoFEM::PairNameFEMethodPtr PairNameFEMethodPtr
Definition: SnesCtx.hpp:40
MoFEMTypes bH
If set to MF_EXIST check if element exist, default MF_EXIST.
Definition: SnesCtx.hpp:34
SnesCtx(Interface &m_field, const std::string &problem_name)
Definition: SnesCtx.hpp:84
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:17
MoFEMErrorCode SnesMoFEMSetBehavior(SNES snes, MoFEMTypes bh)
Set behavior if finite element in sequence does not exist.
Definition: SnesCtx.cpp:205
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:189
BasicMethodsSequence preProcess_Mat
Definition: SnesCtx.hpp:48