v0.15.0
Loading...
Searching...
No Matches
DMMoFEM.hpp
Go to the documentation of this file.
1/** \file DMMoFEM.hpp
2 \brief Discrete manager interface for MoFEM
3
4 This file implements PETSc's DM (Distributed Mesh) interface for MoFEM,
5 providing integration between MoFEM's Problem data structures and PETSc
6 solvers (KSP, SNES, TS, TAO).
7
8 \section dm_problem_connection DM-Problem Integration
9
10 The DM serves as a bridge between MoFEM's core Problem structure and PETSc's
11 algebraic objects (vectors and matrices). Each DM instance contains:
12 - Reference to MoFEM::Interface for mesh and field access
13 - Reference to MoFEM::Problem for DOF management and finite element definitions
14 - Context structures for solver-specific data (KspCtx, SnesCtx, TsCtx)
15
16 \section dm_solver_usage DM Usage with PETSc Solvers
17
18 **KSP (Linear Solvers):**
19 - DM provides matrix and vector creation routines
20 - Integrates with MoFEM finite element assembly for stiffness matrices
21 - Handles DOF mapping between local and global indexing
22
23 **SNES (Nonlinear Solvers):**
24 - DM manages residual and Jacobian evaluation through FEMethod callbacks
25 - Coordinates finite element loops for nonlinear assembly operations
26 - Provides context for Newton-based iteration data
27
28 **TS (Time Steppers):**
29 - DM handles time-dependent mass matrices and residual evaluation
30 - Manages implicit and explicit time integration schemes
31 - Coordinates with MoFEM field history and time-dependent boundary conditions
32
33 \section dm_workflow Typical Workflow
34
35 1. Create DM with DMMoFEMCreateMoFEM() linking to Problem
36 2. Configure solver contexts (KSP/SNES/TS) via DMMoFEMSetXxxCtx()
37 3. Set finite element methods for matrix/vector assembly
38 4. Use standard PETSc solver interfaces with MoFEM-backed DM
39 5. DM automatically handles MoFEM-specific operations during solve
40
41 The DM abstraction enables seamless use of PETSc's solver ecosystem while
42 maintaining access to MoFEM's mesh management and finite element capabilities.
43 */
44
45
46
47#ifndef __DMMOFEM_H
48#define __DMMOFEM_H
49
50#define DM_NO_ELEMENT "DMNONEFE"
51
52namespace MoFEM {
53
54/**
55 * \brief Register MoFEM problem
56 * \ingroup dm
57 *
58 * @param sname Solver name string identifier
59 * @return Error code
60 */
61PetscErrorCode DMRegister_MoFEM(const char sname[]);
62
63/**
64 * \brief Must be called by user to set MoFEM data structures
65 *
66 * \note If problem exist function create DM for it. If you set bit levels,
67 * those bits are to existing bits. Thus if you do not like to change bit ref
68 * level for existing problem, set bits to zero.
69 *
70 * \ingroup dm
71 */
72PetscErrorCode DMMoFEMCreateMoFEM(
73 DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[],
74 const MoFEM::BitRefLevel bit_level,
75 const MoFEM::BitRefLevel bit_mask = MoFEM::BitRefLevel().set());
76
77/**
78 * @brief Duplicate internal data structure
79 *
80 * @param dm
81 * @param dm_duplicate
82 * @return PetscErrorCode
83 */
84PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate);
85
86/**
87 * @brief Swap internal data structure
88 *
89 * @param dm
90 * @param dm_swap
91 * @return PetscErrorCode
92 */
93PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap);
94
95/**
96 * \brief Must be called by user to set Sub DM MoFEM data structures
97 * \ingroup dm
98 *
99 * @param subdm Sub distributed mesh manager to create
100 * @param dm Parent distributed mesh manager
101 * @param problem_name Name of the problem to create sub DM for
102 * @return Error code
103 */
104PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[]);
105
106/**
107 * \brief Get pointer to MoFEM::Interface
108 * @param dm Distributed mesh manager
109 * @param m_field_ptr Pointer to pointer of field interface
110 * @return Error code
111 * \ingroup dm
112 */
113PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr);
114
115/**
116 * \brief Get pointer to problem data structure
117 * \ingroup dm
118 *
119 * @param dm Distributed mesh manager
120 * @param problem_ptr Pointer to pointer of problem data structure
121 * @return Error code
122 */
123PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr);
124
125/**
126 * If this is set to PETSC_TRUE problem is deleted with DM
127 * @param dm the DM object
128 * @param destroy if PETSC_TRUE problem is destroyed
129 * @return error code
130 */
131PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem);
132
133/**
134 * Get if problem will be destroyed with DM
135 * @param dm the DM object
136 * @param destroy return if PETSC_TRUE problem is destroyed
137 * @return error code
138 */
139PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem);
140
141/**
142 * \brief set squared problem
143 * \ingroup dm
144
145 If true, it is assumed that matrix has the same indexing on rows and
146 columns. This reduces interprocessor communication.
147
148 */
149PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem);
150
151/**
152 * \brief get squared problem
153 * \ingroup dm
154
155 If true, it is assumed that matrix has the same indexing on rows and
156 columns. This reduces interprocessor communication.
157
158 */
159PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem);
160
161/**
162 * \brief Resolve shared entities
163 *
164 * @param dm dm
165 * @param fe_name finite element for which shared entities are resolved
166 * @return error code
167 *
168 * \note This function is valid for parallel algebra and serial mesh. It
169 * should be run collectively, i.e. on all processors.
170 *
171 * This allows for tag reduction or tag exchange, f.e.
172 *
173 * \code
174 * CHKERR DMMoFEMResolveSharedFiniteElements(dm,"SHELL_ELEMENT");
175 * Tag th;
176 * CHKERR mField.get_moab().tag_get_handle("ADAPT_ORDER",th);
177 * ParallelComm* pcomm =
178 * ParallelComm::get_pcomm(&mField.get_moab(),MYPCOMM_INDEX);
179 * // CHKERR pcomm->reduce_tags(th,MPI_SUM,prisms);
180 * CHKERR pcomm->exchange_tags(th,prisms);
181 * \endcode
182 *
183 * \ingroup dm
184 */
185PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name);
186
187/**
188 * \brief Get finite elements layout in the problem
189 *
190 * In layout is stored information how many elements are on each processor, for
191 * more information look at petsc documentation
192 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/IS/PetscLayoutCreate.html#PetscLayoutCreate>
193 *
194 * @param dm discrete manager for this problem
195 * @param fe_name finite element name
196 * @param layout pointer to layout, for created layout user takes
197 * responsibility for destroying it.
198 * @return error code
199 *
200 * \ingroup dm
201 */
202PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name,
203 PetscLayout *layout);
204
205/**
206 * \brief add element to dm
207 * \ingroup dm
208 *
209 * \note add_file is a collective, should be executed on all processors.
210 * Otherwise could lead to deadlock.
211 *
212 * @param dm Distributed mesh manager
213 * @param fe_name Finite element name
214 * @return Error code
215 */
216PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name);
217
218/**
219 * \brief add element to dm
220 * \ingroup dm
221 *
222 * \note add_file is a collective, should be executed on all processors.
223 * Otherwise could lead to deadlock.
224 *
225 * @param dm Distributed mesh manager
226 * @param fe_name Vector of finite element names
227 * @return Error code
228 */
229PetscErrorCode DMMoFEMAddElement(DM dm, std::vector<std::string> fe_name);
230
231/**
232 * \brief unset element from dm
233 * \ingroup dm
234 *
235 * @param dm Distributed mesh manager
236 * @param fe_name Finite element name to remove
237 * @return Error code
238 */
239PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name);
240
241/**
242 * \brief set local (or ghosted) vector values on mesh for partition only
243 * \ingroup dm
244
245 * \param l vector
246 * \param mode see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
247 * \param scatter_mode see petsc manual for ScatterMode (The available modes
248 are: SCATTER_FORWARD or SCATTER_REVERSE)
249 *
250 * SCATTER_REVERSE set data to field entities from V vector.
251 * SCATTER_FORWARD set vector V from data field entities
252 *
253 * This code will update ghosted entities on all processors
254 * \code
255 * auto d_vector = createDMVector(dm);
256 * // fill d_vector on owned entities
257 CHKERR DMoFEMMeshToLocalVector(dm, d_vector, INSERT_VALUES,
258 SCATTER_FORWARD);
259 // update ghost values in d_vector
260 CHKERR VecGhostUpdateBegin(d_vector, INSERT_VALUES, SCATTER_FORWARD);
261 CHKERR VecGhostUpdateEnd(d_vector, INSERT_VALUES, SCATTER_FORWARD);
262 // now d_vector has updated ghost values, and following save values of
263 // vector to fields entities
264 CHKERR DMoFEMMeshToLocalVector(dm, d_vector, INSERT_VALUES,
265 SCATTER_REVERSE);
266 * \endcode
267 *
268 * \note you use DMoFEMMeshToGlobalVector as scatter reverse
269 *
270 * For implementations see VecManager::setLocalGhostVector.
271 */
272PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode,
273 ScatterMode scatter_mode);
274
275/**
276 * \brief set ghosted vector values on all existing mesh entities
277 * \ingroup dm
278
279 * \param g vector
280 * \param mode see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
281 * \param scatter_mode see petsc manual for ScatterMode (The available modes
282 are: SCATTER_FORWARD or SCATTER_REVERSE)
283 *
284 * SCATTER_REVERSE set data to field entities from V vector.
285 *
286 * SCATTER_FORWARD set vector V from data field entities
287
288 */
289PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode,
290 ScatterMode scatter_mode);
291
292/**
293 * \brief execute finite element method for each element in dm (problem)
294 * \ingroup dm
295 *
296 * @param dm Distributed mesh manager
297 * @param method Finite element method to execute for preprocessing
298 * @return Error code
299 */
300PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method);
301
302/**
303 * \brief execute finite element method for each element in dm (problem)
304 * \ingroup dm
305 *
306 * @param dm Distributed mesh manager
307 * @param method Finite element method to execute for postprocessing
308 * @return Error code
309 */
310PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method);
311
312/**
313 * \brief Executes FEMethod for finite elements in DM
314 * @param dm MoFEM discrete manager
315 * @param fe_name name of finite element
316 * @param method pointer to \ref MoFEM::FEMethod
317 * @param low_rank lowest rank of processor
318 * @param up_rank upper run of processor
319 * @return Error code
320 * \ingroup dm
321 */
323 DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank,
324 int up_rank, CacheTupleWeakPtr cache_ptr = CacheTupleSharedPtr());
325
326/**
327 * \brief Executes FEMethod for finite elements in DM
328 * @param dm MoFEM discrete manager
329 * @param fe_name name of finite element
330 * @param method pointer to \ref MoFEM::FEMethod
331 * @param low_rank lowest rank of processor
332 * @param up_rank upper run of processor
333 * @return Error code
334 * \ingroup dm
335 */
337 DM dm, const std::string fe_name, boost::shared_ptr<MoFEM::FEMethod> method,
338 int low_rank, int up_rank,
340
341/**
342 * \brief Executes FEMethod for finite elements in DM
343 * @param dm MoFEM discrete manager
344 * @param fe_name name of element
345 * @param method pointer to \ref MOFEM::FEMethod
346 * @return Error code
347 * \ingroup dm
348 */
349PetscErrorCode
350DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method,
352
353/**
354 * \brief Executes FEMethod for finite elements in DM
355 * @param dm MoFEM discrete manager
356 * @param fe_name name of element
357 * @param method pointer to \ref MOFEM::FEMethod
358 * @return Error code
359 * \ingroup dm
360 */
361PetscErrorCode
362DMoFEMLoopFiniteElements(DM dm, const std::string fe_name,
363 boost::shared_ptr<MoFEM::FEMethod> method,
365
366/**
367 * \brief execute method for dofs on field in problem
368 * \ingroup dm
369 */
370PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[],
371 MoFEM::DofMethod *method);
372
373// /**
374// * \brief Set compute operator for KSP solver via sub-matrix and IS
375// *
376// * @param dm DM
377// * @return error code
378// *
379// * \ingroup dm
380// */
381// PetscErrorCode DMMoFEMKSPSetComputeOperatorsViaSubMatrixbByIs(DM dm);
382
383/**
384 * \brief set KSP right hand side evaluation function
385 * \ingroup dm
386 */
387PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[],
388 MoFEM::FEMethod *method,
389 MoFEM::BasicMethod *pre_only,
390 MoFEM::BasicMethod *post_only);
391
392/**
393 * \brief set KSP right hand side evaluation function
394 * \ingroup dm
395 */
396PetscErrorCode
397DMMoFEMKSPSetComputeRHS(DM dm, const std::string fe_name,
398 boost::shared_ptr<MoFEM::FEMethod> method,
399 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
400 boost::shared_ptr<MoFEM::BasicMethod> post_only);
401
402/**
403 * \brief Set KSP operators and push mofem finite element methods
404 *
405 * @param dm DM
406 * @param fe_name finite element name
407 * @param method method on the element (executed for each element in the
408 * problem which given name)
409 * @param pre_only method for pre-process before element method
410 * @param post_only method for post-process after element method
411 * @return error code
412 *
413 * \ingroup dm
414 */
415PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[],
416 MoFEM::FEMethod *method,
417 MoFEM::BasicMethod *pre_only,
418 MoFEM::BasicMethod *post_only);
419
420/**
421 * \brief Set KSP operators and push mofem finite element methods
422 *
423 * @param dm DM
424 * @param fe_name finite element name
425 * @param method method on the element (executed for each element in the
426 * problem which given name)
427 * @param pre_only method for pre-process before element method
428 * @param post_only method for post-process after element method
429 * @return error code
430 *
431 * \ingroup dm
432 */
433PetscErrorCode
434DMMoFEMKSPSetComputeOperators(DM dm, const std::string fe_name,
435 boost::shared_ptr<MoFEM::FEMethod> method,
436 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
437 boost::shared_ptr<MoFEM::BasicMethod> post_only);
438
439/**
440 * \brief set SNES residual evaluation function
441 * \ingroup dm
442 */
443PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[],
444 MoFEM::FEMethod *method,
445 MoFEM::BasicMethod *pre_only,
446 MoFEM::BasicMethod *post_only);
447
448/**
449 * \brief set SNES residual evaluation function
450 * \ingroup dm
451 */
452PetscErrorCode
453DMMoFEMSNESSetFunction(DM dm, const std::string fe_name,
454 boost::shared_ptr<MoFEM::FEMethod> method,
455 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
456 boost::shared_ptr<MoFEM::BasicMethod> post_only);
457
458/**
459 * \brief set SNES Jacobian evaluation function
460 * \ingroup dm
461 */
462PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[],
463 MoFEM::FEMethod *method,
464 MoFEM::BasicMethod *pre_only,
465 MoFEM::BasicMethod *post_only);
466
467/**
468 * \brief set SNES Jacobian evaluation function
469 * \ingroup dm
470 */
471PetscErrorCode
472DMMoFEMSNESSetJacobian(DM dm, const std::string fe_name,
473 boost::shared_ptr<MoFEM::FEMethod> method,
474 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
475 boost::shared_ptr<MoFEM::BasicMethod> post_only);
476
477/**
478 * \brief set TS implicit function evaluation function
479 * \ingroup dm
480 */
481PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[],
482 MoFEM::FEMethod *method,
483 MoFEM::BasicMethod *pre_only,
484 MoFEM::BasicMethod *post_only);
485
486/**
487 * \brief set TS implicit function evaluation function
488 * \ingroup dm
489 */
490PetscErrorCode
491DMMoFEMTSSetIFunction(DM dm, const std::string fe_name,
492 boost::shared_ptr<MoFEM::FEMethod> method,
493 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
494 boost::shared_ptr<MoFEM::BasicMethod> post_only);
495
496/**
497 * \brief set TS Jacobian evaluation function
498 * \ingroup dm
499 */
500PetscErrorCode
501DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name,
502 boost::shared_ptr<MoFEM::FEMethod> method,
503 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
504 boost::shared_ptr<MoFEM::BasicMethod> post_only);
505
506/**
507 * \brief set TS Jacobian evaluation function
508 * \ingroup dm
509 */
510PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const char fe_name[],
511 MoFEM::FEMethod *method,
512 MoFEM::BasicMethod *pre_only,
513 MoFEM::BasicMethod *post_only);
514
515/**
516 * @brief set TS the right hand side function
517 *
518 * <a
519 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetRHSFunction.html#TSSetRHSFunction>See
520 * petsc documentation</a>
521 *
522 * @param dm
523 * @param fe_name
524 * @param method
525 * @param pre_only
526 * @param post_only
527 * @return PetscErrorCode
528 */
529PetscErrorCode
530DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name,
531 boost::shared_ptr<MoFEM::FEMethod> method,
532 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
533 boost::shared_ptr<MoFEM::BasicMethod> post_only);
534
535PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const char fe_name[],
536 MoFEM::FEMethod *method,
537 MoFEM::BasicMethod *pre_only,
538 MoFEM::BasicMethod *post_only);
539
540/**
541 * @brief set TS the right hand side jacobian
542 *
543 * <a
544 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetRHSJacobian.html>See
545 * petsc documentation</a>
546 *
547 * @param dm
548 * @param fe_name
549 * @param method
550 * @param pre_only
551 * @param post_only
552 * @return PetscErrorCode
553 */
554PetscErrorCode
555DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name,
556 boost::shared_ptr<MoFEM::FEMethod> method,
557 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
558 boost::shared_ptr<MoFEM::BasicMethod> post_only);
559
560PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const char fe_name[],
561 MoFEM::FEMethod *method,
562 MoFEM::BasicMethod *pre_only,
563 MoFEM::BasicMethod *post_only);
564
565/**
566 * \brief set TS implicit function evaluation function
567 * \ingroup dm
568 */
569PetscErrorCode
570DMMoFEMTSSetI2Function(DM dm, const std::string fe_name,
571 boost::shared_ptr<MoFEM::FEMethod> method,
572 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
573 boost::shared_ptr<MoFEM::BasicMethod> post_only);
574/**
575 * \brief set TS implicit function evaluation function
576 * \ingroup dm
577 */
578PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const char fe_name[],
579 MoFEM::FEMethod *method,
580 MoFEM::BasicMethod *pre_only,
581 MoFEM::BasicMethod *post_only);
582
583/**
584 * \brief set TS Jacobian evaluation function
585 * \ingroup dm
586 */
587PetscErrorCode
588DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name,
589 boost::shared_ptr<MoFEM::FEMethod> method,
590 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
591 boost::shared_ptr<MoFEM::BasicMethod> post_only);
592/**
593 * \brief set TS Jacobian evaluation function
594 * \ingroup dm
595 */
596PetscErrorCode
597DMMoFEMTSSetI2Jacobian(DM dm, const char fe_name[],
598 MoFEM::FEMethod *method,
599 MoFEM::BasicMethod *pre_only,
600 MoFEM::BasicMethod *post_only);
601
602/**
603 * @brief Set Monitor To TS solver
604 *
605 * <a
606 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSMonitorSet.html>See
607 * PETSc documentation here</a>
608 *
609 * @param dm
610 * @param ts time solver
611 * @param fe_name
612 * @param method
613 * @param pre_only
614 * @param post_only
615 * @return PetscErrorCode
616 */
617PetscErrorCode
618DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name,
619 boost::shared_ptr<MoFEM::FEMethod> method,
620 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
621 boost::shared_ptr<MoFEM::BasicMethod> post_only);
622
623/**
624 * @brief Set Monitor To TS solver
625 *
626 * <a
627 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSMonitorSet.html>See
628 * PETSc documentation here</a>
629 *
630 * @param dm
631 * @param ts time solver
632 * @param fe_name
633 * @param method
634 * @param pre_only
635 * @param post_only
636 * @return PetscErrorCode
637 */
638PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const char fe_name[],
639 MoFEM::FEMethod *method,
640 MoFEM::BasicMethod *pre_only,
641 MoFEM::BasicMethod *post_only);
642
643/**
644 * \brief get MoFEM::KspCtx data structure
645 * \ingroup dm
646 *
647 * @param dm Distributed mesh manager
648 * @param ksp_ctx Pointer to KSP context structure
649 * @return Error code
650 */
651PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx);
652
653/**
654 * \brief get MoFEM::KspCtx data structure
655 * \ingroup dm
656 *
657 * @param dm Distributed mesh manager
658 * @param ksp_ctx Shared pointer to KSP context structure
659 * @return Error code
660 */
661PetscErrorCode
662DMMoFEMGetKspCtx(DM dm, const boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx);
663
664/**
665 * \brief set MoFEM::KspCtx data structure
666 * \ingroup dm
667 *
668 * @param dm Distributed mesh manager
669 * @param ksp_ctx Shared pointer to KSP context structure to set
670 * @return Error code
671 */
672PetscErrorCode DMMoFEMSetKspCtx(DM dm,
673 boost::shared_ptr<MoFEM::KspCtx> ksp_ctx);
674
675/**
676 * \brief get MoFEM::SnesCtx data structure
677 * \ingroup dm
678 *
679 * @param dm Distributed mesh manager
680 * @param snes_ctx Pointer to SNES context structure
681 * @return Error code
682 */
683PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx);
684
685/**
686 * \brief get MoFEM::SnesCtx data structure
687 * \ingroup dm
688 *
689 * @param dm Distributed mesh manager
690 * @param snes_ctx Shared pointer to SNES context structure
691 * @return Error code
692 */
693PetscErrorCode
694DMMoFEMGetSnesCtx(DM dm, const boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx);
695
696/**
697 * \brief Set MoFEM::SnesCtx data structure
698 * \ingroup dm
699
700 */
701PetscErrorCode DMMoFEMSetSnesCtx(DM dm,
702 boost::shared_ptr<MoFEM::SnesCtx> snes_ctx);
703
704/**
705 * \brief get MoFEM::TsCtx data structure
706 * \ingroup dm
707 */
708PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx);
709
710/**
711 * \brief get MoFEM::TsCtx data structure
712 * \ingroup dm
713 */
714PetscErrorCode DMMoFEMGetTsCtx(DM dm,
715 const boost::shared_ptr<MoFEM::TsCtx> &ts_ctx);
716
717/**
718 * \brief Set MoFEM::TsCtx data structure
719 * \ingroup dm
720
721 It take over pointer, do not delete it, DM will destroy pointer
722 when is destroyed.
723
724 */
725PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr<MoFEM::TsCtx> ts_ctx);
726
727/** sets if read mesh is partitioned
728 * \ingroup dm
729 *
730 * @param dm Distributed mesh manager
731 * @param is_partitioned Boolean flag indicating if mesh is partitioned
732 * @return Error code
733 */
734PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned);
735
736/** get if read mesh is partitioned
737 * \ingroup dm
738 *
739 * @param dm Distributed mesh manager
740 * @param is_partitioned Pointer to boolean flag for partition status
741 * @return Error code
742 */
743PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned);
744
745/**
746 * \brief Set operators for MoFEM dm
747 * @param dm Distributed mesh manager
748 * @return Error code
749 * \ingroup dm
750 */
751PetscErrorCode DMSetOperators_MoFEM(DM dm);
752
753/**
754 * \brief Create dm data structure with MoFEM data structure
755 * \ingroup dm
756 *
757 * @param dm Distributed mesh manager to create
758 * @return Error code
759 */
760PetscErrorCode DMCreate_MoFEM(DM dm);
761
762/**
763 * \brief Destroys dm with MoFEM data structure
764 * \ingroup dm
765 *
766 * @param dm Distributed mesh manager to destroy
767 * @return Error code
768 */
769PetscErrorCode DMDestroy_MoFEM(DM dm);
770
771/**
772 * \brief DMShellSetCreateGlobalVector
773 * \ingroup dm
774 *
775 * sets the routine to create a global vector
776 * associated with the shell DM
777 *
778 * @param dm Distributed mesh manager
779 * @param g Pointer to global vector to create
780 * @return Error code
781 */
782PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g);
783
784/**
785 * \brief DMShellSetCreateGlobalVector
786 * \ingroup dm
787 *
788 * sets the routine to create a global vector
789 * associated with the shell DM
790 */
791PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, SmartPetscObj<Vec> &g_ptr);
792
793/**
794 * \brief DMShellSetCreateLocalVector
795 * \ingroup dm
796 *
797 * sets the routine to create a local vector
798 * associated with the shell DM
799 */
800PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l);
801
802/**
803 * DMShellSetCreateMatrix
804 * \ingroup dm
805 *
806 * sets the routine to create a matrix associated with the shell DM
807 */
808PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M);
809
810/**
811 * DMShellSetCreateMatrix
812 * \ingroup dm
813 *
814 * sets the routine to create a matrix associated with the shell DM
815 */
816PetscErrorCode DMCreateMatrix_MoFEM(DM dm, SmartPetscObj<Mat> &M);
817
818/**
819 * Set options for MoFEM DM
820 * \ingroup dm
821 */
822#if PETSC_VERSION_GE(3, 17, 0)
823PetscErrorCode DMSetFromOptions_MoFEM(DM dm,
824 PetscOptionItems *PetscOptionsObject);
825#elif PETSC_VERSION_GE(3, 7, 0)
826PetscErrorCode DMSetFromOptions_MoFEM(PetscOptionItems *PetscOptionsObject,
827 DM dm);
828#elif PETSC_VERSION_GE(3, 5, 3)
829PetscErrorCode DMSetFromOptions_MoFEM(PetscOptions *PetscOptionsObject, DM dm);
830#else
831PetscErrorCode DMSetFromOptions_MoFEM(DM dm);
832#endif
833
834/**
835 * sets up the MoFEM structures inside a DM object
836 * \ingroup dm
837 */
838PetscErrorCode DMSetUp_MoFEM(DM dm);
839
840/**
841 * Sets up the MoFEM structures inside a DM object for sub dm
842 * \ingroup dm
843 *
844 * @param subdm Sub distributed mesh manager to set up
845 * @return Error code
846 */
847PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm);
848
849/**
850 * Add field to sub dm problem on rows
851 * \ingroup dm
852 *
853 * @param dm Distributed mesh manager
854 * @param field_name Name of field to add to rows
855 * @return Error code
856 */
857PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[]);
858
859/** @copydoc DMMoFEMAddSubFieldRow(DM dm, const char field_name[]) */
860PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name);
861
862/**
863 * @brief Add field to sub dm problem on rows
864 * \ingroup dm
865 *
866 * @param dm
867 * @param field_name
868 * @param m
869 * @return PetscErrorCode
870 */
871PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[],
872 boost::shared_ptr<Range> r_ptr);
873
874/** @copydoc DMMoFEMAddSubFieldRow(DM dm, const char field_name[], boost::shared_ptr<Range> r_ptr) */
875PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name,
876 boost::shared_ptr<Range> r_ptr);
877
878/**
879 * Add field to sub dm problem on columns
880 * \ingroup dm
881 */
882PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[]);
883
884/** @copydoc DMMoFEMAddSubFieldCol(DM dm, const char field_name[]) */
885PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name);
886
887/**
888 * @brief Add field to sub dm problem on columns
889 *
890 * @param dm
891 * @param field_name
892 * @param range of entities
893 * @return PetscErrorCode
894 */
895PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[],
896 boost::shared_ptr<Range> r_ptr);
897
898
899/** @copydoc DMMoFEMAddSubFieldCol(DM dm, const char field_name[], boost::shared_ptr<Range> r_ptr); */
900PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name,
901 boost::shared_ptr<Range> r_ptr);
902
903/**
904 * Return true if this DM is sub problem
905 \ingroup dm
906 * @param dm the DM object
907 * @param is_subproblem true if subproblem
908 * @return error code
909 */
910PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm);
911
912/**
913 * \brief get sub problem is
914 * @param dm has to be created with DMMoFEMSetSquareProblem
915 * @param is return is on the row
916 * @return error code
917 *
918 * Returns IS with global indices of the DM used to create SubDM
919 *
920 */
921PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is);
922
923/**
924 * \brief get sub problem is
925 * @param dm has to be created with DMMoFEMSetSquareProblem
926 * @param is return is on the row
927 * @return error code
928 *
929 * Returns IS with global indices of the DM used to create SubDM
930 *
931 */
932PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is);
933
934/**
935 * \brief Add problem to composite DM on row
936 \ingroup dm
937 *
938 * This create block on row with DOFs from problem of given name
939 *
940 * @param dm the DM object
941 * @param prb_name add problem name
942 * @return error code
943 */
944PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[]);
945
946/**
947 * \brief Add problem to composite DM on col
948 * \ingroup dm
949 *
950 * This create block on col with DOFs from problem of given name
951 *
952 * @param dm the DM object
953 * @param prb_name add problem name
954 * @return error code
955 */
956PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[]);
957
958/**
959 * \brief Get if this DM is composite DM
960 * \ingroup dm
961 *
962 * @param dm the DM object
963 * @param is_comp_dm return true if composite problem here
964 * @return error code
965 */
966PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm);
967
968/**
969 * destroy the MoFEM structure
970 * \ingroup dm
971 */
972PetscErrorCode DMDestroy_MoFEM(DM dm);
973
974/**
975 * DMShellSetGlobalToLocal
976 * \ingroup dm
977 *
978 * the routine that begins the global to local scatter
979 */
980PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec);
981
982/**
983 * DMShellSetGlobalToLocal
984 * \ingroup dm
985 *
986 * the routine that begins the global to local scatter
987 */
988PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec);
989
990/**
991 * DMShellSetLocalToGlobal
992 * \ingroup dm
993 *
994 * the routine that begins the local to global scatter
995 */
996PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec);
997
998/**
999 * DMShellSetLocalToGlobal
1000 * \ingroup dm
1001 *
1002 * the routine that ends the local to global scatter
1003 */
1004PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec);
1005
1006/**
1007 * DMShellSetLocalToGlobal
1008 * \ingroup dm
1009 *
1010 * the routine that ends the local to global scatter
1011 */
1012PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec);
1013
1014/**
1015 * Creates a set of IS objects with the global indices of dofs for each field
1016 * @param dm The number of fields (or NULL if not requested)
1017 *
1018 * Output:
1019 * @param numFields The number of fields (or NULL if not requested)
1020 * @param fieldNames The name for each field (or NULL if not requested)
1021 * @param fields The global indices for each field (or NULL if not
1022 requested)
1023 *
1024 * @return error code
1025
1026 * \note The user is responsible for freeing all requested arrays. In
1027 particular,
1028 * every entry of names should be freed with PetscFree(), every entry of fields
1029 * should be destroyed with ISDestroy(), and both arrays should be freed with
1030 * PetscFree().
1031
1032 \ingroup dm
1033
1034 */
1035PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields,
1036 char ***fieldNames, IS **fields);
1037
1038/**
1039 * \brief get field is in the problem
1040 * @param dm the DM objects
1041 * @param rc ROW or COL (no difference is problem is squared)
1042 * @param field_name name of the field
1043 * @param is returned the IS object
1044 * @return error code
1045 *
1046 * \code
1047 * IS is;
1048 * ierr = DMMoFEMGetFieldIS(dm,ROW,"DISP",&is_disp); CHKERRG(ierr);
1049 * \endcode
1050 *
1051 *
1052 \ingroup dm
1053 */
1054PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[],
1055 IS *is);
1056
1057/**
1058 * @brief Set verbosity level
1059 *
1060 * @param dm
1061 * @param verb see VERBOSITY_LEVELS for list of the levels
1062 * @return PetscErrorCode
1063 */
1064PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb);
1065
1066struct BlockStructure;
1067
1068/** @brief Set data for block mat
1069 *
1070 * \note You can reset data by setting nullptr
1071 *
1072 */
1074 boost::shared_ptr<BlockStructure>);
1075
1076/**
1077 * @brief Get data for block mat
1078 *
1079 * @param dm
1080 * @return MoFEMErrorCode
1081 */
1082MoFEMErrorCode DMMoFEMGetBlocMatData(DM dm, boost::shared_ptr<BlockStructure> &);
1083
1084/**
1085 * @brief Create block matrix
1086 \ingroup dm
1087 *
1088 * @param dm
1089 * @param mat
1090 * @return MoFEMErrorCode
1091 */
1092MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat);
1093
1094/**
1095 * @brief Create block matrix
1096 \ingroup dm
1097 *
1098 * @param dm
1099 * @param mat smart pointer
1100 * @return MoFEMErrorCode
1101 */
1103
1104/**
1105 * @brief Set data for nest schur (see specialisation in Schur.hpp)
1106 *
1107 * \note You can reset data by setting nullptr
1108 *
1109 * @tparam T = NestSchurData
1110 * @param dm
1111 * @return MoFEMErrorCode
1112 */
1113template <typename T>
1114MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr<T>);
1115
1116/**
1117 * @brief Create nest schur matrix
1118 *
1119 * @param dm
1120 * @param mat
1121 * @return MoFEMErrorCode
1122 */
1124
1125/**
1126 * @brief Create matrix for hybridised system
1127 *
1128 * @param dm
1129 * @param mat
1130 * @return MoFEMErrorCode
1131 */
1133
1134/**
1135 * \brief PETSc Discrete Manager data structure
1136 *
1137 * This structure should not be accessed or modified by user. Is not
1138 * available from outside MoFEM DM manager. However user can inherit dat
1139 * class and add data for additional functionality.
1140 *
1141 * This is part of implementation for PETSc interface, see more details in
1142 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/index.html>
1143 *
1144 * \ingroup dm
1145 *
1146 */
1147struct DMCtx : public UnknownInterface {
1148
1149 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
1150 UnknownInterface **iface) const;
1151
1152 virtual ~DMCtx() = default;
1153
1154 boost::shared_ptr<KspCtx> kspCtx; ///< data structure KSP
1155 boost::shared_ptr<SnesCtx> snesCtx; ///< data structure SNES
1156 boost::shared_ptr<TsCtx> tsCtx; ///< data structure for TS solver
1157
1158protected:
1159 DMCtx() = default;
1160};
1161
1162/**
1163 * @brief Get the Interface Ptr object
1164 *
1165 * @param dm
1166 * @return auto
1167 */
1168inline auto getInterfacePtr(DM dm) {
1169 MoFEM::Interface *m_field_ptr;
1170 CHK_THROW_MESSAGE(DMoFEMGetInterfacePtr(dm, &m_field_ptr),
1171 "Can not get interface ptr from DM");
1172 return m_field_ptr;
1173};
1174
1175/**
1176 * @brief get problem pointer from DM
1177 *
1178 */
1179inline auto getProblemPtr(DM dm) {
1180 const MoFEM::Problem *problem_ptr;
1181 CHK_THROW_MESSAGE(DMMoFEMGetProblemPtr(dm, &problem_ptr),
1182 "Get cot get problem ptr from DM");
1183 return problem_ptr;
1184};
1185
1186/**
1187 * @brief Get smart matrix from DM
1188 * \ingroup dm
1189 *
1190 */
1191inline auto createDMMatrix(DM dm) {
1194 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1195 return a;
1196};
1197
1198/**
1199 * @brief Get smart hybridised L2 matrix from DM
1200 *
1201 * @param dm
1202 * @return auto
1203 */
1204inline auto createDMHybridisedL2Matrix(DM dm) {
1207 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1208 return a;
1209};
1210
1211inline auto createDMBlockMat(DM dm) {
1212 Mat raw_a;
1213 ierr = DMMoFEMCreateBlockMat(dm, &raw_a);
1214 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1215 return SmartPetscObj<Mat>(raw_a);
1216};
1217
1218inline auto createDMNestSchurMat(DM dm) {
1219 Mat raw_a;
1220 ierr = DMMoFEMCreateNestSchurMat(dm, &raw_a);
1221 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1222 return SmartPetscObj<Mat>(raw_a);
1223};
1224
1225
1226/** @deprecated use createDMMatrix */
1227DEPRECATED inline auto smartCreateDMMatrix(DM dm) { return createDMMatrix(dm); }
1228
1229/**
1230 * @brief Get smart vector from DM
1231 * \ingroup dm
1232 *
1233 */
1234inline auto createDMVector(DM dm) {
1237 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1238 return v;
1239};
1240
1241/** @deprecated use createDMVector*/
1242DEPRECATED inline auto smartCreateDMVector(DM dm) { return createDMVector(dm); }
1243
1244/**
1245 * @brief Get KSP context data structure used by DM
1246 *
1247 */
1248inline auto getDMKspCtx(DM dm) {
1249 boost::shared_ptr<MoFEM::KspCtx> ksp_ctx;
1250 ierr = DMMoFEMGetKspCtx(dm, ksp_ctx);
1251 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1252 return ksp_ctx;
1253};
1254
1255/** @deprecated getDMKspCtx */
1256DEPRECATED inline auto smartGetDMKspCtx(DM dm) { return getDMKspCtx(dm); }
1257
1258/**
1259 * @brief Get SNES context data structure used by DM
1260 *
1261 */
1262inline auto getDMSnesCtx(DM dm) {
1263 boost::shared_ptr<MoFEM::SnesCtx> snes_ctx;
1264 ierr = DMMoFEMGetSnesCtx(dm, snes_ctx);
1265 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1266 return snes_ctx;
1267};
1268
1269/** @deprecated use smartGetDMSnesCtx*/
1270DEPRECATED inline auto smartGetDMSnesCtx(DM dm) { return getDMSnesCtx(dm); }
1271
1272/**
1273 * @brief Get TS context data structure used by DM
1274 *
1275 */
1276inline auto getDMTsCtx(DM dm) {
1277 boost::shared_ptr<MoFEM::TsCtx> ts_ctx;
1279 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1280 return ts_ctx;
1281};
1282
1283/** @deprecated use getDMTsCtx */
1284DEPRECATED inline auto smartGetDMTsCtx(DM dm) { return getDMTsCtx(dm); }
1285
1286/**
1287 * @brief Get sub problem data structure
1288 *
1289 * @param dm
1290 * @return auto
1291 */
1292inline auto getDMSubData(DM dm) {
1293 auto prb_ptr = getProblemPtr(dm);
1294 return prb_ptr->getSubData();
1295};
1296
1297
1298} // namespace MoFEM
1299
1300#endif //__DMMOFEM_H
1301
1302/**
1303 * \defgroup dm Distributed mesh manager
1304 * \brief Implementation of PETSc DM, managing interactions between mesh data
1305 *structures and vectors and matrices
1306 *
1307 * DM objects are used to manage communication between the algebraic structures
1308 *in PETSc (Vec and Mat) and mesh data structures in PDE-based (or other)
1309 * simulations.
1310 *
1311 * DM is abstract interface, here is it particular implementation for MoFEM
1312 *code.
1313 *
1314 * \ingroup mofem
1315 **/
constexpr double a
RowColData
RowColData.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define DEPRECATED
Definition definitions.h:17
PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1391
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1399
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition DMMoFEM.cpp:1113
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition DMMoFEM.cpp:215
PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned)
Definition DMMoFEM.cpp:1124
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES residual evaluation function
Definition DMMoFEM.cpp:708
PetscErrorCode DMMoFEMGetProblemFiniteElementLayout(DM dm, std::string fe_name, PetscLayout *layout)
Get finite elements layout in the problem.
Definition DMMoFEM.cpp:467
PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l)
DMShellSetCreateLocalVector.
Definition DMMoFEM.cpp:1177
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition DMMoFEM.cpp:450
PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set TS implicit function evaluation function
Definition DMMoFEM.cpp:790
PetscErrorCode DMDestroy_MoFEM(DM dm)
Destroys dm with MoFEM data structure.
Definition DMMoFEM.cpp:79
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition DMMoFEM.cpp:114
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:546
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on col.
Definition DMMoFEM.cpp:382
PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1427
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
Definition DMMoFEM.cpp:422
PetscErrorCode DMCreate_MoFEM(DM dm)
Create dm data structure with MoFEM data structure.
Definition DMMoFEM.cpp:71
PetscErrorCode DMMoFEMSetSnesCtx(DM dm, boost::shared_ptr< MoFEM::SnesCtx > snes_ctx)
Set MoFEM::SnesCtx data structure.
Definition DMMoFEM.cpp:1101
PetscErrorCode DMMoFEMSetKspCtx(DM dm, boost::shared_ptr< MoFEM::KspCtx > ksp_ctx)
set MoFEM::KspCtx data structure
Definition DMMoFEM.cpp:1075
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition DMMoFEM.cpp:238
PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm)
Definition DMMoFEM.cpp:322
PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[])
Add problem to composite DM on row.
Definition DMMoFEM.cpp:364
PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set SNES Jacobian evaluation function
Definition DMMoFEM.cpp:749
PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M)
Definition DMMoFEM.cpp:1187
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition DMMoFEM.cpp:1275
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition DMMoFEM.cpp:1132
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx)
get MoFEM::KspCtx data structure
Definition DMMoFEM.cpp:1058
PetscErrorCode DMMoFEMGetSquareProblem(DM dm, PetscBool *square_problem)
get squared problem
Definition DMMoFEM.cpp:480
PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[], MoFEM::DofMethod *method)
execute method for dofs on field in problem
Definition DMMoFEM.cpp:595
PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
set KSP right hand side evaluation function
Definition DMMoFEM.cpp:627
PetscErrorCode DMMoFEMResolveSharedFiniteElements(DM dm, std::string fe_name)
Resolve shared entities.
Definition DMMoFEM.cpp:458
PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr< MoFEM::TsCtx > ts_ctx)
Set MoFEM::TsCtx data structure.
Definition DMMoFEM.cpp:1149
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[], IS *is)
get field is in the problem
Definition DMMoFEM.cpp:1506
PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition DMMoFEM.cpp:843
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat)
Create block matrix.
Definition DMMoFEM.cpp:1543
PetscErrorCode DMoFEMGetInterfacePtr(DM dm, MoFEM::Interface **m_field_ptr)
Get pointer to MoFEM::Interface.
Definition DMMoFEM.cpp:410
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition DMMoFEM.cpp:280
PetscErrorCode DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS Jacobian evaluation function
Definition DMMoFEM.cpp:1007
PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS implicit function evaluation function
Definition DMMoFEM.cpp:965
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1191
PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g)
DMShellSetCreateGlobalVector.
Definition DMMoFEM.cpp:1157
PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx)
get MoFEM::SnesCtx data structure
Definition DMMoFEM.cpp:1084
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode)
set ghosted vector values on all existing mesh entities
Definition DMMoFEM.cpp:525
PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm)
Get if this DM is composite DM.
Definition DMMoFEM.cpp:402
PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm)
Definition DMMoFEM.cpp:1344
PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[], MoFEM::FEMethod *method, MoFEM::BasicMethod *pre_only, MoFEM::BasicMethod *post_only)
Set KSP operators and push mofem finite element methods.
Definition DMMoFEM.cpp:668
PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields, char ***fieldNames, IS **fields)
Definition DMMoFEM.cpp:1461
PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name)
unset element from dm
Definition DMMoFEM.cpp:505
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition DMMoFEM.cpp:1297
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1456
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:557
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:536
PetscErrorCode DMSetOperators_MoFEM(DM dm)
Set operators for MoFEM dm.
Definition DMMoFEM.cpp:49
const double v
phase velocity of light in medium (cm/ns)
MoFEM::TsCtx * ts_ctx
FTensor::Index< 'l', 3 > l
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
DEPRECATED auto smartCreateDMMatrix(DM dm)
Definition DMMoFEM.hpp:1227
MPI_Comm getCommFromPetscObject(PetscObject obj)
Get the Comm From Petsc Object object.
PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
Set Monitor To TS solver.
Definition DMMoFEM.cpp:1046
auto getDMTsCtx(DM dm)
Get TS context data structure used by DM.
Definition DMMoFEM.hpp:1276
DEPRECATED auto smartGetDMKspCtx(DM dm)
Definition DMMoFEM.hpp:1256
PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side function
Definition DMMoFEM.cpp:872
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
Definition DMMoFEM.cpp:434
MoFEMErrorCode DMMoFEMCreateHybridL2Mat(DM dm, SmartPetscObj< Mat > &mat)
Create matrix for hybridised system.
Definition DMMoFEM.cpp:1574
PetscErrorCode DMMoFEMGetDestroyProblem(DM dm, PetscBool *destroy_problem)
Definition DMMoFEM.cpp:442
MoFEMErrorCode DMMoFEMGetBlocMatData(DM dm, boost::shared_ptr< BlockStructure > &)
Get data for block mat.
Definition DMMoFEM.cpp:1534
PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is)
get sub problem is
Definition DMMoFEM.cpp:330
DEPRECATED auto smartCreateDMVector(DM dm)
Definition DMMoFEM.hpp:1242
MoFEMErrorCode DMMoFEMCreateNestSchurMat(DM dm, Mat *mat)
Create nest schur matrix.
Definition DMMoFEM.cpp:1564
PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is)
get sub problem is
Definition DMMoFEM.cpp:347
auto getDMKspCtx(DM dm)
Get KSP context data structure used by DM.
Definition DMMoFEM.hpp:1248
auto createDMHybridisedL2Matrix(DM dm)
Get smart hybridised L2 matrix from DM.
Definition DMMoFEM.hpp:1204
auto getDMSubData(DM dm)
Get sub problem data structure.
Definition DMMoFEM.hpp:1292
MoFEMErrorCode DMMoFEMSetBlocMatData(DM dm, boost::shared_ptr< BlockStructure >)
Set data for block mat.
Definition DMMoFEM.cpp:1525
DEPRECATED auto smartGetDMTsCtx(DM dm)
Definition DMMoFEM.hpp:1284
auto getInterfacePtr(DM dm)
Get the Interface Ptr object.
Definition DMMoFEM.hpp:1168
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
Definition DMMoFEM.cpp:1554
DEPRECATED auto smartGetDMSnesCtx(DM dm)
Definition DMMoFEM.hpp:1270
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1262
PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name, boost::shared_ptr< MoFEM::FEMethod > method, boost::shared_ptr< MoFEM::BasicMethod > pre_only, boost::shared_ptr< MoFEM::BasicMethod > post_only)
set TS the right hand side jacobian
Definition DMMoFEM.cpp:912
PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb)
Set verbosity level.
Definition DMMoFEM.cpp:1517
auto createDMNestSchurMat(DM dm)
Definition DMMoFEM.hpp:1218
PetscErrorCode DMMoFEMSwapDMCtx(DM dm, DM dm_swap)
Swap internal data structure.
Definition DMMoFEM.cpp:196
boost::shared_ptr< CacheTuple > CacheTupleSharedPtr
auto getProblemPtr(DM dm)
get problem pointer from DM
Definition DMMoFEM.hpp:1179
PetscErrorCode DMMoFEMDuplicateDMCtx(DM dm, DM dm_duplicate)
Duplicate internal data structure.
Definition DMMoFEM.cpp:180
boost::weak_ptr< CacheTuple > CacheTupleWeakPtr
auto createDMBlockMat(DM dm)
Definition DMMoFEM.hpp:1211
constexpr auto field_name
constexpr double g
Data structure to exchange data between MoFEM and User Loop Methods.
PETSc Discrete Manager data structure.
Definition DMMoFEM.hpp:1147
virtual ~DMCtx()=default
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
Definition DMMoFEM.hpp:1155
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
Definition DMMoFEM.hpp:1156
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
Definition DMMoFEM.hpp:1154
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition DMMoFEM.cpp:19
DMCtx()=default
Deprecated interface functions.
Data structure for user loop methods on degrees of freedom (DOFs)
Structure for user loop methods on finite elements.
Interface for linear (KSP) solver.
Definition KspCtx.hpp:14
keeps basic data about problem
intrusive_ptr for managing petsc objects
Interface for nonlinear (SNES) solver.
Definition SnesCtx.hpp:15
Interface for Time Stepping (TS) solver.
Definition TsCtx.hpp:17
base class for all interface classes