v0.14.0
PCMGSetUpViaApproxOrders.hpp
Go to the documentation of this file.
1 /** \file PCMGSetUpViaApproxOrders.hpp
2  * \brief header of multi-grid solver for p- adaptivity
3  */
4 
5 #ifndef __PCMGSETUP_VIA_APPROX_ORDERS_HPP__
6 #define __PCMGSETUP_VIA_APPROX_ORDERS_HPP__
7 
8 namespace MoFEM {
9 
10 /**
11  * \brief Set DM ordering
12  *
13  * AO will transform indices from
14  * coarseningIS ordering to ordering used to construct fine matrix.
15  *
16  * @param dm [description]
17  * @param ao [description]
18  * @return [description]
19  */
21 
22 /**
23  * \brief Push back coarsening level to MG via approximation orders
24  *
25  * @param DM discrete manager
26  * @param is Push back IS used for coarsening
27  * @param A Get sub-matrix of A using is (that sets operators for coarsening
28  * levels)
29  * @param subA Returning pointer to created sub matrix
30  * @param subA If true create sub matrix, otherwise in subA has to be valid
31  * pointer to subA
32  * @return Error code
33  *
34  * \ingroup dm
35  */
37  bool create_sub_matrix,
38  bool shell_sub_a);
39 
40 /**
41  * \brief Pop IS form MG via approximation orders
42  * @param DM dm
43  * @return error code
44  *
45  * \ingroup dm
46  */
48 
49 /**
50  * \brief Clear approximation orders
51  * @param DM dm
52  * @return Error code
53  *
54  * \ingroup dm
55  */
57 
58 /**
59  * \brief Register DM for Multi-Grid via approximation orders
60  * @param sname problem/dm registered name
61  * @return error code
62  * \ingroup dm
63  */
65 
66 /**
67  * \brief Create DM data structure for Multi-Grid via approximation orders
68  *
69  * It set data structure and operators needed
70  *
71  * @param dm Discrete manager
72  * @return Error code
73  */
75 
76 /**
77  * @brief Destroy DM
78  *
79  * @param dm
80  * @return * PetscErrorCode
81  */
82 PetscErrorCode DMDestroy_MGViaApproxOrders(DM dm);
83 
84 /**
85  * \brief Create matrix for Multi-Grid via approximation orders
86  *
87  * Not used directly by user
88  *
89  * @param dm Distributed mesh data structure
90  * @param M Matrix
91  * @return Error code
92  * \ingroup dm
93  */
95 
96 /**
97  * \brief Coarsen DM
98  *
99  * Not used directly by user
100  *
101  * @param dm Distributed mesh data structure
102  * @param comm Communicator
103  * @param dmc Coarse distributed mesh data structure
104  * @return Error code
105  *
106  * \ingroup dm
107  */
108 MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc);
109 
110 /**
111  * \brief Create interpolation matrix between data managers dm1 and dm2
112  * @param dm1 Distributed mesh data structure
113  * @param dm2 Distributed mesh data structure
114  * @param mat Pointer to returned interpolation matrix
115  * @param vec Pointer to scaling vector here returned NULL
116  * @return Error code
117  */
119  Vec *vec);
120 
121 /**
122  * \brief Create global vector for DMGViaApproxOrders
123  * @param dm Distributed mesh data structure
124  * @param g returned pointer to vector
125  * @return Error code
126  */
128 
129 struct PCMGSetUpViaApproxOrdersCtx;
130 
131 /**
132  * @brief createPCMGSetUpViaApproxOrdersCtx
133  *
134  * @param dm
135  * @param A
136  * @param use_shell_mat
137  * @return PCMGSetUpViaApproxOrdersCtx
138  */
139 boost::shared_ptr<PCMGSetUpViaApproxOrdersCtx>
140 createPCMGSetUpViaApproxOrdersCtx(DM dm, Mat A, bool use_shell_mat);
141 
142 /**
143  * \brief Function build MG structure
144  * @param pc MG pre-conditioner
145  * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html>
146  * @param ctx MoFEM data structure for MG
147  * @param verb verbosity level
148  * @return error code
149  */
151  PC pc, boost::shared_ptr<PCMGSetUpViaApproxOrdersCtx> ctx, int verb = 0);
152 
153 } // namespace MoFEM
154 
155 #endif //__PCMGSETUP_VIA_APPROX_ORDERS_HPP__
g
constexpr double g
Definition: shallow_wave.cpp:63
MoFEM::createPCMGSetUpViaApproxOrdersCtx
boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > createPCMGSetUpViaApproxOrdersCtx(DM dm, Mat A, bool use_shell_mat)
createPCMGSetUpViaApproxOrdersCtx
Definition: PCMGSetUpViaApproxOrders.cpp:630
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::DMCreateGlobalVector_MGViaApproxOrders
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders(DM dm, Vec *g)
Create global vector for DMGViaApproxOrders.
Definition: PCMGSetUpViaApproxOrders.cpp:484
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::DMDestroy_MGViaApproxOrders
PetscErrorCode DMDestroy_MGViaApproxOrders(DM dm)
Destroy DM.
Definition: PCMGSetUpViaApproxOrders.cpp:335
MoFEM::DMMGViaApproxOrdersSetAO
MoFEMErrorCode DMMGViaApproxOrdersSetAO(DM dm, AO ao)
Set DM ordering.
Definition: PCMGSetUpViaApproxOrders.cpp:219
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::DMRegister_MGViaApproxOrders
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:302
MoFEM::DMMGViaApproxOrdersClearCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS(DM dm)
Clear approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:293
MoFEM::DMMGViaApproxOrdersPopBackCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS(DM dm)
Pop IS form MG via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:281
MoFEM::DMCreateInterpolation_MGViaApproxOrders
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders(DM dm1, DM dm2, Mat *mat, Vec *vec)
Create interpolation matrix between data managers dm1 and dm2.
Definition: PCMGSetUpViaApproxOrders.cpp:391
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::DMCreate_MGViaApproxOrders
MoFEMErrorCode DMCreate_MGViaApproxOrders(DM dm)
Create DM data structure for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:313
MoFEM::DMMGViaApproxOrdersPushBackCoarseningIS
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS(DM dm, IS is, Mat A, bool create_sub_matrix, bool shell_sub_a)
Push back coarsening level to MG via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:227
MoFEM::PCMGSetUpViaApproxOrders
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, boost::shared_ptr< PCMGSetUpViaApproxOrdersCtx > ctx, int verb)
Function build MG structure.
Definition: PCMGSetUpViaApproxOrders.cpp:634
MoFEM::DMCreateMatrix_MGViaApproxOrders
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders(DM dm, Mat *M)
Create matrix for Multi-Grid via approximation orders.
Definition: PCMGSetUpViaApproxOrders.cpp:344
MoFEM::DMCoarsen_MGViaApproxOrders
MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc)
Coarsen DM.
Definition: PCMGSetUpViaApproxOrders.cpp:376