xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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, PetscFormKey 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   const PetscReal   *quadPoints, *quadWeights;
386   PetscInt           qdim, qNc, Nq, q, dE;
387   PetscErrorCode     ierr;
388 
389   PetscFunctionBegin;
390   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
391   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
392   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
393   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
394   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
395   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
396   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
397   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
398   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
399   ierr = PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr);
400   if (!n0 && !n1) PetscFunctionReturn(0);
401   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
402   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
403   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
404   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
405   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
406   if (dsAux) {
407     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
408     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
409     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
410     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
411     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
412     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
413     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);
414   }
415   ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
416   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
417   dE = cgeom->dimEmbed;
418   if (cgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", cgeom->dim, qdim);
419   for (e = 0; e < Ne; ++e) {
420     PetscFEGeom fegeom;
421 
422     fegeom.v = x; /* workspace */
423     ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr);
424     ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dE);CHKERRQ(ierr);
425     for (q = 0; q < Nq; ++q) {
426       PetscReal w;
427       PetscInt  c, d;
428 
429       ierr = PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q*cgeom->dim], &fegeom);CHKERRQ(ierr);
430       w = fegeom.detJ[0]*quadWeights[q];
431       if (debug > 1 && q < cgeom->numPoints) {
432         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
433 #if !defined(PETSC_USE_COMPLEX)
434         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
435 #endif
436       }
437       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d wt %g\n", q, quadWeights[q]);CHKERRQ(ierr);}
438       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
439       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
440       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]);
441       for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w;
442       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]);
443       for (c = 0; c < T[field]->Nc; ++c) for (d = 0; d < dim; ++d) f1[(q*T[field]->Nc+c)*dim+d] *= w;
444     }
445     ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
446     cOffset    += totDim;
447     cOffsetAux += totDimAux;
448   }
449   PetscFunctionReturn(0);
450 }
451 
452 PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
453                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
454 {
455   const PetscInt     debug = 0;
456   const PetscInt     field = key.field;
457   PetscFE            fe;
458   PetscInt           n0,       n1, i;
459   PetscBdPointFunc  *f0_func, *f1_func;
460   PetscQuadrature    quad;
461   PetscTabulation   *Tf, *TfAux = NULL;
462   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
463   const PetscScalar *constants;
464   PetscReal         *x;
465   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
466   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
467   PetscBool          auxOnBd = PETSC_FALSE;
468   const PetscReal   *quadPoints, *quadWeights;
469   PetscInt           qdim, qNc, Nq, q, dE;
470   PetscErrorCode     ierr;
471 
472   PetscFunctionBegin;
473   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
474   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
475   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
476   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
477   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
478   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
479   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
480   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
481   ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr);
482   if (!n0 && !n1) PetscFunctionReturn(0);
483   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
484   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
485   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
486   ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr);
487   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
488   if (dsAux) {
489     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
490     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
491     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
492     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
493     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
494     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
495     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
496     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
497     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
498     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);
499   }
500   NcI = Tf[field]->Nc;
501   ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
502   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
503   dE = fgeom->dimEmbed;
504   /* TODO FIX THIS */
505   fgeom->dim = dim-1;
506   if (fgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", fgeom->dim, qdim);
507   for (e = 0; e < Ne; ++e) {
508     PetscFEGeom    fegeom, cgeom;
509     const PetscInt face = fgeom->face[e][0];
510 
511     fegeom.v = x; /* Workspace */
512     ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr);
513     ierr = PetscArrayzero(f1, Nq*NcI*dE);CHKERRQ(ierr);
514     for (q = 0; q < Nq; ++q) {
515       PetscReal w;
516       PetscInt  c, d;
517 
518       ierr = PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom);CHKERRQ(ierr);
519       ierr = PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom);CHKERRQ(ierr);
520       w = fegeom.detJ[0]*quadWeights[q];
521       if (debug > 1) {
522         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
523           ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
524 #if !defined(PETSC_USE_COMPLEX)
525           ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
526           ierr = DMPrintCellVector(e, "n", dim, fegeom.n);CHKERRQ(ierr);
527 #endif
528         }
529       }
530       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
531       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
532       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*NcI]);
533       for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w;
534       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*NcI*dim]);
535       for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w;
536       if (debug) {
537         ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %D quad point %d\n", e, q);CHKERRQ(ierr);
538         for (c = 0; c < NcI; ++c) {
539           if (n0) {ierr = PetscPrintf(PETSC_COMM_SELF, "  f0[%D] %g\n", c, f0[q*NcI+c]);CHKERRQ(ierr);}
540           if (n1) {
541             for (d = 0; d < dim; ++d) {ierr = PetscPrintf(PETSC_COMM_SELF, "  f1[%D,%D] %g", c, d, f1[(q*NcI + c)*dim + d]);CHKERRQ(ierr);}
542             ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
543           }
544         }
545       }
546     }
547     ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
548     cOffset    += totDim;
549     cOffsetAux += totDimAux;
550   }
551   PetscFunctionReturn(0);
552 }
553 
554 /*
555   BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but
556               all transforms operate in the full space and are square.
557 
558   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
559     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
560     2) We need to assume that the orientation is 0 for both
561     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()
562 */
563 static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
564                                                            const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
565 {
566   const PetscInt     debug = 0;
567   const PetscInt     field = key.field;
568   PetscFE            fe;
569   PetscWeakForm      wf;
570   PetscInt           n0,      n1, i;
571   PetscBdPointFunc  *f0_func, *f1_func;
572   PetscQuadrature    quad;
573   PetscTabulation   *Tf, *TfAux = NULL;
574   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
575   const PetscScalar *constants;
576   PetscReal         *x;
577   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
578   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
579   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
580   const PetscReal   *quadPoints, *quadWeights;
581   PetscInt           qdim, qNc, Nq, q, dE;
582   PetscErrorCode     ierr;
583 
584   PetscFunctionBegin;
585   /* Hybrid discretization is posed directly on faces */
586   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
587   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
588   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
589   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
590   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
591   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
592   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
593   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
594   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
595   ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr);
596   if (!n0 && !n1) PetscFunctionReturn(0);
597   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
598   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
599   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
600   /* NOTE This is a bulk tabulation because the DS is a face discretization */
601   ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr);
602   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
603   if (dsAux) {
604     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
605     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
606     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
607     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
608     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
609     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
610     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
611     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
612     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
613     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);
614   }
615   isCohesiveField = field == Nf-1 ? PETSC_TRUE : PETSC_FALSE;
616   NcI = Tf[field]->Nc;
617   NcS = isCohesiveField ? NcI : 2*NcI;
618   ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
619   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
620   dE = fgeom->dimEmbed;
621   if (fgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", fgeom->dim, qdim);
622   for (e = 0; e < Ne; ++e) {
623     PetscFEGeom    fegeom;
624     const PetscInt face = fgeom->face[e][0];
625 
626     fegeom.v = x; /* Workspace */
627     ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr);
628     ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr);
629     for (q = 0; q < Nq; ++q) {
630       PetscReal w;
631       PetscInt  c, d;
632 
633       ierr = PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom);CHKERRQ(ierr);
634       w = fegeom.detJ[0]*quadWeights[q];
635       if (debug > 1 && q < fgeom->numPoints) {
636         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
637 #if !defined(PETSC_USE_COMPLEX)
638         ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr);
639 #endif
640       }
641       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
642       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
643       ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
644       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);}
645       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]);
646       for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w;
647       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]);
648       for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w;
649     }
650     if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset*2]);}
651     else                 {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset+fOffset*2]);}
652     cOffset    += totDim;
653     cOffsetAux += totDimAux;
654   }
655   PetscFunctionReturn(0);
656 }
657 
658 PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
659                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
660 {
661   const PetscInt     debug      = 0;
662   PetscFE            feI, feJ;
663   PetscWeakForm      wf;
664   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
665   PetscInt           n0, n1, n2, n3, i;
666   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
667   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
668   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
669   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
670   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
671   PetscQuadrature    quad;
672   PetscTabulation   *T, *TAux = NULL;
673   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
674   const PetscScalar *constants;
675   PetscReal         *x;
676   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
677   PetscInt           NcI = 0, NcJ = 0;
678   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
679   PetscInt           dE, Np;
680   PetscBool          isAffine;
681   const PetscReal   *quadPoints, *quadWeights;
682   PetscInt           qNc, Nq, q;
683   PetscErrorCode     ierr;
684 
685   PetscFunctionBegin;
686   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
687   fieldI = key.field / Nf;
688   fieldJ = key.field % Nf;
689   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
690   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
691   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
692   ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr);
693   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
694   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
695   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
696   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
697   switch(jtype) {
698   case PETSCFE_JACOBIAN_DYN: ierr = PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
699   case PETSCFE_JACOBIAN_PRE: ierr = PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
700   case PETSCFE_JACOBIAN:     ierr = PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
701   }
702   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
703   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
704   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
705   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
706   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
707   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
708   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
709   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
710   if (dsAux) {
711     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
712     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
713     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
714     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
715     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
716     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
717     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);
718   }
719   NcI = T[fieldI]->Nc;
720   NcJ = T[fieldJ]->Nc;
721   Np  = cgeom->numPoints;
722   dE  = cgeom->dimEmbed;
723   isAffine = cgeom->isAffine;
724   /* Initialize here in case the function is not defined */
725   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
726   ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
727   ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
728   ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
729   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
730   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
731   for (e = 0; e < Ne; ++e) {
732     PetscFEGeom fegeom;
733 
734     fegeom.dim      = cgeom->dim;
735     fegeom.dimEmbed = cgeom->dimEmbed;
736     if (isAffine) {
737       fegeom.v    = x;
738       fegeom.xi   = cgeom->xi;
739       fegeom.J    = &cgeom->J[e*Np*dE*dE];
740       fegeom.invJ = &cgeom->invJ[e*Np*dE*dE];
741       fegeom.detJ = &cgeom->detJ[e*Np];
742     }
743     for (q = 0; q < Nq; ++q) {
744       PetscReal w;
745       PetscInt  c;
746 
747       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
748       if (isAffine) {
749         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
750       } else {
751         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
752         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
753         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
754         fegeom.detJ = &cgeom->detJ[e*Np+q];
755       }
756       w = fegeom.detJ[0]*quadWeights[q];
757       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
758       if (dsAux)        {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
759       if (n0) {
760         ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
761         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);
762         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
763       }
764       if (n1) {
765         ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
766         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);
767         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
768       }
769       if (n2) {
770         ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
771         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);
772         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
773       }
774       if (n3) {
775         ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
776         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);
777         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
778       }
779 
780       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);
781     }
782     if (debug > 1) {
783       PetscInt fc, f, gc, g;
784 
785       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
786       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
787         for (f = 0; f < T[fieldI]->Nb; ++f) {
788           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
789           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
790             for (g = 0; g < T[fieldJ]->Nb; ++g) {
791               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
792               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
793             }
794           }
795           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
796         }
797       }
798     }
799     cOffset    += totDim;
800     cOffsetAux += totDimAux;
801     eOffset    += PetscSqr(totDim);
802   }
803   PetscFunctionReturn(0);
804 }
805 
806 static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
807                                                        const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
808 {
809   const PetscInt     debug      = 0;
810   PetscFE            feI, feJ;
811   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
812   PetscInt           n0,       n1,       n2,       n3, i;
813   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
814   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
815   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
816   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
817   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
818   PetscQuadrature    quad;
819   PetscTabulation   *T, *TAux = NULL;
820   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
821   const PetscScalar *constants;
822   PetscReal         *x;
823   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
824   PetscInt           NcI = 0, NcJ = 0;
825   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
826   PetscBool          isAffine;
827   const PetscReal   *quadPoints, *quadWeights;
828   PetscInt           qNc, Nq, q, Np, dE;
829   PetscErrorCode     ierr;
830 
831   PetscFunctionBegin;
832   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
833   fieldI = key.field / Nf;
834   fieldJ = key.field % Nf;
835   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
836   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
837   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
838   ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr);
839   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
840   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
841   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
842   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
843   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
844   ierr = PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);
845   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
846   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
847   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
848   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
849   ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr);
850   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
851   if (dsAux) {
852     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
853     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
854     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
855     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
856     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
857     ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);
858   }
859   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
860   Np = fgeom->numPoints;
861   dE = fgeom->dimEmbed;
862   isAffine = fgeom->isAffine;
863   /* Initialize here in case the function is not defined */
864   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
865   ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
866   ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
867   ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
868   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
869   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
870   for (e = 0; e < Ne; ++e) {
871     PetscFEGeom    fegeom, cgeom;
872     const PetscInt face = fgeom->face[e][0];
873     fegeom.n = NULL;
874     fegeom.v = NULL;
875     fegeom.J = NULL;
876     fegeom.detJ = NULL;
877     fegeom.dim      = fgeom->dim;
878     fegeom.dimEmbed = fgeom->dimEmbed;
879     cgeom.dim       = fgeom->dim;
880     cgeom.dimEmbed  = fgeom->dimEmbed;
881     if (isAffine) {
882       fegeom.v    = x;
883       fegeom.xi   = fgeom->xi;
884       fegeom.J    = &fgeom->J[e*Np*dE*dE];
885       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
886       fegeom.detJ = &fgeom->detJ[e*Np];
887       fegeom.n    = &fgeom->n[e*Np*dE];
888 
889       cgeom.J     = &fgeom->suppJ[0][e*Np*dE*dE];
890       cgeom.invJ  = &fgeom->suppInvJ[0][e*Np*dE*dE];
891       cgeom.detJ  = &fgeom->suppDetJ[0][e*Np];
892     }
893     for (q = 0; q < Nq; ++q) {
894       PetscReal w;
895       PetscInt  c;
896 
897       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
898       if (isAffine) {
899         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
900       } else {
901         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
902         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
903         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
904         fegeom.detJ = &fgeom->detJ[e*Np+q];
905         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
906 
907         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
908         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
909         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
910       }
911       w = fegeom.detJ[0]*quadWeights[q];
912       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
913       if (dsAux)        {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
914       if (n0) {
915         ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
916         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, fegeom.n, numConstants, constants, g0);
917         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
918       }
919       if (n1) {
920         ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr);
921         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, fegeom.n, numConstants, constants, g1);
922         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
923       }
924       if (n2) {
925         ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr);
926         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, fegeom.n, numConstants, constants, g2);
927         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
928       }
929       if (n3) {
930         ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr);
931         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, fegeom.n, numConstants, constants, g3);
932         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
933       }
934 
935       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);
936     }
937     if (debug > 1) {
938       PetscInt fc, f, gc, g;
939 
940       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
941       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
942         for (f = 0; f < T[fieldI]->Nb; ++f) {
943           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
944           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
945             for (g = 0; g < T[fieldJ]->Nb; ++g) {
946               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
947               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
948             }
949           }
950           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
951         }
952       }
953     }
954     cOffset    += totDim;
955     cOffsetAux += totDimAux;
956     eOffset    += PetscSqr(totDim);
957   }
958   PetscFunctionReturn(0);
959 }
960 
961 PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
962                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
963 {
964   const PetscInt     debug      = 0;
965   PetscFE            feI, feJ;
966   PetscWeakForm      wf;
967   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
968   PetscInt           n0,       n1,       n2,       n3, i;
969   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
970   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
971   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
972   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
973   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
974   PetscQuadrature    quad;
975   PetscTabulation   *T, *TAux = NULL;
976   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
977   const PetscScalar *constants;
978   PetscReal         *x;
979   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
980   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
981   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
982   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE;
983   const PetscReal   *quadPoints, *quadWeights;
984   PetscInt           qNc, Nq, q, Np, dE;
985   PetscErrorCode     ierr;
986 
987   PetscFunctionBegin;
988   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
989   fieldI = key.field / Nf;
990   fieldJ = key.field % Nf;
991   /* Hybrid discretization is posed directly on faces */
992   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
993   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
994   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
995   ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr);
996   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
997   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
998   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
999   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
1000   switch(jtype) {
1001   case PETSCFE_JACOBIAN_PRE: ierr = PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
1002   case PETSCFE_JACOBIAN:     ierr = PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break;
1003   case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
1004   }
1005   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
1006   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
1007   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
1008   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
1009   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
1010   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
1011   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
1012   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
1013   if (dsAux) {
1014     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
1015     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
1016     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
1017     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
1018     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
1019     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
1020     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
1021     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);}
1022     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);}
1023     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);
1024   }
1025   isCohesiveFieldI = fieldI == Nf-1 ? PETSC_TRUE : PETSC_FALSE;
1026   isCohesiveFieldJ = fieldJ == Nf-1 ? PETSC_TRUE : PETSC_FALSE;
1027   NcI = T[fieldI]->Nc;
1028   NcJ = T[fieldJ]->Nc;
1029   NcS = isCohesiveFieldI ? NcI : 2*NcI;
1030   NcT = isCohesiveFieldJ ? NcJ : 2*NcJ;
1031   Np = fgeom->numPoints;
1032   dE = fgeom->dimEmbed;
1033   isAffine = fgeom->isAffine;
1034   ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr);
1035   ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr);
1036   ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr);
1037   ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr);
1038   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
1039   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
1040   for (e = 0; e < Ne; ++e) {
1041     PetscFEGeom    fegeom;
1042     const PetscInt face = fgeom->face[e][0];
1043 
1044     fegeom.dim      = fgeom->dim;
1045     fegeom.dimEmbed = fgeom->dimEmbed;
1046     if (isAffine) {
1047       fegeom.v    = x;
1048       fegeom.xi   = fgeom->xi;
1049       fegeom.J    = &fgeom->J[e*dE*dE];
1050       fegeom.invJ = &fgeom->invJ[e*dE*dE];
1051       fegeom.detJ = &fgeom->detJ[e];
1052       fegeom.n    = &fgeom->n[e*dE];
1053     }
1054     for (q = 0; q < Nq; ++q) {
1055       PetscReal w;
1056       PetscInt  c;
1057 
1058       if (isAffine) {
1059         /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */
1060         CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
1061       } else {
1062         fegeom.v    = &fegeom.v[(e*Np+q)*dE];
1063         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
1064         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
1065         fegeom.detJ = &fgeom->detJ[e*Np+q];
1066         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
1067       }
1068       w = fegeom.detJ[0]*quadWeights[q];
1069       if (debug > 1 && q < Np) {
1070         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
1071 #if !defined(PETSC_USE_COMPLEX)
1072         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
1073 #endif
1074       }
1075       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
1076       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
1077       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);}
1078       if (n0) {
1079         ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr);
1080         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, fegeom.n, numConstants, constants, g0);
1081         for (c = 0; c < NcS*NcT; ++c) g0[c] *= w;
1082       }
1083       if (n1) {
1084         ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr);
1085         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, fegeom.n, numConstants, constants, g1);
1086         for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w;
1087       }
1088       if (n2) {
1089         ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr);
1090         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, fegeom.n, numConstants, constants, g2);
1091         for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w;
1092       }
1093       if (n3) {
1094         ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr);
1095         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, fegeom.n, numConstants, constants, g3);
1096         for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w;
1097       }
1098 
1099       if (isCohesiveFieldI && isCohesiveFieldJ) {
1100         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);
1101       } else {
1102         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);
1103       }
1104     }
1105     if (debug > 1) {
1106       PetscInt fc, f, gc, g;
1107 
1108       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
1109       for (fc = 0; fc < NcI; ++fc) {
1110         for (f = 0; f < T[fieldI]->Nb; ++f) {
1111           const PetscInt i = offsetI + f*NcI+fc;
1112           for (gc = 0; gc < NcJ; ++gc) {
1113             for (g = 0; g < T[fieldJ]->Nb; ++g) {
1114               const PetscInt j = offsetJ + g*NcJ+gc;
1115               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
1116             }
1117           }
1118           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
1119         }
1120       }
1121     }
1122     cOffset    += totDim;
1123     cOffsetAux += totDimAux;
1124     eOffset    += PetscSqr(totDim);
1125   }
1126   PetscFunctionReturn(0);
1127 }
1128 
1129 static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1130 {
1131   PetscFunctionBegin;
1132   fem->ops->setfromoptions          = NULL;
1133   fem->ops->setup                   = PetscFESetUp_Basic;
1134   fem->ops->view                    = PetscFEView_Basic;
1135   fem->ops->destroy                 = PetscFEDestroy_Basic;
1136   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1137   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
1138   fem->ops->integrate               = PetscFEIntegrate_Basic;
1139   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
1140   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
1141   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
1142   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
1143   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
1144   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
1145   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
1146   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 /*MC
1151   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
1152 
1153   Level: intermediate
1154 
1155 .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
1156 M*/
1157 
1158 PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1159 {
1160   PetscFE_Basic *b;
1161   PetscErrorCode ierr;
1162 
1163   PetscFunctionBegin;
1164   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1165   ierr      = PetscNewLog(fem,&b);CHKERRQ(ierr);
1166   fem->data = b;
1167 
1168   ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr);
1169   PetscFunctionReturn(0);
1170 }
1171