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