v0.14.0
Loading...
Searching...
No Matches
TsCtx.hpp
Go to the documentation of this file.
1/** \file TsCtx.hpp
2 * \brief Context for PETSc Time Stepping
3 */
4
5
6
7#ifndef __TSCTX_HPP__
8#define __TSCTX_HPP__
9
10namespace MoFEM {
11
12/** \brief Interface for Time Stepping (TS) solver
13 * \ingroup mofem_petsc_solvers
14 */
15struct TsCtx {
16
18 moab::Interface &moab;
19
20 std::string problemName;
21 MoFEMTypes bH; ///< If set to MF_EXIST check if element exist
22
26
42
44
45 TsCtx(MoFEM::Interface &m_field, const std::string &problem_name)
46 : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
47 bH(MF_EXIST), zeroMatrix(true) {
48 PetscLogEventRegister("LoopTsIFunction", 0, &MOFEM_EVENT_TsCtxIFunction);
49 PetscLogEventRegister("LoopTsIJacobian", 0, &MOFEM_EVENT_TsCtxIJacobian);
50 PetscLogEventRegister("LoopTsRHSFunction", 0,
52 PetscLogEventRegister("LoopTsRHSJacobian", 0,
54 PetscLogEventRegister("LoopTsMonitor", 0, &MOFEM_EVENT_TsCtxMonitor);
55 PetscLogEventRegister("LoopTsI2Function", 0, &MOFEM_EVENT_TsCtxI2Function);
56 PetscLogEventRegister("LoopTsI2Jacobian", 0, &MOFEM_EVENT_TsCtxI2Jacobian);
57 }
58
59 virtual ~TsCtx() = default;
60
61 /**
62 * @brief Get the loops to do IFunction object
63 *
64 * It is sequence of finite elements used to evaluate the right hand side of
65 * implicit time solver.
66 *
67 * @return FEMethodsSequence&
68 */
70
71 /**
72 * @brief Get the loops to do RHSFunction object
73 *
74 * It is sequence of finite elements used to evaluate the right hand side of
75 * implicit time solver.
76 *
77 * @return FEMethodsSequence&
78 */
80
81 /**
82 * @brief Get the loops to do IJacobian object
83 *
84 * It is sequence of finite elements used to evalite the left hand sie of
85 * implimcit time solver.
86 *
87 * @return FEMethodsSequence&
88 */
90
91 /**
92 * @brief Get the loops to do RHSJacobian object
93 *
94 * It is sequence of finite elements used to evalite the left hand sie of
95 * implimcit time solver.
96 *
97 * @return FEMethodsSequence&
98 */
100
101 /**
102 * @brief Get the loops to do Monitor object
103 *
104 * It is sequence used to monitor solution of time solver.
105 *
106 * @return FEMethodsSequence&
107 */
109
110 /**
111 * @brief Get the preProcess to do IFunction object
112 *
113 * @return BasicMethodsSequence&
114 */
116
117 /**
118 * @brief Get the postProcess to do IFunction object
119 *
120 * @return BasicMethodsSequence&
121 */
124 }
125
126 /**
127 * @brief Get the preProcess to do IJacobian object
128 *
129 * @return BasicMethodsSequence&
130 */
132
133 /**
134 * @brief Get the postProcess to do IJacobian object
135 *
136 * @return BasicMethodsSequence&
137 */
140 }
141
142 /**
143 * @brief Get the preProcess to do Monitor object
144 *
145 * @return BasicMethodsSequence&
146 */
148
149 /**
150 * @brief Get the postProcess to do Monitor object
151 *
152 * @return BasicMethodsSequence&
153 */
155
156 /**
157 * @brief Get the preProcess to do RHSJacobian object
158 *
159 * @return BasicMethodsSequence&
160 */
163 }
164
165 /**
166 * @brief Get the postProcess to do RHSJacobian object
167 *
168 * @return BasicMethodsSequence&
169 */
172 }
173
174 /**
175 * @brief Get the preProcess to do RHSFunction object
176 *
177 * @return BasicMethodsSequence&
178 */
181 }
182
183 /**
184 * @brief Get the postProcess to do RHSFunction object
185 *
186 * @return BasicMethodsSequence&
187 */
190 }
191
192 /**
193 * @brief Clear loops
194 *
195 * @return MoFEMErrorCode
196 */
198
199 friend PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t,
200 Vec F, void *ctx);
201 friend PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec U_t,
202 PetscReal a, Mat A, Mat B, void *ctx);
203 friend PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u,
204 void *ctx);
205 friend PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F,
206 void *ctx);
207 friend PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A,
208 Mat B, void *ctx);
209
210 friend PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec U, Vec U_t,
211 Vec U_tt, Vec F, void *ctx);
212
213 friend PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec U, Vec U_t,
214 Vec U_tt, PetscReal v, PetscReal a,
215 Mat J, Mat P, void *ctx);
216
217private:
225
226 boost::movelib::unique_ptr<bool> vecAssembleSwitch;
227 boost::movelib::unique_ptr<bool> matAssembleSwitch;
228
229};
230
231/**
232 * @brief Set IFunction for TS solver
233 *
234 * <a
235 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetIFunction.html>See
236 * petsc for details</a>
237 *
238 * @param ts
239 * @param t
240 * @param u
241 * @param u_t
242 * @param F
243 * @param ctx
244 * @return PetscErrorCode
245 */
246PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F,
247 void *ctx);
248
249/**
250 * @brief Set function evaluating jacobina in TS solver
251 *
252 * <a
253 * href=https://www.mcs.anl.gov/petsc/petsc-3.1/docs/manualpages/TS/TSSetIJacobian.html>See
254 * PETSc for details</a>
255 *
256 * @param ts
257 * @param t
258 * @param u
259 * @param u_t
260 * @param a
261 * @param A
262 * @param B
263 * @param ctx
264 * @return PetscErrorCode
265 */
266PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a,
267 Mat A, Mat B, void *ctx);
268
269/**
270 * @brief Set monitor for TS solver
271 *
272 * <a
273 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSMonitorSet.html>See
274 * PETSc for details</a>
275 *
276 * @param ts
277 * @param step
278 * @param t
279 * @param u
280 * @param ctx
281 * @return PetscErrorCode
282 */
283PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u,
284 void *ctx);
285
286/**
287 * @brief TS solver function
288 *
289 * <a
290 * href=https://www.mcs.anl.gov/petsc/petsc-3.11/docs/manualpages/TS/TSSetRHSFunction.html#TSSetRHSFunction>See
291 * PETSc for details</a>
292 *
293 * @param ts
294 * @param t
295 * @param u
296 * @param F
297 * @param ctx
298 * @return PetscErrorCode
299 */
300PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx);
301
302/**
303 * @brief TS solver function
304 *
305 * <a
306 * href=https://www.mcs.anl.gov/petsc/petsc-3.11/docs/manualpages/TS/TSSetRHSJacobian.html#TSSetRHSJacobian>See
307 * PETSc for details</a>
308 *
309 * @param ts
310 * @param t
311 * @param u
312 * @param A
313 * @param B
314 * @param ctx
315 * @return PetscErrorCode
316 */
317PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B,
318 void *ctx);
319
320/**
321 * @brief Calculation Jaconian for second order PDE in time
322 *
323 * <a
324 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetI2Jacobian.html>See
325 * PETSc for details</a>
326 *
327 * @param ts
328 * @param t time at step/stage being solved
329 * @param u state vectora
330 * @param u_t time derivative of state vector
331 * @param u_tt second time derivative of state vector
332 * @param a shift for u_t
333 * @param aa shift for u_tt
334 * @param A Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU +
335 * v*dF/dU_t + a*dF/dU_tt
336 * @param B preconditioning matrix for J, may be same as J
337 * @param ctx TsCtx context for matrix evaluation routine
338 * @return PetscErrorCode
339 */
340PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt,
341 PetscReal a, PetscReal aa, Mat A, Mat B,
342 void *ctx);
343
344/**
345 * @brief Calculation the right hand side for second order PDE in time
346 *
347 * <a
348 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetI2Function.html>PETSc
349 * for details</a>
350 *
351 * @param ts
352 * @param t
353 * @param u
354 * @param u_t
355 * @param u_tt
356 * @param F
357 * @param ctx
358 * @return PetscErrorCode
359 */
360PetscErrorCode TsSetI2Function(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt,
361 Vec F, void *ctx);
362
363} // namespace MoFEM
364
365#endif // __TSCTX_HPP__
static Index< 'J', 3 > J
MoFEMTypes
Those types control how functions respond on arguments, f.e. error handling.
Definition: definitions.h:97
@ MF_EXIST
Definition: definitions.h:100
@ F
const double v
phase velocity of light in medium (cm/ns)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec u_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobina in TS solver.
Definition: TsCtx.cpp:131
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:221
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.
Definition: TsCtx.cpp:566
PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition: TsCtx.cpp:281
PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:25
std::deque< BasicMethodPtr > BasicMethodsSequence
Definition: AuxPETSc.hpp:57
PetscErrorCode TsSetI2Jacobian(TS ts, PetscReal t, Vec u, Vec u_t, Vec u_tt, PetscReal a, PetscReal aa, Mat A, Mat B, void *ctx)
Calculation Jaconian for second order PDE in time.
Definition: TsCtx.cpp:469
std::deque< PairNameFEMethodPtr > FEMethodsSequence
Definition: AuxPETSc.hpp:56
PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition: TsCtx.cpp:380
constexpr AssemblyType A
constexpr double t
plate stiffness
Definition: plate.cpp:59
Deprecated interface functions.
Interface for Time Stepping (TS) solver.
Definition: TsCtx.hpp:15
BasicMethodsSequence & getPostProcessRHSFunction()
Get the postProcess to do RHSFunction object.
Definition: TsCtx.hpp:188
PetscLogEvent MOFEM_EVENT_TsCtxI2Jacobian
Definition: TsCtx.hpp:224
MoFEM::BasicMethodsSequence BasicMethodsSequence
Definition: TsCtx.hpp:25
FEMethodsSequence loopsIJacobian
Definition: TsCtx.hpp:27
bool zeroMatrix
Definition: TsCtx.hpp:43
MoFEMTypes bH
If set to MF_EXIST check if element exist.
Definition: TsCtx.hpp:21
moab::Interface & moab
Definition: TsCtx.hpp:18
PetscLogEvent MOFEM_EVENT_TsCtxRHSJacobian
Definition: TsCtx.hpp:219
BasicMethodsSequence & getPostProcessIJacobian()
Get the postProcess to do IJacobian object.
Definition: TsCtx.hpp:138
MoFEM::FEMethodsSequence FEMethodsSequence
Definition: TsCtx.hpp:24
friend PetscErrorCode TsSetIFunction(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
Set IFunction for TS solver.
Definition: TsCtx.cpp:25
FEMethodsSequence loopsRHSFunction
Definition: TsCtx.hpp:31
friend 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 Jaconian for second order PDE in time.
Definition: TsCtx.cpp:469
PetscLogEvent MOFEM_EVENT_TsCtxMonitor
Definition: TsCtx.hpp:222
BasicMethodsSequence postProcessRHSFunction
Definition: TsCtx.hpp:41
FEMethodsSequence & getLoopsMonitor()
Get the loops to do Monitor object.
Definition: TsCtx.hpp:108
BasicMethodsSequence & getPreProcessMonitor()
Get the preProcess to do Monitor object.
Definition: TsCtx.hpp:147
BasicMethodsSequence preProcessMonitor
Definition: TsCtx.hpp:36
BasicMethodsSequence & getPostProcessRHSJacobian()
Get the postProcess to do RHSJacobian object.
Definition: TsCtx.hpp:170
PetscLogEvent MOFEM_EVENT_TsCtxRHSFunction
Definition: TsCtx.hpp:218
BasicMethodsSequence & getPostProcessIFunction()
Get the postProcess to do IFunction object.
Definition: TsCtx.hpp:122
BasicMethodsSequence preProcessIJacobian
Definition: TsCtx.hpp:32
friend PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:221
BasicMethodsSequence & getPreProcessRHSFunction()
Get the preProcess to do RHSFunction object.
Definition: TsCtx.hpp:179
BasicMethodsSequence postProcessIFunction
Definition: TsCtx.hpp:35
friend PetscErrorCode TsSetRHSJacobian(TS ts, PetscReal t, Vec u, Mat A, Mat B, void *ctx)
TS solver function.
Definition: TsCtx.cpp:380
FEMethodsSequence & getLoopsIFunction()
Get the loops to do IFunction object.
Definition: TsCtx.hpp:69
FEMethodsSequence & getLoopsRHSJacobian()
Get the loops to do RHSJacobian object.
Definition: TsCtx.hpp:99
friend PetscErrorCode TsSetIJacobian(TS ts, PetscReal t, Vec u, Vec U_t, PetscReal a, Mat A, Mat B, void *ctx)
Set function evaluating jacobina in TS solver.
Definition: TsCtx.cpp:131
FEMethodsSequence loopsIFunction
Definition: TsCtx.hpp:28
FEMethodsSequence loopsRHSJacobian
Definition: TsCtx.hpp:30
BasicMethodsSequence preProcessIFunction
Definition: TsCtx.hpp:34
FEMethodsSequence loopsMonitor
Definition: TsCtx.hpp:29
virtual ~TsCtx()=default
MoFEM::PairNameFEMethodPtr PairNameFEMethodPtr
Definition: TsCtx.hpp:23
friend 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.
Definition: TsCtx.cpp:566
BasicMethodsSequence postProcessMonitor
Definition: TsCtx.hpp:37
BasicMethodsSequence preProcessRHSJacobian
Definition: TsCtx.hpp:38
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition: TsCtx.hpp:226
BasicMethodsSequence & getPreProcessIFunction()
Get the preProcess to do IFunction object.
Definition: TsCtx.hpp:115
PetscLogEvent MOFEM_EVENT_TsCtxIJacobian
Definition: TsCtx.hpp:221
FEMethodsSequence & getLoopsRHSFunction()
Get the loops to do RHSFunction object.
Definition: TsCtx.hpp:79
FEMethodsSequence & getLoopsIJacobian()
Get the loops to do IJacobian object.
Definition: TsCtx.hpp:89
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition: TsCtx.hpp:227
MoFEMErrorCode clearLoops()
Clear loops.
Definition: TsCtx.cpp:5
BasicMethodsSequence & getPreProcessRHSJacobian()
Get the preProcess to do RHSJacobian object.
Definition: TsCtx.hpp:161
std::string problemName
Definition: TsCtx.hpp:20
TsCtx(MoFEM::Interface &m_field, const std::string &problem_name)
Definition: TsCtx.hpp:45
MoFEM::Interface & mField
Definition: TsCtx.hpp:17
friend PetscErrorCode TsSetRHSFunction(TS ts, PetscReal t, Vec u, Vec F, void *ctx)
TS solver function.
Definition: TsCtx.cpp:281
BasicMethodsSequence preProcessRHSFunction
Definition: TsCtx.hpp:39
BasicMethodsSequence postProcessIJacobian
Definition: TsCtx.hpp:33
BasicMethodsSequence & getPostProcessMonitor()
Get the postProcess to do Monitor object.
Definition: TsCtx.hpp:154
PetscLogEvent MOFEM_EVENT_TsCtxIFunction
Definition: TsCtx.hpp:220
BasicMethodsSequence & getPreProcessIJacobian()
Get the preProcess to do IJacobian object.
Definition: TsCtx.hpp:131
BasicMethodsSequence postProcessRHSJacobian
Definition: TsCtx.hpp:40
PetscLogEvent MOFEM_EVENT_TsCtxI2Function
Definition: TsCtx.hpp:223