xref: /libCEED/examples/fluids/src/petsc_ops.c (revision a267acd12487b2c23311fe9e8f7b0c3e959135b1)
1 // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include "../include/petsc_ops.h"
9 
10 #include <ceed.h>
11 #include <petscdm.h>
12 
13 #include "../navierstokes.h"
14 
15 // @brief Get information about a DM's local vector
16 PetscErrorCode DMGetLocalVectorInfo(DM dm, PetscInt *local_size, PetscInt *global_size, VecType *vec_type) {
17   Vec V_loc;
18 
19   PetscFunctionBeginUser;
20   PetscCall(DMGetLocalVector(dm, &V_loc));
21   if (local_size) PetscCall(VecGetLocalSize(V_loc, local_size));
22   if (global_size) PetscCall(VecGetSize(V_loc, global_size));
23   if (vec_type) PetscCall(VecGetType(V_loc, vec_type));
24   PetscCall(DMRestoreLocalVector(dm, &V_loc));
25   PetscFunctionReturn(PETSC_SUCCESS);
26 }
27 
28 // @brief Get information about a DM's global vector
29 PetscErrorCode DMGetGlobalVectorInfo(DM dm, PetscInt *local_size, PetscInt *global_size, VecType *vec_type) {
30   Vec V;
31 
32   PetscFunctionBeginUser;
33   PetscCall(DMGetGlobalVector(dm, &V));
34   if (local_size) PetscCall(VecGetLocalSize(V, local_size));
35   if (global_size) PetscCall(VecGetSize(V, global_size));
36   if (vec_type) PetscCall(VecGetType(V, vec_type));
37   PetscCall(DMRestoreGlobalVector(dm, &V));
38   PetscFunctionReturn(PETSC_SUCCESS);
39 }
40 
41 /**
42  * @brief Create OperatorApplyContext struct for applying FEM operator in a PETSc context
43  *
44  * All passed in objects are reference copied and may be destroyed if desired (with the exception of `CEED_VECTOR_NONE`).
45  * Resulting context should be destroyed with `OperatorApplyContextDestroy()`.
46  *
47  * @param[in]  dm_x     `DM` associated with the operator active input. May be `NULL`
48  * @param[in]  dm_y     `DM` associated with the operator active output. May be `NULL`
49  * @param[in]  ceed     `Ceed` object
50  * @param[in]  op_apply `CeedOperator` representing the local action of the FEM operator
51  * @param[in]  x_ceed   `CeedVector` for operator active input. May be `CEED_VECTOR_NONE` or `NULL`. If `NULL`, `CeedVector` will be automatically
52  *                      generated.
53  * @param[in]  y_ceed   `CeedVector` for operator active output. May be `CEED_VECTOR_NONE` or `NULL`. If `NULL`, `CeedVector` will be automatically
54  *                      generated.
55  * @param[in]  X_loc    Local `Vec` for operator active input. If `NULL`, vector will be obtained if needed at ApplyCeedOperator time.
56  * @param[in]  Y_loc    Local `Vec` for operator active output. If `NULL`, vector will be obtained if needed at ApplyCeedOperator time.
57  * @param[out] ctx      Struct containing all data necessary for applying the operator
58  */
59 PetscErrorCode OperatorApplyContextCreate(DM dm_x, DM dm_y, Ceed ceed, CeedOperator op_apply, CeedVector x_ceed, CeedVector y_ceed, Vec X_loc,
60                                           Vec Y_loc, OperatorApplyContext *ctx) {
61   CeedSize x_size, y_size;
62 
63   PetscFunctionBeginUser;
64   PetscCallCeed(ceed, CeedOperatorGetActiveVectorLengths(op_apply, &x_size, &y_size));
65   {  // Verify sizes
66     PetscInt X_size, Y_size, dm_X_size, dm_Y_size;
67     CeedSize x_ceed_size, y_ceed_size;
68     if (dm_x) PetscCall(DMGetLocalVectorInfo(dm_x, &dm_X_size, NULL, NULL));
69     if (dm_y) PetscCall(DMGetLocalVectorInfo(dm_y, &dm_Y_size, NULL, NULL));
70     if (X_loc) {
71       PetscCall(VecGetLocalSize(X_loc, &X_size));
72       PetscCheck(X_size == x_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
73                  "X_loc (%" PetscInt_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", X_size, x_size);
74       if (dm_x)
75         PetscCheck(X_size == dm_X_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
76                    "X_loc size (%" PetscInt_FMT ") does not match dm_x local vector size (%" PetscInt_FMT ")", X_size, dm_X_size);
77     }
78     if (Y_loc) {
79       PetscCall(VecGetLocalSize(Y_loc, &Y_size));
80       PetscCheck(Y_size == y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
81                  "Y_loc (%" PetscInt_FMT ") not correct size for CeedOperator active output size (%" CeedSize_FMT ")", Y_size, y_size);
82       if (dm_y)
83         PetscCheck(Y_size == dm_Y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
84                    "Y_loc size (%" PetscInt_FMT ") does not match dm_y local vector size (%" PetscInt_FMT ")", Y_size, dm_Y_size);
85     }
86     if (x_ceed && x_ceed != CEED_VECTOR_NONE) {
87       PetscCallCeed(ceed, CeedVectorGetLength(x_ceed, &x_ceed_size));
88       PetscCheck(x_size >= 0 ? x_ceed_size == x_size : true, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
89                  "x_ceed (%" CeedSize_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", x_ceed_size, x_size);
90       if (dm_x)
91         PetscCheck(x_ceed_size == dm_X_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
92                    "x_ceed size (%" CeedSize_FMT ") does not match dm_x local vector size (%" PetscInt_FMT ")", x_ceed_size, dm_X_size);
93     }
94     if (y_ceed && y_ceed != CEED_VECTOR_NONE) {
95       PetscCallCeed(ceed, CeedVectorGetLength(y_ceed, &y_ceed_size));
96       PetscCheck(y_ceed_size == y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
97                  "y_ceed (%" CeedSize_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", y_ceed_size, y_size);
98       if (dm_y)
99         PetscCheck(y_ceed_size == dm_Y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
100                    "y_ceed size (%" CeedSize_FMT ") does not match dm_y local vector size (%" PetscInt_FMT ")", y_ceed_size, dm_Y_size);
101     }
102   }
103 
104   PetscCall(PetscNew(ctx));
105 
106   // Copy PETSc Objects
107   if (dm_x) PetscCall(PetscObjectReference((PetscObject)dm_x));
108   (*ctx)->dm_x = dm_x;
109   if (dm_y) PetscCall(PetscObjectReference((PetscObject)dm_y));
110   (*ctx)->dm_y = dm_y;
111 
112   if (X_loc) PetscCall(PetscObjectReference((PetscObject)X_loc));
113   (*ctx)->X_loc = X_loc;
114   if (Y_loc) PetscCall(PetscObjectReference((PetscObject)Y_loc));
115   (*ctx)->Y_loc = Y_loc;
116 
117   // Copy libCEED objects
118   if (x_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(x_ceed, &(*ctx)->x_ceed));
119   else PetscCallCeed(ceed, CeedVectorCreate(ceed, x_size, &(*ctx)->x_ceed));
120 
121   if (y_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(y_ceed, &(*ctx)->y_ceed));
122   else PetscCallCeed(ceed, CeedVectorCreate(ceed, y_size, &(*ctx)->y_ceed));
123 
124   PetscCallCeed(ceed, CeedOperatorReferenceCopy(op_apply, &(*ctx)->op));
125   PetscCallCeed(ceed, CeedReferenceCopy(ceed, &(*ctx)->ceed));
126   PetscFunctionReturn(PETSC_SUCCESS);
127 }
128 
129 /**
130  * @brief Destroy OperatorApplyContext struct
131  *
132  * @param[in,out] ctx Context to destroy
133  */
134 PetscErrorCode OperatorApplyContextDestroy(OperatorApplyContext ctx) {
135   PetscFunctionBeginUser;
136   if (!ctx) PetscFunctionReturn(PETSC_SUCCESS);
137   Ceed ceed = ctx->ceed;
138 
139   // Destroy PETSc Objects
140   PetscCall(DMDestroy(&ctx->dm_x));
141   PetscCall(DMDestroy(&ctx->dm_y));
142   PetscCall(VecDestroy(&ctx->X_loc));
143   PetscCall(VecDestroy(&ctx->Y_loc));
144 
145   // Destroy libCEED Objects
146   PetscCallCeed(ceed, CeedVectorDestroy(&ctx->x_ceed));
147   PetscCallCeed(ceed, CeedVectorDestroy(&ctx->y_ceed));
148   PetscCallCeed(ceed, CeedOperatorDestroy(&ctx->op));
149   PetscCallCeed(ceed, CeedDestroy(&ctx->ceed));
150 
151   PetscCall(PetscFree(ctx));
152   PetscFunctionReturn(PETSC_SUCCESS);
153 }
154 
155 //@brief Return VecType for the given DM
156 VecType DMReturnVecType(DM dm) {
157   VecType vec_type;
158   DMGetVecType(dm, &vec_type);
159   return vec_type;
160 }
161 
162 /**
163  * @brief Create local PETSc Vecs for CeedOperator's active input/outputs
164  *
165  * This is primarily used for when the active input/ouput vector does not correspond to a `DM` object, and thus `DMCreateLocalVector` or
166  * `DMGetLocalVector` are not applicable.
167  * For example, if statitics are being store at quadrature points, a `DM`-created `Vec` will not have the
168  * correct size.
169  *
170  * @param[in]  op       Operator to make the Vecs for
171  * @param[in]  vec_type `VecType` for the new Vecs
172  * @param[in]  comm     `MPI_Comm` for the new Vecs
173  * @param[out] input    Vec for CeedOperator active input
174  * @param[out] output   Vec for CeedOperator active output
175  */
176 PetscErrorCode CeedOperatorCreateLocalVecs(CeedOperator op, VecType vec_type, MPI_Comm comm, Vec *input, Vec *output) {
177   CeedSize input_size, output_size;
178   Ceed     ceed;
179 
180   PetscFunctionBeginUser;
181   PetscCall(CeedOperatorGetCeed(op, &ceed));
182   PetscCallCeed(ceed, CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
183   if (input) {
184     PetscCall(VecCreate(comm, input));
185     PetscCall(VecSetType(*input, vec_type));
186     PetscCall(VecSetSizes(*input, input_size, input_size));
187   }
188   if (output) {
189     PetscCall(VecCreate(comm, output));
190     PetscCall(VecSetType(*output, vec_type));
191     PetscCall(VecSetSizes(*output, output_size, output_size));
192   }
193   PetscFunctionReturn(PETSC_SUCCESS);
194 }
195 
196 /**
197  * @brief Apply FEM Operator defined by `OperatorApplyContext` to various input and output vectors
198  *
199  * @param[in]     X             Input global `Vec`, maybe `NULL`
200  * @param[in]     X_loc         Input local `Vec`, maybe `NULL`
201  * @param[in]     x_ceed        Input `CeedVector`, maybe `CEED_VECTOR_NONE`
202  * @param[in,out] y_ceed        Output `CeedVector`, maybe `CEED_VECTOR_NONE`
203  * @param[in,out] Y_loc         Output local `Vec`, maybe `NULL`
204  * @param[in,out] Y             Output global `Vec`, maybe `NULL`
205  * @param[in]     ctx           Context for the operator apply
206  * @param[in]     use_apply_add Whether to use `CeedOperatorApply` or `CeedOperatorApplyAdd`
207  */
208 PetscErrorCode ApplyCeedOperator_Core(Vec X, Vec X_loc, CeedVector x_ceed, CeedVector y_ceed, Vec Y_loc, Vec Y, OperatorApplyContext ctx,
209                                       bool use_apply_add) {
210   PetscMemType x_mem_type, y_mem_type;
211   Ceed         ceed = ctx->ceed;
212 
213   PetscFunctionBeginUser;
214   if (X) PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc));
215   if (X_loc) PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, x_ceed));
216 
217   if (Y_loc) PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, y_ceed));
218 
219   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, X, Y, 0, 0));
220   PetscCall(PetscLogGpuTimeBegin());
221   if (use_apply_add) PetscCallCeed(ceed, CeedOperatorApplyAdd(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE));
222   else PetscCallCeed(ceed, CeedOperatorApply(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE));
223   PetscCall(PetscLogGpuTimeEnd());
224   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, X, Y, 0, 0));
225 
226   if (X_loc) PetscCall(VecReadCeedToPetsc(ctx->x_ceed, x_mem_type, X_loc));
227 
228   if (Y_loc) PetscCall(VecCeedToPetsc(ctx->y_ceed, y_mem_type, Y_loc));
229   if (Y) PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y));
230   PetscFunctionReturn(PETSC_SUCCESS);
231 };
232 
233 PetscErrorCode ApplyCeedOperatorGlobalToGlobal(Vec X, Vec Y, OperatorApplyContext ctx) {
234   Vec X_loc = ctx->X_loc, Y_loc = ctx->Y_loc;
235 
236   PetscFunctionBeginUser;
237   PetscCall(VecZeroEntries(Y));
238 
239   // Get local vectors (if needed)
240   if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
241   if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
242 
243   PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false));
244 
245   // Restore local vector (if needed)
246   if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
247   if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
248   PetscFunctionReturn(PETSC_SUCCESS);
249 }
250 
251 PetscErrorCode ApplyCeedOperatorLocalToGlobal(Vec X_loc, Vec Y, OperatorApplyContext ctx) {
252   Vec Y_loc = ctx->Y_loc;
253 
254   PetscFunctionBeginUser;
255   PetscCall(VecZeroEntries(Y));
256 
257   // Get local vectors (if needed)
258   if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
259 
260   PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false));
261 
262   // Restore local vectors (if needed)
263   if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
264   PetscFunctionReturn(PETSC_SUCCESS);
265 }
266 
267 PetscErrorCode ApplyCeedOperatorGlobalToLocal(Vec X, Vec Y_loc, OperatorApplyContext ctx) {
268   Vec X_loc = ctx->X_loc;
269 
270   PetscFunctionBeginUser;
271   // Get local vectors (if needed)
272   if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
273 
274   PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false));
275 
276   // Restore local vector (if needed)
277   if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
278   PetscFunctionReturn(PETSC_SUCCESS);
279 }
280 
281 PetscErrorCode ApplyCeedOperatorLocalToLocal(Vec X_loc, Vec Y_loc, OperatorApplyContext ctx) {
282   PetscFunctionBeginUser;
283   PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false));
284   PetscFunctionReturn(PETSC_SUCCESS);
285 }
286 
287 PetscErrorCode ApplyAddCeedOperatorLocalToLocal(Vec X_loc, Vec Y_loc, OperatorApplyContext ctx) {
288   PetscFunctionBeginUser;
289   PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, true));
290   PetscFunctionReturn(PETSC_SUCCESS);
291 }
292 
293 /**
294  * @brief Return Mats for KSP solve
295  *
296  * Uses command-line flag with `ksp`'s prefix to determine if mat_ceed should be used directly or whether it should be assembled.
297  * If `Amat` is to be assembled, then `Pmat` is equal to `Amat`.
298  *
299  * If `Amat` uses `mat_ceed`, then `Pmat` is either assembled or uses `mat_ceed` based on the preconditioner choice in `ksp`.
300  *
301  * @param[in]  ksp      `KSP` object for used for solving
302  * @param[in]  mat_ceed `MATCEED` for the linear operator
303  * @param[in]  assemble Whether to assemble `Amat` and `Pmat` if they are not `mat_ceed`
304  * @param[out] Amat     `Mat` to be used for the solver `Amat`
305  * @param[out] Pmat     `Mat` to be used for the solver `Pmat`
306  */
307 PetscErrorCode CreateSolveOperatorsFromMatCeed(KSP ksp, Mat mat_ceed, PetscBool assemble, Mat *Amat, Mat *Pmat) {
308   PetscBool use_matceed_pmat, assemble_amat = PETSC_FALSE;
309   MatType   mat_ceed_inner_type;
310 
311   PetscFunctionBeginUser;
312   PetscCall(MatCeedGetInnerMatType(mat_ceed, &mat_ceed_inner_type));
313   {  // Determine if Amat should be MATCEED or assembled
314     const char *ksp_prefix = NULL;
315 
316     PetscCall(KSPGetOptionsPrefix(ksp, &ksp_prefix));
317     PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), ksp_prefix, "", NULL);
318     PetscCall(PetscOptionsBool("-matceed_assemble_amat", "Assemble the A matrix for KSP solve", NULL, assemble_amat, &assemble_amat, NULL));
319     PetscOptionsEnd();
320   }
321 
322   if (assemble_amat) {
323     PetscCall(MatConvert(mat_ceed, mat_ceed_inner_type, MAT_INITIAL_MATRIX, Amat));
324     if (assemble) PetscCall(MatCeedAssembleCOO(mat_ceed, *Amat));
325 
326     PetscCall(PetscObjectReference((PetscObject)*Amat));
327     *Pmat = *Amat;
328     PetscFunctionReturn(PETSC_SUCCESS);
329   } else {
330     PetscCall(PetscObjectReference((PetscObject)mat_ceed));
331     *Amat = mat_ceed;
332   }
333 
334   {  // Determine if Pmat should be MATCEED or assembled
335     PC     pc;
336     PCType pc_type;
337 
338     PetscCall(KSPGetPC(ksp, &pc));
339     PetscCall(PCGetType(pc, &pc_type));
340     PetscCall(PetscStrcmpAny(pc_type, &use_matceed_pmat, PCJACOBI, PCVPBJACOBI, PCPBJACOBI, ""));
341   }
342 
343   if (use_matceed_pmat) {
344     PetscCall(PetscObjectReference((PetscObject)mat_ceed));
345     *Pmat = mat_ceed;
346   } else {
347     PetscCall(MatConvert(mat_ceed, mat_ceed_inner_type, MAT_INITIAL_MATRIX, Pmat));
348     if (assemble) PetscCall(MatCeedAssembleCOO(mat_ceed, *Pmat));
349   }
350   PetscFunctionReturn(PETSC_SUCCESS);
351 }
352 
353 /**
354  * @brief Runs KSPSetFromOptions and sets Operators based on mat_ceed
355  *
356  * See CreateSolveOperatorsFromMatCeed for details on how the KSPSolve operators are set.
357  *
358  * @param[in] ksp      `KSP` of the solve
359  * @param[in] mat_ceed `MatCeed` linear operator to solve for
360  */
361 PetscErrorCode KSPSetFromOptions_WithMatCeed(KSP ksp, Mat mat_ceed) {
362   Mat Amat, Pmat;
363 
364   PetscFunctionBeginUser;
365   PetscCall(KSPSetFromOptions(ksp));
366   PetscCall(CreateSolveOperatorsFromMatCeed(ksp, mat_ceed, PETSC_TRUE, &Amat, &Pmat));
367   PetscCall(KSPSetOperators(ksp, Amat, Pmat));
368   PetscCall(MatDestroy(&Amat));
369   PetscCall(MatDestroy(&Pmat));
370   PetscFunctionReturn(PETSC_SUCCESS);
371 }
372 // -----------------------------------------------------------------------------
373