xref: /libCEED/examples/fluids/src/mat-ceed.c (revision ed094490f53e580908aa80e9fe815a6fd76d7526)
1c8564c30SJames Wright /// @file
25037d55dSJames Wright /// MatCEED implementation
3c8564c30SJames Wright 
4c8564c30SJames Wright #include <ceed.h>
5c8564c30SJames Wright #include <ceed/backend.h>
6c8564c30SJames Wright #include <mat-ceed-impl.h>
7c8564c30SJames Wright #include <mat-ceed.h>
85037d55dSJames Wright #include <petsc-ceed-utils.h>
95037d55dSJames Wright #include <petsc-ceed.h>
10243afec9SJames Wright #include <petscdm.h>
11243afec9SJames Wright #include <petscmat.h>
12243afec9SJames Wright #include <stdbool.h>
13243afec9SJames Wright #include <stdio.h>
14c8564c30SJames Wright #include <stdlib.h>
15c8564c30SJames Wright #include <string.h>
16c8564c30SJames Wright 
17c8564c30SJames Wright PetscClassId  MATCEED_CLASSID;
1843327b86SJames Wright PetscLogEvent MATCEED_MULT, MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE, MATCEED_MULT_TRANSPOSE_CEEDOP, MATCEED_ASSEMBLE_DIAGONAL,
1943327b86SJames Wright     MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, MATCEED_SETUP_PBDIAGONAL, MATCEED_SETUP_PBDIAGONAL_CEEDOP, MATCEED_ASSEMBLE_PBDIAGONAL,
2043327b86SJames Wright     MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, MATCEED_SETUP_FULL, MATCEED_SETUP_FULL_CEEDOP, MATCEED_ASSEMBLE_FULL, MATCEED_ASSEMBLE_FULL_CEEDOP;
21c8564c30SJames Wright 
22c8564c30SJames Wright /**
23c8564c30SJames Wright   @brief Register MATCEED log events.
24c8564c30SJames Wright 
25c8564c30SJames Wright   Not collective across MPI processes.
26c8564c30SJames Wright 
27c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
28c8564c30SJames Wright **/
29c8564c30SJames Wright static PetscErrorCode MatCeedRegisterLogEvents() {
305037d55dSJames Wright   static PetscBool registered = PETSC_FALSE;
31c8564c30SJames Wright 
32c8564c30SJames Wright   PetscFunctionBeginUser;
33c8564c30SJames Wright   if (registered) PetscFunctionReturn(PETSC_SUCCESS);
3443327b86SJames Wright   PetscCall(PetscClassIdRegister("MatCEED", &MATCEED_CLASSID));
3543327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMul", MATCEED_CLASSID, &MATCEED_MULT));
3643327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMulCeed", MATCEED_CLASSID, &MATCEED_MULT_CEEDOP));
3743327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMulT", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE));
3843327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMulTCeed", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE_CEEDOP));
3943327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmDiag", MATCEED_CLASSID, &MATCEED_ASSEMBLE_DIAGONAL));
4043327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSU", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL));
4143327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL_CEEDOP));
4243327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBD", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL));
4343327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBDCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP));
4443327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmSU", MATCEED_CLASSID, &MATCEED_SETUP_FULL));
4543327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_FULL_CEEDOP));
4643327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsm", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL));
4743327b86SJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL_CEEDOP));
485037d55dSJames Wright   registered = PETSC_TRUE;
49c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
50c8564c30SJames Wright }
51c8564c30SJames Wright 
52c8564c30SJames Wright /**
53c8564c30SJames Wright   @brief Assemble the point block diagonal of a `MATCEED` into a `MATAIJ` or similar.
54c8564c30SJames Wright          The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`.
55c8564c30SJames Wright          The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail.
56c8564c30SJames Wright 
57c8564c30SJames Wright   Collective across MPI processes.
58c8564c30SJames Wright 
59c8564c30SJames Wright   @param[in]      mat_ceed  `MATCEED` to assemble
60c8564c30SJames Wright   @param[in,out]  mat_coo   `MATAIJ` or similar to assemble into
61c8564c30SJames Wright 
62c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
63c8564c30SJames Wright **/
64c8564c30SJames Wright static PetscErrorCode MatCeedAssemblePointBlockDiagonalCOO(Mat mat_ceed, Mat mat_coo) {
65c8564c30SJames Wright   MatCeedContext ctx;
66c8564c30SJames Wright 
67c8564c30SJames Wright   PetscFunctionBeginUser;
68c8564c30SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
69c8564c30SJames Wright 
70c8564c30SJames Wright   // Check if COO pattern set
71c8564c30SJames Wright   {
72c8564c30SJames Wright     PetscInt index = -1;
73c8564c30SJames Wright 
74c8564c30SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) {
75c8564c30SJames Wright       if (ctx->mats_assembled_pbd[i] == mat_coo) index = i;
76c8564c30SJames Wright     }
77c8564c30SJames Wright     if (index == -1) {
78c8564c30SJames Wright       PetscInt     *rows_petsc = NULL, *cols_petsc = NULL;
79c8564c30SJames Wright       CeedInt      *rows_ceed, *cols_ceed;
80c8564c30SJames Wright       PetscCount    num_entries;
81c8564c30SJames Wright       PetscLogStage stage_amg_setup;
82c8564c30SJames Wright 
83c8564c30SJames Wright       // -- Assemble sparsity pattern if mat hasn't been assembled before
8443327b86SJames Wright       PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup));
85c8564c30SJames Wright       if (stage_amg_setup == -1) {
8643327b86SJames Wright         PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup));
87c8564c30SJames Wright       }
88c8564c30SJames Wright       PetscCall(PetscLogStagePush(stage_amg_setup));
8943327b86SJames Wright       PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
9043327b86SJames Wright       PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
9154831c5fSJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed));
9243327b86SJames Wright       PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
93d0593705SJames Wright       PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc));
94d0593705SJames Wright       PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc));
95c8564c30SJames Wright       PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc));
96c8564c30SJames Wright       free(rows_petsc);
97c8564c30SJames Wright       free(cols_petsc);
9854831c5fSJames Wright       if (!ctx->coo_values_pbd) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_pbd));
99c8564c30SJames Wright       PetscCall(PetscRealloc(++ctx->num_mats_assembled_pbd * sizeof(Mat), &ctx->mats_assembled_pbd));
100c8564c30SJames Wright       ctx->mats_assembled_pbd[ctx->num_mats_assembled_pbd - 1] = mat_coo;
10143327b86SJames Wright       PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
102c8564c30SJames Wright       PetscCall(PetscLogStagePop());
103c8564c30SJames Wright     }
104c8564c30SJames Wright   }
105c8564c30SJames Wright 
106c8564c30SJames Wright   // Assemble mat_ceed
10743327b86SJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
108c8564c30SJames Wright   PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY));
109c8564c30SJames Wright   {
110c8564c30SJames Wright     const CeedScalar *values;
111c8564c30SJames Wright     MatType           mat_type;
112c8564c30SJames Wright     CeedMemType       mem_type = CEED_MEM_HOST;
113c8564c30SJames Wright     PetscBool         is_spd, is_spd_known;
114c8564c30SJames Wright 
115c8564c30SJames Wright     PetscCall(MatGetType(mat_coo, &mat_type));
116c8564c30SJames Wright     if (strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
117c8564c30SJames Wright     else if (strstr(mat_type, "kokkos")) mem_type = CEED_MEM_DEVICE;
118c8564c30SJames Wright     else mem_type = CEED_MEM_HOST;
119c8564c30SJames Wright 
12043327b86SJames Wright     PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
12154831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(ctx->op_mult, ctx->coo_values_pbd, CEED_REQUEST_IMMEDIATE));
12243327b86SJames Wright     PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
12354831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_pbd, mem_type, &values));
124c8564c30SJames Wright     PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES));
125c8564c30SJames Wright     PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd));
126c8564c30SJames Wright     if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd));
12754831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_pbd, &values));
128c8564c30SJames Wright   }
129c8564c30SJames Wright   PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY));
13043327b86SJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
131c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
132c8564c30SJames Wright }
133c8564c30SJames Wright 
134c8564c30SJames Wright /**
135c8564c30SJames Wright   @brief Assemble inner `Mat` for diagonal `PC` operations
136c8564c30SJames Wright 
137c8564c30SJames Wright   Collective across MPI processes.
138c8564c30SJames Wright 
139c8564c30SJames Wright   @param[in]   mat_ceed      `MATCEED` to invert
140c8564c30SJames Wright   @param[in]   use_ceed_pbd  Boolean flag to use libCEED PBD assembly
141c8564c30SJames Wright   @param[out]  mat_inner     Inner `Mat` for diagonal operations
142c8564c30SJames Wright 
143c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
144c8564c30SJames Wright **/
145c8564c30SJames Wright static PetscErrorCode MatCeedAssembleInnerBlockDiagonalMat(Mat mat_ceed, PetscBool use_ceed_pbd, Mat *mat_inner) {
146c8564c30SJames Wright   MatCeedContext ctx;
147c8564c30SJames Wright 
148c8564c30SJames Wright   PetscFunctionBeginUser;
149c8564c30SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
150c8564c30SJames Wright   if (use_ceed_pbd) {
151c8564c30SJames Wright     // Check if COO pattern set
1525037d55dSJames Wright     if (!ctx->mat_assembled_pbd_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_pbd_internal));
153c8564c30SJames Wright 
154c8564c30SJames Wright     // Assemble mat_assembled_full_internal
155c8564c30SJames Wright     PetscCall(MatCeedAssemblePointBlockDiagonalCOO(mat_ceed, ctx->mat_assembled_pbd_internal));
156c8564c30SJames Wright     if (mat_inner) *mat_inner = ctx->mat_assembled_pbd_internal;
157c8564c30SJames Wright   } else {
158c8564c30SJames Wright     // Check if COO pattern set
1595037d55dSJames Wright     if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal));
160c8564c30SJames Wright 
161c8564c30SJames Wright     // Assemble mat_assembled_full_internal
162c8564c30SJames Wright     PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal));
163c8564c30SJames Wright     if (mat_inner) *mat_inner = ctx->mat_assembled_full_internal;
164c8564c30SJames Wright   }
165c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
166c8564c30SJames Wright }
167c8564c30SJames Wright 
168c8564c30SJames Wright /**
169243afec9SJames Wright   @brief Get `MATCEED` variable block diagonal for Jacobi.
170243afec9SJames Wright 
171243afec9SJames Wright   Collective across MPI processes.
172243afec9SJames Wright 
173243afec9SJames Wright   @param[in]   mat_ceed   `MATCEED` to invert
174243afec9SJames Wright   @param[out]  mat_vblock  The variable diagonal block matrix
175243afec9SJames Wright 
176243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
177243afec9SJames Wright **/
178243afec9SJames Wright static PetscErrorCode MatGetVariableBlockDiagonal_Ceed(Mat mat_ceed, Mat *mat_vblock) {
179243afec9SJames Wright   MatCeedContext ctx;
180243afec9SJames Wright 
181243afec9SJames Wright   PetscFunctionBeginUser;
182243afec9SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
183243afec9SJames Wright 
184243afec9SJames Wright   // Assemble inner mat if needed
185243afec9SJames Wright   PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_vpbd_valid, mat_vblock));
186243afec9SJames Wright   PetscCall(PetscObjectReference((PetscObject)*mat_vblock));
187243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
188243afec9SJames Wright }
189243afec9SJames Wright 
190243afec9SJames Wright /**
191243afec9SJames Wright   @brief Get `MATCEED` block diagonal for Jacobi.
192243afec9SJames Wright 
193243afec9SJames Wright   Collective across MPI processes.
194243afec9SJames Wright 
195243afec9SJames Wright   @param[in]   mat_ceed   `MATCEED` to invert
196243afec9SJames Wright   @param[out]  mat_block  The variable diagonal block matrix
197243afec9SJames Wright 
198243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
199243afec9SJames Wright **/
200243afec9SJames Wright static PetscErrorCode MatGetBlockDiagonal_Ceed(Mat mat_ceed, Mat *mat_block) {
201243afec9SJames Wright   MatCeedContext ctx;
202243afec9SJames Wright 
203243afec9SJames Wright   PetscFunctionBeginUser;
204243afec9SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
205243afec9SJames Wright 
206243afec9SJames Wright   // Assemble inner mat if needed
207243afec9SJames Wright   PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_pbd_valid, mat_block));
208243afec9SJames Wright   PetscCall(PetscObjectReference((PetscObject)*mat_block));
209243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
210243afec9SJames Wright }
211243afec9SJames Wright 
212243afec9SJames Wright /**
213243afec9SJames Wright   @brief Get on-process diagonal block of `MATCEED`
214243afec9SJames Wright 
215243afec9SJames Wright   This is a block per-process of the diagonal of the global matrix.
216243afec9SJames Wright   This is NOT the diagonal blocks associated with the block size of the matrix (i.e. `MatSetBlockSize()` has no effect on this function).
217c8564c30SJames Wright 
218c8564c30SJames Wright   Collective across MPI processes.
219c8564c30SJames Wright 
220c8564c30SJames Wright   @param[in]   mat_ceed   `MATCEED` to invert
221c8564c30SJames Wright   @param[out]  mat_block  The diagonal block matrix
222c8564c30SJames Wright 
223c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
224c8564c30SJames Wright **/
225c8564c30SJames Wright static PetscErrorCode MatGetDiagonalBlock_Ceed(Mat mat_ceed, Mat *mat_block) {
226c8564c30SJames Wright   MatCeedContext ctx;
227c8564c30SJames Wright 
228c8564c30SJames Wright   PetscFunctionBeginUser;
229c8564c30SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
230c8564c30SJames Wright 
231243afec9SJames Wright   // Check if COO pattern set
232243afec9SJames Wright   if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal));
233c8564c30SJames Wright 
234243afec9SJames Wright   // Assemble mat_assembled_full_internal
235243afec9SJames Wright   PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal));
236c8564c30SJames Wright 
237243afec9SJames Wright   // Get diagonal block
238243afec9SJames Wright   PetscCall(MatGetDiagonalBlock(ctx->mat_assembled_full_internal, mat_block));
239c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
240c8564c30SJames Wright }
241c8564c30SJames Wright 
2423933d9a0SJames Wright /**
2433933d9a0SJames Wright   @brief View `MATCEED`.
2443933d9a0SJames Wright 
2453933d9a0SJames Wright   Collective across MPI processes.
2463933d9a0SJames Wright 
2473933d9a0SJames Wright   @param[in]   mat_ceed  `MATCEED` to view
2483933d9a0SJames Wright   @param[in]   viewer    The visualization context
2493933d9a0SJames Wright 
2503933d9a0SJames Wright   @return An error code: 0 - success, otherwise - failure
2513933d9a0SJames Wright **/
2523933d9a0SJames Wright static PetscErrorCode MatView_Ceed(Mat mat_ceed, PetscViewer viewer) {
2533933d9a0SJames Wright   PetscBool         is_ascii;
2543933d9a0SJames Wright   PetscViewerFormat format;
255243afec9SJames Wright   PetscMPIInt       size, rank;
2563933d9a0SJames Wright   MatCeedContext    ctx;
2573933d9a0SJames Wright 
2583933d9a0SJames Wright   PetscFunctionBeginUser;
2593933d9a0SJames Wright   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2603933d9a0SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
2613933d9a0SJames Wright   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat_ceed), &viewer));
2623933d9a0SJames Wright 
2633933d9a0SJames Wright   PetscCall(PetscViewerGetFormat(viewer, &format));
2643933d9a0SJames Wright   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat_ceed), &size));
2653933d9a0SJames Wright   if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS);
2663933d9a0SJames Wright 
267243afec9SJames Wright   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat_ceed), &rank));
268243afec9SJames Wright   if (rank != 0) PetscFunctionReturn(PETSC_SUCCESS);
269243afec9SJames Wright 
2703933d9a0SJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
2713933d9a0SJames Wright   {
272243afec9SJames Wright     PetscBool is_detailed     = format == PETSC_VIEWER_ASCII_INFO_DETAIL;
273243afec9SJames Wright     char      rank_string[16] = {'\0'};
2743933d9a0SJames Wright     FILE     *file;
2753933d9a0SJames Wright 
276243afec9SJames Wright     PetscCall(PetscViewerASCIIPrintf(viewer, "MatCEED:\n"));
277243afec9SJames Wright     PetscCall(PetscViewerASCIIPushTab(viewer));  // MatCEED
278243afec9SJames Wright     PetscCall(PetscViewerASCIIPrintf(viewer, "Default COO MatType: %s\n", ctx->coo_mat_type));
279243afec9SJames Wright     PetscCall(PetscSNPrintf(rank_string, 16, "on Rank %d", rank));
280243afec9SJames Wright     PetscCall(PetscViewerASCIIPrintf(viewer, "CeedOperator Apply %s:\n", is_detailed ? rank_string : "Summary"));
281243afec9SJames Wright     PetscCall(PetscViewerASCIIPrintf(viewer, "libCEED PB Diagonal Assembly: %s\n", ctx->is_ceed_pbd_valid ? "True" : "False"));
282243afec9SJames Wright     PetscCall(PetscViewerASCIIPrintf(viewer, "libCEED VPB Diagonal Assembly: %s\n", ctx->is_ceed_vpbd_valid ? "True" : "False"));
2833933d9a0SJames Wright     PetscCall(PetscViewerASCIIGetPointer(viewer, &file));
284243afec9SJames Wright     PetscCall(PetscViewerASCIIPushTab(viewer));  // CeedOperator
285243afec9SJames Wright     if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult, file));
286243afec9SJames Wright     else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult, file));
287243afec9SJames Wright     PetscCall(PetscViewerASCIIPopTab(viewer));  // CeedOperator
2883933d9a0SJames Wright     if (ctx->op_mult_transpose) {
289243afec9SJames Wright       PetscCall(PetscViewerASCIIPrintf(viewer, "CeedOperator ApplyTranspose %s:\n", is_detailed ? rank_string : "Summary"));
290243afec9SJames Wright       PetscCall(PetscViewerASCIIPushTab(viewer));  // CeedOperator
291243afec9SJames Wright       if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult_transpose, file));
292243afec9SJames Wright       else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult_transpose, file));
293243afec9SJames Wright       PetscCall(PetscViewerASCIIPopTab(viewer));  // CeedOperator
2943933d9a0SJames Wright     }
295243afec9SJames Wright     PetscCall(PetscViewerASCIIPopTab(viewer));  // MatCEED
2963933d9a0SJames Wright   }
2973933d9a0SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2983933d9a0SJames Wright }
2993933d9a0SJames Wright 
300c8564c30SJames Wright // -----------------------------------------------------------------------------
301c8564c30SJames Wright // MatCeed
302c8564c30SJames Wright // -----------------------------------------------------------------------------
303c8564c30SJames Wright 
304c8564c30SJames Wright /**
305c8564c30SJames Wright   @brief Create PETSc `Mat` from libCEED operators.
306c8564c30SJames Wright 
307c8564c30SJames Wright   Collective across MPI processes.
308c8564c30SJames Wright 
309c8564c30SJames Wright   @param[in]   dm_x                      Input `DM`
310c8564c30SJames Wright   @param[in]   dm_y                      Output `DM`
311c8564c30SJames Wright   @param[in]   op_mult                   `CeedOperator` for forward evaluation
312c8564c30SJames Wright   @param[in]   op_mult_transpose         `CeedOperator` for transpose evaluation
313c8564c30SJames Wright   @param[out]  mat                        New MatCeed
314c8564c30SJames Wright 
315c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
316c8564c30SJames Wright **/
317243afec9SJames Wright PetscErrorCode MatCreateCeed(DM dm_x, DM dm_y, CeedOperator op_mult, CeedOperator op_mult_transpose, Mat *mat) {
318c8564c30SJames Wright   PetscInt       X_l_size, X_g_size, Y_l_size, Y_g_size;
319c8564c30SJames Wright   VecType        vec_type;
320c8564c30SJames Wright   MatCeedContext ctx;
321c8564c30SJames Wright 
322c8564c30SJames Wright   PetscFunctionBeginUser;
323c8564c30SJames Wright   PetscCall(MatCeedRegisterLogEvents());
324c8564c30SJames Wright 
325c8564c30SJames Wright   // Collect context data
326c8564c30SJames Wright   PetscCall(DMGetVecType(dm_x, &vec_type));
327c8564c30SJames Wright   {
328c8564c30SJames Wright     Vec X;
329c8564c30SJames Wright 
330c8564c30SJames Wright     PetscCall(DMGetGlobalVector(dm_x, &X));
331c8564c30SJames Wright     PetscCall(VecGetSize(X, &X_g_size));
332c8564c30SJames Wright     PetscCall(VecGetLocalSize(X, &X_l_size));
333c8564c30SJames Wright     PetscCall(DMRestoreGlobalVector(dm_x, &X));
334c8564c30SJames Wright   }
335c8564c30SJames Wright   if (dm_y) {
336c8564c30SJames Wright     Vec Y;
337c8564c30SJames Wright 
338c8564c30SJames Wright     PetscCall(DMGetGlobalVector(dm_y, &Y));
339c8564c30SJames Wright     PetscCall(VecGetSize(Y, &Y_g_size));
340c8564c30SJames Wright     PetscCall(VecGetLocalSize(Y, &Y_l_size));
341c8564c30SJames Wright     PetscCall(DMRestoreGlobalVector(dm_y, &Y));
342c8564c30SJames Wright   } else {
343c8564c30SJames Wright     dm_y     = dm_x;
344c8564c30SJames Wright     Y_g_size = X_g_size;
345c8564c30SJames Wright     Y_l_size = X_l_size;
346c8564c30SJames Wright   }
3475037d55dSJames Wright 
348c8564c30SJames Wright   // Create context
349c8564c30SJames Wright   {
350c8564c30SJames Wright     Vec X_loc, Y_loc_transpose = NULL;
351c8564c30SJames Wright 
352c8564c30SJames Wright     PetscCall(DMCreateLocalVector(dm_x, &X_loc));
353c8564c30SJames Wright     PetscCall(VecZeroEntries(X_loc));
354c8564c30SJames Wright     if (op_mult_transpose) {
355c8564c30SJames Wright       PetscCall(DMCreateLocalVector(dm_y, &Y_loc_transpose));
356c8564c30SJames Wright       PetscCall(VecZeroEntries(Y_loc_transpose));
357c8564c30SJames Wright     }
35843327b86SJames Wright     PetscCall(MatCeedContextCreate(dm_x, dm_y, X_loc, Y_loc_transpose, op_mult, op_mult_transpose, MATCEED_MULT, MATCEED_MULT_TRANSPOSE,
35943327b86SJames Wright                                    MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE_CEEDOP, &ctx));
360c8564c30SJames Wright     PetscCall(VecDestroy(&X_loc));
361c8564c30SJames Wright     PetscCall(VecDestroy(&Y_loc_transpose));
362c8564c30SJames Wright   }
363c8564c30SJames Wright 
364c8564c30SJames Wright   // Create mat
365c8564c30SJames Wright   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm_x), Y_l_size, X_l_size, Y_g_size, X_g_size, ctx, mat));
366c8564c30SJames Wright   PetscCall(PetscObjectChangeTypeName((PetscObject)*mat, MATCEED));
367c8564c30SJames Wright   // -- Set block and variable block sizes
368c8564c30SJames Wright   if (dm_x == dm_y) {
369c8564c30SJames Wright     MatType dm_mat_type, dm_mat_type_copy;
370c8564c30SJames Wright     Mat     temp_mat;
371c8564c30SJames Wright 
372c8564c30SJames Wright     PetscCall(DMGetMatType(dm_x, &dm_mat_type));
373c8564c30SJames Wright     PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy));
374c8564c30SJames Wright     PetscCall(DMSetMatType(dm_x, MATAIJ));
375c8564c30SJames Wright     PetscCall(DMCreateMatrix(dm_x, &temp_mat));
376c8564c30SJames Wright     PetscCall(DMSetMatType(dm_x, dm_mat_type_copy));
377c8564c30SJames Wright     PetscCall(PetscFree(dm_mat_type_copy));
378c8564c30SJames Wright 
379c8564c30SJames Wright     {
380c8564c30SJames Wright       PetscInt        block_size, num_blocks, max_vblock_size = PETSC_INT_MAX;
381c8564c30SJames Wright       const PetscInt *vblock_sizes;
382c8564c30SJames Wright 
383c8564c30SJames Wright       // -- Get block sizes
384c8564c30SJames Wright       PetscCall(MatGetBlockSize(temp_mat, &block_size));
385c8564c30SJames Wright       PetscCall(MatGetVariableBlockSizes(temp_mat, &num_blocks, &vblock_sizes));
386c8564c30SJames Wright       {
387c8564c30SJames Wright         PetscInt local_min_max[2] = {0}, global_min_max[2] = {0, PETSC_INT_MAX};
388c8564c30SJames Wright 
389c8564c30SJames Wright         for (PetscInt i = 0; i < num_blocks; i++) local_min_max[1] = PetscMax(local_min_max[1], vblock_sizes[i]);
390c8564c30SJames Wright         PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_min_max, global_min_max));
391c8564c30SJames Wright         max_vblock_size = global_min_max[1];
392c8564c30SJames Wright       }
393c8564c30SJames Wright 
394c8564c30SJames Wright       // -- Copy block sizes
395c8564c30SJames Wright       if (block_size > 1) PetscCall(MatSetBlockSize(*mat, block_size));
396c8564c30SJames Wright       if (num_blocks) PetscCall(MatSetVariableBlockSizes(*mat, num_blocks, (PetscInt *)vblock_sizes));
397c8564c30SJames Wright 
398c8564c30SJames Wright       // -- Check libCEED compatibility
399c8564c30SJames Wright       {
400c8564c30SJames Wright         bool is_composite;
401c8564c30SJames Wright 
402c8564c30SJames Wright         ctx->is_ceed_pbd_valid  = PETSC_TRUE;
403c8564c30SJames Wright         ctx->is_ceed_vpbd_valid = PETSC_TRUE;
40454831c5fSJames Wright         PetscCallCeed(ctx->ceed, CeedOperatorIsComposite(op_mult, &is_composite));
405c8564c30SJames Wright         if (is_composite) {
406c8564c30SJames Wright           CeedInt       num_sub_operators;
407c8564c30SJames Wright           CeedOperator *sub_operators;
408c8564c30SJames Wright 
409*ed094490SJeremy L Thompson           PetscCallCeed(ctx->ceed, CeedOperatorCompositeGetNumSub(op_mult, &num_sub_operators));
410*ed094490SJeremy L Thompson           PetscCallCeed(ctx->ceed, CeedOperatorCompositeGetSubList(op_mult, &sub_operators));
411c8564c30SJames Wright           for (CeedInt i = 0; i < num_sub_operators; i++) {
412c8564c30SJames Wright             CeedInt                  num_bases, num_comp;
413c8564c30SJames Wright             CeedBasis               *active_bases;
414c8564c30SJames Wright             CeedOperatorAssemblyData assembly_data;
415c8564c30SJames Wright 
41654831c5fSJames Wright             PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(sub_operators[i], &assembly_data));
41754831c5fSJames Wright             PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL));
41854831c5fSJames Wright             PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp));
419c8564c30SJames Wright             if (num_bases > 1) {
420c8564c30SJames Wright               ctx->is_ceed_pbd_valid  = PETSC_FALSE;
421c8564c30SJames Wright               ctx->is_ceed_vpbd_valid = PETSC_FALSE;
422c8564c30SJames Wright             }
423c8564c30SJames Wright             if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE;
424c8564c30SJames Wright             if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE;
425c8564c30SJames Wright           }
426c8564c30SJames Wright         } else {
427c8564c30SJames Wright           // LCOV_EXCL_START
428c8564c30SJames Wright           CeedInt                  num_bases, num_comp;
429c8564c30SJames Wright           CeedBasis               *active_bases;
430c8564c30SJames Wright           CeedOperatorAssemblyData assembly_data;
431c8564c30SJames Wright 
43254831c5fSJames Wright           PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(op_mult, &assembly_data));
43354831c5fSJames Wright           PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL));
43454831c5fSJames Wright           PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp));
435c8564c30SJames Wright           if (num_bases > 1) {
436c8564c30SJames Wright             ctx->is_ceed_pbd_valid  = PETSC_FALSE;
437c8564c30SJames Wright             ctx->is_ceed_vpbd_valid = PETSC_FALSE;
438c8564c30SJames Wright           }
439c8564c30SJames Wright           if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE;
440c8564c30SJames Wright           if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE;
441c8564c30SJames Wright           // LCOV_EXCL_STOP
442c8564c30SJames Wright         }
443c8564c30SJames Wright         {
444c8564c30SJames Wright           PetscInt local_is_valid[2], global_is_valid[2];
445c8564c30SJames Wright 
446c8564c30SJames Wright           local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_pbd_valid;
447c8564c30SJames Wright           PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid));
448c8564c30SJames Wright           ctx->is_ceed_pbd_valid = global_is_valid[0];
449c8564c30SJames Wright           local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_vpbd_valid;
450c8564c30SJames Wright           PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid));
451c8564c30SJames Wright           ctx->is_ceed_vpbd_valid = global_is_valid[0];
452c8564c30SJames Wright         }
453c8564c30SJames Wright       }
454c8564c30SJames Wright     }
455c8564c30SJames Wright     PetscCall(MatDestroy(&temp_mat));
456c8564c30SJames Wright   }
457c8564c30SJames Wright   // -- Set internal mat type
458c8564c30SJames Wright   {
459c8564c30SJames Wright     VecType vec_type;
4605037d55dSJames Wright     MatType coo_mat_type;
461c8564c30SJames Wright 
462c8564c30SJames Wright     PetscCall(VecGetType(ctx->X_loc, &vec_type));
4635037d55dSJames Wright     if (strstr(vec_type, VECCUDA)) coo_mat_type = MATAIJCUSPARSE;
4645037d55dSJames Wright     else if (strstr(vec_type, VECKOKKOS)) coo_mat_type = MATAIJKOKKOS;
4655037d55dSJames Wright     else coo_mat_type = MATAIJ;
4665037d55dSJames Wright     PetscCall(PetscStrallocpy(coo_mat_type, &ctx->coo_mat_type));
467c8564c30SJames Wright   }
468c8564c30SJames Wright   // -- Set mat operations
469243afec9SJames Wright   PetscCall(MatShellSetContextDestroy(*mat, (PetscCtxDestroyFn *)MatCeedContextDestroy));
4703933d9a0SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_VIEW, (void (*)(void))MatView_Ceed));
471c8564c30SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_MULT, (void (*)(void))MatMult_Ceed));
472c8564c30SJames Wright   if (op_mult_transpose) PetscCall(MatShellSetOperation(*mat, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed));
473c8564c30SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed));
474c8564c30SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed));
475243afec9SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_GET_BLOCK_DIAGONAL, (void (*)(void))MatGetBlockDiagonal_Ceed));
476243afec9SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_GET_VBLOCK_DIAGONAL, (void (*)(void))MatGetVariableBlockDiagonal_Ceed));
477c8564c30SJames Wright   PetscCall(MatShellSetVecType(*mat, vec_type));
478c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
479c8564c30SJames Wright }
480c8564c30SJames Wright 
481c8564c30SJames Wright /**
482c8564c30SJames Wright   @brief Copy `MATCEED` into a compatible `Mat` with type `MatShell` or `MATCEED`.
483c8564c30SJames Wright 
484c8564c30SJames Wright   Collective across MPI processes.
485c8564c30SJames Wright 
486c8564c30SJames Wright   @param[in]   mat_ceed   `MATCEED` to copy from
487c8564c30SJames Wright   @param[out]  mat_other  `MatShell` or `MATCEED` to copy into
488c8564c30SJames Wright 
489c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
490c8564c30SJames Wright **/
491c8564c30SJames Wright PetscErrorCode MatCeedCopy(Mat mat_ceed, Mat mat_other) {
492c8564c30SJames Wright   PetscFunctionBeginUser;
493c8564c30SJames Wright   PetscCall(MatCeedRegisterLogEvents());
494c8564c30SJames Wright 
495c8564c30SJames Wright   // Check type compatibility
496c8564c30SJames Wright   {
4975037d55dSJames Wright     PetscBool is_matceed = PETSC_FALSE, is_matshell = PETSC_FALSE;
498c8564c30SJames Wright     MatType   mat_type_ceed, mat_type_other;
499c8564c30SJames Wright 
500c8564c30SJames Wright     PetscCall(MatGetType(mat_ceed, &mat_type_ceed));
5015037d55dSJames Wright     PetscCall(PetscStrcmp(mat_type_ceed, MATCEED, &is_matceed));
5025037d55dSJames Wright     PetscCheck(is_matceed, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_ceed must have type " MATCEED);
5035037d55dSJames Wright     PetscCall(MatGetType(mat_other, &mat_type_other));
5045037d55dSJames Wright     PetscCall(PetscStrcmp(mat_type_other, MATCEED, &is_matceed));
5055037d55dSJames Wright     PetscCall(PetscStrcmp(mat_type_other, MATSHELL, &is_matceed));
5065037d55dSJames Wright     PetscCheck(is_matceed || is_matshell, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_other must have type " MATCEED " or " MATSHELL);
507c8564c30SJames Wright   }
508c8564c30SJames Wright 
509c8564c30SJames Wright   // Check dimension compatibility
510c8564c30SJames Wright   {
511c8564c30SJames Wright     PetscInt X_l_ceed_size, X_g_ceed_size, Y_l_ceed_size, Y_g_ceed_size, X_l_other_size, X_g_other_size, Y_l_other_size, Y_g_other_size;
512c8564c30SJames Wright 
513c8564c30SJames Wright     PetscCall(MatGetSize(mat_ceed, &Y_g_ceed_size, &X_g_ceed_size));
514c8564c30SJames Wright     PetscCall(MatGetLocalSize(mat_ceed, &Y_l_ceed_size, &X_l_ceed_size));
515c8564c30SJames Wright     PetscCall(MatGetSize(mat_ceed, &Y_g_other_size, &X_g_other_size));
516c8564c30SJames Wright     PetscCall(MatGetLocalSize(mat_ceed, &Y_l_other_size, &X_l_other_size));
517c8564c30SJames Wright     PetscCheck((Y_g_ceed_size == Y_g_other_size) && (X_g_ceed_size == X_g_other_size) && (Y_l_ceed_size == Y_l_other_size) &&
518c8564c30SJames Wright                    (X_l_ceed_size == X_l_other_size),
519c8564c30SJames Wright                PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ,
520c8564c30SJames Wright                "mat_ceed and mat_other must have compatible sizes; found mat_ceed (Global: %" PetscInt_FMT ", %" PetscInt_FMT
521c8564c30SJames Wright                "; Local: %" PetscInt_FMT ", %" PetscInt_FMT ") mat_other (Global: %" PetscInt_FMT ", %" PetscInt_FMT "; Local: %" PetscInt_FMT
522c8564c30SJames Wright                ", %" PetscInt_FMT ")",
523c8564c30SJames Wright                Y_g_ceed_size, X_g_ceed_size, Y_l_ceed_size, X_l_ceed_size, Y_g_other_size, X_g_other_size, Y_l_other_size, X_l_other_size);
524c8564c30SJames Wright   }
525c8564c30SJames Wright 
526c8564c30SJames Wright   // Convert
527c8564c30SJames Wright   {
528c8564c30SJames Wright     VecType        vec_type;
529c8564c30SJames Wright     MatCeedContext ctx;
530c8564c30SJames Wright 
531c8564c30SJames Wright     PetscCall(PetscObjectChangeTypeName((PetscObject)mat_other, MATCEED));
532c8564c30SJames Wright     PetscCall(MatShellGetContext(mat_ceed, &ctx));
533c8564c30SJames Wright     PetscCall(MatCeedContextReference(ctx));
534c8564c30SJames Wright     PetscCall(MatShellSetContext(mat_other, ctx));
535243afec9SJames Wright     PetscCall(MatShellSetContextDestroy(mat_other, (PetscCtxDestroyFn *)MatCeedContextDestroy));
5363933d9a0SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_VIEW, (void (*)(void))MatView_Ceed));
537c8564c30SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_MULT, (void (*)(void))MatMult_Ceed));
538c8564c30SJames Wright     if (ctx->op_mult_transpose) PetscCall(MatShellSetOperation(mat_other, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed));
539c8564c30SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed));
540c8564c30SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed));
541243afec9SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_GET_BLOCK_DIAGONAL, (void (*)(void))MatGetBlockDiagonal_Ceed));
542243afec9SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_GET_VBLOCK_DIAGONAL, (void (*)(void))MatGetVariableBlockDiagonal_Ceed));
543c8564c30SJames Wright     {
544c8564c30SJames Wright       PetscInt block_size;
545c8564c30SJames Wright 
546c8564c30SJames Wright       PetscCall(MatGetBlockSize(mat_ceed, &block_size));
547c8564c30SJames Wright       if (block_size > 1) PetscCall(MatSetBlockSize(mat_other, block_size));
548c8564c30SJames Wright     }
549c8564c30SJames Wright     {
550c8564c30SJames Wright       PetscInt        num_blocks;
551c8564c30SJames Wright       const PetscInt *block_sizes;
552c8564c30SJames Wright 
553c8564c30SJames Wright       PetscCall(MatGetVariableBlockSizes(mat_ceed, &num_blocks, &block_sizes));
554c8564c30SJames Wright       if (num_blocks) PetscCall(MatSetVariableBlockSizes(mat_other, num_blocks, (PetscInt *)block_sizes));
555c8564c30SJames Wright     }
556c8564c30SJames Wright     PetscCall(DMGetVecType(ctx->dm_x, &vec_type));
557c8564c30SJames Wright     PetscCall(MatShellSetVecType(mat_other, vec_type));
558c8564c30SJames Wright   }
559c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
560c8564c30SJames Wright }
561c8564c30SJames Wright 
562c8564c30SJames Wright /**
563243afec9SJames Wright   @brief Mark `CeedQFunction` data as updated and the `CeedQFunction` as requiring re-assembly for a `MatCEED`.
564243afec9SJames Wright 
565243afec9SJames Wright   Collective across MPI processes.
566243afec9SJames Wright 
567243afec9SJames Wright   @param[in]   mat_ceed       `MATCEED`
568243afec9SJames Wright   @param[out]  update_needed  Boolean flag indicating `CeedQFunction` update needed
569243afec9SJames Wright 
570243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
571243afec9SJames Wright **/
572243afec9SJames Wright PetscErrorCode MatCeedSetAssemblyDataUpdateNeeded(Mat mat_ceed, PetscBool update_needed) {
573243afec9SJames Wright   MatCeedContext ctx;
574243afec9SJames Wright 
575243afec9SJames Wright   PetscFunctionBeginUser;
576243afec9SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
577243afec9SJames Wright   PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult, update_needed));
578243afec9SJames Wright   if (ctx->op_mult_transpose) {
579243afec9SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult_transpose, update_needed));
580243afec9SJames Wright   }
581243afec9SJames Wright   if (update_needed) {
582243afec9SJames Wright     PetscCall(MatAssemblyBegin(mat_ceed, MAT_FINAL_ASSEMBLY));
583243afec9SJames Wright     PetscCall(MatAssemblyEnd(mat_ceed, MAT_FINAL_ASSEMBLY));
584243afec9SJames Wright   }
585243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
586243afec9SJames Wright }
587243afec9SJames Wright 
588243afec9SJames Wright /**
5895037d55dSJames Wright   @brief Setup a `Mat` with the same COO pattern as a `MatCEED`.
5905037d55dSJames Wright 
5915037d55dSJames Wright   Collective across MPI processes.
5925037d55dSJames Wright 
5935037d55dSJames Wright   @param[in]   mat_ceed  `MATCEED`
5945037d55dSJames Wright   @param[out]  mat_coo   Sparse `Mat` with same COO pattern
5955037d55dSJames Wright 
5965037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
5975037d55dSJames Wright **/
5985037d55dSJames Wright PetscErrorCode MatCeedCreateMatCOO(Mat mat_ceed, Mat *mat_coo) {
5995037d55dSJames Wright   MatCeedContext ctx;
6005037d55dSJames Wright 
6015037d55dSJames Wright   PetscFunctionBeginUser;
6025037d55dSJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
6035037d55dSJames Wright 
6045037d55dSJames Wright   PetscCheck(ctx->dm_x == ctx->dm_y, PetscObjectComm((PetscObject)mat_ceed), PETSC_ERR_SUP, "COO assembly only supported for MATCEED on a single DM");
6055037d55dSJames Wright 
6065037d55dSJames Wright   // Check cl mat type
6075037d55dSJames Wright   {
6085037d55dSJames Wright     PetscBool is_coo_mat_type_cl = PETSC_FALSE;
6095037d55dSJames Wright     char      coo_mat_type_cl[64];
6105037d55dSJames Wright 
6115037d55dSJames Wright     // Check for specific CL coo mat type for this Mat
6125037d55dSJames Wright     {
6135037d55dSJames Wright       const char *mat_ceed_prefix = NULL;
6145037d55dSJames Wright 
6155037d55dSJames Wright       PetscCall(MatGetOptionsPrefix(mat_ceed, &mat_ceed_prefix));
6165037d55dSJames Wright       PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), mat_ceed_prefix, "", NULL);
6175037d55dSJames Wright       PetscCall(PetscOptionsFList("-ceed_coo_mat_type", "Default MATCEED COO assembly MatType", NULL, MatList, coo_mat_type_cl, coo_mat_type_cl,
6185037d55dSJames Wright                                   sizeof(coo_mat_type_cl), &is_coo_mat_type_cl));
6195037d55dSJames Wright       PetscOptionsEnd();
6205037d55dSJames Wright       if (is_coo_mat_type_cl) {
6215037d55dSJames Wright         PetscCall(PetscFree(ctx->coo_mat_type));
6225037d55dSJames Wright         PetscCall(PetscStrallocpy(coo_mat_type_cl, &ctx->coo_mat_type));
6235037d55dSJames Wright       }
6245037d55dSJames Wright     }
6255037d55dSJames Wright   }
6265037d55dSJames Wright 
6275037d55dSJames Wright   // Create sparse matrix
6285037d55dSJames Wright   {
6295037d55dSJames Wright     MatType dm_mat_type, dm_mat_type_copy;
6305037d55dSJames Wright 
6315037d55dSJames Wright     PetscCall(DMGetMatType(ctx->dm_x, &dm_mat_type));
6325037d55dSJames Wright     PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy));
6335037d55dSJames Wright     PetscCall(DMSetMatType(ctx->dm_x, ctx->coo_mat_type));
6345037d55dSJames Wright     PetscCall(DMCreateMatrix(ctx->dm_x, mat_coo));
6355037d55dSJames Wright     PetscCall(DMSetMatType(ctx->dm_x, dm_mat_type_copy));
6365037d55dSJames Wright     PetscCall(PetscFree(dm_mat_type_copy));
6375037d55dSJames Wright   }
6385037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6395037d55dSJames Wright }
6405037d55dSJames Wright 
6415037d55dSJames Wright /**
6425037d55dSJames Wright   @brief Setup the COO preallocation `MATCEED` into a `MATAIJ` or similar.
643c8564c30SJames Wright          The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail.
644c8564c30SJames Wright 
645c8564c30SJames Wright   Collective across MPI processes.
646c8564c30SJames Wright 
647c8564c30SJames Wright   @param[in]      mat_ceed  `MATCEED` to assemble
648c8564c30SJames Wright   @param[in,out]  mat_coo   `MATAIJ` or similar to assemble into
649c8564c30SJames Wright 
650c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
651c8564c30SJames Wright **/
6525037d55dSJames Wright PetscErrorCode MatCeedSetPreallocationCOO(Mat mat_ceed, Mat mat_coo) {
653c8564c30SJames Wright   MatCeedContext ctx;
654c8564c30SJames Wright 
655c8564c30SJames Wright   PetscFunctionBeginUser;
656c8564c30SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
657c8564c30SJames Wright 
658c8564c30SJames Wright   {
659c8564c30SJames Wright     PetscInt     *rows_petsc = NULL, *cols_petsc = NULL;
660c8564c30SJames Wright     CeedInt      *rows_ceed, *cols_ceed;
661c8564c30SJames Wright     PetscCount    num_entries;
662c8564c30SJames Wright     PetscLogStage stage_amg_setup;
663c8564c30SJames Wright 
664c8564c30SJames Wright     // -- Assemble sparsity pattern if mat hasn't been assembled before
66543327b86SJames Wright     PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup));
666c8564c30SJames Wright     if (stage_amg_setup == -1) {
66743327b86SJames Wright       PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup));
668c8564c30SJames Wright     }
669c8564c30SJames Wright     PetscCall(PetscLogStagePush(stage_amg_setup));
67043327b86SJames Wright     PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL));
67143327b86SJames Wright     PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
67254831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed));
67343327b86SJames Wright     PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
674d0593705SJames Wright     PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc));
675d0593705SJames Wright     PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc));
676c8564c30SJames Wright     PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc));
677c8564c30SJames Wright     free(rows_petsc);
678c8564c30SJames Wright     free(cols_petsc);
67954831c5fSJames Wright     if (!ctx->coo_values_full) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_full));
680c8564c30SJames Wright     PetscCall(PetscRealloc(++ctx->num_mats_assembled_full * sizeof(Mat), &ctx->mats_assembled_full));
681c8564c30SJames Wright     ctx->mats_assembled_full[ctx->num_mats_assembled_full - 1] = mat_coo;
68243327b86SJames Wright     PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL));
683c8564c30SJames Wright     PetscCall(PetscLogStagePop());
684c8564c30SJames Wright   }
6855037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6865037d55dSJames Wright }
6875037d55dSJames Wright 
6885037d55dSJames Wright /**
6895037d55dSJames Wright   @brief Assemble a `MATCEED` into a `MATAIJ` or similar.
6905037d55dSJames Wright          The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`.
6915037d55dSJames Wright          The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail.
6925037d55dSJames Wright 
6935037d55dSJames Wright   Collective across MPI processes.
6945037d55dSJames Wright 
6955037d55dSJames Wright   @param[in]      mat_ceed  `MATCEED` to assemble
6965037d55dSJames Wright   @param[in,out]  mat_coo   `MATAIJ` or similar to assemble into
6975037d55dSJames Wright 
6985037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
6995037d55dSJames Wright **/
7005037d55dSJames Wright PetscErrorCode MatCeedAssembleCOO(Mat mat_ceed, Mat mat_coo) {
7015037d55dSJames Wright   MatCeedContext ctx;
7025037d55dSJames Wright 
7035037d55dSJames Wright   PetscFunctionBeginUser;
7045037d55dSJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
7055037d55dSJames Wright 
7065037d55dSJames Wright   // Set COO pattern if needed
7075037d55dSJames Wright   {
7085037d55dSJames Wright     CeedInt index = -1;
7095037d55dSJames Wright 
7105037d55dSJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) {
7115037d55dSJames Wright       if (ctx->mats_assembled_full[i] == mat_coo) index = i;
7125037d55dSJames Wright     }
7135037d55dSJames Wright     if (index == -1) PetscCall(MatCeedSetPreallocationCOO(mat_ceed, mat_coo));
714c8564c30SJames Wright   }
715c8564c30SJames Wright 
716c8564c30SJames Wright   // Assemble mat_ceed
71743327b86SJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL));
718c8564c30SJames Wright   PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY));
719c8564c30SJames Wright   {
720c8564c30SJames Wright     const CeedScalar *values;
721c8564c30SJames Wright     MatType           mat_type;
722c8564c30SJames Wright     CeedMemType       mem_type = CEED_MEM_HOST;
723c8564c30SJames Wright     PetscBool         is_spd, is_spd_known;
724c8564c30SJames Wright 
725c8564c30SJames Wright     PetscCall(MatGetType(mat_coo, &mat_type));
726c8564c30SJames Wright     if (strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
727c8564c30SJames Wright     else if (strstr(mat_type, "kokkos")) mem_type = CEED_MEM_DEVICE;
728c8564c30SJames Wright     else mem_type = CEED_MEM_HOST;
729c8564c30SJames Wright 
73043327b86SJames Wright     PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
73154831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemble(ctx->op_mult, ctx->coo_values_full));
73243327b86SJames Wright     PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
73354831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_full, mem_type, &values));
734c8564c30SJames Wright     PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES));
735c8564c30SJames Wright     PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd));
736c8564c30SJames Wright     if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd));
73754831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_full, &values));
738c8564c30SJames Wright   }
739c8564c30SJames Wright   PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY));
74043327b86SJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL));
741c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
742c8564c30SJames Wright }
743c8564c30SJames Wright 
744c8564c30SJames Wright /**
7455037d55dSJames Wright   @brief Set the current value of a context field for a `MatCEED`.
7465037d55dSJames Wright 
7475037d55dSJames Wright   Not collective across MPI processes.
7485037d55dSJames Wright 
7495037d55dSJames Wright   @param[in,out]  mat    `MatCEED`
7505037d55dSJames Wright   @param[in]      name   Name of the context field
7515037d55dSJames Wright   @param[in]      value  New context field value
7525037d55dSJames Wright 
7535037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
7545037d55dSJames Wright **/
7555037d55dSJames Wright PetscErrorCode MatCeedSetContextDouble(Mat mat, const char *name, double value) {
7565037d55dSJames Wright   PetscBool      was_updated = PETSC_FALSE;
7575037d55dSJames Wright   MatCeedContext ctx;
7585037d55dSJames Wright 
7595037d55dSJames Wright   PetscFunctionBeginUser;
7605037d55dSJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
7615037d55dSJames Wright   {
7625037d55dSJames Wright     CeedContextFieldLabel label = NULL;
7635037d55dSJames Wright 
7645037d55dSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult, name, &label));
7655037d55dSJames Wright     if (label) {
766243afec9SJames Wright       double set_value = 2 * value + 1.0;
767243afec9SJames Wright 
768243afec9SJames Wright       PetscCall(MatCeedGetContextDouble(mat, name, &set_value));
769243afec9SJames Wright       if (set_value != value) {
7705037d55dSJames Wright         PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult, label, &value));
7715037d55dSJames Wright         was_updated = PETSC_TRUE;
7725037d55dSJames Wright       }
773243afec9SJames Wright     }
7745037d55dSJames Wright     if (ctx->op_mult_transpose) {
7755037d55dSJames Wright       label = NULL;
7765037d55dSJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult_transpose, name, &label));
7775037d55dSJames Wright       if (label) {
778243afec9SJames Wright         double set_value = 2 * value + 1.0;
779243afec9SJames Wright 
780243afec9SJames Wright         PetscCall(MatCeedGetContextDouble(mat, name, &set_value));
781243afec9SJames Wright         if (set_value != value) {
7825037d55dSJames Wright           PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult_transpose, label, &value));
7835037d55dSJames Wright           was_updated = PETSC_TRUE;
7845037d55dSJames Wright         }
7855037d55dSJames Wright       }
7865037d55dSJames Wright     }
787243afec9SJames Wright   }
7885037d55dSJames Wright   if (was_updated) {
7895037d55dSJames Wright     PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
7905037d55dSJames Wright     PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
7915037d55dSJames Wright   }
7925037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
7935037d55dSJames Wright }
7945037d55dSJames Wright 
7955037d55dSJames Wright /**
7965037d55dSJames Wright   @brief Get the current value of a context field for a `MatCEED`.
7975037d55dSJames Wright 
7985037d55dSJames Wright   Not collective across MPI processes.
7995037d55dSJames Wright 
8005037d55dSJames Wright   @param[in]   mat    `MatCEED`
8015037d55dSJames Wright   @param[in]   name   Name of the context field
8025037d55dSJames Wright   @param[out]  value  Current context field value
8035037d55dSJames Wright 
8045037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
8055037d55dSJames Wright **/
8065037d55dSJames Wright PetscErrorCode MatCeedGetContextDouble(Mat mat, const char *name, double *value) {
8075037d55dSJames Wright   MatCeedContext ctx;
8085037d55dSJames Wright 
8095037d55dSJames Wright   PetscFunctionBeginUser;
8105037d55dSJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
8115037d55dSJames Wright   {
8125037d55dSJames Wright     CeedContextFieldLabel label = NULL;
8135037d55dSJames Wright     CeedOperator          op    = ctx->op_mult;
8145037d55dSJames Wright 
8155037d55dSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label));
8165037d55dSJames Wright     if (!label && ctx->op_mult_transpose) {
8175037d55dSJames Wright       op = ctx->op_mult_transpose;
8185037d55dSJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label));
8195037d55dSJames Wright     }
8205037d55dSJames Wright     if (label) {
8215037d55dSJames Wright       PetscSizeT    num_values;
8225037d55dSJames Wright       const double *values_ceed;
8235037d55dSJames Wright 
8245037d55dSJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorGetContextDoubleRead(op, label, &num_values, &values_ceed));
8255037d55dSJames Wright       *value = values_ceed[0];
8265037d55dSJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorRestoreContextDoubleRead(op, label, &values_ceed));
8275037d55dSJames Wright     }
8285037d55dSJames Wright   }
8295037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
8305037d55dSJames Wright }
8315037d55dSJames Wright 
8325037d55dSJames Wright /**
833243afec9SJames Wright   @brief Set the current `PetscReal` value of a context field for a `MatCEED`.
834243afec9SJames Wright 
835243afec9SJames Wright   Not collective across MPI processes.
836243afec9SJames Wright 
837243afec9SJames Wright   @param[in,out]  mat    `MatCEED`
838243afec9SJames Wright   @param[in]      name   Name of the context field
839243afec9SJames Wright   @param[in]      value  New context field value
840243afec9SJames Wright 
841243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
842243afec9SJames Wright **/
843243afec9SJames Wright PetscErrorCode MatCeedSetContextReal(Mat mat, const char *name, PetscReal value) {
844243afec9SJames Wright   double value_double = value;
845243afec9SJames Wright 
846243afec9SJames Wright   PetscFunctionBeginUser;
847243afec9SJames Wright   PetscCall(MatCeedSetContextDouble(mat, name, value_double));
848243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
849243afec9SJames Wright }
850243afec9SJames Wright 
851243afec9SJames Wright /**
852f965f5c6SJames Wright   @brief Get the current `PetscReal` value of a context field for a `MatCEED`.
853f965f5c6SJames Wright 
854f965f5c6SJames Wright   Not collective across MPI processes.
855f965f5c6SJames Wright 
856f965f5c6SJames Wright   @param[in]   mat    `MatCEED`
857f965f5c6SJames Wright   @param[in]   name   Name of the context field
858f965f5c6SJames Wright   @param[out]  value  Current context field value
859f965f5c6SJames Wright 
860f965f5c6SJames Wright   @return An error code: 0 - success, otherwise - failure
861f965f5c6SJames Wright **/
862f965f5c6SJames Wright PetscErrorCode MatCeedGetContextReal(Mat mat, const char *name, PetscReal *value) {
863243afec9SJames Wright   double value_double = 0.0;
864f965f5c6SJames Wright 
865f965f5c6SJames Wright   PetscFunctionBeginUser;
866f965f5c6SJames Wright   PetscCall(MatCeedGetContextDouble(mat, name, &value_double));
867f965f5c6SJames Wright   *value = value_double;
868f965f5c6SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
869f965f5c6SJames Wright }
870f965f5c6SJames Wright 
871f965f5c6SJames Wright /**
872243afec9SJames Wright   @brief Set the current time for a `MatCEED`.
873243afec9SJames Wright 
874243afec9SJames Wright   Not collective across MPI processes.
875243afec9SJames Wright 
876243afec9SJames Wright   @param[in,out]  mat   `MatCEED`
877243afec9SJames Wright   @param[in]      time  Current time
878243afec9SJames Wright 
879243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
880243afec9SJames Wright **/
881243afec9SJames Wright PetscErrorCode MatCeedSetTime(Mat mat, PetscReal time) {
882243afec9SJames Wright   PetscFunctionBeginUser;
883243afec9SJames Wright   {
884243afec9SJames Wright     double time_ceed = time;
885243afec9SJames Wright 
886243afec9SJames Wright     PetscCall(MatCeedSetContextDouble(mat, "time", time_ceed));
887243afec9SJames Wright   }
888243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
889243afec9SJames Wright }
890243afec9SJames Wright 
891243afec9SJames Wright /**
892243afec9SJames Wright   @brief Get the current time for a `MatCEED`.
893243afec9SJames Wright 
894243afec9SJames Wright   Not collective across MPI processes.
895243afec9SJames Wright 
896243afec9SJames Wright   @param[in]   mat   `MatCEED`
897243afec9SJames Wright   @param[out]  time  Current time, or -1.0 if the boundary evaluator has no time field
898243afec9SJames Wright 
899243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
900243afec9SJames Wright **/
901243afec9SJames Wright PetscErrorCode MatCeedGetTime(Mat mat, PetscReal *time) {
902243afec9SJames Wright   PetscFunctionBeginUser;
903243afec9SJames Wright   *time = -1.0;
904243afec9SJames Wright   {
905243afec9SJames Wright     double time_ceed = -1.0;
906243afec9SJames Wright 
907243afec9SJames Wright     PetscCall(MatCeedGetContextDouble(mat, "time", &time_ceed));
908243afec9SJames Wright     *time = time_ceed;
909243afec9SJames Wright   }
910243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
911243afec9SJames Wright }
912243afec9SJames Wright 
913243afec9SJames Wright /**
914243afec9SJames Wright   @brief Set the current time step for a `MatCEED`.
915243afec9SJames Wright 
916243afec9SJames Wright   Not collective across MPI processes.
917243afec9SJames Wright 
918243afec9SJames Wright   @param[in,out]  mat  `MatCEED`
919243afec9SJames Wright   @param[in]      dt   Current time step
920243afec9SJames Wright 
921243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
922243afec9SJames Wright **/
923243afec9SJames Wright PetscErrorCode MatCeedSetDt(Mat mat, PetscReal dt) {
924243afec9SJames Wright   PetscFunctionBeginUser;
925243afec9SJames Wright   {
926243afec9SJames Wright     double dt_ceed = dt;
927243afec9SJames Wright 
928243afec9SJames Wright     PetscCall(MatCeedSetContextDouble(mat, "dt", dt_ceed));
929243afec9SJames Wright   }
930243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
931243afec9SJames Wright }
932243afec9SJames Wright 
933243afec9SJames Wright /**
934243afec9SJames Wright   @brief Set the Jacobian shifts for a `MatCEED`.
935243afec9SJames Wright 
936243afec9SJames Wright   Not collective across MPI processes.
937243afec9SJames Wright 
938243afec9SJames Wright   @param[in,out]  mat      `MatCEED`
939243afec9SJames Wright   @param[in]      shift_v  Velocity shift
940243afec9SJames Wright   @param[in]      shift_a  Acceleration shift
941243afec9SJames Wright 
942243afec9SJames Wright   @return An error code: 0 - success, otherwise - failure
943243afec9SJames Wright **/
944243afec9SJames Wright PetscErrorCode MatCeedSetShifts(Mat mat, PetscReal shift_v, PetscReal shift_a) {
945243afec9SJames Wright   PetscFunctionBeginUser;
946243afec9SJames Wright   {
947243afec9SJames Wright     double shift_v_ceed = shift_v;
948243afec9SJames Wright 
949243afec9SJames Wright     PetscCall(MatCeedSetContextDouble(mat, "shift v", shift_v_ceed));
950243afec9SJames Wright   }
951243afec9SJames Wright   if (shift_a) {
952243afec9SJames Wright     double shift_a_ceed = shift_a;
953243afec9SJames Wright 
954243afec9SJames Wright     PetscCall(MatCeedSetContextDouble(mat, "shift a", shift_a_ceed));
955243afec9SJames Wright   }
956243afec9SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
957243afec9SJames Wright }
958243afec9SJames Wright 
959243afec9SJames Wright /**
960c8564c30SJames Wright   @brief Set user context for a `MATCEED`.
961c8564c30SJames Wright 
962c8564c30SJames Wright   Collective across MPI processes.
963c8564c30SJames Wright 
964c8564c30SJames Wright   @param[in,out]  mat  `MATCEED`
965c8564c30SJames Wright   @param[in]      f    The context destroy function, or NULL
966c8564c30SJames Wright   @param[in]      ctx  User context, or NULL to unset
967c8564c30SJames Wright 
968c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
969c8564c30SJames Wright **/
970243afec9SJames Wright PetscErrorCode MatCeedSetContext(Mat mat, PetscCtxDestroyFn f, void *ctx) {
971c8564c30SJames Wright   PetscContainer user_ctx = NULL;
972c8564c30SJames Wright 
973c8564c30SJames Wright   PetscFunctionBeginUser;
974c8564c30SJames Wright   if (ctx) {
975c8564c30SJames Wright     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)mat), &user_ctx));
976c8564c30SJames Wright     PetscCall(PetscContainerSetPointer(user_ctx, ctx));
977243afec9SJames Wright     PetscCall(PetscContainerSetCtxDestroy(user_ctx, f));
978c8564c30SJames Wright   }
979c8564c30SJames Wright   PetscCall(PetscObjectCompose((PetscObject)mat, "MatCeed user context", (PetscObject)user_ctx));
980c8564c30SJames Wright   PetscCall(PetscContainerDestroy(&user_ctx));
981c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
982c8564c30SJames Wright }
983c8564c30SJames Wright 
984c8564c30SJames Wright /**
985c8564c30SJames Wright   @brief Retrieve the user context for a `MATCEED`.
986c8564c30SJames Wright 
987c8564c30SJames Wright   Collective across MPI processes.
988c8564c30SJames Wright 
989c8564c30SJames Wright   @param[in,out]  mat  `MATCEED`
990c8564c30SJames Wright   @param[in]      ctx  User context
991c8564c30SJames Wright 
992c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
993c8564c30SJames Wright **/
994c8564c30SJames Wright PetscErrorCode MatCeedGetContext(Mat mat, void *ctx) {
995c8564c30SJames Wright   PetscContainer user_ctx;
996c8564c30SJames Wright 
997c8564c30SJames Wright   PetscFunctionBeginUser;
998c8564c30SJames Wright   PetscCall(PetscObjectQuery((PetscObject)mat, "MatCeed user context", (PetscObject *)&user_ctx));
999c8564c30SJames Wright   if (user_ctx) PetscCall(PetscContainerGetPointer(user_ctx, (void **)ctx));
1000c8564c30SJames Wright   else *(void **)ctx = NULL;
1001c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1002c8564c30SJames Wright }
1003c8564c30SJames Wright /**
10045037d55dSJames Wright   @brief Set a user defined matrix operation for a `MATCEED` matrix.
10055037d55dSJames Wright 
10065037d55dSJames Wright   Within each user-defined routine, the user should call `MatCeedGetContext()` to obtain the user-defined context that was set by
10075037d55dSJames Wright `MatCeedSetContext()`.
1008c8564c30SJames Wright 
1009c8564c30SJames Wright   Collective across MPI processes.
1010c8564c30SJames Wright 
1011c8564c30SJames Wright   @param[in,out]  mat  `MATCEED`
10125037d55dSJames Wright   @param[in]      op   Name of the `MatOperation`
10135037d55dSJames Wright   @param[in]      g    Function that provides the operation
1014c8564c30SJames Wright 
1015c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1016c8564c30SJames Wright **/
10175037d55dSJames Wright PetscErrorCode MatCeedSetOperation(Mat mat, MatOperation op, void (*g)(void)) {
10185037d55dSJames Wright   PetscFunctionBeginUser;
10195037d55dSJames Wright   PetscCall(MatShellSetOperation(mat, op, g));
10205037d55dSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
10215037d55dSJames Wright }
10225037d55dSJames Wright 
10235037d55dSJames Wright /**
10245037d55dSJames Wright   @brief Sets the default COO matrix type as a string from the `MATCEED`.
10255037d55dSJames Wright 
10265037d55dSJames Wright   Collective across MPI processes.
10275037d55dSJames Wright 
10285037d55dSJames Wright   @param[in,out]  mat   `MATCEED`
10295037d55dSJames Wright   @param[in]      type  COO `MatType` to set
10305037d55dSJames Wright 
10315037d55dSJames Wright   @return An error code: 0 - success, otherwise - failure
10325037d55dSJames Wright **/
10335037d55dSJames Wright PetscErrorCode MatCeedSetCOOMatType(Mat mat, MatType type) {
1034c8564c30SJames Wright   MatCeedContext ctx;
1035c8564c30SJames Wright 
1036c8564c30SJames Wright   PetscFunctionBeginUser;
1037c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1038c8564c30SJames Wright   // Check if same
1039c8564c30SJames Wright   {
1040c8564c30SJames Wright     size_t    len_old, len_new;
1041c8564c30SJames Wright     PetscBool is_same = PETSC_FALSE;
1042c8564c30SJames Wright 
10435037d55dSJames Wright     PetscCall(PetscStrlen(ctx->coo_mat_type, &len_old));
1044c8564c30SJames Wright     PetscCall(PetscStrlen(type, &len_new));
10455037d55dSJames Wright     if (len_old == len_new) PetscCall(PetscStrncmp(ctx->coo_mat_type, type, len_old, &is_same));
1046c8564c30SJames Wright     if (is_same) PetscFunctionReturn(PETSC_SUCCESS);
1047c8564c30SJames Wright   }
1048c8564c30SJames Wright   // Clean up old mats in different format
1049c8564c30SJames Wright   // LCOV_EXCL_START
1050c8564c30SJames Wright   if (ctx->mat_assembled_full_internal) {
1051c8564c30SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) {
1052c8564c30SJames Wright       if (ctx->mats_assembled_full[i] == ctx->mat_assembled_full_internal) {
1053c8564c30SJames Wright         for (PetscInt j = i + 1; j < ctx->num_mats_assembled_full; j++) {
1054c8564c30SJames Wright           ctx->mats_assembled_full[j - 1] = ctx->mats_assembled_full[j];
1055c8564c30SJames Wright         }
1056c8564c30SJames Wright         ctx->num_mats_assembled_full--;
1057c8564c30SJames Wright         // Note: we'll realloc this array again, so no need to shrink the allocation
1058c8564c30SJames Wright         PetscCall(MatDestroy(&ctx->mat_assembled_full_internal));
1059c8564c30SJames Wright       }
1060c8564c30SJames Wright     }
1061c8564c30SJames Wright   }
1062c8564c30SJames Wright   if (ctx->mat_assembled_pbd_internal) {
1063c8564c30SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) {
1064c8564c30SJames Wright       if (ctx->mats_assembled_pbd[i] == ctx->mat_assembled_pbd_internal) {
1065c8564c30SJames Wright         for (PetscInt j = i + 1; j < ctx->num_mats_assembled_pbd; j++) {
1066c8564c30SJames Wright           ctx->mats_assembled_pbd[j - 1] = ctx->mats_assembled_pbd[j];
1067c8564c30SJames Wright         }
1068c8564c30SJames Wright         // Note: we'll realloc this array again, so no need to shrink the allocation
1069c8564c30SJames Wright         ctx->num_mats_assembled_pbd--;
1070c8564c30SJames Wright         PetscCall(MatDestroy(&ctx->mat_assembled_pbd_internal));
1071c8564c30SJames Wright       }
1072c8564c30SJames Wright     }
1073c8564c30SJames Wright   }
10745037d55dSJames Wright   PetscCall(PetscFree(ctx->coo_mat_type));
10755037d55dSJames Wright   PetscCall(PetscStrallocpy(type, &ctx->coo_mat_type));
1076c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1077c8564c30SJames Wright   // LCOV_EXCL_STOP
1078c8564c30SJames Wright }
1079c8564c30SJames Wright 
1080c8564c30SJames Wright /**
10815037d55dSJames Wright   @brief Gets the default COO matrix type as a string from the `MATCEED`.
1082c8564c30SJames Wright 
1083c8564c30SJames Wright   Collective across MPI processes.
1084c8564c30SJames Wright 
1085c8564c30SJames Wright   @param[in,out]  mat   `MATCEED`
10865037d55dSJames Wright   @param[in]      type  COO `MatType`
1087c8564c30SJames Wright 
1088c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1089c8564c30SJames Wright **/
10905037d55dSJames Wright PetscErrorCode MatCeedGetCOOMatType(Mat mat, MatType *type) {
1091c8564c30SJames Wright   MatCeedContext ctx;
1092c8564c30SJames Wright 
1093c8564c30SJames Wright   PetscFunctionBeginUser;
1094c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
10955037d55dSJames Wright   *type = ctx->coo_mat_type;
1096c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1097c8564c30SJames Wright }
1098c8564c30SJames Wright 
1099c8564c30SJames Wright /**
1100c8564c30SJames Wright   @brief Set input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
1101c8564c30SJames Wright 
1102c8564c30SJames Wright   Not collective across MPI processes.
1103c8564c30SJames Wright 
1104c8564c30SJames Wright   @param[in,out]  mat              `MATCEED`
1105c8564c30SJames Wright   @param[in]      X_loc            Input PETSc local vector, or NULL
1106c8564c30SJames Wright   @param[in]      Y_loc_transpose  Input PETSc local vector for transpose operation, or NULL
1107c8564c30SJames Wright 
1108c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1109c8564c30SJames Wright **/
1110c8564c30SJames Wright PetscErrorCode MatCeedSetLocalVectors(Mat mat, Vec X_loc, Vec Y_loc_transpose) {
1111c8564c30SJames Wright   MatCeedContext ctx;
1112c8564c30SJames Wright 
1113c8564c30SJames Wright   PetscFunctionBeginUser;
1114c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1115c8564c30SJames Wright   if (X_loc) {
1116c8564c30SJames Wright     PetscInt len_old, len_new;
1117c8564c30SJames Wright 
1118c8564c30SJames Wright     PetscCall(VecGetSize(ctx->X_loc, &len_old));
1119c8564c30SJames Wright     PetscCall(VecGetSize(X_loc, &len_new));
1120c8564c30SJames Wright     PetscCheck(len_old == len_new, PETSC_COMM_SELF, PETSC_ERR_LIB, "new X_loc length %" PetscInt_FMT " should match old X_loc length %" PetscInt_FMT,
1121c8564c30SJames Wright                len_new, len_old);
11225037d55dSJames Wright     PetscCall(VecReferenceCopy(X_loc, &ctx->X_loc));
1123c8564c30SJames Wright   }
1124c8564c30SJames Wright   if (Y_loc_transpose) {
1125c8564c30SJames Wright     PetscInt len_old, len_new;
1126c8564c30SJames Wright 
1127c8564c30SJames Wright     PetscCall(VecGetSize(ctx->Y_loc_transpose, &len_old));
1128c8564c30SJames Wright     PetscCall(VecGetSize(Y_loc_transpose, &len_new));
1129c8564c30SJames Wright     PetscCheck(len_old == len_new, PETSC_COMM_SELF, PETSC_ERR_LIB,
1130c8564c30SJames Wright                "new Y_loc_transpose length %" PetscInt_FMT " should match old Y_loc_transpose length %" PetscInt_FMT, len_new, len_old);
11315037d55dSJames Wright     PetscCall(VecReferenceCopy(Y_loc_transpose, &ctx->Y_loc_transpose));
1132c8564c30SJames Wright   }
1133c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1134c8564c30SJames Wright }
1135c8564c30SJames Wright 
1136c8564c30SJames Wright /**
1137c8564c30SJames Wright   @brief Get input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
1138c8564c30SJames Wright 
1139c8564c30SJames Wright   Not collective across MPI processes.
1140c8564c30SJames Wright 
1141c8564c30SJames Wright   @param[in,out]  mat              `MATCEED`
1142c8564c30SJames Wright   @param[out]     X_loc            Input PETSc local vector, or NULL
1143c8564c30SJames Wright   @param[out]     Y_loc_transpose  Input PETSc local vector for transpose operation, or NULL
1144c8564c30SJames Wright 
1145c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1146c8564c30SJames Wright **/
1147c8564c30SJames Wright PetscErrorCode MatCeedGetLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) {
1148c8564c30SJames Wright   MatCeedContext ctx;
1149c8564c30SJames Wright 
1150c8564c30SJames Wright   PetscFunctionBeginUser;
1151c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1152c8564c30SJames Wright   if (X_loc) {
11535037d55dSJames Wright     *X_loc = NULL;
11545037d55dSJames Wright     PetscCall(VecReferenceCopy(ctx->X_loc, X_loc));
1155c8564c30SJames Wright   }
1156c8564c30SJames Wright   if (Y_loc_transpose) {
11575037d55dSJames Wright     *Y_loc_transpose = NULL;
11585037d55dSJames Wright     PetscCall(VecReferenceCopy(ctx->Y_loc_transpose, Y_loc_transpose));
1159c8564c30SJames Wright   }
1160c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1161c8564c30SJames Wright }
1162c8564c30SJames Wright 
1163c8564c30SJames Wright /**
1164c8564c30SJames Wright   @brief Restore input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
1165c8564c30SJames Wright 
1166c8564c30SJames Wright   Not collective across MPI processes.
1167c8564c30SJames Wright 
1168c8564c30SJames Wright   @param[in,out]  mat              MatCeed
1169c8564c30SJames Wright   @param[out]     X_loc            Input PETSc local vector, or NULL
1170c8564c30SJames Wright   @param[out]     Y_loc_transpose  Input PETSc local vector for transpose operation, or NULL
1171c8564c30SJames Wright 
1172c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1173c8564c30SJames Wright **/
1174c8564c30SJames Wright PetscErrorCode MatCeedRestoreLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) {
1175c8564c30SJames Wright   PetscFunctionBeginUser;
1176c8564c30SJames Wright   if (X_loc) PetscCall(VecDestroy(X_loc));
1177c8564c30SJames Wright   if (Y_loc_transpose) PetscCall(VecDestroy(Y_loc_transpose));
1178c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1179c8564c30SJames Wright }
1180c8564c30SJames Wright 
1181c8564c30SJames Wright /**
1182c8564c30SJames Wright   @brief Get libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
1183c8564c30SJames Wright 
1184c8564c30SJames Wright   Not collective across MPI processes.
1185c8564c30SJames Wright 
1186c8564c30SJames Wright   @param[in,out]  mat                MatCeed
1187c8564c30SJames Wright   @param[out]     op_mult            libCEED `CeedOperator` for `MatMult()`, or NULL
1188c8564c30SJames Wright   @param[out]     op_mult_transpose  libCEED `CeedOperator` for `MatMultTranspose()`, or NULL
1189c8564c30SJames Wright 
1190c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1191c8564c30SJames Wright **/
1192c8564c30SJames Wright PetscErrorCode MatCeedGetCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) {
1193c8564c30SJames Wright   MatCeedContext ctx;
1194c8564c30SJames Wright 
1195c8564c30SJames Wright   PetscFunctionBeginUser;
1196c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1197c8564c30SJames Wright   if (op_mult) {
1198c8564c30SJames Wright     *op_mult = NULL;
119954831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult, op_mult));
1200c8564c30SJames Wright   }
1201c8564c30SJames Wright   if (op_mult_transpose) {
1202c8564c30SJames Wright     *op_mult_transpose = NULL;
120354831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult_transpose, op_mult_transpose));
1204c8564c30SJames Wright   }
1205c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1206c8564c30SJames Wright }
1207c8564c30SJames Wright 
1208c8564c30SJames Wright /**
1209c8564c30SJames Wright   @brief Restore libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
1210c8564c30SJames Wright 
1211c8564c30SJames Wright   Not collective across MPI processes.
1212c8564c30SJames Wright 
1213c8564c30SJames Wright   @param[in,out]  mat                MatCeed
1214c8564c30SJames Wright   @param[out]     op_mult            libCEED `CeedOperator` for `MatMult()`, or NULL
1215c8564c30SJames Wright   @param[out]     op_mult_transpose  libCEED `CeedOperator` for `MatMultTranspose()`, or NULL
1216c8564c30SJames Wright 
1217c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1218c8564c30SJames Wright **/
1219c8564c30SJames Wright PetscErrorCode MatCeedRestoreCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) {
1220c8564c30SJames Wright   MatCeedContext ctx;
1221c8564c30SJames Wright 
1222c8564c30SJames Wright   PetscFunctionBeginUser;
1223c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
122454831c5fSJames Wright   if (op_mult) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult));
122554831c5fSJames Wright   if (op_mult_transpose) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult_transpose));
1226c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1227c8564c30SJames Wright }
1228c8564c30SJames Wright 
1229c8564c30SJames Wright /**
1230c8564c30SJames Wright   @brief Set `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
1231c8564c30SJames Wright 
1232c8564c30SJames Wright   Not collective across MPI processes.
1233c8564c30SJames Wright 
1234c8564c30SJames Wright   @param[in,out]  mat                       MatCeed
1235c8564c30SJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward evaluation, or NULL
1236c8564c30SJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose evaluation, or NULL
1237c8564c30SJames Wright 
1238c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1239c8564c30SJames Wright **/
1240c8564c30SJames Wright PetscErrorCode MatCeedSetLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) {
1241c8564c30SJames Wright   MatCeedContext ctx;
1242c8564c30SJames Wright 
1243c8564c30SJames Wright   PetscFunctionBeginUser;
1244c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1245c8564c30SJames Wright   if (log_event_mult) ctx->log_event_mult = log_event_mult;
1246c8564c30SJames Wright   if (log_event_mult_transpose) ctx->log_event_mult_transpose = log_event_mult_transpose;
1247c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1248c8564c30SJames Wright }
1249c8564c30SJames Wright 
1250c8564c30SJames Wright /**
1251c8564c30SJames Wright   @brief Get `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
1252c8564c30SJames Wright 
1253c8564c30SJames Wright   Not collective across MPI processes.
1254c8564c30SJames Wright 
1255c8564c30SJames Wright   @param[in,out]  mat                       MatCeed
1256c8564c30SJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward evaluation, or NULL
1257c8564c30SJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose evaluation, or NULL
1258c8564c30SJames Wright 
1259c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1260c8564c30SJames Wright **/
1261c8564c30SJames Wright PetscErrorCode MatCeedGetLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) {
1262c8564c30SJames Wright   MatCeedContext ctx;
1263c8564c30SJames Wright 
1264c8564c30SJames Wright   PetscFunctionBeginUser;
1265c8564c30SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1266c8564c30SJames Wright   if (log_event_mult) *log_event_mult = ctx->log_event_mult;
1267c8564c30SJames Wright   if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_mult_transpose;
1268c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1269c8564c30SJames Wright }
1270c8564c30SJames Wright 
127143327b86SJames Wright /**
127243327b86SJames Wright   @brief Set `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
127343327b86SJames Wright 
127443327b86SJames Wright   Not collective across MPI processes.
127543327b86SJames Wright 
127643327b86SJames Wright   @param[in,out]  mat                       MatCeed
127743327b86SJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL
127843327b86SJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL
127943327b86SJames Wright 
128043327b86SJames Wright   @return An error code: 0 - success, otherwise - failure
128143327b86SJames Wright **/
128243327b86SJames Wright PetscErrorCode MatCeedSetCeedOperatorLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) {
128343327b86SJames Wright   MatCeedContext ctx;
128443327b86SJames Wright 
128543327b86SJames Wright   PetscFunctionBeginUser;
128643327b86SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
128743327b86SJames Wright   if (log_event_mult) ctx->log_event_ceed_mult = log_event_mult;
128843327b86SJames Wright   if (log_event_mult_transpose) ctx->log_event_ceed_mult_transpose = log_event_mult_transpose;
128943327b86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
129043327b86SJames Wright }
129143327b86SJames Wright 
129243327b86SJames Wright /**
129343327b86SJames Wright   @brief Get `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
129443327b86SJames Wright 
129543327b86SJames Wright   Not collective across MPI processes.
129643327b86SJames Wright 
129743327b86SJames Wright   @param[in,out]  mat                       MatCeed
129843327b86SJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL
129943327b86SJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL
130043327b86SJames Wright 
130143327b86SJames Wright   @return An error code: 0 - success, otherwise - failure
130243327b86SJames Wright **/
130343327b86SJames Wright PetscErrorCode MatCeedGetCeedOperatorLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) {
130443327b86SJames Wright   MatCeedContext ctx;
130543327b86SJames Wright 
130643327b86SJames Wright   PetscFunctionBeginUser;
130743327b86SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
130843327b86SJames Wright   if (log_event_mult) *log_event_mult = ctx->log_event_ceed_mult;
130943327b86SJames Wright   if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_ceed_mult_transpose;
131043327b86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
131143327b86SJames Wright }
131243327b86SJames Wright 
1313c8564c30SJames Wright // -----------------------------------------------------------------------------
1314c8564c30SJames Wright // Operator context data
1315c8564c30SJames Wright // -----------------------------------------------------------------------------
1316c8564c30SJames Wright 
1317c8564c30SJames Wright /**
1318c8564c30SJames Wright   @brief Setup context data for operator application.
1319c8564c30SJames Wright 
1320c8564c30SJames Wright   Collective across MPI processes.
1321c8564c30SJames Wright 
1322c8564c30SJames Wright   @param[in]   dm_x                           Input `DM`
1323c8564c30SJames Wright   @param[in]   dm_y                           Output `DM`
1324c8564c30SJames Wright   @param[in]   X_loc                          Input PETSc local vector, or NULL
1325c8564c30SJames Wright   @param[in]   Y_loc_transpose                Input PETSc local vector for transpose operation, or NULL
1326c8564c30SJames Wright   @param[in]   op_mult                        `CeedOperator` for forward evaluation
1327c8564c30SJames Wright   @param[in]   op_mult_transpose              `CeedOperator` for transpose evaluation
1328c8564c30SJames Wright   @param[in]   log_event_mult                 `PetscLogEvent` for forward evaluation
1329c8564c30SJames Wright   @param[in]   log_event_mult_transpose       `PetscLogEvent` for transpose evaluation
133043327b86SJames Wright   @param[in]   log_event_ceed_mult            `PetscLogEvent` for forward `CeedOperator` evaluation
133143327b86SJames Wright   @param[in]   log_event_ceed_mult_transpose  `PetscLogEvent` for transpose `CeedOperator` evaluation
1332c8564c30SJames Wright   @param[out]  ctx                            Context data for operator evaluation
1333c8564c30SJames Wright 
1334c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1335c8564c30SJames Wright **/
1336c8564c30SJames Wright PetscErrorCode MatCeedContextCreate(DM dm_x, DM dm_y, Vec X_loc, Vec Y_loc_transpose, CeedOperator op_mult, CeedOperator op_mult_transpose,
133743327b86SJames Wright                                     PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose, PetscLogEvent log_event_ceed_mult,
133843327b86SJames Wright                                     PetscLogEvent log_event_ceed_mult_transpose, MatCeedContext *ctx) {
1339c8564c30SJames Wright   CeedSize x_loc_len, y_loc_len;
1340c8564c30SJames Wright 
1341c8564c30SJames Wright   PetscFunctionBeginUser;
1342c8564c30SJames Wright 
1343c8564c30SJames Wright   // Allocate
1344c8564c30SJames Wright   PetscCall(PetscNew(ctx));
1345c8564c30SJames Wright   (*ctx)->ref_count = 1;
1346c8564c30SJames Wright 
1347c8564c30SJames Wright   // Logging
1348c8564c30SJames Wright   (*ctx)->log_event_mult                = log_event_mult;
1349c8564c30SJames Wright   (*ctx)->log_event_mult_transpose      = log_event_mult_transpose;
135043327b86SJames Wright   (*ctx)->log_event_ceed_mult           = log_event_ceed_mult;
135143327b86SJames Wright   (*ctx)->log_event_ceed_mult_transpose = log_event_ceed_mult_transpose;
1352c8564c30SJames Wright 
1353c8564c30SJames Wright   // PETSc objects
13545037d55dSJames Wright   PetscCall(DMReferenceCopy(dm_x, &(*ctx)->dm_x));
13555037d55dSJames Wright   PetscCall(DMReferenceCopy(dm_y, &(*ctx)->dm_y));
13565037d55dSJames Wright   if (X_loc) PetscCall(VecReferenceCopy(X_loc, &(*ctx)->X_loc));
13575037d55dSJames Wright   if (Y_loc_transpose) PetscCall(VecReferenceCopy(Y_loc_transpose, &(*ctx)->Y_loc_transpose));
1358c8564c30SJames Wright 
1359c8564c30SJames Wright   // Memtype
1360c8564c30SJames Wright   {
1361c8564c30SJames Wright     const PetscScalar *x;
1362c8564c30SJames Wright     Vec                X;
1363c8564c30SJames Wright 
1364c8564c30SJames Wright     PetscCall(DMGetLocalVector(dm_x, &X));
1365c8564c30SJames Wright     PetscCall(VecGetArrayReadAndMemType(X, &x, &(*ctx)->mem_type));
1366c8564c30SJames Wright     PetscCall(VecRestoreArrayReadAndMemType(X, &x));
1367c8564c30SJames Wright     PetscCall(DMRestoreLocalVector(dm_x, &X));
1368c8564c30SJames Wright   }
1369c8564c30SJames Wright 
1370c8564c30SJames Wright   // libCEED objects
1371c8564c30SJames Wright   PetscCheck(CeedOperatorGetCeed(op_mult, &(*ctx)->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB,
1372c8564c30SJames Wright              "retrieving Ceed context object failed");
137354831c5fSJames Wright   PetscCallCeed((*ctx)->ceed, CeedOperatorGetActiveVectorLengths(op_mult, &x_loc_len, &y_loc_len));
137454831c5fSJames Wright   PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult, &(*ctx)->op_mult));
137554831c5fSJames Wright   if (op_mult_transpose) PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult_transpose, &(*ctx)->op_mult_transpose));
137654831c5fSJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, x_loc_len, &(*ctx)->x_loc));
137754831c5fSJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, y_loc_len, &(*ctx)->y_loc));
1378c8564c30SJames Wright 
1379c8564c30SJames Wright   // Flop counting
1380c8564c30SJames Wright   {
1381c8564c30SJames Wright     CeedSize ceed_flops_estimate = 0;
1382c8564c30SJames Wright 
138354831c5fSJames Wright     PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult, &ceed_flops_estimate));
1384c8564c30SJames Wright     (*ctx)->flops_mult = ceed_flops_estimate;
1385c8564c30SJames Wright     if (op_mult_transpose) {
138654831c5fSJames Wright       PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult_transpose, &ceed_flops_estimate));
1387c8564c30SJames Wright       (*ctx)->flops_mult_transpose = ceed_flops_estimate;
1388c8564c30SJames Wright     }
1389c8564c30SJames Wright   }
1390c8564c30SJames Wright 
1391c8564c30SJames Wright   // Check sizes
1392c8564c30SJames Wright   if (x_loc_len > 0 || y_loc_len > 0) {
1393c8564c30SJames Wright     CeedSize ctx_x_loc_len, ctx_y_loc_len;
1394c8564c30SJames Wright     PetscInt X_loc_len, dm_x_loc_len, Y_loc_len, dm_y_loc_len;
1395c8564c30SJames Wright     Vec      dm_X_loc, dm_Y_loc;
1396c8564c30SJames Wright 
1397c8564c30SJames Wright     // -- Input
1398c8564c30SJames Wright     PetscCall(DMGetLocalVector(dm_x, &dm_X_loc));
1399c8564c30SJames Wright     PetscCall(VecGetLocalSize(dm_X_loc, &dm_x_loc_len));
1400c8564c30SJames Wright     PetscCall(DMRestoreLocalVector(dm_x, &dm_X_loc));
140154831c5fSJames Wright     PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->x_loc, &ctx_x_loc_len));
1402c86b9822SJames Wright     if (X_loc) {
1403c86b9822SJames Wright       PetscCall(VecGetLocalSize(X_loc, &X_loc_len));
1404c86b9822SJames Wright       PetscCheck(X_loc_len == dm_x_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB,
1405c86b9822SJames Wright                  "X_loc (%" PetscInt_FMT ") must match dm_x (%" PetscInt_FMT ") dimensions", X_loc_len, dm_x_loc_len);
1406c86b9822SJames Wright     }
1407c86b9822SJames Wright     PetscCheck(x_loc_len == dm_x_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, "op (%" CeedSize_FMT ") must match dm_x (%" PetscInt_FMT ") dimensions",
1408c86b9822SJames Wright                x_loc_len, dm_x_loc_len);
1409c86b9822SJames Wright     PetscCheck(x_loc_len == ctx_x_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, "x_loc (%" CeedSize_FMT ") must match op dimensions (%" CeedSize_FMT ")",
1410c86b9822SJames Wright                x_loc_len, ctx_x_loc_len);
1411c8564c30SJames Wright 
1412c8564c30SJames Wright     // -- Output
1413c8564c30SJames Wright     PetscCall(DMGetLocalVector(dm_y, &dm_Y_loc));
1414c8564c30SJames Wright     PetscCall(VecGetLocalSize(dm_Y_loc, &dm_y_loc_len));
1415c8564c30SJames Wright     PetscCall(DMRestoreLocalVector(dm_y, &dm_Y_loc));
141654831c5fSJames Wright     PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->y_loc, &ctx_y_loc_len));
1417c86b9822SJames Wright     PetscCheck(ctx_y_loc_len == dm_y_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, "op (%" CeedSize_FMT ") must match dm_y (%" PetscInt_FMT ") dimensions",
1418c86b9822SJames Wright                ctx_y_loc_len, dm_y_loc_len);
1419c8564c30SJames Wright 
1420c8564c30SJames Wright     // -- Transpose
1421c8564c30SJames Wright     if (Y_loc_transpose) {
1422c8564c30SJames Wright       PetscCall(VecGetLocalSize(Y_loc_transpose, &Y_loc_len));
1423c86b9822SJames Wright       PetscCheck(Y_loc_len == dm_y_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB,
1424c86b9822SJames Wright                  "Y_loc_transpose (%" PetscInt_FMT ") must match dm_y (%" PetscInt_FMT ") dimensions", Y_loc_len, dm_y_loc_len);
1425c8564c30SJames Wright     }
1426c8564c30SJames Wright   }
1427c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1428c8564c30SJames Wright }
1429c8564c30SJames Wright 
1430c8564c30SJames Wright /**
1431c8564c30SJames Wright   @brief Increment reference counter for `MATCEED` context.
1432c8564c30SJames Wright 
1433c8564c30SJames Wright   Not collective across MPI processes.
1434c8564c30SJames Wright 
1435c8564c30SJames Wright   @param[in,out]  ctx  Context data
1436c8564c30SJames Wright 
1437c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1438c8564c30SJames Wright **/
1439c8564c30SJames Wright PetscErrorCode MatCeedContextReference(MatCeedContext ctx) {
1440c8564c30SJames Wright   PetscFunctionBeginUser;
1441c8564c30SJames Wright   ctx->ref_count++;
1442c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1443c8564c30SJames Wright }
1444c8564c30SJames Wright 
1445c8564c30SJames Wright /**
1446c8564c30SJames Wright   @brief Copy reference for `MATCEED`.
1447c8564c30SJames Wright          Note: If `ctx_copy` is non-null, it is assumed to be a valid pointer to a `MatCeedContext`.
1448c8564c30SJames Wright 
1449c8564c30SJames Wright   Not collective across MPI processes.
1450c8564c30SJames Wright 
1451c8564c30SJames Wright   @param[in]   ctx       Context data
1452c8564c30SJames Wright   @param[out]  ctx_copy  Copy of pointer to context data
1453c8564c30SJames Wright 
1454c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1455c8564c30SJames Wright **/
1456c8564c30SJames Wright PetscErrorCode MatCeedContextReferenceCopy(MatCeedContext ctx, MatCeedContext *ctx_copy) {
1457c8564c30SJames Wright   PetscFunctionBeginUser;
1458c8564c30SJames Wright   PetscCall(MatCeedContextReference(ctx));
1459243afec9SJames Wright   PetscCall(MatCeedContextDestroy(ctx_copy));
1460c8564c30SJames Wright   *ctx_copy = ctx;
1461c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1462c8564c30SJames Wright }
1463c8564c30SJames Wright 
1464c8564c30SJames Wright /**
1465c8564c30SJames Wright   @brief Destroy context data for operator application.
1466c8564c30SJames Wright 
1467c8564c30SJames Wright   Collective across MPI processes.
1468c8564c30SJames Wright 
1469c8564c30SJames Wright   @param[in,out]  ctx  Context data for operator evaluation
1470c8564c30SJames Wright 
1471c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1472c8564c30SJames Wright **/
1473243afec9SJames Wright PetscErrorCode MatCeedContextDestroy(MatCeedContext *ctx) {
1474c8564c30SJames Wright   PetscFunctionBeginUser;
1475243afec9SJames Wright   if (!ctx || --(*ctx)->ref_count > 0) PetscFunctionReturn(PETSC_SUCCESS);
1476c8564c30SJames Wright 
1477c8564c30SJames Wright   // PETSc objects
1478243afec9SJames Wright   PetscCall(DMDestroy(&(*ctx)->dm_x));
1479243afec9SJames Wright   PetscCall(DMDestroy(&(*ctx)->dm_y));
1480243afec9SJames Wright   PetscCall(VecDestroy(&(*ctx)->X_loc));
1481243afec9SJames Wright   PetscCall(VecDestroy(&(*ctx)->Y_loc_transpose));
1482243afec9SJames Wright   PetscCall(MatDestroy(&(*ctx)->mat_assembled_full_internal));
1483243afec9SJames Wright   PetscCall(MatDestroy(&(*ctx)->mat_assembled_pbd_internal));
1484243afec9SJames Wright   PetscCall(PetscFree((*ctx)->coo_mat_type));
1485243afec9SJames Wright   PetscCall(PetscFree((*ctx)->mats_assembled_full));
1486243afec9SJames Wright   PetscCall(PetscFree((*ctx)->mats_assembled_pbd));
1487c8564c30SJames Wright 
1488c8564c30SJames Wright   // libCEED objects
1489243afec9SJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->x_loc));
1490243afec9SJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->y_loc));
1491243afec9SJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->coo_values_full));
1492243afec9SJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->coo_values_pbd));
1493243afec9SJames Wright   PetscCallCeed((*ctx)->ceed, CeedOperatorDestroy(&(*ctx)->op_mult));
1494243afec9SJames Wright   PetscCallCeed((*ctx)->ceed, CeedOperatorDestroy(&(*ctx)->op_mult_transpose));
1495243afec9SJames Wright   PetscCheck(CeedDestroy(&(*ctx)->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "destroying libCEED context object failed");
1496c8564c30SJames Wright 
1497c8564c30SJames Wright   // Deallocate
1498243afec9SJames Wright   (*ctx)->is_destroyed = PETSC_TRUE;  // Flag as destroyed in case someone has stale ref
1499243afec9SJames Wright   PetscCall(PetscFree(*ctx));
1500c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1501c8564c30SJames Wright }
1502c8564c30SJames Wright 
1503c8564c30SJames Wright /**
1504c8564c30SJames Wright   @brief Compute the diagonal of an operator via libCEED.
1505c8564c30SJames Wright 
1506c8564c30SJames Wright   Collective across MPI processes.
1507c8564c30SJames Wright 
1508c8564c30SJames Wright   @param[in]   A  `MATCEED`
1509c8564c30SJames Wright   @param[out]  D  Vector holding operator diagonal
1510c8564c30SJames Wright 
1511c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1512c8564c30SJames Wright **/
1513c8564c30SJames Wright PetscErrorCode MatGetDiagonal_Ceed(Mat A, Vec D) {
1514c8564c30SJames Wright   PetscMemType   mem_type;
1515c8564c30SJames Wright   Vec            D_loc;
1516c8564c30SJames Wright   MatCeedContext ctx;
1517c8564c30SJames Wright 
1518c8564c30SJames Wright   PetscFunctionBeginUser;
1519c8564c30SJames Wright   PetscCall(MatShellGetContext(A, &ctx));
1520c8564c30SJames Wright 
1521c8564c30SJames Wright   // Place PETSc vector in libCEED vector
152243327b86SJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL));
1523c8564c30SJames Wright   PetscCall(DMGetLocalVector(ctx->dm_x, &D_loc));
1524d0593705SJames Wright   PetscCall(VecPetscToCeed(D_loc, &mem_type, ctx->x_loc));
1525c8564c30SJames Wright 
1526c8564c30SJames Wright   // Compute Diagonal
152743327b86SJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL));
152854831c5fSJames Wright   PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleDiagonal(ctx->op_mult, ctx->x_loc, CEED_REQUEST_IMMEDIATE));
152943327b86SJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL));
1530c8564c30SJames Wright 
1531c8564c30SJames Wright   // Restore PETSc vector
1532d0593705SJames Wright   PetscCall(VecCeedToPetsc(ctx->x_loc, mem_type, D_loc));
1533c8564c30SJames Wright 
1534c8564c30SJames Wright   // Local-to-Global
1535c8564c30SJames Wright   PetscCall(VecZeroEntries(D));
1536c8564c30SJames Wright   PetscCall(DMLocalToGlobal(ctx->dm_x, D_loc, ADD_VALUES, D));
1537c8564c30SJames Wright   PetscCall(DMRestoreLocalVector(ctx->dm_x, &D_loc));
153843327b86SJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL));
1539c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1540c8564c30SJames Wright }
1541c8564c30SJames Wright 
1542c8564c30SJames Wright /**
1543c8564c30SJames Wright   @brief Compute `A X = Y` for a `MATCEED`.
1544c8564c30SJames Wright 
1545c8564c30SJames Wright   Collective across MPI processes.
1546c8564c30SJames Wright 
1547c8564c30SJames Wright   @param[in]   A  `MATCEED`
1548c8564c30SJames Wright   @param[in]   X  Input PETSc vector
1549c8564c30SJames Wright   @param[out]  Y  Output PETSc vector
1550c8564c30SJames Wright 
1551c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1552c8564c30SJames Wright **/
1553c8564c30SJames Wright PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) {
1554c8564c30SJames Wright   MatCeedContext ctx;
1555c8564c30SJames Wright 
1556c8564c30SJames Wright   PetscFunctionBeginUser;
1557c8564c30SJames Wright   PetscCall(MatShellGetContext(A, &ctx));
155843327b86SJames Wright   PetscCall(PetscLogEventBegin(ctx->log_event_mult, A, X, Y, NULL));
1559c8564c30SJames Wright 
1560c8564c30SJames Wright   {
1561c8564c30SJames Wright     PetscMemType x_mem_type, y_mem_type;
1562c8564c30SJames Wright     Vec          X_loc = ctx->X_loc, Y_loc;
1563c8564c30SJames Wright 
1564c8564c30SJames Wright     // Get local vectors
1565c8564c30SJames Wright     if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
1566c8564c30SJames Wright     PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
1567c8564c30SJames Wright 
1568c8564c30SJames Wright     // Global-to-local
1569c8564c30SJames Wright     PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc));
1570c8564c30SJames Wright 
1571c8564c30SJames Wright     // Setup libCEED vectors
1572d0593705SJames Wright     PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, ctx->x_loc));
1573c8564c30SJames Wright     PetscCall(VecZeroEntries(Y_loc));
1574d0593705SJames Wright     PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc));
1575c8564c30SJames Wright 
1576c8564c30SJames Wright     // Apply libCEED operator
157743327b86SJames Wright     PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult, A, X, Y, NULL));
1578243afec9SJames Wright     PetscCall(PetscLogGpuTimeBegin());
157954831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult, ctx->x_loc, ctx->y_loc, CEED_REQUEST_IMMEDIATE));
1580c8564c30SJames Wright     PetscCall(PetscLogGpuTimeEnd());
1581243afec9SJames Wright     PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult, A, X, Y, NULL));
1582c8564c30SJames Wright 
1583c8564c30SJames Wright     // Restore PETSc vectors
1584d0593705SJames Wright     PetscCall(VecReadCeedToPetsc(ctx->x_loc, x_mem_type, X_loc));
1585d0593705SJames Wright     PetscCall(VecCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc));
1586c8564c30SJames Wright 
1587c8564c30SJames Wright     // Local-to-global
1588c8564c30SJames Wright     PetscCall(VecZeroEntries(Y));
1589c8564c30SJames Wright     PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y));
1590c8564c30SJames Wright 
1591c8564c30SJames Wright     // Restore local vectors, as needed
1592c8564c30SJames Wright     if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
1593c8564c30SJames Wright     PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
1594c8564c30SJames Wright   }
1595c8564c30SJames Wright 
1596c8564c30SJames Wright   // Log flops
1597c8564c30SJames Wright   if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult));
1598c8564c30SJames Wright   else PetscCall(PetscLogFlops(ctx->flops_mult));
159943327b86SJames Wright   PetscCall(PetscLogEventEnd(ctx->log_event_mult, A, X, Y, NULL));
1600c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1601c8564c30SJames Wright }
1602c8564c30SJames Wright 
1603c8564c30SJames Wright /**
1604c8564c30SJames Wright   @brief Compute `A^T Y = X` for a `MATCEED`.
1605c8564c30SJames Wright 
1606c8564c30SJames Wright   Collective across MPI processes.
1607c8564c30SJames Wright 
1608c8564c30SJames Wright   @param[in]   A  `MATCEED`
1609c8564c30SJames Wright   @param[in]   Y  Input PETSc vector
1610c8564c30SJames Wright   @param[out]  X  Output PETSc vector
1611c8564c30SJames Wright 
1612c8564c30SJames Wright   @return An error code: 0 - success, otherwise - failure
1613c8564c30SJames Wright **/
1614c8564c30SJames Wright PetscErrorCode MatMultTranspose_Ceed(Mat A, Vec Y, Vec X) {
1615c8564c30SJames Wright   MatCeedContext ctx;
1616c8564c30SJames Wright 
1617c8564c30SJames Wright   PetscFunctionBeginUser;
1618c8564c30SJames Wright   PetscCall(MatShellGetContext(A, &ctx));
161943327b86SJames Wright   PetscCall(PetscLogEventBegin(ctx->log_event_mult_transpose, A, Y, X, NULL));
1620c8564c30SJames Wright 
1621c8564c30SJames Wright   {
1622c8564c30SJames Wright     PetscMemType x_mem_type, y_mem_type;
1623c8564c30SJames Wright     Vec          X_loc, Y_loc = ctx->Y_loc_transpose;
1624c8564c30SJames Wright 
1625c8564c30SJames Wright     // Get local vectors
1626c8564c30SJames Wright     if (!ctx->Y_loc_transpose) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
1627c8564c30SJames Wright     PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
1628c8564c30SJames Wright 
1629c8564c30SJames Wright     // Global-to-local
1630c8564c30SJames Wright     PetscCall(DMGlobalToLocal(ctx->dm_y, Y, INSERT_VALUES, Y_loc));
1631c8564c30SJames Wright 
1632c8564c30SJames Wright     // Setup libCEED vectors
1633d0593705SJames Wright     PetscCall(VecReadPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc));
1634c8564c30SJames Wright     PetscCall(VecZeroEntries(X_loc));
1635d0593705SJames Wright     PetscCall(VecPetscToCeed(X_loc, &x_mem_type, ctx->x_loc));
1636c8564c30SJames Wright 
1637c8564c30SJames Wright     // Apply libCEED operator
163843327b86SJames Wright     PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL));
1639243afec9SJames Wright     PetscCall(PetscLogGpuTimeBegin());
164054831c5fSJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult_transpose, ctx->y_loc, ctx->x_loc, CEED_REQUEST_IMMEDIATE));
1641c8564c30SJames Wright     PetscCall(PetscLogGpuTimeEnd());
1642243afec9SJames Wright     PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL));
1643c8564c30SJames Wright 
1644c8564c30SJames Wright     // Restore PETSc vectors
1645d0593705SJames Wright     PetscCall(VecReadCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc));
1646d0593705SJames Wright     PetscCall(VecCeedToPetsc(ctx->x_loc, x_mem_type, X_loc));
1647c8564c30SJames Wright 
1648c8564c30SJames Wright     // Local-to-global
1649c8564c30SJames Wright     PetscCall(VecZeroEntries(X));
1650c8564c30SJames Wright     PetscCall(DMLocalToGlobal(ctx->dm_x, X_loc, ADD_VALUES, X));
1651c8564c30SJames Wright 
1652c8564c30SJames Wright     // Restore local vectors, as needed
1653c8564c30SJames Wright     if (!ctx->Y_loc_transpose) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
1654c8564c30SJames Wright     PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
1655c8564c30SJames Wright   }
1656c8564c30SJames Wright 
1657c8564c30SJames Wright   // Log flops
1658c8564c30SJames Wright   if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult_transpose));
1659c8564c30SJames Wright   else PetscCall(PetscLogFlops(ctx->flops_mult_transpose));
166043327b86SJames Wright   PetscCall(PetscLogEventEnd(ctx->log_event_mult_transpose, A, Y, X, NULL));
1661c8564c30SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1662c8564c30SJames Wright }
1663