v0.13.2
Loading...
Searching...
No Matches
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
6
7#ifndef __PCMGSETUP_VIA_APPROX_ORDERS_HPP__
8#define __PCMGSETUP_VIA_APPROX_ORDERS_HPP__
10 Mat A;
11 Vec X, F;
12 IS iS;
13 VecScatter sCat;
14 PCMGSubMatrixCtx(Mat a, IS is);
15 virtual ~PCMGSubMatrixCtx();
16};
17
18/**
19 * \brief Structure for DM for multi-grid via approximation orders
20 * \ingroup dm
21 */
23
24 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
25 MoFEM::UnknownInterface **iface) const;
26
29
30 MoFEMErrorCode destroyCoarseningIS();
31
32 AO aO;
33 std::vector<IS> coarseningIS; ///< Coarsening IS
34 std::vector<Mat> kspOperators; ///< Get KSP operators
35 boost::ptr_vector<PCMGSubMatrixCtx>
36 shellMatrixCtxPtr; ///< Shell sub-matrix context
37};
38
39/**
40 * Get DM Ctx
41 */
42MoFEMErrorCode DMMGViaApproxOrdersGetCtx(DM dm, DMMGViaApproxOrdersCtx **ctx);
43
44/**
45 * \brief Set DM ordering
46 *
47 * IS can be given is some other ordering, AO will transform indices from
48 * coarseningIS ordering to ordering used to construct fine matrix.
49 *
50 * @param dm [description]
51 * @param ao [description]
52 * @return [description]
53 */
54MoFEMErrorCode DMMGViaApproxOrdersSetAO(DM dm, AO ao);
55
56/**
57 * \brief Gets size of coarseningIS in internal data struture
58 * DMMGViaApproxOrdersCtx
59 * @param dm DM
60 * @param size size of coarseningIS
61 * @return Error code
62 */
63MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize(DM dm, int *size);
64
65/**
66 * \brief Push back coarsening level to MG via approximation orders
67 *
68 * @param DM discrete manager
69 * @param is Push back IS used for coarsening
70 * @param A Get sub-matrix of A using is (that sets operators for coarsening
71 * levels)
72 * @param subA Returning pointer to created sub matrix
73 * @param subA If true create sub matrix, otherwise in subA has to be valid
74 * pointer to subA
75 * @return Error code
76 *
77 * \ingroup dm
78 */
79MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS(DM, IS is, Mat A,
80 Mat *subA,
81 bool create_sub_matrix,
82 bool shell_sub_a);
83
84/**
85 * \brief Pop is form MG via approximation orders
86 * @param DM dm
87 * @param is pop back IS
88 * @return error code
89 *
90 * \ingroup dm
91 */
93
94/**
95 * \brief Clear approximation orders
96 * @param DM dm
97 * @return Error code
98 *
99 * \ingroup dm
100 */
101MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS(DM);
102
103/**
104 * \brief Replace coarsening IS in DM via approximation orders
105 * @param dm dm
106 * @param is_vec Pointer to vector of is
107 * @param nb_elems Number of elements
108 * @param A Fine matrix
109 * @return Error code
110 */
111MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS(DM dm, IS *is_vec,
112 int nb_elems, Mat A,
113 int verb = 0);
114
115/**
116 * \brief Get context for DM via approximation orders
117 * @param dm the DM object
118 * @param ctx data context
119 * @return error code
120 */
121MoFEMErrorCode DMMGViaApproxOrdersGetCtx(DM dm,
122 const DMMGViaApproxOrdersCtx **ctx);
123
124/**
125 * \brief Register DM for Multi-Grid via approximation orders
126 * @param sname problem/dm registered name
127 * @return error code
128 * \ingroup dm
129 */
130MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[]);
131
132/**
133 * \brief Create DM data structure for Multi-Grid via approximation orders
134 *
135 * It set data structure and operators needed
136 *
137 * @param dm Discrete manager
138 * @return Error code
139 */
140MoFEMErrorCode DMCreate_MGViaApproxOrders(DM dm);
141
142/**
143 * \brief Create matrix for Multi-Grid via approximation orders
144 *
145 * Not used directly by user
146 *
147 * @param dm Distributed mesh data structure
148 * @param M Matrix
149 * @return Error code
150 * \ingroup dm
151 */
152MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders(DM dm, Mat *M);
153
154/**
155 * \brief Coarsen DM
156 *
157 * Not used directly by user
158 *
159 * @param dm Distributed mesh data structure
160 * @param comm Communicator
161 * @param dmc Coarse distributed mesh data structure
162 * @return Error code
163 *
164 * \ingroup dm
165 */
166MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc);
167
168/**
169 * \brief Create interpolation matrix between data managers dm1 and dm2
170 * @param dm1 Distributed mesh data structure
171 * @param dm2 Distributed mesh data structure
172 * @param mat Pointer to returned interpolation matrix
173 * @param vec Pointer to scaling vector here returned NULL
174 * @return Error code
175 */
176MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders(DM dm1, DM dm2, Mat *mat,
177 Vec *vec);
178
179/**
180 * \brief Create global vector for DMGViaApproxOrders
181 * @param dm Distributed mesh data structure
182 * @param g returned pointer to vector
183 * @return Error code
184 */
185MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders(DM dm, Vec *g);
186
187/**
188 * \brief Set data structures of MG pre-conditioner via approximation orders
189 */
191
192 // Interface *mFieldPtr; ///< MoFEM interface
193 // string problemName; ///< Problem name
194
195 DM dM; ///< Distributed mesh manager
196 Mat A; ///< Matrix at fine level
197
198 PCMGSetUpViaApproxOrdersCtx(DM dm, Mat a, bool shell_sub_a)
199 : dM(dm), A(a), nbLevels(2), coarseOrder(2), orderAtLastLevel(1000),
200 shellSubA(shell_sub_a), verboseLevel(0) {}
201
202 virtual ~PCMGSetUpViaApproxOrdersCtx() = default;
203
204 int nbLevels; ///< number of multi-grid levels
205 int coarseOrder; ///< approximation order of coarse level
206 int orderAtLastLevel; ///< set maximal evaluated order
207
210
211 /**
212 * \brief get options from line command
213 * @return error code
214 */
215 virtual MoFEMErrorCode getOptions();
216
217 /**
218 * \brief Set IS for levels
219 * @param kk level
220 * @param is pointer to IS
221 * @return error code
222 */
223 virtual MoFEMErrorCode createIsAtLevel(int kk, IS *is);
224
225 /**
226 * \brief Destroy IS if internally created
227 * @param kk level
228 * @param is pointer to is
229 * @return error code
230 */
231 virtual MoFEMErrorCode destroyIsAtLevel(int kk, IS *is);
232
233 /**
234 * \brief Set up data structures for MG
235 * @param pc MG pre-conditioner
236 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html>
237 * @param verb verbosity level
238 * @return error code
239 */
240 virtual MoFEMErrorCode buildProlongationOperator(bool use_mat_a,
241 int verb = 0);
242
243 // DEPRECATED virtual MoFEMErrorCode buildProlongationOperator(PC pc,int verb
244 // = 0) {
245 // return buildProlongationOperator(false,verb);
246 // }
247};
248
249/**
250 * \brief Function build MG structure
251 * @param pc MG pre-conditioner
252 * <http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCMG.html>
253 * @param ctx MoFEM data structure for MG
254 * @param verb verbosity level
255 * @return error code
256 */
258 int verb = 0);
259
260#endif //__PCMGSETUP_VIA_APPROX_ORDERS_HPP__
static Index< 'M', 3 > M
MoFEMErrorCode DMCreateGlobalVector_MGViaApproxOrders(DM dm, Vec *g)
Create global vector for DMGViaApproxOrders.
MoFEMErrorCode DMMGViaApproxOrdersSetAO(DM dm, AO ao)
Set DM ordering.
MoFEMErrorCode DMCreate_MGViaApproxOrders(DM dm)
Create DM data structure for Multi-Grid via approximation orders.
MoFEMErrorCode DMMGViaApproxOrdersGetCtx(DM dm, DMMGViaApproxOrdersCtx **ctx)
MoFEMErrorCode DMMGViaApproxOrdersReplaceCoarseningIS(DM dm, IS *is_vec, int nb_elems, Mat A, int verb=0)
Replace coarsening IS in DM via approximation orders.
MoFEMErrorCode PCMGSetUpViaApproxOrders(PC pc, PCMGSetUpViaApproxOrdersCtx *ctx, int verb=0)
Function build MG structure.
MoFEMErrorCode DMMGViaApproxOrdersGetCoarseningISSize(DM dm, int *size)
Gets size of coarseningIS in internal data struture DMMGViaApproxOrdersCtx.
MoFEMErrorCode DMCreateInterpolation_MGViaApproxOrders(DM dm1, DM dm2, Mat *mat, Vec *vec)
Create interpolation matrix between data managers dm1 and dm2.
constexpr double a
MoFEMErrorCode DMMGViaApproxOrdersPushBackCoarseningIS(DM, IS is, Mat A, Mat *subA, bool create_sub_matrix, bool shell_sub_a)
Push back coarsening level to MG via approximation orders.
MoFEMErrorCode DMCreateMatrix_MGViaApproxOrders(DM dm, Mat *M)
Create matrix for Multi-Grid via approximation orders.
MoFEMErrorCode DMCoarsen_MGViaApproxOrders(DM dm, MPI_Comm comm, DM *dmc)
Coarsen DM.
MoFEMErrorCode DMMGViaApproxOrdersPopBackCoarseningIS(DM)
Pop is form MG via approximation orders.
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
MoFEMErrorCode DMMGViaApproxOrdersClearCoarseningIS(DM)
Clear approximation orders.
constexpr AssemblyType A
constexpr double g
Structure for DM for multi-grid via approximation orders.
std::vector< Mat > kspOperators
Get KSP operators.
boost::ptr_vector< PCMGSubMatrixCtx > shellMatrixCtxPtr
Shell sub-matrix context.
std::vector< IS > coarseningIS
Coarsening IS.
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, MoFEM::UnknownInterface **iface) const
PETSc Discrete Manager data structure.
Definition: DMMoFEM.hpp:924
base class for all interface classes
Set data structures of MG pre-conditioner via approximation orders.
PCMGSetUpViaApproxOrdersCtx(DM dm, Mat a, bool shell_sub_a)
virtual MoFEMErrorCode createIsAtLevel(int kk, IS *is)
Set IS for levels.
virtual MoFEMErrorCode getOptions()
get options from line command
virtual MoFEMErrorCode destroyIsAtLevel(int kk, IS *is)
Destroy IS if internally created.
DM dM
Distributed mesh manager.
virtual ~PCMGSetUpViaApproxOrdersCtx()=default
int coarseOrder
approximation order of coarse level
int orderAtLastLevel
set maximal evaluated order
int nbLevels
number of multi-grid levels
virtual MoFEMErrorCode buildProlongationOperator(bool use_mat_a, int verb=0)
Set up data structures for MG.