xref: /libCEED/examples/solids/elasticity.h (revision 3d576824e8d990e1f48c6609089904bee9170514)
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 <ceed.h>
21 #include <petsc.h>
22 #include <petscdmplex.h>
23 #include <petscfe.h>
24 #include <petscksp.h>
25 #include <stdbool.h>
26 #include <string.h>
27 
28 #if PETSC_VERSION_LT(3,14,0)
29 #  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)
30 #endif
31 
32 #ifndef PHYSICS_STRUCT
33 #define PHYSICS_STRUCT
34 typedef struct Physics_private *Physics;
35 struct Physics_private {
36   CeedScalar   nu;      // Poisson's ratio
37   CeedScalar   E;       // Young's Modulus
38 };
39 #endif
40 
41 // -----------------------------------------------------------------------------
42 // Command Line Options
43 // -----------------------------------------------------------------------------
44 // Problem options
45 typedef enum {
46   ELAS_LIN = 0, ELAS_HYPER_SS = 1, ELAS_HYPER_FS = 2
47 } problemType;
48 static const char *const problemTypes[] = {"linElas",
49                                            "hyperSS",
50                                            "hyperFS",
51                                            "problemType","ELAS_",0
52                                           };
53 static const char *const problemTypesForDisp[] = {"Linear elasticity",
54                                                   "Hyper elasticity small strain",
55                                                   "Hyper elasticity finite strain"
56                                                  };
57 
58 // Forcing function options
59 typedef enum {
60   FORCE_NONE = 0, FORCE_CONST = 1, FORCE_MMS = 2
61 } forcingType;
62 static const char *const forcingTypes[] = {"none",
63                                            "constant",
64                                            "mms",
65                                            "forcingType","FORCE_",0
66                                           };
67 static const char *const forcingTypesForDisp[] = {"None",
68                                                   "Constant",
69                                                   "Manufactured solution"
70                                                  };
71 
72 // Multigrid options
73 typedef enum {
74   MULTIGRID_LOGARITHMIC = 0, MULTIGRID_UNIFORM = 1, MULTIGRID_NONE = 2
75 } multigridType;
76 static const char *const multigridTypes [] = {"logarithmic",
77                                               "uniform",
78                                               "none",
79                                               "multigridType","MULTIGRID",0
80                                              };
81 static const char *const multigridTypesForDisp[] = {"P-multigrid, logarithmic coarsening",
82                                                     "P-multigrind, uniform coarsening",
83                                                     "No multigrid"
84                                                    };
85 
86 typedef PetscErrorCode BCFunc(PetscInt, PetscReal, const PetscReal *, PetscInt,
87                               PetscScalar *, void *);
88 // Note: These variables should be updated if additional boundary conditions
89 //         are added to boundary.c.
90 BCFunc BCMMS, BCZero, BCClamp;
91 
92 // -----------------------------------------------------------------------------
93 // Structs
94 // -----------------------------------------------------------------------------
95 // Units
96 typedef struct Units_private *Units;
97 struct Units_private {
98   // Fundamental units
99   PetscScalar meter;
100   PetscScalar kilogram;
101   PetscScalar second;
102   // Derived unit
103   PetscScalar Pascal;
104 };
105 
106 // Application context from user command line options
107 typedef struct AppCtx_private *AppCtx;
108 struct AppCtx_private {
109   char          ceedResource[PETSC_MAX_PATH_LEN];     // libCEED backend
110   char          meshFile[PETSC_MAX_PATH_LEN];         // exodusII mesh file
111   PetscBool     testMode;
112   PetscBool     viewSoln;
113   PetscBool     viewFinalSoln;
114   problemType   problemChoice;
115   forcingType   forcingChoice;
116   multigridType multigridChoice;
117   PetscScalar   nuSmoother;
118   PetscInt      degree;
119   PetscInt      qextra;
120   PetscInt      numLevels;
121   PetscInt      *levelDegrees;
122   PetscInt      numIncrements;                        // Number of steps
123   PetscInt      bcClampCount;
124   PetscInt      bcClampFaces[16];
125   PetscScalar   bcClampMax[16][7];
126   PetscInt      bcTractionCount;
127   PetscInt      bcTractionFaces[16];
128   PetscScalar   bcTractionVector[16][3];
129   PetscScalar   forcingVector[3];
130 };
131 
132 // Problem specific data
133 // *INDENT-OFF*
134 typedef struct {
135   CeedInt           qdatasize;
136   CeedQFunctionUser setupgeo, apply, jacob, energy, diagnostic;
137   const char        *setupgeofname, *applyfname, *jacobfname, *energyfname,
138                     *diagnosticfname;
139   CeedQuadMode      qmode;
140 } problemData;
141 // *INDENT-ON*
142 
143 // Data specific to each problem option
144 extern problemData problemOptions[3];
145 
146 // Forcing function data
147 typedef struct {
148   CeedQFunctionUser setupforcing;
149   const char        *setupforcingfname;
150 } forcingData;
151 
152 extern forcingData forcingOptions[3];
153 
154 // Data for PETSc Matshell
155 typedef struct UserMult_private *UserMult;
156 struct UserMult_private {
157   MPI_Comm        comm;
158   DM              dm;
159   Vec             Xloc, Yloc, NBCs;
160   CeedVector      Xceed, Yceed;
161   CeedOperator    op;
162   CeedQFunction   qf;
163   Ceed            ceed;
164   PetscScalar     loadIncrement;
165   CeedQFunctionContext ctxPhys, ctxPhysSmoother;
166 };
167 
168 // Data for Jacobian setup routine
169 typedef struct FormJacobCtx_private *FormJacobCtx;
170 struct FormJacobCtx_private {
171   UserMult     *jacobCtx;
172   PetscInt     numLevels;
173   SNES         snesCoarse;
174   Mat          *jacobMat, jacobMatCoarse;
175   Vec          Ucoarse;
176 };
177 
178 // Data for PETSc Prolongation/Restriction Matshell
179 typedef struct UserMultProlongRestr_private *UserMultProlongRestr;
180 struct UserMultProlongRestr_private {
181   MPI_Comm     comm;
182   DM           dmC, dmF;
183   Vec          locVecC, locVecF;
184   CeedVector   ceedVecC, ceedVecF;
185   CeedOperator opProlong, opRestrict;
186   Ceed         ceed;
187 };
188 
189 // libCEED data struct for level
190 typedef struct CeedData_private *CeedData;
191 struct CeedData_private {
192   Ceed                ceed;
193   CeedBasis           basisx, basisu, basisCtoF, basisEnergy, basisDiagnostic;
194   CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi,
195                       ErestrictGradui, ErestrictEnergy, ErestrictDiagnostic,
196                       ErestrictqdDiagnostici;
197   CeedQFunction       qfApply, qfJacob, qfEnergy, qfDiagnostic;
198   CeedOperator        opApply, opJacob, opRestrict, opProlong, opEnergy,
199                       opDiagnostic;
200   CeedVector          qdata, qdataDiagnostic, gradu, xceed, yceed, truesoln;
201 };
202 
203 // Translate PetscMemType to CeedMemType
204 static inline CeedMemType MemTypeP2C(PetscMemType mtype) {
205   return PetscMemTypeDevice(mtype) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
206 }
207 
208 // -----------------------------------------------------------------------------
209 // Process command line options
210 // -----------------------------------------------------------------------------
211 // Process general command line options
212 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx appCtx);
213 
214 // Process physics options
215 PetscErrorCode ProcessPhysics(MPI_Comm comm, Physics phys, Units units);
216 
217 // -----------------------------------------------------------------------------
218 // Setup DM
219 // -----------------------------------------------------------------------------
220 PetscErrorCode CreateBCLabel(DM dm, const char name[]);
221 
222 // Create FE by degree
223 PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc,
224                                      PetscBool isSimplex, const char prefix[],
225                                      PetscInt order, PetscFE *fem);
226 
227 // Read mesh and distribute DM in parallel
228 PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx appCtx, DM *dm);
229 
230 // Setup DM with FE space of appropriate degree
231 PetscErrorCode SetupDMByDegree(DM dm, AppCtx appCtx, PetscInt order,
232                                PetscBool boundary, PetscInt ncompu);
233 
234 // -----------------------------------------------------------------------------
235 // libCEED Functions
236 // -----------------------------------------------------------------------------
237 // Destroy libCEED objects
238 PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data);
239 
240 // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
241 PetscInt Involute(PetscInt i);
242 
243 // Utility function to create local CEED restriction from DMPlex
244 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
245     CeedInt height, DMLabel domainLabel, CeedInt value,
246     CeedElemRestriction *Erestrict);
247 
248 // Utility function to get Ceed Restriction for each domain
249 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
250                                        DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
251                                        CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
252                                        CeedElemRestriction *restrictqdi);
253 
254 // Set up libCEED for a given degree
255 PetscErrorCode SetupLibceedFineLevel(DM dm, DM dmEnergy, DM dmDiagnostic,
256                                      Ceed ceed, AppCtx appCtx,
257                                      CeedQFunctionContext physCtx,
258                                      CeedData *data, PetscInt fineLevel,
259                                      PetscInt ncompu, PetscInt Ugsz,
260                                      PetscInt Ulocsz, CeedVector forceCeed,
261                                      CeedVector neumannCeed);
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