v0.15.5
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,
275
276/**
277 * \brief set ghosted vector values on all existing mesh entities
278 * \ingroup dm
279
280 * \param g vector
281 * \param mode see petsc manual for VecSetValue (ADD_VALUES or INSERT_VALUES)
282 * \param scatter_mode see petsc manual for ScatterMode (The available modes
283 are: SCATTER_FORWARD or SCATTER_REVERSE)
284 *
285 * SCATTER_REVERSE set data to field entities from V vector.
286 *
287 * SCATTER_FORWARD set vector V from data field entities
288
289 */
290PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode,
291 ScatterMode scatter_mode,
293
294/**
295 * \brief execute finite element method for each element in dm (problem)
296 * \ingroup dm
297 *
298 * @param dm Distributed mesh manager
299 * @param method Finite element method to execute for preprocessing
300 * @return Error code
301 */
302PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method);
303
304/**
305 * \brief execute finite element method for each element in dm (problem)
306 * \ingroup dm
307 *
308 * @param dm Distributed mesh manager
309 * @param method Finite element method to execute for postprocessing
310 * @return Error code
311 */
312PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method);
313
314/**
315 * \brief Executes FEMethod for finite elements in DM
316 * @param dm MoFEM discrete manager
317 * @param fe_name name of finite element
318 * @param method pointer to \ref MoFEM::FEMethod
319 * @param low_rank lowest rank of processor
320 * @param up_rank upper run of processor
321 * @return Error code
322 * \ingroup dm
323 */
325 DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank,
326 int up_rank, CacheTupleWeakPtr cache_ptr = CacheTupleSharedPtr());
327
328/**
329 * \brief Executes FEMethod for finite elements in DM
330 * @param dm MoFEM discrete manager
331 * @param fe_name name of finite element
332 * @param method pointer to \ref MoFEM::FEMethod
333 * @param low_rank lowest rank of processor
334 * @param up_rank upper run of processor
335 * @return Error code
336 * \ingroup dm
337 */
339 DM dm, const std::string fe_name, boost::shared_ptr<MoFEM::FEMethod> method,
340 int low_rank, int up_rank,
342
343/**
344 * \brief Executes FEMethod for finite elements in DM
345 * @param dm MoFEM discrete manager
346 * @param fe_name name of element
347 * @param method pointer to \ref MOFEM::FEMethod
348 * @return Error code
349 * \ingroup dm
350 */
351PetscErrorCode
352DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method,
354
355/**
356 * \brief Executes FEMethod for finite elements in DM
357 * @param dm MoFEM discrete manager
358 * @param fe_name name of element
359 * @param method pointer to \ref MOFEM::FEMethod
360 * @return Error code
361 * \ingroup dm
362 */
363PetscErrorCode
364DMoFEMLoopFiniteElements(DM dm, const std::string fe_name,
365 boost::shared_ptr<MoFEM::FEMethod> method,
367
368/**
369 * \brief execute method for dofs on field in problem
370 * \ingroup dm
371 */
372PetscErrorCode DMoFEMLoopDofs(DM dm, const char field_name[],
373 MoFEM::DofMethod *method);
374
375// /**
376// * \brief Set compute operator for KSP solver via sub-matrix and IS
377// *
378// * @param dm DM
379// * @return error code
380// *
381// * \ingroup dm
382// */
383// PetscErrorCode DMMoFEMKSPSetComputeOperatorsViaSubMatrixbByIs(DM dm);
384
385/**
386 * \brief set KSP right hand side evaluation function
387 * \ingroup dm
388 */
389PetscErrorCode DMMoFEMKSPSetComputeRHS(DM dm, const char fe_name[],
390 MoFEM::FEMethod *method,
391 MoFEM::BasicMethod *pre_only,
392 MoFEM::BasicMethod *post_only);
393
394/**
395 * \brief set KSP right hand side evaluation function
396 * \ingroup dm
397 */
398PetscErrorCode
399DMMoFEMKSPSetComputeRHS(DM dm, const std::string fe_name,
400 boost::shared_ptr<MoFEM::FEMethod> method,
401 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
402 boost::shared_ptr<MoFEM::BasicMethod> post_only);
403
404/**
405 * \brief Set KSP operators and push mofem finite element methods
406 *
407 * @param dm DM
408 * @param fe_name finite element name
409 * @param method method on the element (executed for each element in the
410 * problem which given name)
411 * @param pre_only method for pre-process before element method
412 * @param post_only method for post-process after element method
413 * @return error code
414 *
415 * \ingroup dm
416 */
417PetscErrorCode DMMoFEMKSPSetComputeOperators(DM dm, const char fe_name[],
418 MoFEM::FEMethod *method,
419 MoFEM::BasicMethod *pre_only,
420 MoFEM::BasicMethod *post_only);
421
422/**
423 * \brief Set KSP operators and push mofem finite element methods
424 *
425 * @param dm DM
426 * @param fe_name finite element name
427 * @param method method on the element (executed for each element in the
428 * problem which given name)
429 * @param pre_only method for pre-process before element method
430 * @param post_only method for post-process after element method
431 * @return error code
432 *
433 * \ingroup dm
434 */
435PetscErrorCode
436DMMoFEMKSPSetComputeOperators(DM dm, const std::string fe_name,
437 boost::shared_ptr<MoFEM::FEMethod> method,
438 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
439 boost::shared_ptr<MoFEM::BasicMethod> post_only);
440
441/**
442 * \brief set SNES residual evaluation function
443 * \ingroup dm
444 */
445PetscErrorCode DMMoFEMSNESSetFunction(DM dm, const char fe_name[],
446 MoFEM::FEMethod *method,
447 MoFEM::BasicMethod *pre_only,
448 MoFEM::BasicMethod *post_only);
449
450/**
451 * \brief set SNES residual evaluation function
452 * \ingroup dm
453 */
454PetscErrorCode
455DMMoFEMSNESSetFunction(DM dm, const std::string fe_name,
456 boost::shared_ptr<MoFEM::FEMethod> method,
457 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
458 boost::shared_ptr<MoFEM::BasicMethod> post_only);
459
460/**
461 * \brief set SNES Jacobian evaluation function
462 * \ingroup dm
463 */
464PetscErrorCode DMMoFEMSNESSetJacobian(DM dm, const char fe_name[],
465 MoFEM::FEMethod *method,
466 MoFEM::BasicMethod *pre_only,
467 MoFEM::BasicMethod *post_only);
468
469/**
470 * \brief set SNES Jacobian evaluation function
471 * \ingroup dm
472 */
473PetscErrorCode
474DMMoFEMSNESSetJacobian(DM dm, const std::string fe_name,
475 boost::shared_ptr<MoFEM::FEMethod> method,
476 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
477 boost::shared_ptr<MoFEM::BasicMethod> post_only);
478
479/**
480 * \brief set TS implicit function evaluation function
481 * \ingroup dm
482 */
483PetscErrorCode DMMoFEMTSSetIFunction(DM dm, const char fe_name[],
484 MoFEM::FEMethod *method,
485 MoFEM::BasicMethod *pre_only,
486 MoFEM::BasicMethod *post_only);
487
488/**
489 * \brief set TS implicit function evaluation function
490 * \ingroup dm
491 */
492PetscErrorCode
493DMMoFEMTSSetIFunction(DM dm, const std::string fe_name,
494 boost::shared_ptr<MoFEM::FEMethod> method,
495 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
496 boost::shared_ptr<MoFEM::BasicMethod> post_only);
497
498/**
499 * \brief set TS Jacobian evaluation function
500 * \ingroup dm
501 */
502PetscErrorCode
503DMMoFEMTSSetIJacobian(DM dm, const std::string fe_name,
504 boost::shared_ptr<MoFEM::FEMethod> method,
505 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
506 boost::shared_ptr<MoFEM::BasicMethod> post_only);
507
508/**
509 * \brief set TS Jacobian evaluation function
510 * \ingroup dm
511 */
512PetscErrorCode DMMoFEMTSSetIJacobian(DM dm, const char fe_name[],
513 MoFEM::FEMethod *method,
514 MoFEM::BasicMethod *pre_only,
515 MoFEM::BasicMethod *post_only);
516
517/**
518 * @brief set TS the right hand side function
519 *
520 * <a
521 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetRHSFunction.html#TSSetRHSFunction>See
522 * petsc documentation</a>
523 *
524 * @param dm
525 * @param fe_name
526 * @param method
527 * @param pre_only
528 * @param post_only
529 * @return PetscErrorCode
530 */
531PetscErrorCode
532DMMoFEMTSSetRHSFunction(DM dm, const std::string fe_name,
533 boost::shared_ptr<MoFEM::FEMethod> method,
534 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
535 boost::shared_ptr<MoFEM::BasicMethod> post_only);
536
537PetscErrorCode DMMoFEMTSSetRHSFunction(DM dm, const char fe_name[],
538 MoFEM::FEMethod *method,
539 MoFEM::BasicMethod *pre_only,
540 MoFEM::BasicMethod *post_only);
541
542/**
543 * @brief set TS the right hand side jacobian
544 *
545 * <a
546 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSSetRHSJacobian.html>See
547 * petsc documentation</a>
548 *
549 * @param dm
550 * @param fe_name
551 * @param method
552 * @param pre_only
553 * @param post_only
554 * @return PetscErrorCode
555 */
556PetscErrorCode
557DMMoFEMTSSetRHSJacobian(DM dm, const std::string fe_name,
558 boost::shared_ptr<MoFEM::FEMethod> method,
559 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
560 boost::shared_ptr<MoFEM::BasicMethod> post_only);
561
562PetscErrorCode DMMoFEMTSSetRHSJacobian(DM dm, const char fe_name[],
563 MoFEM::FEMethod *method,
564 MoFEM::BasicMethod *pre_only,
565 MoFEM::BasicMethod *post_only);
566
567/**
568 * \brief set TS implicit function evaluation function
569 * \ingroup dm
570 */
571PetscErrorCode
572DMMoFEMTSSetI2Function(DM dm, const std::string fe_name,
573 boost::shared_ptr<MoFEM::FEMethod> method,
574 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
575 boost::shared_ptr<MoFEM::BasicMethod> post_only);
576/**
577 * \brief set TS implicit function evaluation function
578 * \ingroup dm
579 */
580PetscErrorCode DMMoFEMTSSetI2Function(DM dm, const char fe_name[],
581 MoFEM::FEMethod *method,
582 MoFEM::BasicMethod *pre_only,
583 MoFEM::BasicMethod *post_only);
584
585/**
586 * \brief set TS Jacobian evaluation function
587 * \ingroup dm
588 */
589PetscErrorCode
590DMMoFEMTSSetI2Jacobian(DM dm, const std::string fe_name,
591 boost::shared_ptr<MoFEM::FEMethod> method,
592 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
593 boost::shared_ptr<MoFEM::BasicMethod> post_only);
594/**
595 * \brief set TS Jacobian evaluation function
596 * \ingroup dm
597 */
598PetscErrorCode
599DMMoFEMTSSetI2Jacobian(DM dm, const char fe_name[],
600 MoFEM::FEMethod *method,
601 MoFEM::BasicMethod *pre_only,
602 MoFEM::BasicMethod *post_only);
603
604/**
605 * @brief Set Monitor To TS solver
606 *
607 * <a
608 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSMonitorSet.html>See
609 * PETSc documentation here</a>
610 *
611 * @param dm
612 * @param ts time solver
613 * @param fe_name
614 * @param method
615 * @param pre_only
616 * @param post_only
617 * @return PetscErrorCode
618 */
619PetscErrorCode
620DMMoFEMTSSetMonitor(DM dm, TS ts, const std::string fe_name,
621 boost::shared_ptr<MoFEM::FEMethod> method,
622 boost::shared_ptr<MoFEM::BasicMethod> pre_only,
623 boost::shared_ptr<MoFEM::BasicMethod> post_only);
624
625/**
626 * @brief Set Monitor To TS solver
627 *
628 * <a
629 * href=https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/TS/TSMonitorSet.html>See
630 * PETSc documentation here</a>
631 *
632 * @param dm
633 * @param ts time solver
634 * @param fe_name
635 * @param method
636 * @param pre_only
637 * @param post_only
638 * @return PetscErrorCode
639 */
640PetscErrorCode DMMoFEMTSSetMonitor(DM dm, TS ts, const char fe_name[],
641 MoFEM::FEMethod *method,
642 MoFEM::BasicMethod *pre_only,
643 MoFEM::BasicMethod *post_only);
644
645/**
646 * \brief get MoFEM::KspCtx data structure
647 * \ingroup dm
648 *
649 * @param dm Distributed mesh manager
650 * @param ksp_ctx Pointer to KSP context structure
651 * @return Error code
652 */
653PetscErrorCode DMMoFEMGetKspCtx(DM dm, MoFEM::KspCtx **ksp_ctx);
654
655/**
656 * \brief get MoFEM::KspCtx data structure
657 * \ingroup dm
658 *
659 * @param dm Distributed mesh manager
660 * @param ksp_ctx Shared pointer to KSP context structure
661 * @return Error code
662 */
663PetscErrorCode
664DMMoFEMGetKspCtx(DM dm, const boost::shared_ptr<MoFEM::KspCtx> &ksp_ctx);
665
666/**
667 * \brief set MoFEM::KspCtx data structure
668 * \ingroup dm
669 *
670 * @param dm Distributed mesh manager
671 * @param ksp_ctx Shared pointer to KSP context structure to set
672 * @return Error code
673 */
674PetscErrorCode DMMoFEMSetKspCtx(DM dm,
675 boost::shared_ptr<MoFEM::KspCtx> ksp_ctx);
676
677/**
678 * \brief get MoFEM::SnesCtx data structure
679 * \ingroup dm
680 *
681 * @param dm Distributed mesh manager
682 * @param snes_ctx Pointer to SNES context structure
683 * @return Error code
684 */
685PetscErrorCode DMMoFEMGetSnesCtx(DM dm, MoFEM::SnesCtx **snes_ctx);
686
687/**
688 * \brief get MoFEM::SnesCtx data structure
689 * \ingroup dm
690 *
691 * @param dm Distributed mesh manager
692 * @param snes_ctx Shared pointer to SNES context structure
693 * @return Error code
694 */
695PetscErrorCode
696DMMoFEMGetSnesCtx(DM dm, const boost::shared_ptr<MoFEM::SnesCtx> &snes_ctx);
697
698/**
699 * \brief Set MoFEM::SnesCtx data structure
700 * \ingroup dm
701
702 */
703PetscErrorCode DMMoFEMSetSnesCtx(DM dm,
704 boost::shared_ptr<MoFEM::SnesCtx> snes_ctx);
705
706/**
707 * \brief get MoFEM::TsCtx data structure
708 * \ingroup dm
709 */
710PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx);
711
712/**
713 * \brief get MoFEM::TsCtx data structure
714 * \ingroup dm
715 */
716PetscErrorCode DMMoFEMGetTsCtx(DM dm,
717 const boost::shared_ptr<MoFEM::TsCtx> &ts_ctx);
718
719/**
720 * \brief Set MoFEM::TsCtx data structure
721 * \ingroup dm
722
723 It take over pointer, do not delete it, DM will destroy pointer
724 when is destroyed.
725
726 */
727PetscErrorCode DMMoFEMSetTsCtx(DM dm, boost::shared_ptr<MoFEM::TsCtx> ts_ctx);
728
729/** sets if read mesh is partitioned
730 * \ingroup dm
731 *
732 * @param dm Distributed mesh manager
733 * @param is_partitioned Boolean flag indicating if mesh is partitioned
734 * @return Error code
735 */
736PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned);
737
738/** get if read mesh is partitioned
739 * \ingroup dm
740 *
741 * @param dm Distributed mesh manager
742 * @param is_partitioned Pointer to boolean flag for partition status
743 * @return Error code
744 */
745PetscErrorCode DMMoFEMGetIsPartitioned(DM dm, PetscBool *is_partitioned);
746
747/**
748 * \brief Set operators for MoFEM dm
749 * @param dm Distributed mesh manager
750 * @return Error code
751 * \ingroup dm
752 */
753PetscErrorCode DMSetOperators_MoFEM(DM dm);
754
755/**
756 * \brief Create dm data structure with MoFEM data structure
757 * \ingroup dm
758 *
759 * @param dm Distributed mesh manager to create
760 * @return Error code
761 */
762PetscErrorCode DMCreate_MoFEM(DM dm);
763
764/**
765 * \brief Destroys dm with MoFEM data structure
766 * \ingroup dm
767 *
768 * @param dm Distributed mesh manager to destroy
769 * @return Error code
770 */
771PetscErrorCode DMDestroy_MoFEM(DM dm);
772
773/**
774 * \brief DMShellSetCreateGlobalVector
775 * \ingroup dm
776 *
777 * sets the routine to create a global vector
778 * associated with the shell DM
779 *
780 * @param dm Distributed mesh manager
781 * @param g Pointer to global vector to create
782 * @return Error code
783 */
784PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, Vec *g);
785
786/**
787 * \brief DMShellSetCreateGlobalVector
788 * \ingroup dm
789 *
790 * sets the routine to create a global vector
791 * associated with the shell DM
792 */
793PetscErrorCode DMCreateGlobalVector_MoFEM(DM dm, SmartPetscObj<Vec> &g_ptr,
795
796/**
797 * \brief DMShellSetCreateLocalVector
798 * \ingroup dm
799 *
800 * sets the routine to create a local vector
801 * associated with the shell DM
802 */
803PetscErrorCode DMCreateLocalVector_MoFEM(DM dm, Vec *l);
804
805/**
806 * DMShellSetCreateMatrix
807 * \ingroup dm
808 *
809 * sets the routine to create a matrix associated with the shell DM
810 */
811PetscErrorCode DMCreateMatrix_MoFEM(DM dm, Mat *M);
812
813/**
814 * DMShellSetCreateMatrix
815 * \ingroup dm
816 *
817 * sets the routine to create a matrix associated with the shell DM
818 */
819PetscErrorCode DMCreateMatrix_MoFEM(DM dm, SmartPetscObj<Mat> &M);
820
821/**
822 * Set options for MoFEM DM
823 * \ingroup dm
824 */
825#if PETSC_VERSION_GE(3, 17, 0)
826PetscErrorCode DMSetFromOptions_MoFEM(DM dm,
827 PetscOptionItems *PetscOptionsObject);
828#elif PETSC_VERSION_GE(3, 7, 0)
829PetscErrorCode DMSetFromOptions_MoFEM(PetscOptionItems *PetscOptionsObject,
830 DM dm);
831#elif PETSC_VERSION_GE(3, 5, 3)
832PetscErrorCode DMSetFromOptions_MoFEM(PetscOptions *PetscOptionsObject, DM dm);
833#else
834PetscErrorCode DMSetFromOptions_MoFEM(DM dm);
835#endif
836
837/**
838 * sets up the MoFEM structures inside a DM object
839 * \ingroup dm
840 */
841PetscErrorCode DMSetUp_MoFEM(DM dm);
842
843/**
844 * Sets up the MoFEM structures inside a DM object for sub dm
845 * \ingroup dm
846 *
847 * @param subdm Sub distributed mesh manager to set up
848 * @return Error code
849 */
850PetscErrorCode DMSubDMSetUp_MoFEM(DM subdm);
851
852/**
853 * Add field to sub dm problem on rows
854 * \ingroup dm
855 *
856 * @param dm Distributed mesh manager
857 * @param field_name Name of field to add to rows
858 * @return Error code
859 */
860PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[]);
861
862/** @copydoc DMMoFEMAddSubFieldRow(DM dm, const char field_name[]) */
863PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name);
864
865/**
866 * @brief Add field to sub dm problem on rows
867 * \ingroup dm
868 *
869 * @param dm
870 * @param field_name
871 * @param m
872 * @return PetscErrorCode
873 */
874PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[],
875 boost::shared_ptr<Range> r_ptr);
876
877/** @copydoc DMMoFEMAddSubFieldRow(DM dm, const char field_name[], boost::shared_ptr<Range> r_ptr) */
878PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, std::string field_name,
879 boost::shared_ptr<Range> r_ptr);
880
881/**
882 * Add field to sub dm problem on columns
883 * \ingroup dm
884 */
885PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[]);
886
887/** @copydoc DMMoFEMAddSubFieldCol(DM dm, const char field_name[]) */
888PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name);
889
890/**
891 * @brief Add field to sub dm problem on columns
892 *
893 * @param dm
894 * @param field_name
895 * @param range of entities
896 * @return PetscErrorCode
897 */
898PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[],
899 boost::shared_ptr<Range> r_ptr);
900
901
902/** @copydoc DMMoFEMAddSubFieldCol(DM dm, const char field_name[], boost::shared_ptr<Range> r_ptr); */
903PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, std::string field_name,
904 boost::shared_ptr<Range> r_ptr);
905
906/**
907 * Return true if this DM is sub problem
908 \ingroup dm
909 * @param dm the DM object
910 * @param is_subproblem true if subproblem
911 * @return error code
912 */
913PetscErrorCode DMMoFEMGetIsSubDM(DM dm, PetscBool *is_sub_dm);
914
915/**
916 * \brief get sub problem is
917 * @param dm has to be created with DMMoFEMSetSquareProblem
918 * @param is return is on the row
919 * @return error code
920 *
921 * Returns IS with global indices of the DM used to create SubDM
922 *
923 */
924PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is);
925
926/**
927 * \brief get sub problem is
928 * @param dm has to be created with DMMoFEMSetSquareProblem
929 * @param is return is on the row
930 * @return error code
931 *
932 * Returns IS with global indices of the DM used to create SubDM
933 *
934 */
935PetscErrorCode DMMoFEMGetSubColIS(DM dm, IS *is);
936
937/**
938 * \brief Add problem to composite DM on row
939 \ingroup dm
940 *
941 * This create block on row with DOFs from problem of given name
942 *
943 * @param dm the DM object
944 * @param prb_name add problem name
945 * @return error code
946 */
947PetscErrorCode DMMoFEMAddRowCompositeProblem(DM dm, const char prb_name[]);
948
949/**
950 * \brief Add problem to composite DM on col
951 * \ingroup dm
952 *
953 * This create block on col with DOFs from problem of given name
954 *
955 * @param dm the DM object
956 * @param prb_name add problem name
957 * @return error code
958 */
959PetscErrorCode DMMoFEMAddColCompositeProblem(DM dm, const char prb_name[]);
960
961/**
962 * \brief Get if this DM is composite DM
963 * \ingroup dm
964 *
965 * @param dm the DM object
966 * @param is_comp_dm return true if composite problem here
967 * @return error code
968 */
969PetscErrorCode DMMoFEMGetIsCompDM(DM dm, PetscBool *is_comp_dm);
970
971/**
972 * destroy the MoFEM structure
973 * \ingroup dm
974 */
975PetscErrorCode DMDestroy_MoFEM(DM dm);
976
977/**
978 * DMShellSetGlobalToLocal
979 * \ingroup dm
980 *
981 * the routine that begins the global to local scatter
982 */
983PetscErrorCode DMGlobalToLocalBegin_MoFEM(DM dm, Vec, InsertMode, Vec);
984
985/**
986 * DMShellSetGlobalToLocal
987 * \ingroup dm
988 *
989 * the routine that begins the global to local scatter
990 */
991PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec);
992
993/**
994 * DMShellSetLocalToGlobal
995 * \ingroup dm
996 *
997 * the routine that begins the local to global scatter
998 */
999PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec);
1000
1001/**
1002 * DMShellSetLocalToGlobal
1003 * \ingroup dm
1004 *
1005 * the routine that ends the local to global scatter
1006 */
1007PetscErrorCode DMLocalToGlobalBegin_MoFEM(DM, Vec, InsertMode, Vec);
1008
1009/**
1010 * DMShellSetLocalToGlobal
1011 * \ingroup dm
1012 *
1013 * the routine that ends the local to global scatter
1014 */
1015PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec);
1016
1017/**
1018 * Creates a set of IS objects with the global indices of dofs for each field
1019 * @param dm The number of fields (or NULL if not requested)
1020 *
1021 * Output:
1022 * @param numFields The number of fields (or NULL if not requested)
1023 * @param fieldNames The name for each field (or NULL if not requested)
1024 * @param fields The global indices for each field (or NULL if not
1025 requested)
1026 *
1027 * @return error code
1028
1029 * \note The user is responsible for freeing all requested arrays. In
1030 particular,
1031 * every entry of names should be freed with PetscFree(), every entry of fields
1032 * should be destroyed with ISDestroy(), and both arrays should be freed with
1033 * PetscFree().
1034
1035 \ingroup dm
1036
1037 */
1038PetscErrorCode DMCreateFieldIS_MoFEM(DM dm, PetscInt *numFields,
1039 char ***fieldNames, IS **fields);
1040
1041/**
1042 * \brief get field is in the problem
1043 * @param dm the DM objects
1044 * @param rc ROW or COL (no difference is problem is squared)
1045 * @param field_name name of the field
1046 * @param is returned the IS object
1047 * @return error code
1048 *
1049 * \code
1050 * IS is;
1051 * ierr = DMMoFEMGetFieldIS(dm,ROW,"DISP",&is_disp); CHKERRG(ierr);
1052 * \endcode
1053 *
1054 *
1055 \ingroup dm
1056 */
1057PetscErrorCode DMMoFEMGetFieldIS(DM dm, RowColData rc, const char field_name[],
1058 IS *is);
1059
1060/**
1061 * @brief Set verbosity level
1062 *
1063 * @param dm
1064 * @param verb see VERBOSITY_LEVELS for list of the levels
1065 * @return PetscErrorCode
1066 */
1067PetscErrorCode DMMoFEMSetVerbosity(DM dm, const int verb);
1068
1069struct BlockStructure;
1070
1071/** @brief Set data for block mat
1072 *
1073 * \note You can reset data by setting nullptr
1074 *
1075 */
1077 boost::shared_ptr<BlockStructure>);
1078
1079/**
1080 * @brief Get data for block mat
1081 *
1082 * @param dm
1083 * @return MoFEMErrorCode
1084 */
1085MoFEMErrorCode DMMoFEMGetBlocMatData(DM dm, boost::shared_ptr<BlockStructure> &);
1086
1087/**
1088 * @brief Create block matrix
1089 \ingroup dm
1090 *
1091 * @param dm
1092 * @param mat
1093 * @return MoFEMErrorCode
1094 */
1095MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat);
1096
1097/**
1098 * @brief Create block matrix
1099 \ingroup dm
1100 *
1101 * @param dm
1102 * @param mat smart pointer
1103 * @return MoFEMErrorCode
1104 */
1106
1107/**
1108 * @brief Set data for nest schur (see specialisation in Schur.hpp)
1109 *
1110 * \note You can reset data by setting nullptr
1111 *
1112 * @tparam T = NestSchurData
1113 * @param dm
1114 * @return MoFEMErrorCode
1115 */
1116template <typename T>
1117MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr<T>);
1118
1119/**
1120 * @brief Create nest schur matrix
1121 *
1122 * @param dm
1123 * @param mat
1124 * @return MoFEMErrorCode
1125 */
1127
1128/**
1129 * @brief Create matrix for hybridised system
1130 *
1131 * @param dm
1132 * @param mat
1133 * @return MoFEMErrorCode
1134 */
1136
1137/**
1138 * \brief PETSc Discrete Manager data structure
1139 *
1140 * This structure should not be accessed or modified by user. Is not
1141 * available from outside MoFEM DM manager. However user can inherit dat
1142 * class and add data for additional functionality.
1143 *
1144 * This is part of implementation for PETSc interface, see more details in
1145 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/DM/index.html>
1146 *
1147 * \ingroup dm
1148 *
1149 */
1150struct DMCtx : public UnknownInterface {
1151
1152 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
1153 UnknownInterface **iface) const;
1154
1155 virtual ~DMCtx() = default;
1156
1157 boost::shared_ptr<KspCtx> kspCtx; ///< data structure KSP
1158 boost::shared_ptr<SnesCtx> snesCtx; ///< data structure SNES
1159 boost::shared_ptr<TsCtx> tsCtx; ///< data structure for TS solver
1160
1161protected:
1162 DMCtx() = default;
1163};
1164
1165/**
1166 * @brief Get the Interface Ptr object
1167 *
1168 * @param dm
1169 * @return auto
1170 */
1171inline auto getInterfacePtr(DM dm) {
1172 MoFEM::Interface *m_field_ptr;
1173 CHK_THROW_MESSAGE(DMoFEMGetInterfacePtr(dm, &m_field_ptr),
1174 "Can not get interface ptr from DM");
1175 return m_field_ptr;
1176};
1177
1178/**
1179 * @brief get problem pointer from DM
1180 *
1181 */
1182inline auto getProblemPtr(DM dm) {
1183 const MoFEM::Problem *problem_ptr;
1184 CHK_THROW_MESSAGE(DMMoFEMGetProblemPtr(dm, &problem_ptr),
1185 "Get cot get problem ptr from DM");
1186 return problem_ptr;
1187};
1188
1189/**
1190 * @brief Get smart matrix from DM
1191 * \ingroup dm
1192 *
1193 */
1194inline auto createDMMatrix(DM dm) {
1197 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1198 return a;
1199};
1200
1201/**
1202 * @brief Get smart hybridised L2 matrix from DM
1203 *
1204 * @param dm
1205 * @return auto
1206 */
1207inline auto createDMHybridisedL2Matrix(DM dm) {
1210 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1211 return a;
1212};
1213
1214inline auto createDMBlockMat(DM dm) {
1215 Mat raw_a;
1216 ierr = DMMoFEMCreateBlockMat(dm, &raw_a);
1217 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1218 return SmartPetscObj<Mat>(raw_a);
1219};
1220
1221inline auto createDMNestSchurMat(DM dm) {
1222 Mat raw_a;
1223 ierr = DMMoFEMCreateNestSchurMat(dm, &raw_a);
1224 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1225 return SmartPetscObj<Mat>(raw_a);
1226};
1227
1228
1229/** @deprecated use createDMMatrix */
1230DEPRECATED inline auto smartCreateDMMatrix(DM dm) { return createDMMatrix(dm); }
1231
1232/**
1233 * @brief Get smart vector from DM
1234 * \ingroup dm
1235 *
1236 */
1240 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1241 return v;
1242};
1243
1244/** @deprecated use createDMVector*/
1245DEPRECATED inline auto smartCreateDMVector(DM dm) { return createDMVector(dm); }
1246
1247/**
1248 * @brief Get KSP context data structure used by DM
1249 *
1250 */
1251inline auto getDMKspCtx(DM dm) {
1252 boost::shared_ptr<MoFEM::KspCtx> ksp_ctx;
1253 ierr = DMMoFEMGetKspCtx(dm, ksp_ctx);
1254 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1255 return ksp_ctx;
1256};
1257
1258/** @deprecated getDMKspCtx */
1259DEPRECATED inline auto smartGetDMKspCtx(DM dm) { return getDMKspCtx(dm); }
1260
1261/**
1262 * @brief Get SNES context data structure used by DM
1263 *
1264 */
1265inline auto getDMSnesCtx(DM dm) {
1266 boost::shared_ptr<MoFEM::SnesCtx> snes_ctx;
1267 ierr = DMMoFEMGetSnesCtx(dm, snes_ctx);
1268 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1269 return snes_ctx;
1270};
1271
1272/** @deprecated use smartGetDMSnesCtx*/
1273DEPRECATED inline auto smartGetDMSnesCtx(DM dm) { return getDMSnesCtx(dm); }
1274
1275/**
1276 * @brief Get TS context data structure used by DM
1277 *
1278 */
1279inline auto getDMTsCtx(DM dm) {
1280 boost::shared_ptr<MoFEM::TsCtx> ts_ctx;
1282 CHKERRABORT(getCommFromPetscObject(reinterpret_cast<PetscObject>(dm)), ierr);
1283 return ts_ctx;
1284};
1285
1286/** @deprecated use getDMTsCtx */
1287DEPRECATED inline auto smartGetDMTsCtx(DM dm) { return getDMTsCtx(dm); }
1288
1289/**
1290 * @brief Get sub problem data structure
1291 *
1292 * @param dm
1293 * @return auto
1294 */
1295inline auto getDMSubData(DM dm) {
1296 auto prb_ptr = getProblemPtr(dm);
1297 return prb_ptr->getSubData();
1298};
1299
1300
1301} // namespace MoFEM
1302
1303#endif //__DMMOFEM_H
1304
1305/**
1306 * \defgroup dm Distributed mesh manager
1307 * \brief Implementation of PETSc DM, managing interactions between mesh data
1308 *structures and vectors and matrices
1309 *
1310 * DM objects are used to manage communication between the algebraic structures
1311 *in PETSc (Vec and Mat) and mesh data structures in PDE-based (or other)
1312 * simulations.
1313 *
1314 * DM is abstract interface, here is it particular implementation for MoFEM
1315 *code.
1316 *
1317 * \ingroup mofem
1318 **/
constexpr double a
RowColData
RowColData.
@ COL
#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:1392
PetscErrorCode DMoFEMMeshToGlobalVector(DM dm, Vec g, InsertMode mode, ScatterMode scatter_mode, RowColData rc=RowColData::COL)
set ghosted vector values on all existing mesh entities
Definition DMMoFEM.cpp:525
PetscErrorCode DMGlobalToLocalEnd_MoFEM(DM dm, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1400
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:1178
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 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:1428
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:1188
PetscErrorCode DMSetFromOptions_MoFEM(DM dm)
Definition DMMoFEM.cpp:1276
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition DMMoFEM.cpp:1132
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode, RowColData rc=RowColData::COL)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
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:1507
auto createDMVector(DM dm, RowColData rc=RowColData::COL)
Get smart vector from DM.
Definition DMMoFEM.hpp:1237
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
MoFEMErrorCode DMMoFEMCreateBlockMat(DM dm, Mat *mat)
Create block matrix.
Definition DMMoFEM.cpp:1544
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:1194
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 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:1345
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:1462
PetscErrorCode DMMoFEMUnSetElement(DM dm, std::string fe_name)
unset element from dm
Definition DMMoFEM.cpp:505
PetscErrorCode DMSetUp_MoFEM(DM dm)
Definition DMMoFEM.cpp:1298
PetscErrorCode DMLocalToGlobalEnd_MoFEM(DM, Vec, InsertMode, Vec)
Definition DMMoFEM.cpp:1457
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:1230
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:1279
DEPRECATED auto smartGetDMKspCtx(DM dm)
Definition DMMoFEM.hpp:1259
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:1575
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:1535
PetscErrorCode DMMoFEMGetSubRowIS(DM dm, IS *is)
get sub problem is
Definition DMMoFEM.cpp:330
DEPRECATED auto smartCreateDMVector(DM dm)
Definition DMMoFEM.hpp:1245
MoFEMErrorCode DMMoFEMCreateNestSchurMat(DM dm, Mat *mat)
Create nest schur matrix.
Definition DMMoFEM.cpp:1565
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:1251
auto createDMHybridisedL2Matrix(DM dm)
Get smart hybridised L2 matrix from DM.
Definition DMMoFEM.hpp:1207
auto getDMSubData(DM dm)
Get sub problem data structure.
Definition DMMoFEM.hpp:1295
MoFEMErrorCode DMMoFEMSetBlocMatData(DM dm, boost::shared_ptr< BlockStructure >)
Set data for block mat.
Definition DMMoFEM.cpp:1526
DEPRECATED auto smartGetDMTsCtx(DM dm)
Definition DMMoFEM.hpp:1287
auto getInterfacePtr(DM dm)
Get the Interface Ptr object.
Definition DMMoFEM.hpp:1171
MoFEMErrorCode DMMoFEMSetNestSchurData(DM dm, boost::shared_ptr< NestSchurData >)
Definition DMMoFEM.cpp:1555
DEPRECATED auto smartGetDMSnesCtx(DM dm)
Definition DMMoFEM.hpp:1273
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1265
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:1518
auto createDMNestSchurMat(DM dm)
Definition DMMoFEM.hpp:1221
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:1182
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:1214
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:1150
virtual ~DMCtx()=default
boost::shared_ptr< SnesCtx > snesCtx
data structure SNES
Definition DMMoFEM.hpp:1158
boost::shared_ptr< TsCtx > tsCtx
data structure for TS solver
Definition DMMoFEM.hpp:1159
boost::shared_ptr< KspCtx > kspCtx
data structure KSP
Definition DMMoFEM.hpp:1157
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