158600ac3SJames Wright /// @file 240d80af1SJames Wright /// MatCEED implementation 358600ac3SJames Wright 458600ac3SJames Wright #include <ceed.h> 558600ac3SJames Wright #include <ceed/backend.h> 658600ac3SJames Wright #include <mat-ceed-impl.h> 758600ac3SJames Wright #include <mat-ceed.h> 840d80af1SJames Wright #include <petsc-ceed-utils.h> 940d80af1SJames Wright #include <petsc-ceed.h> 10*000d2032SJeremy L Thompson #include <petscdm.h> 11*000d2032SJeremy L Thompson #include <petscmat.h> 12*000d2032SJeremy L Thompson #include <stdbool.h> 13*000d2032SJeremy L Thompson #include <stdio.h> 1458600ac3SJames Wright #include <stdlib.h> 1558600ac3SJames Wright #include <string.h> 16*000d2032SJeremy L Thompson #include <petsc/private/petscimpl.h> 1758600ac3SJames Wright 1858600ac3SJames Wright PetscClassId MATCEED_CLASSID; 19c63b910fSJames Wright PetscLogEvent MATCEED_MULT, MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE, MATCEED_MULT_TRANSPOSE_CEEDOP, MATCEED_ASSEMBLE_DIAGONAL, 20c63b910fSJames Wright MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, MATCEED_SETUP_PBDIAGONAL, MATCEED_SETUP_PBDIAGONAL_CEEDOP, MATCEED_ASSEMBLE_PBDIAGONAL, 21c63b910fSJames Wright MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, MATCEED_SETUP_FULL, MATCEED_SETUP_FULL_CEEDOP, MATCEED_ASSEMBLE_FULL, MATCEED_ASSEMBLE_FULL_CEEDOP; 2258600ac3SJames Wright 2358600ac3SJames Wright /** 2458600ac3SJames Wright @brief Register MATCEED log events. 2558600ac3SJames Wright 2658600ac3SJames Wright Not collective across MPI processes. 2758600ac3SJames Wright 2858600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 2958600ac3SJames Wright **/ 3058600ac3SJames Wright static PetscErrorCode MatCeedRegisterLogEvents() { 3140d80af1SJames Wright static PetscBool registered = PETSC_FALSE; 3258600ac3SJames Wright 3358600ac3SJames Wright PetscFunctionBeginUser; 3458600ac3SJames Wright if (registered) PetscFunctionReturn(PETSC_SUCCESS); 35c63b910fSJames Wright PetscCall(PetscClassIdRegister("MatCEED", &MATCEED_CLASSID)); 36c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMul", MATCEED_CLASSID, &MATCEED_MULT)); 37c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMulCeed", MATCEED_CLASSID, &MATCEED_MULT_CEEDOP)); 38c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMulT", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE)); 39c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMulTCeed", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE_CEEDOP)); 40c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmDiag", MATCEED_CLASSID, &MATCEED_ASSEMBLE_DIAGONAL)); 41c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSU", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL)); 42c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL_CEEDOP)); 43c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBD", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL)); 44c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBDCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP)); 45c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmSU", MATCEED_CLASSID, &MATCEED_SETUP_FULL)); 46c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_FULL_CEEDOP)); 47c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsm", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL)); 48c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL_CEEDOP)); 4940d80af1SJames Wright registered = PETSC_TRUE; 5058600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5158600ac3SJames Wright } 5258600ac3SJames Wright 5358600ac3SJames Wright /** 5458600ac3SJames Wright @brief Assemble the point block diagonal of a `MATCEED` into a `MATAIJ` or similar. 5558600ac3SJames Wright The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`. 5658600ac3SJames Wright The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail. 5758600ac3SJames Wright 5858600ac3SJames Wright Collective across MPI processes. 5958600ac3SJames Wright 6058600ac3SJames Wright @param[in] mat_ceed `MATCEED` to assemble 6158600ac3SJames Wright @param[in,out] mat_coo `MATAIJ` or similar to assemble into 6258600ac3SJames Wright 6358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 6458600ac3SJames Wright **/ 6558600ac3SJames Wright static PetscErrorCode MatCeedAssemblePointBlockDiagonalCOO(Mat mat_ceed, Mat mat_coo) { 6658600ac3SJames Wright MatCeedContext ctx; 6758600ac3SJames Wright 6858600ac3SJames Wright PetscFunctionBeginUser; 6958600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 7058600ac3SJames Wright 7158600ac3SJames Wright // Check if COO pattern set 7258600ac3SJames Wright { 7358600ac3SJames Wright PetscInt index = -1; 7458600ac3SJames Wright 7558600ac3SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) { 7658600ac3SJames Wright if (ctx->mats_assembled_pbd[i] == mat_coo) index = i; 7758600ac3SJames Wright } 7858600ac3SJames Wright if (index == -1) { 7958600ac3SJames Wright PetscInt *rows_petsc = NULL, *cols_petsc = NULL; 8058600ac3SJames Wright CeedInt *rows_ceed, *cols_ceed; 8158600ac3SJames Wright PetscCount num_entries; 8258600ac3SJames Wright PetscLogStage stage_amg_setup; 8358600ac3SJames Wright 8458600ac3SJames Wright // -- Assemble sparsity pattern if mat hasn't been assembled before 85c63b910fSJames Wright PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup)); 8658600ac3SJames Wright if (stage_amg_setup == -1) { 87c63b910fSJames Wright PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup)); 8858600ac3SJames Wright } 8958600ac3SJames Wright PetscCall(PetscLogStagePush(stage_amg_setup)); 90c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 91c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 9250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed)); 93c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 94a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc)); 95a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc)); 9658600ac3SJames Wright PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc)); 9758600ac3SJames Wright free(rows_petsc); 9858600ac3SJames Wright free(cols_petsc); 9950f50432SJames Wright if (!ctx->coo_values_pbd) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_pbd)); 10058600ac3SJames Wright PetscCall(PetscRealloc(++ctx->num_mats_assembled_pbd * sizeof(Mat), &ctx->mats_assembled_pbd)); 10158600ac3SJames Wright ctx->mats_assembled_pbd[ctx->num_mats_assembled_pbd - 1] = mat_coo; 102c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 10358600ac3SJames Wright PetscCall(PetscLogStagePop()); 10458600ac3SJames Wright } 10558600ac3SJames Wright } 10658600ac3SJames Wright 10758600ac3SJames Wright // Assemble mat_ceed 108c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 10958600ac3SJames Wright PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY)); 11058600ac3SJames Wright { 11158600ac3SJames Wright const CeedScalar *values; 11258600ac3SJames Wright MatType mat_type; 11358600ac3SJames Wright CeedMemType mem_type = CEED_MEM_HOST; 11458600ac3SJames Wright PetscBool is_spd, is_spd_known; 11558600ac3SJames Wright 11658600ac3SJames Wright PetscCall(MatGetType(mat_coo, &mat_type)); 11758600ac3SJames Wright if (strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 11858600ac3SJames Wright else if (strstr(mat_type, "kokkos")) mem_type = CEED_MEM_DEVICE; 11958600ac3SJames Wright else mem_type = CEED_MEM_HOST; 12058600ac3SJames Wright 121c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 12250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(ctx->op_mult, ctx->coo_values_pbd, CEED_REQUEST_IMMEDIATE)); 123c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 12450f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_pbd, mem_type, &values)); 12558600ac3SJames Wright PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES)); 12658600ac3SJames Wright PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd)); 12758600ac3SJames Wright if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd)); 12850f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_pbd, &values)); 12958600ac3SJames Wright } 13058600ac3SJames Wright PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY)); 131c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 13258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 13358600ac3SJames Wright } 13458600ac3SJames Wright 13558600ac3SJames Wright /** 13658600ac3SJames Wright @brief Assemble inner `Mat` for diagonal `PC` operations 13758600ac3SJames Wright 13858600ac3SJames Wright Collective across MPI processes. 13958600ac3SJames Wright 14058600ac3SJames Wright @param[in] mat_ceed `MATCEED` to invert 14158600ac3SJames Wright @param[in] use_ceed_pbd Boolean flag to use libCEED PBD assembly 14258600ac3SJames Wright @param[out] mat_inner Inner `Mat` for diagonal operations 14358600ac3SJames Wright 14458600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 14558600ac3SJames Wright **/ 14658600ac3SJames Wright static PetscErrorCode MatCeedAssembleInnerBlockDiagonalMat(Mat mat_ceed, PetscBool use_ceed_pbd, Mat *mat_inner) { 14758600ac3SJames Wright MatCeedContext ctx; 14858600ac3SJames Wright 14958600ac3SJames Wright PetscFunctionBeginUser; 15058600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 15158600ac3SJames Wright if (use_ceed_pbd) { 15258600ac3SJames Wright // Check if COO pattern set 15340d80af1SJames Wright if (!ctx->mat_assembled_pbd_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_pbd_internal)); 15458600ac3SJames Wright 15558600ac3SJames Wright // Assemble mat_assembled_full_internal 15658600ac3SJames Wright PetscCall(MatCeedAssemblePointBlockDiagonalCOO(mat_ceed, ctx->mat_assembled_pbd_internal)); 15758600ac3SJames Wright if (mat_inner) *mat_inner = ctx->mat_assembled_pbd_internal; 15858600ac3SJames Wright } else { 15958600ac3SJames Wright // Check if COO pattern set 16040d80af1SJames Wright if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal)); 16158600ac3SJames Wright 16258600ac3SJames Wright // Assemble mat_assembled_full_internal 16358600ac3SJames Wright PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal)); 16458600ac3SJames Wright if (mat_inner) *mat_inner = ctx->mat_assembled_full_internal; 16558600ac3SJames Wright } 16658600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 16758600ac3SJames Wright } 16858600ac3SJames Wright 16958600ac3SJames Wright /** 17058600ac3SJames Wright @brief Get `MATCEED` diagonal block for Jacobi. 17158600ac3SJames Wright 17258600ac3SJames Wright Collective across MPI processes. 17358600ac3SJames Wright 17458600ac3SJames Wright @param[in] mat_ceed `MATCEED` to invert 17558600ac3SJames Wright @param[out] mat_block The diagonal block matrix 17658600ac3SJames Wright 17758600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 17858600ac3SJames Wright **/ 17958600ac3SJames Wright static PetscErrorCode MatGetDiagonalBlock_Ceed(Mat mat_ceed, Mat *mat_block) { 18058600ac3SJames Wright Mat mat_inner = NULL; 18158600ac3SJames Wright MatCeedContext ctx; 18258600ac3SJames Wright 18358600ac3SJames Wright PetscFunctionBeginUser; 18458600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 18558600ac3SJames Wright 18658600ac3SJames Wright // Assemble inner mat if needed 18758600ac3SJames Wright PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_pbd_valid, &mat_inner)); 18858600ac3SJames Wright 18958600ac3SJames Wright // Get block diagonal 19058600ac3SJames Wright PetscCall(MatGetDiagonalBlock(mat_inner, mat_block)); 19158600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 19258600ac3SJames Wright } 19358600ac3SJames Wright 19458600ac3SJames Wright /** 19558600ac3SJames Wright @brief Invert `MATCEED` diagonal block for Jacobi. 19658600ac3SJames Wright 19758600ac3SJames Wright Collective across MPI processes. 19858600ac3SJames Wright 19958600ac3SJames Wright @param[in] mat_ceed `MATCEED` to invert 20058600ac3SJames Wright @param[out] values The block inverses in column major order 20158600ac3SJames Wright 20258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 20358600ac3SJames Wright **/ 20458600ac3SJames Wright static PetscErrorCode MatInvertBlockDiagonal_Ceed(Mat mat_ceed, const PetscScalar **values) { 20558600ac3SJames Wright Mat mat_inner = NULL; 20658600ac3SJames Wright MatCeedContext ctx; 20758600ac3SJames Wright 20858600ac3SJames Wright PetscFunctionBeginUser; 20958600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 21058600ac3SJames Wright 21158600ac3SJames Wright // Assemble inner mat if needed 21258600ac3SJames Wright PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_pbd_valid, &mat_inner)); 21358600ac3SJames Wright 21458600ac3SJames Wright // Invert PB diagonal 21558600ac3SJames Wright PetscCall(MatInvertBlockDiagonal(mat_inner, values)); 21658600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 21758600ac3SJames Wright } 21858600ac3SJames Wright 21958600ac3SJames Wright /** 22058600ac3SJames Wright @brief Invert `MATCEED` variable diagonal block for Jacobi. 22158600ac3SJames Wright 22258600ac3SJames Wright Collective across MPI processes. 22358600ac3SJames Wright 22458600ac3SJames Wright @param[in] mat_ceed `MATCEED` to invert 22558600ac3SJames Wright @param[in] num_blocks The number of blocks on the process 22658600ac3SJames Wright @param[in] block_sizes The size of each block on the process 22758600ac3SJames Wright @param[out] values The block inverses in column major order 22858600ac3SJames Wright 22958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 23058600ac3SJames Wright **/ 23158600ac3SJames Wright static PetscErrorCode MatInvertVariableBlockDiagonal_Ceed(Mat mat_ceed, PetscInt num_blocks, const PetscInt *block_sizes, PetscScalar *values) { 23258600ac3SJames Wright Mat mat_inner = NULL; 23358600ac3SJames Wright MatCeedContext ctx; 23458600ac3SJames Wright 23558600ac3SJames Wright PetscFunctionBeginUser; 23658600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 23758600ac3SJames Wright 23858600ac3SJames Wright // Assemble inner mat if needed 23958600ac3SJames Wright PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_vpbd_valid, &mat_inner)); 24058600ac3SJames Wright 24158600ac3SJames Wright // Invert PB diagonal 24258600ac3SJames Wright PetscCall(MatInvertVariableBlockDiagonal(mat_inner, num_blocks, block_sizes, values)); 24358600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 24458600ac3SJames Wright } 24558600ac3SJames Wright 246e90c2ceeSJames Wright /** 247e90c2ceeSJames Wright @brief View `MATCEED`. 248e90c2ceeSJames Wright 249e90c2ceeSJames Wright Collective across MPI processes. 250e90c2ceeSJames Wright 251e90c2ceeSJames Wright @param[in] mat_ceed `MATCEED` to view 252e90c2ceeSJames Wright @param[in] viewer The visualization context 253e90c2ceeSJames Wright 254e90c2ceeSJames Wright @return An error code: 0 - success, otherwise - failure 255e90c2ceeSJames Wright **/ 256e90c2ceeSJames Wright static PetscErrorCode MatView_Ceed(Mat mat_ceed, PetscViewer viewer) { 257e90c2ceeSJames Wright PetscBool is_ascii; 258e90c2ceeSJames Wright PetscViewerFormat format; 259*000d2032SJeremy L Thompson PetscMPIInt size, rank; 260e90c2ceeSJames Wright MatCeedContext ctx; 261e90c2ceeSJames Wright 262e90c2ceeSJames Wright PetscFunctionBeginUser; 263e90c2ceeSJames Wright PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 264e90c2ceeSJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 265e90c2ceeSJames Wright if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat_ceed), &viewer)); 266e90c2ceeSJames Wright 267e90c2ceeSJames Wright PetscCall(PetscViewerGetFormat(viewer, &format)); 268e90c2ceeSJames Wright PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat_ceed), &size)); 269e90c2ceeSJames Wright if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS); 270e90c2ceeSJames Wright 271*000d2032SJeremy L Thompson PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat_ceed), &rank)); 272*000d2032SJeremy L Thompson if (rank != 0) PetscFunctionReturn(PETSC_SUCCESS); 273*000d2032SJeremy L Thompson 274e90c2ceeSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii)); 275e90c2ceeSJames Wright { 276*000d2032SJeremy L Thompson PetscBool is_detailed = format == PETSC_VIEWER_ASCII_INFO_DETAIL; 277*000d2032SJeremy L Thompson char rank_string[16] = {'\0'}; 278e90c2ceeSJames Wright FILE *file; 279e90c2ceeSJames Wright 28040d80af1SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "MatCEED:\n Default COO MatType:%s\n", ctx->coo_mat_type)); 281e90c2ceeSJames Wright PetscCall(PetscViewerASCIIGetPointer(viewer, &file)); 282*000d2032SJeremy L Thompson PetscCall(PetscSNPrintf(rank_string, 16, "on Rank %d", rank)); 283*000d2032SJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(viewer, " CeedOperator Apply %s:\n", is_detailed ? rank_string : "Summary")); 284*000d2032SJeremy L Thompson if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult, file)); 285*000d2032SJeremy L Thompson else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult, file)); 286e90c2ceeSJames Wright if (ctx->op_mult_transpose) { 287*000d2032SJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(viewer, " CeedOperator ApplyTranspose %s:\n", is_detailed ? rank_string : "Summary")); 288*000d2032SJeremy L Thompson if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult_transpose, file)); 289*000d2032SJeremy L Thompson else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult_transpose, file)); 290e90c2ceeSJames Wright } 291e90c2ceeSJames Wright } 292e90c2ceeSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 293e90c2ceeSJames Wright } 294e90c2ceeSJames Wright 29558600ac3SJames Wright // ----------------------------------------------------------------------------- 29658600ac3SJames Wright // MatCeed 29758600ac3SJames Wright // ----------------------------------------------------------------------------- 29858600ac3SJames Wright 29958600ac3SJames Wright /** 30058600ac3SJames Wright @brief Create PETSc `Mat` from libCEED operators. 30158600ac3SJames Wright 30258600ac3SJames Wright Collective across MPI processes. 30358600ac3SJames Wright 30458600ac3SJames Wright @param[in] dm_x Input `DM` 30558600ac3SJames Wright @param[in] dm_y Output `DM` 30658600ac3SJames Wright @param[in] op_mult `CeedOperator` for forward evaluation 30758600ac3SJames Wright @param[in] op_mult_transpose `CeedOperator` for transpose evaluation 30858600ac3SJames Wright @param[out] mat New MatCeed 30958600ac3SJames Wright 31058600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 31158600ac3SJames Wright **/ 312*000d2032SJeremy L Thompson PetscErrorCode MatCreateCeed(DM dm_x, DM dm_y, CeedOperator op_mult, CeedOperator op_mult_transpose, Mat *mat) { 31358600ac3SJames Wright PetscInt X_l_size, X_g_size, Y_l_size, Y_g_size; 31458600ac3SJames Wright VecType vec_type; 31558600ac3SJames Wright MatCeedContext ctx; 31658600ac3SJames Wright 31758600ac3SJames Wright PetscFunctionBeginUser; 31858600ac3SJames Wright PetscCall(MatCeedRegisterLogEvents()); 31958600ac3SJames Wright 32058600ac3SJames Wright // Collect context data 32158600ac3SJames Wright PetscCall(DMGetVecType(dm_x, &vec_type)); 32258600ac3SJames Wright { 32358600ac3SJames Wright Vec X; 32458600ac3SJames Wright 32558600ac3SJames Wright PetscCall(DMGetGlobalVector(dm_x, &X)); 32658600ac3SJames Wright PetscCall(VecGetSize(X, &X_g_size)); 32758600ac3SJames Wright PetscCall(VecGetLocalSize(X, &X_l_size)); 32858600ac3SJames Wright PetscCall(DMRestoreGlobalVector(dm_x, &X)); 32958600ac3SJames Wright } 33058600ac3SJames Wright if (dm_y) { 33158600ac3SJames Wright Vec Y; 33258600ac3SJames Wright 33358600ac3SJames Wright PetscCall(DMGetGlobalVector(dm_y, &Y)); 33458600ac3SJames Wright PetscCall(VecGetSize(Y, &Y_g_size)); 33558600ac3SJames Wright PetscCall(VecGetLocalSize(Y, &Y_l_size)); 33658600ac3SJames Wright PetscCall(DMRestoreGlobalVector(dm_y, &Y)); 33758600ac3SJames Wright } else { 33858600ac3SJames Wright dm_y = dm_x; 33958600ac3SJames Wright Y_g_size = X_g_size; 34058600ac3SJames Wright Y_l_size = X_l_size; 34158600ac3SJames Wright } 34240d80af1SJames Wright 34358600ac3SJames Wright // Create context 34458600ac3SJames Wright { 34558600ac3SJames Wright Vec X_loc, Y_loc_transpose = NULL; 34658600ac3SJames Wright 34758600ac3SJames Wright PetscCall(DMCreateLocalVector(dm_x, &X_loc)); 34858600ac3SJames Wright PetscCall(VecZeroEntries(X_loc)); 34958600ac3SJames Wright if (op_mult_transpose) { 35058600ac3SJames Wright PetscCall(DMCreateLocalVector(dm_y, &Y_loc_transpose)); 35158600ac3SJames Wright PetscCall(VecZeroEntries(Y_loc_transpose)); 35258600ac3SJames Wright } 353c63b910fSJames Wright PetscCall(MatCeedContextCreate(dm_x, dm_y, X_loc, Y_loc_transpose, op_mult, op_mult_transpose, MATCEED_MULT, MATCEED_MULT_TRANSPOSE, 354c63b910fSJames Wright MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE_CEEDOP, &ctx)); 35558600ac3SJames Wright PetscCall(VecDestroy(&X_loc)); 35658600ac3SJames Wright PetscCall(VecDestroy(&Y_loc_transpose)); 35758600ac3SJames Wright } 35858600ac3SJames Wright 35958600ac3SJames Wright // Create mat 36058600ac3SJames Wright PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm_x), Y_l_size, X_l_size, Y_g_size, X_g_size, ctx, mat)); 36158600ac3SJames Wright PetscCall(PetscObjectChangeTypeName((PetscObject)*mat, MATCEED)); 36258600ac3SJames Wright // -- Set block and variable block sizes 36358600ac3SJames Wright if (dm_x == dm_y) { 36458600ac3SJames Wright MatType dm_mat_type, dm_mat_type_copy; 36558600ac3SJames Wright Mat temp_mat; 36658600ac3SJames Wright 36758600ac3SJames Wright PetscCall(DMGetMatType(dm_x, &dm_mat_type)); 36858600ac3SJames Wright PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy)); 36958600ac3SJames Wright PetscCall(DMSetMatType(dm_x, MATAIJ)); 37058600ac3SJames Wright PetscCall(DMCreateMatrix(dm_x, &temp_mat)); 37158600ac3SJames Wright PetscCall(DMSetMatType(dm_x, dm_mat_type_copy)); 37258600ac3SJames Wright PetscCall(PetscFree(dm_mat_type_copy)); 37358600ac3SJames Wright 37458600ac3SJames Wright { 37558600ac3SJames Wright PetscInt block_size, num_blocks, max_vblock_size = PETSC_INT_MAX; 37658600ac3SJames Wright const PetscInt *vblock_sizes; 37758600ac3SJames Wright 37858600ac3SJames Wright // -- Get block sizes 37958600ac3SJames Wright PetscCall(MatGetBlockSize(temp_mat, &block_size)); 38058600ac3SJames Wright PetscCall(MatGetVariableBlockSizes(temp_mat, &num_blocks, &vblock_sizes)); 38158600ac3SJames Wright { 38258600ac3SJames Wright PetscInt local_min_max[2] = {0}, global_min_max[2] = {0, PETSC_INT_MAX}; 38358600ac3SJames Wright 38458600ac3SJames Wright for (PetscInt i = 0; i < num_blocks; i++) local_min_max[1] = PetscMax(local_min_max[1], vblock_sizes[i]); 38558600ac3SJames Wright PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_min_max, global_min_max)); 38658600ac3SJames Wright max_vblock_size = global_min_max[1]; 38758600ac3SJames Wright } 38858600ac3SJames Wright 38958600ac3SJames Wright // -- Copy block sizes 39058600ac3SJames Wright if (block_size > 1) PetscCall(MatSetBlockSize(*mat, block_size)); 39158600ac3SJames Wright if (num_blocks) PetscCall(MatSetVariableBlockSizes(*mat, num_blocks, (PetscInt *)vblock_sizes)); 39258600ac3SJames Wright 39358600ac3SJames Wright // -- Check libCEED compatibility 39458600ac3SJames Wright { 39558600ac3SJames Wright bool is_composite; 39658600ac3SJames Wright 39758600ac3SJames Wright ctx->is_ceed_pbd_valid = PETSC_TRUE; 39858600ac3SJames Wright ctx->is_ceed_vpbd_valid = PETSC_TRUE; 39950f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorIsComposite(op_mult, &is_composite)); 40058600ac3SJames Wright if (is_composite) { 40158600ac3SJames Wright CeedInt num_sub_operators; 40258600ac3SJames Wright CeedOperator *sub_operators; 40358600ac3SJames Wright 40450f50432SJames Wright PetscCallCeed(ctx->ceed, CeedCompositeOperatorGetNumSub(op_mult, &num_sub_operators)); 40550f50432SJames Wright PetscCallCeed(ctx->ceed, CeedCompositeOperatorGetSubList(op_mult, &sub_operators)); 40658600ac3SJames Wright for (CeedInt i = 0; i < num_sub_operators; i++) { 40758600ac3SJames Wright CeedInt num_bases, num_comp; 40858600ac3SJames Wright CeedBasis *active_bases; 40958600ac3SJames Wright CeedOperatorAssemblyData assembly_data; 41058600ac3SJames Wright 41150f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(sub_operators[i], &assembly_data)); 41250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL)); 41350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp)); 41458600ac3SJames Wright if (num_bases > 1) { 41558600ac3SJames Wright ctx->is_ceed_pbd_valid = PETSC_FALSE; 41658600ac3SJames Wright ctx->is_ceed_vpbd_valid = PETSC_FALSE; 41758600ac3SJames Wright } 41858600ac3SJames Wright if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE; 41958600ac3SJames Wright if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE; 42058600ac3SJames Wright } 42158600ac3SJames Wright } else { 42258600ac3SJames Wright // LCOV_EXCL_START 42358600ac3SJames Wright CeedInt num_bases, num_comp; 42458600ac3SJames Wright CeedBasis *active_bases; 42558600ac3SJames Wright CeedOperatorAssemblyData assembly_data; 42658600ac3SJames Wright 42750f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(op_mult, &assembly_data)); 42850f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL)); 42950f50432SJames Wright PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp)); 43058600ac3SJames Wright if (num_bases > 1) { 43158600ac3SJames Wright ctx->is_ceed_pbd_valid = PETSC_FALSE; 43258600ac3SJames Wright ctx->is_ceed_vpbd_valid = PETSC_FALSE; 43358600ac3SJames Wright } 43458600ac3SJames Wright if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE; 43558600ac3SJames Wright if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE; 43658600ac3SJames Wright // LCOV_EXCL_STOP 43758600ac3SJames Wright } 43858600ac3SJames Wright { 43958600ac3SJames Wright PetscInt local_is_valid[2], global_is_valid[2]; 44058600ac3SJames Wright 44158600ac3SJames Wright local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_pbd_valid; 44258600ac3SJames Wright PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid)); 44358600ac3SJames Wright ctx->is_ceed_pbd_valid = global_is_valid[0]; 44458600ac3SJames Wright local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_vpbd_valid; 44558600ac3SJames Wright PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid)); 44658600ac3SJames Wright ctx->is_ceed_vpbd_valid = global_is_valid[0]; 44758600ac3SJames Wright } 44858600ac3SJames Wright } 44958600ac3SJames Wright } 45058600ac3SJames Wright PetscCall(MatDestroy(&temp_mat)); 45158600ac3SJames Wright } 45258600ac3SJames Wright // -- Set internal mat type 45358600ac3SJames Wright { 45458600ac3SJames Wright VecType vec_type; 45540d80af1SJames Wright MatType coo_mat_type; 45658600ac3SJames Wright 45758600ac3SJames Wright PetscCall(VecGetType(ctx->X_loc, &vec_type)); 45840d80af1SJames Wright if (strstr(vec_type, VECCUDA)) coo_mat_type = MATAIJCUSPARSE; 45940d80af1SJames Wright else if (strstr(vec_type, VECKOKKOS)) coo_mat_type = MATAIJKOKKOS; 46040d80af1SJames Wright else coo_mat_type = MATAIJ; 46140d80af1SJames Wright PetscCall(PetscStrallocpy(coo_mat_type, &ctx->coo_mat_type)); 46258600ac3SJames Wright } 46358600ac3SJames Wright // -- Set mat operations 46458600ac3SJames Wright PetscCall(MatShellSetContextDestroy(*mat, (PetscErrorCode(*)(void *))MatCeedContextDestroy)); 465e90c2ceeSJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_VIEW, (void (*)(void))MatView_Ceed)); 46658600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 46758600ac3SJames Wright if (op_mult_transpose) PetscCall(MatShellSetOperation(*mat, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed)); 46858600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed)); 46958600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed)); 47058600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_INVERT_BLOCK_DIAGONAL, (void (*)(void))MatInvertBlockDiagonal_Ceed)); 47158600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_INVERT_VBLOCK_DIAGONAL, (void (*)(void))MatInvertVariableBlockDiagonal_Ceed)); 47258600ac3SJames Wright PetscCall(MatShellSetVecType(*mat, vec_type)); 47358600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 47458600ac3SJames Wright } 47558600ac3SJames Wright 47658600ac3SJames Wright /** 47758600ac3SJames Wright @brief Copy `MATCEED` into a compatible `Mat` with type `MatShell` or `MATCEED`. 47858600ac3SJames Wright 47958600ac3SJames Wright Collective across MPI processes. 48058600ac3SJames Wright 48158600ac3SJames Wright @param[in] mat_ceed `MATCEED` to copy from 48258600ac3SJames Wright @param[out] mat_other `MatShell` or `MATCEED` to copy into 48358600ac3SJames Wright 48458600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 48558600ac3SJames Wright **/ 48658600ac3SJames Wright PetscErrorCode MatCeedCopy(Mat mat_ceed, Mat mat_other) { 48758600ac3SJames Wright PetscFunctionBeginUser; 48858600ac3SJames Wright PetscCall(MatCeedRegisterLogEvents()); 48958600ac3SJames Wright 49058600ac3SJames Wright // Check type compatibility 49158600ac3SJames Wright { 49240d80af1SJames Wright PetscBool is_matceed = PETSC_FALSE, is_matshell = PETSC_FALSE; 49358600ac3SJames Wright MatType mat_type_ceed, mat_type_other; 49458600ac3SJames Wright 49558600ac3SJames Wright PetscCall(MatGetType(mat_ceed, &mat_type_ceed)); 49640d80af1SJames Wright PetscCall(PetscStrcmp(mat_type_ceed, MATCEED, &is_matceed)); 49740d80af1SJames Wright PetscCheck(is_matceed, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_ceed must have type " MATCEED); 49840d80af1SJames Wright PetscCall(MatGetType(mat_other, &mat_type_other)); 49940d80af1SJames Wright PetscCall(PetscStrcmp(mat_type_other, MATCEED, &is_matceed)); 50040d80af1SJames Wright PetscCall(PetscStrcmp(mat_type_other, MATSHELL, &is_matceed)); 50140d80af1SJames Wright PetscCheck(is_matceed || is_matshell, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_other must have type " MATCEED " or " MATSHELL); 50258600ac3SJames Wright } 50358600ac3SJames Wright 50458600ac3SJames Wright // Check dimension compatibility 50558600ac3SJames Wright { 50658600ac3SJames 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; 50758600ac3SJames Wright 50858600ac3SJames Wright PetscCall(MatGetSize(mat_ceed, &Y_g_ceed_size, &X_g_ceed_size)); 50958600ac3SJames Wright PetscCall(MatGetLocalSize(mat_ceed, &Y_l_ceed_size, &X_l_ceed_size)); 51058600ac3SJames Wright PetscCall(MatGetSize(mat_ceed, &Y_g_other_size, &X_g_other_size)); 51158600ac3SJames Wright PetscCall(MatGetLocalSize(mat_ceed, &Y_l_other_size, &X_l_other_size)); 51258600ac3SJames 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) && 51358600ac3SJames Wright (X_l_ceed_size == X_l_other_size), 51458600ac3SJames Wright PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, 51558600ac3SJames Wright "mat_ceed and mat_other must have compatible sizes; found mat_ceed (Global: %" PetscInt_FMT ", %" PetscInt_FMT 51658600ac3SJames Wright "; Local: %" PetscInt_FMT ", %" PetscInt_FMT ") mat_other (Global: %" PetscInt_FMT ", %" PetscInt_FMT "; Local: %" PetscInt_FMT 51758600ac3SJames Wright ", %" PetscInt_FMT ")", 51858600ac3SJames 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); 51958600ac3SJames Wright } 52058600ac3SJames Wright 52158600ac3SJames Wright // Convert 52258600ac3SJames Wright { 52358600ac3SJames Wright VecType vec_type; 52458600ac3SJames Wright MatCeedContext ctx; 52558600ac3SJames Wright 52658600ac3SJames Wright PetscCall(PetscObjectChangeTypeName((PetscObject)mat_other, MATCEED)); 52758600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 52858600ac3SJames Wright PetscCall(MatCeedContextReference(ctx)); 52958600ac3SJames Wright PetscCall(MatShellSetContext(mat_other, ctx)); 53058600ac3SJames Wright PetscCall(MatShellSetContextDestroy(mat_other, (PetscErrorCode(*)(void *))MatCeedContextDestroy)); 531e90c2ceeSJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_VIEW, (void (*)(void))MatView_Ceed)); 53258600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 53358600ac3SJames Wright if (ctx->op_mult_transpose) PetscCall(MatShellSetOperation(mat_other, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed)); 53458600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed)); 53558600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed)); 53658600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_INVERT_BLOCK_DIAGONAL, (void (*)(void))MatInvertBlockDiagonal_Ceed)); 53758600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_INVERT_VBLOCK_DIAGONAL, (void (*)(void))MatInvertVariableBlockDiagonal_Ceed)); 53858600ac3SJames Wright { 53958600ac3SJames Wright PetscInt block_size; 54058600ac3SJames Wright 54158600ac3SJames Wright PetscCall(MatGetBlockSize(mat_ceed, &block_size)); 54258600ac3SJames Wright if (block_size > 1) PetscCall(MatSetBlockSize(mat_other, block_size)); 54358600ac3SJames Wright } 54458600ac3SJames Wright { 54558600ac3SJames Wright PetscInt num_blocks; 54658600ac3SJames Wright const PetscInt *block_sizes; 54758600ac3SJames Wright 54858600ac3SJames Wright PetscCall(MatGetVariableBlockSizes(mat_ceed, &num_blocks, &block_sizes)); 54958600ac3SJames Wright if (num_blocks) PetscCall(MatSetVariableBlockSizes(mat_other, num_blocks, (PetscInt *)block_sizes)); 55058600ac3SJames Wright } 55158600ac3SJames Wright PetscCall(DMGetVecType(ctx->dm_x, &vec_type)); 55258600ac3SJames Wright PetscCall(MatShellSetVecType(mat_other, vec_type)); 55358600ac3SJames Wright } 55458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 55558600ac3SJames Wright } 55658600ac3SJames Wright 55758600ac3SJames Wright /** 558*000d2032SJeremy L Thompson @brief Mark Mark `CeedQFunction` data as updated and the `CeedQFunction` as requiring re-assembly for a `MatCEED`. 559*000d2032SJeremy L Thompson 560*000d2032SJeremy L Thompson Collective across MPI processes. 561*000d2032SJeremy L Thompson 562*000d2032SJeremy L Thompson @param[in] mat_ceed `MATCEED` 563*000d2032SJeremy L Thompson @param[out] update_needed Boolean flag indicating `CeedQFunction` update needed 564*000d2032SJeremy L Thompson 565*000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 566*000d2032SJeremy L Thompson **/ 567*000d2032SJeremy L Thompson PetscErrorCode MatCeedSetAssemblyDataUpdateNeeded(Mat mat_ceed, PetscBool update_needed) { 568*000d2032SJeremy L Thompson MatCeedContext ctx; 569*000d2032SJeremy L Thompson 570*000d2032SJeremy L Thompson PetscFunctionBeginUser; 571*000d2032SJeremy L Thompson PetscCall(MatShellGetContext(mat_ceed, &ctx)); 572*000d2032SJeremy L Thompson PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult, update_needed)); 573*000d2032SJeremy L Thompson if (ctx->op_mult_transpose) { 574*000d2032SJeremy L Thompson PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult_transpose, update_needed)); 575*000d2032SJeremy L Thompson } 576*000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 577*000d2032SJeremy L Thompson } 578*000d2032SJeremy L Thompson 579*000d2032SJeremy L Thompson /** 58040d80af1SJames Wright @brief Setup a `Mat` with the same COO pattern as a `MatCEED`. 58140d80af1SJames Wright 58240d80af1SJames Wright Collective across MPI processes. 58340d80af1SJames Wright 58440d80af1SJames Wright @param[in] mat_ceed `MATCEED` 58540d80af1SJames Wright @param[out] mat_coo Sparse `Mat` with same COO pattern 58640d80af1SJames Wright 58740d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 58840d80af1SJames Wright **/ 58940d80af1SJames Wright PetscErrorCode MatCeedCreateMatCOO(Mat mat_ceed, Mat *mat_coo) { 59040d80af1SJames Wright MatCeedContext ctx; 59140d80af1SJames Wright 59240d80af1SJames Wright PetscFunctionBeginUser; 59340d80af1SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 59440d80af1SJames Wright 59540d80af1SJames 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"); 59640d80af1SJames Wright 59740d80af1SJames Wright // Check cl mat type 59840d80af1SJames Wright { 59940d80af1SJames Wright PetscBool is_coo_mat_type_cl = PETSC_FALSE; 60040d80af1SJames Wright char coo_mat_type_cl[64]; 60140d80af1SJames Wright 60240d80af1SJames Wright // Check for specific CL coo mat type for this Mat 60340d80af1SJames Wright { 60440d80af1SJames Wright const char *mat_ceed_prefix = NULL; 60540d80af1SJames Wright 60640d80af1SJames Wright PetscCall(MatGetOptionsPrefix(mat_ceed, &mat_ceed_prefix)); 60740d80af1SJames Wright PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), mat_ceed_prefix, "", NULL); 60840d80af1SJames Wright PetscCall(PetscOptionsFList("-ceed_coo_mat_type", "Default MATCEED COO assembly MatType", NULL, MatList, coo_mat_type_cl, coo_mat_type_cl, 60940d80af1SJames Wright sizeof(coo_mat_type_cl), &is_coo_mat_type_cl)); 61040d80af1SJames Wright PetscOptionsEnd(); 61140d80af1SJames Wright if (is_coo_mat_type_cl) { 61240d80af1SJames Wright PetscCall(PetscFree(ctx->coo_mat_type)); 61340d80af1SJames Wright PetscCall(PetscStrallocpy(coo_mat_type_cl, &ctx->coo_mat_type)); 61440d80af1SJames Wright } 61540d80af1SJames Wright } 61640d80af1SJames Wright } 61740d80af1SJames Wright 61840d80af1SJames Wright // Create sparse matrix 61940d80af1SJames Wright { 62040d80af1SJames Wright MatType dm_mat_type, dm_mat_type_copy; 62140d80af1SJames Wright 62240d80af1SJames Wright PetscCall(DMGetMatType(ctx->dm_x, &dm_mat_type)); 62340d80af1SJames Wright PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy)); 62440d80af1SJames Wright PetscCall(DMSetMatType(ctx->dm_x, ctx->coo_mat_type)); 62540d80af1SJames Wright PetscCall(DMCreateMatrix(ctx->dm_x, mat_coo)); 62640d80af1SJames Wright PetscCall(DMSetMatType(ctx->dm_x, dm_mat_type_copy)); 62740d80af1SJames Wright PetscCall(PetscFree(dm_mat_type_copy)); 62840d80af1SJames Wright } 62940d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63040d80af1SJames Wright } 63140d80af1SJames Wright 63240d80af1SJames Wright /** 63340d80af1SJames Wright @brief Setup the COO preallocation `MATCEED` into a `MATAIJ` or similar. 63458600ac3SJames Wright The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail. 63558600ac3SJames Wright 63658600ac3SJames Wright Collective across MPI processes. 63758600ac3SJames Wright 63858600ac3SJames Wright @param[in] mat_ceed `MATCEED` to assemble 63958600ac3SJames Wright @param[in,out] mat_coo `MATAIJ` or similar to assemble into 64058600ac3SJames Wright 64158600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 64258600ac3SJames Wright **/ 64340d80af1SJames Wright PetscErrorCode MatCeedSetPreallocationCOO(Mat mat_ceed, Mat mat_coo) { 64458600ac3SJames Wright MatCeedContext ctx; 64558600ac3SJames Wright 64658600ac3SJames Wright PetscFunctionBeginUser; 64758600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 64858600ac3SJames Wright 64958600ac3SJames Wright { 65058600ac3SJames Wright PetscInt *rows_petsc = NULL, *cols_petsc = NULL; 65158600ac3SJames Wright CeedInt *rows_ceed, *cols_ceed; 65258600ac3SJames Wright PetscCount num_entries; 65358600ac3SJames Wright PetscLogStage stage_amg_setup; 65458600ac3SJames Wright 65558600ac3SJames Wright // -- Assemble sparsity pattern if mat hasn't been assembled before 656c63b910fSJames Wright PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup)); 65758600ac3SJames Wright if (stage_amg_setup == -1) { 658c63b910fSJames Wright PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup)); 65958600ac3SJames Wright } 66058600ac3SJames Wright PetscCall(PetscLogStagePush(stage_amg_setup)); 661c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL)); 662c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 66350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed)); 664c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 665a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc)); 666a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc)); 66758600ac3SJames Wright PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc)); 66858600ac3SJames Wright free(rows_petsc); 66958600ac3SJames Wright free(cols_petsc); 67050f50432SJames Wright if (!ctx->coo_values_full) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_full)); 67158600ac3SJames Wright PetscCall(PetscRealloc(++ctx->num_mats_assembled_full * sizeof(Mat), &ctx->mats_assembled_full)); 67258600ac3SJames Wright ctx->mats_assembled_full[ctx->num_mats_assembled_full - 1] = mat_coo; 673c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL)); 67458600ac3SJames Wright PetscCall(PetscLogStagePop()); 67558600ac3SJames Wright } 67640d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 67740d80af1SJames Wright } 67840d80af1SJames Wright 67940d80af1SJames Wright /** 68040d80af1SJames Wright @brief Assemble a `MATCEED` into a `MATAIJ` or similar. 68140d80af1SJames Wright The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`. 68240d80af1SJames Wright The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail. 68340d80af1SJames Wright 68440d80af1SJames Wright Collective across MPI processes. 68540d80af1SJames Wright 68640d80af1SJames Wright @param[in] mat_ceed `MATCEED` to assemble 68740d80af1SJames Wright @param[in,out] mat_coo `MATAIJ` or similar to assemble into 68840d80af1SJames Wright 68940d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 69040d80af1SJames Wright **/ 69140d80af1SJames Wright PetscErrorCode MatCeedAssembleCOO(Mat mat_ceed, Mat mat_coo) { 69240d80af1SJames Wright MatCeedContext ctx; 69340d80af1SJames Wright 69440d80af1SJames Wright PetscFunctionBeginUser; 69540d80af1SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 69640d80af1SJames Wright 69740d80af1SJames Wright // Set COO pattern if needed 69840d80af1SJames Wright { 69940d80af1SJames Wright CeedInt index = -1; 70040d80af1SJames Wright 70140d80af1SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) { 70240d80af1SJames Wright if (ctx->mats_assembled_full[i] == mat_coo) index = i; 70340d80af1SJames Wright } 70440d80af1SJames Wright if (index == -1) PetscCall(MatCeedSetPreallocationCOO(mat_ceed, mat_coo)); 70558600ac3SJames Wright } 70658600ac3SJames Wright 70758600ac3SJames Wright // Assemble mat_ceed 708c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL)); 70958600ac3SJames Wright PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY)); 71058600ac3SJames Wright { 71158600ac3SJames Wright const CeedScalar *values; 71258600ac3SJames Wright MatType mat_type; 71358600ac3SJames Wright CeedMemType mem_type = CEED_MEM_HOST; 71458600ac3SJames Wright PetscBool is_spd, is_spd_known; 71558600ac3SJames Wright 71658600ac3SJames Wright PetscCall(MatGetType(mat_coo, &mat_type)); 71758600ac3SJames Wright if (strstr(mat_type, "cusparse")) mem_type = CEED_MEM_DEVICE; 71858600ac3SJames Wright else if (strstr(mat_type, "kokkos")) mem_type = CEED_MEM_DEVICE; 71958600ac3SJames Wright else mem_type = CEED_MEM_HOST; 72058600ac3SJames Wright 721c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 72250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemble(ctx->op_mult, ctx->coo_values_full)); 723c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 72450f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_full, mem_type, &values)); 72558600ac3SJames Wright PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES)); 72658600ac3SJames Wright PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd)); 72758600ac3SJames Wright if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd)); 72850f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_full, &values)); 72958600ac3SJames Wright } 73058600ac3SJames Wright PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY)); 731c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL)); 73258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 73358600ac3SJames Wright } 73458600ac3SJames Wright 73558600ac3SJames Wright /** 73640d80af1SJames Wright @brief Set the current value of a context field for a `MatCEED`. 73740d80af1SJames Wright 73840d80af1SJames Wright Not collective across MPI processes. 73940d80af1SJames Wright 74040d80af1SJames Wright @param[in,out] mat `MatCEED` 74140d80af1SJames Wright @param[in] name Name of the context field 74240d80af1SJames Wright @param[in] value New context field value 74340d80af1SJames Wright 74440d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 74540d80af1SJames Wright **/ 74640d80af1SJames Wright PetscErrorCode MatCeedSetContextDouble(Mat mat, const char *name, double value) { 74740d80af1SJames Wright PetscBool was_updated = PETSC_FALSE; 74840d80af1SJames Wright MatCeedContext ctx; 74940d80af1SJames Wright 75040d80af1SJames Wright PetscFunctionBeginUser; 75140d80af1SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 75240d80af1SJames Wright { 75340d80af1SJames Wright CeedContextFieldLabel label = NULL; 75440d80af1SJames Wright 75540d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult, name, &label)); 75640d80af1SJames Wright if (label) { 757*000d2032SJeremy L Thompson double set_value = 2 * value + 1.0; 758*000d2032SJeremy L Thompson 759*000d2032SJeremy L Thompson PetscCall(MatCeedGetContextDouble(mat, name, &set_value)); 760*000d2032SJeremy L Thompson if (set_value != value) { 76140d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult, label, &value)); 76240d80af1SJames Wright was_updated = PETSC_TRUE; 76340d80af1SJames Wright } 764*000d2032SJeremy L Thompson } 76540d80af1SJames Wright if (ctx->op_mult_transpose) { 76640d80af1SJames Wright label = NULL; 76740d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult_transpose, name, &label)); 76840d80af1SJames Wright if (label) { 769*000d2032SJeremy L Thompson double set_value = 2 * value + 1.0; 770*000d2032SJeremy L Thompson 771*000d2032SJeremy L Thompson PetscCall(MatCeedGetContextDouble(mat, name, &set_value)); 772*000d2032SJeremy L Thompson if (set_value != value) { 77340d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult_transpose, label, &value)); 77440d80af1SJames Wright was_updated = PETSC_TRUE; 77540d80af1SJames Wright } 77640d80af1SJames Wright } 77740d80af1SJames Wright } 778*000d2032SJeremy L Thompson } 77940d80af1SJames Wright if (was_updated) { 78040d80af1SJames Wright PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 78140d80af1SJames Wright PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 78240d80af1SJames Wright } 78340d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 78440d80af1SJames Wright } 78540d80af1SJames Wright 78640d80af1SJames Wright /** 78740d80af1SJames Wright @brief Get the current value of a context field for a `MatCEED`. 78840d80af1SJames Wright 78940d80af1SJames Wright Not collective across MPI processes. 79040d80af1SJames Wright 79140d80af1SJames Wright @param[in] mat `MatCEED` 79240d80af1SJames Wright @param[in] name Name of the context field 79340d80af1SJames Wright @param[out] value Current context field value 79440d80af1SJames Wright 79540d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 79640d80af1SJames Wright **/ 79740d80af1SJames Wright PetscErrorCode MatCeedGetContextDouble(Mat mat, const char *name, double *value) { 79840d80af1SJames Wright MatCeedContext ctx; 79940d80af1SJames Wright 80040d80af1SJames Wright PetscFunctionBeginUser; 80140d80af1SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 80240d80af1SJames Wright { 80340d80af1SJames Wright CeedContextFieldLabel label = NULL; 80440d80af1SJames Wright CeedOperator op = ctx->op_mult; 80540d80af1SJames Wright 80640d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label)); 80740d80af1SJames Wright if (!label && ctx->op_mult_transpose) { 80840d80af1SJames Wright op = ctx->op_mult_transpose; 80940d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label)); 81040d80af1SJames Wright } 81140d80af1SJames Wright if (label) { 81240d80af1SJames Wright PetscSizeT num_values; 81340d80af1SJames Wright const double *values_ceed; 81440d80af1SJames Wright 81540d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextDoubleRead(op, label, &num_values, &values_ceed)); 81640d80af1SJames Wright *value = values_ceed[0]; 81740d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorRestoreContextDoubleRead(op, label, &values_ceed)); 81840d80af1SJames Wright } 81940d80af1SJames Wright } 82040d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 82140d80af1SJames Wright } 82240d80af1SJames Wright 82340d80af1SJames Wright /** 824*000d2032SJeremy L Thompson @brief Set the current `PetscReal` value of a context field for a `MatCEED`. 825*000d2032SJeremy L Thompson 826*000d2032SJeremy L Thompson Not collective across MPI processes. 827*000d2032SJeremy L Thompson 828*000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 829*000d2032SJeremy L Thompson @param[in] name Name of the context field 830*000d2032SJeremy L Thompson @param[in] value New context field value 831*000d2032SJeremy L Thompson 832*000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 833*000d2032SJeremy L Thompson **/ 834*000d2032SJeremy L Thompson PetscErrorCode MatCeedSetContextReal(Mat mat, const char *name, PetscReal value) { 835*000d2032SJeremy L Thompson double value_double = value; 836*000d2032SJeremy L Thompson 837*000d2032SJeremy L Thompson PetscFunctionBeginUser; 838*000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, name, value_double)); 839*000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 840*000d2032SJeremy L Thompson } 841*000d2032SJeremy L Thompson 842*000d2032SJeremy L Thompson /** 84351bb547fSJames Wright @brief Get the current `PetscReal` value of a context field for a `MatCEED`. 84451bb547fSJames Wright 84551bb547fSJames Wright Not collective across MPI processes. 84651bb547fSJames Wright 84751bb547fSJames Wright @param[in] mat `MatCEED` 84851bb547fSJames Wright @param[in] name Name of the context field 84951bb547fSJames Wright @param[out] value Current context field value 85051bb547fSJames Wright 85151bb547fSJames Wright @return An error code: 0 - success, otherwise - failure 85251bb547fSJames Wright **/ 85351bb547fSJames Wright PetscErrorCode MatCeedGetContextReal(Mat mat, const char *name, PetscReal *value) { 85487d3884fSJeremy L Thompson double value_double = 0.0; 85551bb547fSJames Wright 85651bb547fSJames Wright PetscFunctionBeginUser; 85751bb547fSJames Wright PetscCall(MatCeedGetContextDouble(mat, name, &value_double)); 85851bb547fSJames Wright *value = value_double; 85951bb547fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 86051bb547fSJames Wright } 86151bb547fSJames Wright 86251bb547fSJames Wright /** 863*000d2032SJeremy L Thompson @brief Set the current time for a `MatCEED`. 864*000d2032SJeremy L Thompson 865*000d2032SJeremy L Thompson Not collective across MPI processes. 866*000d2032SJeremy L Thompson 867*000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 868*000d2032SJeremy L Thompson @param[in] time Current time 869*000d2032SJeremy L Thompson 870*000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 871*000d2032SJeremy L Thompson **/ 872*000d2032SJeremy L Thompson PetscErrorCode MatCeedSetTime(Mat mat, PetscReal time) { 873*000d2032SJeremy L Thompson PetscFunctionBeginUser; 874*000d2032SJeremy L Thompson { 875*000d2032SJeremy L Thompson double time_ceed = time; 876*000d2032SJeremy L Thompson 877*000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "time", time_ceed)); 878*000d2032SJeremy L Thompson } 879*000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 880*000d2032SJeremy L Thompson } 881*000d2032SJeremy L Thompson 882*000d2032SJeremy L Thompson /** 883*000d2032SJeremy L Thompson @brief Get the current time for a `MatCEED`. 884*000d2032SJeremy L Thompson 885*000d2032SJeremy L Thompson Not collective across MPI processes. 886*000d2032SJeremy L Thompson 887*000d2032SJeremy L Thompson @param[in] mat `MatCEED` 888*000d2032SJeremy L Thompson @param[out] time Current time, or -1.0 if the boundary evaluator has no time field 889*000d2032SJeremy L Thompson 890*000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 891*000d2032SJeremy L Thompson **/ 892*000d2032SJeremy L Thompson PetscErrorCode MatCeedGetTime(Mat mat, PetscReal *time) { 893*000d2032SJeremy L Thompson PetscFunctionBeginUser; 894*000d2032SJeremy L Thompson *time = -1.0; 895*000d2032SJeremy L Thompson { 896*000d2032SJeremy L Thompson double time_ceed = -1.0; 897*000d2032SJeremy L Thompson 898*000d2032SJeremy L Thompson PetscCall(MatCeedGetContextDouble(mat, "time", &time_ceed)); 899*000d2032SJeremy L Thompson *time = time_ceed; 900*000d2032SJeremy L Thompson } 901*000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 902*000d2032SJeremy L Thompson } 903*000d2032SJeremy L Thompson 904*000d2032SJeremy L Thompson /** 905*000d2032SJeremy L Thompson @brief Set the current time step for a `MatCEED`. 906*000d2032SJeremy L Thompson 907*000d2032SJeremy L Thompson Not collective across MPI processes. 908*000d2032SJeremy L Thompson 909*000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 910*000d2032SJeremy L Thompson @param[in] dt Current time step 911*000d2032SJeremy L Thompson 912*000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 913*000d2032SJeremy L Thompson **/ 914*000d2032SJeremy L Thompson PetscErrorCode MatCeedSetDt(Mat mat, PetscReal dt) { 915*000d2032SJeremy L Thompson PetscFunctionBeginUser; 916*000d2032SJeremy L Thompson { 917*000d2032SJeremy L Thompson double dt_ceed = dt; 918*000d2032SJeremy L Thompson 919*000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "dt", dt_ceed)); 920*000d2032SJeremy L Thompson } 921*000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 922*000d2032SJeremy L Thompson } 923*000d2032SJeremy L Thompson 924*000d2032SJeremy L Thompson /** 925*000d2032SJeremy L Thompson @brief Set the Jacobian shifts for a `MatCEED`. 926*000d2032SJeremy L Thompson 927*000d2032SJeremy L Thompson Not collective across MPI processes. 928*000d2032SJeremy L Thompson 929*000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 930*000d2032SJeremy L Thompson @param[in] shift_v Velocity shift 931*000d2032SJeremy L Thompson @param[in] shift_a Acceleration shift 932*000d2032SJeremy L Thompson 933*000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 934*000d2032SJeremy L Thompson **/ 935*000d2032SJeremy L Thompson PetscErrorCode MatCeedSetShifts(Mat mat, PetscReal shift_v, PetscReal shift_a) { 936*000d2032SJeremy L Thompson PetscFunctionBeginUser; 937*000d2032SJeremy L Thompson { 938*000d2032SJeremy L Thompson double shift_v_ceed = shift_v; 939*000d2032SJeremy L Thompson 940*000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "shift v", shift_v_ceed)); 941*000d2032SJeremy L Thompson } 942*000d2032SJeremy L Thompson if (shift_a) { 943*000d2032SJeremy L Thompson double shift_a_ceed = shift_a; 944*000d2032SJeremy L Thompson 945*000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "shift a", shift_a_ceed)); 946*000d2032SJeremy L Thompson } 947*000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 948*000d2032SJeremy L Thompson } 949*000d2032SJeremy L Thompson 950*000d2032SJeremy L Thompson /** 95158600ac3SJames Wright @brief Set user context for a `MATCEED`. 95258600ac3SJames Wright 95358600ac3SJames Wright Collective across MPI processes. 95458600ac3SJames Wright 95558600ac3SJames Wright @param[in,out] mat `MATCEED` 95658600ac3SJames Wright @param[in] f The context destroy function, or NULL 95758600ac3SJames Wright @param[in] ctx User context, or NULL to unset 95858600ac3SJames Wright 95958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 96058600ac3SJames Wright **/ 96158600ac3SJames Wright PetscErrorCode MatCeedSetContext(Mat mat, PetscErrorCode (*f)(void *), void *ctx) { 96258600ac3SJames Wright PetscContainer user_ctx = NULL; 96358600ac3SJames Wright 96458600ac3SJames Wright PetscFunctionBeginUser; 96558600ac3SJames Wright if (ctx) { 96658600ac3SJames Wright PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)mat), &user_ctx)); 96758600ac3SJames Wright PetscCall(PetscContainerSetPointer(user_ctx, ctx)); 96858600ac3SJames Wright PetscCall(PetscContainerSetUserDestroy(user_ctx, f)); 96958600ac3SJames Wright } 97058600ac3SJames Wright PetscCall(PetscObjectCompose((PetscObject)mat, "MatCeed user context", (PetscObject)user_ctx)); 97158600ac3SJames Wright PetscCall(PetscContainerDestroy(&user_ctx)); 97258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 97358600ac3SJames Wright } 97458600ac3SJames Wright 97558600ac3SJames Wright /** 97658600ac3SJames Wright @brief Retrieve the user context for a `MATCEED`. 97758600ac3SJames Wright 97858600ac3SJames Wright Collective across MPI processes. 97958600ac3SJames Wright 98058600ac3SJames Wright @param[in,out] mat `MATCEED` 98158600ac3SJames Wright @param[in] ctx User context 98258600ac3SJames Wright 98358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 98458600ac3SJames Wright **/ 98558600ac3SJames Wright PetscErrorCode MatCeedGetContext(Mat mat, void *ctx) { 98658600ac3SJames Wright PetscContainer user_ctx; 98758600ac3SJames Wright 98858600ac3SJames Wright PetscFunctionBeginUser; 98958600ac3SJames Wright PetscCall(PetscObjectQuery((PetscObject)mat, "MatCeed user context", (PetscObject *)&user_ctx)); 99058600ac3SJames Wright if (user_ctx) PetscCall(PetscContainerGetPointer(user_ctx, (void **)ctx)); 99158600ac3SJames Wright else *(void **)ctx = NULL; 99258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 99358600ac3SJames Wright } 99458600ac3SJames Wright /** 99540d80af1SJames Wright @brief Set a user defined matrix operation for a `MATCEED` matrix. 99640d80af1SJames Wright 99740d80af1SJames Wright Within each user-defined routine, the user should call `MatCeedGetContext()` to obtain the user-defined context that was set by 99840d80af1SJames Wright `MatCeedSetContext()`. 99958600ac3SJames Wright 100058600ac3SJames Wright Collective across MPI processes. 100158600ac3SJames Wright 100258600ac3SJames Wright @param[in,out] mat `MATCEED` 100340d80af1SJames Wright @param[in] op Name of the `MatOperation` 100440d80af1SJames Wright @param[in] g Function that provides the operation 100558600ac3SJames Wright 100658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 100758600ac3SJames Wright **/ 100840d80af1SJames Wright PetscErrorCode MatCeedSetOperation(Mat mat, MatOperation op, void (*g)(void)) { 100940d80af1SJames Wright PetscFunctionBeginUser; 101040d80af1SJames Wright PetscCall(MatShellSetOperation(mat, op, g)); 101140d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 101240d80af1SJames Wright } 101340d80af1SJames Wright 101440d80af1SJames Wright /** 101540d80af1SJames Wright @brief Sets the default COO matrix type as a string from the `MATCEED`. 101640d80af1SJames Wright 101740d80af1SJames Wright Collective across MPI processes. 101840d80af1SJames Wright 101940d80af1SJames Wright @param[in,out] mat `MATCEED` 102040d80af1SJames Wright @param[in] type COO `MatType` to set 102140d80af1SJames Wright 102240d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 102340d80af1SJames Wright **/ 102440d80af1SJames Wright PetscErrorCode MatCeedSetCOOMatType(Mat mat, MatType type) { 102558600ac3SJames Wright MatCeedContext ctx; 102658600ac3SJames Wright 102758600ac3SJames Wright PetscFunctionBeginUser; 102858600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 102958600ac3SJames Wright // Check if same 103058600ac3SJames Wright { 103158600ac3SJames Wright size_t len_old, len_new; 103258600ac3SJames Wright PetscBool is_same = PETSC_FALSE; 103358600ac3SJames Wright 103440d80af1SJames Wright PetscCall(PetscStrlen(ctx->coo_mat_type, &len_old)); 103558600ac3SJames Wright PetscCall(PetscStrlen(type, &len_new)); 103640d80af1SJames Wright if (len_old == len_new) PetscCall(PetscStrncmp(ctx->coo_mat_type, type, len_old, &is_same)); 103758600ac3SJames Wright if (is_same) PetscFunctionReturn(PETSC_SUCCESS); 103858600ac3SJames Wright } 103958600ac3SJames Wright // Clean up old mats in different format 104058600ac3SJames Wright // LCOV_EXCL_START 104158600ac3SJames Wright if (ctx->mat_assembled_full_internal) { 104258600ac3SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) { 104358600ac3SJames Wright if (ctx->mats_assembled_full[i] == ctx->mat_assembled_full_internal) { 104458600ac3SJames Wright for (PetscInt j = i + 1; j < ctx->num_mats_assembled_full; j++) { 104558600ac3SJames Wright ctx->mats_assembled_full[j - 1] = ctx->mats_assembled_full[j]; 104658600ac3SJames Wright } 104758600ac3SJames Wright ctx->num_mats_assembled_full--; 104858600ac3SJames Wright // Note: we'll realloc this array again, so no need to shrink the allocation 104958600ac3SJames Wright PetscCall(MatDestroy(&ctx->mat_assembled_full_internal)); 105058600ac3SJames Wright } 105158600ac3SJames Wright } 105258600ac3SJames Wright } 105358600ac3SJames Wright if (ctx->mat_assembled_pbd_internal) { 105458600ac3SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) { 105558600ac3SJames Wright if (ctx->mats_assembled_pbd[i] == ctx->mat_assembled_pbd_internal) { 105658600ac3SJames Wright for (PetscInt j = i + 1; j < ctx->num_mats_assembled_pbd; j++) { 105758600ac3SJames Wright ctx->mats_assembled_pbd[j - 1] = ctx->mats_assembled_pbd[j]; 105858600ac3SJames Wright } 105958600ac3SJames Wright // Note: we'll realloc this array again, so no need to shrink the allocation 106058600ac3SJames Wright ctx->num_mats_assembled_pbd--; 106158600ac3SJames Wright PetscCall(MatDestroy(&ctx->mat_assembled_pbd_internal)); 106258600ac3SJames Wright } 106358600ac3SJames Wright } 106458600ac3SJames Wright } 106540d80af1SJames Wright PetscCall(PetscFree(ctx->coo_mat_type)); 106640d80af1SJames Wright PetscCall(PetscStrallocpy(type, &ctx->coo_mat_type)); 106758600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 106858600ac3SJames Wright // LCOV_EXCL_STOP 106958600ac3SJames Wright } 107058600ac3SJames Wright 107158600ac3SJames Wright /** 107240d80af1SJames Wright @brief Gets the default COO matrix type as a string from the `MATCEED`. 107358600ac3SJames Wright 107458600ac3SJames Wright Collective across MPI processes. 107558600ac3SJames Wright 107658600ac3SJames Wright @param[in,out] mat `MATCEED` 107740d80af1SJames Wright @param[in] type COO `MatType` 107858600ac3SJames Wright 107958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 108058600ac3SJames Wright **/ 108140d80af1SJames Wright PetscErrorCode MatCeedGetCOOMatType(Mat mat, MatType *type) { 108258600ac3SJames Wright MatCeedContext ctx; 108358600ac3SJames Wright 108458600ac3SJames Wright PetscFunctionBeginUser; 108558600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 108640d80af1SJames Wright *type = ctx->coo_mat_type; 108758600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 108858600ac3SJames Wright } 108958600ac3SJames Wright 109058600ac3SJames Wright /** 109158600ac3SJames Wright @brief Set input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 109258600ac3SJames Wright 109358600ac3SJames Wright Not collective across MPI processes. 109458600ac3SJames Wright 109558600ac3SJames Wright @param[in,out] mat `MATCEED` 109658600ac3SJames Wright @param[in] X_loc Input PETSc local vector, or NULL 109758600ac3SJames Wright @param[in] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 109858600ac3SJames Wright 109958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 110058600ac3SJames Wright **/ 110158600ac3SJames Wright PetscErrorCode MatCeedSetLocalVectors(Mat mat, Vec X_loc, Vec Y_loc_transpose) { 110258600ac3SJames Wright MatCeedContext ctx; 110358600ac3SJames Wright 110458600ac3SJames Wright PetscFunctionBeginUser; 110558600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 110658600ac3SJames Wright if (X_loc) { 110758600ac3SJames Wright PetscInt len_old, len_new; 110858600ac3SJames Wright 110958600ac3SJames Wright PetscCall(VecGetSize(ctx->X_loc, &len_old)); 111058600ac3SJames Wright PetscCall(VecGetSize(X_loc, &len_new)); 111158600ac3SJames 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, 111258600ac3SJames Wright len_new, len_old); 111340d80af1SJames Wright PetscCall(VecReferenceCopy(X_loc, &ctx->X_loc)); 111458600ac3SJames Wright } 111558600ac3SJames Wright if (Y_loc_transpose) { 111658600ac3SJames Wright PetscInt len_old, len_new; 111758600ac3SJames Wright 111858600ac3SJames Wright PetscCall(VecGetSize(ctx->Y_loc_transpose, &len_old)); 111958600ac3SJames Wright PetscCall(VecGetSize(Y_loc_transpose, &len_new)); 112058600ac3SJames Wright PetscCheck(len_old == len_new, PETSC_COMM_SELF, PETSC_ERR_LIB, 112158600ac3SJames Wright "new Y_loc_transpose length %" PetscInt_FMT " should match old Y_loc_transpose length %" PetscInt_FMT, len_new, len_old); 112240d80af1SJames Wright PetscCall(VecReferenceCopy(Y_loc_transpose, &ctx->Y_loc_transpose)); 112358600ac3SJames Wright } 112458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 112558600ac3SJames Wright } 112658600ac3SJames Wright 112758600ac3SJames Wright /** 112858600ac3SJames Wright @brief Get input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 112958600ac3SJames Wright 113058600ac3SJames Wright Not collective across MPI processes. 113158600ac3SJames Wright 113258600ac3SJames Wright @param[in,out] mat `MATCEED` 113358600ac3SJames Wright @param[out] X_loc Input PETSc local vector, or NULL 113458600ac3SJames Wright @param[out] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 113558600ac3SJames Wright 113658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 113758600ac3SJames Wright **/ 113858600ac3SJames Wright PetscErrorCode MatCeedGetLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) { 113958600ac3SJames Wright MatCeedContext ctx; 114058600ac3SJames Wright 114158600ac3SJames Wright PetscFunctionBeginUser; 114258600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 114358600ac3SJames Wright if (X_loc) { 114440d80af1SJames Wright *X_loc = NULL; 114540d80af1SJames Wright PetscCall(VecReferenceCopy(ctx->X_loc, X_loc)); 114658600ac3SJames Wright } 114758600ac3SJames Wright if (Y_loc_transpose) { 114840d80af1SJames Wright *Y_loc_transpose = NULL; 114940d80af1SJames Wright PetscCall(VecReferenceCopy(ctx->Y_loc_transpose, Y_loc_transpose)); 115058600ac3SJames Wright } 115158600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 115258600ac3SJames Wright } 115358600ac3SJames Wright 115458600ac3SJames Wright /** 115558600ac3SJames Wright @brief Restore input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 115658600ac3SJames Wright 115758600ac3SJames Wright Not collective across MPI processes. 115858600ac3SJames Wright 115958600ac3SJames Wright @param[in,out] mat MatCeed 116058600ac3SJames Wright @param[out] X_loc Input PETSc local vector, or NULL 116158600ac3SJames Wright @param[out] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 116258600ac3SJames Wright 116358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 116458600ac3SJames Wright **/ 116558600ac3SJames Wright PetscErrorCode MatCeedRestoreLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) { 116658600ac3SJames Wright PetscFunctionBeginUser; 116758600ac3SJames Wright if (X_loc) PetscCall(VecDestroy(X_loc)); 116858600ac3SJames Wright if (Y_loc_transpose) PetscCall(VecDestroy(Y_loc_transpose)); 116958600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 117058600ac3SJames Wright } 117158600ac3SJames Wright 117258600ac3SJames Wright /** 117358600ac3SJames Wright @brief Get libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 117458600ac3SJames Wright 117558600ac3SJames Wright Not collective across MPI processes. 117658600ac3SJames Wright 117758600ac3SJames Wright @param[in,out] mat MatCeed 117858600ac3SJames Wright @param[out] op_mult libCEED `CeedOperator` for `MatMult()`, or NULL 117958600ac3SJames Wright @param[out] op_mult_transpose libCEED `CeedOperator` for `MatMultTranspose()`, or NULL 118058600ac3SJames Wright 118158600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 118258600ac3SJames Wright **/ 118358600ac3SJames Wright PetscErrorCode MatCeedGetCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) { 118458600ac3SJames Wright MatCeedContext ctx; 118558600ac3SJames Wright 118658600ac3SJames Wright PetscFunctionBeginUser; 118758600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 118858600ac3SJames Wright if (op_mult) { 118958600ac3SJames Wright *op_mult = NULL; 119050f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult, op_mult)); 119158600ac3SJames Wright } 119258600ac3SJames Wright if (op_mult_transpose) { 119358600ac3SJames Wright *op_mult_transpose = NULL; 119450f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult_transpose, op_mult_transpose)); 119558600ac3SJames Wright } 119658600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 119758600ac3SJames Wright } 119858600ac3SJames Wright 119958600ac3SJames Wright /** 120058600ac3SJames Wright @brief Restore libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 120158600ac3SJames Wright 120258600ac3SJames Wright Not collective across MPI processes. 120358600ac3SJames Wright 120458600ac3SJames Wright @param[in,out] mat MatCeed 120558600ac3SJames Wright @param[out] op_mult libCEED `CeedOperator` for `MatMult()`, or NULL 120658600ac3SJames Wright @param[out] op_mult_transpose libCEED `CeedOperator` for `MatMultTranspose()`, or NULL 120758600ac3SJames Wright 120858600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 120958600ac3SJames Wright **/ 121058600ac3SJames Wright PetscErrorCode MatCeedRestoreCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) { 121158600ac3SJames Wright MatCeedContext ctx; 121258600ac3SJames Wright 121358600ac3SJames Wright PetscFunctionBeginUser; 121458600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 121550f50432SJames Wright if (op_mult) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult)); 121650f50432SJames Wright if (op_mult_transpose) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult_transpose)); 121758600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 121858600ac3SJames Wright } 121958600ac3SJames Wright 122058600ac3SJames Wright /** 122158600ac3SJames Wright @brief Set `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 122258600ac3SJames Wright 122358600ac3SJames Wright Not collective across MPI processes. 122458600ac3SJames Wright 122558600ac3SJames Wright @param[in,out] mat MatCeed 122658600ac3SJames Wright @param[out] log_event_mult `PetscLogEvent` for forward evaluation, or NULL 122758600ac3SJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose evaluation, or NULL 122858600ac3SJames Wright 122958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 123058600ac3SJames Wright **/ 123158600ac3SJames Wright PetscErrorCode MatCeedSetLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) { 123258600ac3SJames Wright MatCeedContext ctx; 123358600ac3SJames Wright 123458600ac3SJames Wright PetscFunctionBeginUser; 123558600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 123658600ac3SJames Wright if (log_event_mult) ctx->log_event_mult = log_event_mult; 123758600ac3SJames Wright if (log_event_mult_transpose) ctx->log_event_mult_transpose = log_event_mult_transpose; 123858600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 123958600ac3SJames Wright } 124058600ac3SJames Wright 124158600ac3SJames Wright /** 124258600ac3SJames Wright @brief Get `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 124358600ac3SJames Wright 124458600ac3SJames Wright Not collective across MPI processes. 124558600ac3SJames Wright 124658600ac3SJames Wright @param[in,out] mat MatCeed 124758600ac3SJames Wright @param[out] log_event_mult `PetscLogEvent` for forward evaluation, or NULL 124858600ac3SJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose evaluation, or NULL 124958600ac3SJames Wright 125058600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 125158600ac3SJames Wright **/ 125258600ac3SJames Wright PetscErrorCode MatCeedGetLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) { 125358600ac3SJames Wright MatCeedContext ctx; 125458600ac3SJames Wright 125558600ac3SJames Wright PetscFunctionBeginUser; 125658600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 125758600ac3SJames Wright if (log_event_mult) *log_event_mult = ctx->log_event_mult; 125858600ac3SJames Wright if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_mult_transpose; 125958600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 126058600ac3SJames Wright } 126158600ac3SJames Wright 1262c63b910fSJames Wright /** 1263c63b910fSJames Wright @brief Set `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 1264c63b910fSJames Wright 1265c63b910fSJames Wright Not collective across MPI processes. 1266c63b910fSJames Wright 1267c63b910fSJames Wright @param[in,out] mat MatCeed 1268c63b910fSJames Wright @param[out] log_event_mult `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL 1269c63b910fSJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL 1270c63b910fSJames Wright 1271c63b910fSJames Wright @return An error code: 0 - success, otherwise - failure 1272c63b910fSJames Wright **/ 1273c63b910fSJames Wright PetscErrorCode MatCeedSetCeedOperatorLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) { 1274c63b910fSJames Wright MatCeedContext ctx; 1275c63b910fSJames Wright 1276c63b910fSJames Wright PetscFunctionBeginUser; 1277c63b910fSJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 1278c63b910fSJames Wright if (log_event_mult) ctx->log_event_ceed_mult = log_event_mult; 1279c63b910fSJames Wright if (log_event_mult_transpose) ctx->log_event_ceed_mult_transpose = log_event_mult_transpose; 1280c63b910fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1281c63b910fSJames Wright } 1282c63b910fSJames Wright 1283c63b910fSJames Wright /** 1284c63b910fSJames Wright @brief Get `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 1285c63b910fSJames Wright 1286c63b910fSJames Wright Not collective across MPI processes. 1287c63b910fSJames Wright 1288c63b910fSJames Wright @param[in,out] mat MatCeed 1289c63b910fSJames Wright @param[out] log_event_mult `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL 1290c63b910fSJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL 1291c63b910fSJames Wright 1292c63b910fSJames Wright @return An error code: 0 - success, otherwise - failure 1293c63b910fSJames Wright **/ 1294c63b910fSJames Wright PetscErrorCode MatCeedGetCeedOperatorLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) { 1295c63b910fSJames Wright MatCeedContext ctx; 1296c63b910fSJames Wright 1297c63b910fSJames Wright PetscFunctionBeginUser; 1298c63b910fSJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 1299c63b910fSJames Wright if (log_event_mult) *log_event_mult = ctx->log_event_ceed_mult; 1300c63b910fSJames Wright if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_ceed_mult_transpose; 1301c63b910fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1302c63b910fSJames Wright } 1303c63b910fSJames Wright 130458600ac3SJames Wright // ----------------------------------------------------------------------------- 130558600ac3SJames Wright // Operator context data 130658600ac3SJames Wright // ----------------------------------------------------------------------------- 130758600ac3SJames Wright 130858600ac3SJames Wright /** 130958600ac3SJames Wright @brief Setup context data for operator application. 131058600ac3SJames Wright 131158600ac3SJames Wright Collective across MPI processes. 131258600ac3SJames Wright 131358600ac3SJames Wright @param[in] dm_x Input `DM` 131458600ac3SJames Wright @param[in] dm_y Output `DM` 131558600ac3SJames Wright @param[in] X_loc Input PETSc local vector, or NULL 131658600ac3SJames Wright @param[in] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 131758600ac3SJames Wright @param[in] op_mult `CeedOperator` for forward evaluation 131858600ac3SJames Wright @param[in] op_mult_transpose `CeedOperator` for transpose evaluation 131958600ac3SJames Wright @param[in] log_event_mult `PetscLogEvent` for forward evaluation 132058600ac3SJames Wright @param[in] log_event_mult_transpose `PetscLogEvent` for transpose evaluation 1321c63b910fSJames Wright @param[in] log_event_ceed_mult `PetscLogEvent` for forward `CeedOperator` evaluation 1322c63b910fSJames Wright @param[in] log_event_ceed_mult_transpose `PetscLogEvent` for transpose `CeedOperator` evaluation 132358600ac3SJames Wright @param[out] ctx Context data for operator evaluation 132458600ac3SJames Wright 132558600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 132658600ac3SJames Wright **/ 132758600ac3SJames Wright PetscErrorCode MatCeedContextCreate(DM dm_x, DM dm_y, Vec X_loc, Vec Y_loc_transpose, CeedOperator op_mult, CeedOperator op_mult_transpose, 1328c63b910fSJames Wright PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose, PetscLogEvent log_event_ceed_mult, 1329c63b910fSJames Wright PetscLogEvent log_event_ceed_mult_transpose, MatCeedContext *ctx) { 133058600ac3SJames Wright CeedSize x_loc_len, y_loc_len; 133158600ac3SJames Wright 133258600ac3SJames Wright PetscFunctionBeginUser; 133358600ac3SJames Wright 133458600ac3SJames Wright // Allocate 133558600ac3SJames Wright PetscCall(PetscNew(ctx)); 133658600ac3SJames Wright (*ctx)->ref_count = 1; 133758600ac3SJames Wright 133858600ac3SJames Wright // Logging 133958600ac3SJames Wright (*ctx)->log_event_mult = log_event_mult; 134058600ac3SJames Wright (*ctx)->log_event_mult_transpose = log_event_mult_transpose; 1341c63b910fSJames Wright (*ctx)->log_event_ceed_mult = log_event_ceed_mult; 1342c63b910fSJames Wright (*ctx)->log_event_ceed_mult_transpose = log_event_ceed_mult_transpose; 134358600ac3SJames Wright 134458600ac3SJames Wright // PETSc objects 134540d80af1SJames Wright PetscCall(DMReferenceCopy(dm_x, &(*ctx)->dm_x)); 134640d80af1SJames Wright PetscCall(DMReferenceCopy(dm_y, &(*ctx)->dm_y)); 134740d80af1SJames Wright if (X_loc) PetscCall(VecReferenceCopy(X_loc, &(*ctx)->X_loc)); 134840d80af1SJames Wright if (Y_loc_transpose) PetscCall(VecReferenceCopy(Y_loc_transpose, &(*ctx)->Y_loc_transpose)); 134958600ac3SJames Wright 135058600ac3SJames Wright // Memtype 135158600ac3SJames Wright { 135258600ac3SJames Wright const PetscScalar *x; 135358600ac3SJames Wright Vec X; 135458600ac3SJames Wright 135558600ac3SJames Wright PetscCall(DMGetLocalVector(dm_x, &X)); 135658600ac3SJames Wright PetscCall(VecGetArrayReadAndMemType(X, &x, &(*ctx)->mem_type)); 135758600ac3SJames Wright PetscCall(VecRestoreArrayReadAndMemType(X, &x)); 135858600ac3SJames Wright PetscCall(DMRestoreLocalVector(dm_x, &X)); 135958600ac3SJames Wright } 136058600ac3SJames Wright 136158600ac3SJames Wright // libCEED objects 136258600ac3SJames Wright PetscCheck(CeedOperatorGetCeed(op_mult, &(*ctx)->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, 136358600ac3SJames Wright "retrieving Ceed context object failed"); 136450f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedReference((*ctx)->ceed)); 136550f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorGetActiveVectorLengths(op_mult, &x_loc_len, &y_loc_len)); 136650f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult, &(*ctx)->op_mult)); 136750f50432SJames Wright if (op_mult_transpose) PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult_transpose, &(*ctx)->op_mult_transpose)); 136850f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, x_loc_len, &(*ctx)->x_loc)); 136950f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, y_loc_len, &(*ctx)->y_loc)); 137058600ac3SJames Wright 137158600ac3SJames Wright // Flop counting 137258600ac3SJames Wright { 137358600ac3SJames Wright CeedSize ceed_flops_estimate = 0; 137458600ac3SJames Wright 137550f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult, &ceed_flops_estimate)); 137658600ac3SJames Wright (*ctx)->flops_mult = ceed_flops_estimate; 137758600ac3SJames Wright if (op_mult_transpose) { 137850f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult_transpose, &ceed_flops_estimate)); 137958600ac3SJames Wright (*ctx)->flops_mult_transpose = ceed_flops_estimate; 138058600ac3SJames Wright } 138158600ac3SJames Wright } 138258600ac3SJames Wright 138358600ac3SJames Wright // Check sizes 138458600ac3SJames Wright if (x_loc_len > 0 || y_loc_len > 0) { 138558600ac3SJames Wright CeedSize ctx_x_loc_len, ctx_y_loc_len; 138658600ac3SJames Wright PetscInt X_loc_len, dm_x_loc_len, Y_loc_len, dm_y_loc_len; 138758600ac3SJames Wright Vec dm_X_loc, dm_Y_loc; 138858600ac3SJames Wright 138958600ac3SJames Wright // -- Input 139058600ac3SJames Wright PetscCall(DMGetLocalVector(dm_x, &dm_X_loc)); 139158600ac3SJames Wright PetscCall(VecGetLocalSize(dm_X_loc, &dm_x_loc_len)); 139258600ac3SJames Wright PetscCall(DMRestoreLocalVector(dm_x, &dm_X_loc)); 139350f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->x_loc, &ctx_x_loc_len)); 13944c17272bSJames Wright if (X_loc) { 13954c17272bSJames Wright PetscCall(VecGetLocalSize(X_loc, &X_loc_len)); 13964c17272bSJames Wright PetscCheck(X_loc_len == dm_x_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, 13974c17272bSJames Wright "X_loc (%" PetscInt_FMT ") must match dm_x (%" PetscInt_FMT ") dimensions", X_loc_len, dm_x_loc_len); 13984c17272bSJames Wright } 13994c17272bSJames 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", 14004c17272bSJames Wright x_loc_len, dm_x_loc_len); 14014c17272bSJames 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 ")", 14024c17272bSJames Wright x_loc_len, ctx_x_loc_len); 140358600ac3SJames Wright 140458600ac3SJames Wright // -- Output 140558600ac3SJames Wright PetscCall(DMGetLocalVector(dm_y, &dm_Y_loc)); 140658600ac3SJames Wright PetscCall(VecGetLocalSize(dm_Y_loc, &dm_y_loc_len)); 140758600ac3SJames Wright PetscCall(DMRestoreLocalVector(dm_y, &dm_Y_loc)); 140850f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->y_loc, &ctx_y_loc_len)); 14094c17272bSJames 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", 14104c17272bSJames Wright ctx_y_loc_len, dm_y_loc_len); 141158600ac3SJames Wright 141258600ac3SJames Wright // -- Transpose 141358600ac3SJames Wright if (Y_loc_transpose) { 141458600ac3SJames Wright PetscCall(VecGetLocalSize(Y_loc_transpose, &Y_loc_len)); 14154c17272bSJames Wright PetscCheck(Y_loc_len == dm_y_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, 14164c17272bSJames Wright "Y_loc_transpose (%" PetscInt_FMT ") must match dm_y (%" PetscInt_FMT ") dimensions", Y_loc_len, dm_y_loc_len); 141758600ac3SJames Wright } 141858600ac3SJames Wright } 141958600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 142058600ac3SJames Wright } 142158600ac3SJames Wright 142258600ac3SJames Wright /** 142358600ac3SJames Wright @brief Increment reference counter for `MATCEED` context. 142458600ac3SJames Wright 142558600ac3SJames Wright Not collective across MPI processes. 142658600ac3SJames Wright 142758600ac3SJames Wright @param[in,out] ctx Context data 142858600ac3SJames Wright 142958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 143058600ac3SJames Wright **/ 143158600ac3SJames Wright PetscErrorCode MatCeedContextReference(MatCeedContext ctx) { 143258600ac3SJames Wright PetscFunctionBeginUser; 143358600ac3SJames Wright ctx->ref_count++; 143458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 143558600ac3SJames Wright } 143658600ac3SJames Wright 143758600ac3SJames Wright /** 143858600ac3SJames Wright @brief Copy reference for `MATCEED`. 143958600ac3SJames Wright Note: If `ctx_copy` is non-null, it is assumed to be a valid pointer to a `MatCeedContext`. 144058600ac3SJames Wright 144158600ac3SJames Wright Not collective across MPI processes. 144258600ac3SJames Wright 144358600ac3SJames Wright @param[in] ctx Context data 144458600ac3SJames Wright @param[out] ctx_copy Copy of pointer to context data 144558600ac3SJames Wright 144658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 144758600ac3SJames Wright **/ 144858600ac3SJames Wright PetscErrorCode MatCeedContextReferenceCopy(MatCeedContext ctx, MatCeedContext *ctx_copy) { 144958600ac3SJames Wright PetscFunctionBeginUser; 145058600ac3SJames Wright PetscCall(MatCeedContextReference(ctx)); 145158600ac3SJames Wright PetscCall(MatCeedContextDestroy(*ctx_copy)); 145258600ac3SJames Wright *ctx_copy = ctx; 145358600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 145458600ac3SJames Wright } 145558600ac3SJames Wright 145658600ac3SJames Wright /** 145758600ac3SJames Wright @brief Destroy context data for operator application. 145858600ac3SJames Wright 145958600ac3SJames Wright Collective across MPI processes. 146058600ac3SJames Wright 146158600ac3SJames Wright @param[in,out] ctx Context data for operator evaluation 146258600ac3SJames Wright 146358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 146458600ac3SJames Wright **/ 146558600ac3SJames Wright PetscErrorCode MatCeedContextDestroy(MatCeedContext ctx) { 146658600ac3SJames Wright PetscFunctionBeginUser; 146758600ac3SJames Wright if (!ctx || --ctx->ref_count > 0) PetscFunctionReturn(PETSC_SUCCESS); 146858600ac3SJames Wright 146958600ac3SJames Wright // PETSc objects 147058600ac3SJames Wright PetscCall(DMDestroy(&ctx->dm_x)); 147158600ac3SJames Wright PetscCall(DMDestroy(&ctx->dm_y)); 147258600ac3SJames Wright PetscCall(VecDestroy(&ctx->X_loc)); 147358600ac3SJames Wright PetscCall(VecDestroy(&ctx->Y_loc_transpose)); 147458600ac3SJames Wright PetscCall(MatDestroy(&ctx->mat_assembled_full_internal)); 147558600ac3SJames Wright PetscCall(MatDestroy(&ctx->mat_assembled_pbd_internal)); 147640d80af1SJames Wright PetscCall(PetscFree(ctx->coo_mat_type)); 147758600ac3SJames Wright PetscCall(PetscFree(ctx->mats_assembled_full)); 147858600ac3SJames Wright PetscCall(PetscFree(ctx->mats_assembled_pbd)); 147958600ac3SJames Wright 148058600ac3SJames Wright // libCEED objects 148150f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->x_loc)); 148250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->y_loc)); 148350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->coo_values_full)); 148450f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorDestroy(&ctx->coo_values_pbd)); 148550f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorDestroy(&ctx->op_mult)); 148650f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorDestroy(&ctx->op_mult_transpose)); 148758600ac3SJames Wright PetscCheck(CeedDestroy(&ctx->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "destroying libCEED context object failed"); 148858600ac3SJames Wright 148958600ac3SJames Wright // Deallocate 149058600ac3SJames Wright ctx->is_destroyed = PETSC_TRUE; // Flag as destroyed in case someone has stale ref 149158600ac3SJames Wright PetscCall(PetscFree(ctx)); 149258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 149358600ac3SJames Wright } 149458600ac3SJames Wright 149558600ac3SJames Wright /** 149658600ac3SJames Wright @brief Compute the diagonal of an operator via libCEED. 149758600ac3SJames Wright 149858600ac3SJames Wright Collective across MPI processes. 149958600ac3SJames Wright 150058600ac3SJames Wright @param[in] A `MATCEED` 150158600ac3SJames Wright @param[out] D Vector holding operator diagonal 150258600ac3SJames Wright 150358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 150458600ac3SJames Wright **/ 150558600ac3SJames Wright PetscErrorCode MatGetDiagonal_Ceed(Mat A, Vec D) { 150658600ac3SJames Wright PetscMemType mem_type; 150758600ac3SJames Wright Vec D_loc; 150858600ac3SJames Wright MatCeedContext ctx; 150958600ac3SJames Wright 151058600ac3SJames Wright PetscFunctionBeginUser; 151158600ac3SJames Wright PetscCall(MatShellGetContext(A, &ctx)); 151258600ac3SJames Wright 151358600ac3SJames Wright // Place PETSc vector in libCEED vector 1514c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL)); 151558600ac3SJames Wright PetscCall(DMGetLocalVector(ctx->dm_x, &D_loc)); 1516a7dac1d5SJames Wright PetscCall(VecPetscToCeed(D_loc, &mem_type, ctx->x_loc)); 151758600ac3SJames Wright 151858600ac3SJames Wright // Compute Diagonal 1519c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL)); 152050f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleDiagonal(ctx->op_mult, ctx->x_loc, CEED_REQUEST_IMMEDIATE)); 1521c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL)); 152258600ac3SJames Wright 152358600ac3SJames Wright // Restore PETSc vector 1524a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(ctx->x_loc, mem_type, D_loc)); 152558600ac3SJames Wright 152658600ac3SJames Wright // Local-to-Global 152758600ac3SJames Wright PetscCall(VecZeroEntries(D)); 152858600ac3SJames Wright PetscCall(DMLocalToGlobal(ctx->dm_x, D_loc, ADD_VALUES, D)); 152958600ac3SJames Wright PetscCall(DMRestoreLocalVector(ctx->dm_x, &D_loc)); 1530c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL)); 153158600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 153258600ac3SJames Wright } 153358600ac3SJames Wright 153458600ac3SJames Wright /** 153558600ac3SJames Wright @brief Compute `A X = Y` for a `MATCEED`. 153658600ac3SJames Wright 153758600ac3SJames Wright Collective across MPI processes. 153858600ac3SJames Wright 153958600ac3SJames Wright @param[in] A `MATCEED` 154058600ac3SJames Wright @param[in] X Input PETSc vector 154158600ac3SJames Wright @param[out] Y Output PETSc vector 154258600ac3SJames Wright 154358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 154458600ac3SJames Wright **/ 154558600ac3SJames Wright PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) { 154658600ac3SJames Wright MatCeedContext ctx; 154758600ac3SJames Wright 154858600ac3SJames Wright PetscFunctionBeginUser; 154958600ac3SJames Wright PetscCall(MatShellGetContext(A, &ctx)); 1550c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_mult, A, X, Y, NULL)); 155158600ac3SJames Wright 155258600ac3SJames Wright { 155358600ac3SJames Wright PetscMemType x_mem_type, y_mem_type; 155458600ac3SJames Wright Vec X_loc = ctx->X_loc, Y_loc; 155558600ac3SJames Wright 155658600ac3SJames Wright // Get local vectors 155758600ac3SJames Wright if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); 155858600ac3SJames Wright PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); 155958600ac3SJames Wright 156058600ac3SJames Wright // Global-to-local 156158600ac3SJames Wright PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc)); 156258600ac3SJames Wright 156358600ac3SJames Wright // Setup libCEED vectors 1564a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, ctx->x_loc)); 156558600ac3SJames Wright PetscCall(VecZeroEntries(Y_loc)); 1566a7dac1d5SJames Wright PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc)); 156758600ac3SJames Wright 156858600ac3SJames Wright // Apply libCEED operator 156958600ac3SJames Wright PetscCall(PetscLogGpuTimeBegin()); 1570c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult, A, X, Y, NULL)); 157150f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult, ctx->x_loc, ctx->y_loc, CEED_REQUEST_IMMEDIATE)); 1572c63b910fSJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult, A, X, Y, NULL)); 157358600ac3SJames Wright PetscCall(PetscLogGpuTimeEnd()); 157458600ac3SJames Wright 157558600ac3SJames Wright // Restore PETSc vectors 1576a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(ctx->x_loc, x_mem_type, X_loc)); 1577a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc)); 157858600ac3SJames Wright 157958600ac3SJames Wright // Local-to-global 158058600ac3SJames Wright PetscCall(VecZeroEntries(Y)); 158158600ac3SJames Wright PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y)); 158258600ac3SJames Wright 158358600ac3SJames Wright // Restore local vectors, as needed 158458600ac3SJames Wright if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); 158558600ac3SJames Wright PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); 158658600ac3SJames Wright } 158758600ac3SJames Wright 158858600ac3SJames Wright // Log flops 158958600ac3SJames Wright if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult)); 159058600ac3SJames Wright else PetscCall(PetscLogFlops(ctx->flops_mult)); 1591c63b910fSJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_mult, A, X, Y, NULL)); 159258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 159358600ac3SJames Wright } 159458600ac3SJames Wright 159558600ac3SJames Wright /** 159658600ac3SJames Wright @brief Compute `A^T Y = X` for a `MATCEED`. 159758600ac3SJames Wright 159858600ac3SJames Wright Collective across MPI processes. 159958600ac3SJames Wright 160058600ac3SJames Wright @param[in] A `MATCEED` 160158600ac3SJames Wright @param[in] Y Input PETSc vector 160258600ac3SJames Wright @param[out] X Output PETSc vector 160358600ac3SJames Wright 160458600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 160558600ac3SJames Wright **/ 160658600ac3SJames Wright PetscErrorCode MatMultTranspose_Ceed(Mat A, Vec Y, Vec X) { 160758600ac3SJames Wright MatCeedContext ctx; 160858600ac3SJames Wright 160958600ac3SJames Wright PetscFunctionBeginUser; 161058600ac3SJames Wright PetscCall(MatShellGetContext(A, &ctx)); 1611c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_mult_transpose, A, Y, X, NULL)); 161258600ac3SJames Wright 161358600ac3SJames Wright { 161458600ac3SJames Wright PetscMemType x_mem_type, y_mem_type; 161558600ac3SJames Wright Vec X_loc, Y_loc = ctx->Y_loc_transpose; 161658600ac3SJames Wright 161758600ac3SJames Wright // Get local vectors 161858600ac3SJames Wright if (!ctx->Y_loc_transpose) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); 161958600ac3SJames Wright PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); 162058600ac3SJames Wright 162158600ac3SJames Wright // Global-to-local 162258600ac3SJames Wright PetscCall(DMGlobalToLocal(ctx->dm_y, Y, INSERT_VALUES, Y_loc)); 162358600ac3SJames Wright 162458600ac3SJames Wright // Setup libCEED vectors 1625a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc)); 162658600ac3SJames Wright PetscCall(VecZeroEntries(X_loc)); 1627a7dac1d5SJames Wright PetscCall(VecPetscToCeed(X_loc, &x_mem_type, ctx->x_loc)); 162858600ac3SJames Wright 162958600ac3SJames Wright // Apply libCEED operator 163058600ac3SJames Wright PetscCall(PetscLogGpuTimeBegin()); 1631c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL)); 163250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult_transpose, ctx->y_loc, ctx->x_loc, CEED_REQUEST_IMMEDIATE)); 1633c63b910fSJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL)); 163458600ac3SJames Wright PetscCall(PetscLogGpuTimeEnd()); 163558600ac3SJames Wright 163658600ac3SJames Wright // Restore PETSc vectors 1637a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc)); 1638a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(ctx->x_loc, x_mem_type, X_loc)); 163958600ac3SJames Wright 164058600ac3SJames Wright // Local-to-global 164158600ac3SJames Wright PetscCall(VecZeroEntries(X)); 164258600ac3SJames Wright PetscCall(DMLocalToGlobal(ctx->dm_x, X_loc, ADD_VALUES, X)); 164358600ac3SJames Wright 164458600ac3SJames Wright // Restore local vectors, as needed 164558600ac3SJames Wright if (!ctx->Y_loc_transpose) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); 164658600ac3SJames Wright PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); 164758600ac3SJames Wright } 164858600ac3SJames Wright 164958600ac3SJames Wright // Log flops 165058600ac3SJames Wright if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult_transpose)); 165158600ac3SJames Wright else PetscCall(PetscLogFlops(ctx->flops_mult_transpose)); 1652c63b910fSJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_mult_transpose, A, Y, X, NULL)); 165358600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 165458600ac3SJames Wright } 1655