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