xref: /honee/src/mat-ceed.c (revision c1bdbf00d8adb48041bb92d026716db605362268)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
358600ac3SJames Wright /// @file
440d80af1SJames Wright /// MatCEED implementation
558600ac3SJames Wright 
658600ac3SJames Wright #include <ceed.h>
758600ac3SJames Wright #include <ceed/backend.h>
858600ac3SJames Wright #include <mat-ceed-impl.h>
958600ac3SJames Wright #include <mat-ceed.h>
1040d80af1SJames Wright #include <petsc-ceed-utils.h>
1140d80af1SJames Wright #include <petsc-ceed.h>
12000d2032SJeremy L Thompson #include <petscdm.h>
13000d2032SJeremy L Thompson #include <petscmat.h>
14000d2032SJeremy L Thompson #include <stdbool.h>
15000d2032SJeremy L Thompson #include <stdio.h>
1658600ac3SJames Wright #include <stdlib.h>
1758600ac3SJames Wright #include <string.h>
1858600ac3SJames Wright 
1958600ac3SJames Wright PetscClassId  MATCEED_CLASSID;
20c63b910fSJames Wright PetscLogEvent MATCEED_MULT, MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE, MATCEED_MULT_TRANSPOSE_CEEDOP, MATCEED_ASSEMBLE_DIAGONAL,
21c63b910fSJames Wright     MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, MATCEED_SETUP_PBDIAGONAL, MATCEED_SETUP_PBDIAGONAL_CEEDOP, MATCEED_ASSEMBLE_PBDIAGONAL,
22c63b910fSJames Wright     MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, MATCEED_SETUP_FULL, MATCEED_SETUP_FULL_CEEDOP, MATCEED_ASSEMBLE_FULL, MATCEED_ASSEMBLE_FULL_CEEDOP;
2358600ac3SJames Wright 
2458600ac3SJames Wright /**
2558600ac3SJames Wright   @brief Register MATCEED log events.
2658600ac3SJames Wright 
2758600ac3SJames Wright   Not collective across MPI processes.
2858600ac3SJames Wright 
2958600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
3058600ac3SJames Wright **/
3158600ac3SJames Wright static PetscErrorCode MatCeedRegisterLogEvents() {
3240d80af1SJames Wright   static PetscBool registered = PETSC_FALSE;
3358600ac3SJames Wright 
3458600ac3SJames Wright   PetscFunctionBeginUser;
3558600ac3SJames Wright   if (registered) PetscFunctionReturn(PETSC_SUCCESS);
36c63b910fSJames Wright   PetscCall(PetscClassIdRegister("MatCEED", &MATCEED_CLASSID));
37c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMul", MATCEED_CLASSID, &MATCEED_MULT));
38c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMulCeed", MATCEED_CLASSID, &MATCEED_MULT_CEEDOP));
39c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMulT", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE));
40c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDMulTCeed", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE_CEEDOP));
41c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmDiag", MATCEED_CLASSID, &MATCEED_ASSEMBLE_DIAGONAL));
42c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSU", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL));
43c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL_CEEDOP));
44c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBD", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL));
45c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmPBDCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP));
46c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmSU", MATCEED_CLASSID, &MATCEED_SETUP_FULL));
47c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_FULL_CEEDOP));
48c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsm", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL));
49c63b910fSJames Wright   PetscCall(PetscLogEventRegister("MatCEEDAsmCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL_CEEDOP));
5040d80af1SJames Wright   registered = PETSC_TRUE;
5158600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5258600ac3SJames Wright }
5358600ac3SJames Wright 
5458600ac3SJames Wright /**
5558600ac3SJames Wright   @brief Assemble the point block diagonal of a `MATCEED` into a `MATAIJ` or similar.
5658600ac3SJames Wright          The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`.
5758600ac3SJames Wright          The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail.
5858600ac3SJames Wright 
5958600ac3SJames Wright   Collective across MPI processes.
6058600ac3SJames Wright 
6158600ac3SJames Wright   @param[in]      mat_ceed  `MATCEED` to assemble
6258600ac3SJames Wright   @param[in,out]  mat_coo   `MATAIJ` or similar to assemble into
6358600ac3SJames Wright 
6458600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
6558600ac3SJames Wright **/
6658600ac3SJames Wright static PetscErrorCode MatCeedAssemblePointBlockDiagonalCOO(Mat mat_ceed, Mat mat_coo) {
6758600ac3SJames Wright   MatCeedContext ctx;
6858600ac3SJames Wright 
6958600ac3SJames Wright   PetscFunctionBeginUser;
7058600ac3SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
7158600ac3SJames Wright 
7258600ac3SJames Wright   // Check if COO pattern set
7358600ac3SJames Wright   {
7458600ac3SJames Wright     PetscInt index = -1;
7558600ac3SJames Wright 
7658600ac3SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) {
7758600ac3SJames Wright       if (ctx->mats_assembled_pbd[i] == mat_coo) index = i;
7858600ac3SJames Wright     }
7958600ac3SJames Wright     if (index == -1) {
8058600ac3SJames Wright       PetscInt     *rows_petsc = NULL, *cols_petsc = NULL;
8158600ac3SJames Wright       CeedInt      *rows_ceed, *cols_ceed;
8258600ac3SJames Wright       PetscCount    num_entries;
8358600ac3SJames Wright       PetscLogStage stage_amg_setup;
8458600ac3SJames Wright 
8558600ac3SJames Wright       // -- Assemble sparsity pattern if mat hasn't been assembled before
86c63b910fSJames Wright       PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup));
8758600ac3SJames Wright       if (stage_amg_setup == -1) {
88c63b910fSJames Wright         PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup));
8958600ac3SJames Wright       }
9058600ac3SJames Wright       PetscCall(PetscLogStagePush(stage_amg_setup));
91c63b910fSJames Wright       PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
92c63b910fSJames Wright       PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
9350f50432SJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed));
94c63b910fSJames Wright       PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
95a7dac1d5SJames Wright       PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc));
96a7dac1d5SJames Wright       PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc));
9758600ac3SJames Wright       PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc));
9858600ac3SJames Wright       free(rows_petsc);
9958600ac3SJames Wright       free(cols_petsc);
10050f50432SJames Wright       if (!ctx->coo_values_pbd) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_pbd));
10158600ac3SJames Wright       PetscCall(PetscRealloc(++ctx->num_mats_assembled_pbd * sizeof(Mat), &ctx->mats_assembled_pbd));
10258600ac3SJames Wright       ctx->mats_assembled_pbd[ctx->num_mats_assembled_pbd - 1] = mat_coo;
103c63b910fSJames Wright       PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
10458600ac3SJames Wright       PetscCall(PetscLogStagePop());
10558600ac3SJames Wright     }
10658600ac3SJames Wright   }
10758600ac3SJames Wright 
10858600ac3SJames Wright   // Assemble mat_ceed
109c63b910fSJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
11058600ac3SJames Wright   PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY));
11158600ac3SJames Wright   {
11258600ac3SJames Wright     const CeedScalar *values;
11358600ac3SJames Wright     MatType           mat_type;
11458600ac3SJames Wright     CeedMemType       mem_type = CEED_MEM_HOST;
11558600ac3SJames Wright     PetscBool         is_spd, is_spd_known;
11658600ac3SJames Wright 
11758600ac3SJames Wright     PetscCall(MatGetType(mat_coo, &mat_type));
11858600ac3SJames Wright     if (strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
11958600ac3SJames Wright     else if (strstr(mat_type, "kokkos")) mem_type = CEED_MEM_DEVICE;
12058600ac3SJames Wright     else mem_type = CEED_MEM_HOST;
12158600ac3SJames Wright 
122c63b910fSJames Wright     PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
12350f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(ctx->op_mult, ctx->coo_values_pbd, CEED_REQUEST_IMMEDIATE));
124c63b910fSJames Wright     PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
12550f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_pbd, mem_type, &values));
12658600ac3SJames Wright     PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES));
12758600ac3SJames Wright     PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd));
12858600ac3SJames Wright     if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd));
12950f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_pbd, &values));
13058600ac3SJames Wright   }
13158600ac3SJames Wright   PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY));
132c63b910fSJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL));
13358600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
13458600ac3SJames Wright }
13558600ac3SJames Wright 
13658600ac3SJames Wright /**
13758600ac3SJames Wright   @brief Assemble inner `Mat` for diagonal `PC` operations
13858600ac3SJames Wright 
13958600ac3SJames Wright   Collective across MPI processes.
14058600ac3SJames Wright 
14158600ac3SJames Wright   @param[in]   mat_ceed      `MATCEED` to invert
14258600ac3SJames Wright   @param[in]   use_ceed_pbd  Boolean flag to use libCEED PBD assembly
14358600ac3SJames Wright   @param[out]  mat_inner     Inner `Mat` for diagonal operations
14458600ac3SJames Wright 
14558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
14658600ac3SJames Wright **/
14758600ac3SJames Wright static PetscErrorCode MatCeedAssembleInnerBlockDiagonalMat(Mat mat_ceed, PetscBool use_ceed_pbd, Mat *mat_inner) {
14858600ac3SJames Wright   MatCeedContext ctx;
14958600ac3SJames Wright 
15058600ac3SJames Wright   PetscFunctionBeginUser;
15158600ac3SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
15258600ac3SJames Wright   if (use_ceed_pbd) {
15358600ac3SJames Wright     // Check if COO pattern set
15440d80af1SJames Wright     if (!ctx->mat_assembled_pbd_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_pbd_internal));
15558600ac3SJames Wright 
15658600ac3SJames Wright     // Assemble mat_assembled_full_internal
15758600ac3SJames Wright     PetscCall(MatCeedAssemblePointBlockDiagonalCOO(mat_ceed, ctx->mat_assembled_pbd_internal));
15858600ac3SJames Wright     if (mat_inner) *mat_inner = ctx->mat_assembled_pbd_internal;
15958600ac3SJames Wright   } else {
16058600ac3SJames Wright     // Check if COO pattern set
16140d80af1SJames Wright     if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal));
16258600ac3SJames Wright 
16358600ac3SJames Wright     // Assemble mat_assembled_full_internal
16458600ac3SJames Wright     PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal));
16558600ac3SJames Wright     if (mat_inner) *mat_inner = ctx->mat_assembled_full_internal;
16658600ac3SJames Wright   }
16758600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
16858600ac3SJames Wright }
16958600ac3SJames Wright 
17058600ac3SJames Wright /**
171*c1bdbf00SJames Wright   @brief Get on-process diagonal block of `MATCEED`
172*c1bdbf00SJames Wright 
173*c1bdbf00SJames Wright   This is a block per-process of the diagonal of the global matrix.
174*c1bdbf00SJames Wright   This is NOT the diagonal blocks associated with the block size of the matrix (i.e. `MatSetBlockSize()` has no effect on this function).
17558600ac3SJames Wright 
17658600ac3SJames Wright   Collective across MPI processes.
17758600ac3SJames Wright 
17858600ac3SJames Wright   @param[in]   mat_ceed   `MATCEED` to invert
17958600ac3SJames Wright   @param[out]  mat_block  The diagonal block matrix
18058600ac3SJames Wright 
18158600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
18258600ac3SJames Wright **/
18358600ac3SJames Wright static PetscErrorCode MatGetDiagonalBlock_Ceed(Mat mat_ceed, Mat *mat_block) {
18458600ac3SJames Wright   MatCeedContext ctx;
18558600ac3SJames Wright 
18658600ac3SJames Wright   PetscFunctionBeginUser;
18758600ac3SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
18858600ac3SJames Wright 
189*c1bdbf00SJames Wright   // Check if COO pattern set
190*c1bdbf00SJames Wright   if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal));
19158600ac3SJames Wright 
192*c1bdbf00SJames Wright   // Assemble mat_assembled_full_internal
193*c1bdbf00SJames Wright   PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal));
194*c1bdbf00SJames Wright 
195*c1bdbf00SJames Wright   // Get diagonal block
196*c1bdbf00SJames Wright   PetscCall(MatGetDiagonalBlock(ctx->mat_assembled_full_internal, mat_block));
19758600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19858600ac3SJames Wright }
19958600ac3SJames Wright 
20058600ac3SJames Wright /**
20158600ac3SJames Wright   @brief Invert `MATCEED` diagonal block for Jacobi.
20258600ac3SJames Wright 
20358600ac3SJames Wright   Collective across MPI processes.
20458600ac3SJames Wright 
20558600ac3SJames Wright   @param[in]   mat_ceed  `MATCEED` to invert
20658600ac3SJames Wright   @param[out]  values    The block inverses in column major order
20758600ac3SJames Wright 
20858600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
20958600ac3SJames Wright **/
21058600ac3SJames Wright static PetscErrorCode MatInvertBlockDiagonal_Ceed(Mat mat_ceed, const PetscScalar **values) {
21158600ac3SJames Wright   Mat            mat_inner = NULL;
21258600ac3SJames Wright   MatCeedContext ctx;
21358600ac3SJames Wright 
21458600ac3SJames Wright   PetscFunctionBeginUser;
21558600ac3SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
21658600ac3SJames Wright 
21758600ac3SJames Wright   // Assemble inner mat if needed
21858600ac3SJames Wright   PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_pbd_valid, &mat_inner));
21958600ac3SJames Wright 
22058600ac3SJames Wright   // Invert PB diagonal
22158600ac3SJames Wright   PetscCall(MatInvertBlockDiagonal(mat_inner, values));
22258600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
22358600ac3SJames Wright }
22458600ac3SJames Wright 
22558600ac3SJames Wright /**
22658600ac3SJames Wright   @brief Invert `MATCEED` variable diagonal block for Jacobi.
22758600ac3SJames Wright 
22858600ac3SJames Wright   Collective across MPI processes.
22958600ac3SJames Wright 
23058600ac3SJames Wright   @param[in]   mat_ceed     `MATCEED` to invert
23158600ac3SJames Wright   @param[in]   num_blocks   The number of blocks on the process
23258600ac3SJames Wright   @param[in]   block_sizes  The size of each block on the process
23358600ac3SJames Wright   @param[out]  values       The block inverses in column major order
23458600ac3SJames Wright 
23558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
23658600ac3SJames Wright **/
23758600ac3SJames Wright static PetscErrorCode MatInvertVariableBlockDiagonal_Ceed(Mat mat_ceed, PetscInt num_blocks, const PetscInt *block_sizes, PetscScalar *values) {
23858600ac3SJames Wright   Mat            mat_inner = NULL;
23958600ac3SJames Wright   MatCeedContext ctx;
24058600ac3SJames Wright 
24158600ac3SJames Wright   PetscFunctionBeginUser;
24258600ac3SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
24358600ac3SJames Wright 
24458600ac3SJames Wright   // Assemble inner mat if needed
24558600ac3SJames Wright   PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_vpbd_valid, &mat_inner));
24658600ac3SJames Wright 
24758600ac3SJames Wright   // Invert PB diagonal
24858600ac3SJames Wright   PetscCall(MatInvertVariableBlockDiagonal(mat_inner, num_blocks, block_sizes, values));
24958600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
25058600ac3SJames Wright }
25158600ac3SJames Wright 
252e90c2ceeSJames Wright /**
253e90c2ceeSJames Wright   @brief View `MATCEED`.
254e90c2ceeSJames Wright 
255e90c2ceeSJames Wright   Collective across MPI processes.
256e90c2ceeSJames Wright 
257e90c2ceeSJames Wright   @param[in]   mat_ceed  `MATCEED` to view
258e90c2ceeSJames Wright   @param[in]   viewer    The visualization context
259e90c2ceeSJames Wright 
260e90c2ceeSJames Wright   @return An error code: 0 - success, otherwise - failure
261e90c2ceeSJames Wright **/
262e90c2ceeSJames Wright static PetscErrorCode MatView_Ceed(Mat mat_ceed, PetscViewer viewer) {
263e90c2ceeSJames Wright   PetscBool         is_ascii;
264e90c2ceeSJames Wright   PetscViewerFormat format;
265000d2032SJeremy L Thompson   PetscMPIInt       size, rank;
266e90c2ceeSJames Wright   MatCeedContext    ctx;
267e90c2ceeSJames Wright 
268e90c2ceeSJames Wright   PetscFunctionBeginUser;
269e90c2ceeSJames Wright   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
270e90c2ceeSJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
271e90c2ceeSJames Wright   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat_ceed), &viewer));
272e90c2ceeSJames Wright 
273e90c2ceeSJames Wright   PetscCall(PetscViewerGetFormat(viewer, &format));
274e90c2ceeSJames Wright   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat_ceed), &size));
275e90c2ceeSJames Wright   if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS);
276e90c2ceeSJames Wright 
277000d2032SJeremy L Thompson   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat_ceed), &rank));
278000d2032SJeremy L Thompson   if (rank != 0) PetscFunctionReturn(PETSC_SUCCESS);
279000d2032SJeremy L Thompson 
280e90c2ceeSJames Wright   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii));
281e90c2ceeSJames Wright   {
282000d2032SJeremy L Thompson     PetscBool is_detailed     = format == PETSC_VIEWER_ASCII_INFO_DETAIL;
283000d2032SJeremy L Thompson     char      rank_string[16] = {'\0'};
284e90c2ceeSJames Wright     FILE     *file;
285e90c2ceeSJames Wright 
28640d80af1SJames Wright     PetscCall(PetscViewerASCIIPrintf(viewer, "MatCEED:\n  Default COO MatType:%s\n", ctx->coo_mat_type));
287537ec908SJames Wright     PetscCall(PetscViewerASCIIPushTab(viewer));  // MatCEED
288000d2032SJeremy L Thompson     PetscCall(PetscSNPrintf(rank_string, 16, "on Rank %d", rank));
289000d2032SJeremy L Thompson     PetscCall(PetscViewerASCIIPrintf(viewer, "CeedOperator Apply %s:\n", is_detailed ? rank_string : "Summary"));
290537ec908SJames Wright     PetscCall(PetscViewerASCIIGetPointer(viewer, &file));
291537ec908SJames Wright     PetscCall(PetscViewerASCIIPushTab(viewer));  // CeedOperator
292000d2032SJeremy L Thompson     if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult, file));
293000d2032SJeremy L Thompson     else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult, file));
294537ec908SJames Wright     PetscCall(PetscViewerASCIIPopTab(viewer));  // CeedOperator
295e90c2ceeSJames Wright     if (ctx->op_mult_transpose) {
296000d2032SJeremy L Thompson       PetscCall(PetscViewerASCIIPrintf(viewer, "CeedOperator ApplyTranspose %s:\n", is_detailed ? rank_string : "Summary"));
297537ec908SJames Wright       PetscCall(PetscViewerASCIIPushTab(viewer));  // CeedOperator
298000d2032SJeremy L Thompson       if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult_transpose, file));
299000d2032SJeremy L Thompson       else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult_transpose, file));
300537ec908SJames Wright       PetscCall(PetscViewerASCIIPopTab(viewer));  // CeedOperator
301e90c2ceeSJames Wright     }
302537ec908SJames Wright     PetscCall(PetscViewerASCIIPopTab(viewer));  // MatCEED
303e90c2ceeSJames Wright   }
304e90c2ceeSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
305e90c2ceeSJames Wright }
306e90c2ceeSJames Wright 
30758600ac3SJames Wright // -----------------------------------------------------------------------------
30858600ac3SJames Wright // MatCeed
30958600ac3SJames Wright // -----------------------------------------------------------------------------
31058600ac3SJames Wright 
31158600ac3SJames Wright /**
31258600ac3SJames Wright   @brief Create PETSc `Mat` from libCEED operators.
31358600ac3SJames Wright 
31458600ac3SJames Wright   Collective across MPI processes.
31558600ac3SJames Wright 
31658600ac3SJames Wright   @param[in]   dm_x                      Input `DM`
31758600ac3SJames Wright   @param[in]   dm_y                      Output `DM`
31858600ac3SJames Wright   @param[in]   op_mult                   `CeedOperator` for forward evaluation
31958600ac3SJames Wright   @param[in]   op_mult_transpose         `CeedOperator` for transpose evaluation
32058600ac3SJames Wright   @param[out]  mat                        New MatCeed
32158600ac3SJames Wright 
32258600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
32358600ac3SJames Wright **/
324000d2032SJeremy L Thompson PetscErrorCode MatCreateCeed(DM dm_x, DM dm_y, CeedOperator op_mult, CeedOperator op_mult_transpose, Mat *mat) {
32558600ac3SJames Wright   PetscInt       X_l_size, X_g_size, Y_l_size, Y_g_size;
32658600ac3SJames Wright   VecType        vec_type;
32758600ac3SJames Wright   MatCeedContext ctx;
32858600ac3SJames Wright 
32958600ac3SJames Wright   PetscFunctionBeginUser;
33058600ac3SJames Wright   PetscCall(MatCeedRegisterLogEvents());
33158600ac3SJames Wright 
33258600ac3SJames Wright   // Collect context data
33358600ac3SJames Wright   PetscCall(DMGetVecType(dm_x, &vec_type));
33458600ac3SJames Wright   {
33558600ac3SJames Wright     Vec X;
33658600ac3SJames Wright 
33758600ac3SJames Wright     PetscCall(DMGetGlobalVector(dm_x, &X));
33858600ac3SJames Wright     PetscCall(VecGetSize(X, &X_g_size));
33958600ac3SJames Wright     PetscCall(VecGetLocalSize(X, &X_l_size));
34058600ac3SJames Wright     PetscCall(DMRestoreGlobalVector(dm_x, &X));
34158600ac3SJames Wright   }
34258600ac3SJames Wright   if (dm_y) {
34358600ac3SJames Wright     Vec Y;
34458600ac3SJames Wright 
34558600ac3SJames Wright     PetscCall(DMGetGlobalVector(dm_y, &Y));
34658600ac3SJames Wright     PetscCall(VecGetSize(Y, &Y_g_size));
34758600ac3SJames Wright     PetscCall(VecGetLocalSize(Y, &Y_l_size));
34858600ac3SJames Wright     PetscCall(DMRestoreGlobalVector(dm_y, &Y));
34958600ac3SJames Wright   } else {
35058600ac3SJames Wright     dm_y     = dm_x;
35158600ac3SJames Wright     Y_g_size = X_g_size;
35258600ac3SJames Wright     Y_l_size = X_l_size;
35358600ac3SJames Wright   }
35440d80af1SJames Wright 
35558600ac3SJames Wright   // Create context
35658600ac3SJames Wright   {
35758600ac3SJames Wright     Vec X_loc, Y_loc_transpose = NULL;
35858600ac3SJames Wright 
35958600ac3SJames Wright     PetscCall(DMCreateLocalVector(dm_x, &X_loc));
36058600ac3SJames Wright     PetscCall(VecZeroEntries(X_loc));
36158600ac3SJames Wright     if (op_mult_transpose) {
36258600ac3SJames Wright       PetscCall(DMCreateLocalVector(dm_y, &Y_loc_transpose));
36358600ac3SJames Wright       PetscCall(VecZeroEntries(Y_loc_transpose));
36458600ac3SJames Wright     }
365c63b910fSJames Wright     PetscCall(MatCeedContextCreate(dm_x, dm_y, X_loc, Y_loc_transpose, op_mult, op_mult_transpose, MATCEED_MULT, MATCEED_MULT_TRANSPOSE,
366c63b910fSJames Wright                                    MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE_CEEDOP, &ctx));
36758600ac3SJames Wright     PetscCall(VecDestroy(&X_loc));
36858600ac3SJames Wright     PetscCall(VecDestroy(&Y_loc_transpose));
36958600ac3SJames Wright   }
37058600ac3SJames Wright 
37158600ac3SJames Wright   // Create mat
37258600ac3SJames Wright   PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm_x), Y_l_size, X_l_size, Y_g_size, X_g_size, ctx, mat));
37358600ac3SJames Wright   PetscCall(PetscObjectChangeTypeName((PetscObject)*mat, MATCEED));
37458600ac3SJames Wright   // -- Set block and variable block sizes
37558600ac3SJames Wright   if (dm_x == dm_y) {
37658600ac3SJames Wright     MatType dm_mat_type, dm_mat_type_copy;
37758600ac3SJames Wright     Mat     temp_mat;
37858600ac3SJames Wright 
37958600ac3SJames Wright     PetscCall(DMGetMatType(dm_x, &dm_mat_type));
38058600ac3SJames Wright     PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy));
38158600ac3SJames Wright     PetscCall(DMSetMatType(dm_x, MATAIJ));
38258600ac3SJames Wright     PetscCall(DMCreateMatrix(dm_x, &temp_mat));
38358600ac3SJames Wright     PetscCall(DMSetMatType(dm_x, dm_mat_type_copy));
38458600ac3SJames Wright     PetscCall(PetscFree(dm_mat_type_copy));
38558600ac3SJames Wright 
38658600ac3SJames Wright     {
38758600ac3SJames Wright       PetscInt        block_size, num_blocks, max_vblock_size = PETSC_INT_MAX;
38858600ac3SJames Wright       const PetscInt *vblock_sizes;
38958600ac3SJames Wright 
39058600ac3SJames Wright       // -- Get block sizes
39158600ac3SJames Wright       PetscCall(MatGetBlockSize(temp_mat, &block_size));
39258600ac3SJames Wright       PetscCall(MatGetVariableBlockSizes(temp_mat, &num_blocks, &vblock_sizes));
39358600ac3SJames Wright       {
39458600ac3SJames Wright         PetscInt local_min_max[2] = {0}, global_min_max[2] = {0, PETSC_INT_MAX};
39558600ac3SJames Wright 
39658600ac3SJames Wright         for (PetscInt i = 0; i < num_blocks; i++) local_min_max[1] = PetscMax(local_min_max[1], vblock_sizes[i]);
39758600ac3SJames Wright         PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_min_max, global_min_max));
39858600ac3SJames Wright         max_vblock_size = global_min_max[1];
39958600ac3SJames Wright       }
40058600ac3SJames Wright 
40158600ac3SJames Wright       // -- Copy block sizes
40258600ac3SJames Wright       if (block_size > 1) PetscCall(MatSetBlockSize(*mat, block_size));
40358600ac3SJames Wright       if (num_blocks) PetscCall(MatSetVariableBlockSizes(*mat, num_blocks, (PetscInt *)vblock_sizes));
40458600ac3SJames Wright 
40558600ac3SJames Wright       // -- Check libCEED compatibility
40658600ac3SJames Wright       {
40758600ac3SJames Wright         bool is_composite;
40858600ac3SJames Wright 
40958600ac3SJames Wright         ctx->is_ceed_pbd_valid  = PETSC_TRUE;
41058600ac3SJames Wright         ctx->is_ceed_vpbd_valid = PETSC_TRUE;
41150f50432SJames Wright         PetscCallCeed(ctx->ceed, CeedOperatorIsComposite(op_mult, &is_composite));
41258600ac3SJames Wright         if (is_composite) {
41358600ac3SJames Wright           CeedInt       num_sub_operators;
41458600ac3SJames Wright           CeedOperator *sub_operators;
41558600ac3SJames Wright 
41650f50432SJames Wright           PetscCallCeed(ctx->ceed, CeedCompositeOperatorGetNumSub(op_mult, &num_sub_operators));
41750f50432SJames Wright           PetscCallCeed(ctx->ceed, CeedCompositeOperatorGetSubList(op_mult, &sub_operators));
41858600ac3SJames Wright           for (CeedInt i = 0; i < num_sub_operators; i++) {
41958600ac3SJames Wright             CeedInt                  num_bases, num_comp;
42058600ac3SJames Wright             CeedBasis               *active_bases;
42158600ac3SJames Wright             CeedOperatorAssemblyData assembly_data;
42258600ac3SJames Wright 
42350f50432SJames Wright             PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(sub_operators[i], &assembly_data));
42450f50432SJames Wright             PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL));
42550f50432SJames Wright             PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp));
42658600ac3SJames Wright             if (num_bases > 1) {
42758600ac3SJames Wright               ctx->is_ceed_pbd_valid  = PETSC_FALSE;
42858600ac3SJames Wright               ctx->is_ceed_vpbd_valid = PETSC_FALSE;
42958600ac3SJames Wright             }
43058600ac3SJames Wright             if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE;
43158600ac3SJames Wright             if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE;
43258600ac3SJames Wright           }
43358600ac3SJames Wright         } else {
43458600ac3SJames Wright           // LCOV_EXCL_START
43558600ac3SJames Wright           CeedInt                  num_bases, num_comp;
43658600ac3SJames Wright           CeedBasis               *active_bases;
43758600ac3SJames Wright           CeedOperatorAssemblyData assembly_data;
43858600ac3SJames Wright 
43950f50432SJames Wright           PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(op_mult, &assembly_data));
44050f50432SJames Wright           PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL));
44150f50432SJames Wright           PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp));
44258600ac3SJames Wright           if (num_bases > 1) {
44358600ac3SJames Wright             ctx->is_ceed_pbd_valid  = PETSC_FALSE;
44458600ac3SJames Wright             ctx->is_ceed_vpbd_valid = PETSC_FALSE;
44558600ac3SJames Wright           }
44658600ac3SJames Wright           if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE;
44758600ac3SJames Wright           if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE;
44858600ac3SJames Wright           // LCOV_EXCL_STOP
44958600ac3SJames Wright         }
45058600ac3SJames Wright         {
45158600ac3SJames Wright           PetscInt local_is_valid[2], global_is_valid[2];
45258600ac3SJames Wright 
45358600ac3SJames Wright           local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_pbd_valid;
45458600ac3SJames Wright           PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid));
45558600ac3SJames Wright           ctx->is_ceed_pbd_valid = global_is_valid[0];
45658600ac3SJames Wright           local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_vpbd_valid;
45758600ac3SJames Wright           PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid));
45858600ac3SJames Wright           ctx->is_ceed_vpbd_valid = global_is_valid[0];
45958600ac3SJames Wright         }
46058600ac3SJames Wright       }
46158600ac3SJames Wright     }
46258600ac3SJames Wright     PetscCall(MatDestroy(&temp_mat));
46358600ac3SJames Wright   }
46458600ac3SJames Wright   // -- Set internal mat type
46558600ac3SJames Wright   {
46658600ac3SJames Wright     VecType vec_type;
46740d80af1SJames Wright     MatType coo_mat_type;
46858600ac3SJames Wright 
46958600ac3SJames Wright     PetscCall(VecGetType(ctx->X_loc, &vec_type));
47040d80af1SJames Wright     if (strstr(vec_type, VECCUDA)) coo_mat_type = MATAIJCUSPARSE;
47140d80af1SJames Wright     else if (strstr(vec_type, VECKOKKOS)) coo_mat_type = MATAIJKOKKOS;
47240d80af1SJames Wright     else coo_mat_type = MATAIJ;
47340d80af1SJames Wright     PetscCall(PetscStrallocpy(coo_mat_type, &ctx->coo_mat_type));
47458600ac3SJames Wright   }
47558600ac3SJames Wright   // -- Set mat operations
47658600ac3SJames Wright   PetscCall(MatShellSetContextDestroy(*mat, (PetscErrorCode(*)(void *))MatCeedContextDestroy));
477e90c2ceeSJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_VIEW, (void (*)(void))MatView_Ceed));
47858600ac3SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_MULT, (void (*)(void))MatMult_Ceed));
47958600ac3SJames Wright   if (op_mult_transpose) PetscCall(MatShellSetOperation(*mat, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed));
48058600ac3SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed));
48158600ac3SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed));
48258600ac3SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_INVERT_BLOCK_DIAGONAL, (void (*)(void))MatInvertBlockDiagonal_Ceed));
48358600ac3SJames Wright   PetscCall(MatShellSetOperation(*mat, MATOP_INVERT_VBLOCK_DIAGONAL, (void (*)(void))MatInvertVariableBlockDiagonal_Ceed));
48458600ac3SJames Wright   PetscCall(MatShellSetVecType(*mat, vec_type));
48558600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
48658600ac3SJames Wright }
48758600ac3SJames Wright 
48858600ac3SJames Wright /**
48958600ac3SJames Wright   @brief Copy `MATCEED` into a compatible `Mat` with type `MatShell` or `MATCEED`.
49058600ac3SJames Wright 
49158600ac3SJames Wright   Collective across MPI processes.
49258600ac3SJames Wright 
49358600ac3SJames Wright   @param[in]   mat_ceed   `MATCEED` to copy from
49458600ac3SJames Wright   @param[out]  mat_other  `MatShell` or `MATCEED` to copy into
49558600ac3SJames Wright 
49658600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
49758600ac3SJames Wright **/
49858600ac3SJames Wright PetscErrorCode MatCeedCopy(Mat mat_ceed, Mat mat_other) {
49958600ac3SJames Wright   PetscFunctionBeginUser;
50058600ac3SJames Wright   PetscCall(MatCeedRegisterLogEvents());
50158600ac3SJames Wright 
50258600ac3SJames Wright   // Check type compatibility
50358600ac3SJames Wright   {
50440d80af1SJames Wright     PetscBool is_matceed = PETSC_FALSE, is_matshell = PETSC_FALSE;
50558600ac3SJames Wright     MatType   mat_type_ceed, mat_type_other;
50658600ac3SJames Wright 
50758600ac3SJames Wright     PetscCall(MatGetType(mat_ceed, &mat_type_ceed));
50840d80af1SJames Wright     PetscCall(PetscStrcmp(mat_type_ceed, MATCEED, &is_matceed));
50940d80af1SJames Wright     PetscCheck(is_matceed, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_ceed must have type " MATCEED);
51040d80af1SJames Wright     PetscCall(MatGetType(mat_other, &mat_type_other));
51140d80af1SJames Wright     PetscCall(PetscStrcmp(mat_type_other, MATCEED, &is_matceed));
51240d80af1SJames Wright     PetscCall(PetscStrcmp(mat_type_other, MATSHELL, &is_matceed));
51340d80af1SJames Wright     PetscCheck(is_matceed || is_matshell, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_other must have type " MATCEED " or " MATSHELL);
51458600ac3SJames Wright   }
51558600ac3SJames Wright 
51658600ac3SJames Wright   // Check dimension compatibility
51758600ac3SJames Wright   {
51858600ac3SJames 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;
51958600ac3SJames Wright 
52058600ac3SJames Wright     PetscCall(MatGetSize(mat_ceed, &Y_g_ceed_size, &X_g_ceed_size));
52158600ac3SJames Wright     PetscCall(MatGetLocalSize(mat_ceed, &Y_l_ceed_size, &X_l_ceed_size));
52258600ac3SJames Wright     PetscCall(MatGetSize(mat_ceed, &Y_g_other_size, &X_g_other_size));
52358600ac3SJames Wright     PetscCall(MatGetLocalSize(mat_ceed, &Y_l_other_size, &X_l_other_size));
52458600ac3SJames 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) &&
52558600ac3SJames Wright                    (X_l_ceed_size == X_l_other_size),
52658600ac3SJames Wright                PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ,
52758600ac3SJames Wright                "mat_ceed and mat_other must have compatible sizes; found mat_ceed (Global: %" PetscInt_FMT ", %" PetscInt_FMT
52858600ac3SJames Wright                "; Local: %" PetscInt_FMT ", %" PetscInt_FMT ") mat_other (Global: %" PetscInt_FMT ", %" PetscInt_FMT "; Local: %" PetscInt_FMT
52958600ac3SJames Wright                ", %" PetscInt_FMT ")",
53058600ac3SJames 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);
53158600ac3SJames Wright   }
53258600ac3SJames Wright 
53358600ac3SJames Wright   // Convert
53458600ac3SJames Wright   {
53558600ac3SJames Wright     VecType        vec_type;
53658600ac3SJames Wright     MatCeedContext ctx;
53758600ac3SJames Wright 
53858600ac3SJames Wright     PetscCall(PetscObjectChangeTypeName((PetscObject)mat_other, MATCEED));
53958600ac3SJames Wright     PetscCall(MatShellGetContext(mat_ceed, &ctx));
54058600ac3SJames Wright     PetscCall(MatCeedContextReference(ctx));
54158600ac3SJames Wright     PetscCall(MatShellSetContext(mat_other, ctx));
54258600ac3SJames Wright     PetscCall(MatShellSetContextDestroy(mat_other, (PetscErrorCode(*)(void *))MatCeedContextDestroy));
543e90c2ceeSJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_VIEW, (void (*)(void))MatView_Ceed));
54458600ac3SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_MULT, (void (*)(void))MatMult_Ceed));
54558600ac3SJames Wright     if (ctx->op_mult_transpose) PetscCall(MatShellSetOperation(mat_other, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed));
54658600ac3SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed));
54758600ac3SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed));
54858600ac3SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_INVERT_BLOCK_DIAGONAL, (void (*)(void))MatInvertBlockDiagonal_Ceed));
54958600ac3SJames Wright     PetscCall(MatShellSetOperation(mat_other, MATOP_INVERT_VBLOCK_DIAGONAL, (void (*)(void))MatInvertVariableBlockDiagonal_Ceed));
55058600ac3SJames Wright     {
55158600ac3SJames Wright       PetscInt block_size;
55258600ac3SJames Wright 
55358600ac3SJames Wright       PetscCall(MatGetBlockSize(mat_ceed, &block_size));
55458600ac3SJames Wright       if (block_size > 1) PetscCall(MatSetBlockSize(mat_other, block_size));
55558600ac3SJames Wright     }
55658600ac3SJames Wright     {
55758600ac3SJames Wright       PetscInt        num_blocks;
55858600ac3SJames Wright       const PetscInt *block_sizes;
55958600ac3SJames Wright 
56058600ac3SJames Wright       PetscCall(MatGetVariableBlockSizes(mat_ceed, &num_blocks, &block_sizes));
56158600ac3SJames Wright       if (num_blocks) PetscCall(MatSetVariableBlockSizes(mat_other, num_blocks, (PetscInt *)block_sizes));
56258600ac3SJames Wright     }
56358600ac3SJames Wright     PetscCall(DMGetVecType(ctx->dm_x, &vec_type));
56458600ac3SJames Wright     PetscCall(MatShellSetVecType(mat_other, vec_type));
56558600ac3SJames Wright   }
56658600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
56758600ac3SJames Wright }
56858600ac3SJames Wright 
56958600ac3SJames Wright /**
570ba95ebb4SJeremy L Thompson   @brief Mark `CeedQFunction` data as updated and the `CeedQFunction` as requiring re-assembly for a `MatCEED`.
571000d2032SJeremy L Thompson 
572000d2032SJeremy L Thompson   Collective across MPI processes.
573000d2032SJeremy L Thompson 
574000d2032SJeremy L Thompson   @param[in]   mat_ceed       `MATCEED`
575000d2032SJeremy L Thompson   @param[out]  update_needed  Boolean flag indicating `CeedQFunction` update needed
576000d2032SJeremy L Thompson 
577000d2032SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
578000d2032SJeremy L Thompson **/
579000d2032SJeremy L Thompson PetscErrorCode MatCeedSetAssemblyDataUpdateNeeded(Mat mat_ceed, PetscBool update_needed) {
580000d2032SJeremy L Thompson   MatCeedContext ctx;
581000d2032SJeremy L Thompson 
582000d2032SJeremy L Thompson   PetscFunctionBeginUser;
583000d2032SJeremy L Thompson   PetscCall(MatShellGetContext(mat_ceed, &ctx));
584000d2032SJeremy L Thompson   PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult, update_needed));
585000d2032SJeremy L Thompson   if (ctx->op_mult_transpose) {
586000d2032SJeremy L Thompson     PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult_transpose, update_needed));
587000d2032SJeremy L Thompson   }
588397c0187SJeremy L Thompson   if (update_needed) {
589397c0187SJeremy L Thompson     PetscCall(MatAssemblyBegin(mat_ceed, MAT_FINAL_ASSEMBLY));
590397c0187SJeremy L Thompson     PetscCall(MatAssemblyEnd(mat_ceed, MAT_FINAL_ASSEMBLY));
591397c0187SJeremy L Thompson   }
592000d2032SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
593000d2032SJeremy L Thompson }
594000d2032SJeremy L Thompson 
595000d2032SJeremy L Thompson /**
59640d80af1SJames Wright   @brief Setup a `Mat` with the same COO pattern as a `MatCEED`.
59740d80af1SJames Wright 
59840d80af1SJames Wright   Collective across MPI processes.
59940d80af1SJames Wright 
60040d80af1SJames Wright   @param[in]   mat_ceed  `MATCEED`
60140d80af1SJames Wright   @param[out]  mat_coo   Sparse `Mat` with same COO pattern
60240d80af1SJames Wright 
60340d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
60440d80af1SJames Wright **/
60540d80af1SJames Wright PetscErrorCode MatCeedCreateMatCOO(Mat mat_ceed, Mat *mat_coo) {
60640d80af1SJames Wright   MatCeedContext ctx;
60740d80af1SJames Wright 
60840d80af1SJames Wright   PetscFunctionBeginUser;
60940d80af1SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
61040d80af1SJames Wright 
61140d80af1SJames 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");
61240d80af1SJames Wright 
61340d80af1SJames Wright   // Check cl mat type
61440d80af1SJames Wright   {
61540d80af1SJames Wright     PetscBool is_coo_mat_type_cl = PETSC_FALSE;
61640d80af1SJames Wright     char      coo_mat_type_cl[64];
61740d80af1SJames Wright 
61840d80af1SJames Wright     // Check for specific CL coo mat type for this Mat
61940d80af1SJames Wright     {
62040d80af1SJames Wright       const char *mat_ceed_prefix = NULL;
62140d80af1SJames Wright 
62240d80af1SJames Wright       PetscCall(MatGetOptionsPrefix(mat_ceed, &mat_ceed_prefix));
62340d80af1SJames Wright       PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), mat_ceed_prefix, "", NULL);
62440d80af1SJames Wright       PetscCall(PetscOptionsFList("-ceed_coo_mat_type", "Default MATCEED COO assembly MatType", NULL, MatList, coo_mat_type_cl, coo_mat_type_cl,
62540d80af1SJames Wright                                   sizeof(coo_mat_type_cl), &is_coo_mat_type_cl));
62640d80af1SJames Wright       PetscOptionsEnd();
62740d80af1SJames Wright       if (is_coo_mat_type_cl) {
62840d80af1SJames Wright         PetscCall(PetscFree(ctx->coo_mat_type));
62940d80af1SJames Wright         PetscCall(PetscStrallocpy(coo_mat_type_cl, &ctx->coo_mat_type));
63040d80af1SJames Wright       }
63140d80af1SJames Wright     }
63240d80af1SJames Wright   }
63340d80af1SJames Wright 
63440d80af1SJames Wright   // Create sparse matrix
63540d80af1SJames Wright   {
63640d80af1SJames Wright     MatType dm_mat_type, dm_mat_type_copy;
63740d80af1SJames Wright 
63840d80af1SJames Wright     PetscCall(DMGetMatType(ctx->dm_x, &dm_mat_type));
63940d80af1SJames Wright     PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy));
64040d80af1SJames Wright     PetscCall(DMSetMatType(ctx->dm_x, ctx->coo_mat_type));
64140d80af1SJames Wright     PetscCall(DMCreateMatrix(ctx->dm_x, mat_coo));
64240d80af1SJames Wright     PetscCall(DMSetMatType(ctx->dm_x, dm_mat_type_copy));
64340d80af1SJames Wright     PetscCall(PetscFree(dm_mat_type_copy));
64440d80af1SJames Wright   }
64540d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
64640d80af1SJames Wright }
64740d80af1SJames Wright 
64840d80af1SJames Wright /**
64940d80af1SJames Wright   @brief Setup the COO preallocation `MATCEED` into a `MATAIJ` or similar.
65058600ac3SJames Wright          The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail.
65158600ac3SJames Wright 
65258600ac3SJames Wright   Collective across MPI processes.
65358600ac3SJames Wright 
65458600ac3SJames Wright   @param[in]      mat_ceed  `MATCEED` to assemble
65558600ac3SJames Wright   @param[in,out]  mat_coo   `MATAIJ` or similar to assemble into
65658600ac3SJames Wright 
65758600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
65858600ac3SJames Wright **/
65940d80af1SJames Wright PetscErrorCode MatCeedSetPreallocationCOO(Mat mat_ceed, Mat mat_coo) {
66058600ac3SJames Wright   MatCeedContext ctx;
66158600ac3SJames Wright 
66258600ac3SJames Wright   PetscFunctionBeginUser;
66358600ac3SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
66458600ac3SJames Wright 
66558600ac3SJames Wright   {
66658600ac3SJames Wright     PetscInt     *rows_petsc = NULL, *cols_petsc = NULL;
66758600ac3SJames Wright     CeedInt      *rows_ceed, *cols_ceed;
66858600ac3SJames Wright     PetscCount    num_entries;
66958600ac3SJames Wright     PetscLogStage stage_amg_setup;
67058600ac3SJames Wright 
67158600ac3SJames Wright     // -- Assemble sparsity pattern if mat hasn't been assembled before
672c63b910fSJames Wright     PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup));
67358600ac3SJames Wright     if (stage_amg_setup == -1) {
674c63b910fSJames Wright       PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup));
67558600ac3SJames Wright     }
67658600ac3SJames Wright     PetscCall(PetscLogStagePush(stage_amg_setup));
677c63b910fSJames Wright     PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL));
678c63b910fSJames Wright     PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
67950f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed));
680c63b910fSJames Wright     PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
681a7dac1d5SJames Wright     PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc));
682a7dac1d5SJames Wright     PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc));
68358600ac3SJames Wright     PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc));
68458600ac3SJames Wright     free(rows_petsc);
68558600ac3SJames Wright     free(cols_petsc);
68650f50432SJames Wright     if (!ctx->coo_values_full) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_full));
68758600ac3SJames Wright     PetscCall(PetscRealloc(++ctx->num_mats_assembled_full * sizeof(Mat), &ctx->mats_assembled_full));
68858600ac3SJames Wright     ctx->mats_assembled_full[ctx->num_mats_assembled_full - 1] = mat_coo;
689c63b910fSJames Wright     PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL));
69058600ac3SJames Wright     PetscCall(PetscLogStagePop());
69158600ac3SJames Wright   }
69240d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
69340d80af1SJames Wright }
69440d80af1SJames Wright 
69540d80af1SJames Wright /**
69640d80af1SJames Wright   @brief Assemble a `MATCEED` into a `MATAIJ` or similar.
69740d80af1SJames Wright          The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`.
69840d80af1SJames Wright          The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail.
69940d80af1SJames Wright 
70040d80af1SJames Wright   Collective across MPI processes.
70140d80af1SJames Wright 
70240d80af1SJames Wright   @param[in]      mat_ceed  `MATCEED` to assemble
70340d80af1SJames Wright   @param[in,out]  mat_coo   `MATAIJ` or similar to assemble into
70440d80af1SJames Wright 
70540d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
70640d80af1SJames Wright **/
70740d80af1SJames Wright PetscErrorCode MatCeedAssembleCOO(Mat mat_ceed, Mat mat_coo) {
70840d80af1SJames Wright   MatCeedContext ctx;
70940d80af1SJames Wright 
71040d80af1SJames Wright   PetscFunctionBeginUser;
71140d80af1SJames Wright   PetscCall(MatShellGetContext(mat_ceed, &ctx));
71240d80af1SJames Wright 
71340d80af1SJames Wright   // Set COO pattern if needed
71440d80af1SJames Wright   {
71540d80af1SJames Wright     CeedInt index = -1;
71640d80af1SJames Wright 
71740d80af1SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) {
71840d80af1SJames Wright       if (ctx->mats_assembled_full[i] == mat_coo) index = i;
71940d80af1SJames Wright     }
72040d80af1SJames Wright     if (index == -1) PetscCall(MatCeedSetPreallocationCOO(mat_ceed, mat_coo));
72158600ac3SJames Wright   }
72258600ac3SJames Wright 
72358600ac3SJames Wright   // Assemble mat_ceed
724c63b910fSJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL));
72558600ac3SJames Wright   PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY));
72658600ac3SJames Wright   {
72758600ac3SJames Wright     const CeedScalar *values;
72858600ac3SJames Wright     MatType           mat_type;
72958600ac3SJames Wright     CeedMemType       mem_type = CEED_MEM_HOST;
73058600ac3SJames Wright     PetscBool         is_spd, is_spd_known;
73158600ac3SJames Wright 
73258600ac3SJames Wright     PetscCall(MatGetType(mat_coo, &mat_type));
73358600ac3SJames Wright     if (strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE;
73458600ac3SJames Wright     else if (strstr(mat_type, "kokkos")) mem_type = CEED_MEM_DEVICE;
73558600ac3SJames Wright     else mem_type = CEED_MEM_HOST;
73658600ac3SJames Wright 
737c63b910fSJames Wright     PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
73850f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemble(ctx->op_mult, ctx->coo_values_full));
739c63b910fSJames Wright     PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL));
74050f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_full, mem_type, &values));
74158600ac3SJames Wright     PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES));
74258600ac3SJames Wright     PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd));
74358600ac3SJames Wright     if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd));
74450f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_full, &values));
74558600ac3SJames Wright   }
74658600ac3SJames Wright   PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY));
747c63b910fSJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL));
74858600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
74958600ac3SJames Wright }
75058600ac3SJames Wright 
75158600ac3SJames Wright /**
75240d80af1SJames Wright   @brief Set the current value of a context field for a `MatCEED`.
75340d80af1SJames Wright 
75440d80af1SJames Wright   Not collective across MPI processes.
75540d80af1SJames Wright 
75640d80af1SJames Wright   @param[in,out]  mat    `MatCEED`
75740d80af1SJames Wright   @param[in]      name   Name of the context field
75840d80af1SJames Wright   @param[in]      value  New context field value
75940d80af1SJames Wright 
76040d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
76140d80af1SJames Wright **/
76240d80af1SJames Wright PetscErrorCode MatCeedSetContextDouble(Mat mat, const char *name, double value) {
76340d80af1SJames Wright   PetscBool      was_updated = PETSC_FALSE;
76440d80af1SJames Wright   MatCeedContext ctx;
76540d80af1SJames Wright 
76640d80af1SJames Wright   PetscFunctionBeginUser;
76740d80af1SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
76840d80af1SJames Wright   {
76940d80af1SJames Wright     CeedContextFieldLabel label = NULL;
77040d80af1SJames Wright 
77140d80af1SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult, name, &label));
77240d80af1SJames Wright     if (label) {
773000d2032SJeremy L Thompson       double set_value = 2 * value + 1.0;
774000d2032SJeremy L Thompson 
775000d2032SJeremy L Thompson       PetscCall(MatCeedGetContextDouble(mat, name, &set_value));
776000d2032SJeremy L Thompson       if (set_value != value) {
77740d80af1SJames Wright         PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult, label, &value));
77840d80af1SJames Wright         was_updated = PETSC_TRUE;
77940d80af1SJames Wright       }
780000d2032SJeremy L Thompson     }
78140d80af1SJames Wright     if (ctx->op_mult_transpose) {
78240d80af1SJames Wright       label = NULL;
78340d80af1SJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult_transpose, name, &label));
78440d80af1SJames Wright       if (label) {
785000d2032SJeremy L Thompson         double set_value = 2 * value + 1.0;
786000d2032SJeremy L Thompson 
787000d2032SJeremy L Thompson         PetscCall(MatCeedGetContextDouble(mat, name, &set_value));
788000d2032SJeremy L Thompson         if (set_value != value) {
78940d80af1SJames Wright           PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult_transpose, label, &value));
79040d80af1SJames Wright           was_updated = PETSC_TRUE;
79140d80af1SJames Wright         }
79240d80af1SJames Wright       }
79340d80af1SJames Wright     }
794000d2032SJeremy L Thompson   }
79540d80af1SJames Wright   if (was_updated) {
79640d80af1SJames Wright     PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
79740d80af1SJames Wright     PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
79840d80af1SJames Wright   }
79940d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
80040d80af1SJames Wright }
80140d80af1SJames Wright 
80240d80af1SJames Wright /**
80340d80af1SJames Wright   @brief Get the current value of a context field for a `MatCEED`.
80440d80af1SJames Wright 
80540d80af1SJames Wright   Not collective across MPI processes.
80640d80af1SJames Wright 
80740d80af1SJames Wright   @param[in]   mat    `MatCEED`
80840d80af1SJames Wright   @param[in]   name   Name of the context field
80940d80af1SJames Wright   @param[out]  value  Current context field value
81040d80af1SJames Wright 
81140d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
81240d80af1SJames Wright **/
81340d80af1SJames Wright PetscErrorCode MatCeedGetContextDouble(Mat mat, const char *name, double *value) {
81440d80af1SJames Wright   MatCeedContext ctx;
81540d80af1SJames Wright 
81640d80af1SJames Wright   PetscFunctionBeginUser;
81740d80af1SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
81840d80af1SJames Wright   {
81940d80af1SJames Wright     CeedContextFieldLabel label = NULL;
82040d80af1SJames Wright     CeedOperator          op    = ctx->op_mult;
82140d80af1SJames Wright 
82240d80af1SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label));
82340d80af1SJames Wright     if (!label && ctx->op_mult_transpose) {
82440d80af1SJames Wright       op = ctx->op_mult_transpose;
82540d80af1SJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label));
82640d80af1SJames Wright     }
82740d80af1SJames Wright     if (label) {
82840d80af1SJames Wright       PetscSizeT    num_values;
82940d80af1SJames Wright       const double *values_ceed;
83040d80af1SJames Wright 
83140d80af1SJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorGetContextDoubleRead(op, label, &num_values, &values_ceed));
83240d80af1SJames Wright       *value = values_ceed[0];
83340d80af1SJames Wright       PetscCallCeed(ctx->ceed, CeedOperatorRestoreContextDoubleRead(op, label, &values_ceed));
83440d80af1SJames Wright     }
83540d80af1SJames Wright   }
83640d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
83740d80af1SJames Wright }
83840d80af1SJames Wright 
83940d80af1SJames Wright /**
840000d2032SJeremy L Thompson   @brief Set the current `PetscReal` value of a context field for a `MatCEED`.
841000d2032SJeremy L Thompson 
842000d2032SJeremy L Thompson   Not collective across MPI processes.
843000d2032SJeremy L Thompson 
844000d2032SJeremy L Thompson   @param[in,out]  mat    `MatCEED`
845000d2032SJeremy L Thompson   @param[in]      name   Name of the context field
846000d2032SJeremy L Thompson   @param[in]      value  New context field value
847000d2032SJeremy L Thompson 
848000d2032SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
849000d2032SJeremy L Thompson **/
850000d2032SJeremy L Thompson PetscErrorCode MatCeedSetContextReal(Mat mat, const char *name, PetscReal value) {
851000d2032SJeremy L Thompson   double value_double = value;
852000d2032SJeremy L Thompson 
853000d2032SJeremy L Thompson   PetscFunctionBeginUser;
854000d2032SJeremy L Thompson   PetscCall(MatCeedSetContextDouble(mat, name, value_double));
855000d2032SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
856000d2032SJeremy L Thompson }
857000d2032SJeremy L Thompson 
858000d2032SJeremy L Thompson /**
85951bb547fSJames Wright   @brief Get the current `PetscReal` value of a context field for a `MatCEED`.
86051bb547fSJames Wright 
86151bb547fSJames Wright   Not collective across MPI processes.
86251bb547fSJames Wright 
86351bb547fSJames Wright   @param[in]   mat    `MatCEED`
86451bb547fSJames Wright   @param[in]   name   Name of the context field
86551bb547fSJames Wright   @param[out]  value  Current context field value
86651bb547fSJames Wright 
86751bb547fSJames Wright   @return An error code: 0 - success, otherwise - failure
86851bb547fSJames Wright **/
86951bb547fSJames Wright PetscErrorCode MatCeedGetContextReal(Mat mat, const char *name, PetscReal *value) {
87087d3884fSJeremy L Thompson   double value_double = 0.0;
87151bb547fSJames Wright 
87251bb547fSJames Wright   PetscFunctionBeginUser;
87351bb547fSJames Wright   PetscCall(MatCeedGetContextDouble(mat, name, &value_double));
87451bb547fSJames Wright   *value = value_double;
87551bb547fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
87651bb547fSJames Wright }
87751bb547fSJames Wright 
87851bb547fSJames Wright /**
879000d2032SJeremy L Thompson   @brief Set the current time for a `MatCEED`.
880000d2032SJeremy L Thompson 
881000d2032SJeremy L Thompson   Not collective across MPI processes.
882000d2032SJeremy L Thompson 
883000d2032SJeremy L Thompson   @param[in,out]  mat   `MatCEED`
884000d2032SJeremy L Thompson   @param[in]      time  Current time
885000d2032SJeremy L Thompson 
886000d2032SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
887000d2032SJeremy L Thompson **/
888000d2032SJeremy L Thompson PetscErrorCode MatCeedSetTime(Mat mat, PetscReal time) {
889000d2032SJeremy L Thompson   PetscFunctionBeginUser;
890000d2032SJeremy L Thompson   {
891000d2032SJeremy L Thompson     double time_ceed = time;
892000d2032SJeremy L Thompson 
893000d2032SJeremy L Thompson     PetscCall(MatCeedSetContextDouble(mat, "time", time_ceed));
894000d2032SJeremy L Thompson   }
895000d2032SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
896000d2032SJeremy L Thompson }
897000d2032SJeremy L Thompson 
898000d2032SJeremy L Thompson /**
899000d2032SJeremy L Thompson   @brief Get the current time for a `MatCEED`.
900000d2032SJeremy L Thompson 
901000d2032SJeremy L Thompson   Not collective across MPI processes.
902000d2032SJeremy L Thompson 
903000d2032SJeremy L Thompson   @param[in]   mat   `MatCEED`
904000d2032SJeremy L Thompson   @param[out]  time  Current time, or -1.0 if the boundary evaluator has no time field
905000d2032SJeremy L Thompson 
906000d2032SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
907000d2032SJeremy L Thompson **/
908000d2032SJeremy L Thompson PetscErrorCode MatCeedGetTime(Mat mat, PetscReal *time) {
909000d2032SJeremy L Thompson   PetscFunctionBeginUser;
910000d2032SJeremy L Thompson   *time = -1.0;
911000d2032SJeremy L Thompson   {
912000d2032SJeremy L Thompson     double time_ceed = -1.0;
913000d2032SJeremy L Thompson 
914000d2032SJeremy L Thompson     PetscCall(MatCeedGetContextDouble(mat, "time", &time_ceed));
915000d2032SJeremy L Thompson     *time = time_ceed;
916000d2032SJeremy L Thompson   }
917000d2032SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
918000d2032SJeremy L Thompson }
919000d2032SJeremy L Thompson 
920000d2032SJeremy L Thompson /**
921000d2032SJeremy L Thompson   @brief Set the current time step for a `MatCEED`.
922000d2032SJeremy L Thompson 
923000d2032SJeremy L Thompson   Not collective across MPI processes.
924000d2032SJeremy L Thompson 
925000d2032SJeremy L Thompson   @param[in,out]  mat  `MatCEED`
926000d2032SJeremy L Thompson   @param[in]      dt   Current time step
927000d2032SJeremy L Thompson 
928000d2032SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
929000d2032SJeremy L Thompson **/
930000d2032SJeremy L Thompson PetscErrorCode MatCeedSetDt(Mat mat, PetscReal dt) {
931000d2032SJeremy L Thompson   PetscFunctionBeginUser;
932000d2032SJeremy L Thompson   {
933000d2032SJeremy L Thompson     double dt_ceed = dt;
934000d2032SJeremy L Thompson 
935000d2032SJeremy L Thompson     PetscCall(MatCeedSetContextDouble(mat, "dt", dt_ceed));
936000d2032SJeremy L Thompson   }
937000d2032SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
938000d2032SJeremy L Thompson }
939000d2032SJeremy L Thompson 
940000d2032SJeremy L Thompson /**
941000d2032SJeremy L Thompson   @brief Set the Jacobian shifts for a `MatCEED`.
942000d2032SJeremy L Thompson 
943000d2032SJeremy L Thompson   Not collective across MPI processes.
944000d2032SJeremy L Thompson 
945000d2032SJeremy L Thompson   @param[in,out]  mat      `MatCEED`
946000d2032SJeremy L Thompson   @param[in]      shift_v  Velocity shift
947000d2032SJeremy L Thompson   @param[in]      shift_a  Acceleration shift
948000d2032SJeremy L Thompson 
949000d2032SJeremy L Thompson   @return An error code: 0 - success, otherwise - failure
950000d2032SJeremy L Thompson **/
951000d2032SJeremy L Thompson PetscErrorCode MatCeedSetShifts(Mat mat, PetscReal shift_v, PetscReal shift_a) {
952000d2032SJeremy L Thompson   PetscFunctionBeginUser;
953000d2032SJeremy L Thompson   {
954000d2032SJeremy L Thompson     double shift_v_ceed = shift_v;
955000d2032SJeremy L Thompson 
956000d2032SJeremy L Thompson     PetscCall(MatCeedSetContextDouble(mat, "shift v", shift_v_ceed));
957000d2032SJeremy L Thompson   }
958000d2032SJeremy L Thompson   if (shift_a) {
959000d2032SJeremy L Thompson     double shift_a_ceed = shift_a;
960000d2032SJeremy L Thompson 
961000d2032SJeremy L Thompson     PetscCall(MatCeedSetContextDouble(mat, "shift a", shift_a_ceed));
962000d2032SJeremy L Thompson   }
963000d2032SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
964000d2032SJeremy L Thompson }
965000d2032SJeremy L Thompson 
966000d2032SJeremy L Thompson /**
96758600ac3SJames Wright   @brief Set user context for a `MATCEED`.
96858600ac3SJames Wright 
96958600ac3SJames Wright   Collective across MPI processes.
97058600ac3SJames Wright 
97158600ac3SJames Wright   @param[in,out]  mat  `MATCEED`
97258600ac3SJames Wright   @param[in]      f    The context destroy function, or NULL
97358600ac3SJames Wright   @param[in]      ctx  User context, or NULL to unset
97458600ac3SJames Wright 
97558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
97658600ac3SJames Wright **/
97758600ac3SJames Wright PetscErrorCode MatCeedSetContext(Mat mat, PetscErrorCode (*f)(void *), void *ctx) {
97858600ac3SJames Wright   PetscContainer user_ctx = NULL;
97958600ac3SJames Wright 
98058600ac3SJames Wright   PetscFunctionBeginUser;
98158600ac3SJames Wright   if (ctx) {
98258600ac3SJames Wright     PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)mat), &user_ctx));
98358600ac3SJames Wright     PetscCall(PetscContainerSetPointer(user_ctx, ctx));
98458600ac3SJames Wright     PetscCall(PetscContainerSetUserDestroy(user_ctx, f));
98558600ac3SJames Wright   }
98658600ac3SJames Wright   PetscCall(PetscObjectCompose((PetscObject)mat, "MatCeed user context", (PetscObject)user_ctx));
98758600ac3SJames Wright   PetscCall(PetscContainerDestroy(&user_ctx));
98858600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
98958600ac3SJames Wright }
99058600ac3SJames Wright 
99158600ac3SJames Wright /**
99258600ac3SJames Wright   @brief Retrieve the user context for a `MATCEED`.
99358600ac3SJames Wright 
99458600ac3SJames Wright   Collective across MPI processes.
99558600ac3SJames Wright 
99658600ac3SJames Wright   @param[in,out]  mat  `MATCEED`
99758600ac3SJames Wright   @param[in]      ctx  User context
99858600ac3SJames Wright 
99958600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
100058600ac3SJames Wright **/
100158600ac3SJames Wright PetscErrorCode MatCeedGetContext(Mat mat, void *ctx) {
100258600ac3SJames Wright   PetscContainer user_ctx;
100358600ac3SJames Wright 
100458600ac3SJames Wright   PetscFunctionBeginUser;
100558600ac3SJames Wright   PetscCall(PetscObjectQuery((PetscObject)mat, "MatCeed user context", (PetscObject *)&user_ctx));
100658600ac3SJames Wright   if (user_ctx) PetscCall(PetscContainerGetPointer(user_ctx, (void **)ctx));
100758600ac3SJames Wright   else *(void **)ctx = NULL;
100858600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
100958600ac3SJames Wright }
101058600ac3SJames Wright /**
101140d80af1SJames Wright   @brief Set a user defined matrix operation for a `MATCEED` matrix.
101240d80af1SJames Wright 
101340d80af1SJames Wright   Within each user-defined routine, the user should call `MatCeedGetContext()` to obtain the user-defined context that was set by
101440d80af1SJames Wright `MatCeedSetContext()`.
101558600ac3SJames Wright 
101658600ac3SJames Wright   Collective across MPI processes.
101758600ac3SJames Wright 
101858600ac3SJames Wright   @param[in,out]  mat  `MATCEED`
101940d80af1SJames Wright   @param[in]      op   Name of the `MatOperation`
102040d80af1SJames Wright   @param[in]      g    Function that provides the operation
102158600ac3SJames Wright 
102258600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
102358600ac3SJames Wright **/
102440d80af1SJames Wright PetscErrorCode MatCeedSetOperation(Mat mat, MatOperation op, void (*g)(void)) {
102540d80af1SJames Wright   PetscFunctionBeginUser;
102640d80af1SJames Wright   PetscCall(MatShellSetOperation(mat, op, g));
102740d80af1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
102840d80af1SJames Wright }
102940d80af1SJames Wright 
103040d80af1SJames Wright /**
103140d80af1SJames Wright   @brief Sets the default COO matrix type as a string from the `MATCEED`.
103240d80af1SJames Wright 
103340d80af1SJames Wright   Collective across MPI processes.
103440d80af1SJames Wright 
103540d80af1SJames Wright   @param[in,out]  mat   `MATCEED`
103640d80af1SJames Wright   @param[in]      type  COO `MatType` to set
103740d80af1SJames Wright 
103840d80af1SJames Wright   @return An error code: 0 - success, otherwise - failure
103940d80af1SJames Wright **/
104040d80af1SJames Wright PetscErrorCode MatCeedSetCOOMatType(Mat mat, MatType type) {
104158600ac3SJames Wright   MatCeedContext ctx;
104258600ac3SJames Wright 
104358600ac3SJames Wright   PetscFunctionBeginUser;
104458600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
104558600ac3SJames Wright   // Check if same
104658600ac3SJames Wright   {
104758600ac3SJames Wright     size_t    len_old, len_new;
104858600ac3SJames Wright     PetscBool is_same = PETSC_FALSE;
104958600ac3SJames Wright 
105040d80af1SJames Wright     PetscCall(PetscStrlen(ctx->coo_mat_type, &len_old));
105158600ac3SJames Wright     PetscCall(PetscStrlen(type, &len_new));
105240d80af1SJames Wright     if (len_old == len_new) PetscCall(PetscStrncmp(ctx->coo_mat_type, type, len_old, &is_same));
105358600ac3SJames Wright     if (is_same) PetscFunctionReturn(PETSC_SUCCESS);
105458600ac3SJames Wright   }
105558600ac3SJames Wright   // Clean up old mats in different format
105658600ac3SJames Wright   // LCOV_EXCL_START
105758600ac3SJames Wright   if (ctx->mat_assembled_full_internal) {
105858600ac3SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) {
105958600ac3SJames Wright       if (ctx->mats_assembled_full[i] == ctx->mat_assembled_full_internal) {
106058600ac3SJames Wright         for (PetscInt j = i + 1; j < ctx->num_mats_assembled_full; j++) {
106158600ac3SJames Wright           ctx->mats_assembled_full[j - 1] = ctx->mats_assembled_full[j];
106258600ac3SJames Wright         }
106358600ac3SJames Wright         ctx->num_mats_assembled_full--;
106458600ac3SJames Wright         // Note: we'll realloc this array again, so no need to shrink the allocation
106558600ac3SJames Wright         PetscCall(MatDestroy(&ctx->mat_assembled_full_internal));
106658600ac3SJames Wright       }
106758600ac3SJames Wright     }
106858600ac3SJames Wright   }
106958600ac3SJames Wright   if (ctx->mat_assembled_pbd_internal) {
107058600ac3SJames Wright     for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) {
107158600ac3SJames Wright       if (ctx->mats_assembled_pbd[i] == ctx->mat_assembled_pbd_internal) {
107258600ac3SJames Wright         for (PetscInt j = i + 1; j < ctx->num_mats_assembled_pbd; j++) {
107358600ac3SJames Wright           ctx->mats_assembled_pbd[j - 1] = ctx->mats_assembled_pbd[j];
107458600ac3SJames Wright         }
107558600ac3SJames Wright         // Note: we'll realloc this array again, so no need to shrink the allocation
107658600ac3SJames Wright         ctx->num_mats_assembled_pbd--;
107758600ac3SJames Wright         PetscCall(MatDestroy(&ctx->mat_assembled_pbd_internal));
107858600ac3SJames Wright       }
107958600ac3SJames Wright     }
108058600ac3SJames Wright   }
108140d80af1SJames Wright   PetscCall(PetscFree(ctx->coo_mat_type));
108240d80af1SJames Wright   PetscCall(PetscStrallocpy(type, &ctx->coo_mat_type));
108358600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
108458600ac3SJames Wright   // LCOV_EXCL_STOP
108558600ac3SJames Wright }
108658600ac3SJames Wright 
108758600ac3SJames Wright /**
108840d80af1SJames Wright   @brief Gets the default COO matrix type as a string from the `MATCEED`.
108958600ac3SJames Wright 
109058600ac3SJames Wright   Collective across MPI processes.
109158600ac3SJames Wright 
109258600ac3SJames Wright   @param[in,out]  mat   `MATCEED`
109340d80af1SJames Wright   @param[in]      type  COO `MatType`
109458600ac3SJames Wright 
109558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
109658600ac3SJames Wright **/
109740d80af1SJames Wright PetscErrorCode MatCeedGetCOOMatType(Mat mat, MatType *type) {
109858600ac3SJames Wright   MatCeedContext ctx;
109958600ac3SJames Wright 
110058600ac3SJames Wright   PetscFunctionBeginUser;
110158600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
110240d80af1SJames Wright   *type = ctx->coo_mat_type;
110358600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
110458600ac3SJames Wright }
110558600ac3SJames Wright 
110658600ac3SJames Wright /**
110758600ac3SJames Wright   @brief Set input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
110858600ac3SJames Wright 
110958600ac3SJames Wright   Not collective across MPI processes.
111058600ac3SJames Wright 
111158600ac3SJames Wright   @param[in,out]  mat              `MATCEED`
111258600ac3SJames Wright   @param[in]      X_loc            Input PETSc local vector, or NULL
111358600ac3SJames Wright   @param[in]      Y_loc_transpose  Input PETSc local vector for transpose operation, or NULL
111458600ac3SJames Wright 
111558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
111658600ac3SJames Wright **/
111758600ac3SJames Wright PetscErrorCode MatCeedSetLocalVectors(Mat mat, Vec X_loc, Vec Y_loc_transpose) {
111858600ac3SJames Wright   MatCeedContext ctx;
111958600ac3SJames Wright 
112058600ac3SJames Wright   PetscFunctionBeginUser;
112158600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
112258600ac3SJames Wright   if (X_loc) {
112358600ac3SJames Wright     PetscInt len_old, len_new;
112458600ac3SJames Wright 
112558600ac3SJames Wright     PetscCall(VecGetSize(ctx->X_loc, &len_old));
112658600ac3SJames Wright     PetscCall(VecGetSize(X_loc, &len_new));
112758600ac3SJames 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,
112858600ac3SJames Wright                len_new, len_old);
112940d80af1SJames Wright     PetscCall(VecReferenceCopy(X_loc, &ctx->X_loc));
113058600ac3SJames Wright   }
113158600ac3SJames Wright   if (Y_loc_transpose) {
113258600ac3SJames Wright     PetscInt len_old, len_new;
113358600ac3SJames Wright 
113458600ac3SJames Wright     PetscCall(VecGetSize(ctx->Y_loc_transpose, &len_old));
113558600ac3SJames Wright     PetscCall(VecGetSize(Y_loc_transpose, &len_new));
113658600ac3SJames Wright     PetscCheck(len_old == len_new, PETSC_COMM_SELF, PETSC_ERR_LIB,
113758600ac3SJames Wright                "new Y_loc_transpose length %" PetscInt_FMT " should match old Y_loc_transpose length %" PetscInt_FMT, len_new, len_old);
113840d80af1SJames Wright     PetscCall(VecReferenceCopy(Y_loc_transpose, &ctx->Y_loc_transpose));
113958600ac3SJames Wright   }
114058600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
114158600ac3SJames Wright }
114258600ac3SJames Wright 
114358600ac3SJames Wright /**
114458600ac3SJames Wright   @brief Get input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
114558600ac3SJames Wright 
114658600ac3SJames Wright   Not collective across MPI processes.
114758600ac3SJames Wright 
114858600ac3SJames Wright   @param[in,out]  mat              `MATCEED`
114958600ac3SJames Wright   @param[out]     X_loc            Input PETSc local vector, or NULL
115058600ac3SJames Wright   @param[out]     Y_loc_transpose  Input PETSc local vector for transpose operation, or NULL
115158600ac3SJames Wright 
115258600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
115358600ac3SJames Wright **/
115458600ac3SJames Wright PetscErrorCode MatCeedGetLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) {
115558600ac3SJames Wright   MatCeedContext ctx;
115658600ac3SJames Wright 
115758600ac3SJames Wright   PetscFunctionBeginUser;
115858600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
115958600ac3SJames Wright   if (X_loc) {
116040d80af1SJames Wright     *X_loc = NULL;
116140d80af1SJames Wright     PetscCall(VecReferenceCopy(ctx->X_loc, X_loc));
116258600ac3SJames Wright   }
116358600ac3SJames Wright   if (Y_loc_transpose) {
116440d80af1SJames Wright     *Y_loc_transpose = NULL;
116540d80af1SJames Wright     PetscCall(VecReferenceCopy(ctx->Y_loc_transpose, Y_loc_transpose));
116658600ac3SJames Wright   }
116758600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
116858600ac3SJames Wright }
116958600ac3SJames Wright 
117058600ac3SJames Wright /**
117158600ac3SJames Wright   @brief Restore input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
117258600ac3SJames Wright 
117358600ac3SJames Wright   Not collective across MPI processes.
117458600ac3SJames Wright 
117558600ac3SJames Wright   @param[in,out]  mat              MatCeed
117658600ac3SJames Wright   @param[out]     X_loc            Input PETSc local vector, or NULL
117758600ac3SJames Wright   @param[out]     Y_loc_transpose  Input PETSc local vector for transpose operation, or NULL
117858600ac3SJames Wright 
117958600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
118058600ac3SJames Wright **/
118158600ac3SJames Wright PetscErrorCode MatCeedRestoreLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) {
118258600ac3SJames Wright   PetscFunctionBeginUser;
118358600ac3SJames Wright   if (X_loc) PetscCall(VecDestroy(X_loc));
118458600ac3SJames Wright   if (Y_loc_transpose) PetscCall(VecDestroy(Y_loc_transpose));
118558600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
118658600ac3SJames Wright }
118758600ac3SJames Wright 
118858600ac3SJames Wright /**
118958600ac3SJames Wright   @brief Get libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
119058600ac3SJames Wright 
119158600ac3SJames Wright   Not collective across MPI processes.
119258600ac3SJames Wright 
119358600ac3SJames Wright   @param[in,out]  mat                MatCeed
119458600ac3SJames Wright   @param[out]     op_mult            libCEED `CeedOperator` for `MatMult()`, or NULL
119558600ac3SJames Wright   @param[out]     op_mult_transpose  libCEED `CeedOperator` for `MatMultTranspose()`, or NULL
119658600ac3SJames Wright 
119758600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
119858600ac3SJames Wright **/
119958600ac3SJames Wright PetscErrorCode MatCeedGetCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) {
120058600ac3SJames Wright   MatCeedContext ctx;
120158600ac3SJames Wright 
120258600ac3SJames Wright   PetscFunctionBeginUser;
120358600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
120458600ac3SJames Wright   if (op_mult) {
120558600ac3SJames Wright     *op_mult = NULL;
120650f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult, op_mult));
120758600ac3SJames Wright   }
120858600ac3SJames Wright   if (op_mult_transpose) {
120958600ac3SJames Wright     *op_mult_transpose = NULL;
121050f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult_transpose, op_mult_transpose));
121158600ac3SJames Wright   }
121258600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
121358600ac3SJames Wright }
121458600ac3SJames Wright 
121558600ac3SJames Wright /**
121658600ac3SJames Wright   @brief Restore libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations.
121758600ac3SJames Wright 
121858600ac3SJames Wright   Not collective across MPI processes.
121958600ac3SJames Wright 
122058600ac3SJames Wright   @param[in,out]  mat                MatCeed
122158600ac3SJames Wright   @param[out]     op_mult            libCEED `CeedOperator` for `MatMult()`, or NULL
122258600ac3SJames Wright   @param[out]     op_mult_transpose  libCEED `CeedOperator` for `MatMultTranspose()`, or NULL
122358600ac3SJames Wright 
122458600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
122558600ac3SJames Wright **/
122658600ac3SJames Wright PetscErrorCode MatCeedRestoreCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) {
122758600ac3SJames Wright   MatCeedContext ctx;
122858600ac3SJames Wright 
122958600ac3SJames Wright   PetscFunctionBeginUser;
123058600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
123150f50432SJames Wright   if (op_mult) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult));
123250f50432SJames Wright   if (op_mult_transpose) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult_transpose));
123358600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
123458600ac3SJames Wright }
123558600ac3SJames Wright 
123658600ac3SJames Wright /**
123758600ac3SJames Wright   @brief Set `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
123858600ac3SJames Wright 
123958600ac3SJames Wright   Not collective across MPI processes.
124058600ac3SJames Wright 
124158600ac3SJames Wright   @param[in,out]  mat                       MatCeed
124258600ac3SJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward evaluation, or NULL
124358600ac3SJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose evaluation, or NULL
124458600ac3SJames Wright 
124558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
124658600ac3SJames Wright **/
124758600ac3SJames Wright PetscErrorCode MatCeedSetLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) {
124858600ac3SJames Wright   MatCeedContext ctx;
124958600ac3SJames Wright 
125058600ac3SJames Wright   PetscFunctionBeginUser;
125158600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
125258600ac3SJames Wright   if (log_event_mult) ctx->log_event_mult = log_event_mult;
125358600ac3SJames Wright   if (log_event_mult_transpose) ctx->log_event_mult_transpose = log_event_mult_transpose;
125458600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
125558600ac3SJames Wright }
125658600ac3SJames Wright 
125758600ac3SJames Wright /**
125858600ac3SJames Wright   @brief Get `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
125958600ac3SJames Wright 
126058600ac3SJames Wright   Not collective across MPI processes.
126158600ac3SJames Wright 
126258600ac3SJames Wright   @param[in,out]  mat                       MatCeed
126358600ac3SJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward evaluation, or NULL
126458600ac3SJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose evaluation, or NULL
126558600ac3SJames Wright 
126658600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
126758600ac3SJames Wright **/
126858600ac3SJames Wright PetscErrorCode MatCeedGetLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) {
126958600ac3SJames Wright   MatCeedContext ctx;
127058600ac3SJames Wright 
127158600ac3SJames Wright   PetscFunctionBeginUser;
127258600ac3SJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
127358600ac3SJames Wright   if (log_event_mult) *log_event_mult = ctx->log_event_mult;
127458600ac3SJames Wright   if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_mult_transpose;
127558600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
127658600ac3SJames Wright }
127758600ac3SJames Wright 
1278c63b910fSJames Wright /**
1279c63b910fSJames Wright   @brief Set `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
1280c63b910fSJames Wright 
1281c63b910fSJames Wright   Not collective across MPI processes.
1282c63b910fSJames Wright 
1283c63b910fSJames Wright   @param[in,out]  mat                       MatCeed
1284c63b910fSJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL
1285c63b910fSJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL
1286c63b910fSJames Wright 
1287c63b910fSJames Wright   @return An error code: 0 - success, otherwise - failure
1288c63b910fSJames Wright **/
1289c63b910fSJames Wright PetscErrorCode MatCeedSetCeedOperatorLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) {
1290c63b910fSJames Wright   MatCeedContext ctx;
1291c63b910fSJames Wright 
1292c63b910fSJames Wright   PetscFunctionBeginUser;
1293c63b910fSJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1294c63b910fSJames Wright   if (log_event_mult) ctx->log_event_ceed_mult = log_event_mult;
1295c63b910fSJames Wright   if (log_event_mult_transpose) ctx->log_event_ceed_mult_transpose = log_event_mult_transpose;
1296c63b910fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1297c63b910fSJames Wright }
1298c63b910fSJames Wright 
1299c63b910fSJames Wright /**
1300c63b910fSJames Wright   @brief Get `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators.
1301c63b910fSJames Wright 
1302c63b910fSJames Wright   Not collective across MPI processes.
1303c63b910fSJames Wright 
1304c63b910fSJames Wright   @param[in,out]  mat                       MatCeed
1305c63b910fSJames Wright   @param[out]     log_event_mult            `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL
1306c63b910fSJames Wright   @param[out]     log_event_mult_transpose  `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL
1307c63b910fSJames Wright 
1308c63b910fSJames Wright   @return An error code: 0 - success, otherwise - failure
1309c63b910fSJames Wright **/
1310c63b910fSJames Wright PetscErrorCode MatCeedGetCeedOperatorLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) {
1311c63b910fSJames Wright   MatCeedContext ctx;
1312c63b910fSJames Wright 
1313c63b910fSJames Wright   PetscFunctionBeginUser;
1314c63b910fSJames Wright   PetscCall(MatShellGetContext(mat, &ctx));
1315c63b910fSJames Wright   if (log_event_mult) *log_event_mult = ctx->log_event_ceed_mult;
1316c63b910fSJames Wright   if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_ceed_mult_transpose;
1317c63b910fSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1318c63b910fSJames Wright }
1319c63b910fSJames Wright 
132058600ac3SJames Wright // -----------------------------------------------------------------------------
132158600ac3SJames Wright // Operator context data
132258600ac3SJames Wright // -----------------------------------------------------------------------------
132358600ac3SJames Wright 
132458600ac3SJames Wright /**
132558600ac3SJames Wright   @brief Setup context data for operator application.
132658600ac3SJames Wright 
132758600ac3SJames Wright   Collective across MPI processes.
132858600ac3SJames Wright 
132958600ac3SJames Wright   @param[in]   dm_x                           Input `DM`
133058600ac3SJames Wright   @param[in]   dm_y                           Output `DM`
133158600ac3SJames Wright   @param[in]   X_loc                          Input PETSc local vector, or NULL
133258600ac3SJames Wright   @param[in]   Y_loc_transpose                Input PETSc local vector for transpose operation, or NULL
133358600ac3SJames Wright   @param[in]   op_mult                        `CeedOperator` for forward evaluation
133458600ac3SJames Wright   @param[in]   op_mult_transpose              `CeedOperator` for transpose evaluation
133558600ac3SJames Wright   @param[in]   log_event_mult                 `PetscLogEvent` for forward evaluation
133658600ac3SJames Wright   @param[in]   log_event_mult_transpose       `PetscLogEvent` for transpose evaluation
1337c63b910fSJames Wright   @param[in]   log_event_ceed_mult            `PetscLogEvent` for forward `CeedOperator` evaluation
1338c63b910fSJames Wright   @param[in]   log_event_ceed_mult_transpose  `PetscLogEvent` for transpose `CeedOperator` evaluation
133958600ac3SJames Wright   @param[out]  ctx                            Context data for operator evaluation
134058600ac3SJames Wright 
134158600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
134258600ac3SJames Wright **/
134358600ac3SJames Wright PetscErrorCode MatCeedContextCreate(DM dm_x, DM dm_y, Vec X_loc, Vec Y_loc_transpose, CeedOperator op_mult, CeedOperator op_mult_transpose,
1344c63b910fSJames Wright                                     PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose, PetscLogEvent log_event_ceed_mult,
1345c63b910fSJames Wright                                     PetscLogEvent log_event_ceed_mult_transpose, MatCeedContext *ctx) {
134658600ac3SJames Wright   CeedSize x_loc_len, y_loc_len;
134758600ac3SJames Wright 
134858600ac3SJames Wright   PetscFunctionBeginUser;
134958600ac3SJames Wright 
135058600ac3SJames Wright   // Allocate
135158600ac3SJames Wright   PetscCall(PetscNew(ctx));
135258600ac3SJames Wright   (*ctx)->ref_count = 1;
135358600ac3SJames Wright 
135458600ac3SJames Wright   // Logging
135558600ac3SJames Wright   (*ctx)->log_event_mult                = log_event_mult;
135658600ac3SJames Wright   (*ctx)->log_event_mult_transpose      = log_event_mult_transpose;
1357c63b910fSJames Wright   (*ctx)->log_event_ceed_mult           = log_event_ceed_mult;
1358c63b910fSJames Wright   (*ctx)->log_event_ceed_mult_transpose = log_event_ceed_mult_transpose;
135958600ac3SJames Wright 
136058600ac3SJames Wright   // PETSc objects
136140d80af1SJames Wright   PetscCall(DMReferenceCopy(dm_x, &(*ctx)->dm_x));
136240d80af1SJames Wright   PetscCall(DMReferenceCopy(dm_y, &(*ctx)->dm_y));
136340d80af1SJames Wright   if (X_loc) PetscCall(VecReferenceCopy(X_loc, &(*ctx)->X_loc));
136440d80af1SJames Wright   if (Y_loc_transpose) PetscCall(VecReferenceCopy(Y_loc_transpose, &(*ctx)->Y_loc_transpose));
136558600ac3SJames Wright 
136658600ac3SJames Wright   // Memtype
136758600ac3SJames Wright   {
136858600ac3SJames Wright     const PetscScalar *x;
136958600ac3SJames Wright     Vec                X;
137058600ac3SJames Wright 
137158600ac3SJames Wright     PetscCall(DMGetLocalVector(dm_x, &X));
137258600ac3SJames Wright     PetscCall(VecGetArrayReadAndMemType(X, &x, &(*ctx)->mem_type));
137358600ac3SJames Wright     PetscCall(VecRestoreArrayReadAndMemType(X, &x));
137458600ac3SJames Wright     PetscCall(DMRestoreLocalVector(dm_x, &X));
137558600ac3SJames Wright   }
137658600ac3SJames Wright 
137758600ac3SJames Wright   // libCEED objects
137858600ac3SJames Wright   PetscCheck(CeedOperatorGetCeed(op_mult, &(*ctx)->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB,
137958600ac3SJames Wright              "retrieving Ceed context object failed");
138050f50432SJames Wright   PetscCallCeed((*ctx)->ceed, CeedReference((*ctx)->ceed));
138150f50432SJames Wright   PetscCallCeed((*ctx)->ceed, CeedOperatorGetActiveVectorLengths(op_mult, &x_loc_len, &y_loc_len));
138250f50432SJames Wright   PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult, &(*ctx)->op_mult));
138350f50432SJames Wright   if (op_mult_transpose) PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult_transpose, &(*ctx)->op_mult_transpose));
138450f50432SJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, x_loc_len, &(*ctx)->x_loc));
138550f50432SJames Wright   PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, y_loc_len, &(*ctx)->y_loc));
138658600ac3SJames Wright 
138758600ac3SJames Wright   // Flop counting
138858600ac3SJames Wright   {
138958600ac3SJames Wright     CeedSize ceed_flops_estimate = 0;
139058600ac3SJames Wright 
139150f50432SJames Wright     PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult, &ceed_flops_estimate));
139258600ac3SJames Wright     (*ctx)->flops_mult = ceed_flops_estimate;
139358600ac3SJames Wright     if (op_mult_transpose) {
139450f50432SJames Wright       PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult_transpose, &ceed_flops_estimate));
139558600ac3SJames Wright       (*ctx)->flops_mult_transpose = ceed_flops_estimate;
139658600ac3SJames Wright     }
139758600ac3SJames Wright   }
139858600ac3SJames Wright 
139958600ac3SJames Wright   // Check sizes
140058600ac3SJames Wright   if (x_loc_len > 0 || y_loc_len > 0) {
140158600ac3SJames Wright     CeedSize ctx_x_loc_len, ctx_y_loc_len;
140258600ac3SJames Wright     PetscInt X_loc_len, dm_x_loc_len, Y_loc_len, dm_y_loc_len;
140358600ac3SJames Wright     Vec      dm_X_loc, dm_Y_loc;
140458600ac3SJames Wright 
140558600ac3SJames Wright     // -- Input
140658600ac3SJames Wright     PetscCall(DMGetLocalVector(dm_x, &dm_X_loc));
140758600ac3SJames Wright     PetscCall(VecGetLocalSize(dm_X_loc, &dm_x_loc_len));
140858600ac3SJames Wright     PetscCall(DMRestoreLocalVector(dm_x, &dm_X_loc));
140950f50432SJames Wright     PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->x_loc, &ctx_x_loc_len));
14104c17272bSJames Wright     if (X_loc) {
14114c17272bSJames Wright       PetscCall(VecGetLocalSize(X_loc, &X_loc_len));
14124c17272bSJames Wright       PetscCheck(X_loc_len == dm_x_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB,
14134c17272bSJames Wright                  "X_loc (%" PetscInt_FMT ") must match dm_x (%" PetscInt_FMT ") dimensions", X_loc_len, dm_x_loc_len);
14144c17272bSJames Wright     }
14154c17272bSJames 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",
14164c17272bSJames Wright                x_loc_len, dm_x_loc_len);
14174c17272bSJames 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 ")",
14184c17272bSJames Wright                x_loc_len, ctx_x_loc_len);
141958600ac3SJames Wright 
142058600ac3SJames Wright     // -- Output
142158600ac3SJames Wright     PetscCall(DMGetLocalVector(dm_y, &dm_Y_loc));
142258600ac3SJames Wright     PetscCall(VecGetLocalSize(dm_Y_loc, &dm_y_loc_len));
142358600ac3SJames Wright     PetscCall(DMRestoreLocalVector(dm_y, &dm_Y_loc));
142450f50432SJames Wright     PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->y_loc, &ctx_y_loc_len));
14254c17272bSJames 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",
14264c17272bSJames Wright                ctx_y_loc_len, dm_y_loc_len);
142758600ac3SJames Wright 
142858600ac3SJames Wright     // -- Transpose
142958600ac3SJames Wright     if (Y_loc_transpose) {
143058600ac3SJames Wright       PetscCall(VecGetLocalSize(Y_loc_transpose, &Y_loc_len));
14314c17272bSJames Wright       PetscCheck(Y_loc_len == dm_y_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB,
14324c17272bSJames Wright                  "Y_loc_transpose (%" PetscInt_FMT ") must match dm_y (%" PetscInt_FMT ") dimensions", Y_loc_len, dm_y_loc_len);
143358600ac3SJames Wright     }
143458600ac3SJames Wright   }
143558600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
143658600ac3SJames Wright }
143758600ac3SJames Wright 
143858600ac3SJames Wright /**
143958600ac3SJames Wright   @brief Increment reference counter for `MATCEED` context.
144058600ac3SJames Wright 
144158600ac3SJames Wright   Not collective across MPI processes.
144258600ac3SJames Wright 
144358600ac3SJames Wright   @param[in,out]  ctx  Context data
144458600ac3SJames Wright 
144558600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
144658600ac3SJames Wright **/
144758600ac3SJames Wright PetscErrorCode MatCeedContextReference(MatCeedContext ctx) {
144858600ac3SJames Wright   PetscFunctionBeginUser;
144958600ac3SJames Wright   ctx->ref_count++;
145058600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
145158600ac3SJames Wright }
145258600ac3SJames Wright 
145358600ac3SJames Wright /**
145458600ac3SJames Wright   @brief Copy reference for `MATCEED`.
145558600ac3SJames Wright          Note: If `ctx_copy` is non-null, it is assumed to be a valid pointer to a `MatCeedContext`.
145658600ac3SJames Wright 
145758600ac3SJames Wright   Not collective across MPI processes.
145858600ac3SJames Wright 
145958600ac3SJames Wright   @param[in]   ctx       Context data
146058600ac3SJames Wright   @param[out]  ctx_copy  Copy of pointer to context data
146158600ac3SJames Wright 
146258600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
146358600ac3SJames Wright **/
146458600ac3SJames Wright PetscErrorCode MatCeedContextReferenceCopy(MatCeedContext ctx, MatCeedContext *ctx_copy) {
146558600ac3SJames Wright   PetscFunctionBeginUser;
146658600ac3SJames Wright   PetscCall(MatCeedContextReference(ctx));
146758600ac3SJames Wright   PetscCall(MatCeedContextDestroy(*ctx_copy));
146858600ac3SJames Wright   *ctx_copy = ctx;
146958600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
147058600ac3SJames Wright }
147158600ac3SJames Wright 
147258600ac3SJames Wright /**
147358600ac3SJames Wright   @brief Destroy context data for operator application.
147458600ac3SJames Wright 
147558600ac3SJames Wright   Collective across MPI processes.
147658600ac3SJames Wright 
147758600ac3SJames Wright   @param[in,out]  ctx  Context data for operator evaluation
147858600ac3SJames Wright 
147958600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
148058600ac3SJames Wright **/
148158600ac3SJames Wright PetscErrorCode MatCeedContextDestroy(MatCeedContext ctx) {
148258600ac3SJames Wright   PetscFunctionBeginUser;
148358600ac3SJames Wright   if (!ctx || --ctx->ref_count > 0) PetscFunctionReturn(PETSC_SUCCESS);
148458600ac3SJames Wright 
148558600ac3SJames Wright   // PETSc objects
148658600ac3SJames Wright   PetscCall(DMDestroy(&ctx->dm_x));
148758600ac3SJames Wright   PetscCall(DMDestroy(&ctx->dm_y));
148858600ac3SJames Wright   PetscCall(VecDestroy(&ctx->X_loc));
148958600ac3SJames Wright   PetscCall(VecDestroy(&ctx->Y_loc_transpose));
149058600ac3SJames Wright   PetscCall(MatDestroy(&ctx->mat_assembled_full_internal));
149158600ac3SJames Wright   PetscCall(MatDestroy(&ctx->mat_assembled_pbd_internal));
149240d80af1SJames Wright   PetscCall(PetscFree(ctx->coo_mat_type));
149358600ac3SJames Wright   PetscCall(PetscFree(ctx->mats_assembled_full));
149458600ac3SJames Wright   PetscCall(PetscFree(ctx->mats_assembled_pbd));
149558600ac3SJames Wright 
149658600ac3SJames Wright   // libCEED objects
149750f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->x_loc));
149850f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->y_loc));
149950f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->coo_values_full));
150050f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->coo_values_pbd));
150150f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedOperatorDestroy(&ctx->op_mult));
150250f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedOperatorDestroy(&ctx->op_mult_transpose));
150358600ac3SJames Wright   PetscCheck(CeedDestroy(&ctx->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "destroying libCEED context object failed");
150458600ac3SJames Wright 
150558600ac3SJames Wright   // Deallocate
150658600ac3SJames Wright   ctx->is_destroyed = PETSC_TRUE;  // Flag as destroyed in case someone has stale ref
150758600ac3SJames Wright   PetscCall(PetscFree(ctx));
150858600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
150958600ac3SJames Wright }
151058600ac3SJames Wright 
151158600ac3SJames Wright /**
151258600ac3SJames Wright   @brief Compute the diagonal of an operator via libCEED.
151358600ac3SJames Wright 
151458600ac3SJames Wright   Collective across MPI processes.
151558600ac3SJames Wright 
151658600ac3SJames Wright   @param[in]   A  `MATCEED`
151758600ac3SJames Wright   @param[out]  D  Vector holding operator diagonal
151858600ac3SJames Wright 
151958600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
152058600ac3SJames Wright **/
152158600ac3SJames Wright PetscErrorCode MatGetDiagonal_Ceed(Mat A, Vec D) {
152258600ac3SJames Wright   PetscMemType   mem_type;
152358600ac3SJames Wright   Vec            D_loc;
152458600ac3SJames Wright   MatCeedContext ctx;
152558600ac3SJames Wright 
152658600ac3SJames Wright   PetscFunctionBeginUser;
152758600ac3SJames Wright   PetscCall(MatShellGetContext(A, &ctx));
152858600ac3SJames Wright 
152958600ac3SJames Wright   // Place PETSc vector in libCEED vector
1530c63b910fSJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL));
153158600ac3SJames Wright   PetscCall(DMGetLocalVector(ctx->dm_x, &D_loc));
1532a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(D_loc, &mem_type, ctx->x_loc));
153358600ac3SJames Wright 
153458600ac3SJames Wright   // Compute Diagonal
1535c63b910fSJames Wright   PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL));
153650f50432SJames Wright   PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleDiagonal(ctx->op_mult, ctx->x_loc, CEED_REQUEST_IMMEDIATE));
1537c63b910fSJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL));
153858600ac3SJames Wright 
153958600ac3SJames Wright   // Restore PETSc vector
1540a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(ctx->x_loc, mem_type, D_loc));
154158600ac3SJames Wright 
154258600ac3SJames Wright   // Local-to-Global
154358600ac3SJames Wright   PetscCall(VecZeroEntries(D));
154458600ac3SJames Wright   PetscCall(DMLocalToGlobal(ctx->dm_x, D_loc, ADD_VALUES, D));
154558600ac3SJames Wright   PetscCall(DMRestoreLocalVector(ctx->dm_x, &D_loc));
1546c63b910fSJames Wright   PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL));
154758600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
154858600ac3SJames Wright }
154958600ac3SJames Wright 
155058600ac3SJames Wright /**
155158600ac3SJames Wright   @brief Compute `A X = Y` for a `MATCEED`.
155258600ac3SJames Wright 
155358600ac3SJames Wright   Collective across MPI processes.
155458600ac3SJames Wright 
155558600ac3SJames Wright   @param[in]   A  `MATCEED`
155658600ac3SJames Wright   @param[in]   X  Input PETSc vector
155758600ac3SJames Wright   @param[out]  Y  Output PETSc vector
155858600ac3SJames Wright 
155958600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
156058600ac3SJames Wright **/
156158600ac3SJames Wright PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) {
156258600ac3SJames Wright   MatCeedContext ctx;
156358600ac3SJames Wright 
156458600ac3SJames Wright   PetscFunctionBeginUser;
156558600ac3SJames Wright   PetscCall(MatShellGetContext(A, &ctx));
1566c63b910fSJames Wright   PetscCall(PetscLogEventBegin(ctx->log_event_mult, A, X, Y, NULL));
156758600ac3SJames Wright 
156858600ac3SJames Wright   {
156958600ac3SJames Wright     PetscMemType x_mem_type, y_mem_type;
157058600ac3SJames Wright     Vec          X_loc = ctx->X_loc, Y_loc;
157158600ac3SJames Wright 
157258600ac3SJames Wright     // Get local vectors
157358600ac3SJames Wright     if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
157458600ac3SJames Wright     PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
157558600ac3SJames Wright 
157658600ac3SJames Wright     // Global-to-local
157758600ac3SJames Wright     PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc));
157858600ac3SJames Wright 
157958600ac3SJames Wright     // Setup libCEED vectors
1580a7dac1d5SJames Wright     PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, ctx->x_loc));
158158600ac3SJames Wright     PetscCall(VecZeroEntries(Y_loc));
1582a7dac1d5SJames Wright     PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc));
158358600ac3SJames Wright 
158458600ac3SJames Wright     // Apply libCEED operator
1585c63b910fSJames Wright     PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult, A, X, Y, NULL));
1586537ec908SJames Wright     PetscCall(PetscLogGpuTimeBegin());
158750f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult, ctx->x_loc, ctx->y_loc, CEED_REQUEST_IMMEDIATE));
158858600ac3SJames Wright     PetscCall(PetscLogGpuTimeEnd());
1589537ec908SJames Wright     PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult, A, X, Y, NULL));
159058600ac3SJames Wright 
159158600ac3SJames Wright     // Restore PETSc vectors
1592a7dac1d5SJames Wright     PetscCall(VecReadCeedToPetsc(ctx->x_loc, x_mem_type, X_loc));
1593a7dac1d5SJames Wright     PetscCall(VecCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc));
159458600ac3SJames Wright 
159558600ac3SJames Wright     // Local-to-global
159658600ac3SJames Wright     PetscCall(VecZeroEntries(Y));
159758600ac3SJames Wright     PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y));
159858600ac3SJames Wright 
159958600ac3SJames Wright     // Restore local vectors, as needed
160058600ac3SJames Wright     if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
160158600ac3SJames Wright     PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
160258600ac3SJames Wright   }
160358600ac3SJames Wright 
160458600ac3SJames Wright   // Log flops
160558600ac3SJames Wright   if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult));
160658600ac3SJames Wright   else PetscCall(PetscLogFlops(ctx->flops_mult));
1607c63b910fSJames Wright   PetscCall(PetscLogEventEnd(ctx->log_event_mult, A, X, Y, NULL));
160858600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
160958600ac3SJames Wright }
161058600ac3SJames Wright 
161158600ac3SJames Wright /**
161258600ac3SJames Wright   @brief Compute `A^T Y = X` for a `MATCEED`.
161358600ac3SJames Wright 
161458600ac3SJames Wright   Collective across MPI processes.
161558600ac3SJames Wright 
161658600ac3SJames Wright   @param[in]   A  `MATCEED`
161758600ac3SJames Wright   @param[in]   Y  Input PETSc vector
161858600ac3SJames Wright   @param[out]  X  Output PETSc vector
161958600ac3SJames Wright 
162058600ac3SJames Wright   @return An error code: 0 - success, otherwise - failure
162158600ac3SJames Wright **/
162258600ac3SJames Wright PetscErrorCode MatMultTranspose_Ceed(Mat A, Vec Y, Vec X) {
162358600ac3SJames Wright   MatCeedContext ctx;
162458600ac3SJames Wright 
162558600ac3SJames Wright   PetscFunctionBeginUser;
162658600ac3SJames Wright   PetscCall(MatShellGetContext(A, &ctx));
1627c63b910fSJames Wright   PetscCall(PetscLogEventBegin(ctx->log_event_mult_transpose, A, Y, X, NULL));
162858600ac3SJames Wright 
162958600ac3SJames Wright   {
163058600ac3SJames Wright     PetscMemType x_mem_type, y_mem_type;
163158600ac3SJames Wright     Vec          X_loc, Y_loc = ctx->Y_loc_transpose;
163258600ac3SJames Wright 
163358600ac3SJames Wright     // Get local vectors
163458600ac3SJames Wright     if (!ctx->Y_loc_transpose) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc));
163558600ac3SJames Wright     PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc));
163658600ac3SJames Wright 
163758600ac3SJames Wright     // Global-to-local
163858600ac3SJames Wright     PetscCall(DMGlobalToLocal(ctx->dm_y, Y, INSERT_VALUES, Y_loc));
163958600ac3SJames Wright 
164058600ac3SJames Wright     // Setup libCEED vectors
1641a7dac1d5SJames Wright     PetscCall(VecReadPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc));
164258600ac3SJames Wright     PetscCall(VecZeroEntries(X_loc));
1643a7dac1d5SJames Wright     PetscCall(VecPetscToCeed(X_loc, &x_mem_type, ctx->x_loc));
164458600ac3SJames Wright 
164558600ac3SJames Wright     // Apply libCEED operator
1646c63b910fSJames Wright     PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL));
1647537ec908SJames Wright     PetscCall(PetscLogGpuTimeBegin());
164850f50432SJames Wright     PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult_transpose, ctx->y_loc, ctx->x_loc, CEED_REQUEST_IMMEDIATE));
164958600ac3SJames Wright     PetscCall(PetscLogGpuTimeEnd());
1650537ec908SJames Wright     PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL));
165158600ac3SJames Wright 
165258600ac3SJames Wright     // Restore PETSc vectors
1653a7dac1d5SJames Wright     PetscCall(VecReadCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc));
1654a7dac1d5SJames Wright     PetscCall(VecCeedToPetsc(ctx->x_loc, x_mem_type, X_loc));
165558600ac3SJames Wright 
165658600ac3SJames Wright     // Local-to-global
165758600ac3SJames Wright     PetscCall(VecZeroEntries(X));
165858600ac3SJames Wright     PetscCall(DMLocalToGlobal(ctx->dm_x, X_loc, ADD_VALUES, X));
165958600ac3SJames Wright 
166058600ac3SJames Wright     // Restore local vectors, as needed
166158600ac3SJames Wright     if (!ctx->Y_loc_transpose) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc));
166258600ac3SJames Wright     PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc));
166358600ac3SJames Wright   }
166458600ac3SJames Wright 
166558600ac3SJames Wright   // Log flops
166658600ac3SJames Wright   if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult_transpose));
166758600ac3SJames Wright   else PetscCall(PetscLogFlops(ctx->flops_mult_transpose));
1668c63b910fSJames Wright   PetscCall(PetscLogEventEnd(ctx->log_event_mult_transpose, A, Y, X, NULL));
166958600ac3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
167058600ac3SJames Wright }
1671