xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 87510d7d5654fee7c2f87f0a72f2261ac85db6ad)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
5d71ae5a4SJacob Faibussowitsch {
620cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *)fem->data;
720cf1dd8SToby Isaac 
820cf1dd8SToby Isaac   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1120cf1dd8SToby Isaac }
1220cf1dd8SToby Isaac 
13d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
14d71ae5a4SJacob Faibussowitsch {
15d9bac1caSLisandro Dalcin   PetscInt        dim, Nc;
16d9bac1caSLisandro Dalcin   PetscSpace      basis = NULL;
17d9bac1caSLisandro Dalcin   PetscDualSpace  dual  = NULL;
18d9bac1caSLisandro Dalcin   PetscQuadrature quad  = NULL;
1920cf1dd8SToby Isaac 
2020cf1dd8SToby Isaac   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe, &Nc));
239566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &basis));
249566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dual));
259566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(v));
2763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc));
289566063dSJacob Faibussowitsch   if (basis) PetscCall(PetscSpaceView(basis, v));
299566063dSJacob Faibussowitsch   if (dual) PetscCall(PetscDualSpaceView(dual, v));
309566063dSJacob Faibussowitsch   if (quad) PetscCall(PetscQuadratureView(quad, v));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(v));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
36d71ae5a4SJacob Faibussowitsch {
3720cf1dd8SToby Isaac   PetscBool iascii;
3820cf1dd8SToby Isaac 
3920cf1dd8SToby Isaac   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
419566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscFEView_Basic_Ascii(fe, v));
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4320cf1dd8SToby Isaac }
4420cf1dd8SToby Isaac 
4520cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
46d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
47d71ae5a4SJacob Faibussowitsch {
48b9d4cb8dSJed Brown   PetscReal    *work;
4920cf1dd8SToby Isaac   PetscBLASInt *pivots;
5020cf1dd8SToby Isaac   PetscBLASInt  n, info;
5120cf1dd8SToby Isaac   PetscInt      pdim, j;
5220cf1dd8SToby Isaac 
5320cf1dd8SToby Isaac   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pdim * pdim, &fem->invV));
5620cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
5720cf1dd8SToby Isaac     PetscReal       *Bf;
5820cf1dd8SToby Isaac     PetscQuadrature  f;
5920cf1dd8SToby Isaac     const PetscReal *points, *weights;
6020cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
6120cf1dd8SToby Isaac 
629566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f));
639566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf));
659566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
6620cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
6720cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
68b9d4cb8dSJed Brown       fem->invV[j * pdim + k] = 0.0;
6920cf1dd8SToby Isaac 
7020cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
71b9d4cb8dSJed Brown         for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c];
7220cf1dd8SToby Isaac       }
7320cf1dd8SToby Isaac     }
749566063dSJacob Faibussowitsch     PetscCall(PetscFree(Bf));
7520cf1dd8SToby Isaac   }
76ea2bdf6dSBarry Smith 
779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work));
7820cf1dd8SToby Isaac   n = pdim;
79792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
8063a3b9bcSJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info);
81792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
8263a3b9bcSJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)info);
839566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pivots, work));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8520cf1dd8SToby Isaac }
8620cf1dd8SToby Isaac 
87d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
88d71ae5a4SJacob Faibussowitsch {
8920cf1dd8SToby Isaac   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9220cf1dd8SToby Isaac }
9320cf1dd8SToby Isaac 
94b9d4cb8dSJed Brown /* Tensor contraction on the middle index,
95b9d4cb8dSJed Brown  *    C[m,n,p] = A[m,k,p] * B[k,n]
96b9d4cb8dSJed Brown  * where all matrices use C-style ordering.
97b9d4cb8dSJed Brown  */
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C)
99d71ae5a4SJacob Faibussowitsch {
100b9d4cb8dSJed Brown   PetscInt i;
101b9d4cb8dSJed Brown 
102b9d4cb8dSJed Brown   PetscFunctionBegin;
103aa9788aaSMatthew G. Knepley   PetscCheck(n && p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Empty tensor is not allowed %" PetscInt_FMT " %" PetscInt_FMT, n, p);
104b9d4cb8dSJed Brown   for (i = 0; i < m; i++) {
105b9d4cb8dSJed Brown     PetscBLASInt n_, p_, k_, lda, ldb, ldc;
106b9d4cb8dSJed Brown     PetscReal    one = 1, zero = 0;
107b9d4cb8dSJed Brown     /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
108b9d4cb8dSJed Brown      * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
109b9d4cb8dSJed Brown      */
1109566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &n_));
1119566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(p, &p_));
1129566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(k, &k_));
113b9d4cb8dSJed Brown     lda = p_;
114b9d4cb8dSJed Brown     ldb = n_;
115b9d4cb8dSJed Brown     ldc = p_;
116792fecdfSBarry Smith     PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc));
117b9d4cb8dSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2. * m * n * p * k));
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120b9d4cb8dSJed Brown }
121b9d4cb8dSJed Brown 
122d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
123d71ae5a4SJacob Faibussowitsch {
12420cf1dd8SToby Isaac   DM         dm;
12520cf1dd8SToby Isaac   PetscInt   pdim; /* Dimension of FE space P */
12620cf1dd8SToby Isaac   PetscInt   dim;  /* Spatial dimension */
12720cf1dd8SToby Isaac   PetscInt   Nc;   /* Field components */
128ef0bb6c7SMatthew G. Knepley   PetscReal *B    = K >= 0 ? T->T[0] : NULL;
129ef0bb6c7SMatthew G. Knepley   PetscReal *D    = K >= 1 ? T->T[1] : NULL;
130ef0bb6c7SMatthew G. Knepley   PetscReal *H    = K >= 2 ? T->T[2] : NULL;
131ef0bb6c7SMatthew G. Knepley   PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
13220cf1dd8SToby Isaac 
13320cf1dd8SToby Isaac   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
1359566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1369566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
13820cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
1399566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1409566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1419566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1429566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
143b9d4cb8dSJed Brown   /* Translate from prime to nodal basis */
14420cf1dd8SToby Isaac   if (B) {
145b9d4cb8dSJed Brown     /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
1469566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B));
14720cf1dd8SToby Isaac   }
148aa9788aaSMatthew G. Knepley   if (D && dim) {
149b9d4cb8dSJed Brown     /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
1509566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D));
15120cf1dd8SToby Isaac   }
152aa9788aaSMatthew G. Knepley   if (H && dim) {
153b9d4cb8dSJed Brown     /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
1549566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H));
15520cf1dd8SToby Isaac   }
1569566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1579566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1589566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16020cf1dd8SToby Isaac }
16120cf1dd8SToby Isaac 
1622dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
163d71ae5a4SJacob Faibussowitsch {
164b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
1654bee2e38SMatthew G. Knepley   PetscFE            fe;
16620cf1dd8SToby Isaac   PetscPointFunc     obj_func;
16720cf1dd8SToby Isaac   PetscQuadrature    quad;
168ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
1694bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
17020cf1dd8SToby Isaac   const PetscScalar *constants;
171*87510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
172ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
17320cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
17420cf1dd8SToby Isaac   PetscBool          isAffine;
17520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
17620cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
17720cf1dd8SToby Isaac 
17820cf1dd8SToby Isaac   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetObjective(ds, field, &obj_func));
1803ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
1819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
183*87510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
1849566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
1859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
1889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
1899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
1909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
1919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
192*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
1939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1944bee2e38SMatthew G. Knepley   if (dsAux) {
1959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1999566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
2009566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
20163a3b9bcSJacob Faibussowitsch     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);
20220cf1dd8SToby Isaac   }
2039566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
20463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
20520cf1dd8SToby Isaac   Np       = cgeom->numPoints;
20620cf1dd8SToby Isaac   dE       = cgeom->dimEmbed;
20720cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
20820cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2094bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
21020cf1dd8SToby Isaac 
21127f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
21227f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
21320cf1dd8SToby Isaac     if (isAffine) {
2144bee2e38SMatthew G. Knepley       fegeom.v    = x;
2154bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2167132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
2177132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
2187132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
21920cf1dd8SToby Isaac     }
2204bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
221d627b919SMatthew G. Knepley       PetscScalar integrand = 0.;
2224bee2e38SMatthew G. Knepley       PetscReal   w;
2234bee2e38SMatthew G. Knepley 
2244bee2e38SMatthew G. Knepley       if (isAffine) {
2257132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
2264bee2e38SMatthew G. Knepley       } else {
2274bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
2284bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
2294bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
2304bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
2314bee2e38SMatthew G. Knepley       }
232*87510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
2334bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
23420cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
23563a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
2367be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2379566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
23820cf1dd8SToby Isaac #endif
23920cf1dd8SToby Isaac       }
24063a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
2419566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
2429566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
2434bee2e38SMatthew G. Knepley       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);
2444bee2e38SMatthew G. Knepley       integrand *= w;
24520cf1dd8SToby Isaac       integral[e * Nf + field] += integrand;
2469566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[field])));
24720cf1dd8SToby Isaac     }
24820cf1dd8SToby Isaac     cOffset += totDim;
24920cf1dd8SToby Isaac     cOffsetAux += totDimAux;
25020cf1dd8SToby Isaac   }
2513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25220cf1dd8SToby Isaac }
25320cf1dd8SToby Isaac 
2542dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFunc obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
255d71ae5a4SJacob Faibussowitsch {
256b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
2574bee2e38SMatthew G. Knepley   PetscFE            fe;
258afe6d6adSToby Isaac   PetscQuadrature    quad;
259ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2604bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
261afe6d6adSToby Isaac   const PetscScalar *constants;
262*87510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
263ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
264afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
265afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
266afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
267afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
268afe6d6adSToby Isaac 
269afe6d6adSToby Isaac   PetscFunctionBegin;
2703ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
2719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
2729566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
273*87510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
2749566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
2759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
2789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
2799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
2809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
2819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
282*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
2839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
2844bee2e38SMatthew G. Knepley   if (dsAux) {
2859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
2869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
2889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
291afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
2929566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
2939566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
29463a3b9bcSJacob Faibussowitsch     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);
295afe6d6adSToby Isaac   }
2969566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
29763a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
29879ab67a3SMatthew G. Knepley   if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Field: %" PetscInt_FMT " Nface: %" PetscInt_FMT " Nq: %" PetscInt_FMT "\n", field, Ne, Nq));
299afe6d6adSToby Isaac   Np       = fgeom->numPoints;
300afe6d6adSToby Isaac   dE       = fgeom->dimEmbed;
301afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
302afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
3039f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
304afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
305ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
306ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
307ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
308b2deab97SMatthew G. Knepley     fegeom.invJ         = NULL;
309ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
31027f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
31127f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
31227f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
31327f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
3144bee2e38SMatthew G. Knepley     if (isAffine) {
3154bee2e38SMatthew G. Knepley       fegeom.v    = x;
3164bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3177132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
3187132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
3197132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
3207132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
3219f209ee4SMatthew G. Knepley 
3227132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
3237132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
3247132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
3254bee2e38SMatthew G. Knepley     }
326afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
32779ab67a3SMatthew G. Knepley       PetscScalar integrand = 0.;
3284bee2e38SMatthew G. Knepley       PetscReal   w;
329afe6d6adSToby Isaac 
330afe6d6adSToby Isaac       if (isAffine) {
3317132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
332afe6d6adSToby Isaac       } else {
3333fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
3349f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
3359f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
3364bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
3374bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
3389f209ee4SMatthew G. Knepley 
3399f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
3409f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
3419f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
342afe6d6adSToby Isaac       }
343*87510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
3444bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
345afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
34663a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
347afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3489566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
349afe6d6adSToby Isaac #endif
350afe6d6adSToby Isaac       }
35163a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
35279ab67a3SMatthew G. Knepley       if (debug > 3) {
35379ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    x_q ("));
35479ab67a3SMatthew G. Knepley         for (PetscInt d = 0; d < dE; ++d) {
35579ab67a3SMatthew G. Knepley           if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
35679ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.v[d]));
35779ab67a3SMatthew G. Knepley         }
35879ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
35979ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    n_q ("));
36079ab67a3SMatthew G. Knepley         for (PetscInt d = 0; d < dE; ++d) {
36179ab67a3SMatthew G. Knepley           if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
36279ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.n[d]));
36379ab67a3SMatthew G. Knepley         }
36479ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
36579ab67a3SMatthew G. Knepley         for (PetscInt f = 0; f < Nf; ++f) {
36679ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "    u_%" PetscInt_FMT " (", f));
36779ab67a3SMatthew G. Knepley           for (PetscInt c = 0; c < uOff[f + 1] - uOff[f]; ++c) {
36879ab67a3SMatthew G. Knepley             if (c) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
36979ab67a3SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(u[uOff[f] + c])));
37079ab67a3SMatthew G. Knepley           }
37179ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
37279ab67a3SMatthew G. Knepley         }
37379ab67a3SMatthew G. Knepley       }
3749566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3759566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3764bee2e38SMatthew G. Knepley       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);
3774bee2e38SMatthew G. Knepley       integrand *= w;
378afe6d6adSToby Isaac       integral[e * Nf + field] += integrand;
37979ab67a3SMatthew G. Knepley       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g tot: %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
380afe6d6adSToby Isaac     }
381afe6d6adSToby Isaac     cOffset += totDim;
382afe6d6adSToby Isaac     cOffsetAux += totDimAux;
383afe6d6adSToby Isaac   }
3843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
385afe6d6adSToby Isaac }
386afe6d6adSToby Isaac 
387d71ae5a4SJacob Faibussowitsch 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[])
388d71ae5a4SJacob Faibussowitsch {
389b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
3906528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3914bee2e38SMatthew G. Knepley   PetscFE            fe;
3926528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3936528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
3946528b96dSMatthew G. Knepley   PetscPointFunc    *f0_func, *f1_func;
39520cf1dd8SToby Isaac   PetscQuadrature    quad;
396ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3974bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
39820cf1dd8SToby Isaac   const PetscScalar *constants;
399*87510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
400ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
401ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
40220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
4036587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
40420cf1dd8SToby Isaac 
40520cf1dd8SToby Isaac   PetscFunctionBegin;
4069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
4079566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
408*87510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
4099566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
4109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
4169566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
4173ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
4189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
422*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
4239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4244bee2e38SMatthew G. Knepley   if (dsAux) {
4259566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4269566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4279566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4289566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4299566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4309566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
43163a3b9bcSJacob Faibussowitsch     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);
43220cf1dd8SToby Isaac   }
4339566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
43463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
43520cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
43663a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
43720cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4384bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
43920cf1dd8SToby Isaac 
4406587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
4419566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
4429566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
44320cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4444bee2e38SMatthew G. Knepley       PetscReal w;
4454bee2e38SMatthew G. Knepley       PetscInt  c, d;
44620cf1dd8SToby Isaac 
4479566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
448*87510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
4494bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
4506587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
45163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4527be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4539566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
45420cf1dd8SToby Isaac #endif
45520cf1dd8SToby Isaac       }
45616cd844bSPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
4579566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
4586528b96dSMatthew G. Knepley       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]);
459ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
4606528b96dSMatthew G. Knepley       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]);
4619371c9d4SSatish Balay       for (c = 0; c < T[field]->Nc; ++c)
4629371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w;
463b8025e53SMatthew G. Knepley       if (debug) {
464e8e188d2SZach Atkins         // LCOV_EXCL_START
465e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q]));
466e8e188d2SZach Atkins         for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c]));
467e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
468b8025e53SMatthew G. Knepley         if (debug > 2) {
46963a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
47063a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
4719566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
472e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field der %" PetscInt_FMT ":", field));
473e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c])));
474e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
47563a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
47663a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
4779566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
478e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  res der %" PetscInt_FMT ":", field));
479e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc; ++c) {
480e8e188d2SZach Atkins             for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dim + d])));
481b8025e53SMatthew G. Knepley           }
482e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
483e8e188d2SZach Atkins         }
484e8e188d2SZach Atkins         // LCOV_EXCL_STOP
485b8025e53SMatthew G. Knepley       }
48620cf1dd8SToby Isaac     }
4879566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
48820cf1dd8SToby Isaac     cOffset += totDim;
48920cf1dd8SToby Isaac     cOffsetAux += totDimAux;
49020cf1dd8SToby Isaac   }
4913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49220cf1dd8SToby Isaac }
49320cf1dd8SToby Isaac 
494d71ae5a4SJacob Faibussowitsch 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[])
495d71ae5a4SJacob Faibussowitsch {
496b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
49706d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
4984bee2e38SMatthew G. Knepley   PetscFE            fe;
49906d8a0d3SMatthew G. Knepley   PetscInt           n0, n1, i;
50006d8a0d3SMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
50120cf1dd8SToby Isaac   PetscQuadrature    quad;
502ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
5034bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
50420cf1dd8SToby Isaac   const PetscScalar *constants;
505*87510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
506ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
507ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
5086587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
50920cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
5106587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
51120cf1dd8SToby Isaac 
51220cf1dd8SToby Isaac   PetscFunctionBegin;
5139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
5149566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
515*87510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
5169566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
5179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
5209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
5219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
5229566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
5233ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
5249566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
5259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
5269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
5279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
528*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
5299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
5304bee2e38SMatthew G. Knepley   if (dsAux) {
5319566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
5329566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
5339566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5349566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
5359566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
5369566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
5377be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
5389566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
5399566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
54063a3b9bcSJacob Faibussowitsch     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);
54120cf1dd8SToby Isaac   }
542ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
5439566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
54463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
54520cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
5466587ee25SMatthew G. Knepley   /* TODO FIX THIS */
5476587ee25SMatthew G. Knepley   fgeom->dim = dim - 1;
54863a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
54920cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5509f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
55120cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5529f209ee4SMatthew G. Knepley 
5536587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
5549566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcI));
5559566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
55620cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5574bee2e38SMatthew G. Knepley       PetscReal w;
5584bee2e38SMatthew G. Knepley       PetscInt  c, d;
5594bee2e38SMatthew G. Knepley 
5609566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
5619566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
562*87510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
5634bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
56462bd480fSMatthew G. Knepley       if (debug > 1) {
5656587ee25SMatthew G. Knepley         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
56663a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5677be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5689566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5699566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
57020cf1dd8SToby Isaac #endif
57120cf1dd8SToby Isaac         }
57262bd480fSMatthew G. Knepley       }
5738e3a54c0SPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
5749566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
57506d8a0d3SMatthew G. Knepley       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]);
5764bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w;
57706d8a0d3SMatthew G. Knepley       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]);
5789371c9d4SSatish Balay       for (c = 0; c < NcI; ++c)
5799371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w;
58062bd480fSMatthew G. Knepley       if (debug) {
58163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
58262bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
58363a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c])));
58462bd480fSMatthew G. Knepley           if (n1) {
58563a3b9bcSJacob Faibussowitsch             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])));
5869566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
58762bd480fSMatthew G. Knepley           }
58862bd480fSMatthew G. Knepley         }
58962bd480fSMatthew G. Knepley       }
59020cf1dd8SToby Isaac     }
5919566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
59220cf1dd8SToby Isaac     cOffset += totDim;
59320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
59420cf1dd8SToby Isaac   }
5953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59620cf1dd8SToby Isaac }
59720cf1dd8SToby Isaac 
59827f02ce8SMatthew G. Knepley /*
59927f02ce8SMatthew G. Knepley   BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but
60027f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
60127f02ce8SMatthew G. Knepley 
60227f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
60327f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
60427f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
60527f02ce8SMatthew G. Knepley     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()
60627f02ce8SMatthew G. Knepley */
6072dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
608d71ae5a4SJacob Faibussowitsch {
609b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
6106528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
61127f02ce8SMatthew G. Knepley   PetscFE            fe;
6126528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6136528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
6146528b96dSMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
61527f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
6160e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
61707218a29SMatthew G. Knepley   PetscTabulation   *Tf, *TfIn, *TfAux = NULL;
61827f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
61927f02ce8SMatthew G. Knepley   const PetscScalar *constants;
62027f02ce8SMatthew G. Knepley   PetscReal         *x;
621665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
62207218a29SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
6236587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
62427f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
6256587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
62627f02ce8SMatthew G. Knepley 
62727f02ce8SMatthew G. Knepley   PetscFunctionBegin;
62827f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
6299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
6309566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
6319566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
6329566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6339566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
63407218a29SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
635429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets
63607218a29SMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x));
6379566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
6389566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
6399566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
6403ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
6419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
6429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
6439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
64427f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
6459566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
64607218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
647*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
6489566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
64927f02ce8SMatthew G. Knepley   if (dsAux) {
6509566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
6519566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6529566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6539566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6549566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6559566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
65601907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
6579566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
6589566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
65963a3b9bcSJacob Faibussowitsch     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);
66027f02ce8SMatthew G. Knepley   }
6619566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
662665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
663c2b7495fSMatthew G. Knepley   NcS = NcI;
6640abb75b6SMatthew G. Knepley   if (!isCohesiveField && s == 2) {
6650abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
6660abb75b6SMatthew G. Knepley     NcS *= 2;
6670abb75b6SMatthew G. Knepley   }
6689566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
6690e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
67063a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
67127f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
67263a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
67327f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
67427f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
6750e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
6760e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
6774e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
67827f02ce8SMatthew G. Knepley 
6796587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6809566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcS));
6819566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
68227f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
6830e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
68427f02ce8SMatthew G. Knepley       PetscReal w;
68527f02ce8SMatthew G. Knepley       PetscInt  c, d;
68627f02ce8SMatthew G. Knepley 
6874e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0]));
6884e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1]));
68907218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
69027f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
6916587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
69263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
69327f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
6949566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
69527f02ce8SMatthew G. Knepley #endif
69627f02ce8SMatthew G. Knepley       }
697a4158a15SMatthew G. Knepley       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0]));
69827f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
6998e3a54c0SPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], PetscSafePointerPlusOffset(coefficients_t, cOffsetIn), u, u_x, u_t));
70007218a29SMatthew G. Knepley       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
7016528b96dSMatthew G. Knepley       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]);
70227f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
7039ee2af8cSMatthew G. Knepley       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]);
7049371c9d4SSatish Balay       for (c = 0; c < NcS; ++c)
7059371c9d4SSatish Balay         for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
70627f02ce8SMatthew G. Knepley     }
7079371c9d4SSatish Balay     if (isCohesiveField) {
7083ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
7099371c9d4SSatish Balay     } else {
7103ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
7119371c9d4SSatish Balay     }
71227f02ce8SMatthew G. Knepley     cOffset += totDim;
71307218a29SMatthew G. Knepley     cOffsetIn += totDimIn;
71427f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
71527f02ce8SMatthew G. Knepley   }
7163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71727f02ce8SMatthew G. Knepley }
71827f02ce8SMatthew G. Knepley 
719d71ae5a4SJacob Faibussowitsch 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[])
720d71ae5a4SJacob Faibussowitsch {
721b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
7224bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
7236528b96dSMatthew G. Knepley   PetscWeakForm      wf;
7246528b96dSMatthew G. Knepley   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
7256528b96dSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
72620cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
72720cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
72820cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
72920cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
73020cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
73120cf1dd8SToby Isaac   PetscQuadrature    quad;
732ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
7334bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
73420cf1dd8SToby Isaac   const PetscScalar *constants;
735*87510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
736ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
737ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
7386528b96dSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
73920cf1dd8SToby Isaac   PetscInt           dE, Np;
74020cf1dd8SToby Isaac   PetscBool          isAffine;
74120cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
74220cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
74320cf1dd8SToby Isaac 
74420cf1dd8SToby Isaac   PetscFunctionBegin;
7459566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
7466528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
7476528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
7489566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
7499566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
7509566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
751*87510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
7529566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
7539566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
7549566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
7559566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
7569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
75720cf1dd8SToby Isaac   switch (jtype) {
758d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
759d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
760d71ae5a4SJacob Faibussowitsch     break;
761d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
762d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
763d71ae5a4SJacob Faibussowitsch     break;
764d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
765d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
766d71ae5a4SJacob Faibussowitsch     break;
76720cf1dd8SToby Isaac   }
7683ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
7699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
7709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
7719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
7729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
7739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
7749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
775a102dd69SStefano Zampini   PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
7769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
7774bee2e38SMatthew G. Knepley   if (dsAux) {
7789566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7799566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
7809566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
7819566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
7829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
7839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
78463a3b9bcSJacob Faibussowitsch     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);
78520cf1dd8SToby Isaac   }
78627f02ce8SMatthew G. Knepley   NcI      = T[fieldI]->Nc;
78727f02ce8SMatthew G. Knepley   NcJ      = T[fieldJ]->Nc;
7884bee2e38SMatthew G. Knepley   Np       = cgeom->numPoints;
7894bee2e38SMatthew G. Knepley   dE       = cgeom->dimEmbed;
7904bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
79127f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
7929566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
7939566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7949566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7959566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7969566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
79763a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
7984bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
7994bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
8004bee2e38SMatthew G. Knepley 
80127f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
80227f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
8034bee2e38SMatthew G. Knepley     if (isAffine) {
8044bee2e38SMatthew G. Knepley       fegeom.v    = x;
8054bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
8067132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
8077132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
8087132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
8094bee2e38SMatthew G. Knepley     }
81020cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
81120cf1dd8SToby Isaac       PetscReal w;
8124bee2e38SMatthew G. Knepley       PetscInt  c;
81320cf1dd8SToby Isaac 
81420cf1dd8SToby Isaac       if (isAffine) {
8157132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
81620cf1dd8SToby Isaac       } else {
8174bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
8184bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
8194bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
8204bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
82120cf1dd8SToby Isaac       }
822*87510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
8239566063dSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0]));
8244bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
82516cd844bSPierre Jolivet       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
8269566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
827ea672e62SMatthew G. Knepley       if (n0) {
8289566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
8296528b96dSMatthew G. Knepley         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);
83020cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
83120cf1dd8SToby Isaac       }
832ea672e62SMatthew G. Knepley       if (n1) {
8339566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
8346528b96dSMatthew G. Knepley         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);
8354bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
83620cf1dd8SToby Isaac       }
837ea672e62SMatthew G. Knepley       if (n2) {
8389566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
8396528b96dSMatthew G. Knepley         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);
8404bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
84120cf1dd8SToby Isaac       }
842ea672e62SMatthew G. Knepley       if (n3) {
8439566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8446528b96dSMatthew G. Knepley         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);
8454bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
84620cf1dd8SToby Isaac       }
84720cf1dd8SToby Isaac 
8489566063dSJacob Faibussowitsch       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));
84920cf1dd8SToby Isaac     }
85020cf1dd8SToby Isaac     if (debug > 1) {
85120cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
85220cf1dd8SToby Isaac 
85363a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
854ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
855ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
856ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
857ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
858ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
859ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
86063a3b9bcSJacob Faibussowitsch               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])));
86120cf1dd8SToby Isaac             }
86220cf1dd8SToby Isaac           }
8639566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
86420cf1dd8SToby Isaac         }
86520cf1dd8SToby Isaac       }
86620cf1dd8SToby Isaac     }
86720cf1dd8SToby Isaac     cOffset += totDim;
86820cf1dd8SToby Isaac     cOffsetAux += totDimAux;
86920cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
87020cf1dd8SToby Isaac   }
8713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87220cf1dd8SToby Isaac }
87320cf1dd8SToby Isaac 
874e3d591f2SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, 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[])
875d71ae5a4SJacob Faibussowitsch {
876b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
8774bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
87845480ffeSMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
87945480ffeSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
88020cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
88120cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
88220cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
88320cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
88420cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
88520cf1dd8SToby Isaac   PetscQuadrature    quad;
886ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
8874bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
88820cf1dd8SToby Isaac   const PetscScalar *constants;
889*87510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
890ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
891ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
89245480ffeSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
89320cf1dd8SToby Isaac   PetscBool          isAffine;
89420cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
89520cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
89620cf1dd8SToby Isaac 
89720cf1dd8SToby Isaac   PetscFunctionBegin;
8989566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
89945480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
90045480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
9019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
9029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
9039566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
904*87510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
9059566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
9069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
9079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
9089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
9099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
9109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
911e3d591f2SMatthew G. Knepley   switch (jtype) {
912e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN_PRE:
913e3d591f2SMatthew G. Knepley     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
914e3d591f2SMatthew G. Knepley     break;
915e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN:
9169566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
917e3d591f2SMatthew G. Knepley     break;
918e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN:
919e3d591f2SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()");
920e3d591f2SMatthew G. Knepley   }
9213ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
9229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
9239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
9249566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
9259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
926*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
9279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
9284bee2e38SMatthew G. Knepley   if (dsAux) {
9299566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
9309566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
9319566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
9329566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
9339566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
9349566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
93520cf1dd8SToby Isaac   }
936ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
93720cf1dd8SToby Isaac   Np       = fgeom->numPoints;
93820cf1dd8SToby Isaac   dE       = fgeom->dimEmbed;
93920cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
94027f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
9419566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
9429566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
9439566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
9449566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
9459566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
94663a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
94720cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
9489f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
94920cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
950ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
951ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
952ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
953ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
95427f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
95527f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
95627f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
95727f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
9584bee2e38SMatthew G. Knepley     if (isAffine) {
9594bee2e38SMatthew G. Knepley       fegeom.v    = x;
9604bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
9617132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
9627132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
9637132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
9647132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
9659f209ee4SMatthew G. Knepley 
9667132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
9677132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
9687132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
9694bee2e38SMatthew G. Knepley     }
97020cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
97120cf1dd8SToby Isaac       PetscReal w;
9724bee2e38SMatthew G. Knepley       PetscInt  c;
97320cf1dd8SToby Isaac 
97463a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
97520cf1dd8SToby Isaac       if (isAffine) {
9767132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
97720cf1dd8SToby Isaac       } else {
9783fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
9799f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
9809f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
9814bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
9824bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
9839f209ee4SMatthew G. Knepley 
9849f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
9859f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
9869f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
98720cf1dd8SToby Isaac       }
988*87510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
9894bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
9909566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
9919566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
992ea672e62SMatthew G. Knepley       if (n0) {
9939566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
99445480ffeSMatthew G. Knepley         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);
99520cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
99620cf1dd8SToby Isaac       }
997ea672e62SMatthew G. Knepley       if (n1) {
9989566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
99945480ffeSMatthew G. Knepley         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);
10004bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
100120cf1dd8SToby Isaac       }
1002ea672e62SMatthew G. Knepley       if (n2) {
10039566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
100445480ffeSMatthew G. Knepley         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);
10054bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
100620cf1dd8SToby Isaac       }
1007ea672e62SMatthew G. Knepley       if (n3) {
10089566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
100945480ffeSMatthew G. Knepley         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);
10104bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
101120cf1dd8SToby Isaac       }
101220cf1dd8SToby Isaac 
10139566063dSJacob Faibussowitsch       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));
101420cf1dd8SToby Isaac     }
101520cf1dd8SToby Isaac     if (debug > 1) {
101620cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
101720cf1dd8SToby Isaac 
101863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
1019ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
1020ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
1021ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
1022ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
1023ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
1024ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
102563a3b9bcSJacob Faibussowitsch               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])));
102620cf1dd8SToby Isaac             }
102720cf1dd8SToby Isaac           }
10289566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
102920cf1dd8SToby Isaac         }
103020cf1dd8SToby Isaac       }
103120cf1dd8SToby Isaac     }
103220cf1dd8SToby Isaac     cOffset += totDim;
103320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
103420cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
103520cf1dd8SToby Isaac   }
10363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103720cf1dd8SToby Isaac }
103820cf1dd8SToby Isaac 
10392dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscDS dsIn, 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[])
1040d71ae5a4SJacob Faibussowitsch {
1041b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
104227f02ce8SMatthew G. Knepley   PetscFE            feI, feJ;
1043148442b3SMatthew G. Knepley   PetscWeakForm      wf;
1044148442b3SMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
1045148442b3SMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
104627f02ce8SMatthew G. Knepley   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
104727f02ce8SMatthew G. Knepley   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
104827f02ce8SMatthew G. Knepley   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
104927f02ce8SMatthew G. Knepley   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
105027f02ce8SMatthew G. Knepley   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
1051665f567fSMatthew G. Knepley   PetscQuadrature    quad;
10520e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
105307218a29SMatthew G. Knepley   PetscTabulation   *T, *TfIn, *TAux = NULL;
105427f02ce8SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
105527f02ce8SMatthew G. Knepley   const PetscScalar *constants;
105627f02ce8SMatthew G. Knepley   PetscReal         *x;
1057665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
1058665f567fSMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
105945480ffeSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
106007218a29SMatthew G. Knepley   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE;
106127f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
10620502970dSMatthew G. Knepley   PetscInt           qNc, Nq, q;
106327f02ce8SMatthew G. Knepley 
106427f02ce8SMatthew G. Knepley   PetscFunctionBegin;
10659566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
106645480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
106745480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
106827f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
10699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
10709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
10719566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
10729566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
10739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1074429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets
10759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
10769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
107727f02ce8SMatthew G. Knepley   switch (jtype) {
1078d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
1079d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1080d71ae5a4SJacob Faibussowitsch     break;
1081d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
1082d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1083d71ae5a4SJacob Faibussowitsch     break;
1084d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
1085d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
108627f02ce8SMatthew G. Knepley   }
10873ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
10889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
10899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
10909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
10919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
109207218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
10939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
10949566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
1095*87510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
10969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
109727f02ce8SMatthew G. Knepley   if (dsAux) {
10989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
10999566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
11009566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
11019566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
11029566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
11039566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
110401907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
11059566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
11069566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
110763a3b9bcSJacob Faibussowitsch     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);
110827f02ce8SMatthew G. Knepley   }
11099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
11109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1111665f567fSMatthew G. Knepley   NcI = T[fieldI]->Nc;
1112665f567fSMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
111327f02ce8SMatthew G. Knepley   NcS = isCohesiveFieldI ? NcI : 2 * NcI;
111427f02ce8SMatthew G. Knepley   NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ;
11150abb75b6SMatthew G. Knepley   if (!isCohesiveFieldI && s == 2) {
11160abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
11170abb75b6SMatthew G. Knepley     NcS *= 2;
11180abb75b6SMatthew G. Knepley   }
11190abb75b6SMatthew G. Knepley   if (!isCohesiveFieldJ && s == 2) {
11200abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
11210abb75b6SMatthew G. Knepley     NcT *= 2;
11220abb75b6SMatthew G. Knepley   }
11230502970dSMatthew G. Knepley   // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though
11240502970dSMatthew G. Knepley   // the coordinates are in dE dimensions
11259566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS * NcT));
11260502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
11270502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
11280502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
11299566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
11300e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
113163a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
113227f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
113327f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
11340e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
11350e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
11364e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
113727f02ce8SMatthew G. Knepley 
113807218a29SMatthew G. Knepley     fegeom.v = x; /* Workspace */
113927f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
11400e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
114127f02ce8SMatthew G. Knepley       PetscReal w;
114227f02ce8SMatthew G. Knepley       PetscInt  c;
114327f02ce8SMatthew G. Knepley 
11444e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0]));
11454e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1]));
114607218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
114727f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
114807218a29SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
114963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
115027f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
11519566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
115227f02ce8SMatthew G. Knepley #endif
115327f02ce8SMatthew G. Knepley       }
115463a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
11558e3a54c0SPierre Jolivet       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
115607218a29SMatthew G. Knepley       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1157ea672e62SMatthew G. Knepley       if (n0) {
11589566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS * NcT));
1159148442b3SMatthew G. Knepley         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);
116027f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT; ++c) g0[c] *= w;
116127f02ce8SMatthew G. Knepley       }
1162ea672e62SMatthew G. Knepley       if (n1) {
11630502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
1164148442b3SMatthew G. Knepley         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);
11650502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w;
116627f02ce8SMatthew G. Knepley       }
1167ea672e62SMatthew G. Knepley       if (n2) {
11680502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
1169148442b3SMatthew G. Knepley         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);
11700502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w;
117127f02ce8SMatthew G. Knepley       }
1172ea672e62SMatthew G. Knepley       if (n3) {
11730502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
1174148442b3SMatthew G. Knepley         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);
11750502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w;
117627f02ce8SMatthew G. Knepley       }
117727f02ce8SMatthew G. Knepley 
11785fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
11795fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
11809566063dSJacob Faibussowitsch           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));
118127f02ce8SMatthew G. Knepley         } else {
11820abb75b6SMatthew G. Knepley           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
11830abb75b6SMatthew G. Knepley           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 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));
11840abb75b6SMatthew G. Knepley         }
11850abb75b6SMatthew G. Knepley       } else {
11860abb75b6SMatthew G. Knepley         if (s == 2) {
11870abb75b6SMatthew G. Knepley           if (isCohesiveFieldJ) {
11880abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
11890abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 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));
11900abb75b6SMatthew G. Knepley           } else {
11910abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
11920abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 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));
11930abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 2], &g1[NcI * NcJ * dim * 2], &g2[NcI * NcJ * dim * 2], &g3[NcI * NcJ * dim * dim * 2], eOffset, totDim, offsetI, offsetJ, elemMat));
11940abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 3], &g1[NcI * NcJ * dim * 3], &g2[NcI * NcJ * dim * 3], &g3[NcI * NcJ * dim * dim * 3], eOffset, totDim, offsetI, offsetJ, elemMat));
11955fedec97SMatthew G. Knepley           }
11969371c9d4SSatish Balay         } else
11970abb75b6SMatthew G. Knepley           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
11980abb75b6SMatthew G. Knepley       }
119927f02ce8SMatthew G. Knepley     }
120027f02ce8SMatthew G. Knepley     if (debug > 1) {
12014e913f38SMatthew G. Knepley       const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb));
12024e913f38SMatthew G. Knepley       const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb));
12034e913f38SMatthew G. Knepley       const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb));
12044e913f38SMatthew G. Knepley       const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb));
12054e913f38SMatthew G. Knepley       PetscInt       f, g;
120627f02ce8SMatthew G. Knepley 
12074e913f38SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT " s %s totDim %" PetscInt_FMT " offsets (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", fieldI, fieldJ, s ? (s > 1 ? "Coh" : "Pos") : "Neg", totDim, eOffset, offsetI, offsetJ));
12084e913f38SMatthew G. Knepley       for (f = fS; f < fE; ++f) {
12094e913f38SMatthew G. Knepley         const PetscInt i = offsetI + f;
12104e913f38SMatthew G. Knepley         for (g = gS; g < gE; ++g) {
12114e913f38SMatthew G. Knepley           const PetscInt j = offsetJ + g;
1212e978a55eSPierre Jolivet           PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, f, i, g, j);
12134e913f38SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "    elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f / NcI, f % NcI, g / NcJ, g % NcJ, (double)PetscRealPart(elemMat[eOffset + i * totDim + j])));
121427f02ce8SMatthew G. Knepley         }
12159566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
121627f02ce8SMatthew G. Knepley       }
121727f02ce8SMatthew G. Knepley     }
121827f02ce8SMatthew G. Knepley     cOffset += totDim;
121927f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
122027f02ce8SMatthew G. Knepley     eOffset += PetscSqr(totDim);
122127f02ce8SMatthew G. Knepley   }
12223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122327f02ce8SMatthew G. Knepley }
122427f02ce8SMatthew G. Knepley 
1225d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1226d71ae5a4SJacob Faibussowitsch {
122720cf1dd8SToby Isaac   PetscFunctionBegin;
122820cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
122920cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
123020cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
123120cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
123220cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1233ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
123420cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1235afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
123620cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
123720cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
123827f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
123920cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
124020cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
124120cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
124227f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
12433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124420cf1dd8SToby Isaac }
124520cf1dd8SToby Isaac 
124620cf1dd8SToby Isaac /*MC
1247dce8aebaSBarry Smith   PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization
124820cf1dd8SToby Isaac 
124920cf1dd8SToby Isaac   Level: intermediate
125020cf1dd8SToby Isaac 
1251dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
125220cf1dd8SToby Isaac M*/
125320cf1dd8SToby Isaac 
1254d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1255d71ae5a4SJacob Faibussowitsch {
125620cf1dd8SToby Isaac   PetscFE_Basic *b;
125720cf1dd8SToby Isaac 
125820cf1dd8SToby Isaac   PetscFunctionBegin;
125920cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
12604dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
126120cf1dd8SToby Isaac   fem->data = b;
126220cf1dd8SToby Isaac 
12639566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
12643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126520cf1dd8SToby Isaac }
1266