xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 0abb75b6235d24c021c3aa34ce9bcf81918e8852)
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 {
16420cf1dd8SToby Isaac   const PetscInt     debug = 0;
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;
17120cf1dd8SToby Isaac   PetscReal         *x;
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));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
1849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
1879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
1889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
1899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
1909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1924bee2e38SMatthew G. Knepley   if (dsAux) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
1989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
19963a3b9bcSJacob 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);
20020cf1dd8SToby Isaac   }
2019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
20263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
20320cf1dd8SToby Isaac   Np       = cgeom->numPoints;
20420cf1dd8SToby Isaac   dE       = cgeom->dimEmbed;
20520cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
20620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2074bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
20820cf1dd8SToby Isaac 
20927f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
21027f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
21120cf1dd8SToby Isaac     if (isAffine) {
2124bee2e38SMatthew G. Knepley       fegeom.v    = x;
2134bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2147132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
2157132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
2167132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
21720cf1dd8SToby Isaac     }
2184bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
2194bee2e38SMatthew G. Knepley       PetscScalar integrand;
2204bee2e38SMatthew G. Knepley       PetscReal   w;
2214bee2e38SMatthew G. Knepley 
2224bee2e38SMatthew G. Knepley       if (isAffine) {
2237132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
2244bee2e38SMatthew G. Knepley       } else {
2254bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
2264bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
2274bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
2284bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
2294bee2e38SMatthew G. Knepley       }
2304bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
23120cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
23263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
2337be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2349566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
23520cf1dd8SToby Isaac #endif
23620cf1dd8SToby Isaac       }
23763a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
2389566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
2399566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
2404bee2e38SMatthew 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);
2414bee2e38SMatthew G. Knepley       integrand *= w;
24220cf1dd8SToby Isaac       integral[e * Nf + field] += integrand;
2439566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[field])));
24420cf1dd8SToby Isaac     }
24520cf1dd8SToby Isaac     cOffset += totDim;
24620cf1dd8SToby Isaac     cOffsetAux += totDimAux;
24720cf1dd8SToby Isaac   }
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24920cf1dd8SToby Isaac }
25020cf1dd8SToby Isaac 
2512dce792eSToby 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[])
252d71ae5a4SJacob Faibussowitsch {
253afe6d6adSToby Isaac   const PetscInt     debug = 0;
2544bee2e38SMatthew G. Knepley   PetscFE            fe;
255afe6d6adSToby Isaac   PetscQuadrature    quad;
256ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2574bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
258afe6d6adSToby Isaac   const PetscScalar *constants;
259afe6d6adSToby Isaac   PetscReal         *x;
260ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
261afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
262afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
263afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
264afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
265afe6d6adSToby Isaac 
266afe6d6adSToby Isaac   PetscFunctionBegin;
2673ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
2689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
2719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
2749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
2759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
2769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
2779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
2789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
2794bee2e38SMatthew G. Knepley   if (dsAux) {
2809566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
2819566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
2839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
286afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
2879566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
2889566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
28963a3b9bcSJacob 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);
290afe6d6adSToby Isaac   }
2919566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
29263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
293afe6d6adSToby Isaac   Np       = fgeom->numPoints;
294afe6d6adSToby Isaac   dE       = fgeom->dimEmbed;
295afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
296afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
2979f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
298afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
299ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
300ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
301ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
302ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
30327f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
30427f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
30527f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
30627f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
3074bee2e38SMatthew G. Knepley     if (isAffine) {
3084bee2e38SMatthew G. Knepley       fegeom.v    = x;
3094bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3107132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
3117132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
3127132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
3137132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
3149f209ee4SMatthew G. Knepley 
3157132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
3167132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
3177132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
3184bee2e38SMatthew G. Knepley     }
319afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
320afe6d6adSToby Isaac       PetscScalar integrand;
3214bee2e38SMatthew G. Knepley       PetscReal   w;
322afe6d6adSToby Isaac 
323afe6d6adSToby Isaac       if (isAffine) {
3247132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
325afe6d6adSToby Isaac       } else {
3263fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
3279f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
3289f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
3294bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
3304bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
3319f209ee4SMatthew G. Knepley 
3329f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
3339f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
3349f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
335afe6d6adSToby Isaac       }
3364bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
337afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
33863a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
339afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3409566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
341afe6d6adSToby Isaac #endif
342afe6d6adSToby Isaac       }
34363a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
3449566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3459566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3464bee2e38SMatthew 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);
3474bee2e38SMatthew G. Knepley       integrand *= w;
348afe6d6adSToby Isaac       integral[e * Nf + field] += integrand;
3499566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
350afe6d6adSToby Isaac     }
351afe6d6adSToby Isaac     cOffset += totDim;
352afe6d6adSToby Isaac     cOffsetAux += totDimAux;
353afe6d6adSToby Isaac   }
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355afe6d6adSToby Isaac }
356afe6d6adSToby Isaac 
357d71ae5a4SJacob 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[])
358d71ae5a4SJacob Faibussowitsch {
35920cf1dd8SToby Isaac   const PetscInt     debug = 0;
3606528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3614bee2e38SMatthew G. Knepley   PetscFE            fe;
3626528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3636528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
3646528b96dSMatthew G. Knepley   PetscPointFunc    *f0_func, *f1_func;
36520cf1dd8SToby Isaac   PetscQuadrature    quad;
366ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3674bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
36820cf1dd8SToby Isaac   const PetscScalar *constants;
36920cf1dd8SToby Isaac   PetscReal         *x;
370ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
371ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
37220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3736587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac   PetscFunctionBegin;
3769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
3779566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
3789566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
3799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
3819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
3849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
3859566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
3863ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
3879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
3889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
3899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
3909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
3919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3924bee2e38SMatthew G. Knepley   if (dsAux) {
3939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
3959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
3989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
39963a3b9bcSJacob 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);
40020cf1dd8SToby Isaac   }
4019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
40263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
40320cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
40463a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
40520cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4064bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
40720cf1dd8SToby Isaac 
4086587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
4099566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
4109566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
41120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4124bee2e38SMatthew G. Knepley       PetscReal w;
4134bee2e38SMatthew G. Knepley       PetscInt  c, d;
41420cf1dd8SToby Isaac 
4159566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
4164bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
4176587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
41863a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4197be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4209566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
42120cf1dd8SToby Isaac #endif
42220cf1dd8SToby Isaac       }
4239566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
4249566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
4256528b96dSMatthew 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]);
426ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
4276528b96dSMatthew 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]);
4289371c9d4SSatish Balay       for (c = 0; c < T[field]->Nc; ++c)
4299371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w;
430b8025e53SMatthew G. Knepley       if (debug) {
431e8e188d2SZach Atkins         // LCOV_EXCL_START
432e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q]));
433e8e188d2SZach Atkins         for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c]));
434e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
435b8025e53SMatthew G. Knepley         if (debug > 2) {
43663a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
43763a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
4389566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
439e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field der %" PetscInt_FMT ":", field));
440e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c])));
441e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
44263a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
44363a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
4449566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
445e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  res der %" PetscInt_FMT ":", field));
446e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc; ++c) {
447e8e188d2SZach Atkins             for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dim + d])));
448b8025e53SMatthew G. Knepley           }
449e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
450e8e188d2SZach Atkins         }
451e8e188d2SZach Atkins         // LCOV_EXCL_STOP
452b8025e53SMatthew G. Knepley       }
45320cf1dd8SToby Isaac     }
4549566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
45520cf1dd8SToby Isaac     cOffset += totDim;
45620cf1dd8SToby Isaac     cOffsetAux += totDimAux;
45720cf1dd8SToby Isaac   }
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45920cf1dd8SToby Isaac }
46020cf1dd8SToby Isaac 
461d71ae5a4SJacob 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[])
462d71ae5a4SJacob Faibussowitsch {
46320cf1dd8SToby Isaac   const PetscInt     debug = 0;
46406d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
4654bee2e38SMatthew G. Knepley   PetscFE            fe;
46606d8a0d3SMatthew G. Knepley   PetscInt           n0, n1, i;
46706d8a0d3SMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
46820cf1dd8SToby Isaac   PetscQuadrature    quad;
469ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
4704bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
47120cf1dd8SToby Isaac   const PetscScalar *constants;
47220cf1dd8SToby Isaac   PetscReal         *x;
473ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
474ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
4756587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
47620cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
4776587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
47820cf1dd8SToby Isaac 
47920cf1dd8SToby Isaac   PetscFunctionBegin;
4809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
4819566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
4829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
4839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4889566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
4893ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
4909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
4949566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4954bee2e38SMatthew G. Knepley   if (dsAux) {
4969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
4979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4999566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
5009566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
5019566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
5027be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
5039566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
5049566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
50563a3b9bcSJacob 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);
50620cf1dd8SToby Isaac   }
507ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
5089566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
50963a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
51020cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
5116587ee25SMatthew G. Knepley   /* TODO FIX THIS */
5126587ee25SMatthew G. Knepley   fgeom->dim = dim - 1;
51363a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
51420cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5159f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
51620cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5179f209ee4SMatthew G. Knepley 
5186587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
5199566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcI));
5209566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
52120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5224bee2e38SMatthew G. Knepley       PetscReal w;
5234bee2e38SMatthew G. Knepley       PetscInt  c, d;
5244bee2e38SMatthew G. Knepley 
5259566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
5269566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
5274bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
52862bd480fSMatthew G. Knepley       if (debug > 1) {
5296587ee25SMatthew G. Knepley         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
53063a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5317be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5329566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5339566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
53420cf1dd8SToby Isaac #endif
53520cf1dd8SToby Isaac         }
53662bd480fSMatthew G. Knepley       }
5379566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
5389566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
53906d8a0d3SMatthew 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]);
5404bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w;
54106d8a0d3SMatthew 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]);
5429371c9d4SSatish Balay       for (c = 0; c < NcI; ++c)
5439371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w;
54462bd480fSMatthew G. Knepley       if (debug) {
54563a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
54662bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
54763a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c])));
54862bd480fSMatthew G. Knepley           if (n1) {
54963a3b9bcSJacob 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])));
5509566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
55162bd480fSMatthew G. Knepley           }
55262bd480fSMatthew G. Knepley         }
55362bd480fSMatthew G. Knepley       }
55420cf1dd8SToby Isaac     }
5559566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
55620cf1dd8SToby Isaac     cOffset += totDim;
55720cf1dd8SToby Isaac     cOffsetAux += totDimAux;
55820cf1dd8SToby Isaac   }
5593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56020cf1dd8SToby Isaac }
56120cf1dd8SToby Isaac 
56227f02ce8SMatthew G. Knepley /*
56327f02ce8SMatthew 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
56427f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
56527f02ce8SMatthew G. Knepley 
56627f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
56727f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
56827f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
56927f02ce8SMatthew 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()
57027f02ce8SMatthew G. Knepley */
5712dce792eSToby 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[])
572d71ae5a4SJacob Faibussowitsch {
57327f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
5746528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
57527f02ce8SMatthew G. Knepley   PetscFE            fe;
5766528b96dSMatthew G. Knepley   PetscWeakForm      wf;
5776528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
5786528b96dSMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
57927f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
5800e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
58107218a29SMatthew G. Knepley   PetscTabulation   *Tf, *TfIn, *TfAux = NULL;
58227f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
58327f02ce8SMatthew G. Knepley   const PetscScalar *constants;
58427f02ce8SMatthew G. Knepley   PetscReal         *x;
585665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
58607218a29SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
5876587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
58827f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
5896587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
59027f02ce8SMatthew G. Knepley 
59127f02ce8SMatthew G. Knepley   PetscFunctionBegin;
59227f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
5939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
5949566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
5959566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
5969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5979566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
59807218a29SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
599429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets
60007218a29SMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x));
6019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
6029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
6039566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
6043ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
6059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
6069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
6079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
60827f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
6099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
61007218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
6119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
61227f02ce8SMatthew G. Knepley   if (dsAux) {
6139566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
6149566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6159566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6169566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6179566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
61901907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
6209566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
6219566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
62263a3b9bcSJacob 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);
62327f02ce8SMatthew G. Knepley   }
6249566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
625665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
626c2b7495fSMatthew G. Knepley   NcS = NcI;
627*0abb75b6SMatthew G. Knepley   if (!isCohesiveField && s == 2) {
628*0abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
629*0abb75b6SMatthew G. Knepley     NcS *= 2;
630*0abb75b6SMatthew G. Knepley   }
6319566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
6320e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
63363a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
63427f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
63563a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
63627f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
63727f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
6380e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
6390e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
6404e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
64127f02ce8SMatthew G. Knepley 
6426587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6439566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcS));
6449566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
64527f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
6460e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
64727f02ce8SMatthew G. Knepley       PetscReal w;
64827f02ce8SMatthew G. Knepley       PetscInt  c, d;
64927f02ce8SMatthew G. Knepley 
6504e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0]));
6514e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1]));
65207218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
65327f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
6546587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
65563a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
65627f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
6579566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
65827f02ce8SMatthew G. Knepley #endif
65927f02ce8SMatthew G. Knepley       }
660a4158a15SMatthew 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]));
66127f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
66207218a29SMatthew G. Knepley       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], &coefficients_t[cOffsetIn], u, u_x, u_t));
66307218a29SMatthew 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));
6646528b96dSMatthew 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]);
66527f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
6669ee2af8cSMatthew 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]);
6679371c9d4SSatish Balay       for (c = 0; c < NcS; ++c)
6689371c9d4SSatish Balay         for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
66927f02ce8SMatthew G. Knepley     }
6709371c9d4SSatish Balay     if (isCohesiveField) {
6713ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
6729371c9d4SSatish Balay     } else {
6733ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
6749371c9d4SSatish Balay     }
67527f02ce8SMatthew G. Knepley     cOffset += totDim;
67607218a29SMatthew G. Knepley     cOffsetIn += totDimIn;
67727f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
67827f02ce8SMatthew G. Knepley   }
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
68027f02ce8SMatthew G. Knepley }
68127f02ce8SMatthew G. Knepley 
682d71ae5a4SJacob 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[])
683d71ae5a4SJacob Faibussowitsch {
68420cf1dd8SToby Isaac   const PetscInt     debug = 0;
6854bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
6866528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6876528b96dSMatthew G. Knepley   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
6886528b96dSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
68920cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
69020cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
69120cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
69220cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
69320cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
69420cf1dd8SToby Isaac   PetscQuadrature    quad;
695ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
6964bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
69720cf1dd8SToby Isaac   const PetscScalar *constants;
69820cf1dd8SToby Isaac   PetscReal         *x;
699ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
700ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
7016528b96dSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
70220cf1dd8SToby Isaac   PetscInt           dE, Np;
70320cf1dd8SToby Isaac   PetscBool          isAffine;
70420cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
70520cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
70620cf1dd8SToby Isaac 
70720cf1dd8SToby Isaac   PetscFunctionBegin;
7089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
7096528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
7106528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
7119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
7129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
7139566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
7149566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
7159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
7169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
7179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
7189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
71920cf1dd8SToby Isaac   switch (jtype) {
720d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
721d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
722d71ae5a4SJacob Faibussowitsch     break;
723d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
724d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
725d71ae5a4SJacob Faibussowitsch     break;
726d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
727d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
728d71ae5a4SJacob Faibussowitsch     break;
72920cf1dd8SToby Isaac   }
7303ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
7319566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
7329566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
7339566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
7349566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
7359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
7369566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
7379566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
7384bee2e38SMatthew G. Knepley   if (dsAux) {
7399566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7409566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
7419566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
7429566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
7439566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
7449566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
74563a3b9bcSJacob 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);
74620cf1dd8SToby Isaac   }
74727f02ce8SMatthew G. Knepley   NcI      = T[fieldI]->Nc;
74827f02ce8SMatthew G. Knepley   NcJ      = T[fieldJ]->Nc;
7494bee2e38SMatthew G. Knepley   Np       = cgeom->numPoints;
7504bee2e38SMatthew G. Knepley   dE       = cgeom->dimEmbed;
7514bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
75227f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
7539566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
7549566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7559566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7569566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7579566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
75863a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
7594bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
7604bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
7614bee2e38SMatthew G. Knepley 
76227f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
76327f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
7644bee2e38SMatthew G. Knepley     if (isAffine) {
7654bee2e38SMatthew G. Knepley       fegeom.v    = x;
7664bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
7677132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
7687132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
7697132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
7704bee2e38SMatthew G. Knepley     }
77120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
77220cf1dd8SToby Isaac       PetscReal w;
7734bee2e38SMatthew G. Knepley       PetscInt  c;
77420cf1dd8SToby Isaac 
77520cf1dd8SToby Isaac       if (isAffine) {
7767132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
77720cf1dd8SToby Isaac       } else {
7784bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
7794bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
7804bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
7814bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
78220cf1dd8SToby Isaac       }
7839566063dSJacob 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]));
7844bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
7859566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
7869566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
787ea672e62SMatthew G. Knepley       if (n0) {
7889566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
7896528b96dSMatthew 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);
79020cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
79120cf1dd8SToby Isaac       }
792ea672e62SMatthew G. Knepley       if (n1) {
7939566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7946528b96dSMatthew 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);
7954bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
79620cf1dd8SToby Isaac       }
797ea672e62SMatthew G. Knepley       if (n2) {
7989566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7996528b96dSMatthew 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);
8004bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
80120cf1dd8SToby Isaac       }
802ea672e62SMatthew G. Knepley       if (n3) {
8039566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8046528b96dSMatthew 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);
8054bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
80620cf1dd8SToby Isaac       }
80720cf1dd8SToby Isaac 
8089566063dSJacob 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));
80920cf1dd8SToby Isaac     }
81020cf1dd8SToby Isaac     if (debug > 1) {
81120cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
81220cf1dd8SToby Isaac 
81363a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
814ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
815ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
816ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
817ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
818ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
819ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
82063a3b9bcSJacob 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])));
82120cf1dd8SToby Isaac             }
82220cf1dd8SToby Isaac           }
8239566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
82420cf1dd8SToby Isaac         }
82520cf1dd8SToby Isaac       }
82620cf1dd8SToby Isaac     }
82720cf1dd8SToby Isaac     cOffset += totDim;
82820cf1dd8SToby Isaac     cOffsetAux += totDimAux;
82920cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
83020cf1dd8SToby Isaac   }
8313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83220cf1dd8SToby Isaac }
83320cf1dd8SToby Isaac 
8342dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
835d71ae5a4SJacob Faibussowitsch {
83620cf1dd8SToby Isaac   const PetscInt     debug = 0;
8374bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
83845480ffeSMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
83945480ffeSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
84020cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
84120cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
84220cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
84320cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
84420cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
84520cf1dd8SToby Isaac   PetscQuadrature    quad;
846ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
8474bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
84820cf1dd8SToby Isaac   const PetscScalar *constants;
84920cf1dd8SToby Isaac   PetscReal         *x;
850ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
851ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
85245480ffeSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
85320cf1dd8SToby Isaac   PetscBool          isAffine;
85420cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
85520cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
85620cf1dd8SToby Isaac 
85720cf1dd8SToby Isaac   PetscFunctionBegin;
8589566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
85945480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
86045480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
8619566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
8629566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
8639566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
8649566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
8659566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
8669566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
8679566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
8689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
8699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
8709566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
8713ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
8729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
8739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
8749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
8759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
8769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
8774bee2e38SMatthew G. Knepley   if (dsAux) {
8789566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8799566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
8809566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
8819566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
8829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
8839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
88420cf1dd8SToby Isaac   }
885ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
88620cf1dd8SToby Isaac   Np       = fgeom->numPoints;
88720cf1dd8SToby Isaac   dE       = fgeom->dimEmbed;
88820cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
88927f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
8909566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
8919566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
8929566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
8939566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8949566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
89563a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
89620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
8979f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
89820cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
899ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
900ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
901ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
902ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
90327f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
90427f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
90527f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
90627f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
9074bee2e38SMatthew G. Knepley     if (isAffine) {
9084bee2e38SMatthew G. Knepley       fegeom.v    = x;
9094bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
9107132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
9117132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
9127132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
9137132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
9149f209ee4SMatthew G. Knepley 
9157132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
9167132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
9177132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
9184bee2e38SMatthew G. Knepley     }
91920cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
92020cf1dd8SToby Isaac       PetscReal w;
9214bee2e38SMatthew G. Knepley       PetscInt  c;
92220cf1dd8SToby Isaac 
92363a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
92420cf1dd8SToby Isaac       if (isAffine) {
9257132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
92620cf1dd8SToby Isaac       } else {
9273fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
9289f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
9299f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
9304bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
9314bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
9329f209ee4SMatthew G. Knepley 
9339f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
9349f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
9359f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
93620cf1dd8SToby Isaac       }
9374bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
9389566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
9399566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
940ea672e62SMatthew G. Knepley       if (n0) {
9419566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
94245480ffeSMatthew 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);
94320cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
94420cf1dd8SToby Isaac       }
945ea672e62SMatthew G. Knepley       if (n1) {
9469566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
94745480ffeSMatthew 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);
9484bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
94920cf1dd8SToby Isaac       }
950ea672e62SMatthew G. Knepley       if (n2) {
9519566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
95245480ffeSMatthew 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);
9534bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
95420cf1dd8SToby Isaac       }
955ea672e62SMatthew G. Knepley       if (n3) {
9569566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
95745480ffeSMatthew 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);
9584bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
95920cf1dd8SToby Isaac       }
96020cf1dd8SToby Isaac 
9619566063dSJacob 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));
96220cf1dd8SToby Isaac     }
96320cf1dd8SToby Isaac     if (debug > 1) {
96420cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
96520cf1dd8SToby Isaac 
96663a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
967ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
968ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
969ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
970ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
971ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
972ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
97363a3b9bcSJacob 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])));
97420cf1dd8SToby Isaac             }
97520cf1dd8SToby Isaac           }
9769566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
97720cf1dd8SToby Isaac         }
97820cf1dd8SToby Isaac       }
97920cf1dd8SToby Isaac     }
98020cf1dd8SToby Isaac     cOffset += totDim;
98120cf1dd8SToby Isaac     cOffsetAux += totDimAux;
98220cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
98320cf1dd8SToby Isaac   }
9843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98520cf1dd8SToby Isaac }
98620cf1dd8SToby Isaac 
9872dce792eSToby 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[])
988d71ae5a4SJacob Faibussowitsch {
98927f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
99027f02ce8SMatthew G. Knepley   PetscFE            feI, feJ;
991148442b3SMatthew G. Knepley   PetscWeakForm      wf;
992148442b3SMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
993148442b3SMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
99427f02ce8SMatthew G. Knepley   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
99527f02ce8SMatthew G. Knepley   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
99627f02ce8SMatthew G. Knepley   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
99727f02ce8SMatthew G. Knepley   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
99827f02ce8SMatthew G. Knepley   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
999665f567fSMatthew G. Knepley   PetscQuadrature    quad;
10000e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
100107218a29SMatthew G. Knepley   PetscTabulation   *T, *TfIn, *TAux = NULL;
100227f02ce8SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
100327f02ce8SMatthew G. Knepley   const PetscScalar *constants;
100427f02ce8SMatthew G. Knepley   PetscReal         *x;
1005665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
1006665f567fSMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
100745480ffeSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
100807218a29SMatthew G. Knepley   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE;
100927f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
10100502970dSMatthew G. Knepley   PetscInt           qNc, Nq, q;
101127f02ce8SMatthew G. Knepley 
101227f02ce8SMatthew G. Knepley   PetscFunctionBegin;
10139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
101445480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
101545480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
101627f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
10179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
10189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
10199566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
10209566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
10219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1022429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets
10239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
10249566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
102527f02ce8SMatthew G. Knepley   switch (jtype) {
1026d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
1027d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1028d71ae5a4SJacob Faibussowitsch     break;
1029d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
1030d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1031d71ae5a4SJacob Faibussowitsch     break;
1032d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
1033d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
103427f02ce8SMatthew G. Knepley   }
10353ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
10369566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
10379566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
10389566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
10399566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
104007218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
10419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
10429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
10439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
104427f02ce8SMatthew G. Knepley   if (dsAux) {
10459566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
10469566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
10479566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
10489566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
10499566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
10509566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
105101907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
10529566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
10539566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
105463a3b9bcSJacob 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);
105527f02ce8SMatthew G. Knepley   }
10569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
10579566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1058665f567fSMatthew G. Knepley   NcI = T[fieldI]->Nc;
1059665f567fSMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
106027f02ce8SMatthew G. Knepley   NcS = isCohesiveFieldI ? NcI : 2 * NcI;
106127f02ce8SMatthew G. Knepley   NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ;
1062*0abb75b6SMatthew G. Knepley   if (!isCohesiveFieldI && s == 2) {
1063*0abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
1064*0abb75b6SMatthew G. Knepley     NcS *= 2;
1065*0abb75b6SMatthew G. Knepley   }
1066*0abb75b6SMatthew G. Knepley   if (!isCohesiveFieldJ && s == 2) {
1067*0abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
1068*0abb75b6SMatthew G. Knepley     NcT *= 2;
1069*0abb75b6SMatthew G. Knepley   }
10700502970dSMatthew G. Knepley   // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though
10710502970dSMatthew G. Knepley   // the coordinates are in dE dimensions
10729566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS * NcT));
10730502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
10740502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
10750502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
10769566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
10770e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
107863a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
107927f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
108027f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
10810e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
10820e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
10834e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
108427f02ce8SMatthew G. Knepley 
108507218a29SMatthew G. Knepley     fegeom.v = x; /* Workspace */
108627f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
10870e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
108827f02ce8SMatthew G. Knepley       PetscReal w;
108927f02ce8SMatthew G. Knepley       PetscInt  c;
109027f02ce8SMatthew G. Knepley 
10914e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0]));
10924e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1]));
109307218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
109427f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
109507218a29SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
109663a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
109727f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
10989566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
109927f02ce8SMatthew G. Knepley #endif
110027f02ce8SMatthew G. Knepley       }
110163a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
110207218a29SMatthew G. Knepley       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
110307218a29SMatthew 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));
1104ea672e62SMatthew G. Knepley       if (n0) {
11059566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS * NcT));
1106148442b3SMatthew 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);
110727f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT; ++c) g0[c] *= w;
110827f02ce8SMatthew G. Knepley       }
1109ea672e62SMatthew G. Knepley       if (n1) {
11100502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
1111148442b3SMatthew 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);
11120502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w;
111327f02ce8SMatthew G. Knepley       }
1114ea672e62SMatthew G. Knepley       if (n2) {
11150502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
1116148442b3SMatthew 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);
11170502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w;
111827f02ce8SMatthew G. Knepley       }
1119ea672e62SMatthew G. Knepley       if (n3) {
11200502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
1121148442b3SMatthew 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);
11220502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w;
112327f02ce8SMatthew G. Knepley       }
112427f02ce8SMatthew G. Knepley 
11255fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
11265fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
11279566063dSJacob 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));
112827f02ce8SMatthew G. Knepley         } else {
1129*0abb75b6SMatthew 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));
1130*0abb75b6SMatthew 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));
1131*0abb75b6SMatthew G. Knepley         }
1132*0abb75b6SMatthew G. Knepley       } else {
1133*0abb75b6SMatthew G. Knepley         if (s == 2) {
1134*0abb75b6SMatthew G. Knepley           if (isCohesiveFieldJ) {
1135*0abb75b6SMatthew 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));
1136*0abb75b6SMatthew 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));
1137*0abb75b6SMatthew G. Knepley           } else {
1138*0abb75b6SMatthew 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));
1139*0abb75b6SMatthew 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));
1140*0abb75b6SMatthew 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));
1141*0abb75b6SMatthew 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));
11425fedec97SMatthew G. Knepley           }
11439371c9d4SSatish Balay         } else
1144*0abb75b6SMatthew 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));
1145*0abb75b6SMatthew G. Knepley       }
114627f02ce8SMatthew G. Knepley     }
114727f02ce8SMatthew G. Knepley     if (debug > 1) {
11484e913f38SMatthew G. Knepley       const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb));
11494e913f38SMatthew G. Knepley       const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb));
11504e913f38SMatthew G. Knepley       const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb));
11514e913f38SMatthew G. Knepley       const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb));
11524e913f38SMatthew G. Knepley       PetscInt       f, g;
115327f02ce8SMatthew G. Knepley 
11544e913f38SMatthew 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));
11554e913f38SMatthew G. Knepley       for (f = fS; f < fE; ++f) {
11564e913f38SMatthew G. Knepley         const PetscInt i = offsetI + f;
11574e913f38SMatthew G. Knepley         for (g = gS; g < gE; ++g) {
11584e913f38SMatthew G. Knepley           const PetscInt j = offsetJ + g;
11594e913f38SMatthew G. Knepley           PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", f, i, g, j);
11604e913f38SMatthew 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])));
116127f02ce8SMatthew G. Knepley         }
11629566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
116327f02ce8SMatthew G. Knepley       }
116427f02ce8SMatthew G. Knepley     }
116527f02ce8SMatthew G. Knepley     cOffset += totDim;
116627f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
116727f02ce8SMatthew G. Knepley     eOffset += PetscSqr(totDim);
116827f02ce8SMatthew G. Knepley   }
11693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117027f02ce8SMatthew G. Knepley }
117127f02ce8SMatthew G. Knepley 
1172d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1173d71ae5a4SJacob Faibussowitsch {
117420cf1dd8SToby Isaac   PetscFunctionBegin;
117520cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
117620cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
117720cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
117820cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
117920cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1180ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
118120cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1182afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
118320cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
118420cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
118527f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
118620cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
118720cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
118820cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
118927f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
11903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119120cf1dd8SToby Isaac }
119220cf1dd8SToby Isaac 
119320cf1dd8SToby Isaac /*MC
1194dce8aebaSBarry Smith   PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization
119520cf1dd8SToby Isaac 
119620cf1dd8SToby Isaac   Level: intermediate
119720cf1dd8SToby Isaac 
1198dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
119920cf1dd8SToby Isaac M*/
120020cf1dd8SToby Isaac 
1201d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1202d71ae5a4SJacob Faibussowitsch {
120320cf1dd8SToby Isaac   PetscFE_Basic *b;
120420cf1dd8SToby Isaac 
120520cf1dd8SToby Isaac   PetscFunctionBegin;
120620cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
12074dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
120820cf1dd8SToby Isaac   fem->data = b;
120920cf1dd8SToby Isaac 
12109566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
12113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121220cf1dd8SToby Isaac }
1213