1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 358600ac3SJames Wright /// @file 440d80af1SJames Wright /// MatCEED implementation 558600ac3SJames Wright 658600ac3SJames Wright #include <ceed.h> 758600ac3SJames Wright #include <ceed/backend.h> 858600ac3SJames Wright #include <mat-ceed-impl.h> 958600ac3SJames Wright #include <mat-ceed.h> 1040d80af1SJames Wright #include <petsc-ceed-utils.h> 1140d80af1SJames Wright #include <petsc-ceed.h> 12000d2032SJeremy L Thompson #include <petscdm.h> 13000d2032SJeremy L Thompson #include <petscmat.h> 14000d2032SJeremy L Thompson #include <stdbool.h> 15000d2032SJeremy L Thompson #include <stdio.h> 1658600ac3SJames Wright #include <stdlib.h> 1758600ac3SJames Wright #include <string.h> 1858600ac3SJames Wright 1958600ac3SJames Wright PetscClassId MATCEED_CLASSID; 20c63b910fSJames Wright PetscLogEvent MATCEED_MULT, MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE, MATCEED_MULT_TRANSPOSE_CEEDOP, MATCEED_ASSEMBLE_DIAGONAL, 21c63b910fSJames Wright MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, MATCEED_SETUP_PBDIAGONAL, MATCEED_SETUP_PBDIAGONAL_CEEDOP, MATCEED_ASSEMBLE_PBDIAGONAL, 22c63b910fSJames Wright MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, MATCEED_SETUP_FULL, MATCEED_SETUP_FULL_CEEDOP, MATCEED_ASSEMBLE_FULL, MATCEED_ASSEMBLE_FULL_CEEDOP; 2358600ac3SJames Wright 2458600ac3SJames Wright /** 2558600ac3SJames Wright @brief Register MATCEED log events. 2658600ac3SJames Wright 2758600ac3SJames Wright Not collective across MPI processes. 2858600ac3SJames Wright 2958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 3058600ac3SJames Wright **/ 3158600ac3SJames Wright static PetscErrorCode MatCeedRegisterLogEvents() { 3240d80af1SJames Wright static PetscBool registered = PETSC_FALSE; 3358600ac3SJames Wright 3458600ac3SJames Wright PetscFunctionBeginUser; 3558600ac3SJames Wright if (registered) PetscFunctionReturn(PETSC_SUCCESS); 36c63b910fSJames Wright PetscCall(PetscClassIdRegister("MatCEED", &MATCEED_CLASSID)); 37c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMul", MATCEED_CLASSID, &MATCEED_MULT)); 38c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMulCeed", MATCEED_CLASSID, &MATCEED_MULT_CEEDOP)); 39c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMulT", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE)); 40c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDMulTCeed", MATCEED_CLASSID, &MATCEED_MULT_TRANSPOSE_CEEDOP)); 41c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmDiag", MATCEED_CLASSID, &MATCEED_ASSEMBLE_DIAGONAL)); 42c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSU", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL)); 43c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBDSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_PBDIAGONAL_CEEDOP)); 44c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBD", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL)); 45c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmPBDCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP)); 46c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmSU", MATCEED_CLASSID, &MATCEED_SETUP_FULL)); 47c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmSUCeed", MATCEED_CLASSID, &MATCEED_SETUP_FULL_CEEDOP)); 48c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsm", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL)); 49c63b910fSJames Wright PetscCall(PetscLogEventRegister("MatCEEDAsmCeed", MATCEED_CLASSID, &MATCEED_ASSEMBLE_FULL_CEEDOP)); 5040d80af1SJames Wright registered = PETSC_TRUE; 5158600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5258600ac3SJames Wright } 5358600ac3SJames Wright 5458600ac3SJames Wright /** 5558600ac3SJames Wright @brief Assemble the point block diagonal of a `MATCEED` into a `MATAIJ` or similar. 5658600ac3SJames Wright The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`. 5758600ac3SJames Wright The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail. 5858600ac3SJames Wright 5958600ac3SJames Wright Collective across MPI processes. 6058600ac3SJames Wright 6158600ac3SJames Wright @param[in] mat_ceed `MATCEED` to assemble 6258600ac3SJames Wright @param[in,out] mat_coo `MATAIJ` or similar to assemble into 6358600ac3SJames Wright 6458600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 6558600ac3SJames Wright **/ 6658600ac3SJames Wright static PetscErrorCode MatCeedAssemblePointBlockDiagonalCOO(Mat mat_ceed, Mat mat_coo) { 6758600ac3SJames Wright MatCeedContext ctx; 6858600ac3SJames Wright 6958600ac3SJames Wright PetscFunctionBeginUser; 7058600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 7158600ac3SJames Wright 7258600ac3SJames Wright // Check if COO pattern set 7358600ac3SJames Wright { 7458600ac3SJames Wright PetscInt index = -1; 7558600ac3SJames Wright 7658600ac3SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) { 7758600ac3SJames Wright if (ctx->mats_assembled_pbd[i] == mat_coo) index = i; 7858600ac3SJames Wright } 7958600ac3SJames Wright if (index == -1) { 8058600ac3SJames Wright PetscInt *rows_petsc = NULL, *cols_petsc = NULL; 8158600ac3SJames Wright CeedInt *rows_ceed, *cols_ceed; 8258600ac3SJames Wright PetscCount num_entries; 8358600ac3SJames Wright PetscLogStage stage_amg_setup; 8458600ac3SJames Wright 8558600ac3SJames Wright // -- Assemble sparsity pattern if mat hasn't been assembled before 86c63b910fSJames Wright PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup)); 8758600ac3SJames Wright if (stage_amg_setup == -1) { 88c63b910fSJames Wright PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup)); 8958600ac3SJames Wright } 9058600ac3SJames Wright PetscCall(PetscLogStagePush(stage_amg_setup)); 91c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 92c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 9350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonalSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed)); 94c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 95a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc)); 96a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc)); 9758600ac3SJames Wright PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc)); 9858600ac3SJames Wright free(rows_petsc); 9958600ac3SJames Wright free(cols_petsc); 10050f50432SJames Wright if (!ctx->coo_values_pbd) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_pbd)); 10158600ac3SJames Wright PetscCall(PetscRealloc(++ctx->num_mats_assembled_pbd * sizeof(Mat), &ctx->mats_assembled_pbd)); 10258600ac3SJames Wright ctx->mats_assembled_pbd[ctx->num_mats_assembled_pbd - 1] = mat_coo; 103c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 10458600ac3SJames Wright PetscCall(PetscLogStagePop()); 10558600ac3SJames Wright } 10658600ac3SJames Wright } 10758600ac3SJames Wright 10858600ac3SJames Wright // Assemble mat_ceed 109c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 11058600ac3SJames Wright PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY)); 11158600ac3SJames Wright { 11258600ac3SJames Wright const CeedScalar *values; 113ed5c6999SJames Wright PetscMemType mem_type_petsc; 11458600ac3SJames Wright CeedMemType mem_type = CEED_MEM_HOST; 11558600ac3SJames Wright PetscBool is_spd, is_spd_known; 11658600ac3SJames Wright 117*2da92326SJames Wright PetscCall(MatGetMemTypeFromVecType(mat_coo, &mem_type_petsc)); 118ed5c6999SJames Wright mem_type = MemTypePetscToCeed(mem_type_petsc); 119c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 12050f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemblePointBlockDiagonal(ctx->op_mult, ctx->coo_values_pbd, CEED_REQUEST_IMMEDIATE)); 121c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 12250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_pbd, mem_type, &values)); 12358600ac3SJames Wright PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES)); 12458600ac3SJames Wright PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd)); 12558600ac3SJames Wright if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd)); 12650f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_pbd, &values)); 12758600ac3SJames Wright } 12858600ac3SJames Wright PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY)); 129c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_PBDIAGONAL, mat_ceed, mat_coo, NULL, NULL)); 13058600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 13158600ac3SJames Wright } 13258600ac3SJames Wright 13358600ac3SJames Wright /** 13458600ac3SJames Wright @brief Assemble inner `Mat` for diagonal `PC` operations 13558600ac3SJames Wright 13658600ac3SJames Wright Collective across MPI processes. 13758600ac3SJames Wright 13858600ac3SJames Wright @param[in] mat_ceed `MATCEED` to invert 13958600ac3SJames Wright @param[in] use_ceed_pbd Boolean flag to use libCEED PBD assembly 14058600ac3SJames Wright @param[out] mat_inner Inner `Mat` for diagonal operations 14158600ac3SJames Wright 14258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 14358600ac3SJames Wright **/ 14458600ac3SJames Wright static PetscErrorCode MatCeedAssembleInnerBlockDiagonalMat(Mat mat_ceed, PetscBool use_ceed_pbd, Mat *mat_inner) { 14558600ac3SJames Wright MatCeedContext ctx; 14658600ac3SJames Wright 14758600ac3SJames Wright PetscFunctionBeginUser; 14858600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 14958600ac3SJames Wright if (use_ceed_pbd) { 15058600ac3SJames Wright // Check if COO pattern set 15140d80af1SJames Wright if (!ctx->mat_assembled_pbd_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_pbd_internal)); 15258600ac3SJames Wright 15358600ac3SJames Wright // Assemble mat_assembled_full_internal 15458600ac3SJames Wright PetscCall(MatCeedAssemblePointBlockDiagonalCOO(mat_ceed, ctx->mat_assembled_pbd_internal)); 15558600ac3SJames Wright if (mat_inner) *mat_inner = ctx->mat_assembled_pbd_internal; 15658600ac3SJames Wright } else { 15758600ac3SJames Wright // Check if COO pattern set 15840d80af1SJames Wright if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal)); 15958600ac3SJames Wright 16058600ac3SJames Wright // Assemble mat_assembled_full_internal 16158600ac3SJames Wright PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal)); 16258600ac3SJames Wright if (mat_inner) *mat_inner = ctx->mat_assembled_full_internal; 16358600ac3SJames Wright } 16458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 16558600ac3SJames Wright } 16658600ac3SJames Wright 16758600ac3SJames Wright /** 168e8ff1987SJames Wright @brief Get `MATCEED` variable block diagonal for Jacobi. 169e8ff1987SJames Wright 170e8ff1987SJames Wright Collective across MPI processes. 171e8ff1987SJames Wright 172e8ff1987SJames Wright @param[in] mat_ceed `MATCEED` to invert 173e8ff1987SJames Wright @param[out] mat_vblock The variable diagonal block matrix 174e8ff1987SJames Wright 175e8ff1987SJames Wright @return An error code: 0 - success, otherwise - failure 176e8ff1987SJames Wright **/ 177e8ff1987SJames Wright static PetscErrorCode MatGetVariableBlockDiagonal_Ceed(Mat mat_ceed, Mat *mat_vblock) { 178e8ff1987SJames Wright MatCeedContext ctx; 179e8ff1987SJames Wright 180e8ff1987SJames Wright PetscFunctionBeginUser; 181e8ff1987SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 182e8ff1987SJames Wright 183e8ff1987SJames Wright // Assemble inner mat if needed 184e8ff1987SJames Wright PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_vpbd_valid, mat_vblock)); 185e8ff1987SJames Wright PetscCall(PetscObjectReference((PetscObject)*mat_vblock)); 186e8ff1987SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 187e8ff1987SJames Wright } 188e8ff1987SJames Wright 189e8ff1987SJames Wright /** 190e8ff1987SJames Wright @brief Get `MATCEED` block diagonal for Jacobi. 191e8ff1987SJames Wright 192e8ff1987SJames Wright Collective across MPI processes. 193e8ff1987SJames Wright 194e8ff1987SJames Wright @param[in] mat_ceed `MATCEED` to invert 195e8ff1987SJames Wright @param[out] mat_block The variable diagonal block matrix 196e8ff1987SJames Wright 197e8ff1987SJames Wright @return An error code: 0 - success, otherwise - failure 198e8ff1987SJames Wright **/ 199e8ff1987SJames Wright static PetscErrorCode MatGetBlockDiagonal_Ceed(Mat mat_ceed, Mat *mat_block) { 200e8ff1987SJames Wright MatCeedContext ctx; 201e8ff1987SJames Wright 202e8ff1987SJames Wright PetscFunctionBeginUser; 203e8ff1987SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 204e8ff1987SJames Wright 205e8ff1987SJames Wright // Assemble inner mat if needed 206e8ff1987SJames Wright PetscCall(MatCeedAssembleInnerBlockDiagonalMat(mat_ceed, ctx->is_ceed_pbd_valid, mat_block)); 207e8ff1987SJames Wright PetscCall(PetscObjectReference((PetscObject)*mat_block)); 208e8ff1987SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 209e8ff1987SJames Wright } 210e8ff1987SJames Wright 211e8ff1987SJames Wright /** 212c1bdbf00SJames Wright @brief Get on-process diagonal block of `MATCEED` 213c1bdbf00SJames Wright 214c1bdbf00SJames Wright This is a block per-process of the diagonal of the global matrix. 215c1bdbf00SJames Wright This is NOT the diagonal blocks associated with the block size of the matrix (i.e. `MatSetBlockSize()` has no effect on this function). 21658600ac3SJames Wright 21758600ac3SJames Wright Collective across MPI processes. 21858600ac3SJames Wright 21958600ac3SJames Wright @param[in] mat_ceed `MATCEED` to invert 22058600ac3SJames Wright @param[out] mat_block The diagonal block matrix 22158600ac3SJames Wright 22258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 22358600ac3SJames Wright **/ 22458600ac3SJames Wright static PetscErrorCode MatGetDiagonalBlock_Ceed(Mat mat_ceed, Mat *mat_block) { 22558600ac3SJames Wright MatCeedContext ctx; 22658600ac3SJames Wright 22758600ac3SJames Wright PetscFunctionBeginUser; 22858600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 22958600ac3SJames Wright 230c1bdbf00SJames Wright // Check if COO pattern set 231c1bdbf00SJames Wright if (!ctx->mat_assembled_full_internal) PetscCall(MatCeedCreateMatCOO(mat_ceed, &ctx->mat_assembled_full_internal)); 23258600ac3SJames Wright 233c1bdbf00SJames Wright // Assemble mat_assembled_full_internal 234c1bdbf00SJames Wright PetscCall(MatCeedAssembleCOO(mat_ceed, ctx->mat_assembled_full_internal)); 235c1bdbf00SJames Wright 236c1bdbf00SJames Wright // Get diagonal block 237c1bdbf00SJames Wright PetscCall(MatGetDiagonalBlock(ctx->mat_assembled_full_internal, mat_block)); 23858600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 23958600ac3SJames Wright } 24058600ac3SJames Wright 24158600ac3SJames Wright /** 242e90c2ceeSJames Wright @brief View `MATCEED`. 243e90c2ceeSJames Wright 244e90c2ceeSJames Wright Collective across MPI processes. 245e90c2ceeSJames Wright 246e90c2ceeSJames Wright @param[in] mat_ceed `MATCEED` to view 247e90c2ceeSJames Wright @param[in] viewer The visualization context 248e90c2ceeSJames Wright 249e90c2ceeSJames Wright @return An error code: 0 - success, otherwise - failure 250e90c2ceeSJames Wright **/ 251e90c2ceeSJames Wright static PetscErrorCode MatView_Ceed(Mat mat_ceed, PetscViewer viewer) { 252e90c2ceeSJames Wright PetscBool is_ascii; 253e90c2ceeSJames Wright PetscViewerFormat format; 254000d2032SJeremy L Thompson PetscMPIInt size, rank; 255e90c2ceeSJames Wright MatCeedContext ctx; 256e90c2ceeSJames Wright 257e90c2ceeSJames Wright PetscFunctionBeginUser; 258e90c2ceeSJames Wright PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 259e90c2ceeSJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 260e90c2ceeSJames Wright if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat_ceed), &viewer)); 261e90c2ceeSJames Wright 262e90c2ceeSJames Wright PetscCall(PetscViewerGetFormat(viewer, &format)); 263e90c2ceeSJames Wright PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat_ceed), &size)); 264e90c2ceeSJames Wright if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS); 265e90c2ceeSJames Wright 266000d2032SJeremy L Thompson PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat_ceed), &rank)); 267000d2032SJeremy L Thompson if (rank != 0) PetscFunctionReturn(PETSC_SUCCESS); 268000d2032SJeremy L Thompson 269e90c2ceeSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &is_ascii)); 270e90c2ceeSJames Wright { 271000d2032SJeremy L Thompson PetscBool is_detailed = format == PETSC_VIEWER_ASCII_INFO_DETAIL; 272000d2032SJeremy L Thompson char rank_string[16] = {'\0'}; 273e90c2ceeSJames Wright FILE *file; 274e90c2ceeSJames Wright 275dfdd2b92SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "MatCEED:\n")); 276537ec908SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // MatCEED 277dfdd2b92SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "Default COO MatType: %s\n", ctx->coo_mat_type)); 278000d2032SJeremy L Thompson PetscCall(PetscSNPrintf(rank_string, 16, "on Rank %d", rank)); 279000d2032SJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(viewer, "CeedOperator Apply %s:\n", is_detailed ? rank_string : "Summary")); 280e8ff1987SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "libCEED PB Diagonal Assembly: %s\n", ctx->is_ceed_pbd_valid ? "True" : "False")); 281e8ff1987SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "libCEED VPB Diagonal Assembly: %s\n", ctx->is_ceed_vpbd_valid ? "True" : "False")); 282537ec908SJames Wright PetscCall(PetscViewerASCIIGetPointer(viewer, &file)); 283537ec908SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // CeedOperator 284000d2032SJeremy L Thompson if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult, file)); 285000d2032SJeremy L Thompson else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult, file)); 286537ec908SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // CeedOperator 287e90c2ceeSJames Wright if (ctx->op_mult_transpose) { 288000d2032SJeremy L Thompson PetscCall(PetscViewerASCIIPrintf(viewer, "CeedOperator ApplyTranspose %s:\n", is_detailed ? rank_string : "Summary")); 289537ec908SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // CeedOperator 290000d2032SJeremy L Thompson if (is_detailed) PetscCallCeed(ctx->ceed, CeedOperatorView(ctx->op_mult_transpose, file)); 291000d2032SJeremy L Thompson else PetscCallCeed(ctx->ceed, CeedOperatorViewTerse(ctx->op_mult_transpose, file)); 292537ec908SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // CeedOperator 293e90c2ceeSJames Wright } 294537ec908SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // MatCEED 295e90c2ceeSJames Wright } 296e90c2ceeSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 297e90c2ceeSJames Wright } 298e90c2ceeSJames Wright 29958600ac3SJames Wright // ----------------------------------------------------------------------------- 30058600ac3SJames Wright // MatCeed 30158600ac3SJames Wright // ----------------------------------------------------------------------------- 30258600ac3SJames Wright 30358600ac3SJames Wright /** 30458600ac3SJames Wright @brief Create PETSc `Mat` from libCEED operators. 30558600ac3SJames Wright 30658600ac3SJames Wright Collective across MPI processes. 30758600ac3SJames Wright 30858600ac3SJames Wright @param[in] dm_x Input `DM` 30958600ac3SJames Wright @param[in] dm_y Output `DM` 31058600ac3SJames Wright @param[in] op_mult `CeedOperator` for forward evaluation 31158600ac3SJames Wright @param[in] op_mult_transpose `CeedOperator` for transpose evaluation 31258600ac3SJames Wright @param[out] mat New MatCeed 31358600ac3SJames Wright 31458600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 31558600ac3SJames Wright **/ 316000d2032SJeremy L Thompson PetscErrorCode MatCreateCeed(DM dm_x, DM dm_y, CeedOperator op_mult, CeedOperator op_mult_transpose, Mat *mat) { 31758600ac3SJames Wright PetscInt X_l_size, X_g_size, Y_l_size, Y_g_size; 31858600ac3SJames Wright VecType vec_type; 31958600ac3SJames Wright MatCeedContext ctx; 32058600ac3SJames Wright 32158600ac3SJames Wright PetscFunctionBeginUser; 32258600ac3SJames Wright PetscCall(MatCeedRegisterLogEvents()); 32358600ac3SJames Wright 32458600ac3SJames Wright // Collect context data 32558600ac3SJames Wright PetscCall(DMGetVecType(dm_x, &vec_type)); 32658600ac3SJames Wright { 32758600ac3SJames Wright Vec X; 32858600ac3SJames Wright 32958600ac3SJames Wright PetscCall(DMGetGlobalVector(dm_x, &X)); 33058600ac3SJames Wright PetscCall(VecGetSize(X, &X_g_size)); 33158600ac3SJames Wright PetscCall(VecGetLocalSize(X, &X_l_size)); 33258600ac3SJames Wright PetscCall(DMRestoreGlobalVector(dm_x, &X)); 33358600ac3SJames Wright } 33458600ac3SJames Wright if (dm_y) { 33558600ac3SJames Wright Vec Y; 33658600ac3SJames Wright 33758600ac3SJames Wright PetscCall(DMGetGlobalVector(dm_y, &Y)); 33858600ac3SJames Wright PetscCall(VecGetSize(Y, &Y_g_size)); 33958600ac3SJames Wright PetscCall(VecGetLocalSize(Y, &Y_l_size)); 34058600ac3SJames Wright PetscCall(DMRestoreGlobalVector(dm_y, &Y)); 34158600ac3SJames Wright } else { 34258600ac3SJames Wright dm_y = dm_x; 34358600ac3SJames Wright Y_g_size = X_g_size; 34458600ac3SJames Wright Y_l_size = X_l_size; 34558600ac3SJames Wright } 34640d80af1SJames Wright 34758600ac3SJames Wright // Create context 34858600ac3SJames Wright { 34958600ac3SJames Wright Vec X_loc, Y_loc_transpose = NULL; 35058600ac3SJames Wright 35158600ac3SJames Wright PetscCall(DMCreateLocalVector(dm_x, &X_loc)); 35258600ac3SJames Wright PetscCall(VecZeroEntries(X_loc)); 35358600ac3SJames Wright if (op_mult_transpose) { 35458600ac3SJames Wright PetscCall(DMCreateLocalVector(dm_y, &Y_loc_transpose)); 35558600ac3SJames Wright PetscCall(VecZeroEntries(Y_loc_transpose)); 35658600ac3SJames Wright } 357c63b910fSJames Wright PetscCall(MatCeedContextCreate(dm_x, dm_y, X_loc, Y_loc_transpose, op_mult, op_mult_transpose, MATCEED_MULT, MATCEED_MULT_TRANSPOSE, 358c63b910fSJames Wright MATCEED_MULT_CEEDOP, MATCEED_MULT_TRANSPOSE_CEEDOP, &ctx)); 35958600ac3SJames Wright PetscCall(VecDestroy(&X_loc)); 36058600ac3SJames Wright PetscCall(VecDestroy(&Y_loc_transpose)); 36158600ac3SJames Wright } 36258600ac3SJames Wright 36358600ac3SJames Wright // Create mat 36458600ac3SJames Wright PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm_x), Y_l_size, X_l_size, Y_g_size, X_g_size, ctx, mat)); 36558600ac3SJames Wright PetscCall(PetscObjectChangeTypeName((PetscObject)*mat, MATCEED)); 36658600ac3SJames Wright // -- Set block and variable block sizes 36758600ac3SJames Wright if (dm_x == dm_y) { 36858600ac3SJames Wright MatType dm_mat_type, dm_mat_type_copy; 36958600ac3SJames Wright Mat temp_mat; 37058600ac3SJames Wright 37158600ac3SJames Wright PetscCall(DMGetMatType(dm_x, &dm_mat_type)); 37258600ac3SJames Wright PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy)); 37358600ac3SJames Wright PetscCall(DMSetMatType(dm_x, MATAIJ)); 37458600ac3SJames Wright PetscCall(DMCreateMatrix(dm_x, &temp_mat)); 37558600ac3SJames Wright PetscCall(DMSetMatType(dm_x, dm_mat_type_copy)); 37658600ac3SJames Wright PetscCall(PetscFree(dm_mat_type_copy)); 37758600ac3SJames Wright 37858600ac3SJames Wright { 37958600ac3SJames Wright PetscInt block_size, num_blocks, max_vblock_size = PETSC_INT_MAX; 38058600ac3SJames Wright const PetscInt *vblock_sizes; 38158600ac3SJames Wright 38258600ac3SJames Wright // -- Get block sizes 38358600ac3SJames Wright PetscCall(MatGetBlockSize(temp_mat, &block_size)); 38458600ac3SJames Wright PetscCall(MatGetVariableBlockSizes(temp_mat, &num_blocks, &vblock_sizes)); 38558600ac3SJames Wright { 38658600ac3SJames Wright PetscInt local_min_max[2] = {0}, global_min_max[2] = {0, PETSC_INT_MAX}; 38758600ac3SJames Wright 38858600ac3SJames Wright for (PetscInt i = 0; i < num_blocks; i++) local_min_max[1] = PetscMax(local_min_max[1], vblock_sizes[i]); 38958600ac3SJames Wright PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_min_max, global_min_max)); 39058600ac3SJames Wright max_vblock_size = global_min_max[1]; 39158600ac3SJames Wright } 39258600ac3SJames Wright 39358600ac3SJames Wright // -- Copy block sizes 39458600ac3SJames Wright if (block_size > 1) PetscCall(MatSetBlockSize(*mat, block_size)); 39558600ac3SJames Wright if (num_blocks) PetscCall(MatSetVariableBlockSizes(*mat, num_blocks, (PetscInt *)vblock_sizes)); 39658600ac3SJames Wright 39758600ac3SJames Wright // -- Check libCEED compatibility 39858600ac3SJames Wright { 39958600ac3SJames Wright bool is_composite; 40058600ac3SJames Wright 40158600ac3SJames Wright ctx->is_ceed_pbd_valid = PETSC_TRUE; 40258600ac3SJames Wright ctx->is_ceed_vpbd_valid = PETSC_TRUE; 40350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorIsComposite(op_mult, &is_composite)); 40458600ac3SJames Wright if (is_composite) { 40558600ac3SJames Wright CeedInt num_sub_operators; 40658600ac3SJames Wright CeedOperator *sub_operators; 40758600ac3SJames Wright 40850f50432SJames Wright PetscCallCeed(ctx->ceed, CeedCompositeOperatorGetNumSub(op_mult, &num_sub_operators)); 40950f50432SJames Wright PetscCallCeed(ctx->ceed, CeedCompositeOperatorGetSubList(op_mult, &sub_operators)); 41058600ac3SJames Wright for (CeedInt i = 0; i < num_sub_operators; i++) { 41158600ac3SJames Wright CeedInt num_bases, num_comp; 41258600ac3SJames Wright CeedBasis *active_bases; 41358600ac3SJames Wright CeedOperatorAssemblyData assembly_data; 41458600ac3SJames Wright 41550f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(sub_operators[i], &assembly_data)); 41650f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL)); 41750f50432SJames Wright PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp)); 41858600ac3SJames Wright if (num_bases > 1) { 41958600ac3SJames Wright ctx->is_ceed_pbd_valid = PETSC_FALSE; 42058600ac3SJames Wright ctx->is_ceed_vpbd_valid = PETSC_FALSE; 42158600ac3SJames Wright } 42258600ac3SJames Wright if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE; 42358600ac3SJames Wright if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE; 42458600ac3SJames Wright } 42558600ac3SJames Wright } else { 42658600ac3SJames Wright // LCOV_EXCL_START 42758600ac3SJames Wright CeedInt num_bases, num_comp; 42858600ac3SJames Wright CeedBasis *active_bases; 42958600ac3SJames Wright CeedOperatorAssemblyData assembly_data; 43058600ac3SJames Wright 43150f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetOperatorAssemblyData(op_mult, &assembly_data)); 43250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorAssemblyDataGetBases(assembly_data, &num_bases, &active_bases, NULL, NULL, NULL, NULL)); 43350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedBasisGetNumComponents(active_bases[0], &num_comp)); 43458600ac3SJames Wright if (num_bases > 1) { 43558600ac3SJames Wright ctx->is_ceed_pbd_valid = PETSC_FALSE; 43658600ac3SJames Wright ctx->is_ceed_vpbd_valid = PETSC_FALSE; 43758600ac3SJames Wright } 43858600ac3SJames Wright if (num_comp != block_size) ctx->is_ceed_pbd_valid = PETSC_FALSE; 43958600ac3SJames Wright if (num_comp < max_vblock_size) ctx->is_ceed_vpbd_valid = PETSC_FALSE; 44058600ac3SJames Wright // LCOV_EXCL_STOP 44158600ac3SJames Wright } 44258600ac3SJames Wright { 44358600ac3SJames Wright PetscInt local_is_valid[2], global_is_valid[2]; 44458600ac3SJames Wright 44558600ac3SJames Wright local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_pbd_valid; 44658600ac3SJames Wright PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid)); 44758600ac3SJames Wright ctx->is_ceed_pbd_valid = global_is_valid[0]; 44858600ac3SJames Wright local_is_valid[0] = local_is_valid[1] = ctx->is_ceed_vpbd_valid; 44958600ac3SJames Wright PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)dm_x), local_is_valid, global_is_valid)); 45058600ac3SJames Wright ctx->is_ceed_vpbd_valid = global_is_valid[0]; 45158600ac3SJames Wright } 45258600ac3SJames Wright } 45358600ac3SJames Wright } 45458600ac3SJames Wright PetscCall(MatDestroy(&temp_mat)); 45558600ac3SJames Wright } 45658600ac3SJames Wright // -- Set internal mat type 45758600ac3SJames Wright { 45858600ac3SJames Wright VecType vec_type; 45940d80af1SJames Wright MatType coo_mat_type; 46058600ac3SJames Wright 46158600ac3SJames Wright PetscCall(VecGetType(ctx->X_loc, &vec_type)); 462ed5c6999SJames Wright PetscCall(DefaultMatTypeFromVecType(vec_type, &coo_mat_type)); 46340d80af1SJames Wright PetscCall(PetscStrallocpy(coo_mat_type, &ctx->coo_mat_type)); 46458600ac3SJames Wright } 46558600ac3SJames Wright // -- Set mat operations 46667aa9f91SJames Wright PetscCall(MatShellSetContextDestroy(*mat, (PetscCtxDestroyFn *)MatCeedContextDestroy)); 467e90c2ceeSJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_VIEW, (void (*)(void))MatView_Ceed)); 46858600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 46958600ac3SJames Wright if (op_mult_transpose) PetscCall(MatShellSetOperation(*mat, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed)); 47058600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed)); 47158600ac3SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed)); 472e8ff1987SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_GET_BLOCK_DIAGONAL, (void (*)(void))MatGetBlockDiagonal_Ceed)); 473e8ff1987SJames Wright PetscCall(MatShellSetOperation(*mat, MATOP_GET_VBLOCK_DIAGONAL, (void (*)(void))MatGetVariableBlockDiagonal_Ceed)); 47458600ac3SJames Wright PetscCall(MatShellSetVecType(*mat, vec_type)); 47558600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 47658600ac3SJames Wright } 47758600ac3SJames Wright 47858600ac3SJames Wright /** 47958600ac3SJames Wright @brief Copy `MATCEED` into a compatible `Mat` with type `MatShell` or `MATCEED`. 48058600ac3SJames Wright 48158600ac3SJames Wright Collective across MPI processes. 48258600ac3SJames Wright 48358600ac3SJames Wright @param[in] mat_ceed `MATCEED` to copy from 48458600ac3SJames Wright @param[out] mat_other `MatShell` or `MATCEED` to copy into 48558600ac3SJames Wright 48658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 48758600ac3SJames Wright **/ 48858600ac3SJames Wright PetscErrorCode MatCeedCopy(Mat mat_ceed, Mat mat_other) { 48958600ac3SJames Wright PetscFunctionBeginUser; 49058600ac3SJames Wright PetscCall(MatCeedRegisterLogEvents()); 49158600ac3SJames Wright 49258600ac3SJames Wright // Check type compatibility 49358600ac3SJames Wright { 49440d80af1SJames Wright PetscBool is_matceed = PETSC_FALSE, is_matshell = PETSC_FALSE; 49558600ac3SJames Wright MatType mat_type_ceed, mat_type_other; 49658600ac3SJames Wright 49758600ac3SJames Wright PetscCall(MatGetType(mat_ceed, &mat_type_ceed)); 49840d80af1SJames Wright PetscCall(PetscStrcmp(mat_type_ceed, MATCEED, &is_matceed)); 49940d80af1SJames Wright PetscCheck(is_matceed, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_ceed must have type " MATCEED); 50040d80af1SJames Wright PetscCall(MatGetType(mat_other, &mat_type_other)); 50140d80af1SJames Wright PetscCall(PetscStrcmp(mat_type_other, MATCEED, &is_matceed)); 50240d80af1SJames Wright PetscCall(PetscStrcmp(mat_type_other, MATSHELL, &is_matceed)); 50340d80af1SJames Wright PetscCheck(is_matceed || is_matshell, PETSC_COMM_SELF, PETSC_ERR_LIB, "mat_other must have type " MATCEED " or " MATSHELL); 50458600ac3SJames Wright } 50558600ac3SJames Wright 50658600ac3SJames Wright // Check dimension compatibility 50758600ac3SJames Wright { 50858600ac3SJames 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; 50958600ac3SJames Wright 51058600ac3SJames Wright PetscCall(MatGetSize(mat_ceed, &Y_g_ceed_size, &X_g_ceed_size)); 51158600ac3SJames Wright PetscCall(MatGetLocalSize(mat_ceed, &Y_l_ceed_size, &X_l_ceed_size)); 51258600ac3SJames Wright PetscCall(MatGetSize(mat_ceed, &Y_g_other_size, &X_g_other_size)); 51358600ac3SJames Wright PetscCall(MatGetLocalSize(mat_ceed, &Y_l_other_size, &X_l_other_size)); 51458600ac3SJames 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) && 51558600ac3SJames Wright (X_l_ceed_size == X_l_other_size), 51658600ac3SJames Wright PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, 51758600ac3SJames Wright "mat_ceed and mat_other must have compatible sizes; found mat_ceed (Global: %" PetscInt_FMT ", %" PetscInt_FMT 51858600ac3SJames Wright "; Local: %" PetscInt_FMT ", %" PetscInt_FMT ") mat_other (Global: %" PetscInt_FMT ", %" PetscInt_FMT "; Local: %" PetscInt_FMT 51958600ac3SJames Wright ", %" PetscInt_FMT ")", 52058600ac3SJames 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); 52158600ac3SJames Wright } 52258600ac3SJames Wright 52358600ac3SJames Wright // Convert 52458600ac3SJames Wright { 52558600ac3SJames Wright VecType vec_type; 52658600ac3SJames Wright MatCeedContext ctx; 52758600ac3SJames Wright 52858600ac3SJames Wright PetscCall(PetscObjectChangeTypeName((PetscObject)mat_other, MATCEED)); 52958600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 53058600ac3SJames Wright PetscCall(MatCeedContextReference(ctx)); 53158600ac3SJames Wright PetscCall(MatShellSetContext(mat_other, ctx)); 53267aa9f91SJames Wright PetscCall(MatShellSetContextDestroy(mat_other, (PetscCtxDestroyFn *)MatCeedContextDestroy)); 533e90c2ceeSJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_VIEW, (void (*)(void))MatView_Ceed)); 53458600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 53558600ac3SJames Wright if (ctx->op_mult_transpose) PetscCall(MatShellSetOperation(mat_other, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_Ceed)); 53658600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiagonal_Ceed)); 53758600ac3SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_GET_DIAGONAL_BLOCK, (void (*)(void))MatGetDiagonalBlock_Ceed)); 538e8ff1987SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_GET_BLOCK_DIAGONAL, (void (*)(void))MatGetBlockDiagonal_Ceed)); 539e8ff1987SJames Wright PetscCall(MatShellSetOperation(mat_other, MATOP_GET_VBLOCK_DIAGONAL, (void (*)(void))MatGetVariableBlockDiagonal_Ceed)); 54058600ac3SJames Wright { 54158600ac3SJames Wright PetscInt block_size; 54258600ac3SJames Wright 54358600ac3SJames Wright PetscCall(MatGetBlockSize(mat_ceed, &block_size)); 54458600ac3SJames Wright if (block_size > 1) PetscCall(MatSetBlockSize(mat_other, block_size)); 54558600ac3SJames Wright } 54658600ac3SJames Wright { 54758600ac3SJames Wright PetscInt num_blocks; 54858600ac3SJames Wright const PetscInt *block_sizes; 54958600ac3SJames Wright 55058600ac3SJames Wright PetscCall(MatGetVariableBlockSizes(mat_ceed, &num_blocks, &block_sizes)); 55158600ac3SJames Wright if (num_blocks) PetscCall(MatSetVariableBlockSizes(mat_other, num_blocks, (PetscInt *)block_sizes)); 55258600ac3SJames Wright } 55358600ac3SJames Wright PetscCall(DMGetVecType(ctx->dm_x, &vec_type)); 55458600ac3SJames Wright PetscCall(MatShellSetVecType(mat_other, vec_type)); 55558600ac3SJames Wright } 55658600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 55758600ac3SJames Wright } 55858600ac3SJames Wright 55958600ac3SJames Wright /** 560ba95ebb4SJeremy L Thompson @brief Mark `CeedQFunction` data as updated and the `CeedQFunction` as requiring re-assembly for a `MatCEED`. 561000d2032SJeremy L Thompson 562000d2032SJeremy L Thompson Collective across MPI processes. 563000d2032SJeremy L Thompson 564000d2032SJeremy L Thompson @param[in] mat_ceed `MATCEED` 565000d2032SJeremy L Thompson @param[out] update_needed Boolean flag indicating `CeedQFunction` update needed 566000d2032SJeremy L Thompson 567000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 568000d2032SJeremy L Thompson **/ 569000d2032SJeremy L Thompson PetscErrorCode MatCeedSetAssemblyDataUpdateNeeded(Mat mat_ceed, PetscBool update_needed) { 570000d2032SJeremy L Thompson MatCeedContext ctx; 571000d2032SJeremy L Thompson 572000d2032SJeremy L Thompson PetscFunctionBeginUser; 573000d2032SJeremy L Thompson PetscCall(MatShellGetContext(mat_ceed, &ctx)); 574000d2032SJeremy L Thompson PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult, update_needed)); 575000d2032SJeremy L Thompson if (ctx->op_mult_transpose) { 576000d2032SJeremy L Thompson PetscCallCeed(ctx->ceed, CeedOperatorSetQFunctionAssemblyDataUpdateNeeded(ctx->op_mult_transpose, update_needed)); 577000d2032SJeremy L Thompson } 578397c0187SJeremy L Thompson if (update_needed) { 579397c0187SJeremy L Thompson PetscCall(MatAssemblyBegin(mat_ceed, MAT_FINAL_ASSEMBLY)); 580397c0187SJeremy L Thompson PetscCall(MatAssemblyEnd(mat_ceed, MAT_FINAL_ASSEMBLY)); 581397c0187SJeremy L Thompson } 582000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 583000d2032SJeremy L Thompson } 584000d2032SJeremy L Thompson 585000d2032SJeremy L Thompson /** 58640d80af1SJames Wright @brief Setup a `Mat` with the same COO pattern as a `MatCEED`. 58740d80af1SJames Wright 58840d80af1SJames Wright Collective across MPI processes. 58940d80af1SJames Wright 59040d80af1SJames Wright @param[in] mat_ceed `MATCEED` 59140d80af1SJames Wright @param[out] mat_coo Sparse `Mat` with same COO pattern 59240d80af1SJames Wright 59340d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 59440d80af1SJames Wright **/ 59540d80af1SJames Wright PetscErrorCode MatCeedCreateMatCOO(Mat mat_ceed, Mat *mat_coo) { 59640d80af1SJames Wright MatCeedContext ctx; 59740d80af1SJames Wright 59840d80af1SJames Wright PetscFunctionBeginUser; 59940d80af1SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 60040d80af1SJames Wright 60140d80af1SJames 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"); 60240d80af1SJames Wright 60340d80af1SJames Wright // Check cl mat type 60440d80af1SJames Wright { 60540d80af1SJames Wright PetscBool is_coo_mat_type_cl = PETSC_FALSE; 606ed5c6999SJames Wright char coo_mat_type_cl[PETSC_MAX_OPTION_NAME]; 60740d80af1SJames Wright 60840d80af1SJames Wright // Check for specific CL coo mat type for this Mat 60940d80af1SJames Wright { 61040d80af1SJames Wright const char *mat_ceed_prefix = NULL; 61140d80af1SJames Wright 61240d80af1SJames Wright PetscCall(MatGetOptionsPrefix(mat_ceed, &mat_ceed_prefix)); 61340d80af1SJames Wright PetscOptionsBegin(PetscObjectComm((PetscObject)mat_ceed), mat_ceed_prefix, "", NULL); 61440d80af1SJames Wright PetscCall(PetscOptionsFList("-ceed_coo_mat_type", "Default MATCEED COO assembly MatType", NULL, MatList, coo_mat_type_cl, coo_mat_type_cl, 61540d80af1SJames Wright sizeof(coo_mat_type_cl), &is_coo_mat_type_cl)); 61640d80af1SJames Wright PetscOptionsEnd(); 61740d80af1SJames Wright if (is_coo_mat_type_cl) { 61840d80af1SJames Wright PetscCall(PetscFree(ctx->coo_mat_type)); 61940d80af1SJames Wright PetscCall(PetscStrallocpy(coo_mat_type_cl, &ctx->coo_mat_type)); 62040d80af1SJames Wright } 62140d80af1SJames Wright } 62240d80af1SJames Wright } 62340d80af1SJames Wright 62440d80af1SJames Wright // Create sparse matrix 62540d80af1SJames Wright { 62640d80af1SJames Wright MatType dm_mat_type, dm_mat_type_copy; 62740d80af1SJames Wright 62840d80af1SJames Wright PetscCall(DMGetMatType(ctx->dm_x, &dm_mat_type)); 62940d80af1SJames Wright PetscCall(PetscStrallocpy(dm_mat_type, (char **)&dm_mat_type_copy)); 63040d80af1SJames Wright PetscCall(DMSetMatType(ctx->dm_x, ctx->coo_mat_type)); 63140d80af1SJames Wright PetscCall(DMCreateMatrix(ctx->dm_x, mat_coo)); 63240d80af1SJames Wright PetscCall(DMSetMatType(ctx->dm_x, dm_mat_type_copy)); 63340d80af1SJames Wright PetscCall(PetscFree(dm_mat_type_copy)); 63440d80af1SJames Wright } 63540d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63640d80af1SJames Wright } 63740d80af1SJames Wright 63840d80af1SJames Wright /** 63940d80af1SJames Wright @brief Setup the COO preallocation `MATCEED` into a `MATAIJ` or similar. 64058600ac3SJames Wright The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail. 64158600ac3SJames Wright 64258600ac3SJames Wright Collective across MPI processes. 64358600ac3SJames Wright 64458600ac3SJames Wright @param[in] mat_ceed `MATCEED` to assemble 64558600ac3SJames Wright @param[in,out] mat_coo `MATAIJ` or similar to assemble into 64658600ac3SJames Wright 64758600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 64858600ac3SJames Wright **/ 64940d80af1SJames Wright PetscErrorCode MatCeedSetPreallocationCOO(Mat mat_ceed, Mat mat_coo) { 65058600ac3SJames Wright MatCeedContext ctx; 65158600ac3SJames Wright 65258600ac3SJames Wright PetscFunctionBeginUser; 65358600ac3SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 65458600ac3SJames Wright 65558600ac3SJames Wright { 65658600ac3SJames Wright PetscInt *rows_petsc = NULL, *cols_petsc = NULL; 65758600ac3SJames Wright CeedInt *rows_ceed, *cols_ceed; 65858600ac3SJames Wright PetscCount num_entries; 65958600ac3SJames Wright PetscLogStage stage_amg_setup; 66058600ac3SJames Wright 66158600ac3SJames Wright // -- Assemble sparsity pattern if mat hasn't been assembled before 662c63b910fSJames Wright PetscCall(PetscLogStageGetId("MatCEED Asm Setup", &stage_amg_setup)); 66358600ac3SJames Wright if (stage_amg_setup == -1) { 664c63b910fSJames Wright PetscCall(PetscLogStageRegister("MatCEED Asm Setup", &stage_amg_setup)); 66558600ac3SJames Wright } 66658600ac3SJames Wright PetscCall(PetscLogStagePush(stage_amg_setup)); 667c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL)); 668c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 66950f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleSymbolic(ctx->op_mult, &num_entries, &rows_ceed, &cols_ceed)); 670c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 671a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc)); 672a7dac1d5SJames Wright PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc)); 67358600ac3SJames Wright PetscCall(MatSetPreallocationCOOLocal(mat_coo, num_entries, rows_petsc, cols_petsc)); 67458600ac3SJames Wright free(rows_petsc); 67558600ac3SJames Wright free(cols_petsc); 67650f50432SJames Wright if (!ctx->coo_values_full) PetscCallCeed(ctx->ceed, CeedVectorCreate(ctx->ceed, num_entries, &ctx->coo_values_full)); 67758600ac3SJames Wright PetscCall(PetscRealloc(++ctx->num_mats_assembled_full * sizeof(Mat), &ctx->mats_assembled_full)); 67858600ac3SJames Wright ctx->mats_assembled_full[ctx->num_mats_assembled_full - 1] = mat_coo; 679c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_SETUP_FULL, mat_ceed, mat_coo, NULL, NULL)); 68058600ac3SJames Wright PetscCall(PetscLogStagePop()); 68158600ac3SJames Wright } 68240d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 68340d80af1SJames Wright } 68440d80af1SJames Wright 68540d80af1SJames Wright /** 68640d80af1SJames Wright @brief Assemble a `MATCEED` into a `MATAIJ` or similar. 68740d80af1SJames Wright The `mat_coo` preallocation is set to match the sparsity pattern of `mat_ceed`. 68840d80af1SJames Wright The caller is responsible for assuring the global and local sizes are compatible, otherwise this function will fail. 68940d80af1SJames Wright 69040d80af1SJames Wright Collective across MPI processes. 69140d80af1SJames Wright 69240d80af1SJames Wright @param[in] mat_ceed `MATCEED` to assemble 69340d80af1SJames Wright @param[in,out] mat_coo `MATAIJ` or similar to assemble into 69440d80af1SJames Wright 69540d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 69640d80af1SJames Wright **/ 69740d80af1SJames Wright PetscErrorCode MatCeedAssembleCOO(Mat mat_ceed, Mat mat_coo) { 69840d80af1SJames Wright MatCeedContext ctx; 69940d80af1SJames Wright 70040d80af1SJames Wright PetscFunctionBeginUser; 70140d80af1SJames Wright PetscCall(MatShellGetContext(mat_ceed, &ctx)); 70240d80af1SJames Wright 70340d80af1SJames Wright // Set COO pattern if needed 70440d80af1SJames Wright { 70540d80af1SJames Wright CeedInt index = -1; 70640d80af1SJames Wright 70740d80af1SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) { 70840d80af1SJames Wright if (ctx->mats_assembled_full[i] == mat_coo) index = i; 70940d80af1SJames Wright } 71040d80af1SJames Wright if (index == -1) PetscCall(MatCeedSetPreallocationCOO(mat_ceed, mat_coo)); 71158600ac3SJames Wright } 71258600ac3SJames Wright 71358600ac3SJames Wright // Assemble mat_ceed 714c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL)); 71558600ac3SJames Wright PetscCall(MatAssemblyBegin(mat_coo, MAT_FINAL_ASSEMBLY)); 71658600ac3SJames Wright { 71758600ac3SJames Wright const CeedScalar *values; 71858600ac3SJames Wright CeedMemType mem_type = CEED_MEM_HOST; 71958600ac3SJames Wright PetscBool is_spd, is_spd_known; 720ed5c6999SJames Wright PetscMemType mem_type_petsc; 72158600ac3SJames Wright 722*2da92326SJames Wright PetscCall(MatGetMemTypeFromVecType(mat_coo, &mem_type_petsc)); 723ed5c6999SJames Wright mem_type = MemTypePetscToCeed(mem_type_petsc); 724c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 72550f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssemble(ctx->op_mult, ctx->coo_values_full)); 726c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL_CEEDOP, mat_ceed, mat_coo, NULL, NULL)); 72750f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorGetArrayRead(ctx->coo_values_full, mem_type, &values)); 72858600ac3SJames Wright PetscCall(MatSetValuesCOO(mat_coo, values, INSERT_VALUES)); 72958600ac3SJames Wright PetscCall(MatIsSPDKnown(mat_ceed, &is_spd_known, &is_spd)); 73058600ac3SJames Wright if (is_spd_known) PetscCall(MatSetOption(mat_coo, MAT_SPD, is_spd)); 73150f50432SJames Wright PetscCallCeed(ctx->ceed, CeedVectorRestoreArrayRead(ctx->coo_values_full, &values)); 73258600ac3SJames Wright } 73358600ac3SJames Wright PetscCall(MatAssemblyEnd(mat_coo, MAT_FINAL_ASSEMBLY)); 734c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_FULL, mat_ceed, mat_coo, NULL, NULL)); 73558600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 73658600ac3SJames Wright } 73758600ac3SJames Wright 73858600ac3SJames Wright /** 73940d80af1SJames Wright @brief Set the current value of a context field for a `MatCEED`. 74040d80af1SJames Wright 74140d80af1SJames Wright Not collective across MPI processes. 74240d80af1SJames Wright 74340d80af1SJames Wright @param[in,out] mat `MatCEED` 74440d80af1SJames Wright @param[in] name Name of the context field 74540d80af1SJames Wright @param[in] value New context field value 74640d80af1SJames Wright 74740d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 74840d80af1SJames Wright **/ 74940d80af1SJames Wright PetscErrorCode MatCeedSetContextDouble(Mat mat, const char *name, double value) { 75040d80af1SJames Wright PetscBool was_updated = PETSC_FALSE; 75140d80af1SJames Wright MatCeedContext ctx; 75240d80af1SJames Wright 75340d80af1SJames Wright PetscFunctionBeginUser; 75440d80af1SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 75540d80af1SJames Wright { 75640d80af1SJames Wright CeedContextFieldLabel label = NULL; 75740d80af1SJames Wright 75840d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult, name, &label)); 75940d80af1SJames Wright if (label) { 760000d2032SJeremy L Thompson double set_value = 2 * value + 1.0; 761000d2032SJeremy L Thompson 762000d2032SJeremy L Thompson PetscCall(MatCeedGetContextDouble(mat, name, &set_value)); 763000d2032SJeremy L Thompson if (set_value != value) { 76440d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult, label, &value)); 76540d80af1SJames Wright was_updated = PETSC_TRUE; 76640d80af1SJames Wright } 767000d2032SJeremy L Thompson } 76840d80af1SJames Wright if (ctx->op_mult_transpose) { 76940d80af1SJames Wright label = NULL; 77040d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(ctx->op_mult_transpose, name, &label)); 77140d80af1SJames Wright if (label) { 772000d2032SJeremy L Thompson double set_value = 2 * value + 1.0; 773000d2032SJeremy L Thompson 774000d2032SJeremy L Thompson PetscCall(MatCeedGetContextDouble(mat, name, &set_value)); 775000d2032SJeremy L Thompson if (set_value != value) { 77640d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorSetContextDouble(ctx->op_mult_transpose, label, &value)); 77740d80af1SJames Wright was_updated = PETSC_TRUE; 77840d80af1SJames Wright } 77940d80af1SJames Wright } 78040d80af1SJames Wright } 781000d2032SJeremy L Thompson } 78240d80af1SJames Wright if (was_updated) { 78340d80af1SJames Wright PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY)); 78440d80af1SJames Wright PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY)); 78540d80af1SJames Wright } 78640d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 78740d80af1SJames Wright } 78840d80af1SJames Wright 78940d80af1SJames Wright /** 79040d80af1SJames Wright @brief Get the current value of a context field for a `MatCEED`. 79140d80af1SJames Wright 79240d80af1SJames Wright Not collective across MPI processes. 79340d80af1SJames Wright 79440d80af1SJames Wright @param[in] mat `MatCEED` 79540d80af1SJames Wright @param[in] name Name of the context field 79640d80af1SJames Wright @param[out] value Current context field value 79740d80af1SJames Wright 79840d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 79940d80af1SJames Wright **/ 80040d80af1SJames Wright PetscErrorCode MatCeedGetContextDouble(Mat mat, const char *name, double *value) { 80140d80af1SJames Wright MatCeedContext ctx; 80240d80af1SJames Wright 80340d80af1SJames Wright PetscFunctionBeginUser; 80440d80af1SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 80540d80af1SJames Wright { 80640d80af1SJames Wright CeedContextFieldLabel label = NULL; 80740d80af1SJames Wright CeedOperator op = ctx->op_mult; 80840d80af1SJames Wright 80940d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label)); 81040d80af1SJames Wright if (!label && ctx->op_mult_transpose) { 81140d80af1SJames Wright op = ctx->op_mult_transpose; 81240d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextFieldLabel(op, name, &label)); 81340d80af1SJames Wright } 81440d80af1SJames Wright if (label) { 81540d80af1SJames Wright PetscSizeT num_values; 81640d80af1SJames Wright const double *values_ceed; 81740d80af1SJames Wright 81840d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorGetContextDoubleRead(op, label, &num_values, &values_ceed)); 81940d80af1SJames Wright *value = values_ceed[0]; 82040d80af1SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorRestoreContextDoubleRead(op, label, &values_ceed)); 82140d80af1SJames Wright } 82240d80af1SJames Wright } 82340d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 82440d80af1SJames Wright } 82540d80af1SJames Wright 82640d80af1SJames Wright /** 827000d2032SJeremy L Thompson @brief Set the current `PetscReal` value of a context field for a `MatCEED`. 828000d2032SJeremy L Thompson 829000d2032SJeremy L Thompson Not collective across MPI processes. 830000d2032SJeremy L Thompson 831000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 832000d2032SJeremy L Thompson @param[in] name Name of the context field 833000d2032SJeremy L Thompson @param[in] value New context field value 834000d2032SJeremy L Thompson 835000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 836000d2032SJeremy L Thompson **/ 837000d2032SJeremy L Thompson PetscErrorCode MatCeedSetContextReal(Mat mat, const char *name, PetscReal value) { 838000d2032SJeremy L Thompson double value_double = value; 839000d2032SJeremy L Thompson 840000d2032SJeremy L Thompson PetscFunctionBeginUser; 841000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, name, value_double)); 842000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 843000d2032SJeremy L Thompson } 844000d2032SJeremy L Thompson 845000d2032SJeremy L Thompson /** 84651bb547fSJames Wright @brief Get the current `PetscReal` value of a context field for a `MatCEED`. 84751bb547fSJames Wright 84851bb547fSJames Wright Not collective across MPI processes. 84951bb547fSJames Wright 85051bb547fSJames Wright @param[in] mat `MatCEED` 85151bb547fSJames Wright @param[in] name Name of the context field 85251bb547fSJames Wright @param[out] value Current context field value 85351bb547fSJames Wright 85451bb547fSJames Wright @return An error code: 0 - success, otherwise - failure 85551bb547fSJames Wright **/ 85651bb547fSJames Wright PetscErrorCode MatCeedGetContextReal(Mat mat, const char *name, PetscReal *value) { 85787d3884fSJeremy L Thompson double value_double = 0.0; 85851bb547fSJames Wright 85951bb547fSJames Wright PetscFunctionBeginUser; 86051bb547fSJames Wright PetscCall(MatCeedGetContextDouble(mat, name, &value_double)); 86151bb547fSJames Wright *value = value_double; 86251bb547fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 86351bb547fSJames Wright } 86451bb547fSJames Wright 86551bb547fSJames Wright /** 866000d2032SJeremy L Thompson @brief Set the current time for a `MatCEED`. 867000d2032SJeremy L Thompson 868000d2032SJeremy L Thompson Not collective across MPI processes. 869000d2032SJeremy L Thompson 870000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 871000d2032SJeremy L Thompson @param[in] time Current time 872000d2032SJeremy L Thompson 873000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 874000d2032SJeremy L Thompson **/ 875000d2032SJeremy L Thompson PetscErrorCode MatCeedSetTime(Mat mat, PetscReal time) { 876000d2032SJeremy L Thompson PetscFunctionBeginUser; 877000d2032SJeremy L Thompson { 878000d2032SJeremy L Thompson double time_ceed = time; 879000d2032SJeremy L Thompson 880000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "time", time_ceed)); 881000d2032SJeremy L Thompson } 882000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 883000d2032SJeremy L Thompson } 884000d2032SJeremy L Thompson 885000d2032SJeremy L Thompson /** 886000d2032SJeremy L Thompson @brief Get the current time for a `MatCEED`. 887000d2032SJeremy L Thompson 888000d2032SJeremy L Thompson Not collective across MPI processes. 889000d2032SJeremy L Thompson 890000d2032SJeremy L Thompson @param[in] mat `MatCEED` 891000d2032SJeremy L Thompson @param[out] time Current time, or -1.0 if the boundary evaluator has no time field 892000d2032SJeremy L Thompson 893000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 894000d2032SJeremy L Thompson **/ 895000d2032SJeremy L Thompson PetscErrorCode MatCeedGetTime(Mat mat, PetscReal *time) { 896000d2032SJeremy L Thompson PetscFunctionBeginUser; 897000d2032SJeremy L Thompson *time = -1.0; 898000d2032SJeremy L Thompson { 899000d2032SJeremy L Thompson double time_ceed = -1.0; 900000d2032SJeremy L Thompson 901000d2032SJeremy L Thompson PetscCall(MatCeedGetContextDouble(mat, "time", &time_ceed)); 902000d2032SJeremy L Thompson *time = time_ceed; 903000d2032SJeremy L Thompson } 904000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 905000d2032SJeremy L Thompson } 906000d2032SJeremy L Thompson 907000d2032SJeremy L Thompson /** 908000d2032SJeremy L Thompson @brief Set the current time step for a `MatCEED`. 909000d2032SJeremy L Thompson 910000d2032SJeremy L Thompson Not collective across MPI processes. 911000d2032SJeremy L Thompson 912000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 913000d2032SJeremy L Thompson @param[in] dt Current time step 914000d2032SJeremy L Thompson 915000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 916000d2032SJeremy L Thompson **/ 917000d2032SJeremy L Thompson PetscErrorCode MatCeedSetDt(Mat mat, PetscReal dt) { 918000d2032SJeremy L Thompson PetscFunctionBeginUser; 919000d2032SJeremy L Thompson { 920000d2032SJeremy L Thompson double dt_ceed = dt; 921000d2032SJeremy L Thompson 922000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "dt", dt_ceed)); 923000d2032SJeremy L Thompson } 924000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 925000d2032SJeremy L Thompson } 926000d2032SJeremy L Thompson 927000d2032SJeremy L Thompson /** 928000d2032SJeremy L Thompson @brief Set the Jacobian shifts for a `MatCEED`. 929000d2032SJeremy L Thompson 930000d2032SJeremy L Thompson Not collective across MPI processes. 931000d2032SJeremy L Thompson 932000d2032SJeremy L Thompson @param[in,out] mat `MatCEED` 933000d2032SJeremy L Thompson @param[in] shift_v Velocity shift 934000d2032SJeremy L Thompson @param[in] shift_a Acceleration shift 935000d2032SJeremy L Thompson 936000d2032SJeremy L Thompson @return An error code: 0 - success, otherwise - failure 937000d2032SJeremy L Thompson **/ 938000d2032SJeremy L Thompson PetscErrorCode MatCeedSetShifts(Mat mat, PetscReal shift_v, PetscReal shift_a) { 939000d2032SJeremy L Thompson PetscFunctionBeginUser; 940000d2032SJeremy L Thompson { 941000d2032SJeremy L Thompson double shift_v_ceed = shift_v; 942000d2032SJeremy L Thompson 943000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "shift v", shift_v_ceed)); 944000d2032SJeremy L Thompson } 945000d2032SJeremy L Thompson if (shift_a) { 946000d2032SJeremy L Thompson double shift_a_ceed = shift_a; 947000d2032SJeremy L Thompson 948000d2032SJeremy L Thompson PetscCall(MatCeedSetContextDouble(mat, "shift a", shift_a_ceed)); 949000d2032SJeremy L Thompson } 950000d2032SJeremy L Thompson PetscFunctionReturn(PETSC_SUCCESS); 951000d2032SJeremy L Thompson } 952000d2032SJeremy L Thompson 953000d2032SJeremy L Thompson /** 95458600ac3SJames Wright @brief Set user context for a `MATCEED`. 95558600ac3SJames Wright 95658600ac3SJames Wright Collective across MPI processes. 95758600ac3SJames Wright 95858600ac3SJames Wright @param[in,out] mat `MATCEED` 95958600ac3SJames Wright @param[in] f The context destroy function, or NULL 96058600ac3SJames Wright @param[in] ctx User context, or NULL to unset 96158600ac3SJames Wright 96258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 96358600ac3SJames Wright **/ 96467aa9f91SJames Wright PetscErrorCode MatCeedSetContext(Mat mat, PetscCtxDestroyFn f, void *ctx) { 96558600ac3SJames Wright PetscContainer user_ctx = NULL; 96658600ac3SJames Wright 96758600ac3SJames Wright PetscFunctionBeginUser; 96858600ac3SJames Wright if (ctx) { 96958600ac3SJames Wright PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)mat), &user_ctx)); 97058600ac3SJames Wright PetscCall(PetscContainerSetPointer(user_ctx, ctx)); 97167aa9f91SJames Wright PetscCall(PetscContainerSetCtxDestroy(user_ctx, f)); 97258600ac3SJames Wright } 97358600ac3SJames Wright PetscCall(PetscObjectCompose((PetscObject)mat, "MatCeed user context", (PetscObject)user_ctx)); 97458600ac3SJames Wright PetscCall(PetscContainerDestroy(&user_ctx)); 97558600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 97658600ac3SJames Wright } 97758600ac3SJames Wright 97858600ac3SJames Wright /** 97958600ac3SJames Wright @brief Retrieve the user context for a `MATCEED`. 98058600ac3SJames Wright 98158600ac3SJames Wright Collective across MPI processes. 98258600ac3SJames Wright 98358600ac3SJames Wright @param[in,out] mat `MATCEED` 98458600ac3SJames Wright @param[in] ctx User context 98558600ac3SJames Wright 98658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 98758600ac3SJames Wright **/ 98858600ac3SJames Wright PetscErrorCode MatCeedGetContext(Mat mat, void *ctx) { 98958600ac3SJames Wright PetscContainer user_ctx; 99058600ac3SJames Wright 99158600ac3SJames Wright PetscFunctionBeginUser; 99258600ac3SJames Wright PetscCall(PetscObjectQuery((PetscObject)mat, "MatCeed user context", (PetscObject *)&user_ctx)); 99358600ac3SJames Wright if (user_ctx) PetscCall(PetscContainerGetPointer(user_ctx, (void **)ctx)); 99458600ac3SJames Wright else *(void **)ctx = NULL; 99558600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 99658600ac3SJames Wright } 99758600ac3SJames Wright /** 99840d80af1SJames Wright @brief Set a user defined matrix operation for a `MATCEED` matrix. 99940d80af1SJames Wright 100040d80af1SJames Wright Within each user-defined routine, the user should call `MatCeedGetContext()` to obtain the user-defined context that was set by 100140d80af1SJames Wright `MatCeedSetContext()`. 100258600ac3SJames Wright 100358600ac3SJames Wright Collective across MPI processes. 100458600ac3SJames Wright 100558600ac3SJames Wright @param[in,out] mat `MATCEED` 100640d80af1SJames Wright @param[in] op Name of the `MatOperation` 100740d80af1SJames Wright @param[in] g Function that provides the operation 100858600ac3SJames Wright 100958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 101058600ac3SJames Wright **/ 101140d80af1SJames Wright PetscErrorCode MatCeedSetOperation(Mat mat, MatOperation op, void (*g)(void)) { 101240d80af1SJames Wright PetscFunctionBeginUser; 101340d80af1SJames Wright PetscCall(MatShellSetOperation(mat, op, g)); 101440d80af1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 101540d80af1SJames Wright } 101640d80af1SJames Wright 101740d80af1SJames Wright /** 101840d80af1SJames Wright @brief Sets the default COO matrix type as a string from the `MATCEED`. 101940d80af1SJames Wright 102040d80af1SJames Wright Collective across MPI processes. 102140d80af1SJames Wright 102240d80af1SJames Wright @param[in,out] mat `MATCEED` 102340d80af1SJames Wright @param[in] type COO `MatType` to set 102440d80af1SJames Wright 102540d80af1SJames Wright @return An error code: 0 - success, otherwise - failure 102640d80af1SJames Wright **/ 102740d80af1SJames Wright PetscErrorCode MatCeedSetCOOMatType(Mat mat, MatType type) { 102858600ac3SJames Wright MatCeedContext ctx; 102958600ac3SJames Wright 103058600ac3SJames Wright PetscFunctionBeginUser; 103158600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 103258600ac3SJames Wright // Check if same 103358600ac3SJames Wright { 103458600ac3SJames Wright size_t len_old, len_new; 103558600ac3SJames Wright PetscBool is_same = PETSC_FALSE; 103658600ac3SJames Wright 103740d80af1SJames Wright PetscCall(PetscStrlen(ctx->coo_mat_type, &len_old)); 103858600ac3SJames Wright PetscCall(PetscStrlen(type, &len_new)); 103940d80af1SJames Wright if (len_old == len_new) PetscCall(PetscStrncmp(ctx->coo_mat_type, type, len_old, &is_same)); 104058600ac3SJames Wright if (is_same) PetscFunctionReturn(PETSC_SUCCESS); 104158600ac3SJames Wright } 104258600ac3SJames Wright // Clean up old mats in different format 104358600ac3SJames Wright // LCOV_EXCL_START 104458600ac3SJames Wright if (ctx->mat_assembled_full_internal) { 104558600ac3SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_full; i++) { 104658600ac3SJames Wright if (ctx->mats_assembled_full[i] == ctx->mat_assembled_full_internal) { 104758600ac3SJames Wright for (PetscInt j = i + 1; j < ctx->num_mats_assembled_full; j++) { 104858600ac3SJames Wright ctx->mats_assembled_full[j - 1] = ctx->mats_assembled_full[j]; 104958600ac3SJames Wright } 105058600ac3SJames Wright ctx->num_mats_assembled_full--; 105158600ac3SJames Wright // Note: we'll realloc this array again, so no need to shrink the allocation 105258600ac3SJames Wright PetscCall(MatDestroy(&ctx->mat_assembled_full_internal)); 105358600ac3SJames Wright } 105458600ac3SJames Wright } 105558600ac3SJames Wright } 105658600ac3SJames Wright if (ctx->mat_assembled_pbd_internal) { 105758600ac3SJames Wright for (PetscInt i = 0; i < ctx->num_mats_assembled_pbd; i++) { 105858600ac3SJames Wright if (ctx->mats_assembled_pbd[i] == ctx->mat_assembled_pbd_internal) { 105958600ac3SJames Wright for (PetscInt j = i + 1; j < ctx->num_mats_assembled_pbd; j++) { 106058600ac3SJames Wright ctx->mats_assembled_pbd[j - 1] = ctx->mats_assembled_pbd[j]; 106158600ac3SJames Wright } 106258600ac3SJames Wright // Note: we'll realloc this array again, so no need to shrink the allocation 106358600ac3SJames Wright ctx->num_mats_assembled_pbd--; 106458600ac3SJames Wright PetscCall(MatDestroy(&ctx->mat_assembled_pbd_internal)); 106558600ac3SJames Wright } 106658600ac3SJames Wright } 106758600ac3SJames Wright } 106840d80af1SJames Wright PetscCall(PetscFree(ctx->coo_mat_type)); 106940d80af1SJames Wright PetscCall(PetscStrallocpy(type, &ctx->coo_mat_type)); 107058600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 107158600ac3SJames Wright // LCOV_EXCL_STOP 107258600ac3SJames Wright } 107358600ac3SJames Wright 107458600ac3SJames Wright /** 107540d80af1SJames Wright @brief Gets the default COO matrix type as a string from the `MATCEED`. 107658600ac3SJames Wright 107758600ac3SJames Wright Collective across MPI processes. 107858600ac3SJames Wright 107958600ac3SJames Wright @param[in,out] mat `MATCEED` 108040d80af1SJames Wright @param[in] type COO `MatType` 108158600ac3SJames Wright 108258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 108358600ac3SJames Wright **/ 108440d80af1SJames Wright PetscErrorCode MatCeedGetCOOMatType(Mat mat, MatType *type) { 108558600ac3SJames Wright MatCeedContext ctx; 108658600ac3SJames Wright 108758600ac3SJames Wright PetscFunctionBeginUser; 108858600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 108940d80af1SJames Wright *type = ctx->coo_mat_type; 109058600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 109158600ac3SJames Wright } 109258600ac3SJames Wright 109358600ac3SJames Wright /** 109458600ac3SJames Wright @brief Set input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 109558600ac3SJames Wright 109658600ac3SJames Wright Not collective across MPI processes. 109758600ac3SJames Wright 109858600ac3SJames Wright @param[in,out] mat `MATCEED` 109958600ac3SJames Wright @param[in] X_loc Input PETSc local vector, or NULL 110058600ac3SJames Wright @param[in] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 110158600ac3SJames Wright 110258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 110358600ac3SJames Wright **/ 110458600ac3SJames Wright PetscErrorCode MatCeedSetLocalVectors(Mat mat, Vec X_loc, Vec Y_loc_transpose) { 110558600ac3SJames Wright MatCeedContext ctx; 110658600ac3SJames Wright 110758600ac3SJames Wright PetscFunctionBeginUser; 110858600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 110958600ac3SJames Wright if (X_loc) { 111058600ac3SJames Wright PetscInt len_old, len_new; 111158600ac3SJames Wright 111258600ac3SJames Wright PetscCall(VecGetSize(ctx->X_loc, &len_old)); 111358600ac3SJames Wright PetscCall(VecGetSize(X_loc, &len_new)); 111458600ac3SJames 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, 111558600ac3SJames Wright len_new, len_old); 111640d80af1SJames Wright PetscCall(VecReferenceCopy(X_loc, &ctx->X_loc)); 111758600ac3SJames Wright } 111858600ac3SJames Wright if (Y_loc_transpose) { 111958600ac3SJames Wright PetscInt len_old, len_new; 112058600ac3SJames Wright 112158600ac3SJames Wright PetscCall(VecGetSize(ctx->Y_loc_transpose, &len_old)); 112258600ac3SJames Wright PetscCall(VecGetSize(Y_loc_transpose, &len_new)); 112358600ac3SJames Wright PetscCheck(len_old == len_new, PETSC_COMM_SELF, PETSC_ERR_LIB, 112458600ac3SJames Wright "new Y_loc_transpose length %" PetscInt_FMT " should match old Y_loc_transpose length %" PetscInt_FMT, len_new, len_old); 112540d80af1SJames Wright PetscCall(VecReferenceCopy(Y_loc_transpose, &ctx->Y_loc_transpose)); 112658600ac3SJames Wright } 112758600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 112858600ac3SJames Wright } 112958600ac3SJames Wright 113058600ac3SJames Wright /** 113158600ac3SJames Wright @brief Get input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 113258600ac3SJames Wright 113358600ac3SJames Wright Not collective across MPI processes. 113458600ac3SJames Wright 113558600ac3SJames Wright @param[in,out] mat `MATCEED` 113658600ac3SJames Wright @param[out] X_loc Input PETSc local vector, or NULL 113758600ac3SJames Wright @param[out] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 113858600ac3SJames Wright 113958600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 114058600ac3SJames Wright **/ 114158600ac3SJames Wright PetscErrorCode MatCeedGetLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) { 114258600ac3SJames Wright MatCeedContext ctx; 114358600ac3SJames Wright 114458600ac3SJames Wright PetscFunctionBeginUser; 114558600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 114658600ac3SJames Wright if (X_loc) { 114740d80af1SJames Wright *X_loc = NULL; 114840d80af1SJames Wright PetscCall(VecReferenceCopy(ctx->X_loc, X_loc)); 114958600ac3SJames Wright } 115058600ac3SJames Wright if (Y_loc_transpose) { 115140d80af1SJames Wright *Y_loc_transpose = NULL; 115240d80af1SJames Wright PetscCall(VecReferenceCopy(ctx->Y_loc_transpose, Y_loc_transpose)); 115358600ac3SJames Wright } 115458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 115558600ac3SJames Wright } 115658600ac3SJames Wright 115758600ac3SJames Wright /** 115858600ac3SJames Wright @brief Restore input local vectors for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 115958600ac3SJames Wright 116058600ac3SJames Wright Not collective across MPI processes. 116158600ac3SJames Wright 116258600ac3SJames Wright @param[in,out] mat MatCeed 116358600ac3SJames Wright @param[out] X_loc Input PETSc local vector, or NULL 116458600ac3SJames Wright @param[out] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 116558600ac3SJames Wright 116658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 116758600ac3SJames Wright **/ 116858600ac3SJames Wright PetscErrorCode MatCeedRestoreLocalVectors(Mat mat, Vec *X_loc, Vec *Y_loc_transpose) { 116958600ac3SJames Wright PetscFunctionBeginUser; 117058600ac3SJames Wright if (X_loc) PetscCall(VecDestroy(X_loc)); 117158600ac3SJames Wright if (Y_loc_transpose) PetscCall(VecDestroy(Y_loc_transpose)); 117258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 117358600ac3SJames Wright } 117458600ac3SJames Wright 117558600ac3SJames Wright /** 117658600ac3SJames Wright @brief Get libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 117758600ac3SJames Wright 117858600ac3SJames Wright Not collective across MPI processes. 117958600ac3SJames Wright 118058600ac3SJames Wright @param[in,out] mat MatCeed 118158600ac3SJames Wright @param[out] op_mult libCEED `CeedOperator` for `MatMult()`, or NULL 118258600ac3SJames Wright @param[out] op_mult_transpose libCEED `CeedOperator` for `MatMultTranspose()`, or NULL 118358600ac3SJames Wright 118458600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 118558600ac3SJames Wright **/ 118658600ac3SJames Wright PetscErrorCode MatCeedGetCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) { 118758600ac3SJames Wright MatCeedContext ctx; 118858600ac3SJames Wright 118958600ac3SJames Wright PetscFunctionBeginUser; 119058600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 119158600ac3SJames Wright if (op_mult) { 119258600ac3SJames Wright *op_mult = NULL; 119350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult, op_mult)); 119458600ac3SJames Wright } 119558600ac3SJames Wright if (op_mult_transpose) { 119658600ac3SJames Wright *op_mult_transpose = NULL; 119750f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorReferenceCopy(ctx->op_mult_transpose, op_mult_transpose)); 119858600ac3SJames Wright } 119958600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 120058600ac3SJames Wright } 120158600ac3SJames Wright 120258600ac3SJames Wright /** 120358600ac3SJames Wright @brief Restore libCEED `CeedOperator` for `MATCEED` `MatMult()` and `MatMultTranspose()` operations. 120458600ac3SJames Wright 120558600ac3SJames Wright Not collective across MPI processes. 120658600ac3SJames Wright 120758600ac3SJames Wright @param[in,out] mat MatCeed 120858600ac3SJames Wright @param[out] op_mult libCEED `CeedOperator` for `MatMult()`, or NULL 120958600ac3SJames Wright @param[out] op_mult_transpose libCEED `CeedOperator` for `MatMultTranspose()`, or NULL 121058600ac3SJames Wright 121158600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 121258600ac3SJames Wright **/ 121358600ac3SJames Wright PetscErrorCode MatCeedRestoreCeedOperators(Mat mat, CeedOperator *op_mult, CeedOperator *op_mult_transpose) { 121458600ac3SJames Wright MatCeedContext ctx; 121558600ac3SJames Wright 121658600ac3SJames Wright PetscFunctionBeginUser; 121758600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 121850f50432SJames Wright if (op_mult) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult)); 121950f50432SJames Wright if (op_mult_transpose) PetscCallCeed(ctx->ceed, CeedOperatorDestroy(op_mult_transpose)); 122058600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 122158600ac3SJames Wright } 122258600ac3SJames Wright 122358600ac3SJames Wright /** 122458600ac3SJames Wright @brief Set `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 122558600ac3SJames Wright 122658600ac3SJames Wright Not collective across MPI processes. 122758600ac3SJames Wright 122858600ac3SJames Wright @param[in,out] mat MatCeed 122958600ac3SJames Wright @param[out] log_event_mult `PetscLogEvent` for forward evaluation, or NULL 123058600ac3SJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose evaluation, or NULL 123158600ac3SJames Wright 123258600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 123358600ac3SJames Wright **/ 123458600ac3SJames Wright PetscErrorCode MatCeedSetLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) { 123558600ac3SJames Wright MatCeedContext ctx; 123658600ac3SJames Wright 123758600ac3SJames Wright PetscFunctionBeginUser; 123858600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 123958600ac3SJames Wright if (log_event_mult) ctx->log_event_mult = log_event_mult; 124058600ac3SJames Wright if (log_event_mult_transpose) ctx->log_event_mult_transpose = log_event_mult_transpose; 124158600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 124258600ac3SJames Wright } 124358600ac3SJames Wright 124458600ac3SJames Wright /** 124558600ac3SJames Wright @brief Get `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 124658600ac3SJames Wright 124758600ac3SJames Wright Not collective across MPI processes. 124858600ac3SJames Wright 124958600ac3SJames Wright @param[in,out] mat MatCeed 125058600ac3SJames Wright @param[out] log_event_mult `PetscLogEvent` for forward evaluation, or NULL 125158600ac3SJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose evaluation, or NULL 125258600ac3SJames Wright 125358600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 125458600ac3SJames Wright **/ 125558600ac3SJames Wright PetscErrorCode MatCeedGetLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) { 125658600ac3SJames Wright MatCeedContext ctx; 125758600ac3SJames Wright 125858600ac3SJames Wright PetscFunctionBeginUser; 125958600ac3SJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 126058600ac3SJames Wright if (log_event_mult) *log_event_mult = ctx->log_event_mult; 126158600ac3SJames Wright if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_mult_transpose; 126258600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 126358600ac3SJames Wright } 126458600ac3SJames Wright 1265c63b910fSJames Wright /** 1266c63b910fSJames Wright @brief Set `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 1267c63b910fSJames Wright 1268c63b910fSJames Wright Not collective across MPI processes. 1269c63b910fSJames Wright 1270c63b910fSJames Wright @param[in,out] mat MatCeed 1271c63b910fSJames Wright @param[out] log_event_mult `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL 1272c63b910fSJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL 1273c63b910fSJames Wright 1274c63b910fSJames Wright @return An error code: 0 - success, otherwise - failure 1275c63b910fSJames Wright **/ 1276c63b910fSJames Wright PetscErrorCode MatCeedSetCeedOperatorLogEvents(Mat mat, PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose) { 1277c63b910fSJames Wright MatCeedContext ctx; 1278c63b910fSJames Wright 1279c63b910fSJames Wright PetscFunctionBeginUser; 1280c63b910fSJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 1281c63b910fSJames Wright if (log_event_mult) ctx->log_event_ceed_mult = log_event_mult; 1282c63b910fSJames Wright if (log_event_mult_transpose) ctx->log_event_ceed_mult_transpose = log_event_mult_transpose; 1283c63b910fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1284c63b910fSJames Wright } 1285c63b910fSJames Wright 1286c63b910fSJames Wright /** 1287c63b910fSJames Wright @brief Get `CeedOperator` `PetscLogEvent` for `MATCEED` `MatMult()` and `MatMultTranspose()` operators. 1288c63b910fSJames Wright 1289c63b910fSJames Wright Not collective across MPI processes. 1290c63b910fSJames Wright 1291c63b910fSJames Wright @param[in,out] mat MatCeed 1292c63b910fSJames Wright @param[out] log_event_mult `PetscLogEvent` for forward `CeedOperator` evaluation, or NULL 1293c63b910fSJames Wright @param[out] log_event_mult_transpose `PetscLogEvent` for transpose `CeedOperator` evaluation, or NULL 1294c63b910fSJames Wright 1295c63b910fSJames Wright @return An error code: 0 - success, otherwise - failure 1296c63b910fSJames Wright **/ 1297c63b910fSJames Wright PetscErrorCode MatCeedGetCeedOperatorLogEvents(Mat mat, PetscLogEvent *log_event_mult, PetscLogEvent *log_event_mult_transpose) { 1298c63b910fSJames Wright MatCeedContext ctx; 1299c63b910fSJames Wright 1300c63b910fSJames Wright PetscFunctionBeginUser; 1301c63b910fSJames Wright PetscCall(MatShellGetContext(mat, &ctx)); 1302c63b910fSJames Wright if (log_event_mult) *log_event_mult = ctx->log_event_ceed_mult; 1303c63b910fSJames Wright if (log_event_mult_transpose) *log_event_mult_transpose = ctx->log_event_ceed_mult_transpose; 1304c63b910fSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1305c63b910fSJames Wright } 1306c63b910fSJames Wright 130758600ac3SJames Wright // ----------------------------------------------------------------------------- 130858600ac3SJames Wright // Operator context data 130958600ac3SJames Wright // ----------------------------------------------------------------------------- 131058600ac3SJames Wright 131158600ac3SJames Wright /** 131258600ac3SJames Wright @brief Setup context data for operator application. 131358600ac3SJames Wright 131458600ac3SJames Wright Collective across MPI processes. 131558600ac3SJames Wright 131658600ac3SJames Wright @param[in] dm_x Input `DM` 131758600ac3SJames Wright @param[in] dm_y Output `DM` 131858600ac3SJames Wright @param[in] X_loc Input PETSc local vector, or NULL 131958600ac3SJames Wright @param[in] Y_loc_transpose Input PETSc local vector for transpose operation, or NULL 132058600ac3SJames Wright @param[in] op_mult `CeedOperator` for forward evaluation 132158600ac3SJames Wright @param[in] op_mult_transpose `CeedOperator` for transpose evaluation 132258600ac3SJames Wright @param[in] log_event_mult `PetscLogEvent` for forward evaluation 132358600ac3SJames Wright @param[in] log_event_mult_transpose `PetscLogEvent` for transpose evaluation 1324c63b910fSJames Wright @param[in] log_event_ceed_mult `PetscLogEvent` for forward `CeedOperator` evaluation 1325c63b910fSJames Wright @param[in] log_event_ceed_mult_transpose `PetscLogEvent` for transpose `CeedOperator` evaluation 132658600ac3SJames Wright @param[out] ctx Context data for operator evaluation 132758600ac3SJames Wright 132858600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 132958600ac3SJames Wright **/ 133058600ac3SJames Wright PetscErrorCode MatCeedContextCreate(DM dm_x, DM dm_y, Vec X_loc, Vec Y_loc_transpose, CeedOperator op_mult, CeedOperator op_mult_transpose, 1331c63b910fSJames Wright PetscLogEvent log_event_mult, PetscLogEvent log_event_mult_transpose, PetscLogEvent log_event_ceed_mult, 1332c63b910fSJames Wright PetscLogEvent log_event_ceed_mult_transpose, MatCeedContext *ctx) { 133358600ac3SJames Wright CeedSize x_loc_len, y_loc_len; 133458600ac3SJames Wright 133558600ac3SJames Wright PetscFunctionBeginUser; 133658600ac3SJames Wright 133758600ac3SJames Wright // Allocate 133858600ac3SJames Wright PetscCall(PetscNew(ctx)); 133958600ac3SJames Wright (*ctx)->ref_count = 1; 134058600ac3SJames Wright 134158600ac3SJames Wright // Logging 134258600ac3SJames Wright (*ctx)->log_event_mult = log_event_mult; 134358600ac3SJames Wright (*ctx)->log_event_mult_transpose = log_event_mult_transpose; 1344c63b910fSJames Wright (*ctx)->log_event_ceed_mult = log_event_ceed_mult; 1345c63b910fSJames Wright (*ctx)->log_event_ceed_mult_transpose = log_event_ceed_mult_transpose; 134658600ac3SJames Wright 134758600ac3SJames Wright // PETSc objects 134840d80af1SJames Wright PetscCall(DMReferenceCopy(dm_x, &(*ctx)->dm_x)); 134940d80af1SJames Wright PetscCall(DMReferenceCopy(dm_y, &(*ctx)->dm_y)); 135040d80af1SJames Wright if (X_loc) PetscCall(VecReferenceCopy(X_loc, &(*ctx)->X_loc)); 135140d80af1SJames Wright if (Y_loc_transpose) PetscCall(VecReferenceCopy(Y_loc_transpose, &(*ctx)->Y_loc_transpose)); 135258600ac3SJames Wright 135358600ac3SJames Wright // Memtype 135458600ac3SJames Wright { 135558600ac3SJames Wright const PetscScalar *x; 135658600ac3SJames Wright Vec X; 135758600ac3SJames Wright 135858600ac3SJames Wright PetscCall(DMGetLocalVector(dm_x, &X)); 135958600ac3SJames Wright PetscCall(VecGetArrayReadAndMemType(X, &x, &(*ctx)->mem_type)); 136058600ac3SJames Wright PetscCall(VecRestoreArrayReadAndMemType(X, &x)); 136158600ac3SJames Wright PetscCall(DMRestoreLocalVector(dm_x, &X)); 136258600ac3SJames Wright } 136358600ac3SJames Wright 136458600ac3SJames Wright // libCEED objects 136558600ac3SJames Wright PetscCheck(CeedOperatorGetCeed(op_mult, &(*ctx)->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, 136658600ac3SJames Wright "retrieving Ceed context object failed"); 136750f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorGetActiveVectorLengths(op_mult, &x_loc_len, &y_loc_len)); 136850f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult, &(*ctx)->op_mult)); 136950f50432SJames Wright if (op_mult_transpose) PetscCallCeed((*ctx)->ceed, CeedOperatorReferenceCopy(op_mult_transpose, &(*ctx)->op_mult_transpose)); 137050f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, x_loc_len, &(*ctx)->x_loc)); 137150f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorCreate((*ctx)->ceed, y_loc_len, &(*ctx)->y_loc)); 137258600ac3SJames Wright 137358600ac3SJames Wright // Flop counting 137458600ac3SJames Wright { 137558600ac3SJames Wright CeedSize ceed_flops_estimate = 0; 137658600ac3SJames Wright 137750f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult, &ceed_flops_estimate)); 137858600ac3SJames Wright (*ctx)->flops_mult = ceed_flops_estimate; 137958600ac3SJames Wright if (op_mult_transpose) { 138050f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorGetFlopsEstimate(op_mult_transpose, &ceed_flops_estimate)); 138158600ac3SJames Wright (*ctx)->flops_mult_transpose = ceed_flops_estimate; 138258600ac3SJames Wright } 138358600ac3SJames Wright } 138458600ac3SJames Wright 138558600ac3SJames Wright // Check sizes 138658600ac3SJames Wright if (x_loc_len > 0 || y_loc_len > 0) { 138758600ac3SJames Wright CeedSize ctx_x_loc_len, ctx_y_loc_len; 138858600ac3SJames Wright PetscInt X_loc_len, dm_x_loc_len, Y_loc_len, dm_y_loc_len; 138958600ac3SJames Wright Vec dm_X_loc, dm_Y_loc; 139058600ac3SJames Wright 139158600ac3SJames Wright // -- Input 139258600ac3SJames Wright PetscCall(DMGetLocalVector(dm_x, &dm_X_loc)); 139358600ac3SJames Wright PetscCall(VecGetLocalSize(dm_X_loc, &dm_x_loc_len)); 139458600ac3SJames Wright PetscCall(DMRestoreLocalVector(dm_x, &dm_X_loc)); 139550f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->x_loc, &ctx_x_loc_len)); 13964c17272bSJames Wright if (X_loc) { 13974c17272bSJames Wright PetscCall(VecGetLocalSize(X_loc, &X_loc_len)); 13984c17272bSJames Wright PetscCheck(X_loc_len == dm_x_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, 13994c17272bSJames Wright "X_loc (%" PetscInt_FMT ") must match dm_x (%" PetscInt_FMT ") dimensions", X_loc_len, dm_x_loc_len); 14004c17272bSJames Wright } 14014c17272bSJames 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", 14024c17272bSJames Wright x_loc_len, dm_x_loc_len); 14034c17272bSJames 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 ")", 14044c17272bSJames Wright x_loc_len, ctx_x_loc_len); 140558600ac3SJames Wright 140658600ac3SJames Wright // -- Output 140758600ac3SJames Wright PetscCall(DMGetLocalVector(dm_y, &dm_Y_loc)); 140858600ac3SJames Wright PetscCall(VecGetLocalSize(dm_Y_loc, &dm_y_loc_len)); 140958600ac3SJames Wright PetscCall(DMRestoreLocalVector(dm_y, &dm_Y_loc)); 141050f50432SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorGetLength((*ctx)->y_loc, &ctx_y_loc_len)); 14114c17272bSJames 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", 14124c17272bSJames Wright ctx_y_loc_len, dm_y_loc_len); 141358600ac3SJames Wright 141458600ac3SJames Wright // -- Transpose 141558600ac3SJames Wright if (Y_loc_transpose) { 141658600ac3SJames Wright PetscCall(VecGetLocalSize(Y_loc_transpose, &Y_loc_len)); 14174c17272bSJames Wright PetscCheck(Y_loc_len == dm_y_loc_len, PETSC_COMM_SELF, PETSC_ERR_LIB, 14184c17272bSJames Wright "Y_loc_transpose (%" PetscInt_FMT ") must match dm_y (%" PetscInt_FMT ") dimensions", Y_loc_len, dm_y_loc_len); 141958600ac3SJames Wright } 142058600ac3SJames Wright } 142158600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 142258600ac3SJames Wright } 142358600ac3SJames Wright 142458600ac3SJames Wright /** 142558600ac3SJames Wright @brief Increment reference counter for `MATCEED` context. 142658600ac3SJames Wright 142758600ac3SJames Wright Not collective across MPI processes. 142858600ac3SJames Wright 142958600ac3SJames Wright @param[in,out] ctx Context data 143058600ac3SJames Wright 143158600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 143258600ac3SJames Wright **/ 143358600ac3SJames Wright PetscErrorCode MatCeedContextReference(MatCeedContext ctx) { 143458600ac3SJames Wright PetscFunctionBeginUser; 143558600ac3SJames Wright ctx->ref_count++; 143658600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 143758600ac3SJames Wright } 143858600ac3SJames Wright 143958600ac3SJames Wright /** 144058600ac3SJames Wright @brief Copy reference for `MATCEED`. 144158600ac3SJames Wright Note: If `ctx_copy` is non-null, it is assumed to be a valid pointer to a `MatCeedContext`. 144258600ac3SJames Wright 144358600ac3SJames Wright Not collective across MPI processes. 144458600ac3SJames Wright 144558600ac3SJames Wright @param[in] ctx Context data 144658600ac3SJames Wright @param[out] ctx_copy Copy of pointer to context data 144758600ac3SJames Wright 144858600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 144958600ac3SJames Wright **/ 145058600ac3SJames Wright PetscErrorCode MatCeedContextReferenceCopy(MatCeedContext ctx, MatCeedContext *ctx_copy) { 145158600ac3SJames Wright PetscFunctionBeginUser; 145258600ac3SJames Wright PetscCall(MatCeedContextReference(ctx)); 145367aa9f91SJames Wright PetscCall(MatCeedContextDestroy(ctx_copy)); 145458600ac3SJames Wright *ctx_copy = ctx; 145558600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 145658600ac3SJames Wright } 145758600ac3SJames Wright 145858600ac3SJames Wright /** 145958600ac3SJames Wright @brief Destroy context data for operator application. 146058600ac3SJames Wright 146158600ac3SJames Wright Collective across MPI processes. 146258600ac3SJames Wright 146358600ac3SJames Wright @param[in,out] ctx Context data for operator evaluation 146458600ac3SJames Wright 146558600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 146658600ac3SJames Wright **/ 146767aa9f91SJames Wright PetscErrorCode MatCeedContextDestroy(MatCeedContext *ctx) { 146858600ac3SJames Wright PetscFunctionBeginUser; 146967aa9f91SJames Wright if (!ctx || --(*ctx)->ref_count > 0) PetscFunctionReturn(PETSC_SUCCESS); 147058600ac3SJames Wright 147158600ac3SJames Wright // PETSc objects 147267aa9f91SJames Wright PetscCall(DMDestroy(&(*ctx)->dm_x)); 147367aa9f91SJames Wright PetscCall(DMDestroy(&(*ctx)->dm_y)); 147467aa9f91SJames Wright PetscCall(VecDestroy(&(*ctx)->X_loc)); 147567aa9f91SJames Wright PetscCall(VecDestroy(&(*ctx)->Y_loc_transpose)); 147667aa9f91SJames Wright PetscCall(MatDestroy(&(*ctx)->mat_assembled_full_internal)); 147767aa9f91SJames Wright PetscCall(MatDestroy(&(*ctx)->mat_assembled_pbd_internal)); 147867aa9f91SJames Wright PetscCall(PetscFree((*ctx)->coo_mat_type)); 147967aa9f91SJames Wright PetscCall(PetscFree((*ctx)->mats_assembled_full)); 148067aa9f91SJames Wright PetscCall(PetscFree((*ctx)->mats_assembled_pbd)); 148158600ac3SJames Wright 148258600ac3SJames Wright // libCEED objects 148367aa9f91SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->x_loc)); 148467aa9f91SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->y_loc)); 148567aa9f91SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->coo_values_full)); 148667aa9f91SJames Wright PetscCallCeed((*ctx)->ceed, CeedVectorDestroy(&(*ctx)->coo_values_pbd)); 148767aa9f91SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorDestroy(&(*ctx)->op_mult)); 148867aa9f91SJames Wright PetscCallCeed((*ctx)->ceed, CeedOperatorDestroy(&(*ctx)->op_mult_transpose)); 148967aa9f91SJames Wright PetscCheck(CeedDestroy(&(*ctx)->ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "destroying libCEED context object failed"); 149058600ac3SJames Wright 149158600ac3SJames Wright // Deallocate 149267aa9f91SJames Wright (*ctx)->is_destroyed = PETSC_TRUE; // Flag as destroyed in case someone has stale ref 149367aa9f91SJames Wright PetscCall(PetscFree(*ctx)); 149458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 149558600ac3SJames Wright } 149658600ac3SJames Wright 149758600ac3SJames Wright /** 149858600ac3SJames Wright @brief Compute the diagonal of an operator via libCEED. 149958600ac3SJames Wright 150058600ac3SJames Wright Collective across MPI processes. 150158600ac3SJames Wright 150258600ac3SJames Wright @param[in] A `MATCEED` 150358600ac3SJames Wright @param[out] D Vector holding operator diagonal 150458600ac3SJames Wright 150558600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 150658600ac3SJames Wright **/ 150758600ac3SJames Wright PetscErrorCode MatGetDiagonal_Ceed(Mat A, Vec D) { 150858600ac3SJames Wright PetscMemType mem_type; 150958600ac3SJames Wright Vec D_loc; 151058600ac3SJames Wright MatCeedContext ctx; 151158600ac3SJames Wright 151258600ac3SJames Wright PetscFunctionBeginUser; 151358600ac3SJames Wright PetscCall(MatShellGetContext(A, &ctx)); 151458600ac3SJames Wright 151558600ac3SJames Wright // Place PETSc vector in libCEED vector 1516c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL)); 151758600ac3SJames Wright PetscCall(DMGetLocalVector(ctx->dm_x, &D_loc)); 1518a7dac1d5SJames Wright PetscCall(VecPetscToCeed(D_loc, &mem_type, ctx->x_loc)); 151958600ac3SJames Wright 152058600ac3SJames Wright // Compute Diagonal 1521c63b910fSJames Wright PetscCall(PetscLogEventBegin(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL)); 152250f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorLinearAssembleDiagonal(ctx->op_mult, ctx->x_loc, CEED_REQUEST_IMMEDIATE)); 1523c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL_CEEDOP, A, D, NULL, NULL)); 152458600ac3SJames Wright 152558600ac3SJames Wright // Restore PETSc vector 1526a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(ctx->x_loc, mem_type, D_loc)); 152758600ac3SJames Wright 152858600ac3SJames Wright // Local-to-Global 152958600ac3SJames Wright PetscCall(VecZeroEntries(D)); 153058600ac3SJames Wright PetscCall(DMLocalToGlobal(ctx->dm_x, D_loc, ADD_VALUES, D)); 153158600ac3SJames Wright PetscCall(DMRestoreLocalVector(ctx->dm_x, &D_loc)); 1532c63b910fSJames Wright PetscCall(PetscLogEventEnd(MATCEED_ASSEMBLE_DIAGONAL, A, D, NULL, NULL)); 153358600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 153458600ac3SJames Wright } 153558600ac3SJames Wright 153658600ac3SJames Wright /** 153758600ac3SJames Wright @brief Compute `A X = Y` for a `MATCEED`. 153858600ac3SJames Wright 153958600ac3SJames Wright Collective across MPI processes. 154058600ac3SJames Wright 154158600ac3SJames Wright @param[in] A `MATCEED` 154258600ac3SJames Wright @param[in] X Input PETSc vector 154358600ac3SJames Wright @param[out] Y Output PETSc vector 154458600ac3SJames Wright 154558600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 154658600ac3SJames Wright **/ 154758600ac3SJames Wright PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) { 154858600ac3SJames Wright MatCeedContext ctx; 154958600ac3SJames Wright 155058600ac3SJames Wright PetscFunctionBeginUser; 155158600ac3SJames Wright PetscCall(MatShellGetContext(A, &ctx)); 1552c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_mult, A, X, Y, NULL)); 155358600ac3SJames Wright 155458600ac3SJames Wright { 155558600ac3SJames Wright PetscMemType x_mem_type, y_mem_type; 155658600ac3SJames Wright Vec X_loc = ctx->X_loc, Y_loc; 155758600ac3SJames Wright 155858600ac3SJames Wright // Get local vectors 155958600ac3SJames Wright if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); 156058600ac3SJames Wright PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); 156158600ac3SJames Wright 156258600ac3SJames Wright // Global-to-local 156358600ac3SJames Wright PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc)); 156458600ac3SJames Wright 156558600ac3SJames Wright // Setup libCEED vectors 1566a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(X_loc, &x_mem_type, ctx->x_loc)); 156758600ac3SJames Wright PetscCall(VecZeroEntries(Y_loc)); 1568a7dac1d5SJames Wright PetscCall(VecPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc)); 156958600ac3SJames Wright 157058600ac3SJames Wright // Apply libCEED operator 1571c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult, A, X, Y, NULL)); 1572537ec908SJames Wright PetscCall(PetscLogGpuTimeBegin()); 157350f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult, ctx->x_loc, ctx->y_loc, CEED_REQUEST_IMMEDIATE)); 1574ed5c6999SJames Wright // Log flops 1575ed5c6999SJames Wright if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult)); 1576ed5c6999SJames Wright else PetscCall(PetscLogFlops(ctx->flops_mult)); 157758600ac3SJames Wright PetscCall(PetscLogGpuTimeEnd()); 1578537ec908SJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult, A, X, Y, NULL)); 157958600ac3SJames Wright 158058600ac3SJames Wright // Restore PETSc vectors 1581a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(ctx->x_loc, x_mem_type, X_loc)); 1582a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc)); 158358600ac3SJames Wright 158458600ac3SJames Wright // Local-to-global 158558600ac3SJames Wright PetscCall(VecZeroEntries(Y)); 158658600ac3SJames Wright PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y)); 158758600ac3SJames Wright 158858600ac3SJames Wright // Restore local vectors, as needed 158958600ac3SJames Wright if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); 159058600ac3SJames Wright PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); 159158600ac3SJames Wright } 159258600ac3SJames Wright 1593c63b910fSJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_mult, A, X, Y, NULL)); 159458600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 159558600ac3SJames Wright } 159658600ac3SJames Wright 159758600ac3SJames Wright /** 159858600ac3SJames Wright @brief Compute `A^T Y = X` for a `MATCEED`. 159958600ac3SJames Wright 160058600ac3SJames Wright Collective across MPI processes. 160158600ac3SJames Wright 160258600ac3SJames Wright @param[in] A `MATCEED` 160358600ac3SJames Wright @param[in] Y Input PETSc vector 160458600ac3SJames Wright @param[out] X Output PETSc vector 160558600ac3SJames Wright 160658600ac3SJames Wright @return An error code: 0 - success, otherwise - failure 160758600ac3SJames Wright **/ 160858600ac3SJames Wright PetscErrorCode MatMultTranspose_Ceed(Mat A, Vec Y, Vec X) { 160958600ac3SJames Wright MatCeedContext ctx; 161058600ac3SJames Wright 161158600ac3SJames Wright PetscFunctionBeginUser; 161258600ac3SJames Wright PetscCall(MatShellGetContext(A, &ctx)); 1613c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_mult_transpose, A, Y, X, NULL)); 161458600ac3SJames Wright 161558600ac3SJames Wright { 161658600ac3SJames Wright PetscMemType x_mem_type, y_mem_type; 161758600ac3SJames Wright Vec X_loc, Y_loc = ctx->Y_loc_transpose; 161858600ac3SJames Wright 161958600ac3SJames Wright // Get local vectors 162058600ac3SJames Wright if (!ctx->Y_loc_transpose) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); 162158600ac3SJames Wright PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); 162258600ac3SJames Wright 162358600ac3SJames Wright // Global-to-local 162458600ac3SJames Wright PetscCall(DMGlobalToLocal(ctx->dm_y, Y, INSERT_VALUES, Y_loc)); 162558600ac3SJames Wright 162658600ac3SJames Wright // Setup libCEED vectors 1627a7dac1d5SJames Wright PetscCall(VecReadPetscToCeed(Y_loc, &y_mem_type, ctx->y_loc)); 162858600ac3SJames Wright PetscCall(VecZeroEntries(X_loc)); 1629a7dac1d5SJames Wright PetscCall(VecPetscToCeed(X_loc, &x_mem_type, ctx->x_loc)); 163058600ac3SJames Wright 163158600ac3SJames Wright // Apply libCEED operator 1632c63b910fSJames Wright PetscCall(PetscLogEventBegin(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL)); 1633537ec908SJames Wright PetscCall(PetscLogGpuTimeBegin()); 163450f50432SJames Wright PetscCallCeed(ctx->ceed, CeedOperatorApplyAdd(ctx->op_mult_transpose, ctx->y_loc, ctx->x_loc, CEED_REQUEST_IMMEDIATE)); 1635ed5c6999SJames Wright // Log flops 1636ed5c6999SJames Wright if (PetscMemTypeDevice(ctx->mem_type)) PetscCall(PetscLogGpuFlops(ctx->flops_mult_transpose)); 1637ed5c6999SJames Wright else PetscCall(PetscLogFlops(ctx->flops_mult_transpose)); 163858600ac3SJames Wright PetscCall(PetscLogGpuTimeEnd()); 1639537ec908SJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_ceed_mult_transpose, A, Y, X, NULL)); 164058600ac3SJames Wright 164158600ac3SJames Wright // Restore PETSc vectors 1642a7dac1d5SJames Wright PetscCall(VecReadCeedToPetsc(ctx->y_loc, y_mem_type, Y_loc)); 1643a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(ctx->x_loc, x_mem_type, X_loc)); 164458600ac3SJames Wright 164558600ac3SJames Wright // Local-to-global 164658600ac3SJames Wright PetscCall(VecZeroEntries(X)); 164758600ac3SJames Wright PetscCall(DMLocalToGlobal(ctx->dm_x, X_loc, ADD_VALUES, X)); 164858600ac3SJames Wright 164958600ac3SJames Wright // Restore local vectors, as needed 165058600ac3SJames Wright if (!ctx->Y_loc_transpose) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); 165158600ac3SJames Wright PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); 165258600ac3SJames Wright } 165358600ac3SJames Wright 1654c63b910fSJames Wright PetscCall(PetscLogEventEnd(ctx->log_event_mult_transpose, A, Y, X, NULL)); 165558600ac3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 165658600ac3SJames Wright } 1657