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