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