xref: /honee/src/petsc_ops.c (revision c1c64bfc7dc4f9d635fe2c12e841413ce2125953)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2f7325489SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3f7325489SJames Wright //
4f7325489SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5f7325489SJames Wright //
6f7325489SJames Wright // This file is part of CEED:  http://github.com/ceed
7f7325489SJames Wright 
8f7325489SJames Wright #include "../include/petsc_ops.h"
9f7325489SJames Wright 
10f7325489SJames Wright #include <ceed.h>
11f7325489SJames Wright #include <petscdm.h>
12f7325489SJames Wright 
13f7325489SJames Wright #include "../navierstokes.h"
14f7325489SJames Wright 
1529afcdd1SJames Wright // @brief Get information about a DM's local vector
1629afcdd1SJames Wright PetscErrorCode DMGetLocalVectorInfo(DM dm, PetscInt *local_size, PetscInt *global_size, VecType *vec_type) {
1729afcdd1SJames Wright   Vec V_loc;
1829afcdd1SJames Wright 
1929afcdd1SJames Wright   PetscFunctionBeginUser;
2029afcdd1SJames Wright   PetscCall(DMGetLocalVector(dm, &V_loc));
2129afcdd1SJames Wright   if (local_size) PetscCall(VecGetLocalSize(V_loc, local_size));
2229afcdd1SJames Wright   if (global_size) PetscCall(VecGetSize(V_loc, global_size));
2329afcdd1SJames Wright   if (vec_type) PetscCall(VecGetType(V_loc, vec_type));
2429afcdd1SJames Wright   PetscCall(DMRestoreLocalVector(dm, &V_loc));
25d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2629afcdd1SJames Wright }
2729afcdd1SJames Wright 
2829afcdd1SJames Wright // @brief Get information about a DM's global vector
2929afcdd1SJames Wright PetscErrorCode DMGetGlobalVectorInfo(DM dm, PetscInt *local_size, PetscInt *global_size, VecType *vec_type) {
3029afcdd1SJames Wright   Vec V;
3129afcdd1SJames Wright 
3229afcdd1SJames Wright   PetscFunctionBeginUser;
3329afcdd1SJames Wright   PetscCall(DMGetGlobalVector(dm, &V));
3429afcdd1SJames Wright   if (local_size) PetscCall(VecGetLocalSize(V, local_size));
3529afcdd1SJames Wright   if (global_size) PetscCall(VecGetSize(V, global_size));
3629afcdd1SJames Wright   if (vec_type) PetscCall(VecGetType(V, vec_type));
3729afcdd1SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &V));
38d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3929afcdd1SJames Wright }
4029afcdd1SJames Wright 
41fbcfc4fcSJames Wright /**
42fbcfc4fcSJames Wright  * @brief Create OperatorApplyContext struct for applying FEM operator in a PETSc context
43fbcfc4fcSJames Wright  *
44fbcfc4fcSJames Wright  * All passed in objects are reference copied and may be destroyed if desired (with the exception of `CEED_VECTOR_NONE`).
45fbcfc4fcSJames Wright  * Resulting context should be destroyed with `OperatorApplyContextDestroy()`.
46fbcfc4fcSJames Wright  *
47fbcfc4fcSJames Wright  * @param[in]  dm_x     `DM` associated with the operator active input. May be `NULL`
48fbcfc4fcSJames Wright  * @param[in]  dm_y     `DM` associated with the operator active output. May be `NULL`
49fbcfc4fcSJames Wright  * @param[in]  ceed     `Ceed` object
50fbcfc4fcSJames Wright  * @param[in]  op_apply `CeedOperator` representing the local action of the FEM operator
51fbcfc4fcSJames Wright  * @param[in]  x_ceed   `CeedVector` for operator active input. May be `CEED_VECTOR_NONE` or `NULL`. If `NULL`, `CeedVector` will be automatically
52fbcfc4fcSJames Wright  *                      generated.
53fbcfc4fcSJames Wright  * @param[in]  y_ceed   `CeedVector` for operator active output. May be `CEED_VECTOR_NONE` or `NULL`. If `NULL`, `CeedVector` will be automatically
54fbcfc4fcSJames Wright  *                      generated.
55fbcfc4fcSJames Wright  * @param[in]  X_loc    Local `Vec` for operator active input. If `NULL`, vector will be obtained if needed at ApplyCeedOperator time.
56fbcfc4fcSJames Wright  * @param[in]  Y_loc    Local `Vec` for operator active output. If `NULL`, vector will be obtained if needed at ApplyCeedOperator time.
57fbcfc4fcSJames Wright  * @param[out] ctx      Struct containing all data necessary for applying the operator
58fbcfc4fcSJames Wright  */
59f7325489SJames Wright PetscErrorCode OperatorApplyContextCreate(DM dm_x, DM dm_y, Ceed ceed, CeedOperator op_apply, CeedVector x_ceed, CeedVector y_ceed, Vec X_loc,
60fbcfc4fcSJames Wright                                           Vec Y_loc, OperatorApplyContext *ctx) {
61f7325489SJames Wright   CeedSize x_size, y_size;
62fbcfc4fcSJames Wright 
63fbcfc4fcSJames Wright   PetscFunctionBeginUser;
64b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorGetActiveVectorLengths(op_apply, &x_size, &y_size));
65fbcfc4fcSJames Wright   {  // Verify sizes
6682560a70SJames Wright     PetscInt X_size, Y_size, dm_X_size, dm_Y_size;
6782560a70SJames Wright     CeedSize x_ceed_size, y_ceed_size;
6882560a70SJames Wright     if (dm_x) PetscCall(DMGetLocalVectorInfo(dm_x, &dm_X_size, NULL, NULL));
6982560a70SJames Wright     if (dm_y) PetscCall(DMGetLocalVectorInfo(dm_y, &dm_Y_size, NULL, NULL));
70f7325489SJames Wright     if (X_loc) {
71f7325489SJames Wright       PetscCall(VecGetLocalSize(X_loc, &X_size));
72f7325489SJames Wright       PetscCheck(X_size == x_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
73f7325489SJames Wright                  "X_loc (%" PetscInt_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", X_size, x_size);
7482560a70SJames Wright       if (dm_x)
7582560a70SJames Wright         PetscCheck(X_size == dm_X_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
7682560a70SJames Wright                    "X_loc size (%" PetscInt_FMT ") does not match dm_x local vector size (%" PetscInt_FMT ")", X_size, dm_X_size);
77f7325489SJames Wright     }
78f7325489SJames Wright     if (Y_loc) {
79f7325489SJames Wright       PetscCall(VecGetLocalSize(Y_loc, &Y_size));
80f7325489SJames Wright       PetscCheck(Y_size == y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
81f7325489SJames Wright                  "Y_loc (%" PetscInt_FMT ") not correct size for CeedOperator active output size (%" CeedSize_FMT ")", Y_size, y_size);
8282560a70SJames Wright       if (dm_y)
8382560a70SJames Wright         PetscCheck(Y_size == dm_Y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
8482560a70SJames Wright                    "Y_loc size (%" PetscInt_FMT ") does not match dm_y local vector size (%" PetscInt_FMT ")", Y_size, dm_Y_size);
8582560a70SJames Wright     }
8682560a70SJames Wright     if (x_ceed && x_ceed != CEED_VECTOR_NONE) {
87b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedVectorGetLength(x_ceed, &x_ceed_size));
8882560a70SJames Wright       PetscCheck(x_size >= 0 ? x_ceed_size == x_size : true, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
8982560a70SJames Wright                  "x_ceed (%" CeedSize_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", x_ceed_size, x_size);
9082560a70SJames Wright       if (dm_x)
9182560a70SJames Wright         PetscCheck(x_ceed_size == dm_X_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
9282560a70SJames Wright                    "x_ceed size (%" CeedSize_FMT ") does not match dm_x local vector size (%" PetscInt_FMT ")", x_ceed_size, dm_X_size);
9382560a70SJames Wright     }
9482560a70SJames Wright     if (y_ceed && y_ceed != CEED_VECTOR_NONE) {
95b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedVectorGetLength(y_ceed, &y_ceed_size));
9682560a70SJames Wright       PetscCheck(y_ceed_size == y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
9782560a70SJames Wright                  "y_ceed (%" CeedSize_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", y_ceed_size, y_size);
9882560a70SJames Wright       if (dm_y)
9982560a70SJames Wright         PetscCheck(y_ceed_size == dm_Y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ,
10082560a70SJames Wright                    "y_ceed size (%" CeedSize_FMT ") does not match dm_y local vector size (%" PetscInt_FMT ")", y_ceed_size, dm_Y_size);
101f7325489SJames Wright     }
102f7325489SJames Wright   }
103f7325489SJames Wright 
104fbcfc4fcSJames Wright   PetscCall(PetscNew(ctx));
105f7325489SJames Wright 
106f7325489SJames Wright   // Copy PETSc Objects
107f7325489SJames Wright   if (dm_x) PetscCall(PetscObjectReference((PetscObject)dm_x));
108fbcfc4fcSJames Wright   (*ctx)->dm_x = dm_x;
109f7325489SJames Wright   if (dm_y) PetscCall(PetscObjectReference((PetscObject)dm_y));
110fbcfc4fcSJames Wright   (*ctx)->dm_y = dm_y;
111f7325489SJames Wright 
112f7325489SJames Wright   if (X_loc) PetscCall(PetscObjectReference((PetscObject)X_loc));
113fbcfc4fcSJames Wright   (*ctx)->X_loc = X_loc;
114f7325489SJames Wright   if (Y_loc) PetscCall(PetscObjectReference((PetscObject)Y_loc));
115fbcfc4fcSJames Wright   (*ctx)->Y_loc = Y_loc;
116f7325489SJames Wright 
117f7325489SJames Wright   // Copy libCEED objects
118b4c37c5cSJames Wright   if (x_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(x_ceed, &(*ctx)->x_ceed));
119b4c37c5cSJames Wright   else PetscCallCeed(ceed, CeedVectorCreate(ceed, x_size, &(*ctx)->x_ceed));
120fbcfc4fcSJames Wright 
121b4c37c5cSJames Wright   if (y_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(y_ceed, &(*ctx)->y_ceed));
122b4c37c5cSJames Wright   else PetscCallCeed(ceed, CeedVectorCreate(ceed, y_size, &(*ctx)->y_ceed));
123fbcfc4fcSJames Wright 
124b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorReferenceCopy(op_apply, &(*ctx)->op));
125b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedReferenceCopy(ceed, &(*ctx)->ceed));
126d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
127f7325489SJames Wright }
128f7325489SJames Wright 
129fbcfc4fcSJames Wright /**
130fbcfc4fcSJames Wright  * @brief Destroy OperatorApplyContext struct
131fbcfc4fcSJames Wright  *
132fbcfc4fcSJames Wright  * @param[in,out] ctx Context to destroy
133fbcfc4fcSJames Wright  */
134fbcfc4fcSJames Wright PetscErrorCode OperatorApplyContextDestroy(OperatorApplyContext ctx) {
135f7325489SJames Wright   PetscFunctionBeginUser;
136d949ddfcSJames Wright   if (!ctx) PetscFunctionReturn(PETSC_SUCCESS);
137b4c37c5cSJames Wright   Ceed ceed = ctx->ceed;
138f7325489SJames Wright 
139f7325489SJames Wright   // Destroy PETSc Objects
140fbcfc4fcSJames Wright   PetscCall(DMDestroy(&ctx->dm_x));
141fbcfc4fcSJames Wright   PetscCall(DMDestroy(&ctx->dm_y));
142fbcfc4fcSJames Wright   PetscCall(VecDestroy(&ctx->X_loc));
143fbcfc4fcSJames Wright   PetscCall(VecDestroy(&ctx->Y_loc));
144f7325489SJames Wright 
145f7325489SJames Wright   // Destroy libCEED Objects
146b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&ctx->x_ceed));
147b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&ctx->y_ceed));
148b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&ctx->op));
149b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedDestroy(&ctx->ceed));
150f7325489SJames Wright 
151fbcfc4fcSJames Wright   PetscCall(PetscFree(ctx));
152d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153f7325489SJames Wright }
154f7325489SJames Wright 
1556cbc221eSJames Wright //@brief Return VecType for the given DM
1566cbc221eSJames Wright VecType DMReturnVecType(DM dm) {
1576cbc221eSJames Wright   VecType vec_type;
1586cbc221eSJames Wright   DMGetVecType(dm, &vec_type);
1596cbc221eSJames Wright   return vec_type;
1606cbc221eSJames Wright }
1616cbc221eSJames Wright 
1626cbc221eSJames Wright /**
1636cbc221eSJames Wright  * @brief Create local PETSc Vecs for CeedOperator's active input/outputs
1646cbc221eSJames Wright  *
1656cbc221eSJames Wright  * This is primarily used for when the active input/ouput vector does not correspond to a `DM` object, and thus `DMCreateLocalVector` or
1666cbc221eSJames Wright  * `DMGetLocalVector` are not applicable.
1676cbc221eSJames Wright  * For example, if statitics are being store at quadrature points, a `DM`-created `Vec` will not have the
1686cbc221eSJames Wright  * correct size.
1696cbc221eSJames Wright  *
1706cbc221eSJames Wright  * @param[in]  op       Operator to make the Vecs for
17195327578SJames Wright  * @param[in]  vec_type `VecType` for the new Vecs
17295327578SJames Wright  * @param[in]  comm     `MPI_Comm` for the new Vecs
1736cbc221eSJames Wright  * @param[out] input    Vec for CeedOperator active input
1746cbc221eSJames Wright  * @param[out] output   Vec for CeedOperator active output
1756cbc221eSJames Wright  */
1766cbc221eSJames Wright PetscErrorCode CeedOperatorCreateLocalVecs(CeedOperator op, VecType vec_type, MPI_Comm comm, Vec *input, Vec *output) {
1776cbc221eSJames Wright   CeedSize input_size, output_size;
178b4c37c5cSJames Wright   Ceed     ceed;
1791fe3d3f0SJames Wright   int      comm_size;
1806cbc221eSJames Wright 
1816cbc221eSJames Wright   PetscFunctionBeginUser;
182b4c37c5cSJames Wright   PetscCall(CeedOperatorGetCeed(op, &ceed));
1831fe3d3f0SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
1841fe3d3f0SJames Wright   PetscCheck(comm_size == 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "MPI_Comm must be of size 1, recieved comm of size %d", comm_size);
185b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size));
1866cbc221eSJames Wright   if (input) {
1876cbc221eSJames Wright     PetscCall(VecCreate(comm, input));
1886cbc221eSJames Wright     PetscCall(VecSetType(*input, vec_type));
1896cbc221eSJames Wright     PetscCall(VecSetSizes(*input, input_size, input_size));
1906cbc221eSJames Wright   }
1916cbc221eSJames Wright   if (output) {
1926cbc221eSJames Wright     PetscCall(VecCreate(comm, output));
1936cbc221eSJames Wright     PetscCall(VecSetType(*output, vec_type));
1946cbc221eSJames Wright     PetscCall(VecSetSizes(*output, output_size, output_size));
1956cbc221eSJames Wright   }
196d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1976cbc221eSJames Wright }
1986cbc221eSJames Wright 
199f7325489SJames Wright /**
200f7325489SJames Wright  * @brief Apply FEM Operator defined by `OperatorApplyContext` to various input and output vectors
201f7325489SJames Wright  *
2023fc17c4bSJames Wright  * @param[in]     X             Input global `Vec`, maybe `NULL`
2033fc17c4bSJames Wright  * @param[in]     X_loc         Input local `Vec`, maybe `NULL`
2043fc17c4bSJames Wright  * @param[in]     x_ceed        Input `CeedVector`, maybe `CEED_VECTOR_NONE`
2053fc17c4bSJames Wright  * @param[in,out] y_ceed        Output `CeedVector`, maybe `CEED_VECTOR_NONE`
2063fc17c4bSJames Wright  * @param[in,out] Y_loc         Output local `Vec`, maybe `NULL`
2073fc17c4bSJames Wright  * @param[in,out] Y             Output global `Vec`, maybe `NULL`
2083fc17c4bSJames Wright  * @param[in]     ctx           Context for the operator apply
2093fc17c4bSJames Wright  * @param[in]     use_apply_add Whether to use `CeedOperatorApply` or `CeedOperatorApplyAdd`
210f7325489SJames Wright  */
211f7325489SJames Wright PetscErrorCode ApplyCeedOperator_Core(Vec X, Vec X_loc, CeedVector x_ceed, CeedVector y_ceed, Vec Y_loc, Vec Y, OperatorApplyContext ctx,
212f7325489SJames Wright                                       bool use_apply_add) {
213f7325489SJames Wright   PetscMemType x_mem_type, y_mem_type;
214b4c37c5cSJames Wright   Ceed         ceed = ctx->ceed;
215f7325489SJames Wright 
216f7325489SJames Wright   PetscFunctionBeginUser;
217f7325489SJames Wright   if (X) PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc));
218a7dac1d5SJames Wright   if (X_loc) PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, x_ceed));
219f7325489SJames Wright 
220a7dac1d5SJames Wright   if (Y_loc) PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, y_ceed));
221f7325489SJames Wright 
2227eedc94cSJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_CeedOperatorApply, X, Y, 0, 0));
2237eedc94cSJames Wright   PetscCall(PetscLogGpuTimeBegin());
224b4c37c5cSJames Wright   if (use_apply_add) PetscCallCeed(ceed, CeedOperatorApplyAdd(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE));
225b4c37c5cSJames Wright   else PetscCallCeed(ceed, CeedOperatorApply(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE));
2267eedc94cSJames Wright   PetscCall(PetscLogGpuTimeEnd());
2277eedc94cSJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_CeedOperatorApply, X, Y, 0, 0));
228f7325489SJames Wright 
229a7dac1d5SJames Wright   if (X_loc) PetscCall(VecReadCeedToPetsc(ctx->x_ceed, x_mem_type, X_loc));
230f7325489SJames Wright 
231a7dac1d5SJames Wright   if (Y_loc) PetscCall(VecCeedToPetsc(ctx->y_ceed, y_mem_type, Y_loc));
232f7325489SJames Wright   if (Y) PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y));
233d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
234f7325489SJames Wright };
235f7325489SJames Wright 
236f7325489SJames Wright PetscErrorCode ApplyCeedOperatorGlobalToGlobal(Vec X, Vec Y, OperatorApplyContext ctx) {
237f7325489SJames Wright   Vec X_loc = ctx->X_loc, Y_loc = ctx->Y_loc;
238f7325489SJames Wright 
239f7325489SJames Wright   PetscFunctionBeginUser;
240f7325489SJames Wright   PetscCall(VecZeroEntries(Y));
241f7325489SJames Wright 
242f7325489SJames Wright   // Get local vectors (if needed)
243f7325489SJames Wright   if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
244f7325489SJames Wright   if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
245f7325489SJames Wright 
246f7325489SJames Wright   PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false));
247f7325489SJames Wright 
248f7325489SJames Wright   // Restore local vector (if needed)
249f7325489SJames Wright   if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
250f7325489SJames Wright   if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
251d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
252f7325489SJames Wright }
253f7325489SJames Wright 
254f7325489SJames Wright PetscErrorCode ApplyCeedOperatorLocalToGlobal(Vec X_loc, Vec Y, OperatorApplyContext ctx) {
255f7325489SJames Wright   Vec Y_loc = ctx->Y_loc;
256f7325489SJames Wright 
257f7325489SJames Wright   PetscFunctionBeginUser;
258f7325489SJames Wright   PetscCall(VecZeroEntries(Y));
259f7325489SJames Wright 
260f7325489SJames Wright   // Get local vectors (if needed)
261f7325489SJames Wright   if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
262f7325489SJames Wright 
263f7325489SJames Wright   PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false));
264f7325489SJames Wright 
265f7325489SJames Wright   // Restore local vectors (if needed)
266f7325489SJames Wright   if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
267d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
268f7325489SJames Wright }
269f7325489SJames Wright 
270f7325489SJames Wright PetscErrorCode ApplyCeedOperatorGlobalToLocal(Vec X, Vec Y_loc, OperatorApplyContext ctx) {
271f7325489SJames Wright   Vec X_loc = ctx->X_loc;
272f7325489SJames Wright 
273f7325489SJames Wright   PetscFunctionBeginUser;
274f7325489SJames Wright   // Get local vectors (if needed)
275f7325489SJames Wright   if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
276f7325489SJames Wright 
277f7325489SJames Wright   PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false));
278f7325489SJames Wright 
279f7325489SJames Wright   // Restore local vector (if needed)
280f7325489SJames Wright   if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
281d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
282f7325489SJames Wright }
283f7325489SJames Wright 
2845ce09d65SJames Wright PetscErrorCode ApplyCeedOperatorLocalToLocal(Vec X_loc, Vec Y_loc, OperatorApplyContext ctx) {
2855ce09d65SJames Wright   PetscFunctionBeginUser;
2865ce09d65SJames Wright   PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false));
287d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2885ce09d65SJames Wright }
2895ce09d65SJames Wright 
290f7325489SJames Wright PetscErrorCode ApplyAddCeedOperatorLocalToLocal(Vec X_loc, Vec Y_loc, OperatorApplyContext ctx) {
291f7325489SJames Wright   PetscFunctionBeginUser;
292f7325489SJames Wright   PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, true));
293d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
294f7325489SJames Wright }
295f7325489SJames Wright 
2963170c09fSJames Wright /**
2973170c09fSJames Wright  * @brief Return Mats for KSP solve
2983170c09fSJames Wright  *
2993170c09fSJames Wright  * Uses command-line flag with `ksp`'s prefix to determine if mat_ceed should be used directly or whether it should be assembled.
3003170c09fSJames Wright  * If `Amat` is to be assembled, then `Pmat` is equal to `Amat`.
3013170c09fSJames Wright  *
3023170c09fSJames Wright  * If `Amat` uses `mat_ceed`, then `Pmat` is either assembled or uses `mat_ceed` based on the preconditioner choice in `ksp`.
3033170c09fSJames Wright  *
3043170c09fSJames Wright  * @param[in]  ksp      `KSP` object for used for solving
3053170c09fSJames Wright  * @param[in]  mat_ceed `MATCEED` for the linear operator
3063170c09fSJames Wright  * @param[in]  assemble Whether to assemble `Amat` and `Pmat` if they are not `mat_ceed`
3073170c09fSJames Wright  * @param[out] Amat     `Mat` to be used for the solver `Amat`
3083170c09fSJames Wright  * @param[out] Pmat     `Mat` to be used for the solver `Pmat`
3093170c09fSJames Wright  */
3103170c09fSJames Wright PetscErrorCode CreateSolveOperatorsFromMatCeed(KSP ksp, Mat mat_ceed, PetscBool assemble, Mat *Amat, Mat *Pmat) {
3113170c09fSJames Wright   PetscBool use_matceed_pmat, assemble_amat = PETSC_FALSE;
3123170c09fSJames Wright 
3133170c09fSJames Wright   PetscFunctionBeginUser;
3143170c09fSJames Wright   {  // Determine if Amat should be MATCEED or assembled
3153170c09fSJames Wright     const char *ksp_prefix = NULL;
3163170c09fSJames Wright 
3173170c09fSJames Wright     PetscCall(KSPGetOptionsPrefix(ksp, &ksp_prefix));
3183170c09fSJames Wright     PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), ksp_prefix, "", NULL);
3193170c09fSJames Wright     PetscCall(PetscOptionsBool("-matceed_assemble_amat", "Assemble the A matrix for KSP solve", NULL, assemble_amat, &assemble_amat, NULL));
3203170c09fSJames Wright     PetscOptionsEnd();
3213170c09fSJames Wright   }
3223170c09fSJames Wright 
3233170c09fSJames Wright   if (assemble_amat) {
324*c1c64bfcSJames Wright     PetscCall(MatCeedCreateMatCOO(mat_ceed, Amat));
3253170c09fSJames Wright     if (assemble) PetscCall(MatCeedAssembleCOO(mat_ceed, *Amat));
3263170c09fSJames Wright 
3273170c09fSJames Wright     PetscCall(PetscObjectReference((PetscObject)*Amat));
3283170c09fSJames Wright     *Pmat = *Amat;
3293170c09fSJames Wright     PetscFunctionReturn(PETSC_SUCCESS);
3303170c09fSJames Wright   } else {
3313170c09fSJames Wright     PetscCall(PetscObjectReference((PetscObject)mat_ceed));
3323170c09fSJames Wright     *Amat = mat_ceed;
3333170c09fSJames Wright   }
3343170c09fSJames Wright 
3353170c09fSJames Wright   {  // Determine if Pmat should be MATCEED or assembled
3363170c09fSJames Wright     PC     pc;
3373170c09fSJames Wright     PCType pc_type;
3383170c09fSJames Wright 
3393170c09fSJames Wright     PetscCall(KSPGetPC(ksp, &pc));
3403170c09fSJames Wright     PetscCall(PCGetType(pc, &pc_type));
3413170c09fSJames Wright     PetscCall(PetscStrcmpAny(pc_type, &use_matceed_pmat, PCJACOBI, PCVPBJACOBI, PCPBJACOBI, ""));
3423170c09fSJames Wright   }
3433170c09fSJames Wright 
3443170c09fSJames Wright   if (use_matceed_pmat) {
3453170c09fSJames Wright     PetscCall(PetscObjectReference((PetscObject)mat_ceed));
3463170c09fSJames Wright     *Pmat = mat_ceed;
3473170c09fSJames Wright   } else {
348*c1c64bfcSJames Wright     PetscCall(MatCeedCreateMatCOO(mat_ceed, Pmat));
3493170c09fSJames Wright     if (assemble) PetscCall(MatCeedAssembleCOO(mat_ceed, *Pmat));
3503170c09fSJames Wright   }
3513170c09fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3523170c09fSJames Wright }
3533170c09fSJames Wright 
3543170c09fSJames Wright /**
3553170c09fSJames Wright  * @brief Runs KSPSetFromOptions and sets Operators based on mat_ceed
3563170c09fSJames Wright  *
3573170c09fSJames Wright  * See CreateSolveOperatorsFromMatCeed for details on how the KSPSolve operators are set.
3583170c09fSJames Wright  *
3593170c09fSJames Wright  * @param[in] ksp      `KSP` of the solve
3603170c09fSJames Wright  * @param[in] mat_ceed `MatCeed` linear operator to solve for
3613170c09fSJames Wright  */
3623170c09fSJames Wright PetscErrorCode KSPSetFromOptions_WithMatCeed(KSP ksp, Mat mat_ceed) {
3633170c09fSJames Wright   Mat Amat, Pmat;
3643170c09fSJames Wright 
3653170c09fSJames Wright   PetscFunctionBeginUser;
3663170c09fSJames Wright   PetscCall(KSPSetFromOptions(ksp));
3673170c09fSJames Wright   PetscCall(CreateSolveOperatorsFromMatCeed(ksp, mat_ceed, PETSC_TRUE, &Amat, &Pmat));
3683170c09fSJames Wright   PetscCall(KSPSetOperators(ksp, Amat, Pmat));
3693170c09fSJames Wright   PetscCall(MatDestroy(&Amat));
3703170c09fSJames Wright   PetscCall(MatDestroy(&Pmat));
3713170c09fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3723170c09fSJames Wright }
373f7325489SJames Wright // -----------------------------------------------------------------------------
374