1 // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2 // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3 // reserved. See files LICENSE and NOTICE for details. 4 // 5 // This file is part of CEED, a collection of benchmarks, miniapps, software 6 // libraries and APIs for efficient high-order finite element and spectral 7 // element discretizations for exascale applications. For more information and 8 // source code availability see http://github.com/ceed. 9 // 10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11 // a collaborative effort of two U.S. Department of Energy organizations (Office 12 // of Science and the National Nuclear Security Administration) responsible for 13 // the planning and preparation of a capable exascale ecosystem, including 14 // software, applications, hardware, advanced system engineering and early 15 // testbed platforms, in support of the nation's exascale computing imperative. 16 17 #ifndef setup_h 18 #define setup_h 19 20 #include <stdbool.h> 21 #include <string.h> 22 23 #include <petsc.h> 24 #include <petscdmplex.h> 25 #include <petscksp.h> 26 #include <petscfe.h> 27 28 #include <ceed.h> 29 30 #ifndef PHYSICS_STRUCT 31 #define PHYSICS_STRUCT 32 typedef struct Physics_private *Physics; 33 struct Physics_private { 34 CeedScalar nu; // Poisson's ratio 35 CeedScalar E; // Young's Modulus 36 }; 37 #endif 38 39 // ----------------------------------------------------------------------------- 40 // Command Line Options 41 // ----------------------------------------------------------------------------- 42 // Problem options 43 typedef enum { 44 ELAS_LIN = 0, ELAS_HYPER_SS = 1, ELAS_HYPER_FS = 2 45 } problemType; 46 static const char *const problemTypes[] = {"linElas", 47 "hyperSS", 48 "hyperFS", 49 "problemType","ELAS_",0 50 }; 51 static const char *const problemTypesForDisp[] = {"Linear elasticity", 52 "Hyper elasticity small strain", 53 "Hyper elasticity finite strain" 54 }; 55 56 // Forcing function options 57 typedef enum { 58 FORCE_NONE = 0, FORCE_CONST = 1, FORCE_MMS = 2 59 } forcingType; 60 static const char *const forcingTypes[] = {"none", 61 "constant", 62 "mms", 63 "forcingType","FORCE_",0 64 }; 65 static const char *const forcingTypesForDisp[] = {"None", 66 "Constant", 67 "Manufactured solution" 68 }; 69 70 // Multigrid options 71 typedef enum { 72 MULTIGRID_LOGARITHMIC = 0, MULTIGRID_UNIFORM = 1, MULTIGRID_NONE = 2 73 } multigridType; 74 static const char *const multigridTypes [] = {"logarithmic", 75 "uniform", 76 "none", 77 "multigridType","MULTIGRID",0 78 }; 79 static const char *const multigridTypesForDisp[] = {"P-multigrid, logarithmic coarsening", 80 "P-multigrind, uniform coarsening", 81 "No multigrid" 82 }; 83 84 typedef PetscErrorCode BCFunc(PetscInt, PetscReal, const PetscReal *, PetscInt, 85 PetscScalar *, void *); 86 // Note: These variables should be updated if additional boundary conditions 87 // are added to boundary.c. 88 BCFunc BCMMS, BCZero, BCClamp; 89 90 // MemType Options 91 static const char *const memTypes[] = {"host","device","memType", 92 "CEED_MEM_",0 93 }; 94 95 // ----------------------------------------------------------------------------- 96 // Structs 97 // ----------------------------------------------------------------------------- 98 // Units 99 typedef struct Units_private *Units; 100 struct Units_private { 101 // Fundamental units 102 PetscScalar meter; 103 PetscScalar kilogram; 104 PetscScalar second; 105 // Derived unit 106 PetscScalar Pascal; 107 }; 108 109 // Application context from user command line options 110 typedef struct AppCtx_private *AppCtx; 111 struct AppCtx_private { 112 char ceedResource[PETSC_MAX_PATH_LEN]; // libCEED backend 113 char meshFile[PETSC_MAX_PATH_LEN]; // exodusII mesh file 114 PetscBool testMode; 115 PetscBool viewSoln; 116 PetscBool viewFinalSoln; 117 problemType problemChoice; 118 forcingType forcingChoice; 119 multigridType multigridChoice; 120 PetscScalar nuSmoother; 121 PetscInt degree; 122 PetscInt qextra; 123 PetscInt numLevels; 124 PetscInt *levelDegrees; 125 PetscInt numIncrements; // Number of steps 126 PetscInt bcClampFaces[16]; 127 PetscInt bcClampCount; 128 PetscScalar bcClampMax[16][7]; 129 PetscScalar forcingVector[3]; 130 PetscBool petscHaveCuda, setMemTypeRequest; 131 CeedMemType memTypeRequested; 132 }; 133 134 // Problem specific data 135 // *INDENT-OFF* 136 typedef struct { 137 CeedInt qdatasize; 138 CeedQFunctionUser setupgeo, apply, jacob, energy, diagnostic; 139 const char *setupgeofname, *applyfname, *jacobfname, *energyfname, 140 *diagnosticfname; 141 CeedQuadMode qmode; 142 } problemData; 143 // *INDENT-ON* 144 145 // Data specific to each problem option 146 extern problemData problemOptions[3]; 147 148 // Forcing function data 149 typedef struct { 150 CeedQFunctionUser setupforcing; 151 const char *setupforcingfname; 152 } forcingData; 153 154 extern forcingData forcingOptions[3]; 155 156 // Data for PETSc Matshell 157 typedef struct UserMult_private *UserMult; 158 struct UserMult_private { 159 MPI_Comm comm; 160 DM dm; 161 Vec Xloc, Yloc; 162 CeedVector Xceed, Yceed; 163 CeedOperator op; 164 CeedQFunction qf; 165 Ceed ceed; 166 PetscScalar loadIncrement; 167 Physics phys, physSmoother; 168 CeedMemType memType; 169 int (*VecGetArray)(Vec, PetscScalar **); 170 int (*VecGetArrayRead)(Vec, const PetscScalar **); 171 int (*VecRestoreArray)(Vec, PetscScalar **); 172 int (*VecRestoreArrayRead)(Vec, const PetscScalar **); 173 }; 174 175 // Data for Jacobian setup routine 176 typedef struct FormJacobCtx_private *FormJacobCtx; 177 struct FormJacobCtx_private { 178 UserMult *jacobCtx; 179 PetscInt numLevels; 180 SNES snesCoarse; 181 Mat *jacobMat, jacobMatCoarse; 182 Vec Ucoarse; 183 }; 184 185 // Data for PETSc Prolongation/Restriction Matshell 186 typedef struct UserMultProlongRestr_private *UserMultProlongRestr; 187 struct UserMultProlongRestr_private { 188 MPI_Comm comm; 189 DM dmC, dmF; 190 Vec locVecC, locVecF; 191 CeedVector ceedVecC, ceedVecF; 192 CeedOperator opProlong, opRestrict; 193 Ceed ceed; 194 CeedMemType memType; 195 int (*VecGetArray)(Vec, PetscScalar **); 196 int (*VecGetArrayRead)(Vec, const PetscScalar **); 197 int (*VecRestoreArray)(Vec, PetscScalar **); 198 int (*VecRestoreArrayRead)(Vec, const PetscScalar **); 199 }; 200 201 // libCEED data struct for level 202 typedef struct CeedData_private *CeedData; 203 struct CeedData_private { 204 Ceed ceed; 205 CeedBasis basisx, basisu, basisCtoF, basisEnergy, basisDiagnostic; 206 CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi, 207 ErestrictGradui, ErestrictEnergy, ErestrictDiagnostic, 208 ErestrictqdDiagnostici; 209 CeedQFunction qfApply, qfJacob, qfEnergy, qfDiagnostic; 210 CeedOperator opApply, opJacob, opRestrict, opProlong, opEnergy, 211 opDiagnostic; 212 CeedVector qdata, qdataDiagnostic, gradu, xceed, yceed, truesoln; 213 }; 214 215 // ----------------------------------------------------------------------------- 216 // Process command line options 217 // ----------------------------------------------------------------------------- 218 // Process general command line options 219 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx appCtx); 220 221 // Process physics options 222 PetscErrorCode ProcessPhysics(MPI_Comm comm, Physics phys, Units units); 223 224 // ----------------------------------------------------------------------------- 225 // Setup DM 226 // ----------------------------------------------------------------------------- 227 PetscErrorCode CreateBCLabel(DM dm, const char name[]); 228 229 // Create FE by degree 230 PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc, 231 PetscBool isSimplex, const char prefix[], 232 PetscInt order, PetscFE *fem); 233 234 // Read mesh and distribute DM in parallel 235 PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx appCtx, DM *dm); 236 237 // Setup DM with FE space of appropriate degree 238 PetscErrorCode SetupDMByDegree(DM dm, AppCtx appCtx, PetscInt order, 239 PetscBool boundary, PetscInt ncompu); 240 241 // ----------------------------------------------------------------------------- 242 // libCEED Functions 243 // ----------------------------------------------------------------------------- 244 // Destroy libCEED objects 245 PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data); 246 247 // Get libCEED restriction data from DMPlex 248 PetscErrorCode CreateRestrictionPlex(Ceed ceed, CeedInt P, CeedInt ncomp, 249 CeedElemRestriction *Erestrict, DM dm); 250 251 // Set up libCEED for a given degree 252 PetscErrorCode SetupLibceedFineLevel(DM dm, DM dmEnergy, DM dmDiagnostic, 253 Ceed ceed, AppCtx appCtx, Physics phys, 254 CeedData *data, PetscInt fineLevel, 255 PetscInt ncompu, PetscInt Ugsz, 256 PetscInt Ulocsz, CeedVector forceCeed, 257 CeedQFunction qfRestrict, 258 CeedQFunction qfProlong); 259 260 // Set up libCEED multigrid level for a given degree 261 PetscErrorCode SetupLibceedLevel(DM dm, Ceed ceed, AppCtx appCtx, Physics phys, 262 CeedData *data, PetscInt level, 263 PetscInt ncompu, PetscInt Ugsz, 264 PetscInt Ulocsz, CeedVector fineMult, 265 CeedQFunction qfRestrict, 266 CeedQFunction qfProlong); 267 268 // Setup context data for Jacobian evaluation 269 PetscErrorCode SetupJacobianCtx(MPI_Comm comm, AppCtx appCtx, DM dm, Vec V, 270 Vec Vloc, CeedData ceedData, Ceed ceed, 271 Physics phys, Physics physSmoother, 272 UserMult jacobianCtx); 273 274 // Setup context data for prolongation and restriction operators 275 PetscErrorCode SetupProlongRestrictCtx(MPI_Comm comm, AppCtx appCtx, DM dmC, 276 DM dmF, Vec VF, Vec VlocC, Vec VlocF, 277 CeedData ceedDataC, CeedData ceedDataF, 278 Ceed ceed, 279 UserMultProlongRestr prolongRestrCtx); 280 281 // ----------------------------------------------------------------------------- 282 // Jacobian setup 283 // ----------------------------------------------------------------------------- 284 PetscErrorCode FormJacobian(SNES snes, Vec U, Mat J, Mat Jpre, void *ctx); 285 286 // ----------------------------------------------------------------------------- 287 // Solution output 288 // ----------------------------------------------------------------------------- 289 PetscErrorCode ViewSolution(MPI_Comm comm, Vec U, PetscInt increment, 290 PetscScalar loadIncrement); 291 292 PetscErrorCode ViewDiagnosticQuantities(MPI_Comm comm, DM dmU, 293 UserMult user, Vec U, 294 CeedElemRestriction ErestrictDiagnostic); 295 296 // ----------------------------------------------------------------------------- 297 // libCEED Operators for MatShell 298 // ----------------------------------------------------------------------------- 299 // This function uses libCEED to compute the local action of an operator 300 PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user); 301 302 // This function uses libCEED to compute the non-linear residual 303 PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx); 304 305 // This function uses libCEED to apply the Jacobian for assembly via a SNES 306 PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx); 307 308 // This function uses libCEED to compute the action of the Jacobian 309 PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y); 310 311 // This function uses libCEED to compute the action of the prolongation operator 312 PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y); 313 314 // This function uses libCEED to compute the action of the restriction operator 315 PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y); 316 317 // This function returns the computed diagonal of the operator 318 PetscErrorCode GetDiag_Ceed(Mat A, Vec D); 319 320 // This function calculates the strain energy in the final solution 321 PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user, 322 CeedOperator opEnergy, Vec X, 323 PetscReal *energy); 324 325 // ----------------------------------------------------------------------------- 326 // Boundary Functions 327 // ----------------------------------------------------------------------------- 328 // Note: If additional boundary conditions are added, an update is needed in 329 // elasticity.h for the boundaryOptions variable. 330 331 // BCMMS - boundary function 332 // Values on all points of the mesh is set based on given solution below 333 // for u[0], u[1], u[2] 334 PetscErrorCode BCMMS(PetscInt dim, PetscReal loadIncrement, 335 const PetscReal coords[], PetscInt ncompu, 336 PetscScalar *u, void *ctx); 337 338 // BCClamp - fix boundary values with affine transformation at fraction of load 339 // increment 340 PetscErrorCode BCClamp(PetscInt dim, PetscReal loadIncrement, 341 const PetscReal coords[], PetscInt ncompu, 342 PetscScalar *u, void *ctx); 343 344 #endif //setup_h 345