xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision ee12ae39415b2e672d944cdca066227dadbf8b14)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 #include <petscblaslapack.h>
3 
4 static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
5 {
6   PetscFE_Basic *b = (PetscFE_Basic *) fem->data;
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   ierr = PetscFree(b);CHKERRQ(ierr);
11   PetscFunctionReturn(0);
12 }
13 
14 static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
15 {
16   PetscInt          dim, Nc;
17   PetscSpace        basis = NULL;
18   PetscDualSpace    dual = NULL;
19   PetscQuadrature   quad = NULL;
20   PetscErrorCode    ierr;
21 
22   PetscFunctionBegin;
23   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
24   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
25   ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr);
26   ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr);
27   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
28   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
29   ierr = PetscViewerASCIIPrintf(v, "Basic Finite Element in %D dimensions with %D components\n",dim,Nc);CHKERRQ(ierr);
30   if (basis) {ierr = PetscSpaceView(basis, v);CHKERRQ(ierr);}
31   if (dual)  {ierr = PetscDualSpaceView(dual, v);CHKERRQ(ierr);}
32   if (quad)  {ierr = PetscQuadratureView(quad, v);CHKERRQ(ierr);}
33   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
38 {
39   PetscBool      iascii;
40   PetscErrorCode ierr;
41 
42   PetscFunctionBegin;
43   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
44   if (iascii) {ierr = PetscFEView_Basic_Ascii(fe, v);CHKERRQ(ierr);}
45   PetscFunctionReturn(0);
46 }
47 
48 /* Construct the change of basis from prime basis to nodal basis */
49 PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
50 {
51   PetscReal     *work;
52   PetscBLASInt  *pivots;
53   PetscBLASInt   n, info;
54   PetscInt       pdim, j;
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
59   ierr = PetscMalloc1(pdim*pdim,&fem->invV);CHKERRQ(ierr);
60   for (j = 0; j < pdim; ++j) {
61     PetscReal       *Bf;
62     PetscQuadrature  f;
63     const PetscReal *points, *weights;
64     PetscInt         Nc, Nq, q, k, c;
65 
66     ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr);
67     ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
68     ierr = PetscMalloc1(Nc*Nq*pdim,&Bf);CHKERRQ(ierr);
69     ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr);
70     for (k = 0; k < pdim; ++k) {
71       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
72       fem->invV[j*pdim+k] = 0.0;
73 
74       for (q = 0; q < Nq; ++q) {
75         for (c = 0; c < Nc; ++c) fem->invV[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c];
76       }
77     }
78     ierr = PetscFree(Bf);CHKERRQ(ierr);
79   }
80 
81   ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr);
82   n = pdim;
83   PetscStackCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
84   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
85   PetscStackCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
86   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
87   ierr = PetscFree2(pivots,work);CHKERRQ(ierr);
88   PetscFunctionReturn(0);
89 }
90 
91 PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
92 {
93   PetscErrorCode ierr;
94 
95   PetscFunctionBegin;
96   ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 /* Tensor contraction on the middle index,
101  *    C[m,n,p] = A[m,k,p] * B[k,n]
102  * where all matrices use C-style ordering.
103  */
104 static PetscErrorCode TensorContract_Private(PetscInt m,PetscInt n,PetscInt p,PetscInt k,const PetscReal *A,const PetscReal *B,PetscReal *C)
105 {
106   PetscErrorCode ierr;
107   PetscInt i;
108 
109   PetscFunctionBegin;
110   for (i=0; i<m; i++) {
111     PetscBLASInt n_,p_,k_,lda,ldb,ldc;
112     PetscReal one = 1, zero = 0;
113     /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
114      * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
115      */
116     ierr = PetscBLASIntCast(n,&n_);CHKERRQ(ierr);
117     ierr = PetscBLASIntCast(p,&p_);CHKERRQ(ierr);
118     ierr = PetscBLASIntCast(k,&k_);CHKERRQ(ierr);
119     lda = p_;
120     ldb = n_;
121     ldc = p_;
122     PetscStackCallBLAS("BLASgemm",BLASREALgemm_("N","T",&p_,&n_,&k_,&one,A+i*k*p,&lda,B,&ldb,&zero,C+i*n*p,&ldc));
123   }
124   ierr = PetscLogFlops(2.*m*n*p*k);CHKERRQ(ierr);
125   PetscFunctionReturn(0);
126 }
127 
128 PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
129 {
130   DM               dm;
131   PetscInt         pdim; /* Dimension of FE space P */
132   PetscInt         dim;  /* Spatial dimension */
133   PetscInt         Nc;   /* Field components */
134   PetscReal       *B = K >= 0 ? T->T[0] : NULL;
135   PetscReal       *D = K >= 1 ? T->T[1] : NULL;
136   PetscReal       *H = K >= 2 ? T->T[2] : NULL;
137   PetscReal       *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
138   PetscErrorCode   ierr;
139 
140   PetscFunctionBegin;
141   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
142   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
143   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
144   ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr);
145   /* Evaluate the prime basis functions at all points */
146   if (K >= 0) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
147   if (K >= 1) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
148   if (K >= 2) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
149   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH);CHKERRQ(ierr);
150   /* Translate from prime to nodal basis */
151   if (B) {
152     /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
153     ierr = TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B);CHKERRQ(ierr);
154   }
155   if (D) {
156     /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
157     ierr = TensorContract_Private(npoints, pdim, Nc*dim, pdim, tmpD, fem->invV, D);CHKERRQ(ierr);
158   }
159   if (H) {
160     /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
161     ierr = TensorContract_Private(npoints, pdim, Nc*dim*dim, pdim, tmpH, fem->invV, H);CHKERRQ(ierr);
162   }
163   if (K >= 0) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
164   if (K >= 1) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
165   if (K >= 2) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
166   PetscFunctionReturn(0);
167 }
168 
169 static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
170                                              const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
171 {
172   const PetscInt     debug = 0;
173   PetscFE            fe;
174   PetscPointFunc     obj_func;
175   PetscQuadrature    quad;
176   PetscTabulation   *T, *TAux = NULL;
177   PetscScalar       *u, *u_x, *a, *a_x;
178   const PetscScalar *constants;
179   PetscReal         *x;
180   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
181   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
182   PetscBool          isAffine;
183   const PetscReal   *quadPoints, *quadWeights;
184   PetscInt           qNc, Nq, q;
185   PetscErrorCode     ierr;
186 
187   PetscFunctionBegin;
188   ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr);
189   if (!obj_func) PetscFunctionReturn(0);
190   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
191   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
192   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
193   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
194   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
195   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
196   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
197   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
198   ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr);
199   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
200   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
201   if (dsAux) {
202     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
203     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
204     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
205     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
206     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
207     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
208     if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
209   }
210   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
211   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
212   Np = cgeom->numPoints;
213   dE = cgeom->dimEmbed;
214   isAffine = cgeom->isAffine;
215   for (e = 0; e < Ne; ++e) {
216     PetscFEGeom fegeom;
217 
218     fegeom.dim      = cgeom->dim;
219     fegeom.dimEmbed = cgeom->dimEmbed;
220     if (isAffine) {
221       fegeom.v    = x;
222       fegeom.xi   = cgeom->xi;
223       fegeom.J    = &cgeom->J[e*Np*dE*dE];
224       fegeom.invJ = &cgeom->invJ[e*Np*dE*dE];
225       fegeom.detJ = &cgeom->detJ[e*Np];
226     }
227     for (q = 0; q < Nq; ++q) {
228       PetscScalar integrand;
229       PetscReal   w;
230 
231       if (isAffine) {
232         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
233       } else {
234         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
235         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
236         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
237         fegeom.detJ = &cgeom->detJ[e*Np+q];
238       }
239       w = fegeom.detJ[0]*quadWeights[q];
240       if (debug > 1 && q < Np) {
241         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
242 #if !defined(PETSC_USE_COMPLEX)
243         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
244 #endif
245       }
246       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
247       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr);
248       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
249       obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, numConstants, constants, &integrand);
250       integrand *= w;
251       integral[e*Nf+field] += integrand;
252       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);}
253     }
254     cOffset    += totDim;
255     cOffsetAux += totDimAux;
256   }
257   PetscFunctionReturn(0);
258 }
259 
260 static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field,
261                                                PetscBdPointFunc obj_func,
262                                                PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
263 {
264   const PetscInt     debug = 0;
265   PetscFE            fe;
266   PetscQuadrature    quad;
267   PetscTabulation   *Tf, *TfAux = NULL;
268   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
269   const PetscScalar *constants;
270   PetscReal         *x;
271   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
272   PetscBool          isAffine, auxOnBd;
273   const PetscReal   *quadPoints, *quadWeights;
274   PetscInt           qNc, Nq, q, Np, dE;
275   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
276   PetscErrorCode     ierr;
277 
278   PetscFunctionBegin;
279   if (!obj_func) PetscFunctionReturn(0);
280   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
281   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
282   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
283   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
284   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
285   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
286   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
287   ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr);
288   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
289   ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr);
290   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
291   if (dsAux) {
292     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
293     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
294     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
295     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
296     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
297     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
298     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
299     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
300     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
301     if (Tf[0]->Np != TfAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
302   }
303   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
304   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
305   Np = fgeom->numPoints;
306   dE = fgeom->dimEmbed;
307   isAffine = fgeom->isAffine;
308   for (e = 0; e < Ne; ++e) {
309     PetscFEGeom    fegeom, cgeom;
310     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
311     fegeom.n = NULL;
312     fegeom.v = NULL;
313     fegeom.J = NULL;
314     fegeom.detJ = NULL;
315     fegeom.dim      = fgeom->dim;
316     fegeom.dimEmbed = fgeom->dimEmbed;
317     cgeom.dim       = fgeom->dim;
318     cgeom.dimEmbed  = fgeom->dimEmbed;
319     if (isAffine) {
320       fegeom.v    = x;
321       fegeom.xi   = fgeom->xi;
322       fegeom.J    = &fgeom->J[e*Np*dE*dE];
323       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
324       fegeom.detJ = &fgeom->detJ[e*Np];
325       fegeom.n    = &fgeom->n[e*Np*dE];
326 
327       cgeom.J     = &fgeom->suppJ[0][e*Np*dE*dE];
328       cgeom.invJ  = &fgeom->suppInvJ[0][e*Np*dE*dE];
329       cgeom.detJ  = &fgeom->suppDetJ[0][e*Np];
330     }
331     for (q = 0; q < Nq; ++q) {
332       PetscScalar integrand;
333       PetscReal   w;
334 
335       if (isAffine) {
336         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
337       } else {
338         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
339         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
340         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
341         fegeom.detJ = &fgeom->detJ[e*Np+q];
342         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
343 
344         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
345         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
346         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
347       }
348       w = fegeom.detJ[0]*quadWeights[q];
349       if (debug > 1 && q < Np) {
350         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
351 #ifndef PETSC_USE_COMPLEX
352         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
353 #endif
354       }
355       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
356       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr);
357       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
358       obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, fegeom.n, numConstants, constants, &integrand);
359       integrand *= w;
360       integral[e*Nf+field] += integrand;
361       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);}
362     }
363     cOffset    += totDim;
364     cOffsetAux += totDimAux;
365   }
366   PetscFunctionReturn(0);
367 }
368 
369 PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
370                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
371 {
372   const PetscInt     debug = 0;
373   const PetscInt     field = key.field;
374   PetscFE            fe;
375   PetscWeakForm      wf;
376   PetscInt           n0,      n1, i;
377   PetscPointFunc    *f0_func, *f1_func;
378   PetscQuadrature    quad;
379   PetscTabulation   *T, *TAux = NULL;
380   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
381   const PetscScalar *constants;
382   PetscReal         *x;
383   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
384   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
385   PetscBool          isAffine;
386   const PetscReal   *quadPoints, *quadWeights;
387   PetscInt           qNc, Nq, q, Np, dE;
388   PetscErrorCode     ierr;
389 
390   PetscFunctionBegin;
391   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
392   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
393   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
394   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
395   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
396   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
397   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
398   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
399   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
400   ierr = PetscWeakFormGetResidual(wf, key.label, key.value, key.field, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr);
401   if (!n0 && !n1) PetscFunctionReturn(0);
402   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
403   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
404   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
405   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
406   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
407   if (dsAux) {
408     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
409     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
410     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
411     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
412     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
413     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
414     if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
415   }
416   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
417   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
418   Np = cgeom->numPoints;
419   dE = cgeom->dimEmbed;
420   isAffine = cgeom->isAffine;
421   for (e = 0; e < Ne; ++e) {
422     PetscFEGeom fegeom;
423 
424     fegeom.dim      = cgeom->dim;
425     fegeom.dimEmbed = cgeom->dimEmbed;
426     if (isAffine) {
427       fegeom.v    = x;
428       fegeom.xi   = cgeom->xi;
429       fegeom.J    = &cgeom->J[e*Np*dE*dE];
430       fegeom.invJ = &cgeom->invJ[e*Np*dE*dE];
431       fegeom.detJ = &cgeom->detJ[e*Np];
432     }
433     ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr);
434     ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dE);CHKERRQ(ierr);
435     for (q = 0; q < Nq; ++q) {
436       PetscReal w;
437       PetscInt  c, d;
438 
439       if (isAffine) {
440         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
441       } else {
442         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
443         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
444         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
445         fegeom.detJ = &cgeom->detJ[e*Np+q];
446       }
447       w = fegeom.detJ[0]*quadWeights[q];
448       if (debug > 1 && q < Np) {
449         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
450 #if !defined(PETSC_USE_COMPLEX)
451         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
452 #endif
453       }
454       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
455       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
456       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
457       for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q*T[field]->Nc]);
458       for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w;
459       for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q*T[field]->Nc*dim]);
460       for (c = 0; c < T[field]->Nc; ++c) for (d = 0; d < dim; ++d) f1[(q*T[field]->Nc+c)*dim+d] *= w;
461     }
462     ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
463     cOffset    += totDim;
464     cOffsetAux += totDimAux;
465   }
466   PetscFunctionReturn(0);
467 }
468 
469 PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
470                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
471 {
472   const PetscInt     debug = 0;
473   PetscFE            fe;
474   PetscBdPointFunc   f0_func;
475   PetscBdPointFunc   f1_func;
476   PetscQuadrature    quad;
477   PetscTabulation   *Tf, *TfAux = NULL;
478   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
479   const PetscScalar *constants;
480   PetscReal         *x;
481   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
482   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
483   PetscBool          isAffine, auxOnBd = PETSC_FALSE;
484   const PetscReal   *quadPoints, *quadWeights;
485   PetscInt           qNc, Nq, q, Np, dE;
486   PetscErrorCode     ierr;
487 
488   PetscFunctionBegin;
489   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
490   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
491   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
492   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
493   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
494   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
495   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
496   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
497   ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr);
498   if (!f0_func && !f1_func) PetscFunctionReturn(0);
499   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
500   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
501   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
502   ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr);
503   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
504   if (dsAux) {
505     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
506     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
507     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
508     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
509     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
510     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
511     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
512     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
513     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
514     if (Tf[0]->Np != TfAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
515   }
516   NcI = Tf[field]->Nc;
517   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
518   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
519   Np = fgeom->numPoints;
520   dE = fgeom->dimEmbed;
521   isAffine = fgeom->isAffine;
522   for (e = 0; e < Ne; ++e) {
523     PetscFEGeom    fegeom, cgeom;
524     const PetscInt face = fgeom->face[e][0];
525     fegeom.n = NULL;
526     fegeom.v = NULL;
527     fegeom.J = NULL;
528     fegeom.detJ = NULL;
529     fegeom.dim      = fgeom->dim;
530     fegeom.dimEmbed = fgeom->dimEmbed;
531     cgeom.dim       = fgeom->dim;
532     cgeom.dimEmbed  = fgeom->dimEmbed;
533     if (isAffine) {
534       fegeom.v    = x;
535       fegeom.xi   = fgeom->xi;
536       fegeom.J    = &fgeom->J[e*Np*dE*dE];
537       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
538       fegeom.detJ = &fgeom->detJ[e*Np];
539       fegeom.n    = &fgeom->n[e*Np*dE];
540 
541       cgeom.J     = &fgeom->suppJ[0][e*Np*dE*dE];
542       cgeom.invJ  = &fgeom->suppInvJ[0][e*Np*dE*dE];
543       cgeom.detJ  = &fgeom->suppDetJ[0][e*Np];
544     }
545     ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr);
546     ierr = PetscArrayzero(f1, Nq*NcI*dE);CHKERRQ(ierr);
547     for (q = 0; q < Nq; ++q) {
548       PetscReal w;
549       PetscInt  c, d;
550 
551       if (isAffine) {
552         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
553       } else {
554         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
555         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
556         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
557         fegeom.detJ = &fgeom->detJ[e*Np+q];
558         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
559 
560         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
561         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
562         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
563       }
564       w = fegeom.detJ[0]*quadWeights[q];
565       if (debug > 1 && q < Np) {
566         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
567 #if !defined(PETSC_USE_COMPLEX)
568         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
569 #endif
570       }
571       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
572       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
573       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
574       if (f0_func) {
575         f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q*NcI]);
576         for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w;
577       }
578       if (f1_func) {
579         f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q*NcI*dim]);
580         for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w;
581       }
582     }
583     ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, &cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
584     cOffset    += totDim;
585     cOffsetAux += totDimAux;
586   }
587   PetscFunctionReturn(0);
588 }
589 
590 /*
591   BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but
592               all transforms operate in the full space and are square.
593 
594   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
595     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
596     2) We need to assume that the orientation is 0 for both
597     3) TODO We need to use a non-square Jacobian for the derivative maps, meaning the embedding dimension has to go to EvaluateFieldJets() and UpdateElementVec()
598 */
599 static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
600                                                            const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
601 {
602   const PetscInt     debug = 0;
603   const PetscInt     field = key.field;
604   PetscFE            fe;
605   PetscWeakForm      wf;
606   PetscInt           n0,      n1, i;
607   PetscBdPointFunc  *f0_func, *f1_func;
608   PetscQuadrature    quad;
609   PetscTabulation   *Tf, *TfAux = NULL;
610   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
611   const PetscScalar *constants;
612   PetscReal         *x;
613   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
614   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
615   PetscBool          isCohesiveField, isAffine, auxOnBd = PETSC_FALSE;
616   const PetscReal   *quadPoints, *quadWeights;
617   PetscInt           qNc, Nq, q, Np, dE;
618   PetscErrorCode     ierr;
619 
620   PetscFunctionBegin;
621   /* Hybrid discretization is posed directly on faces */
622   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
623   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
624   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
625   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
626   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
627   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
628   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
629   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
630   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
631   ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr);
632   if (!n0 && !n1) PetscFunctionReturn(0);
633   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
634   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
635   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
636   /* NOTE This is a bulk tabulation because the DS is a face discretization */
637   ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr);
638   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
639   if (dsAux) {
640     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
641     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
642     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
643     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
644     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
645     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
646     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
647     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
648     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
649     if (Tf[0]->Np != TfAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
650   }
651   isCohesiveField = field == Nf-1 ? PETSC_TRUE : PETSC_FALSE;
652   NcI = Tf[field]->Nc;
653   NcS = isCohesiveField ? NcI : 2*NcI;
654   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
655   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
656   Np = fgeom->numPoints;
657   dE = fgeom->dimEmbed;
658   isAffine = fgeom->isAffine;
659   for (e = 0; e < Ne; ++e) {
660     PetscFEGeom    fegeom;
661     const PetscInt face = fgeom->face[e][0];
662 
663     fegeom.dim      = fgeom->dim;
664     fegeom.dimEmbed = fgeom->dimEmbed;
665     if (isAffine) {
666       fegeom.v    = x;
667       fegeom.xi   = fgeom->xi;
668       fegeom.J    = &fgeom->J[e*dE*dE];
669       fegeom.invJ = &fgeom->invJ[e*dE*dE];
670       fegeom.detJ = &fgeom->detJ[e];
671       fegeom.n    = &fgeom->n[e*dE];
672     }
673     ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr);
674     ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr);
675     for (q = 0; q < Nq; ++q) {
676       PetscReal w;
677       PetscInt  c, d;
678 
679       if (isAffine) {
680         CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
681       } else {
682         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
683         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
684         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
685         fegeom.detJ = &fgeom->detJ[e*Np+q];
686         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
687       }
688       w = fegeom.detJ[0]*quadWeights[q];
689       if (debug > 1 && q < Np) {
690         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
691 #if !defined(PETSC_USE_COMPLEX)
692         ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr);
693 #endif
694       }
695       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
696       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
697       ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
698       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
699       for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q*NcS]);
700       for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w;
701       for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q*NcS*dim]);
702       for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w;
703     }
704     if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);}
705     else                 {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);}
706     cOffset    += totDim;
707     cOffsetAux += totDimAux;
708   }
709   PetscFunctionReturn(0);
710 }
711 
712 PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
713                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
714 {
715   const PetscInt     debug      = 0;
716   PetscFE            feI, feJ;
717   PetscWeakForm      wf;
718   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
719   PetscInt           n0, n1, n2, n3, i;
720   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
721   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
722   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
723   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
724   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
725   PetscQuadrature    quad;
726   PetscTabulation   *T, *TAux = NULL;
727   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
728   const PetscScalar *constants;
729   PetscReal         *x;
730   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
731   PetscInt           NcI = 0, NcJ = 0;
732   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
733   PetscInt           dE, Np;
734   PetscBool          isAffine;
735   const PetscReal   *quadPoints, *quadWeights;
736   PetscInt           qNc, Nq, q;
737   PetscErrorCode     ierr;
738 
739   PetscFunctionBegin;
740   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
741   fieldI = key.field / Nf;
742   fieldJ = key.field % Nf;
743   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
744   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
745   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
746   ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr);
747   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
748   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
749   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
750   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
751   switch(jtype) {
752   case PETSCFE_JACOBIAN_DYN: ierr = PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
753   case PETSCFE_JACOBIAN_PRE: ierr = PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
754   case PETSCFE_JACOBIAN:     ierr = PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
755   }
756   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
757   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
758   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
759   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
760   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
761   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
762   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
763   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
764   if (dsAux) {
765     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
766     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
767     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
768     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
769     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
770     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
771     if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
772   }
773   NcI = T[fieldI]->Nc;
774   NcJ = T[fieldJ]->Nc;
775   Np  = cgeom->numPoints;
776   dE  = cgeom->dimEmbed;
777   isAffine = cgeom->isAffine;
778   /* Initialize here in case the function is not defined */
779   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
780   ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
781   ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
782   ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
783   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
784   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
785   for (e = 0; e < Ne; ++e) {
786     PetscFEGeom fegeom;
787 
788     fegeom.dim      = cgeom->dim;
789     fegeom.dimEmbed = cgeom->dimEmbed;
790     if (isAffine) {
791       fegeom.v    = x;
792       fegeom.xi   = cgeom->xi;
793       fegeom.J    = &cgeom->J[e*Np*dE*dE];
794       fegeom.invJ = &cgeom->invJ[e*Np*dE*dE];
795       fegeom.detJ = &cgeom->detJ[e*Np];
796     }
797     for (q = 0; q < Nq; ++q) {
798       PetscReal w;
799       PetscInt  c;
800 
801       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
802       if (isAffine) {
803         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
804       } else {
805         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
806         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
807         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
808         fegeom.detJ = &cgeom->detJ[e*Np+q];
809       }
810       w = fegeom.detJ[0]*quadWeights[q];
811       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
812       if (dsAux)        {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
813       if (g0_func) {
814         ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
815         for (i = 0; i < n0; ++i) g0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g0);
816         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
817       }
818       if (g1_func) {
819         ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
820         for (i = 0; i < n1; ++i) g1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g1);
821         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
822       }
823       if (g2_func) {
824         ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
825         for (i = 0; i < n2; ++i) g2_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g2);
826         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
827       }
828       if (g3_func) {
829         ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
830         for (i = 0; i < n3; ++i) g3_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g3);
831         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
832       }
833 
834       ierr = PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr);
835     }
836     if (debug > 1) {
837       PetscInt fc, f, gc, g;
838 
839       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
840       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
841         for (f = 0; f < T[fieldI]->Nb; ++f) {
842           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
843           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
844             for (g = 0; g < T[fieldJ]->Nb; ++g) {
845               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
846               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
847             }
848           }
849           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
850         }
851       }
852     }
853     cOffset    += totDim;
854     cOffsetAux += totDimAux;
855     eOffset    += PetscSqr(totDim);
856   }
857   PetscFunctionReturn(0);
858 }
859 
860 static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
861                                                        const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
862 {
863   const PetscInt     debug      = 0;
864   PetscFE            feI, feJ;
865   PetscBdPointJac    g0_func, g1_func, g2_func, g3_func;
866   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
867   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
868   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
869   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
870   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
871   PetscQuadrature    quad;
872   PetscTabulation   *T, *TAux = NULL;
873   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
874   const PetscScalar *constants;
875   PetscReal         *x;
876   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
877   PetscInt           NcI = 0, NcJ = 0;
878   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
879   PetscBool          isAffine;
880   const PetscReal   *quadPoints, *quadWeights;
881   PetscInt           qNc, Nq, q, Np, dE;
882   PetscErrorCode     ierr;
883 
884   PetscFunctionBegin;
885   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
886   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
887   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
888   ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr);
889   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
890   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
891   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
892   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
893   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
894   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
895   ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);
896   if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0);
897   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
898   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
899   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
900   ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr);
901   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
902   if (dsAux) {
903     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
904     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
905     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
906     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
907     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
908     ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);
909   }
910   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
911   Np = fgeom->numPoints;
912   dE = fgeom->dimEmbed;
913   isAffine = fgeom->isAffine;
914   /* Initialize here in case the function is not defined */
915   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
916   ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
917   ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
918   ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
919   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
920   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
921   for (e = 0; e < Ne; ++e) {
922     PetscFEGeom    fegeom, cgeom;
923     const PetscInt face = fgeom->face[e][0];
924     fegeom.n = NULL;
925     fegeom.v = NULL;
926     fegeom.J = NULL;
927     fegeom.detJ = NULL;
928     fegeom.dim      = fgeom->dim;
929     fegeom.dimEmbed = fgeom->dimEmbed;
930     cgeom.dim       = fgeom->dim;
931     cgeom.dimEmbed  = fgeom->dimEmbed;
932     if (isAffine) {
933       fegeom.v    = x;
934       fegeom.xi   = fgeom->xi;
935       fegeom.J    = &fgeom->J[e*Np*dE*dE];
936       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
937       fegeom.detJ = &fgeom->detJ[e*Np];
938       fegeom.n    = &fgeom->n[e*Np*dE];
939 
940       cgeom.J     = &fgeom->suppJ[0][e*Np*dE*dE];
941       cgeom.invJ  = &fgeom->suppInvJ[0][e*Np*dE*dE];
942       cgeom.detJ  = &fgeom->suppDetJ[0][e*Np];
943     }
944     for (q = 0; q < Nq; ++q) {
945       PetscReal w;
946       PetscInt  c;
947 
948       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
949       if (isAffine) {
950         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
951       } else {
952         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
953         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
954         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
955         fegeom.detJ = &fgeom->detJ[e*Np+q];
956         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
957 
958         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
959         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
960         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
961       }
962       w = fegeom.detJ[0]*quadWeights[q];
963       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
964       if (dsAux)        {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
965       if (g0_func) {
966         ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
967         g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0);
968         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
969       }
970       if (g1_func) {
971         ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
972         g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1);
973         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
974       }
975       if (g2_func) {
976         ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
977         g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2);
978         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
979       }
980       if (g3_func) {
981         ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
982         g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3);
983         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
984       }
985 
986       ierr = PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr);
987     }
988     if (debug > 1) {
989       PetscInt fc, f, gc, g;
990 
991       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
992       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
993         for (f = 0; f < T[fieldI]->Nb; ++f) {
994           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
995           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
996             for (g = 0; g < T[fieldJ]->Nb; ++g) {
997               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
998               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
999             }
1000           }
1001           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
1002         }
1003       }
1004     }
1005     cOffset    += totDim;
1006     cOffsetAux += totDimAux;
1007     eOffset    += PetscSqr(totDim);
1008   }
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
1013                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1014 {
1015   const PetscInt     debug      = 0;
1016   PetscFE            feI, feJ;
1017   PetscBdPointJac    g0_func, g1_func, g2_func, g3_func;
1018   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
1019   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
1020   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
1021   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
1022   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
1023   PetscQuadrature    quad;
1024   PetscTabulation   *T, *TAux = NULL;
1025   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
1026   const PetscScalar *constants;
1027   PetscReal         *x;
1028   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
1029   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
1030   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
1031   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE;
1032   const PetscReal   *quadPoints, *quadWeights;
1033   PetscInt           qNc, Nq, q, Np, dE;
1034   PetscErrorCode     ierr;
1035 
1036   PetscFunctionBegin;
1037   /* Hybrid discretization is posed directly on faces */
1038   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
1039   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
1040   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
1041   ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr);
1042   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
1043   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
1044   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
1045   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
1046   switch(jtype) {
1047   case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetBdJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
1048   case PETSCFE_JACOBIAN:     ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
1049   case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
1050   }
1051   if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0);
1052   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
1053   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
1054   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
1055   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
1056   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
1057   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
1058   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
1059   if (dsAux) {
1060     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
1061     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
1062     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
1063     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
1064     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
1065     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
1066     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
1067     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);}
1068     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);}
1069     if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
1070   }
1071   isCohesiveFieldI = fieldI == Nf-1 ? PETSC_TRUE : PETSC_FALSE;
1072   isCohesiveFieldJ = fieldJ == Nf-1 ? PETSC_TRUE : PETSC_FALSE;
1073   NcI = T[fieldI]->Nc;
1074   NcJ = T[fieldJ]->Nc;
1075   NcS = isCohesiveFieldI ? NcI : 2*NcI;
1076   NcT = isCohesiveFieldJ ? NcJ : 2*NcJ;
1077   Np = fgeom->numPoints;
1078   dE = fgeom->dimEmbed;
1079   isAffine = fgeom->isAffine;
1080   ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr);
1081   ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr);
1082   ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr);
1083   ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr);
1084   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
1085   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
1086   for (e = 0; e < Ne; ++e) {
1087     PetscFEGeom    fegeom;
1088     const PetscInt face = fgeom->face[e][0];
1089 
1090     fegeom.dim      = fgeom->dim;
1091     fegeom.dimEmbed = fgeom->dimEmbed;
1092     if (isAffine) {
1093       fegeom.v    = x;
1094       fegeom.xi   = fgeom->xi;
1095       fegeom.J    = &fgeom->J[e*dE*dE];
1096       fegeom.invJ = &fgeom->invJ[e*dE*dE];
1097       fegeom.detJ = &fgeom->detJ[e];
1098       fegeom.n    = &fgeom->n[e*dE];
1099     }
1100     for (q = 0; q < Nq; ++q) {
1101       PetscReal w;
1102       PetscInt  c;
1103 
1104       if (isAffine) {
1105         /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */
1106         CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
1107       } else {
1108         fegeom.v    = &fegeom.v[(e*Np+q)*dE];
1109         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
1110         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
1111         fegeom.detJ = &fgeom->detJ[e*Np+q];
1112         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
1113       }
1114       w = fegeom.detJ[0]*quadWeights[q];
1115       if (debug > 1 && q < Np) {
1116         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
1117 #if !defined(PETSC_USE_COMPLEX)
1118         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
1119 #endif
1120       }
1121       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
1122       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
1123       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
1124       if (g0_func) {
1125         ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr);
1126         g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0);
1127         for (c = 0; c < NcS*NcT; ++c) g0[c] *= w;
1128       }
1129       if (g1_func) {
1130         ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr);
1131         g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1);
1132         for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w;
1133       }
1134       if (g2_func) {
1135         ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr);
1136         g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2);
1137         for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w;
1138       }
1139       if (g3_func) {
1140         ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr);
1141         g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3);
1142         for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w;
1143       }
1144 
1145       if (isCohesiveFieldI && isCohesiveFieldJ) {
1146         ierr = PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI*2, offsetJ*2, elemMat);CHKERRQ(ierr);
1147       } else {
1148         ierr = PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI*2, offsetJ*2, elemMat);CHKERRQ(ierr);
1149       }
1150     }
1151     if (debug > 1) {
1152       PetscInt fc, f, gc, g;
1153 
1154       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
1155       for (fc = 0; fc < NcI; ++fc) {
1156         for (f = 0; f < T[fieldI]->Nb; ++f) {
1157           const PetscInt i = offsetI + f*NcI+fc;
1158           for (gc = 0; gc < NcJ; ++gc) {
1159             for (g = 0; g < T[fieldJ]->Nb; ++g) {
1160               const PetscInt j = offsetJ + g*NcJ+gc;
1161               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
1162             }
1163           }
1164           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
1165         }
1166       }
1167     }
1168     cOffset    += totDim;
1169     cOffsetAux += totDimAux;
1170     eOffset    += PetscSqr(totDim);
1171   }
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1176 {
1177   PetscFunctionBegin;
1178   fem->ops->setfromoptions          = NULL;
1179   fem->ops->setup                   = PetscFESetUp_Basic;
1180   fem->ops->view                    = PetscFEView_Basic;
1181   fem->ops->destroy                 = PetscFEDestroy_Basic;
1182   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1183   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
1184   fem->ops->integrate               = PetscFEIntegrate_Basic;
1185   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
1186   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
1187   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
1188   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
1189   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
1190   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
1191   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
1192   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 /*MC
1197   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
1198 
1199   Level: intermediate
1200 
1201 .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1202 M*/
1203 
1204 PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1205 {
1206   PetscFE_Basic *b;
1207   PetscErrorCode ierr;
1208 
1209   PetscFunctionBegin;
1210   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1211   ierr      = PetscNewLog(fem,&b);CHKERRQ(ierr);
1212   fem->data = b;
1213 
1214   ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr);
1215   PetscFunctionReturn(0);
1216 }
1217