Actual source code: petscfeimpl.h

petsc-main 2021-04-20
Report Typos and Errors
  1: #if !defined(PETSCFEIMPL_H)
  2: #define PETSCFEIMPL_H

  4: #include <petscfe.h>
  5: #include <petscds.h>
  6: #include <petsc/private/petscimpl.h>
  7: #include <petsc/private/dmpleximpl.h>

  9: PETSC_EXTERN PetscBool PetscSpaceRegisterAllCalled;
 10: PETSC_EXTERN PetscBool PetscDualSpaceRegisterAllCalled;
 11: PETSC_EXTERN PetscBool PetscFERegisterAllCalled;
 12: PETSC_EXTERN PetscErrorCode PetscSpaceRegisterAll(void);
 13: PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterAll(void);
 14: PETSC_EXTERN PetscErrorCode PetscFERegisterAll(void);

 16: PETSC_EXTERN PetscBool FEcite;
 17: PETSC_EXTERN const char FECitation[];

 19: PETSC_EXTERN PetscLogEvent PETSCDUALSPACE_SetUp, PETSCFE_SetUp;

 21: typedef struct _PetscSpaceOps *PetscSpaceOps;
 22: struct _PetscSpaceOps {
 23:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscSpace);
 24:   PetscErrorCode (*setup)(PetscSpace);
 25:   PetscErrorCode (*view)(PetscSpace,PetscViewer);
 26:   PetscErrorCode (*destroy)(PetscSpace);

 28:   PetscErrorCode (*getdimension)(PetscSpace,PetscInt*);
 29:   PetscErrorCode (*evaluate)(PetscSpace,PetscInt,const PetscReal*,PetscReal*,PetscReal*,PetscReal*);
 30:   PetscErrorCode (*getheightsubspace)(PetscSpace,PetscInt,PetscSpace *);
 31: };

 33: struct _p_PetscSpace {
 34:   PETSCHEADER(struct _PetscSpaceOps);
 35:   void                   *data;          /* Implementation object */
 36:   PetscInt                degree;        /* The approximation order of the space */
 37:   PetscInt                maxDegree;     /* The containing approximation order of the space */
 38:   PetscInt                Nc;            /* The number of components */
 39:   PetscInt                Nv;            /* The number of variables in the space, e.g. x and y */
 40:   PetscInt                dim;           /* The dimension of the space */
 41:   DM                      dm;            /* Shell to use for temp allocation */
 42: };

 44: typedef struct {
 45:   PetscBool                symmetric;   /* Use only symmetric polynomials */
 46:   PetscBool                tensor;      /* Flag for tensor product */
 47:   PetscInt                *degrees;     /* Degrees of single variable which we need to compute */
 48:   PetscSpacePolynomialType ptype;       /* Allows us to make the Hdiv and Hcurl spaces */
 49:   PetscBool                setupCalled;
 50:   PetscSpace              *subspaces;   /* Subspaces for each dimension */
 51: } PetscSpace_Poly;

 53: typedef struct {
 54:   PetscSpace *tensspaces;
 55:   PetscInt    numTensSpaces;
 56:   PetscInt    dim;
 57:   PetscBool   uniform;
 58:   PetscBool   setupCalled;
 59:   PetscSpace *heightsubspaces;    /* Height subspaces */
 60: } PetscSpace_Tensor;

 62: typedef struct {
 63:   PetscSpace *sumspaces;
 64:   PetscInt    numSumSpaces;
 65:   PetscBool   concatenate;
 66:   PetscBool   setupCalled;
 67: } PetscSpace_Sum;

 69: typedef struct {
 70:   PetscQuadrature quad;         /* The points defining the space */
 71: } PetscSpace_Point;

 73: typedef struct _PetscDualSpaceOps *PetscDualSpaceOps;
 74: struct _PetscDualSpaceOps {
 75:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscDualSpace);
 76:   PetscErrorCode (*setup)(PetscDualSpace);
 77:   PetscErrorCode (*view)(PetscDualSpace,PetscViewer);
 78:   PetscErrorCode (*destroy)(PetscDualSpace);

 80:   PetscErrorCode (*duplicate)(PetscDualSpace,PetscDualSpace);
 81:   PetscErrorCode (*createheightsubspace)(PetscDualSpace,PetscInt,PetscDualSpace *);
 82:   PetscErrorCode (*createpointsubspace)(PetscDualSpace,PetscInt,PetscDualSpace *);
 83:   PetscErrorCode (*getsymmetries)(PetscDualSpace,const PetscInt****,const PetscScalar****);
 84:   PetscErrorCode (*apply)(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
 85:   PetscErrorCode (*applyall)(PetscDualSpace, const PetscScalar *, PetscScalar *);
 86:   PetscErrorCode (*applyint)(PetscDualSpace, const PetscScalar *, PetscScalar *);
 87:   PetscErrorCode (*createalldata)(PetscDualSpace, PetscQuadrature *, Mat *);
 88:   PetscErrorCode (*createintdata)(PetscDualSpace, PetscQuadrature *, Mat *);
 89: };

 91: struct _p_PetscDualSpace {
 92:   PETSCHEADER(struct _PetscDualSpaceOps);
 93:   void            *data;       /* Implementation object */
 94:   DM               dm;         /* The integration region K */
 95:   PetscInt         order;      /* The approximation order of the space */
 96:   PetscInt         Nc;         /* The number of components */
 97:   PetscQuadrature *functional; /* The basis of functionals for this space */
 98:   Mat              allMat;
 99:   PetscQuadrature  allNodes;   /* Collects all quadrature points representing functionals in the basis */
100:   Vec              allNodeValues;
101:   Vec              allDofValues;
102:   Mat              intMat;
103:   PetscQuadrature  intNodes;   /* Collects all quadrature points representing functionals in the basis in the interior of the cell */
104:   Vec              intNodeValues;
105:   Vec              intDofValues;
106:   PetscInt         spdim;      /* The dual-space dimension */
107:   PetscInt         spintdim;   /* The dual-space interior dimension */
108:   PetscInt         k;          /* k-simplex corresponding to the dofs in this basis (we always use the 3D complex right now) */
109:   PetscBool        uniform;
110:   PetscBool        setupcalled;
111:   PetscBool        setfromoptionscalled;
112:   PetscSection     pointSection;
113:   PetscDualSpace  *pointSpaces;
114:   PetscDualSpace  *heightSpaces;
115:   PetscInt        *numDof;
116: };

