xref: /libCEED/examples/solids/elasticity.h (revision 5c25879a8ff80373d5be013119a867f98bfe2528)
1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details.
4ccaff030SJeremy L Thompson //
5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and
8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed.
9ccaff030SJeremy L Thompson //
10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16ccaff030SJeremy L Thompson 
17ccaff030SJeremy L Thompson #ifndef setup_h
18ccaff030SJeremy L Thompson #define setup_h
19ccaff030SJeremy L Thompson 
20ccaff030SJeremy L Thompson #include <stdbool.h>
21ccaff030SJeremy L Thompson #include <string.h>
22ccaff030SJeremy L Thompson 
23ccaff030SJeremy L Thompson #include <petsc.h>
24ccaff030SJeremy L Thompson #include <petscdmplex.h>
25ccaff030SJeremy L Thompson #include <petscksp.h>
26ccaff030SJeremy L Thompson #include <petscfe.h>
27ccaff030SJeremy L Thompson 
28ccaff030SJeremy L Thompson #include <ceed.h>
29ccaff030SJeremy L Thompson 
30ccaff030SJeremy L Thompson #ifndef PHYSICS_STRUCT
31ccaff030SJeremy L Thompson #define PHYSICS_STRUCT
32ccaff030SJeremy L Thompson typedef struct Physics_private *Physics;
33ccaff030SJeremy L Thompson struct Physics_private {
34ccaff030SJeremy L Thompson   CeedScalar   nu;      // Poisson's ratio
35ccaff030SJeremy L Thompson   CeedScalar   E;       // Young's Modulus
36ccaff030SJeremy L Thompson };
37ccaff030SJeremy L Thompson #endif
38ccaff030SJeremy L Thompson 
39ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
40ccaff030SJeremy L Thompson // Command Line Options
41ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
42ccaff030SJeremy L Thompson // Problem options
43ccaff030SJeremy L Thompson typedef enum {
44ccaff030SJeremy L Thompson   ELAS_LIN = 0, ELAS_HYPER_SS = 1, ELAS_HYPER_FS = 2
45ccaff030SJeremy L Thompson } problemType;
46ccaff030SJeremy L Thompson static const char *const problemTypes[] = {"linElas",
47ccaff030SJeremy L Thompson                                            "hyperSS",
48ccaff030SJeremy L Thompson                                            "hyperFS",
49ccaff030SJeremy L Thompson                                            "problemType","ELAS_",0
50ccaff030SJeremy L Thompson                                           };
51ccaff030SJeremy L Thompson static const char *const problemTypesForDisp[] = {"Linear elasticity",
52ccaff030SJeremy L Thompson                                                   "Hyper elasticity small strain",
53ccaff030SJeremy L Thompson                                                   "Hyper elasticity finite strain"
54ccaff030SJeremy L Thompson                                                  };
55ccaff030SJeremy L Thompson 
56ccaff030SJeremy L Thompson // Forcing function options
57ccaff030SJeremy L Thompson typedef enum {
58ccaff030SJeremy L Thompson   FORCE_NONE = 0, FORCE_CONST = 1, FORCE_MMS = 2
59ccaff030SJeremy L Thompson } forcingType;
60ccaff030SJeremy L Thompson static const char *const forcingTypes[] = {"none",
61ccaff030SJeremy L Thompson                                            "constant",
62ccaff030SJeremy L Thompson                                            "mms",
63ccaff030SJeremy L Thompson                                            "forcingType","FORCE_",0
64ccaff030SJeremy L Thompson                                           };
65ccaff030SJeremy L Thompson static const char *const forcingTypesForDisp[] = {"None",
66ccaff030SJeremy L Thompson                                                   "Constant",
67ccaff030SJeremy L Thompson                                                   "Manufactured solution"
68ccaff030SJeremy L Thompson                                                  };
69ccaff030SJeremy L Thompson 
70ccaff030SJeremy L Thompson // Multigrid options
71ccaff030SJeremy L Thompson typedef enum {
72ccaff030SJeremy L Thompson   MULTIGRID_LOGARITHMIC = 0, MULTIGRID_UNIFORM = 1, MULTIGRID_NONE = 2
73ccaff030SJeremy L Thompson } multigridType;
74ccaff030SJeremy L Thompson static const char *const multigridTypes [] = {"logarithmic",
75ccaff030SJeremy L Thompson                                               "uniform",
76ccaff030SJeremy L Thompson                                               "none",
77ccaff030SJeremy L Thompson                                               "multigridType","MULTIGRID",0
78ccaff030SJeremy L Thompson                                              };
79ccaff030SJeremy L Thompson static const char *const multigridTypesForDisp[] = {"P-multigrid, logarithmic coarsening",
80ccaff030SJeremy L Thompson                                                     "P-multigrind, uniform coarsening",
81ccaff030SJeremy L Thompson                                                     "No multigrid"
82ccaff030SJeremy L Thompson                                                    };
83ccaff030SJeremy L Thompson 
84ccaff030SJeremy L Thompson typedef PetscErrorCode BCFunc(PetscInt, PetscReal, const PetscReal *, PetscInt,
85ccaff030SJeremy L Thompson                               PetscScalar *, void *);
86ccaff030SJeremy L Thompson // Note: These variables should be updated if additional boundary conditions
87ccaff030SJeremy L Thompson //         are added to boundary.c.
88ccaff030SJeremy L Thompson BCFunc BCMMS, BCZero, BCClamp;
89ccaff030SJeremy L Thompson 
90ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
91ccaff030SJeremy L Thompson // Structs
92ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
93ccaff030SJeremy L Thompson // Units
94ccaff030SJeremy L Thompson typedef struct Units_private *Units;
95ccaff030SJeremy L Thompson struct Units_private {
96ccaff030SJeremy L Thompson   // Fundamental units
97ccaff030SJeremy L Thompson   PetscScalar meter;
98ccaff030SJeremy L Thompson   PetscScalar kilogram;
99ccaff030SJeremy L Thompson   PetscScalar second;
100ccaff030SJeremy L Thompson   // Derived unit
101ccaff030SJeremy L Thompson   PetscScalar Pascal;
102ccaff030SJeremy L Thompson };
103ccaff030SJeremy L Thompson 
104ccaff030SJeremy L Thompson // Application context from user command line options
105ccaff030SJeremy L Thompson typedef struct AppCtx_private *AppCtx;
106ccaff030SJeremy L Thompson struct AppCtx_private {
107ccaff030SJeremy L Thompson   char          ceedResource[PETSC_MAX_PATH_LEN];     // libCEED backend
108ccaff030SJeremy L Thompson   char          ceedResourceFine[PETSC_MAX_PATH_LEN]; // libCEED for fine grid
109ccaff030SJeremy L Thompson   char          meshFile[PETSC_MAX_PATH_LEN];         // exodusII mesh file
110ccaff030SJeremy L Thompson   PetscBool     testMode;
111ccaff030SJeremy L Thompson   PetscBool     viewSoln;
11212b49b89SJeremy L Thompson   PetscBool     viewFinalSoln;
113ccaff030SJeremy L Thompson   problemType   problemChoice;
114ccaff030SJeremy L Thompson   forcingType   forcingChoice;
115ccaff030SJeremy L Thompson   multigridType multigridChoice;
116ccaff030SJeremy L Thompson   PetscInt      degree;
11712b49b89SJeremy L Thompson   PetscInt      qextra;
118ccaff030SJeremy L Thompson   PetscInt      numLevels;
119ccaff030SJeremy L Thompson   PetscInt      *levelDegrees;
120ccaff030SJeremy L Thompson   PetscInt      numIncrements;                        // Number of steps
121d642641fSjeremylt   PetscInt      bcClampFaces[16];
122d642641fSjeremylt   PetscInt      bcClampCount;
123d642641fSjeremylt   PetscScalar   bcClampMax[16][7];
124c83beb56Sjeremylt   PetscScalar   forcingVector[3];
125ccaff030SJeremy L Thompson };
126ccaff030SJeremy L Thompson 
127ccaff030SJeremy L Thompson // Problem specific data
128ccaff030SJeremy L Thompson typedef struct {
129ccaff030SJeremy L Thompson   CeedInt           qdatasize;
130*5c25879aSJeremy L Thompson   CeedQFunctionUser setupgeo, apply, jacob, energy, diagnostic;
131*5c25879aSJeremy L Thompson   const char        *setupgeofname, *applyfname, *jacobfname, *energyfname,
132*5c25879aSJeremy L Thompson                     *diagnosticfname;
133ccaff030SJeremy L Thompson   CeedQuadMode      qmode;
134ccaff030SJeremy L Thompson } problemData;
135ccaff030SJeremy L Thompson 
136ccaff030SJeremy L Thompson // Data specific to each problem option
1379e9355a0SJed Brown extern problemData problemOptions[3];
138ccaff030SJeremy L Thompson 
139ccaff030SJeremy L Thompson // Forcing function data
140ccaff030SJeremy L Thompson typedef struct {
141ccaff030SJeremy L Thompson   CeedQFunctionUser setupforcing;
142ccaff030SJeremy L Thompson   const char        *setupforcingfname;
143ccaff030SJeremy L Thompson } forcingData;
144ccaff030SJeremy L Thompson 
1459e9355a0SJed Brown extern forcingData forcingOptions[3];
146ccaff030SJeremy L Thompson 
147ccaff030SJeremy L Thompson // Data for PETSc Matshell
148ccaff030SJeremy L Thompson typedef struct UserMult_private *UserMult;
149ccaff030SJeremy L Thompson struct UserMult_private {
150ccaff030SJeremy L Thompson   MPI_Comm     comm;
151ccaff030SJeremy L Thompson   DM           dm;
152ccaff030SJeremy L Thompson   Vec          Xloc, Yloc;
153ccaff030SJeremy L Thompson   CeedVector   Xceed, Yceed;
154ccaff030SJeremy L Thompson   CeedOperator op;
155ccaff030SJeremy L Thompson   Ceed         ceed;
156ccaff030SJeremy L Thompson   PetscScalar  loadIncrement;
157ccaff030SJeremy L Thompson };
158ccaff030SJeremy L Thompson 
159ccaff030SJeremy L Thompson // Data for Jacobian setup routine
160ccaff030SJeremy L Thompson typedef struct FormJacobCtx_private *FormJacobCtx;
161ccaff030SJeremy L Thompson struct FormJacobCtx_private {
162ccaff030SJeremy L Thompson   UserMult     *jacobCtx;
163ccaff030SJeremy L Thompson   PetscInt     numLevels;
164ccaff030SJeremy L Thompson   SNES         snesCoarse;
165ccaff030SJeremy L Thompson   Mat          *jacobMat, jacobMatCoarse;
166ccaff030SJeremy L Thompson   Vec          Ucoarse;
167ccaff030SJeremy L Thompson };
168ccaff030SJeremy L Thompson 
169ccaff030SJeremy L Thompson // Data for PETSc Prolongation/Restriction Matshell
170ccaff030SJeremy L Thompson typedef struct UserMultProlongRestr_private *UserMultProlongRestr;
171ccaff030SJeremy L Thompson struct UserMultProlongRestr_private {
172ccaff030SJeremy L Thompson   MPI_Comm     comm;
173ccaff030SJeremy L Thompson   DM           dmC, dmF;
174ccaff030SJeremy L Thompson   Vec          locVecC, locVecF, multVec;
175ccaff030SJeremy L Thompson   CeedVector   ceedVecC, ceedVecF;
176ccaff030SJeremy L Thompson   CeedOperator opProlong, opRestrict;
177ccaff030SJeremy L Thompson   Ceed         ceed;
178ccaff030SJeremy L Thompson };
179ccaff030SJeremy L Thompson 
180ccaff030SJeremy L Thompson // libCEED data struct for level
181ccaff030SJeremy L Thompson typedef struct CeedData_private *CeedData;
182ccaff030SJeremy L Thompson struct CeedData_private {
183ccaff030SJeremy L Thompson   Ceed                ceed;
184*5c25879aSJeremy L Thompson   CeedBasis           basisx, basisu, basisCtoF, basisEnergy, basisDiagnostic;
1852d93065eSjeremylt   CeedElemRestriction Erestrictx, Erestrictu, Erestrictqdi,
186*5c25879aSJeremy L Thompson                       ErestrictGradui, ErestrictEnergy, ErestrictDiagnostic,
187*5c25879aSJeremy L Thompson                       ErestrictqdDiagnostici;
188*5c25879aSJeremy L Thompson   CeedQFunction       qfApply, qfJacob, qfEnergy, qfDiagnostic;
189*5c25879aSJeremy L Thompson   CeedOperator        opApply, opJacob, opRestrict, opProlong, opEnergy,
190*5c25879aSJeremy L Thompson                       opDiagnostic;
191*5c25879aSJeremy L Thompson   CeedVector          qdata, qdataDiagnostic, gradu, xceed, yceed, truesoln;
192ccaff030SJeremy L Thompson };
193ccaff030SJeremy L Thompson 
194ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
195ccaff030SJeremy L Thompson // Process command line options
196ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
197ccaff030SJeremy L Thompson // Process general command line options
198ccaff030SJeremy L Thompson PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx appCtx);
199ccaff030SJeremy L Thompson 
200ccaff030SJeremy L Thompson // Process physics options
201ccaff030SJeremy L Thompson PetscErrorCode ProcessPhysics(MPI_Comm comm, Physics phys, Units units);
202ccaff030SJeremy L Thompson 
203ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
204ccaff030SJeremy L Thompson // Setup DM
205ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
206ccaff030SJeremy L Thompson PetscErrorCode CreateBCLabel(DM dm, const char name[]);
207ccaff030SJeremy L Thompson 
208ccaff030SJeremy L Thompson // Create FE by degree
209ccaff030SJeremy L Thompson PetscErrorCode PetscFECreateByDegree(DM dm, PetscInt dim, PetscInt Nc,
210ccaff030SJeremy L Thompson                                      PetscBool isSimplex, const char prefix[],
211ccaff030SJeremy L Thompson                                      PetscInt order, PetscFE *fem);
212ccaff030SJeremy L Thompson 
213ccaff030SJeremy L Thompson // Read mesh and distribute DM in parallel
214ccaff030SJeremy L Thompson PetscErrorCode CreateDistributedDM(MPI_Comm comm, AppCtx appCtx, DM *dm);
215ccaff030SJeremy L Thompson 
216ccaff030SJeremy L Thompson // Setup DM with FE space of appropriate degree
217ccaff030SJeremy L Thompson PetscErrorCode SetupDMByDegree(DM dm, AppCtx appCtx, PetscInt order,
218a3c02c40SJeremy L Thompson                                PetscBool boundary, PetscInt ncompu);
219ccaff030SJeremy L Thompson 
220ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
221ccaff030SJeremy L Thompson // libCEED Functions
222ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
223ccaff030SJeremy L Thompson // Destroy libCEED objects
224ccaff030SJeremy L Thompson PetscErrorCode CeedDataDestroy(CeedInt level, CeedData data);
225ccaff030SJeremy L Thompson 
226ccaff030SJeremy L Thompson // Get libCEED restriction data from DMPlex
227d979a051Sjeremylt PetscErrorCode CreateRestrictionPlex(Ceed ceed, CeedInt P, CeedInt ncomp,
228ccaff030SJeremy L Thompson                                      CeedElemRestriction *Erestrict, DM dm);
229ccaff030SJeremy L Thompson 
230ccaff030SJeremy L Thompson // Set up libCEED for a given degree
231a3c02c40SJeremy L Thompson PetscErrorCode SetupLibceedFineLevel(DM dm, DM dmEnergy, DM dmDiagnostic,
232a3c02c40SJeremy L Thompson                                      Ceed ceed, AppCtx appCtx, Physics phys,
233a3c02c40SJeremy L Thompson                                      CeedData *data, PetscInt fineLevel,
234a3c02c40SJeremy L Thompson                                      PetscInt ncompu, PetscInt Ugsz,
235a3c02c40SJeremy L Thompson                                      PetscInt Ulocsz, CeedVector forceCeed,
236ccaff030SJeremy L Thompson                                      CeedQFunction qfRestrict,
237ccaff030SJeremy L Thompson                                      CeedQFunction qfProlong);
238ccaff030SJeremy L Thompson 
239ccaff030SJeremy L Thompson // Set up libCEED for a given degree
240ccaff030SJeremy L Thompson PetscErrorCode SetupLibceedLevel(DM dm, Ceed ceed, AppCtx appCtx, Physics phys,
241ccaff030SJeremy L Thompson                                  CeedData *data, PetscInt level,
242ccaff030SJeremy L Thompson                                  PetscInt ncompu, PetscInt Ugsz,
243ccaff030SJeremy L Thompson                                  PetscInt Ulocsz, CeedVector forceCeed,
244ccaff030SJeremy L Thompson                                  CeedQFunction qfRestrict,
245ccaff030SJeremy L Thompson                                  CeedQFunction qfProlong);
246ccaff030SJeremy L Thompson 
247ccaff030SJeremy L Thompson // Setup context data for Jacobian evaluation
248ccaff030SJeremy L Thompson PetscErrorCode SetupJacobianCtx(MPI_Comm comm, AppCtx appCtx, DM dm, Vec V,
249ccaff030SJeremy L Thompson                                 Vec Vloc, CeedData ceedData, Ceed ceed,
250ccaff030SJeremy L Thompson                                 UserMult jacobianCtx);
251ccaff030SJeremy L Thompson 
252ccaff030SJeremy L Thompson // Setup context data for prolongation and restriction operators
253ccaff030SJeremy L Thompson PetscErrorCode SetupProlongRestrictCtx(MPI_Comm comm, DM dmC, DM dmF, Vec VF,
254ccaff030SJeremy L Thompson                                        Vec VlocC, Vec VlocF, CeedData ceedDataC,
255ccaff030SJeremy L Thompson                                        CeedData ceedDataF, Ceed ceed,
256ccaff030SJeremy L Thompson                                        UserMultProlongRestr prolongRestrCtx);
257ccaff030SJeremy L Thompson 
258ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
259ccaff030SJeremy L Thompson // Jacobian setup
260ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
261ccaff030SJeremy L Thompson PetscErrorCode FormJacobian(SNES snes, Vec U, Mat J, Mat Jpre, void *ctx);
262ccaff030SJeremy L Thompson 
263ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
264*5c25879aSJeremy L Thompson // Solution output
265ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
266ccaff030SJeremy L Thompson PetscErrorCode ViewSolution(MPI_Comm comm, Vec U, PetscInt increment,
267ccaff030SJeremy L Thompson                             PetscScalar loadIncrement);
268ccaff030SJeremy L Thompson 
269*5c25879aSJeremy L Thompson PetscErrorCode ViewDiagnosticQuantities(MPI_Comm comm, DM dmU,
270*5c25879aSJeremy L Thompson                                         UserMult user, Vec U,
271*5c25879aSJeremy L Thompson                                         CeedElemRestriction ErestrictDiagnostic);
272*5c25879aSJeremy L Thompson 
273ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
274ccaff030SJeremy L Thompson // libCEED Operators for MatShell
275ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
276ccaff030SJeremy L Thompson // This function uses libCEED to compute the local action of an operator
277ccaff030SJeremy L Thompson PetscErrorCode ApplyLocalCeedOp(Vec X, Vec Y, UserMult user);
278ccaff030SJeremy L Thompson 
279ccaff030SJeremy L Thompson // This function uses libCEED to compute the non-linear residual
280ccaff030SJeremy L Thompson PetscErrorCode FormResidual_Ceed(SNES snes, Vec X, Vec Y, void *ctx);
281ccaff030SJeremy L Thompson 
282ccaff030SJeremy L Thompson // This function uses libCEED to apply the Jacobian for assembly via a SNES
283ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobianCoarse_Ceed(SNES snes, Vec X, Vec Y, void *ctx);
284ccaff030SJeremy L Thompson 
285ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the Jacobian
286ccaff030SJeremy L Thompson PetscErrorCode ApplyJacobian_Ceed(Mat A, Vec X, Vec Y);
287ccaff030SJeremy L Thompson 
288ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the prolongation operator
289ccaff030SJeremy L Thompson PetscErrorCode Prolong_Ceed(Mat A, Vec X, Vec Y);
290ccaff030SJeremy L Thompson 
291ccaff030SJeremy L Thompson // This function uses libCEED to compute the action of the restriction operator
292ccaff030SJeremy L Thompson PetscErrorCode Restrict_Ceed(Mat A, Vec X, Vec Y);
2932d93065eSjeremylt 
294ccaff030SJeremy L Thompson // This function returns the computed diagonal of the operator
295ccaff030SJeremy L Thompson PetscErrorCode GetDiag_Ceed(Mat A, Vec D);
296ccaff030SJeremy L Thompson 
2972d93065eSjeremylt // This function calculates the strain energy in the final solution
298a3c02c40SJeremy L Thompson PetscErrorCode ComputeStrainEnergy(DM dmEnergy, UserMult user,
299a3c02c40SJeremy L Thompson                                    CeedOperator opEnergy, Vec X,
300a3c02c40SJeremy L Thompson                                    PetscReal *energy);
3012d93065eSjeremylt 
302ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
303ccaff030SJeremy L Thompson // Boundary Functions
304ccaff030SJeremy L Thompson // -----------------------------------------------------------------------------
305ccaff030SJeremy L Thompson // Note: If additional boundary conditions are added, an update is needed in
306ccaff030SJeremy L Thompson //         elasticity.h for the boundaryOptions variable.
307ccaff030SJeremy L Thompson 
308ccaff030SJeremy L Thompson // BCMMS - boundary function
309ccaff030SJeremy L Thompson // Values on all points of the mesh is set based on given solution below
310ccaff030SJeremy L Thompson // for u[0], u[1], u[2]
311ccaff030SJeremy L Thompson PetscErrorCode BCMMS(PetscInt dim, PetscReal loadIncrement,
312ccaff030SJeremy L Thompson                      const PetscReal coords[], PetscInt ncompu,
313ccaff030SJeremy L Thompson                      PetscScalar *u, void *ctx);
314ccaff030SJeremy L Thompson 
315d642641fSjeremylt // BCClamp - fix boundary values with affine transformation at fraction of load
316d642641fSjeremylt //   increment
317d642641fSjeremylt PetscErrorCode BCClamp(PetscInt dim, PetscReal loadIncrement,
318ccaff030SJeremy L Thompson                        const PetscReal coords[], PetscInt ncompu,
319ccaff030SJeremy L Thompson                        PetscScalar *u, void *ctx);
320ccaff030SJeremy L Thompson 
321ccaff030SJeremy L Thompson #endif //setup_h
322