1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details. 4ccaff030SJeremy L Thompson // 5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software 6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral 7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and 8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed. 9ccaff030SJeremy L Thompson // 10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office 12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for 13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including 14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early 15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative. 16ccaff030SJeremy L Thompson 17ccaff030SJeremy L Thompson /// @file 18ccaff030SJeremy L Thompson /// Matrix shell operations for solid mechanics example using PETSc 19ccaff030SJeremy L Thompson 20ccaff030SJeremy L Thompson #include "../elasticity.h" 21ccaff030SJeremy L Thompson 22ccaff030SJeremy L Thompson // ----------------------------------------------------------------------------- 23ccaff030SJeremy L Thompson // libCEED Operators for MatShell 24ccaff030SJeremy L Thompson // ----------------------------------------------------------------------------- 25ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator 26ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user) { 27ccaff030SJeremy L Thompson PetscErrorCode ierr; 28ccaff030SJeremy L Thompson PetscScalar *x, *y; 29*d1d35e2fSjeremylt PetscMemType x_mem_type, y_mem_type; 30ccaff030SJeremy L Thompson 31ccaff030SJeremy L Thompson PetscFunctionBeginUser; 32ccaff030SJeremy L Thompson 33ccaff030SJeremy L Thompson // Global-to-local 34*d1d35e2fSjeremylt ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr); 35*d1d35e2fSjeremylt ierr = VecZeroEntries(user->Y_loc); CHKERRQ(ierr); 36ccaff030SJeremy L Thompson 37ccaff030SJeremy L Thompson // Setup CEED vectors 38*d1d35e2fSjeremylt ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, 39*d1d35e2fSjeremylt &x_mem_type); CHKERRQ(ierr); 40*d1d35e2fSjeremylt ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr); 41*d1d35e2fSjeremylt CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); 42*d1d35e2fSjeremylt CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y); 43ccaff030SJeremy L Thompson 44ccaff030SJeremy L Thompson // Apply CEED operator 45*d1d35e2fSjeremylt CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, CEED_REQUEST_IMMEDIATE); 46ccaff030SJeremy L Thompson 47ccaff030SJeremy L Thompson // Restore PETSc vectors 48*d1d35e2fSjeremylt CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL); 49*d1d35e2fSjeremylt CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL); 50*d1d35e2fSjeremylt ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x); 51ccaff030SJeremy L Thompson CHKERRQ(ierr); 52*d1d35e2fSjeremylt ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr); 53ccaff030SJeremy L Thompson 54ccaff030SJeremy L Thompson // Local-to-global 55ccaff030SJeremy L Thompson ierr = VecZeroEntries(Y); CHKERRQ(ierr); 56*d1d35e2fSjeremylt ierr = DMLocalToGlobal(user->dm, user->Y_loc, ADD_VALUES, Y); CHKERRQ(ierr); 57ccaff030SJeremy L Thompson 58ccaff030SJeremy L Thompson PetscFunctionReturn(0); 59ccaff030SJeremy L Thompson }; 60ccaff030SJeremy L Thompson 61ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual 62ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx) { 63ccaff030SJeremy L Thompson PetscErrorCode ierr; 64ccaff030SJeremy L Thompson UserMult user = (UserMult)ctx; 65ccaff030SJeremy L Thompson 66ccaff030SJeremy L Thompson PetscFunctionBeginUser; 67ccaff030SJeremy L Thompson 68ccaff030SJeremy L Thompson // Use computed BCs 69*d1d35e2fSjeremylt ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr); 70*d1d35e2fSjeremylt ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc, 71*d1d35e2fSjeremylt user->load_increment, NULL, NULL, NULL); 72ccaff030SJeremy L Thompson CHKERRQ(ierr); 73ccaff030SJeremy L Thompson 74ccaff030SJeremy L Thompson // libCEED for local action of residual evaluator 75ccaff030SJeremy L Thompson ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr); 76ccaff030SJeremy L Thompson 77fe394131Sjeremylt // Neumann BCs 78*d1d35e2fSjeremylt if (user->neumann_bcs) { 79*d1d35e2fSjeremylt ierr = VecAXPY(Y, -user->load_increment, user->neumann_bcs); CHKERRQ(ierr); 80fe394131Sjeremylt } 81fe394131Sjeremylt 82ccaff030SJeremy L Thompson PetscFunctionReturn(0); 83ccaff030SJeremy L Thompson }; 84ccaff030SJeremy L Thompson 85ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES 86ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx) { 87ccaff030SJeremy L Thompson PetscErrorCode ierr; 88ccaff030SJeremy L Thompson UserMult user = (UserMult)ctx; 89ccaff030SJeremy L Thompson 90ccaff030SJeremy L Thompson PetscFunctionBeginUser; 91ccaff030SJeremy L Thompson 9262e9c006SJeremy L Thompson // Zero boundary values 93*d1d35e2fSjeremylt ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr); 94ccaff030SJeremy L Thompson 95ccaff030SJeremy L Thompson // libCEED for local action of residual evaluator 96ccaff030SJeremy L Thompson ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr); 97ccaff030SJeremy L Thompson 98ccaff030SJeremy L Thompson PetscFunctionReturn(0); 99ccaff030SJeremy L Thompson }; 100ccaff030SJeremy L Thompson 101ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian 102ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y) { 103ccaff030SJeremy L Thompson PetscErrorCode ierr; 104ccaff030SJeremy L Thompson UserMult user; 105ccaff030SJeremy L Thompson 106ccaff030SJeremy L Thompson PetscFunctionBeginUser; 107ccaff030SJeremy L Thompson 108ccaff030SJeremy L Thompson // Zero boundary values 109ccaff030SJeremy L Thompson ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 110*d1d35e2fSjeremylt ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr); 111ccaff030SJeremy L Thompson 112ccaff030SJeremy L Thompson // libCEED for local action of Jacobian 113ccaff030SJeremy L Thompson ierr = ApplyLocalCeedOp(X, Y, user); CHKERRQ(ierr); 114ccaff030SJeremy L Thompson 115ccaff030SJeremy L Thompson PetscFunctionReturn(0); 116ccaff030SJeremy L Thompson }; 117ccaff030SJeremy L Thompson 118ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator 119ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y) { 120ccaff030SJeremy L Thompson PetscErrorCode ierr; 121ccaff030SJeremy L Thompson UserMultProlongRestr user; 122ccaff030SJeremy L Thompson PetscScalar *c, *f; 123*d1d35e2fSjeremylt PetscMemType c_mem_type, f_mem_type; 124ccaff030SJeremy L Thompson 125ccaff030SJeremy L Thompson PetscFunctionBeginUser; 126ccaff030SJeremy L Thompson 127ccaff030SJeremy L Thompson ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 128ccaff030SJeremy L Thompson 129ccaff030SJeremy L Thompson // Global-to-local 130*d1d35e2fSjeremylt ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr); 131*d1d35e2fSjeremylt ierr = DMGlobalToLocal(user->dm_c, X, INSERT_VALUES, user->loc_vec_c); 132ccaff030SJeremy L Thompson CHKERRQ(ierr); 133*d1d35e2fSjeremylt ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr); 134ccaff030SJeremy L Thompson 135ccaff030SJeremy L Thompson // Setup CEED vectors 136*d1d35e2fSjeremylt ierr = VecGetArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c, 137*d1d35e2fSjeremylt &c_mem_type); CHKERRQ(ierr); 138*d1d35e2fSjeremylt ierr = VecGetArrayAndMemType(user->loc_vec_f, &f, &f_mem_type); CHKERRQ(ierr); 139*d1d35e2fSjeremylt CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER, 140*d1d35e2fSjeremylt c); 141*d1d35e2fSjeremylt CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER, 142*d1d35e2fSjeremylt f); 143ccaff030SJeremy L Thompson 144ccaff030SJeremy L Thompson // Apply CEED operator 145*d1d35e2fSjeremylt CeedOperatorApply(user->op_prolong, user->ceed_vec_c, user->ceed_vec_f, 146ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 147ccaff030SJeremy L Thompson 148ccaff030SJeremy L Thompson // Restore PETSc vectors 149*d1d35e2fSjeremylt CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL); 150*d1d35e2fSjeremylt CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL); 151*d1d35e2fSjeremylt ierr = VecRestoreArrayReadAndMemType(user->loc_vec_c, (const PetscScalar **)&c); 152ccaff030SJeremy L Thompson CHKERRQ(ierr); 153*d1d35e2fSjeremylt ierr = VecRestoreArrayAndMemType(user->loc_vec_f, &f); CHKERRQ(ierr); 154ccaff030SJeremy L Thompson 155ccaff030SJeremy L Thompson // Local-to-global 156ccaff030SJeremy L Thompson ierr = VecZeroEntries(Y); CHKERRQ(ierr); 157*d1d35e2fSjeremylt ierr = DMLocalToGlobal(user->dm_f, user->loc_vec_f, ADD_VALUES, Y); 158ccaff030SJeremy L Thompson CHKERRQ(ierr); 159ccaff030SJeremy L Thompson 160ccaff030SJeremy L Thompson PetscFunctionReturn(0); 161ccaff030SJeremy L Thompson } 162ccaff030SJeremy L Thompson 163ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator 164ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y) { 165ccaff030SJeremy L Thompson PetscErrorCode ierr; 166ccaff030SJeremy L Thompson UserMultProlongRestr user; 167ccaff030SJeremy L Thompson PetscScalar *c, *f; 168*d1d35e2fSjeremylt PetscMemType c_mem_type, f_mem_type; 169ccaff030SJeremy L Thompson 170ccaff030SJeremy L Thompson PetscFunctionBeginUser; 171ccaff030SJeremy L Thompson 172ccaff030SJeremy L Thompson ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 173ccaff030SJeremy L Thompson 174ccaff030SJeremy L Thompson // Global-to-local 175*d1d35e2fSjeremylt ierr = VecZeroEntries(user->loc_vec_f); CHKERRQ(ierr); 176*d1d35e2fSjeremylt ierr = DMGlobalToLocal(user->dm_f, X, INSERT_VALUES, user->loc_vec_f); 177ccaff030SJeremy L Thompson CHKERRQ(ierr); 178*d1d35e2fSjeremylt ierr = VecZeroEntries(user->loc_vec_c); CHKERRQ(ierr); 179ccaff030SJeremy L Thompson 180ccaff030SJeremy L Thompson // Setup CEED vectors 181*d1d35e2fSjeremylt ierr = VecGetArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f, 182*d1d35e2fSjeremylt &f_mem_type); CHKERRQ(ierr); 183*d1d35e2fSjeremylt ierr = VecGetArrayAndMemType(user->loc_vec_c, &c, &c_mem_type); CHKERRQ(ierr); 184*d1d35e2fSjeremylt CeedVectorSetArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), CEED_USE_POINTER, 185*d1d35e2fSjeremylt f); 186*d1d35e2fSjeremylt CeedVectorSetArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), CEED_USE_POINTER, 187*d1d35e2fSjeremylt c); 188ccaff030SJeremy L Thompson 189ccaff030SJeremy L Thompson // Apply CEED operator 190*d1d35e2fSjeremylt CeedOperatorApply(user->op_restrict, user->ceed_vec_f, user->ceed_vec_c, 191ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 192ccaff030SJeremy L Thompson 193ccaff030SJeremy L Thompson // Restore PETSc vectors 194*d1d35e2fSjeremylt CeedVectorTakeArray(user->ceed_vec_f, MemTypeP2C(f_mem_type), NULL); 195*d1d35e2fSjeremylt CeedVectorTakeArray(user->ceed_vec_c, MemTypeP2C(c_mem_type), NULL); 196*d1d35e2fSjeremylt ierr = VecRestoreArrayReadAndMemType(user->loc_vec_f, (const PetscScalar **)&f); 197ccaff030SJeremy L Thompson CHKERRQ(ierr); 198*d1d35e2fSjeremylt ierr = VecRestoreArrayAndMemType(user->loc_vec_c, &c); CHKERRQ(ierr); 199ccaff030SJeremy L Thompson 200ccaff030SJeremy L Thompson // Local-to-global 201ccaff030SJeremy L Thompson ierr = VecZeroEntries(Y); CHKERRQ(ierr); 202*d1d35e2fSjeremylt ierr = DMLocalToGlobal(user->dm_c, user->loc_vec_c, ADD_VALUES, Y); 203ccaff030SJeremy L Thompson CHKERRQ(ierr); 204ccaff030SJeremy L Thompson 205ccaff030SJeremy L Thompson PetscFunctionReturn(0); 206ccaff030SJeremy L Thompson }; 2072d93065eSjeremylt 208ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator 209ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D) { 210ccaff030SJeremy L Thompson PetscErrorCode ierr; 211ccaff030SJeremy L Thompson UserMult user; 212ccaff030SJeremy L Thompson 213ccaff030SJeremy L Thompson PetscFunctionBeginUser; 214ccaff030SJeremy L Thompson 215ccaff030SJeremy L Thompson ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 216ccaff030SJeremy L Thompson 217f7b4142eSJeremy L Thompson // -- Set physics context 218*d1d35e2fSjeremylt if (user->ctx_phys_smoother) 219*d1d35e2fSjeremylt CeedQFunctionSetContext(user->qf, user->ctx_phys_smoother); 220f7b4142eSJeremy L Thompson 2212bba3ffaSJeremy L Thompson // Compute Diagonal via libCEED 2222bba3ffaSJeremy L Thompson PetscScalar *x; 223*d1d35e2fSjeremylt PetscMemType x_mem_type; 2242bba3ffaSJeremy L Thompson 2252bba3ffaSJeremy L Thompson // -- Place PETSc vector in libCEED vector 226*d1d35e2fSjeremylt ierr = VecGetArrayAndMemType(user->X_loc, &x, &x_mem_type); CHKERRQ(ierr); 227*d1d35e2fSjeremylt CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); 2282bba3ffaSJeremy L Thompson 2292bba3ffaSJeremy L Thompson // -- Compute Diagonal 230*d1d35e2fSjeremylt CeedOperatorLinearAssembleDiagonal(user->op, user->x_ceed, 231ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 232ccaff030SJeremy L Thompson 233f7b4142eSJeremy L Thompson // -- Reset physics context 234*d1d35e2fSjeremylt if (user->ctx_phys_smoother) 235*d1d35e2fSjeremylt CeedQFunctionSetContext(user->qf, user->ctx_phys); 236f7b4142eSJeremy L Thompson 237ccaff030SJeremy L Thompson // -- Local-to-Global 238*d1d35e2fSjeremylt CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL); 239*d1d35e2fSjeremylt ierr = VecRestoreArrayAndMemType(user->X_loc, &x); CHKERRQ(ierr); 240ccaff030SJeremy L Thompson ierr = VecZeroEntries(D); CHKERRQ(ierr); 241*d1d35e2fSjeremylt ierr = DMLocalToGlobal(user->dm, user->X_loc, ADD_VALUES, D); CHKERRQ(ierr); 242ccaff030SJeremy L Thompson 2432bba3ffaSJeremy L Thompson // Cleanup 244*d1d35e2fSjeremylt ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr); 245ccaff030SJeremy L Thompson 246ccaff030SJeremy L Thompson PetscFunctionReturn(0); 247ccaff030SJeremy L Thompson }; 2482d93065eSjeremylt 2492d93065eSjeremylt // This function calculates the strain energy in the final solution 250a3c02c40SJeremy L Thompson PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user, 251*d1d35e2fSjeremylt CeedOperator op_energy, Vec X, 252a3c02c40SJeremy L Thompson PetscReal *energy) { 2532d93065eSjeremylt PetscErrorCode ierr; 2542d93065eSjeremylt PetscScalar *x; 255*d1d35e2fSjeremylt PetscMemType x_mem_type; 2562d93065eSjeremylt CeedInt length; 2572d93065eSjeremylt 2582d93065eSjeremylt PetscFunctionBeginUser; 2592d93065eSjeremylt 2602d93065eSjeremylt // Global-to-local 261*d1d35e2fSjeremylt ierr = VecZeroEntries(user->X_loc); CHKERRQ(ierr); 262*d1d35e2fSjeremylt ierr = DMGlobalToLocal(user->dm, X, INSERT_VALUES, user->X_loc); CHKERRQ(ierr); 263*d1d35e2fSjeremylt ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, user->X_loc, 264*d1d35e2fSjeremylt user->load_increment, NULL, NULL, NULL); 2655c25879aSJeremy L Thompson CHKERRQ(ierr); 2662d93065eSjeremylt 267a3c02c40SJeremy L Thompson // Setup libCEED input vector 268*d1d35e2fSjeremylt ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, 269*d1d35e2fSjeremylt &x_mem_type); CHKERRQ(ierr); 270*d1d35e2fSjeremylt CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); 2712d93065eSjeremylt 272a3c02c40SJeremy L Thompson // Setup libCEED output vector 273*d1d35e2fSjeremylt Vec E_loc; 274*d1d35e2fSjeremylt CeedVector e_loc; 275*d1d35e2fSjeremylt ierr = DMCreateLocalVector(dmEnergy, &E_loc); CHKERRQ(ierr); 276*d1d35e2fSjeremylt ierr = VecGetSize(E_loc, &length); CHKERRQ(ierr); 277*d1d35e2fSjeremylt ierr = VecDestroy(&E_loc); CHKERRQ(ierr); 278*d1d35e2fSjeremylt CeedVectorCreate(user->ceed, length, &e_loc); 279a3c02c40SJeremy L Thompson 2802d93065eSjeremylt // Apply libCEED operator 281*d1d35e2fSjeremylt CeedOperatorApply(op_energy, user->x_ceed, e_loc, CEED_REQUEST_IMMEDIATE); 2822d93065eSjeremylt 2832d93065eSjeremylt // Restore PETSc vector 284*d1d35e2fSjeremylt CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL); 285*d1d35e2fSjeremylt ierr = VecRestoreArrayRead(user->X_loc, (const PetscScalar **)&x); 2862d93065eSjeremylt CHKERRQ(ierr); 2872d93065eSjeremylt 2882d93065eSjeremylt // Reduce max error 2892d93065eSjeremylt const CeedScalar *e; 290*d1d35e2fSjeremylt CeedVectorGetArrayRead(e_loc, CEED_MEM_HOST, &e); 2912d93065eSjeremylt (*energy) = 0; 2922d93065eSjeremylt for (CeedInt i=0; i<length; i++) 2932d93065eSjeremylt (*energy) += e[i]; 294*d1d35e2fSjeremylt CeedVectorRestoreArrayRead(e_loc, &e); 295*d1d35e2fSjeremylt CeedVectorDestroy(&e_loc); 2962d93065eSjeremylt 2972d93065eSjeremylt ierr = MPI_Allreduce(MPI_IN_PLACE, energy, 1, MPIU_REAL, MPIU_SUM, 2982d93065eSjeremylt user->comm); CHKERRQ(ierr); 2992d93065eSjeremylt 3002d93065eSjeremylt PetscFunctionReturn(0); 3012d93065eSjeremylt }; 302