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