118: typedef struct _n_Petsc1DNodeFamily *Petsc1DNodeFamily;
119: typedef struct _n_PetscLagNodeIndices *PetscLagNodeIndices;

121: PETSC_EXTERN PetscErrorCode PetscLagNodeIndicesGetData_Internal(PetscLagNodeIndices, PetscInt *, PetscInt *, PetscInt *, const PetscInt *[], const PetscReal *[]);
122: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(PetscDualSpace sp, PetscInt ornt, Mat *symMat);

124: typedef struct {
125:   /* these describe the types of dual spaces implemented */
126:   PetscBool         tensorCell;  /* Flag for tensor product cell */
127:   PetscBool         tensorSpace; /* Flag for tensor product space of polynomials, as opposed to a space of maximum degree */
128:   PetscBool         trimmed;     /* Flag for dual space of trimmed polynomial spaces */
129:   PetscBool         continuous;  /* Flag for a continuous basis, as opposed to discontinuous across element boundaries */

131:   PetscBool         interiorOnly; /* To make setup faster for tensor elements, only construct interior dofs in recursive calls */

133:   /* these keep track of symmetries */
134:   PetscInt       ***symperms;
135:   PetscScalar    ***symflips;
136:   PetscInt          numSelfSym;
137:   PetscInt          selfSymOff;
138:   PetscBool         symComputed;

140:   /* these describe different schemes of placing nodes in a simplex, from
141:    * which are derived all dofs in Lagrange dual spaces */
142:   PetscDTNodeType   nodeType;
143:   PetscBool         endNodes;
144:   PetscReal         nodeExponent;
145:   PetscInt          numNodeSkip; /* The number of end nodes from the 1D Node family to skip */
146:   Petsc1DNodeFamily nodeFamily;

148:   PetscInt          numCopies;

150:   PetscBool         useMoments;  /* Use moments for functionals */
151:   PetscInt          momentOrder; /* Order for moment quadrature */

153:   /* these are ways of indexing nodes in a way that makes
154:    * the computation of symmetries programmatic */
155:   PetscLagNodeIndices vertIndices;
156:   PetscLagNodeIndices intNodeIndices;
157:   PetscLagNodeIndices allNodeIndices;
158: } PetscDualSpace_Lag;

