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