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