160: typedef struct {
161:   PetscInt  dim;
162:   PetscInt *numDof;
163: } PetscDualSpace_Simple;

165: typedef struct _PetscFEOps *PetscFEOps;
166: struct _PetscFEOps {
167:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,PetscFE);
168:   PetscErrorCode (*setup)(PetscFE);
169:   PetscErrorCode (*view)(PetscFE,PetscViewer);
170:   PetscErrorCode (*destroy)(PetscFE);
171:   PetscErrorCode (*getdimension)(PetscFE,PetscInt*);
172:   PetscErrorCode (*createtabulation)(PetscFE,PetscInt,const PetscReal*,PetscInt,PetscTabulation);
173:   /* Element integration */
174:   PetscErrorCode (*integrate)(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
175:   PetscErrorCode (*integratebd)(PetscDS, PetscInt, PetscBdPointFunc, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
176:   PetscErrorCode (*integrateresidual)(PetscDS, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
177:   PetscErrorCode (*integratebdresidual)(PetscDS, PetscWeakForm, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
178:   PetscErrorCode (*integratehybridresidual)(PetscDS, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
179:   PetscErrorCode (*integratejacobianaction)(PetscFE, PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
180:   PetscErrorCode (*integratejacobian)(PetscDS, PetscFEJacobianType, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
181:   PetscErrorCode (*integratebdjacobian)(PetscDS, PetscWeakForm, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
182:   PetscErrorCode (*integratehybridjacobian)(PetscDS, PetscFEJacobianType, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
183: };

185: struct _p_PetscFE {
186:   PETSCHEADER(struct _PetscFEOps);
187:   void           *data;                  /* Implementation object */
188:   PetscSpace      basisSpace;            /* The basis space P */
189:   PetscDualSpace  dualSpace;             /* The dual space P' */
190:   PetscInt        numComponents;         /* The number of field components */
191:   PetscQuadrature quadrature;            /* Suitable quadrature on K */
192:   PetscQuadrature faceQuadrature;        /* Suitable face quadrature on \partial K */
193:   PetscFE        *subspaces;             /* Subspaces for each dimension */
194:   PetscReal      *invV;                  /* Change of basis matrix, from prime to nodal basis set */
195:   PetscTabulation T;                     /* Tabulation of basis and derivatives at quadrature points */
196:   PetscTabulation Tf;                    /* Tabulation of basis and derivatives at quadrature points on each face */
197:   PetscTabulation Tc;                    /* Tabulation of basis at face centroids */
198:   PetscInt        blockSize, numBlocks;  /* Blocks are processed concurrently */
199:   PetscInt        batchSize, numBatches; /* A batch is made up of blocks, Batches are processed in serial */
200:   PetscBool       setupcalled;
201: };

203: typedef struct {
204:   PetscInt cellType;
205: } PetscFE_Basic;

207: #ifdef PETSC_HAVE_OPENCL

209: #ifdef __APPLE__
210: #include <OpenCL/cl.h>
211: #else
212: #include <CL/cl.h>
213: #endif

215: typedef struct {
216:   cl_platform_id   pf_id;
217:   cl_device_id     dev_id;
218:   cl_context       ctx_id;
219:   cl_command_queue queue_id;
220:   PetscDataType    realType;
221:   PetscLogEvent    residualEvent;
222:   PetscInt         op; /* ANDY: Stand-in for real equation code generation */
223: } PetscFE_OpenCL;
224: #endif

226: typedef struct {
227:   PetscInt   numSubelements; /* The number of subelements */
228:   PetscReal *v0;             /* The affine transformation for each subelement */
229:   PetscReal *jac, *invjac;
230:   PetscInt  *embedding;      /* Map from subelements dofs to element dofs */
231: } PetscFE_Composite;

233: /* Utility functions */
234: PETSC_STATIC_INLINE void CoordinatesRefToReal(PetscInt dimReal, PetscInt dimRef, const PetscReal xi0[], const PetscReal v0[], const PetscReal J[], const PetscReal xi[], PetscReal x[])
235: {
236:   PetscInt d, e;

238:   for (d = 0; d < dimReal; ++d) {
239:     x[d] = v0[d];
240:     for (e = 0; e < dimRef; ++e) {
241:       x[d] += J[d*dimReal+e]*(xi[e] - xi0[e]);
242:     }
243:   }
244: }

246: PETSC_STATIC_INLINE void CoordinatesRealToRef(PetscInt dimReal, PetscInt dimRef, const PetscReal xi0[], const PetscReal v0[], const PetscReal invJ[], const PetscReal x[], PetscReal xi[])
247: {
248:   PetscInt d, e;

250:   for (d = 0; d < dimRef; ++d) {
251:     xi[d] = xi0[d];
252:     for (e = 0; e < dimReal; ++e) {
253:       xi[d] += invJ[d*dimReal+e]*(x[e] - v0[e]);
254:     }
255:   }
256: }

258: PETSC_STATIC_INLINE PetscErrorCode PetscFEInterpolate_Static(PetscFE fe, const PetscScalar x[], PetscFEGeom *fegeom, PetscInt q, PetscScalar interpolant[])
259: {
260:   PetscTabulation T;
261:   PetscInt        fc, f;
262:   PetscErrorCode  ierr;

265:   PetscFEGetCellTabulation(fe, 0, &T);
266:   {
267:     const PetscReal *basis = T->T[0];
268:     const PetscInt   Nb    = T->Nb;
269:     const PetscInt   Nc    = T->Nc;
270:     for (fc = 0; fc < Nc; ++fc) {
271:       interpolant[fc] = 0.0;
272:       for (f = 0; f < Nb; ++f) {
273:         interpolant[fc] += x[f]*basis[(q*Nb + f)*Nc + fc];
274:       }
275:     }
276:   }
277:   PetscFEPushforward(fe, fegeom, 1, interpolant);
278:   return(0);
279: }

281: PETSC_STATIC_INLINE PetscErrorCode PetscFEInterpolateGradient_Static(PetscFE fe, PetscInt k, const PetscScalar x[], PetscFEGeom *fegeom, PetscInt q, PetscScalar interpolant[])
282: {
283:   PetscTabulation T;
284:   PetscInt        fc, f, d;
285:   PetscErrorCode  ierr;

288:   PetscFEGetCellTabulation(fe, k, &T);
289:   {
290:     const PetscReal *basisDer = T->T[1];
291:     const PetscReal *basisHes = k > 1 ? T->T[2] : NULL;
292:     const PetscInt   Nb       = T->Nb;
293:     const PetscInt   Nc       = T->Nc;
294:     const PetscInt   cdim     = T->cdim;

296:     if (cdim != fegeom->dimEmbed) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Geometry dim %D must match tabulation dim %D", fegeom->dimEmbed, cdim);
297:     for (fc = 0; fc < Nc; ++fc) {
298:       for (d = 0; d < cdim; ++d) interpolant[fc*cdim+d] = 0.0;
299:       for (f = 0; f < Nb; ++f) {
300:         for (d = 0; d < cdim; ++d) {
301:           interpolant[fc*cdim+d] += x[f]*basisDer[((q*Nb + f)*Nc + fc)*cdim + d];
302:         }
303:       }
304:     }
305:     if (k > 1) {
306:       const PetscInt off = Nc*cdim;

308:       for (fc = 0; fc < Nc; ++fc) {
309:         for (d = 0; d < cdim*cdim; ++d) interpolant[off+fc*cdim*cdim+d] = 0.0;
310:         for (f = 0; f < Nb; ++f) {
311:           for (d = 0; d < cdim*cdim; ++d) interpolant[off+fc*cdim+d] += x[f]*basisHes[((q*Nb + f)*Nc + fc)*cdim*cdim + d];
312:         }
313:       }
314:     }
315:   }
316:   PetscFEPushforwardGradient(fe, fegeom, 1, interpolant);
317:   return(0);
318: }

320: PETSC_STATIC_INLINE PetscErrorCode PetscFEInterpolateFieldAndGradient_Static(PetscFE fe, PetscInt k, const PetscScalar x[], PetscFEGeom *fegeom, PetscInt q, PetscScalar interpolant[], PetscScalar interpolantGrad[])
321: {
322:   PetscTabulation T;
323:   PetscInt        fc, f, d;
324:   PetscErrorCode  ierr;

327:   PetscFEGetCellTabulation(fe, k, &T);
328:   {
329:     const PetscReal *basis    = T->T[0];
330:     const PetscReal *basisDer = T->T[1];
331:     const PetscReal *basisHes = k > 1 ? T->T[2] : NULL;
332:     const PetscInt   Nb       = T->Nb;
333:     const PetscInt   Nc       = T->Nc;
334:     const PetscInt   cdim     = T->cdim;

336:     if (cdim != fegeom->dimEmbed) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Geometry dim %D must match tabulation dim %D", fegeom->dimEmbed, cdim);
337:     for (fc = 0; fc < Nc; ++fc) {
338:       interpolant[fc] = 0.0;
339:       for (d = 0; d < cdim; ++d) interpolantGrad[fc*cdim+d] = 0.0;
340:       for (f = 0; f < Nb; ++f) {
341:         interpolant[fc] += x[f]*basis[(q*Nb + f)*Nc + fc];
342:         for (d = 0; d < cdim; ++d) interpolantGrad[fc*cdim+d] += x[f]*basisDer[((q*Nb + f)*Nc + fc)*cdim + d];
343:       }
344:     }
345:     if (k > 1) {
346:       const PetscInt off = Nc*cdim;

348:       for (fc = 0; fc < Nc; ++fc) {
349:         for (d = 0; d < cdim*cdim; ++d) interpolantGrad[off+fc*cdim*cdim+d] = 0.0;
350:         for (f = 0; f < Nb; ++f) {
351:           for (d = 0; d < cdim*cdim; ++d) interpolantGrad[off+fc*cdim+d] += x[f]*basisHes[((q*Nb + f)*Nc + fc)*cdim*cdim + d];
352:         }
353:       }
354:       PetscFEPushforwardHessian(fe, fegeom, 1, &interpolantGrad[off]);
355:     }
356:   }
357:   PetscFEPushforward(fe, fegeom, 1, interpolant);
358:   PetscFEPushforwardGradient(fe, fegeom, 1, interpolantGrad);
359:   return(0);
360: }

362: PETSC_INTERN PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt, PetscInt, PetscInt[]);
363: PETSC_INTERN PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt, PetscInt, PetscInt[]);

365: PETSC_INTERN PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace, PetscSection*);
366: PETSC_INTERN PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace, PetscSection);
367: PETSC_INTERN PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace, PetscInt, PetscInt);

369: PETSC_INTERN PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS, PetscInt, PetscInt, PetscInt, PetscTabulation[], PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscScalar[], PetscScalar[], PetscScalar[]);
370: PETSC_INTERN PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
371: PETSC_INTERN PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE, PetscTabulation, PetscInt, PetscScalar[], PetscScalar[], PetscFEGeom *, PetscScalar[], PetscScalar[], PetscScalar[]);
372: PETSC_INTERN PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE, PetscFE, PetscInt, PetscInt, PetscTabulation, PetscScalar[], PetscScalar[], PetscTabulation, PetscScalar[], PetscScalar[], PetscFEGeom *, const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[]);

