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