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