1e83e87a5Sjeremylt #include "../include/matops.h" 2e83e87a5Sjeremylt #include "../include/petscutils.h" 3e83e87a5Sjeremylt 4e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 5e83e87a5Sjeremylt // This function returns the computed diagonal of the operator 6e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 7e83e87a5Sjeremylt PetscErrorCode MatGetDiag(Mat A, Vec D) { 8e83e87a5Sjeremylt PetscErrorCode ierr; 9*d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx; 10e83e87a5Sjeremylt 11e83e87a5Sjeremylt PetscFunctionBeginUser; 12e83e87a5Sjeremylt 13*d4d45553Srezgarshakeri ierr = MatShellGetContext(A, &op_apply_ctx); CHKERRQ(ierr); 14e83e87a5Sjeremylt 15e83e87a5Sjeremylt // Compute Diagonal via libCEED 16f1c40530SJeremy L Thompson PetscScalar *y; 179b072555Sjeremylt PetscMemType mem_type; 18e83e87a5Sjeremylt 19e83e87a5Sjeremylt // -- Place PETSc vector in libCEED vector 20*d4d45553Srezgarshakeri ierr = VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &mem_type); CHKERRQ(ierr); 21*d4d45553Srezgarshakeri CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(mem_type), 22*d4d45553Srezgarshakeri CEED_USE_POINTER, y); 23e83e87a5Sjeremylt 24e83e87a5Sjeremylt // -- Compute Diagonal 25*d4d45553Srezgarshakeri CeedOperatorLinearAssembleDiagonal(op_apply_ctx->op, op_apply_ctx->y_ceed, 26e83e87a5Sjeremylt CEED_REQUEST_IMMEDIATE); 27e83e87a5Sjeremylt 28e83e87a5Sjeremylt // -- Local-to-Global 29*d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(mem_type), NULL); 30*d4d45553Srezgarshakeri ierr = VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y); CHKERRQ(ierr); 31e83e87a5Sjeremylt ierr = VecZeroEntries(D); CHKERRQ(ierr); 32*d4d45553Srezgarshakeri ierr = DMLocalToGlobal(op_apply_ctx->dm, op_apply_ctx->Y_loc, ADD_VALUES, D); 33*d4d45553Srezgarshakeri CHKERRQ(ierr); 34e83e87a5Sjeremylt 35e83e87a5Sjeremylt PetscFunctionReturn(0); 36e83e87a5Sjeremylt }; 37e83e87a5Sjeremylt 38e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 39e83e87a5Sjeremylt // This function uses libCEED to compute the action of the Laplacian with 40e83e87a5Sjeremylt // Dirichlet boundary conditions 41e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 42*d4d45553Srezgarshakeri PetscErrorCode ApplyLocal_Ceed(Vec X, Vec Y, 43*d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx) { 44e83e87a5Sjeremylt PetscErrorCode ierr; 45e83e87a5Sjeremylt PetscScalar *x, *y; 469b072555Sjeremylt PetscMemType x_mem_type, y_mem_type; 47e83e87a5Sjeremylt 48e83e87a5Sjeremylt PetscFunctionBeginUser; 49e83e87a5Sjeremylt 50e83e87a5Sjeremylt // Global-to-local 51*d4d45553Srezgarshakeri ierr = DMGlobalToLocal(op_apply_ctx->dm, X, INSERT_VALUES, op_apply_ctx->X_loc); 52*d4d45553Srezgarshakeri CHKERRQ(ierr); 53e83e87a5Sjeremylt 54e83e87a5Sjeremylt // Setup libCEED vectors 55*d4d45553Srezgarshakeri ierr = VecGetArrayReadAndMemType(op_apply_ctx->X_loc, (const PetscScalar **)&x, 569b072555Sjeremylt &x_mem_type); CHKERRQ(ierr); 57*d4d45553Srezgarshakeri ierr = VecGetArrayAndMemType(op_apply_ctx->Y_loc, &y, &y_mem_type); 58*d4d45553Srezgarshakeri CHKERRQ(ierr); 59*d4d45553Srezgarshakeri CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), 60*d4d45553Srezgarshakeri CEED_USE_POINTER, x); 61*d4d45553Srezgarshakeri CeedVectorSetArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), 62*d4d45553Srezgarshakeri CEED_USE_POINTER, y); 63e83e87a5Sjeremylt 64e83e87a5Sjeremylt // Apply libCEED operator 65*d4d45553Srezgarshakeri CeedOperatorApply(op_apply_ctx->op, op_apply_ctx->x_ceed, op_apply_ctx->y_ceed, 66*d4d45553Srezgarshakeri CEED_REQUEST_IMMEDIATE); 67e83e87a5Sjeremylt 68e83e87a5Sjeremylt // Restore PETSc vectors 69*d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(x_mem_type), NULL); 70*d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->y_ceed, MemTypeP2C(y_mem_type), NULL); 71*d4d45553Srezgarshakeri ierr = VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, 72*d4d45553Srezgarshakeri (const PetscScalar **)&x); 73e83e87a5Sjeremylt CHKERRQ(ierr); 74*d4d45553Srezgarshakeri ierr = VecRestoreArrayAndMemType(op_apply_ctx->Y_loc, &y); CHKERRQ(ierr); 75e83e87a5Sjeremylt 76e83e87a5Sjeremylt // Local-to-global 77e83e87a5Sjeremylt ierr = VecZeroEntries(Y); CHKERRQ(ierr); 78*d4d45553Srezgarshakeri ierr = DMLocalToGlobal(op_apply_ctx->dm, op_apply_ctx->Y_loc, ADD_VALUES, Y); 79*d4d45553Srezgarshakeri CHKERRQ(ierr); 80e83e87a5Sjeremylt 81e83e87a5Sjeremylt PetscFunctionReturn(0); 82e83e87a5Sjeremylt }; 83e83e87a5Sjeremylt 84e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 85e83e87a5Sjeremylt // This function wraps the libCEED operator for a MatShell 86e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 87e83e87a5Sjeremylt PetscErrorCode MatMult_Ceed(Mat A, Vec X, Vec Y) { 88e83e87a5Sjeremylt PetscErrorCode ierr; 89*d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx; 90e83e87a5Sjeremylt 91e83e87a5Sjeremylt PetscFunctionBeginUser; 92e83e87a5Sjeremylt 93*d4d45553Srezgarshakeri ierr = MatShellGetContext(A, &op_apply_ctx); CHKERRQ(ierr); 94e83e87a5Sjeremylt 95e83e87a5Sjeremylt // libCEED for local action of residual evaluator 96*d4d45553Srezgarshakeri ierr = ApplyLocal_Ceed(X, Y, op_apply_ctx); CHKERRQ(ierr); 97e83e87a5Sjeremylt 98e83e87a5Sjeremylt PetscFunctionReturn(0); 99e83e87a5Sjeremylt }; 100e83e87a5Sjeremylt 101e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 102e83e87a5Sjeremylt // This function wraps the libCEED operator for a SNES residual evaluation 103e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 104e83e87a5Sjeremylt PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) { 105e83e87a5Sjeremylt PetscErrorCode ierr; 106*d4d45553Srezgarshakeri OperatorApplyContext op_apply_ctx = (OperatorApplyContext)ctx; 107e83e87a5Sjeremylt 108e83e87a5Sjeremylt PetscFunctionBeginUser; 109e83e87a5Sjeremylt 110e83e87a5Sjeremylt // libCEED for local action of residual evaluator 111*d4d45553Srezgarshakeri ierr = ApplyLocal_Ceed(X, Y, op_apply_ctx); CHKERRQ(ierr); 112e83e87a5Sjeremylt 113e83e87a5Sjeremylt PetscFunctionReturn(0); 114e83e87a5Sjeremylt }; 115e83e87a5Sjeremylt 116e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 117e83e87a5Sjeremylt // This function uses libCEED to compute the action of the prolongation operator 118e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 119e83e87a5Sjeremylt PetscErrorCode MatMult_Prolong(Mat A, Vec X, Vec Y) { 120e83e87a5Sjeremylt PetscErrorCode ierr; 121*d4d45553Srezgarshakeri ProlongRestrContext pr_restr_ctx; 122e83e87a5Sjeremylt PetscScalar *c, *f; 1239b072555Sjeremylt PetscMemType c_mem_type, f_mem_type; 124e83e87a5Sjeremylt 125e83e87a5Sjeremylt PetscFunctionBeginUser; 126e83e87a5Sjeremylt 127*d4d45553Srezgarshakeri ierr = MatShellGetContext(A, &pr_restr_ctx); CHKERRQ(ierr); 128e83e87a5Sjeremylt 129e83e87a5Sjeremylt // Global-to-local 130*d4d45553Srezgarshakeri ierr = VecZeroEntries(pr_restr_ctx->loc_vec_c); CHKERRQ(ierr); 131*d4d45553Srezgarshakeri ierr = DMGlobalToLocal(pr_restr_ctx->dmc, X, INSERT_VALUES, 132*d4d45553Srezgarshakeri pr_restr_ctx->loc_vec_c); 133e83e87a5Sjeremylt CHKERRQ(ierr); 134e83e87a5Sjeremylt 135e83e87a5Sjeremylt // Setup libCEED vectors 136*d4d45553Srezgarshakeri ierr = VecGetArrayReadAndMemType(pr_restr_ctx->loc_vec_c, 137*d4d45553Srezgarshakeri (const PetscScalar **)&c, 1389b072555Sjeremylt &c_mem_type); CHKERRQ(ierr); 139*d4d45553Srezgarshakeri ierr = VecGetArrayAndMemType(pr_restr_ctx->loc_vec_f, &f, &f_mem_type); 140*d4d45553Srezgarshakeri CHKERRQ(ierr); 141*d4d45553Srezgarshakeri CeedVectorSetArray(pr_restr_ctx->ceed_vec_c, MemTypeP2C(c_mem_type), 142*d4d45553Srezgarshakeri CEED_USE_POINTER, c); 143*d4d45553Srezgarshakeri CeedVectorSetArray(pr_restr_ctx->ceed_vec_f, MemTypeP2C(f_mem_type), 144*d4d45553Srezgarshakeri CEED_USE_POINTER, f); 145e83e87a5Sjeremylt 146e83e87a5Sjeremylt // Apply libCEED operator 147*d4d45553Srezgarshakeri CeedOperatorApply(pr_restr_ctx->op_prolong, pr_restr_ctx->ceed_vec_c, 148*d4d45553Srezgarshakeri pr_restr_ctx->ceed_vec_f, CEED_REQUEST_IMMEDIATE); 149e83e87a5Sjeremylt 150e83e87a5Sjeremylt // Restore PETSc vectors 151*d4d45553Srezgarshakeri CeedVectorTakeArray(pr_restr_ctx->ceed_vec_c, MemTypeP2C(c_mem_type), NULL); 152*d4d45553Srezgarshakeri CeedVectorTakeArray(pr_restr_ctx->ceed_vec_f, MemTypeP2C(f_mem_type), NULL); 153*d4d45553Srezgarshakeri ierr = VecRestoreArrayReadAndMemType(pr_restr_ctx->loc_vec_c, 154*d4d45553Srezgarshakeri (const PetscScalar **)&c); 155e83e87a5Sjeremylt CHKERRQ(ierr); 156*d4d45553Srezgarshakeri ierr = VecRestoreArrayAndMemType(pr_restr_ctx->loc_vec_f, &f); CHKERRQ(ierr); 157e83e87a5Sjeremylt 158e83e87a5Sjeremylt // Multiplicity 159*d4d45553Srezgarshakeri ierr = VecPointwiseMult(pr_restr_ctx->loc_vec_f, pr_restr_ctx->loc_vec_f, 160*d4d45553Srezgarshakeri pr_restr_ctx->mult_vec); 161e83e87a5Sjeremylt 162e83e87a5Sjeremylt // Local-to-global 163e83e87a5Sjeremylt ierr = VecZeroEntries(Y); CHKERRQ(ierr); 164*d4d45553Srezgarshakeri ierr = DMLocalToGlobal(pr_restr_ctx->dmf, pr_restr_ctx->loc_vec_f, ADD_VALUES, 165*d4d45553Srezgarshakeri Y); CHKERRQ(ierr); 166e83e87a5Sjeremylt 167e83e87a5Sjeremylt PetscFunctionReturn(0); 168e83e87a5Sjeremylt }; 169e83e87a5Sjeremylt 170e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 171e83e87a5Sjeremylt // This function uses libCEED to compute the action of the restriction operator 172e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 173e83e87a5Sjeremylt PetscErrorCode MatMult_Restrict(Mat A, Vec X, Vec Y) { 174e83e87a5Sjeremylt PetscErrorCode ierr; 175*d4d45553Srezgarshakeri ProlongRestrContext pr_restr_ctx; 176e83e87a5Sjeremylt PetscScalar *c, *f; 1779b072555Sjeremylt PetscMemType c_mem_type, f_mem_type; 178e83e87a5Sjeremylt 179e83e87a5Sjeremylt PetscFunctionBeginUser; 180e83e87a5Sjeremylt 181*d4d45553Srezgarshakeri ierr = MatShellGetContext(A, &pr_restr_ctx); CHKERRQ(ierr); 182e83e87a5Sjeremylt 183e83e87a5Sjeremylt // Global-to-local 184*d4d45553Srezgarshakeri ierr = VecZeroEntries(pr_restr_ctx->loc_vec_f); CHKERRQ(ierr); 185*d4d45553Srezgarshakeri ierr = DMGlobalToLocal(pr_restr_ctx->dmf, X, INSERT_VALUES, 186*d4d45553Srezgarshakeri pr_restr_ctx->loc_vec_f); 187e83e87a5Sjeremylt CHKERRQ(ierr); 188e83e87a5Sjeremylt 189e83e87a5Sjeremylt // Multiplicity 190*d4d45553Srezgarshakeri ierr = VecPointwiseMult(pr_restr_ctx->loc_vec_f, pr_restr_ctx->loc_vec_f, 191*d4d45553Srezgarshakeri pr_restr_ctx->mult_vec); 192e83e87a5Sjeremylt CHKERRQ(ierr); 193e83e87a5Sjeremylt 194e83e87a5Sjeremylt // Setup libCEED vectors 195*d4d45553Srezgarshakeri ierr = VecGetArrayReadAndMemType(pr_restr_ctx->loc_vec_f, 196*d4d45553Srezgarshakeri (const PetscScalar **)&f, 1979b072555Sjeremylt &f_mem_type); CHKERRQ(ierr); 198*d4d45553Srezgarshakeri ierr = VecGetArrayAndMemType(pr_restr_ctx->loc_vec_c, &c, &c_mem_type); 199*d4d45553Srezgarshakeri CHKERRQ(ierr); 200*d4d45553Srezgarshakeri CeedVectorSetArray(pr_restr_ctx->ceed_vec_f, MemTypeP2C(f_mem_type), 201*d4d45553Srezgarshakeri CEED_USE_POINTER, f); 202*d4d45553Srezgarshakeri CeedVectorSetArray(pr_restr_ctx->ceed_vec_c, MemTypeP2C(c_mem_type), 203*d4d45553Srezgarshakeri CEED_USE_POINTER, c); 204e83e87a5Sjeremylt 205e83e87a5Sjeremylt // Apply CEED operator 206*d4d45553Srezgarshakeri CeedOperatorApply(pr_restr_ctx->op_restrict, pr_restr_ctx->ceed_vec_f, 207*d4d45553Srezgarshakeri pr_restr_ctx->ceed_vec_c, CEED_REQUEST_IMMEDIATE); 208e83e87a5Sjeremylt 209e83e87a5Sjeremylt // Restore PETSc vectors 210*d4d45553Srezgarshakeri CeedVectorTakeArray(pr_restr_ctx->ceed_vec_c, MemTypeP2C(c_mem_type), NULL); 211*d4d45553Srezgarshakeri CeedVectorTakeArray(pr_restr_ctx->ceed_vec_f, MemTypeP2C(f_mem_type), NULL); 212*d4d45553Srezgarshakeri ierr = VecRestoreArrayReadAndMemType(pr_restr_ctx->loc_vec_f, 213*d4d45553Srezgarshakeri (const PetscScalar **)&f); CHKERRQ(ierr); 214*d4d45553Srezgarshakeri ierr = VecRestoreArrayAndMemType(pr_restr_ctx->loc_vec_c, &c); CHKERRQ(ierr); 215e83e87a5Sjeremylt 216e83e87a5Sjeremylt // Local-to-global 217e83e87a5Sjeremylt ierr = VecZeroEntries(Y); CHKERRQ(ierr); 218*d4d45553Srezgarshakeri ierr = DMLocalToGlobal(pr_restr_ctx->dmc, pr_restr_ctx->loc_vec_c, ADD_VALUES, 219*d4d45553Srezgarshakeri Y); 220e83e87a5Sjeremylt CHKERRQ(ierr); 221e83e87a5Sjeremylt 222e83e87a5Sjeremylt PetscFunctionReturn(0); 223e83e87a5Sjeremylt }; 224e83e87a5Sjeremylt 225e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 226e83e87a5Sjeremylt // This function calculates the error in the final solution 227e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 228*d4d45553Srezgarshakeri PetscErrorCode ComputeErrorMax(OperatorApplyContext op_apply_ctx, 229*d4d45553Srezgarshakeri CeedOperator op_error, 230e83e87a5Sjeremylt Vec X, CeedVector target, 2319b072555Sjeremylt PetscScalar *max_error) { 232e83e87a5Sjeremylt PetscErrorCode ierr; 233e83e87a5Sjeremylt PetscScalar *x; 2349b072555Sjeremylt PetscMemType mem_type; 235e83e87a5Sjeremylt CeedVector collocated_error; 2361f9221feSJeremy L Thompson CeedSize length; 237e83e87a5Sjeremylt 238e83e87a5Sjeremylt PetscFunctionBeginUser; 239e83e87a5Sjeremylt CeedVectorGetLength(target, &length); 240*d4d45553Srezgarshakeri CeedVectorCreate(op_apply_ctx->ceed, length, &collocated_error); 241e83e87a5Sjeremylt 242e83e87a5Sjeremylt // Global-to-local 243*d4d45553Srezgarshakeri ierr = DMGlobalToLocal(op_apply_ctx->dm, X, INSERT_VALUES, op_apply_ctx->X_loc); 244*d4d45553Srezgarshakeri CHKERRQ(ierr); 245e83e87a5Sjeremylt 246e83e87a5Sjeremylt // Setup CEED vector 247*d4d45553Srezgarshakeri ierr = VecGetArrayAndMemType(op_apply_ctx->X_loc, &x, &mem_type); CHKERRQ(ierr); 248*d4d45553Srezgarshakeri CeedVectorSetArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), 249*d4d45553Srezgarshakeri CEED_USE_POINTER, x); 250e83e87a5Sjeremylt 251e83e87a5Sjeremylt // Apply CEED operator 252*d4d45553Srezgarshakeri CeedOperatorApply(op_error, op_apply_ctx->x_ceed, collocated_error, 253e83e87a5Sjeremylt CEED_REQUEST_IMMEDIATE); 254e83e87a5Sjeremylt 255e83e87a5Sjeremylt // Restore PETSc vector 256*d4d45553Srezgarshakeri CeedVectorTakeArray(op_apply_ctx->x_ceed, MemTypeP2C(mem_type), NULL); 257*d4d45553Srezgarshakeri ierr = VecRestoreArrayReadAndMemType(op_apply_ctx->X_loc, 258*d4d45553Srezgarshakeri (const PetscScalar **)&x); 259e83e87a5Sjeremylt CHKERRQ(ierr); 260e83e87a5Sjeremylt 261e83e87a5Sjeremylt // Reduce max error 2629b072555Sjeremylt *max_error = 0; 263e83e87a5Sjeremylt const CeedScalar *e; 264e83e87a5Sjeremylt CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e); 265e83e87a5Sjeremylt for (CeedInt i=0; i<length; i++) { 2669b072555Sjeremylt *max_error = PetscMax(*max_error, PetscAbsScalar(e[i])); 267e83e87a5Sjeremylt } 268e83e87a5Sjeremylt CeedVectorRestoreArrayRead(collocated_error, &e); 2699b072555Sjeremylt ierr = MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX, 270*d4d45553Srezgarshakeri op_apply_ctx->comm); 271e83e87a5Sjeremylt CHKERRQ(ierr); 272e83e87a5Sjeremylt 273e83e87a5Sjeremylt // Cleanup 274e83e87a5Sjeremylt CeedVectorDestroy(&collocated_error); 275e83e87a5Sjeremylt 276e83e87a5Sjeremylt PetscFunctionReturn(0); 277e83e87a5Sjeremylt }; 278e83e87a5Sjeremylt 279e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 280