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 // ----------------------------------------------------------------------------- 91 // Structs 92 // ----------------------------------------------------------------------------- 93 // Units 94 typedef struct Units_private *Units; 95 struct Units_private { 96 // Fundamental units 97 PetscScalar meter; 98 PetscScalar kilogram; 99 PetscScalar second; 100 // Derived unit 101 PetscScalar Pascal; 102 }; 103 104 // Application context from user command line options 105 typedef struct AppCtx_private *AppCtx; 106 struct AppCtx_private { 107 char ceedResource[PETSC_MAX_PATH_LEN]; // libCEED backend 108 char ceedResourceFine[PETSC_MAX_PATH_LEN]; // libCEED for fine grid 109 char meshFile[PETSC_MAX_PATH_LEN]; // exodusII mesh file 110 PetscBool testMode; 111 PetscBool viewSoln; 112 PetscBool viewFinalSoln; 113 problemType problemChoice; 114 forcingType forcingChoice; 115 multigridType multigridChoice; 116 PetscInt degree; 117 PetscInt qextra; 118 PetscInt numLevels; 119 PetscInt *levelDegrees; 120 PetscInt numIncrements; // Number of steps 121 PetscInt bcClampFaces[16]; 122 PetscInt bcClampCount; 123 PetscScalar bcClampMax[16][7]; 124 PetscScalar forcingVector[3]; 125 }; 126 127 // Problem specific data 128 typedef struct { 129 CeedInt qdatasize; 130 CeedQFunctionUser setupgeo, apply, jacob, energy, diagnostic; 131 const char *setupgeofname, *applyfname, *jacobfname, *energyfname, 132 *diagnosticfname; 133 CeedQuadMode qmode; 134 } problemData; 135 136 // Data specific to each problem option 137 extern problemData problemOptions[3]; 138 139 // Forcing function data 140 typedef struct { 141 CeedQFunctionUser setupforcing; 142 const char *setupforcingfname; 143 } forcingData; 144 145 extern forcingData forcingOptions[3]; 146 147 // Data for PETSc Matshell 148 typedef struct UserMult_private *UserMult; 149 struct UserMult_private { 150 MPI_Comm comm; 151 DM dm; 152 Vec Xloc, Yloc; 153 CeedVector Xceed, Yceed; 154 CeedOperator op; 155 Ceed ceed; 156 PetscScalar loadIncrement; 157 }; 158 159 // Data for Jacobian setup routine 160 typedef struct FormJacobCtx_private *FormJacobCtx; 161 struct FormJacobCtx_private { 162 UserMult *jacobCtx; 163 PetscInt numLevels; 164 SNES snesCoarse; 165 Mat *jacobMat, jacobMatCoarse; 166 Vec Ucoarse; 167 }; 168 169 // Data for PETSc Prolongation/Restriction Matshell 170 typedef struct UserMultProlongRestr_private *UserMultProlongRestr; 171 struct UserMultProlongRestr_private { 172 MPI_Comm comm; 173 DM dmC, dmF; 174 Vec locVecC, locVecF, multVec; 175 CeedVector ceedVecC, ceedVecF; 176 CeedOperator opProlong, opRestrict; 177 Ceed ceed; 178 }; 179 180 // libCEED data struct for level 181 typedef struct CeedData_private *CeedData; 182 struct CeedData_private { 183 Ceed ceed; 184 CeedBasis basisx, basisu, basisCtoF, basisEnergy, basisDiagnostic; 185 CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi, 186 ErestrictGradui, ErestrictEnergy, ErestrictDiagnostic, 187 ErestrictqdDiagnostici; 188 CeedQFunction qfApply, qfJacob, qfEnergy, qfDiagnostic; 189 CeedOperator opApply, opJacob, opRestrict, opProlong, opEnergy, 190 opDiagnostic; 191 CeedVector qdata, qdataDiagnostic, gradu, xceed, yceed, truesoln; 192 }; 193 194 // ----------------------------------------------------------------------------- 195 // Process command line options 196 // ----------------------------------------------------------------------------- 197 // Process general command line options 198 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx appCtx); 199 200 // Process physics options 201 PetscErrorCode ProcessPhysics(MPI_Comm comm, Physics phys, Units units); 202 203 // ----------------------------------------------------------------------------- 204 // Setup DM 205 // ----------------------------------------------------------------------------- 206 PetscErrorCode CreateBCLabel(DM dm, const char name[]); 207 208 // Create FE by degree 209 PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc, 210 PetscBool isSimplex, const char prefix[], 211 PetscInt order, PetscFE *fem); 212 213 // Read mesh and distribute DM in parallel 214 PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx appCtx, DM *dm); 215 216 // Setup DM with FE space of appropriate degree 217 PetscErrorCode SetupDMByDegree(DM dm, AppCtx appCtx, PetscInt order, 218 PetscBool boundary, PetscInt ncompu); 219 220 // ----------------------------------------------------------------------------- 221 // libCEED Functions 222 // ----------------------------------------------------------------------------- 223 // Destroy libCEED objects 224 PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data); 225 226 // Get libCEED restriction data from DMPlex 227 PetscErrorCode CreateRestrictionPlex(Ceed ceed, CeedInt P, CeedInt ncomp, 228 CeedElemRestriction *Erestrict, DM dm); 229 230 // Set up libCEED for a given degree 231 PetscErrorCode SetupLibceedFineLevel(DM dm, DM dmEnergy, DM dmDiagnostic, 232 Ceed ceed, AppCtx appCtx, Physics phys, 233 CeedData *data, PetscInt fineLevel, 234 PetscInt ncompu, PetscInt Ugsz, 235 PetscInt Ulocsz, CeedVector forceCeed, 236 CeedQFunction qfRestrict, 237 CeedQFunction qfProlong); 238 239 // Set up libCEED for a given degree 240 PetscErrorCode SetupLibceedLevel(DM dm, Ceed ceed, AppCtx appCtx, Physics phys, 241 CeedData *data, PetscInt level, 242 PetscInt ncompu, PetscInt Ugsz, 243 PetscInt Ulocsz, CeedVector forceCeed, 244 CeedQFunction qfRestrict, 245 CeedQFunction qfProlong); 246 247 // Setup context data for Jacobian evaluation 248 PetscErrorCode SetupJacobianCtx(MPI_Comm comm, AppCtx appCtx, DM dm, Vec V, 249 Vec Vloc, CeedData ceedData, Ceed ceed, 250 UserMult jacobianCtx); 251 252 // Setup context data for prolongation and restriction operators 253 PetscErrorCode SetupProlongRestrictCtx(MPI_Comm comm, DM dmC, DM dmF, Vec VF, 254 Vec VlocC, Vec VlocF, CeedData ceedDataC, 255 CeedData ceedDataF, Ceed ceed, 256 UserMultProlongRestr prolongRestrCtx); 257 258 // ----------------------------------------------------------------------------- 259 // Jacobian setup 260 // ----------------------------------------------------------------------------- 261 PetscErrorCode FormJacobian(SNES snes, Vec U, Mat J, Mat Jpre, void *ctx); 262 263 // ----------------------------------------------------------------------------- 264 // Solution output 265 // ----------------------------------------------------------------------------- 266 PetscErrorCode ViewSolution(MPI_Comm comm, Vec U, PetscInt increment, 267 PetscScalar loadIncrement); 268 269 PetscErrorCode ViewDiagnosticQuantities(MPI_Comm comm, DM dmU, 270 UserMult user, Vec U, 271 CeedElemRestriction ErestrictDiagnostic); 272 273 // ----------------------------------------------------------------------------- 274 // libCEED Operators for MatShell 275 // ----------------------------------------------------------------------------- 276 // This function uses libCEED to compute the local action of an operator 277 PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user); 278 279 // This function uses libCEED to compute the non-linear residual 280 PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx); 281 282 // This function uses libCEED to apply the Jacobian for assembly via a SNES 283 PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx); 284 285 // This function uses libCEED to compute the action of the Jacobian 286 PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y); 287 288 // This function uses libCEED to compute the action of the prolongation operator 289 PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y); 290 291 // This function uses libCEED to compute the action of the restriction operator 292 PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y); 293 294 // This function returns the computed diagonal of the operator 295 PetscErrorCode GetDiag_Ceed(Mat A, Vec D); 296 297 // This function calculates the strain energy in the final solution 298 PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user, 299 CeedOperator opEnergy, Vec X, 300 PetscReal *energy); 301 302 // ----------------------------------------------------------------------------- 303 // Boundary Functions 304 // ----------------------------------------------------------------------------- 305 // Note: If additional boundary conditions are added, an update is needed in 306 // elasticity.h for the boundaryOptions variable. 307 308 // BCMMS - boundary function 309 // Values on all points of the mesh is set based on given solution below 310 // for u[0], u[1], u[2] 311 PetscErrorCode BCMMS(PetscInt dim, PetscReal loadIncrement, 312 const PetscReal coords[], PetscInt ncompu, 313 PetscScalar *u, void *ctx); 314 315 // BCClamp - fix boundary values with affine transformation at fraction of load 316 // increment 317 PetscErrorCode BCClamp(PetscInt dim, PetscReal loadIncrement, 318 const PetscReal coords[], PetscInt ncompu, 319 PetscScalar *u, void *ctx); 320 321 #endif //setup_h 322