1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3f7325489SJames Wright 4149fb536SJames Wright #include <petsc_ops.h> 5f7325489SJames Wright 6f7325489SJames Wright #include <ceed.h> 7f7325489SJames Wright #include <petscdm.h> 8f7325489SJames Wright 9149fb536SJames Wright #include <navierstokes.h> 10f7325489SJames Wright 1129afcdd1SJames Wright // @brief Get information about a DM's local vector 1229afcdd1SJames Wright PetscErrorCode DMGetLocalVectorInfo(DM dm, PetscInt *local_size, PetscInt *global_size, VecType *vec_type) { 1329afcdd1SJames Wright Vec V_loc; 1429afcdd1SJames Wright 1529afcdd1SJames Wright PetscFunctionBeginUser; 1629afcdd1SJames Wright PetscCall(DMGetLocalVector(dm, &V_loc)); 1729afcdd1SJames Wright if (local_size) PetscCall(VecGetLocalSize(V_loc, local_size)); 1829afcdd1SJames Wright if (global_size) PetscCall(VecGetSize(V_loc, global_size)); 1929afcdd1SJames Wright if (vec_type) PetscCall(VecGetType(V_loc, vec_type)); 2029afcdd1SJames Wright PetscCall(DMRestoreLocalVector(dm, &V_loc)); 21d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2229afcdd1SJames Wright } 2329afcdd1SJames Wright 2429afcdd1SJames Wright // @brief Get information about a DM's global vector 2529afcdd1SJames Wright PetscErrorCode DMGetGlobalVectorInfo(DM dm, PetscInt *local_size, PetscInt *global_size, VecType *vec_type) { 2629afcdd1SJames Wright Vec V; 2729afcdd1SJames Wright 2829afcdd1SJames Wright PetscFunctionBeginUser; 2929afcdd1SJames Wright PetscCall(DMGetGlobalVector(dm, &V)); 3029afcdd1SJames Wright if (local_size) PetscCall(VecGetLocalSize(V, local_size)); 3129afcdd1SJames Wright if (global_size) PetscCall(VecGetSize(V, global_size)); 3229afcdd1SJames Wright if (vec_type) PetscCall(VecGetType(V, vec_type)); 3329afcdd1SJames Wright PetscCall(DMRestoreGlobalVector(dm, &V)); 34d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3529afcdd1SJames Wright } 3629afcdd1SJames Wright 37fbcfc4fcSJames Wright /** 38fbcfc4fcSJames Wright * @brief Create OperatorApplyContext struct for applying FEM operator in a PETSc context 39fbcfc4fcSJames Wright * 40fbcfc4fcSJames Wright * All passed in objects are reference copied and may be destroyed if desired (with the exception of `CEED_VECTOR_NONE`). 41fbcfc4fcSJames Wright * Resulting context should be destroyed with `OperatorApplyContextDestroy()`. 42fbcfc4fcSJames Wright * 43fbcfc4fcSJames Wright * @param[in] dm_x `DM` associated with the operator active input. May be `NULL` 44fbcfc4fcSJames Wright * @param[in] dm_y `DM` associated with the operator active output. May be `NULL` 45fbcfc4fcSJames Wright * @param[in] ceed `Ceed` object 46fbcfc4fcSJames Wright * @param[in] op_apply `CeedOperator` representing the local action of the FEM operator 47fbcfc4fcSJames Wright * @param[in] x_ceed `CeedVector` for operator active input. May be `CEED_VECTOR_NONE` or `NULL`. If `NULL`, `CeedVector` will be automatically 48fbcfc4fcSJames Wright * generated. 49fbcfc4fcSJames Wright * @param[in] y_ceed `CeedVector` for operator active output. May be `CEED_VECTOR_NONE` or `NULL`. If `NULL`, `CeedVector` will be automatically 50fbcfc4fcSJames Wright * generated. 51fbcfc4fcSJames Wright * @param[in] X_loc Local `Vec` for operator active input. If `NULL`, vector will be obtained if needed at ApplyCeedOperator time. 52fbcfc4fcSJames Wright * @param[in] Y_loc Local `Vec` for operator active output. If `NULL`, vector will be obtained if needed at ApplyCeedOperator time. 53fbcfc4fcSJames Wright * @param[out] ctx Struct containing all data necessary for applying the operator 54fbcfc4fcSJames Wright */ 55f7325489SJames Wright PetscErrorCode OperatorApplyContextCreate(DM dm_x, DM dm_y, Ceed ceed, CeedOperator op_apply, CeedVector x_ceed, CeedVector y_ceed, Vec X_loc, 56fbcfc4fcSJames Wright Vec Y_loc, OperatorApplyContext *ctx) { 57f7325489SJames Wright CeedSize x_size, y_size; 58fbcfc4fcSJames Wright 59fbcfc4fcSJames Wright PetscFunctionBeginUser; 60b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetActiveVectorLengths(op_apply, &x_size, &y_size)); 61fbcfc4fcSJames Wright { // Verify sizes 6282560a70SJames Wright PetscInt X_size, Y_size, dm_X_size, dm_Y_size; 6382560a70SJames Wright CeedSize x_ceed_size, y_ceed_size; 6482560a70SJames Wright if (dm_x) PetscCall(DMGetLocalVectorInfo(dm_x, &dm_X_size, NULL, NULL)); 6582560a70SJames Wright if (dm_y) PetscCall(DMGetLocalVectorInfo(dm_y, &dm_Y_size, NULL, NULL)); 66f7325489SJames Wright if (X_loc) { 67f7325489SJames Wright PetscCall(VecGetLocalSize(X_loc, &X_size)); 68f7325489SJames Wright PetscCheck(X_size == x_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 69f7325489SJames Wright "X_loc (%" PetscInt_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", X_size, x_size); 7082560a70SJames Wright if (dm_x) 7182560a70SJames Wright PetscCheck(X_size == dm_X_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 7282560a70SJames Wright "X_loc size (%" PetscInt_FMT ") does not match dm_x local vector size (%" PetscInt_FMT ")", X_size, dm_X_size); 73f7325489SJames Wright } 74f7325489SJames Wright if (Y_loc) { 75f7325489SJames Wright PetscCall(VecGetLocalSize(Y_loc, &Y_size)); 76f7325489SJames Wright PetscCheck(Y_size == y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 77f7325489SJames Wright "Y_loc (%" PetscInt_FMT ") not correct size for CeedOperator active output size (%" CeedSize_FMT ")", Y_size, y_size); 7882560a70SJames Wright if (dm_y) 7982560a70SJames Wright PetscCheck(Y_size == dm_Y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 8082560a70SJames Wright "Y_loc size (%" PetscInt_FMT ") does not match dm_y local vector size (%" PetscInt_FMT ")", Y_size, dm_Y_size); 8182560a70SJames Wright } 8282560a70SJames Wright if (x_ceed && x_ceed != CEED_VECTOR_NONE) { 83b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorGetLength(x_ceed, &x_ceed_size)); 8482560a70SJames Wright PetscCheck(x_size >= 0 ? x_ceed_size == x_size : true, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 8582560a70SJames Wright "x_ceed (%" CeedSize_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", x_ceed_size, x_size); 8682560a70SJames Wright if (dm_x) 8782560a70SJames Wright PetscCheck(x_ceed_size == dm_X_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 8882560a70SJames Wright "x_ceed size (%" CeedSize_FMT ") does not match dm_x local vector size (%" PetscInt_FMT ")", x_ceed_size, dm_X_size); 8982560a70SJames Wright } 9082560a70SJames Wright if (y_ceed && y_ceed != CEED_VECTOR_NONE) { 91b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorGetLength(y_ceed, &y_ceed_size)); 9282560a70SJames Wright PetscCheck(y_ceed_size == y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 9382560a70SJames Wright "y_ceed (%" CeedSize_FMT ") not correct size for CeedOperator active input size (%" CeedSize_FMT ")", y_ceed_size, y_size); 9482560a70SJames Wright if (dm_y) 9582560a70SJames Wright PetscCheck(y_ceed_size == dm_Y_size, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, 9682560a70SJames Wright "y_ceed size (%" CeedSize_FMT ") does not match dm_y local vector size (%" PetscInt_FMT ")", y_ceed_size, dm_Y_size); 97f7325489SJames Wright } 98f7325489SJames Wright } 99f7325489SJames Wright 100fbcfc4fcSJames Wright PetscCall(PetscNew(ctx)); 101f7325489SJames Wright 102f7325489SJames Wright // Copy PETSc Objects 103f7325489SJames Wright if (dm_x) PetscCall(PetscObjectReference((PetscObject)dm_x)); 104fbcfc4fcSJames Wright (*ctx)->dm_x = dm_x; 105f7325489SJames Wright if (dm_y) PetscCall(PetscObjectReference((PetscObject)dm_y)); 106fbcfc4fcSJames Wright (*ctx)->dm_y = dm_y; 107f7325489SJames Wright 108f7325489SJames Wright if (X_loc) PetscCall(PetscObjectReference((PetscObject)X_loc)); 109fbcfc4fcSJames Wright (*ctx)->X_loc = X_loc; 110f7325489SJames Wright if (Y_loc) PetscCall(PetscObjectReference((PetscObject)Y_loc)); 111fbcfc4fcSJames Wright (*ctx)->Y_loc = Y_loc; 112f7325489SJames Wright 113f7325489SJames Wright // Copy libCEED objects 114b4c37c5cSJames Wright if (x_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(x_ceed, &(*ctx)->x_ceed)); 115b4c37c5cSJames Wright else PetscCallCeed(ceed, CeedVectorCreate(ceed, x_size, &(*ctx)->x_ceed)); 116fbcfc4fcSJames Wright 117b4c37c5cSJames Wright if (y_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(y_ceed, &(*ctx)->y_ceed)); 118b4c37c5cSJames Wright else PetscCallCeed(ceed, CeedVectorCreate(ceed, y_size, &(*ctx)->y_ceed)); 119fbcfc4fcSJames Wright 120b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorReferenceCopy(op_apply, &(*ctx)->op)); 121b4c37c5cSJames Wright PetscCallCeed(ceed, CeedReferenceCopy(ceed, &(*ctx)->ceed)); 122d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 123f7325489SJames Wright } 124f7325489SJames Wright 125fbcfc4fcSJames Wright /** 126fbcfc4fcSJames Wright * @brief Destroy OperatorApplyContext struct 127fbcfc4fcSJames Wright * 128fbcfc4fcSJames Wright * @param[in,out] ctx Context to destroy 129fbcfc4fcSJames Wright */ 130fbcfc4fcSJames Wright PetscErrorCode OperatorApplyContextDestroy(OperatorApplyContext ctx) { 131f7325489SJames Wright PetscFunctionBeginUser; 132d949ddfcSJames Wright if (!ctx) PetscFunctionReturn(PETSC_SUCCESS); 133b4c37c5cSJames Wright Ceed ceed = ctx->ceed; 134f7325489SJames Wright 135f7325489SJames Wright // Destroy PETSc Objects 136fbcfc4fcSJames Wright PetscCall(DMDestroy(&ctx->dm_x)); 137fbcfc4fcSJames Wright PetscCall(DMDestroy(&ctx->dm_y)); 138fbcfc4fcSJames Wright PetscCall(VecDestroy(&ctx->X_loc)); 139fbcfc4fcSJames Wright PetscCall(VecDestroy(&ctx->Y_loc)); 140f7325489SJames Wright 141f7325489SJames Wright // Destroy libCEED Objects 142b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&ctx->x_ceed)); 143b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&ctx->y_ceed)); 144b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&ctx->op)); 145b4c37c5cSJames Wright PetscCallCeed(ceed, CeedDestroy(&ctx->ceed)); 146f7325489SJames Wright 147fbcfc4fcSJames Wright PetscCall(PetscFree(ctx)); 148d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 149f7325489SJames Wright } 150f7325489SJames Wright 1516cbc221eSJames Wright //@brief Return VecType for the given DM 1526cbc221eSJames Wright VecType DMReturnVecType(DM dm) { 1536cbc221eSJames Wright VecType vec_type; 15487d3884fSJeremy L Thompson PetscInt ierr = DMGetVecType(dm, &vec_type); 15587d3884fSJeremy L Thompson return ierr ? NULL : vec_type; 1566cbc221eSJames Wright } 1576cbc221eSJames Wright 1586cbc221eSJames Wright /** 1596cbc221eSJames Wright * @brief Create local PETSc Vecs for CeedOperator's active input/outputs 1606cbc221eSJames Wright * 1616cbc221eSJames Wright * This is primarily used for when the active input/ouput vector does not correspond to a `DM` object, and thus `DMCreateLocalVector` or 1626cbc221eSJames Wright * `DMGetLocalVector` are not applicable. 1636cbc221eSJames Wright * For example, if statitics are being store at quadrature points, a `DM`-created `Vec` will not have the 1646cbc221eSJames Wright * correct size. 1656cbc221eSJames Wright * 1666cbc221eSJames Wright * @param[in] op Operator to make the Vecs for 16795327578SJames Wright * @param[in] vec_type `VecType` for the new Vecs 16895327578SJames Wright * @param[in] comm `MPI_Comm` for the new Vecs 1696cbc221eSJames Wright * @param[out] input Vec for CeedOperator active input 1706cbc221eSJames Wright * @param[out] output Vec for CeedOperator active output 1716cbc221eSJames Wright */ 1726cbc221eSJames Wright PetscErrorCode CeedOperatorCreateLocalVecs(CeedOperator op, VecType vec_type, MPI_Comm comm, Vec *input, Vec *output) { 1736cbc221eSJames Wright CeedSize input_size, output_size; 174b4c37c5cSJames Wright Ceed ceed; 1751fe3d3f0SJames Wright int comm_size; 1766cbc221eSJames Wright 1776cbc221eSJames Wright PetscFunctionBeginUser; 178b4c37c5cSJames Wright PetscCall(CeedOperatorGetCeed(op, &ceed)); 1791fe3d3f0SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 1801fe3d3f0SJames 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); 181b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetActiveVectorLengths(op, &input_size, &output_size)); 1826cbc221eSJames Wright if (input) { 1836cbc221eSJames Wright PetscCall(VecCreate(comm, input)); 1846cbc221eSJames Wright PetscCall(VecSetType(*input, vec_type)); 1856cbc221eSJames Wright PetscCall(VecSetSizes(*input, input_size, input_size)); 1866cbc221eSJames Wright } 1876cbc221eSJames Wright if (output) { 1886cbc221eSJames Wright PetscCall(VecCreate(comm, output)); 1896cbc221eSJames Wright PetscCall(VecSetType(*output, vec_type)); 1906cbc221eSJames Wright PetscCall(VecSetSizes(*output, output_size, output_size)); 1916cbc221eSJames Wright } 192519781aeSJames Wright PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, comm, PETSC_ERR_LIB, "Destroying Ceed object failed"); 193d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1946cbc221eSJames Wright } 1956cbc221eSJames Wright 196f7325489SJames Wright /** 197f7325489SJames Wright * @brief Apply FEM Operator defined by `OperatorApplyContext` to various input and output vectors 198f7325489SJames Wright * 1993fc17c4bSJames Wright * @param[in] X Input global `Vec`, maybe `NULL` 2003fc17c4bSJames Wright * @param[in] X_loc Input local `Vec`, maybe `NULL` 2013fc17c4bSJames Wright * @param[in] x_ceed Input `CeedVector`, maybe `CEED_VECTOR_NONE` 2023fc17c4bSJames Wright * @param[in,out] y_ceed Output `CeedVector`, maybe `CEED_VECTOR_NONE` 2033fc17c4bSJames Wright * @param[in,out] Y_loc Output local `Vec`, maybe `NULL` 2043fc17c4bSJames Wright * @param[in,out] Y Output global `Vec`, maybe `NULL` 2053fc17c4bSJames Wright * @param[in] ctx Context for the operator apply 2063fc17c4bSJames Wright * @param[in] use_apply_add Whether to use `CeedOperatorApply` or `CeedOperatorApplyAdd` 207f7325489SJames Wright */ 208f7325489SJames Wright PetscErrorCode ApplyCeedOperator_Core(Vec X, Vec X_loc, CeedVector x_ceed, CeedVector y_ceed, Vec Y_loc, Vec Y, OperatorApplyContext ctx, 209f7325489SJames Wright bool use_apply_add) { 210f7325489SJames Wright PetscMemType x_mem_type, y_mem_type; 211b4c37c5cSJames Wright Ceed ceed = ctx->ceed; 212f7325489SJames Wright 213f7325489SJames Wright PetscFunctionBeginUser; 214f7325489SJames Wright if (X) PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc)); 215a7dac1d5SJames Wright if (X_loc) PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, x_ceed)); 216f7325489SJames Wright 217a7dac1d5SJames Wright if (Y_loc) PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, y_ceed)); 218f7325489SJames Wright 219*ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_CeedOperatorApply, X, Y, 0, 0)); 2207eedc94cSJames Wright PetscCall(PetscLogGpuTimeBegin()); 221b4c37c5cSJames Wright if (use_apply_add) PetscCallCeed(ceed, CeedOperatorApplyAdd(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE)); 222b4c37c5cSJames Wright else PetscCallCeed(ceed, CeedOperatorApply(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE)); 2237eedc94cSJames Wright PetscCall(PetscLogGpuTimeEnd()); 224*ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_CeedOperatorApply, X, Y, 0, 0)); 225f7325489SJames Wright 226a7dac1d5SJames Wright if (X_loc) PetscCall(VecReadCeedToPetsc(ctx->x_ceed, x_mem_type, X_loc)); 227f7325489SJames Wright 228a7dac1d5SJames Wright if (Y_loc) PetscCall(VecCeedToPetsc(ctx->y_ceed, y_mem_type, Y_loc)); 229f7325489SJames Wright if (Y) PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y)); 230d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 231f7325489SJames Wright }; 232f7325489SJames Wright 233f7325489SJames Wright PetscErrorCode ApplyCeedOperatorGlobalToGlobal(Vec X, Vec Y, OperatorApplyContext ctx) { 234f7325489SJames Wright Vec X_loc = ctx->X_loc, Y_loc = ctx->Y_loc; 235f7325489SJames Wright 236f7325489SJames Wright PetscFunctionBeginUser; 237f7325489SJames Wright PetscCall(VecZeroEntries(Y)); 238f7325489SJames Wright 239f7325489SJames Wright // Get local vectors (if needed) 240f7325489SJames Wright if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); 241f7325489SJames Wright if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); 242f7325489SJames Wright 243f7325489SJames Wright PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false)); 244f7325489SJames Wright 245f7325489SJames Wright // Restore local vector (if needed) 246f7325489SJames Wright if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); 247f7325489SJames Wright if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); 248d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 249f7325489SJames Wright } 250f7325489SJames Wright 251f7325489SJames Wright PetscErrorCode ApplyCeedOperatorLocalToGlobal(Vec X_loc, Vec Y, OperatorApplyContext ctx) { 252f7325489SJames Wright Vec Y_loc = ctx->Y_loc; 253f7325489SJames Wright 254f7325489SJames Wright PetscFunctionBeginUser; 255f7325489SJames Wright PetscCall(VecZeroEntries(Y)); 256f7325489SJames Wright 257f7325489SJames Wright // Get local vectors (if needed) 258f7325489SJames Wright if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); 259f7325489SJames Wright 260f7325489SJames Wright PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false)); 261f7325489SJames Wright 262f7325489SJames Wright // Restore local vectors (if needed) 263f7325489SJames Wright if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); 264d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 265f7325489SJames Wright } 266f7325489SJames Wright 267f7325489SJames Wright PetscErrorCode ApplyCeedOperatorGlobalToLocal(Vec X, Vec Y_loc, OperatorApplyContext ctx) { 268f7325489SJames Wright Vec X_loc = ctx->X_loc; 269f7325489SJames Wright 270f7325489SJames Wright PetscFunctionBeginUser; 271f7325489SJames Wright // Get local vectors (if needed) 272f7325489SJames Wright if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); 273f7325489SJames Wright 274f7325489SJames Wright PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false)); 275f7325489SJames Wright 276f7325489SJames Wright // Restore local vector (if needed) 277f7325489SJames Wright if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); 278d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 279f7325489SJames Wright } 280f7325489SJames Wright 2815ce09d65SJames Wright PetscErrorCode ApplyCeedOperatorLocalToLocal(Vec X_loc, Vec Y_loc, OperatorApplyContext ctx) { 2825ce09d65SJames Wright PetscFunctionBeginUser; 2835ce09d65SJames Wright PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false)); 284d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2855ce09d65SJames Wright } 2865ce09d65SJames Wright 287f7325489SJames Wright PetscErrorCode ApplyAddCeedOperatorLocalToLocal(Vec X_loc, Vec Y_loc, OperatorApplyContext ctx) { 288f7325489SJames Wright PetscFunctionBeginUser; 289f7325489SJames Wright PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, true)); 290d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 291f7325489SJames Wright } 292f7325489SJames Wright 2933170c09fSJames Wright /** 2943170c09fSJames Wright * @brief Return Mats for KSP solve 2953170c09fSJames Wright * 2963170c09fSJames Wright * Uses command-line flag with `ksp`'s prefix to determine if mat_ceed should be used directly or whether it should be assembled. 2973170c09fSJames Wright * If `Amat` is to be assembled, then `Pmat` is equal to `Amat`. 2983170c09fSJames Wright * 2993170c09fSJames Wright * If `Amat` uses `mat_ceed`, then `Pmat` is either assembled or uses `mat_ceed` based on the preconditioner choice in `ksp`. 3003170c09fSJames Wright * 3013170c09fSJames Wright * @param[in] ksp `KSP` object for used for solving 3023170c09fSJames Wright * @param[in] mat_ceed `MATCEED` for the linear operator 3033170c09fSJames Wright * @param[in] assemble Whether to assemble `Amat` and `Pmat` if they are not `mat_ceed` 3043170c09fSJames Wright * @param[out] Amat `Mat` to be used for the solver `Amat` 3053170c09fSJames Wright * @param[out] Pmat `Mat` to be used for the solver `Pmat` 3063170c09fSJames Wright */ 3073170c09fSJames Wright PetscErrorCode CreateSolveOperatorsFromMatCeed(KSP ksp, Mat mat_ceed, PetscBool assemble, Mat *Amat, Mat *Pmat) { 3083170c09fSJames Wright PetscBool use_matceed_pmat, assemble_amat = PETSC_FALSE; 3093170c09fSJames Wright 3103170c09fSJames Wright PetscFunctionBeginUser; 3113170c09fSJames Wright { // Determine if Amat should be MATCEED or assembled 3123170c09fSJames Wright const char *ksp_prefix = NULL; 3133170c09fSJames Wright 3143170c09fSJames Wright PetscCall(KSPGetOptionsPrefix(ksp, &ksp_prefix)); 3153170c09fSJames Wright PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), ksp_prefix, "", NULL); 3163170c09fSJames Wright PetscCall(PetscOptionsBool("-matceed_assemble_amat", "Assemble the A matrix for KSP solve", NULL, assemble_amat, &assemble_amat, NULL)); 3173170c09fSJames Wright PetscOptionsEnd(); 3183170c09fSJames Wright } 3193170c09fSJames Wright 3203170c09fSJames Wright if (assemble_amat) { 321c1c64bfcSJames Wright PetscCall(MatCeedCreateMatCOO(mat_ceed, Amat)); 3223170c09fSJames Wright if (assemble) PetscCall(MatCeedAssembleCOO(mat_ceed, *Amat)); 3233170c09fSJames Wright 3243170c09fSJames Wright PetscCall(PetscObjectReference((PetscObject)*Amat)); 3253170c09fSJames Wright *Pmat = *Amat; 3263170c09fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3273170c09fSJames Wright } else { 3283170c09fSJames Wright PetscCall(PetscObjectReference((PetscObject)mat_ceed)); 3293170c09fSJames Wright *Amat = mat_ceed; 3303170c09fSJames Wright } 3313170c09fSJames Wright 3323170c09fSJames Wright { // Determine if Pmat should be MATCEED or assembled 3333170c09fSJames Wright PC pc; 3343170c09fSJames Wright PCType pc_type; 3353170c09fSJames Wright 3363170c09fSJames Wright PetscCall(KSPGetPC(ksp, &pc)); 3373170c09fSJames Wright PetscCall(PCGetType(pc, &pc_type)); 3389ff51368SJames Wright PetscCall(PetscStrcmpAny(pc_type, &use_matceed_pmat, PCNONE, PCJACOBI, PCVPBJACOBI, PCPBJACOBI, PCBJACOBI, "")); 3393170c09fSJames Wright } 3403170c09fSJames Wright 3413170c09fSJames Wright if (use_matceed_pmat) { 3423170c09fSJames Wright PetscCall(PetscObjectReference((PetscObject)mat_ceed)); 3433170c09fSJames Wright *Pmat = mat_ceed; 3443170c09fSJames Wright } else { 345c1c64bfcSJames Wright PetscCall(MatCeedCreateMatCOO(mat_ceed, Pmat)); 3463170c09fSJames Wright if (assemble) PetscCall(MatCeedAssembleCOO(mat_ceed, *Pmat)); 3473170c09fSJames Wright } 3483170c09fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3493170c09fSJames Wright } 3503170c09fSJames Wright 3513170c09fSJames Wright /** 3523170c09fSJames Wright * @brief Runs KSPSetFromOptions and sets Operators based on mat_ceed 3533170c09fSJames Wright * 3543170c09fSJames Wright * See CreateSolveOperatorsFromMatCeed for details on how the KSPSolve operators are set. 3553170c09fSJames Wright * 3563170c09fSJames Wright * @param[in] ksp `KSP` of the solve 3573170c09fSJames Wright * @param[in] mat_ceed `MatCeed` linear operator to solve for 3583170c09fSJames Wright */ 3593170c09fSJames Wright PetscErrorCode KSPSetFromOptions_WithMatCeed(KSP ksp, Mat mat_ceed) { 3603170c09fSJames Wright Mat Amat, Pmat; 3613170c09fSJames Wright 3623170c09fSJames Wright PetscFunctionBeginUser; 3633170c09fSJames Wright PetscCall(KSPSetFromOptions(ksp)); 3643170c09fSJames Wright PetscCall(CreateSolveOperatorsFromMatCeed(ksp, mat_ceed, PETSC_TRUE, &Amat, &Pmat)); 3653170c09fSJames Wright PetscCall(KSPSetOperators(ksp, Amat, Pmat)); 3663170c09fSJames Wright PetscCall(MatDestroy(&Amat)); 3673170c09fSJames Wright PetscCall(MatDestroy(&Pmat)); 3683170c09fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3693170c09fSJames Wright } 370f7325489SJames Wright // ----------------------------------------------------------------------------- 371