xref: /libCEED/examples/solids/elasticity.h (revision 7be1e82bee12be3372d0edb3771f3d66d6ac97b8)
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 #if PETSC_VERSION_LT(3,14,0)
31 #  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)
32 #endif
33 
34 #ifndef PHYSICS_STRUCT
35 #define PHYSICS_STRUCT
36 typedef struct Physics_private *Physics;
37 struct Physics_private {
38   CeedScalar   nu;      // Poisson's ratio
39   CeedScalar   E;       // Young's Modulus
40 };
41 #endif
42 
43 // -----------------------------------------------------------------------------
44 // Command Line Options
45 // -----------------------------------------------------------------------------
46 // Problem options
47 typedef enum {
48   ELAS_LIN = 0, ELAS_HYPER_SS = 1, ELAS_HYPER_FS = 2
49 } problemType;
50 static const char *const problemTypes[] = {"linElas",
51                                            "hyperSS",
52                                            "hyperFS",
53                                            "problemType","ELAS_",0
54                                           };
55 static const char *const problemTypesForDisp[] = {"Linear elasticity",
56                                                   "Hyper elasticity small strain",
57                                                   "Hyper elasticity finite strain"
58                                                  };
59 
60 // Forcing function options
61 typedef enum {
62   FORCE_NONE = 0, FORCE_CONST = 1, FORCE_MMS = 2
63 } forcingType;
64 static const char *const forcingTypes[] = {"none",
65                                            "constant",
66                                            "mms",
67                                            "forcingType","FORCE_",0
68                                           };
69 static const char *const forcingTypesForDisp[] = {"None",
70                                                   "Constant",
71                                                   "Manufactured solution"
72                                                  };
73 
74 // Multigrid options
75 typedef enum {
76   MULTIGRID_LOGARITHMIC = 0, MULTIGRID_UNIFORM = 1, MULTIGRID_NONE = 2
77 } multigridType;
78 static const char *const multigridTypes [] = {"logarithmic",
79                                               "uniform",
80                                               "none",
81                                               "multigridType","MULTIGRID",0
82                                              };
83 static const char *const multigridTypesForDisp[] = {"P-multigrid, logarithmic coarsening",
84                                                     "P-multigrind, uniform coarsening",
85                                                     "No multigrid"
86                                                    };
87 
88 typedef PetscErrorCode BCFunc(PetscInt, PetscReal, const PetscReal *, PetscInt,
89                               PetscScalar *, void *);
90 // Note: These variables should be updated if additional boundary conditions
91 //         are added to boundary.c.
92 BCFunc BCMMS, BCZero, BCClamp;
93 
94 // -----------------------------------------------------------------------------
95 // Structs
96 // -----------------------------------------------------------------------------
97 // Units
98 typedef struct Units_private *Units;
99 struct Units_private {
100   // Fundamental units
101   PetscScalar meter;
102   PetscScalar kilogram;
103   PetscScalar second;
104   // Derived unit
105   PetscScalar Pascal;
106 };
107 
108 // Application context from user command line options
109 typedef struct AppCtx_private *AppCtx;
110 struct AppCtx_private {
111   char          ceedResource[PETSC_MAX_PATH_LEN];     // libCEED backend
112   char          meshFile[PETSC_MAX_PATH_LEN];         // exodusII mesh file
113   PetscBool     testMode;
114   PetscBool     viewSoln;
115   PetscBool     viewFinalSoln;
116   problemType   problemChoice;
117   forcingType   forcingChoice;
118   multigridType multigridChoice;
119   PetscScalar   nuSmoother;
120   PetscInt      degree;
121   PetscInt      qextra;
122   PetscInt      numLevels;
123   PetscInt      *levelDegrees;
124   PetscInt      numIncrements;                        // Number of steps
125   PetscInt      bcClampCount;
126   PetscInt      bcClampFaces[16];
127   PetscScalar   bcClampMax[16][7];
128   PetscInt      bcTractionCount;
129   PetscInt      bcTractionFaces[16];
130   PetscScalar   bcTractionVector[16][3];
131   PetscScalar   forcingVector[3];
132 };
133 
134 // Problem specific data
135 // *INDENT-OFF*
136 typedef struct {
137   CeedInt           qdatasize;
138   CeedQFunctionUser setupgeo, apply, jacob, energy, diagnostic;
139   const char        *setupgeofname, *applyfname, *jacobfname, *energyfname,
140                     *diagnosticfname;
141   CeedQuadMode      qmode;
142 } problemData;
143 // *INDENT-ON*
144 
145 // Data specific to each problem option
146 extern problemData problemOptions[3];
147 
148 // Forcing function data
149 typedef struct {
150   CeedQFunctionUser setupforcing;
151   const char        *setupforcingfname;
152 } forcingData;
153 
154 extern forcingData forcingOptions[3];
155 
156 // Data for PETSc Matshell
157 typedef struct UserMult_private *UserMult;
158 struct UserMult_private {
159   MPI_Comm        comm;
160   DM              dm;
161   Vec             Xloc, Yloc, NBCs;
162   CeedVector      Xceed, Yceed;
163   CeedOperator    op;
164   CeedQFunction   qf;
165   Ceed            ceed;
166   PetscScalar     loadIncrement;
167   CeedQFunctionContext ctxPhys, ctxPhysSmoother;
168 };
169 
170 // Data for Jacobian setup routine
171 typedef struct FormJacobCtx_private *FormJacobCtx;
172 struct FormJacobCtx_private {
173   UserMult     *jacobCtx;
174   PetscInt     numLevels;
175   SNES         snesCoarse;
176   Mat          *jacobMat, jacobMatCoarse;
177   Vec          Ucoarse;
178 };
179 
180 // Data for PETSc Prolongation/Restriction Matshell
181 typedef struct UserMultProlongRestr_private *UserMultProlongRestr;
182 struct UserMultProlongRestr_private {
183   MPI_Comm     comm;
184   DM           dmC, dmF;
185   Vec          locVecC, locVecF;
186   CeedVector   ceedVecC, ceedVecF;
187   CeedOperator opProlong, opRestrict;
188   Ceed         ceed;
189 };
190 
191 // libCEED data struct for level
192 typedef struct CeedData_private *CeedData;
193 struct CeedData_private {
194   Ceed                ceed;
195   CeedBasis           basisx, basisu, basisCtoF, basisEnergy, basisDiagnostic;
196   CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi,
197                       ErestrictGradui, ErestrictEnergy, ErestrictDiagnostic,
198                       ErestrictqdDiagnostici;
199   CeedQFunction       qfApply, qfJacob, qfEnergy, qfDiagnostic;
200   CeedOperator        opApply, opJacob, opRestrict, opProlong, opEnergy,
201                       opDiagnostic;
202   CeedVector          qdata, qdataDiagnostic, gradu, xceed, yceed, truesoln;
203 };
204 
205 // Translate PetscMemType to CeedMemType
206 static inline CeedMemType MemTypeP2C(PetscMemType mtype) {
207   return PetscMemTypeDevice(mtype) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
208 }
209 
210 // -----------------------------------------------------------------------------
211 // Process command line options
212 // -----------------------------------------------------------------------------
213 // Process general command line options
214 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx appCtx);
215 
216 // Process physics options
217 PetscErrorCode ProcessPhysics(MPI_Comm comm, Physics phys, Units units);
218 
219 // -----------------------------------------------------------------------------
220 // Setup DM
221 // -----------------------------------------------------------------------------
222 PetscErrorCode CreateBCLabel(DM dm, const char name[]);
223 
224 // Create FE by degree
225 PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc,
226                                      PetscBool isSimplex, const char prefix[],
227                                      PetscInt order, PetscFE *fem);
228 
229 // Read mesh and distribute DM in parallel
230 PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx appCtx, DM *dm);
231 
232 // Setup DM with FE space of appropriate degree
233 PetscErrorCode SetupDMByDegree(DM dm, AppCtx appCtx, PetscInt order,
234                                PetscBool boundary, PetscInt ncompu);
235 
236 // -----------------------------------------------------------------------------
237 // libCEED Functions
238 // -----------------------------------------------------------------------------
239 // Destroy libCEED objects
240 PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data);
241 
242 // Utility function - essential BC dofs are encoded in closure indices as -(i+1)
243 PetscInt Involute(PetscInt i);
244 
245 // Utility function to create local CEED restriction from DMPlex
246 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
247     CeedInt height, DMLabel domainLabel, CeedInt value,
248     CeedElemRestriction *Erestrict);
249 
250 // Utility function to get Ceed Restriction for each domain
251 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
252                                        DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
253                                        CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
254                                        CeedElemRestriction *restrictqdi);
255 
256 // Set up libCEED for a given degree
257 PetscErrorCode SetupLibceedFineLevel(DM dm, DM dmEnergy, DM dmDiagnostic,
258                                      Ceed ceed, AppCtx appCtx,
259                                      CeedQFunctionContext physCtx,
260                                      CeedData *data, PetscInt fineLevel,
261                                      PetscInt ncompu, PetscInt Ugsz,
262                                      PetscInt Ulocsz, CeedVector forceCeed,
263                                      CeedVector neumannCeed);
264 
265 // Set up libCEED multigrid level for a given degree
266 PetscErrorCode SetupLibceedLevel(DM dm, Ceed ceed, AppCtx appCtx,
267                                  CeedData *data, PetscInt level,
268                                  PetscInt ncompu, PetscInt Ugsz,
269                                  PetscInt Ulocsz, CeedVector fineMult);
270 
271 // Setup context data for Jacobian evaluation
272 PetscErrorCode SetupJacobianCtx(MPI_Comm comm, AppCtx appCtx, DM dm, Vec V,
273                                 Vec Vloc, CeedData ceedData, Ceed ceed,
274                                 CeedQFunctionContext ctxPhys,
275                                 CeedQFunctionContext ctxPhysSmoother,
276                                 UserMult jacobianCtx);
277 
278 // Setup context data for prolongation and restriction operators
279 PetscErrorCode SetupProlongRestrictCtx(MPI_Comm comm, AppCtx appCtx, DM dmC,
280                                        DM dmF, Vec VF, Vec VlocC, Vec VlocF,
281                                        CeedData ceedDataC, CeedData ceedDataF,
282                                        Ceed ceed,
283                                        UserMultProlongRestr prolongRestrCtx);
284 
285 // -----------------------------------------------------------------------------
286 // Jacobian setup
287 // -----------------------------------------------------------------------------
288 PetscErrorCode FormJacobian(SNES snes, Vec U, Mat J, Mat Jpre, void *ctx);
289 
290 // -----------------------------------------------------------------------------
291 // Solution output
292 // -----------------------------------------------------------------------------
293 PetscErrorCode ViewSolution(MPI_Comm comm, Vec U, PetscInt increment,
294                             PetscScalar loadIncrement);
295 
296 PetscErrorCode ViewDiagnosticQuantities(MPI_Comm comm, DM dmU,
297                                         UserMult user, Vec U,
298                                         CeedElemRestriction ErestrictDiagnostic);
299 
300 // -----------------------------------------------------------------------------
301 // libCEED Operators for MatShell
302 // -----------------------------------------------------------------------------
303 // This function uses libCEED to compute the local action of an operator
304 PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user);
305 
306 // This function uses libCEED to compute the non-linear residual
307 PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx);
308 
309 // This function uses libCEED to apply the Jacobian for assembly via a SNES
310 PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx);
311 
312 // This function uses libCEED to compute the action of the Jacobian
313 PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y);
314 
315 // This function uses libCEED to compute the action of the prolongation operator
316 PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y);
317 
318 // This function uses libCEED to compute the action of the restriction operator
319 PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y);
320 
321 // This function returns the computed diagonal of the operator
322 PetscErrorCode GetDiag_Ceed(Mat A, Vec D);
323 
324 // This function calculates the strain energy in the final solution
325 PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
326                                    CeedOperator opEnergy, Vec X,
327                                    PetscReal *energy);
328 
329 // -----------------------------------------------------------------------------
330 // Boundary Functions
331 // -----------------------------------------------------------------------------
332 // Note: If additional boundary conditions are added, an update is needed in
333 //         elasticity.h for the boundaryOptions variable.
334 
335 // BCMMS - boundary function
336 // Values on all points of the mesh is set based on given solution below
337 // for u[0], u[1], u[2]
338 PetscErrorCode BCMMS(PetscInt dim, PetscReal loadIncrement,
339                      const PetscReal coords[], PetscInt ncompu,
340                      PetscScalar *u, void *ctx);
341 
342 // BCClamp - fix boundary values with affine transformation at fraction of load
343 //   increment
344 PetscErrorCode BCClamp(PetscInt dim, PetscReal loadIncrement,
345                        const PetscReal coords[], PetscInt ncompu,
346                        PetscScalar *u, void *ctx);
347 
348 #endif //setup_h
349