xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 79ab67a34e24ffe79d2fdfccb4206de89a5e4b2d)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
5d71ae5a4SJacob Faibussowitsch {
620cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *)fem->data;
720cf1dd8SToby Isaac 
820cf1dd8SToby Isaac   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1120cf1dd8SToby Isaac }
1220cf1dd8SToby Isaac 
13d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
14d71ae5a4SJacob Faibussowitsch {
15d9bac1caSLisandro Dalcin   PetscInt        dim, Nc;
16d9bac1caSLisandro Dalcin   PetscSpace      basis = NULL;
17d9bac1caSLisandro Dalcin   PetscDualSpace  dual  = NULL;
18d9bac1caSLisandro Dalcin   PetscQuadrature quad  = NULL;
1920cf1dd8SToby Isaac 
2020cf1dd8SToby Isaac   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe, &Nc));
239566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &basis));
249566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dual));
259566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(v));
2763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc));
289566063dSJacob Faibussowitsch   if (basis) PetscCall(PetscSpaceView(basis, v));
299566063dSJacob Faibussowitsch   if (dual) PetscCall(PetscDualSpaceView(dual, v));
309566063dSJacob Faibussowitsch   if (quad) PetscCall(PetscQuadratureView(quad, v));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(v));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
36d71ae5a4SJacob Faibussowitsch {
3720cf1dd8SToby Isaac   PetscBool iascii;
3820cf1dd8SToby Isaac 
3920cf1dd8SToby Isaac   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
419566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscFEView_Basic_Ascii(fe, v));
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4320cf1dd8SToby Isaac }
4420cf1dd8SToby Isaac 
4520cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
46d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
47d71ae5a4SJacob Faibussowitsch {
48b9d4cb8dSJed Brown   PetscReal    *work;
4920cf1dd8SToby Isaac   PetscBLASInt *pivots;
5020cf1dd8SToby Isaac   PetscBLASInt  n, info;
5120cf1dd8SToby Isaac   PetscInt      pdim, j;
5220cf1dd8SToby Isaac 
5320cf1dd8SToby Isaac   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pdim * pdim, &fem->invV));
5620cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
5720cf1dd8SToby Isaac     PetscReal       *Bf;
5820cf1dd8SToby Isaac     PetscQuadrature  f;
5920cf1dd8SToby Isaac     const PetscReal *points, *weights;
6020cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
6120cf1dd8SToby Isaac 
629566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f));
639566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf));
659566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
6620cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
6720cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
68b9d4cb8dSJed Brown       fem->invV[j * pdim + k] = 0.0;
6920cf1dd8SToby Isaac 
7020cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
71b9d4cb8dSJed Brown         for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c];
7220cf1dd8SToby Isaac       }
7320cf1dd8SToby Isaac     }
749566063dSJacob Faibussowitsch     PetscCall(PetscFree(Bf));
7520cf1dd8SToby Isaac   }
76ea2bdf6dSBarry Smith 
779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work));
7820cf1dd8SToby Isaac   n = pdim;
79792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
8063a3b9bcSJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info);
81792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
8263a3b9bcSJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)info);
839566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pivots, work));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8520cf1dd8SToby Isaac }
8620cf1dd8SToby Isaac 
87d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
88d71ae5a4SJacob Faibussowitsch {
8920cf1dd8SToby Isaac   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9220cf1dd8SToby Isaac }
9320cf1dd8SToby Isaac 
94b9d4cb8dSJed Brown /* Tensor contraction on the middle index,
95b9d4cb8dSJed Brown  *    C[m,n,p] = A[m,k,p] * B[k,n]
96b9d4cb8dSJed Brown  * where all matrices use C-style ordering.
97b9d4cb8dSJed Brown  */
98d71ae5a4SJacob Faibussowitsch static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C)
99d71ae5a4SJacob Faibussowitsch {
100b9d4cb8dSJed Brown   PetscInt i;
101b9d4cb8dSJed Brown 
102b9d4cb8dSJed Brown   PetscFunctionBegin;
103aa9788aaSMatthew G. Knepley   PetscCheck(n && p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Empty tensor is not allowed %" PetscInt_FMT " %" PetscInt_FMT, n, p);
104b9d4cb8dSJed Brown   for (i = 0; i < m; i++) {
105b9d4cb8dSJed Brown     PetscBLASInt n_, p_, k_, lda, ldb, ldc;
106b9d4cb8dSJed Brown     PetscReal    one = 1, zero = 0;
107b9d4cb8dSJed Brown     /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
108b9d4cb8dSJed Brown      * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
109b9d4cb8dSJed Brown      */
1109566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &n_));
1119566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(p, &p_));
1129566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(k, &k_));
113b9d4cb8dSJed Brown     lda = p_;
114b9d4cb8dSJed Brown     ldb = n_;
115b9d4cb8dSJed Brown     ldc = p_;
116792fecdfSBarry Smith     PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc));
117b9d4cb8dSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2. * m * n * p * k));
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120b9d4cb8dSJed Brown }
121b9d4cb8dSJed Brown 
122d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
123d71ae5a4SJacob Faibussowitsch {
12420cf1dd8SToby Isaac   DM         dm;
12520cf1dd8SToby Isaac   PetscInt   pdim; /* Dimension of FE space P */
12620cf1dd8SToby Isaac   PetscInt   dim;  /* Spatial dimension */
12720cf1dd8SToby Isaac   PetscInt   Nc;   /* Field components */
128ef0bb6c7SMatthew G. Knepley   PetscReal *B    = K >= 0 ? T->T[0] : NULL;
129ef0bb6c7SMatthew G. Knepley   PetscReal *D    = K >= 1 ? T->T[1] : NULL;
130ef0bb6c7SMatthew G. Knepley   PetscReal *H    = K >= 2 ? T->T[2] : NULL;
131ef0bb6c7SMatthew G. Knepley   PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
13220cf1dd8SToby Isaac 
13320cf1dd8SToby Isaac   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
1359566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1369566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
13820cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
1399566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1409566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1419566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1429566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
143b9d4cb8dSJed Brown   /* Translate from prime to nodal basis */
14420cf1dd8SToby Isaac   if (B) {
145b9d4cb8dSJed Brown     /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
1469566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B));
14720cf1dd8SToby Isaac   }
148aa9788aaSMatthew G. Knepley   if (D && dim) {
149b9d4cb8dSJed Brown     /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
1509566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D));
15120cf1dd8SToby Isaac   }
152aa9788aaSMatthew G. Knepley   if (H && dim) {
153b9d4cb8dSJed Brown     /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
1549566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H));
15520cf1dd8SToby Isaac   }
1569566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1579566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1589566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16020cf1dd8SToby Isaac }
16120cf1dd8SToby Isaac 
1622dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
163d71ae5a4SJacob Faibussowitsch {
164b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
1654bee2e38SMatthew G. Knepley   PetscFE            fe;
16620cf1dd8SToby Isaac   PetscPointFunc     obj_func;
16720cf1dd8SToby Isaac   PetscQuadrature    quad;
168ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
1694bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
17020cf1dd8SToby Isaac   const PetscScalar *constants;
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) {
219d627b919SMatthew G. Knepley       PetscScalar integrand = 0.;
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 {
253b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
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);
293*79ab67a3SMatthew G. Knepley   if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Field: %" PetscInt_FMT " Nface: %" PetscInt_FMT " Nq: %" PetscInt_FMT "\n", field, Ne, Nq));
294afe6d6adSToby Isaac   Np       = fgeom->numPoints;
295afe6d6adSToby Isaac   dE       = fgeom->dimEmbed;
296afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
297afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
2989f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
299afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
300ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
301ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
302ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
303b2deab97SMatthew G. Knepley     fegeom.invJ         = NULL;
304ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
30527f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
30627f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
30727f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
30827f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
3094bee2e38SMatthew G. Knepley     if (isAffine) {
3104bee2e38SMatthew G. Knepley       fegeom.v    = x;
3114bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3127132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
3137132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
3147132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
3157132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
3169f209ee4SMatthew G. Knepley 
3177132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
3187132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
3197132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
3204bee2e38SMatthew G. Knepley     }
321afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
322*79ab67a3SMatthew G. Knepley       PetscScalar integrand = 0.;
3234bee2e38SMatthew G. Knepley       PetscReal   w;
324afe6d6adSToby Isaac 
325afe6d6adSToby Isaac       if (isAffine) {
3267132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
327afe6d6adSToby Isaac       } else {
3283fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
3299f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
3309f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
3314bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
3324bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
3339f209ee4SMatthew G. Knepley 
3349f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
3359f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
3369f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
337afe6d6adSToby Isaac       }
3384bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
339afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
34063a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
341afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3429566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
343afe6d6adSToby Isaac #endif
344afe6d6adSToby Isaac       }
34563a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
346*79ab67a3SMatthew G. Knepley       if (debug > 3) {
347*79ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    x_q ("));
348*79ab67a3SMatthew G. Knepley         for (PetscInt d = 0; d < dE; ++d) {
349*79ab67a3SMatthew G. Knepley           if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
350*79ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.v[d]));
351*79ab67a3SMatthew G. Knepley         }
352*79ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
353*79ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    n_q ("));
354*79ab67a3SMatthew G. Knepley         for (PetscInt d = 0; d < dE; ++d) {
355*79ab67a3SMatthew G. Knepley           if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
356*79ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.n[d]));
357*79ab67a3SMatthew G. Knepley         }
358*79ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
359*79ab67a3SMatthew G. Knepley         for (PetscInt f = 0; f < Nf; ++f) {
360*79ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "    u_%" PetscInt_FMT " (", f));
361*79ab67a3SMatthew G. Knepley           for (PetscInt c = 0; c < uOff[f + 1] - uOff[f]; ++c) {
362*79ab67a3SMatthew G. Knepley             if (c) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
363*79ab67a3SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(u[uOff[f] + c])));
364*79ab67a3SMatthew G. Knepley           }
365*79ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
366*79ab67a3SMatthew G. Knepley         }
367*79ab67a3SMatthew G. Knepley       }
3689566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3699566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3704bee2e38SMatthew 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);
3714bee2e38SMatthew G. Knepley       integrand *= w;
372afe6d6adSToby Isaac       integral[e * Nf + field] += integrand;
373*79ab67a3SMatthew G. Knepley       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g tot: %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
374afe6d6adSToby Isaac     }
375afe6d6adSToby Isaac     cOffset += totDim;
376afe6d6adSToby Isaac     cOffsetAux += totDimAux;
377afe6d6adSToby Isaac   }
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
379afe6d6adSToby Isaac }
380afe6d6adSToby Isaac 
381d71ae5a4SJacob 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[])
382d71ae5a4SJacob Faibussowitsch {
383b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
3846528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3854bee2e38SMatthew G. Knepley   PetscFE            fe;
3866528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3876528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
3886528b96dSMatthew G. Knepley   PetscPointFunc    *f0_func, *f1_func;
38920cf1dd8SToby Isaac   PetscQuadrature    quad;
390ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3914bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
39220cf1dd8SToby Isaac   const PetscScalar *constants;
39320cf1dd8SToby Isaac   PetscReal         *x;
394ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
395ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
39620cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3976587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
39820cf1dd8SToby Isaac 
39920cf1dd8SToby Isaac   PetscFunctionBegin;
4009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
4019566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
4029566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
4039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
4099566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
4103ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
4119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
4159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4164bee2e38SMatthew G. Knepley   if (dsAux) {
4179566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4199566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4209566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4219566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4229566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
42363a3b9bcSJacob 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);
42420cf1dd8SToby Isaac   }
4259566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
42663a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
42720cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
42863a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
42920cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4304bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
43120cf1dd8SToby Isaac 
4326587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
4339566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
4349566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
43520cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4364bee2e38SMatthew G. Knepley       PetscReal w;
4374bee2e38SMatthew G. Knepley       PetscInt  c, d;
43820cf1dd8SToby Isaac 
4399566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
4404bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
4416587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
44263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4437be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4449566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
44520cf1dd8SToby Isaac #endif
44620cf1dd8SToby Isaac       }
44716cd844bSPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
4489566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
4496528b96dSMatthew 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]);
450ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
4516528b96dSMatthew 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]);
4529371c9d4SSatish Balay       for (c = 0; c < T[field]->Nc; ++c)
4539371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w;
454b8025e53SMatthew G. Knepley       if (debug) {
455e8e188d2SZach Atkins         // LCOV_EXCL_START
456e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q]));
457e8e188d2SZach Atkins         for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c]));
458e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
459b8025e53SMatthew G. Knepley         if (debug > 2) {
46063a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
46163a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
4629566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
463e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field der %" PetscInt_FMT ":", field));
464e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c])));
465e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
46663a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
46763a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
4689566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
469e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  res der %" PetscInt_FMT ":", field));
470e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc; ++c) {
471e8e188d2SZach Atkins             for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dim + d])));
472b8025e53SMatthew G. Knepley           }
473e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
474e8e188d2SZach Atkins         }
475e8e188d2SZach Atkins         // LCOV_EXCL_STOP
476b8025e53SMatthew G. Knepley       }
47720cf1dd8SToby Isaac     }
4789566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
47920cf1dd8SToby Isaac     cOffset += totDim;
48020cf1dd8SToby Isaac     cOffsetAux += totDimAux;
48120cf1dd8SToby Isaac   }
4823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48320cf1dd8SToby Isaac }
48420cf1dd8SToby Isaac 
485d71ae5a4SJacob 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[])
486d71ae5a4SJacob Faibussowitsch {
487b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
48806d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
4894bee2e38SMatthew G. Knepley   PetscFE            fe;
49006d8a0d3SMatthew G. Knepley   PetscInt           n0, n1, i;
49106d8a0d3SMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
49220cf1dd8SToby Isaac   PetscQuadrature    quad;
493ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
4944bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
49520cf1dd8SToby Isaac   const PetscScalar *constants;
49620cf1dd8SToby Isaac   PetscReal         *x;
497ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
498ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
4996587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
50020cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
5016587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
50220cf1dd8SToby Isaac 
50320cf1dd8SToby Isaac   PetscFunctionBegin;
5049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
5059566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
5069566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
5079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
5109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
5119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
5129566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
5133ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
5149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
5159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
5169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
5179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
5189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
5194bee2e38SMatthew G. Knepley   if (dsAux) {
5209566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
5219566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
5229566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5239566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
5249566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
5259566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
5267be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
5279566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
5289566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
52963a3b9bcSJacob 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);
53020cf1dd8SToby Isaac   }
531ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
5329566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
53363a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
53420cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
5356587ee25SMatthew G. Knepley   /* TODO FIX THIS */
5366587ee25SMatthew G. Knepley   fgeom->dim = dim - 1;
53763a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
53820cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5399f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
54020cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5419f209ee4SMatthew G. Knepley 
5426587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
5439566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcI));
5449566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
54520cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5464bee2e38SMatthew G. Knepley       PetscReal w;
5474bee2e38SMatthew G. Knepley       PetscInt  c, d;
5484bee2e38SMatthew G. Knepley 
5499566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
5509566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
5514bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
55262bd480fSMatthew G. Knepley       if (debug > 1) {
5536587ee25SMatthew G. Knepley         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
55463a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5557be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5569566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5579566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
55820cf1dd8SToby Isaac #endif
55920cf1dd8SToby Isaac         }
56062bd480fSMatthew G. Knepley       }
5618e3a54c0SPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
5629566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
56306d8a0d3SMatthew 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]);
5644bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w;
56506d8a0d3SMatthew 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]);
5669371c9d4SSatish Balay       for (c = 0; c < NcI; ++c)
5679371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w;
56862bd480fSMatthew G. Knepley       if (debug) {
56963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
57062bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
57163a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c])));
57262bd480fSMatthew G. Knepley           if (n1) {
57363a3b9bcSJacob 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])));
5749566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
57562bd480fSMatthew G. Knepley           }
57662bd480fSMatthew G. Knepley         }
57762bd480fSMatthew G. Knepley       }
57820cf1dd8SToby Isaac     }
5799566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
58020cf1dd8SToby Isaac     cOffset += totDim;
58120cf1dd8SToby Isaac     cOffsetAux += totDimAux;
58220cf1dd8SToby Isaac   }
5833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58420cf1dd8SToby Isaac }
58520cf1dd8SToby Isaac 
58627f02ce8SMatthew G. Knepley /*
58727f02ce8SMatthew 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
58827f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
58927f02ce8SMatthew G. Knepley 
59027f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
59127f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
59227f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
59327f02ce8SMatthew 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()
59427f02ce8SMatthew G. Knepley */
5952dce792eSToby 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[])
596d71ae5a4SJacob Faibussowitsch {
597b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
5986528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
59927f02ce8SMatthew G. Knepley   PetscFE            fe;
6006528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6016528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
6026528b96dSMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
60327f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
6040e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
60507218a29SMatthew G. Knepley   PetscTabulation   *Tf, *TfIn, *TfAux = NULL;
60627f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
60727f02ce8SMatthew G. Knepley   const PetscScalar *constants;
60827f02ce8SMatthew G. Knepley   PetscReal         *x;
609665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
61007218a29SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
6116587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
61227f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
6136587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
61427f02ce8SMatthew G. Knepley 
61527f02ce8SMatthew G. Knepley   PetscFunctionBegin;
61627f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
6179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
6189566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
6199566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
6209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
62207218a29SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
623429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets
62407218a29SMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x));
6259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
6269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
6279566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
6283ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
6299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
6309566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
6319566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
63227f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
6339566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
63407218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
6359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
63627f02ce8SMatthew G. Knepley   if (dsAux) {
6379566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
6389566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6399566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6409566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6419566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6429566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
64301907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
6449566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
6459566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
64663a3b9bcSJacob 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);
64727f02ce8SMatthew G. Knepley   }
6489566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
649665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
650c2b7495fSMatthew G. Knepley   NcS = NcI;
6510abb75b6SMatthew G. Knepley   if (!isCohesiveField && s == 2) {
6520abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
6530abb75b6SMatthew G. Knepley     NcS *= 2;
6540abb75b6SMatthew G. Knepley   }
6559566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
6560e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
65763a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
65827f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
65963a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
66027f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
66127f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
6620e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
6630e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
6644e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
66527f02ce8SMatthew G. Knepley 
6666587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6679566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcS));
6689566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
66927f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
6700e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
67127f02ce8SMatthew G. Knepley       PetscReal w;
67227f02ce8SMatthew G. Knepley       PetscInt  c, d;
67327f02ce8SMatthew G. Knepley 
6744e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0]));
6754e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1]));
67607218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
67727f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
6786587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
67963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
68027f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
6819566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
68227f02ce8SMatthew G. Knepley #endif
68327f02ce8SMatthew G. Knepley       }
684a4158a15SMatthew 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]));
68527f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
6868e3a54c0SPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], PetscSafePointerPlusOffset(coefficients_t, cOffsetIn), u, u_x, u_t));
68707218a29SMatthew 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));
6886528b96dSMatthew 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]);
68927f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
6909ee2af8cSMatthew 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]);
6919371c9d4SSatish Balay       for (c = 0; c < NcS; ++c)
6929371c9d4SSatish Balay         for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
69327f02ce8SMatthew G. Knepley     }
6949371c9d4SSatish Balay     if (isCohesiveField) {
6953ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
6969371c9d4SSatish Balay     } else {
6973ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
6989371c9d4SSatish Balay     }
69927f02ce8SMatthew G. Knepley     cOffset += totDim;
70007218a29SMatthew G. Knepley     cOffsetIn += totDimIn;
70127f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
70227f02ce8SMatthew G. Knepley   }
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70427f02ce8SMatthew G. Knepley }
70527f02ce8SMatthew G. Knepley 
706d71ae5a4SJacob 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[])
707d71ae5a4SJacob Faibussowitsch {
708b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
7094bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
7106528b96dSMatthew G. Knepley   PetscWeakForm      wf;
7116528b96dSMatthew G. Knepley   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
7126528b96dSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
71320cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
71420cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
71520cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
71620cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
71720cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
71820cf1dd8SToby Isaac   PetscQuadrature    quad;
719ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
7204bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
72120cf1dd8SToby Isaac   const PetscScalar *constants;
72220cf1dd8SToby Isaac   PetscReal         *x;
723ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
724ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
7256528b96dSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
72620cf1dd8SToby Isaac   PetscInt           dE, Np;
72720cf1dd8SToby Isaac   PetscBool          isAffine;
72820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
72920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
73020cf1dd8SToby Isaac 
73120cf1dd8SToby Isaac   PetscFunctionBegin;
7329566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
7336528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
7346528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
7359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
7369566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
7379566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
7389566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
7399566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
7409566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
7419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
7429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
74320cf1dd8SToby Isaac   switch (jtype) {
744d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
745d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
746d71ae5a4SJacob Faibussowitsch     break;
747d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
748d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
749d71ae5a4SJacob Faibussowitsch     break;
750d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
751d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
752d71ae5a4SJacob Faibussowitsch     break;
75320cf1dd8SToby Isaac   }
7543ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
7559566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
7569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
7579566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
7589566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
7599566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
7609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
7619566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
7624bee2e38SMatthew G. Knepley   if (dsAux) {
7639566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7649566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
7659566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
7669566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
7679566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
7689566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
76963a3b9bcSJacob 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);
77020cf1dd8SToby Isaac   }
77127f02ce8SMatthew G. Knepley   NcI      = T[fieldI]->Nc;
77227f02ce8SMatthew G. Knepley   NcJ      = T[fieldJ]->Nc;
7734bee2e38SMatthew G. Knepley   Np       = cgeom->numPoints;
7744bee2e38SMatthew G. Knepley   dE       = cgeom->dimEmbed;
7754bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
77627f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
7779566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
7789566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7799566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7809566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7819566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
78263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
7834bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
7844bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
7854bee2e38SMatthew G. Knepley 
78627f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
78727f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
7884bee2e38SMatthew G. Knepley     if (isAffine) {
7894bee2e38SMatthew G. Knepley       fegeom.v    = x;
7904bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
7917132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
7927132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
7937132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
7944bee2e38SMatthew G. Knepley     }
79520cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
79620cf1dd8SToby Isaac       PetscReal w;
7974bee2e38SMatthew G. Knepley       PetscInt  c;
79820cf1dd8SToby Isaac 
79920cf1dd8SToby Isaac       if (isAffine) {
8007132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
80120cf1dd8SToby Isaac       } else {
8024bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
8034bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
8044bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
8054bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
80620cf1dd8SToby Isaac       }
8079566063dSJacob 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]));
8084bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
80916cd844bSPierre Jolivet       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
8109566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
811ea672e62SMatthew G. Knepley       if (n0) {
8129566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
8136528b96dSMatthew 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);
81420cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
81520cf1dd8SToby Isaac       }
816ea672e62SMatthew G. Knepley       if (n1) {
8179566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
8186528b96dSMatthew 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);
8194bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
82020cf1dd8SToby Isaac       }
821ea672e62SMatthew G. Knepley       if (n2) {
8229566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
8236528b96dSMatthew 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);
8244bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
82520cf1dd8SToby Isaac       }
826ea672e62SMatthew G. Knepley       if (n3) {
8279566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8286528b96dSMatthew 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);
8294bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
83020cf1dd8SToby Isaac       }
83120cf1dd8SToby Isaac 
8329566063dSJacob 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));
83320cf1dd8SToby Isaac     }
83420cf1dd8SToby Isaac     if (debug > 1) {
83520cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
83620cf1dd8SToby Isaac 
83763a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
838ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
839ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
840ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
841ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
842ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
843ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
84463a3b9bcSJacob 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])));
84520cf1dd8SToby Isaac             }
84620cf1dd8SToby Isaac           }
8479566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
84820cf1dd8SToby Isaac         }
84920cf1dd8SToby Isaac       }
85020cf1dd8SToby Isaac     }
85120cf1dd8SToby Isaac     cOffset += totDim;
85220cf1dd8SToby Isaac     cOffsetAux += totDimAux;
85320cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
85420cf1dd8SToby Isaac   }
8553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
85620cf1dd8SToby Isaac }
85720cf1dd8SToby Isaac 
858e3d591f2SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
859d71ae5a4SJacob Faibussowitsch {
860b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
8614bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
86245480ffeSMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
86345480ffeSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
86420cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
86520cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
86620cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
86720cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
86820cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
86920cf1dd8SToby Isaac   PetscQuadrature    quad;
870ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
8714bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
87220cf1dd8SToby Isaac   const PetscScalar *constants;
87320cf1dd8SToby Isaac   PetscReal         *x;
874ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
875ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
87645480ffeSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
87720cf1dd8SToby Isaac   PetscBool          isAffine;
87820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
87920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
88020cf1dd8SToby Isaac 
88120cf1dd8SToby Isaac   PetscFunctionBegin;
8829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
88345480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
88445480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
8859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
8869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
8879566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
8889566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
8899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
8909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
8919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
8929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
8939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
894e3d591f2SMatthew G. Knepley   switch (jtype) {
895e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN_PRE:
896e3d591f2SMatthew G. Knepley     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
897e3d591f2SMatthew G. Knepley     break;
898e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN:
8999566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
900e3d591f2SMatthew G. Knepley     break;
901e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN:
902e3d591f2SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()");
903e3d591f2SMatthew G. Knepley   }
9043ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
9059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
9069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
9079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
9089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
9099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
9104bee2e38SMatthew G. Knepley   if (dsAux) {
9119566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
9129566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
9139566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
9149566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
9159566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
9169566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
91720cf1dd8SToby Isaac   }
918ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
91920cf1dd8SToby Isaac   Np       = fgeom->numPoints;
92020cf1dd8SToby Isaac   dE       = fgeom->dimEmbed;
92120cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
92227f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
9239566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
9249566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
9259566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
9269566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
9279566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
92863a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
92920cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
9309f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
93120cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
932ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
933ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
934ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
935ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
93627f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
93727f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
93827f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
93927f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
9404bee2e38SMatthew G. Knepley     if (isAffine) {
9414bee2e38SMatthew G. Knepley       fegeom.v    = x;
9424bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
9437132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
9447132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
9457132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
9467132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
9479f209ee4SMatthew G. Knepley 
9487132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
9497132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
9507132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
9514bee2e38SMatthew G. Knepley     }
95220cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
95320cf1dd8SToby Isaac       PetscReal w;
9544bee2e38SMatthew G. Knepley       PetscInt  c;
95520cf1dd8SToby Isaac 
95663a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
95720cf1dd8SToby Isaac       if (isAffine) {
9587132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
95920cf1dd8SToby Isaac       } else {
9603fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
9619f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
9629f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
9634bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
9644bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
9659f209ee4SMatthew G. Knepley 
9669f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
9679f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
9689f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
96920cf1dd8SToby Isaac       }
9704bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
9719566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
9729566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
973ea672e62SMatthew G. Knepley       if (n0) {
9749566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
97545480ffeSMatthew 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);
97620cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
97720cf1dd8SToby Isaac       }
978ea672e62SMatthew G. Knepley       if (n1) {
9799566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
98045480ffeSMatthew 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);
9814bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
98220cf1dd8SToby Isaac       }
983ea672e62SMatthew G. Knepley       if (n2) {
9849566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
98545480ffeSMatthew 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);
9864bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
98720cf1dd8SToby Isaac       }
988ea672e62SMatthew G. Knepley       if (n3) {
9899566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
99045480ffeSMatthew 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);
9914bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
99220cf1dd8SToby Isaac       }
99320cf1dd8SToby Isaac 
9949566063dSJacob 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));
99520cf1dd8SToby Isaac     }
99620cf1dd8SToby Isaac     if (debug > 1) {
99720cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
99820cf1dd8SToby Isaac 
99963a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
1000ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
1001ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
1002ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
1003ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
1004ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
1005ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
100663a3b9bcSJacob 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])));
100720cf1dd8SToby Isaac             }
100820cf1dd8SToby Isaac           }
10099566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
101020cf1dd8SToby Isaac         }
101120cf1dd8SToby Isaac       }
101220cf1dd8SToby Isaac     }
101320cf1dd8SToby Isaac     cOffset += totDim;
101420cf1dd8SToby Isaac     cOffsetAux += totDimAux;
101520cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
101620cf1dd8SToby Isaac   }
10173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
101820cf1dd8SToby Isaac }
101920cf1dd8SToby Isaac 
10202dce792eSToby 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[])
1021d71ae5a4SJacob Faibussowitsch {
1022b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
102327f02ce8SMatthew G. Knepley   PetscFE            feI, feJ;
1024148442b3SMatthew G. Knepley   PetscWeakForm      wf;
1025148442b3SMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
1026148442b3SMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
102727f02ce8SMatthew G. Knepley   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
102827f02ce8SMatthew G. Knepley   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
102927f02ce8SMatthew G. Knepley   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
103027f02ce8SMatthew G. Knepley   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
103127f02ce8SMatthew G. Knepley   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
1032665f567fSMatthew G. Knepley   PetscQuadrature    quad;
10330e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
103407218a29SMatthew G. Knepley   PetscTabulation   *T, *TfIn, *TAux = NULL;
103527f02ce8SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
103627f02ce8SMatthew G. Knepley   const PetscScalar *constants;
103727f02ce8SMatthew G. Knepley   PetscReal         *x;
1038665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
1039665f567fSMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
104045480ffeSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
104107218a29SMatthew G. Knepley   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE;
104227f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
10430502970dSMatthew G. Knepley   PetscInt           qNc, Nq, q;
104427f02ce8SMatthew G. Knepley 
104527f02ce8SMatthew G. Knepley   PetscFunctionBegin;
10469566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
104745480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
104845480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
104927f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
10509566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
10519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
10529566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
10539566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
10549566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1055429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets
10569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
10579566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
105827f02ce8SMatthew G. Knepley   switch (jtype) {
1059d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
1060d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1061d71ae5a4SJacob Faibussowitsch     break;
1062d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
1063d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1064d71ae5a4SJacob Faibussowitsch     break;
1065d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
1066d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
106727f02ce8SMatthew G. Knepley   }
10683ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
10699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
10709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
10719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
10729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
107307218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
10749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
10759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
10769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
107727f02ce8SMatthew G. Knepley   if (dsAux) {
10789566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
10799566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
10809566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
10819566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
10829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
10839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
108401907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
10859566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
10869566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
108763a3b9bcSJacob 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);
108827f02ce8SMatthew G. Knepley   }
10899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
10909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1091665f567fSMatthew G. Knepley   NcI = T[fieldI]->Nc;
1092665f567fSMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
109327f02ce8SMatthew G. Knepley   NcS = isCohesiveFieldI ? NcI : 2 * NcI;
109427f02ce8SMatthew G. Knepley   NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ;
10950abb75b6SMatthew G. Knepley   if (!isCohesiveFieldI && s == 2) {
10960abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
10970abb75b6SMatthew G. Knepley     NcS *= 2;
10980abb75b6SMatthew G. Knepley   }
10990abb75b6SMatthew G. Knepley   if (!isCohesiveFieldJ && s == 2) {
11000abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
11010abb75b6SMatthew G. Knepley     NcT *= 2;
11020abb75b6SMatthew G. Knepley   }
11030502970dSMatthew G. Knepley   // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though
11040502970dSMatthew G. Knepley   // the coordinates are in dE dimensions
11059566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS * NcT));
11060502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
11070502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
11080502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
11099566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
11100e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
111163a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
111227f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
111327f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
11140e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
11150e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
11164e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
111727f02ce8SMatthew G. Knepley 
111807218a29SMatthew G. Knepley     fegeom.v = x; /* Workspace */
111927f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
11200e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
112127f02ce8SMatthew G. Knepley       PetscReal w;
112227f02ce8SMatthew G. Knepley       PetscInt  c;
112327f02ce8SMatthew G. Knepley 
11244e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0]));
11254e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1]));
112607218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
112727f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
112807218a29SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
112963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
113027f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
11319566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
113227f02ce8SMatthew G. Knepley #endif
113327f02ce8SMatthew G. Knepley       }
113463a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
11358e3a54c0SPierre Jolivet       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
113607218a29SMatthew 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));
1137ea672e62SMatthew G. Knepley       if (n0) {
11389566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS * NcT));
1139148442b3SMatthew 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);
114027f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT; ++c) g0[c] *= w;
114127f02ce8SMatthew G. Knepley       }
1142ea672e62SMatthew G. Knepley       if (n1) {
11430502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
1144148442b3SMatthew 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);
11450502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w;
114627f02ce8SMatthew G. Knepley       }
1147ea672e62SMatthew G. Knepley       if (n2) {
11480502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
1149148442b3SMatthew 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);
11500502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w;
115127f02ce8SMatthew G. Knepley       }
1152ea672e62SMatthew G. Knepley       if (n3) {
11530502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
1154148442b3SMatthew 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);
11550502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w;
115627f02ce8SMatthew G. Knepley       }
115727f02ce8SMatthew G. Knepley 
11585fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
11595fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
11609566063dSJacob 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));
116127f02ce8SMatthew G. Knepley         } else {
11620abb75b6SMatthew 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));
11630abb75b6SMatthew 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));
11640abb75b6SMatthew G. Knepley         }
11650abb75b6SMatthew G. Knepley       } else {
11660abb75b6SMatthew G. Knepley         if (s == 2) {
11670abb75b6SMatthew G. Knepley           if (isCohesiveFieldJ) {
11680abb75b6SMatthew 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));
11690abb75b6SMatthew 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));
11700abb75b6SMatthew G. Knepley           } else {
11710abb75b6SMatthew 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));
11720abb75b6SMatthew 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));
11730abb75b6SMatthew 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));
11740abb75b6SMatthew 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));
11755fedec97SMatthew G. Knepley           }
11769371c9d4SSatish Balay         } else
11770abb75b6SMatthew 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));
11780abb75b6SMatthew G. Knepley       }
117927f02ce8SMatthew G. Knepley     }
118027f02ce8SMatthew G. Knepley     if (debug > 1) {
11814e913f38SMatthew G. Knepley       const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb));
11824e913f38SMatthew G. Knepley       const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb));
11834e913f38SMatthew G. Knepley       const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb));
11844e913f38SMatthew G. Knepley       const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb));
11854e913f38SMatthew G. Knepley       PetscInt       f, g;
118627f02ce8SMatthew G. Knepley 
11874e913f38SMatthew 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));
11884e913f38SMatthew G. Knepley       for (f = fS; f < fE; ++f) {
11894e913f38SMatthew G. Knepley         const PetscInt i = offsetI + f;
11904e913f38SMatthew G. Knepley         for (g = gS; g < gE; ++g) {
11914e913f38SMatthew G. Knepley           const PetscInt j = offsetJ + g;
1192e978a55eSPierre Jolivet           PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, f, i, g, j);
11934e913f38SMatthew 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])));
119427f02ce8SMatthew G. Knepley         }
11959566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
119627f02ce8SMatthew G. Knepley       }
119727f02ce8SMatthew G. Knepley     }
119827f02ce8SMatthew G. Knepley     cOffset += totDim;
119927f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
120027f02ce8SMatthew G. Knepley     eOffset += PetscSqr(totDim);
120127f02ce8SMatthew G. Knepley   }
12023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120327f02ce8SMatthew G. Knepley }
120427f02ce8SMatthew G. Knepley 
1205d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1206d71ae5a4SJacob Faibussowitsch {
120720cf1dd8SToby Isaac   PetscFunctionBegin;
120820cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
120920cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
121020cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
121120cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
121220cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1213ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
121420cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1215afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
121620cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
121720cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
121827f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
121920cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
122020cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
122120cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
122227f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
12233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122420cf1dd8SToby Isaac }
122520cf1dd8SToby Isaac 
122620cf1dd8SToby Isaac /*MC
1227dce8aebaSBarry Smith   PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization
122820cf1dd8SToby Isaac 
122920cf1dd8SToby Isaac   Level: intermediate
123020cf1dd8SToby Isaac 
1231dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
123220cf1dd8SToby Isaac M*/
123320cf1dd8SToby Isaac 
1234d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1235d71ae5a4SJacob Faibussowitsch {
123620cf1dd8SToby Isaac   PetscFE_Basic *b;
123720cf1dd8SToby Isaac 
123820cf1dd8SToby Isaac   PetscFunctionBegin;
123920cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
12404dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
124120cf1dd8SToby Isaac   fem->data = b;
124220cf1dd8SToby Isaac 
12439566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
12443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124520cf1dd8SToby Isaac }
1246