374: PETSC_INTERN PetscErrorCode PetscFEEvaluateFieldJets_Hybrid_Internal(PetscDS, PetscInt, PetscInt, PetscInt, PetscTabulation[], PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscScalar[], PetscScalar[], PetscScalar[]);
375: PETSC_INTERN PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE, PetscTabulation, PetscInt, PetscScalar[], PetscScalar[], PetscFEGeom *, PetscScalar[], PetscScalar[], PetscScalar[]);
376: PETSC_INTERN PetscErrorCode PetscFEUpdateElementMat_Hybrid_Internal(PetscFE, PetscBool, PetscFE, PetscBool, PetscInt, PetscInt, PetscTabulation, PetscScalar[], PetscScalar[], PetscTabulation, PetscScalar[], PetscScalar[], PetscFEGeom *, const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[]);

378: PETSC_EXTERN PetscErrorCode PetscFEGetDimension_Basic(PetscFE, PetscInt *);
379: PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar [], const PetscScalar [], PetscDS, const PetscScalar [], PetscReal, PetscScalar []);
380: PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS, PetscWeakForm, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar [], const PetscScalar [], PetscDS, const PetscScalar [], PetscReal, PetscScalar[]);
381: PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS, PetscFEJacobianType, PetscHashFormKey, PetscInt, PetscFEGeom *, const PetscScalar [], const PetscScalar [], PetscDS, const PetscScalar [], PetscReal, PetscReal, PetscScalar []);
382: #endif