xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 2dce792e531186164765a9583d36d03ffc15e9ea)
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 
162*2dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
163d71ae5a4SJacob Faibussowitsch {
16420cf1dd8SToby Isaac   const PetscInt     debug = 0;
1654bee2e38SMatthew G. Knepley   PetscFE            fe;
16620cf1dd8SToby Isaac   PetscPointFunc     obj_func;
16720cf1dd8SToby Isaac   PetscQuadrature    quad;
168ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
1694bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
17020cf1dd8SToby Isaac   const PetscScalar *constants;
17120cf1dd8SToby Isaac   PetscReal         *x;
172ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
17320cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
17420cf1dd8SToby Isaac   PetscBool          isAffine;
17520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
17620cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
17720cf1dd8SToby Isaac 
17820cf1dd8SToby Isaac   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetObjective(ds, field, &obj_func));
1803ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
1819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
1849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
1879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
1889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
1899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
1909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1924bee2e38SMatthew G. Knepley   if (dsAux) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
1989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
19963a3b9bcSJacob Faibussowitsch     PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
20020cf1dd8SToby Isaac   }
2019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
20263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
20320cf1dd8SToby Isaac   Np       = cgeom->numPoints;
20420cf1dd8SToby Isaac   dE       = cgeom->dimEmbed;
20520cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
20620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2074bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
20820cf1dd8SToby Isaac 
20927f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
21027f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
21120cf1dd8SToby Isaac     if (isAffine) {
2124bee2e38SMatthew G. Knepley       fegeom.v    = x;
2134bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2147132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
2157132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
2167132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
21720cf1dd8SToby Isaac     }
2184bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
2194bee2e38SMatthew G. Knepley       PetscScalar integrand;
2204bee2e38SMatthew G. Knepley       PetscReal   w;
2214bee2e38SMatthew G. Knepley 
2224bee2e38SMatthew G. Knepley       if (isAffine) {
2237132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
2244bee2e38SMatthew G. Knepley       } else {
2254bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
2264bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
2274bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
2284bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
2294bee2e38SMatthew G. Knepley       }
2304bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
23120cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
23263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
2337be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2349566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
23520cf1dd8SToby Isaac #endif
23620cf1dd8SToby Isaac       }
23763a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
2389566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
2399566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
2404bee2e38SMatthew G. Knepley       obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, numConstants, constants, &integrand);
2414bee2e38SMatthew G. Knepley       integrand *= w;
24220cf1dd8SToby Isaac       integral[e * Nf + field] += integrand;
2439566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[field])));
24420cf1dd8SToby Isaac     }
24520cf1dd8SToby Isaac     cOffset += totDim;
24620cf1dd8SToby Isaac     cOffsetAux += totDimAux;
24720cf1dd8SToby Isaac   }
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24920cf1dd8SToby Isaac }
25020cf1dd8SToby Isaac 
251*2dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFunc obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
252d71ae5a4SJacob Faibussowitsch {
253afe6d6adSToby Isaac   const PetscInt     debug = 0;
2544bee2e38SMatthew G. Knepley   PetscFE            fe;
255afe6d6adSToby Isaac   PetscQuadrature    quad;
256ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2574bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
258afe6d6adSToby Isaac   const PetscScalar *constants;
259afe6d6adSToby Isaac   PetscReal         *x;
260ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
261afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
262afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
263afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
264afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
265afe6d6adSToby Isaac 
266afe6d6adSToby Isaac   PetscFunctionBegin;
2673ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
2689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
2719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
2749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
2759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
2769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
2779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
2789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
2794bee2e38SMatthew G. Knepley   if (dsAux) {
2809566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
2819566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
2839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
286afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
2879566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
2889566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
28963a3b9bcSJacob Faibussowitsch     PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
290afe6d6adSToby Isaac   }
2919566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
29263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
293afe6d6adSToby Isaac   Np       = fgeom->numPoints;
294afe6d6adSToby Isaac   dE       = fgeom->dimEmbed;
295afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
296afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
2979f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
298afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
299ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
300ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
301ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
302ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
30327f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
30427f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
30527f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
30627f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
3074bee2e38SMatthew G. Knepley     if (isAffine) {
3084bee2e38SMatthew G. Knepley       fegeom.v    = x;
3094bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3107132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
3117132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
3127132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
3137132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
3149f209ee4SMatthew G. Knepley 
3157132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
3167132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
3177132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
3184bee2e38SMatthew G. Knepley     }
319afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
320afe6d6adSToby Isaac       PetscScalar integrand;
3214bee2e38SMatthew G. Knepley       PetscReal   w;
322afe6d6adSToby Isaac 
323afe6d6adSToby Isaac       if (isAffine) {
3247132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
325afe6d6adSToby Isaac       } else {
3263fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
3279f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
3289f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
3294bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
3304bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
3319f209ee4SMatthew G. Knepley 
3329f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
3339f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
3349f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
335afe6d6adSToby Isaac       }
3364bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
337afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
33863a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
339afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3409566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
341afe6d6adSToby Isaac #endif
342afe6d6adSToby Isaac       }
34363a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
3449566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3459566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3464bee2e38SMatthew G. Knepley       obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, fegeom.n, numConstants, constants, &integrand);
3474bee2e38SMatthew G. Knepley       integrand *= w;
348afe6d6adSToby Isaac       integral[e * Nf + field] += integrand;
3499566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
350afe6d6adSToby Isaac     }
351afe6d6adSToby Isaac     cOffset += totDim;
352afe6d6adSToby Isaac     cOffsetAux += totDimAux;
353afe6d6adSToby Isaac   }
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355afe6d6adSToby Isaac }
356afe6d6adSToby Isaac 
357d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
358d71ae5a4SJacob Faibussowitsch {
35920cf1dd8SToby Isaac   const PetscInt     debug = 0;
3606528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3614bee2e38SMatthew G. Knepley   PetscFE            fe;
3626528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3636528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
3646528b96dSMatthew G. Knepley   PetscPointFunc    *f0_func, *f1_func;
36520cf1dd8SToby Isaac   PetscQuadrature    quad;
366ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3674bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
36820cf1dd8SToby Isaac   const PetscScalar *constants;
36920cf1dd8SToby Isaac   PetscReal         *x;
370ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
371ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
37220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3736587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac   PetscFunctionBegin;
3769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
3779566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
3789566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
3799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
3819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
3849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
3859566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
3863ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
3879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
3889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
3899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
3909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
3919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3924bee2e38SMatthew G. Knepley   if (dsAux) {
3939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
3959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
3989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
39963a3b9bcSJacob Faibussowitsch     PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
40020cf1dd8SToby Isaac   }
4019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
40263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
40320cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
40463a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
40520cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4064bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
40720cf1dd8SToby Isaac 
4086587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
4099566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
4109566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
41120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4124bee2e38SMatthew G. Knepley       PetscReal w;
4134bee2e38SMatthew G. Knepley       PetscInt  c, d;
41420cf1dd8SToby Isaac 
4159566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
4164bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
4176587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
41863a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4197be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4209566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
42120cf1dd8SToby Isaac #endif
42220cf1dd8SToby Isaac       }
4239566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
4249566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
4256528b96dSMatthew G. Knepley       for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q * T[field]->Nc]);
426ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
4276528b96dSMatthew G. Knepley       for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q * T[field]->Nc * dim]);
4289371c9d4SSatish Balay       for (c = 0; c < T[field]->Nc; ++c)
4299371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w;
430b8025e53SMatthew G. Knepley       if (debug) {
43163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g\n", q, (double)quadWeights[q]));
432b8025e53SMatthew G. Knepley         if (debug > 2) {
43363a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
43463a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
4359566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
43663a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
43763a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
4389566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
439b8025e53SMatthew G. Knepley         }
440b8025e53SMatthew G. Knepley       }
44120cf1dd8SToby Isaac     }
4429566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
44320cf1dd8SToby Isaac     cOffset += totDim;
44420cf1dd8SToby Isaac     cOffsetAux += totDimAux;
44520cf1dd8SToby Isaac   }
4463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44720cf1dd8SToby Isaac }
44820cf1dd8SToby Isaac 
449d71ae5a4SJacob 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[])
450d71ae5a4SJacob Faibussowitsch {
45120cf1dd8SToby Isaac   const PetscInt     debug = 0;
45206d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
4534bee2e38SMatthew G. Knepley   PetscFE            fe;
45406d8a0d3SMatthew G. Knepley   PetscInt           n0, n1, i;
45506d8a0d3SMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
45620cf1dd8SToby Isaac   PetscQuadrature    quad;
457ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
4584bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
45920cf1dd8SToby Isaac   const PetscScalar *constants;
46020cf1dd8SToby Isaac   PetscReal         *x;
461ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
462ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
4636587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
46420cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
4656587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
46620cf1dd8SToby Isaac 
46720cf1dd8SToby Isaac   PetscFunctionBegin;
4689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
4699566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
4709566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
4719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4769566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
4773ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
4789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
4829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4834bee2e38SMatthew G. Knepley   if (dsAux) {
4849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
4859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4907be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
4919566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
4929566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
49363a3b9bcSJacob 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);
49420cf1dd8SToby Isaac   }
495ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
4969566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
49763a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
49820cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
4996587ee25SMatthew G. Knepley   /* TODO FIX THIS */
5006587ee25SMatthew G. Knepley   fgeom->dim = dim - 1;
50163a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
50220cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5039f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
50420cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5059f209ee4SMatthew G. Knepley 
5066587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
5079566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcI));
5089566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
50920cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5104bee2e38SMatthew G. Knepley       PetscReal w;
5114bee2e38SMatthew G. Knepley       PetscInt  c, d;
5124bee2e38SMatthew G. Knepley 
5139566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
5149566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
5154bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
51662bd480fSMatthew G. Knepley       if (debug > 1) {
5176587ee25SMatthew G. Knepley         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
51863a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5197be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5209566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5219566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
52220cf1dd8SToby Isaac #endif
52320cf1dd8SToby Isaac         }
52462bd480fSMatthew G. Knepley       }
5259566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
5269566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
52706d8a0d3SMatthew 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]);
5284bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w;
52906d8a0d3SMatthew 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]);
5309371c9d4SSatish Balay       for (c = 0; c < NcI; ++c)
5319371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w;
53262bd480fSMatthew G. Knepley       if (debug) {
53363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
53462bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
53563a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c])));
53662bd480fSMatthew G. Knepley           if (n1) {
53763a3b9bcSJacob 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])));
5389566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
53962bd480fSMatthew G. Knepley           }
54062bd480fSMatthew G. Knepley         }
54162bd480fSMatthew G. Knepley       }
54220cf1dd8SToby Isaac     }
5439566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
54420cf1dd8SToby Isaac     cOffset += totDim;
54520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
54620cf1dd8SToby Isaac   }
5473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54820cf1dd8SToby Isaac }
54920cf1dd8SToby Isaac 
55027f02ce8SMatthew G. Knepley /*
55127f02ce8SMatthew 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
55227f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
55327f02ce8SMatthew G. Knepley 
55427f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
55527f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
55627f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
55727f02ce8SMatthew 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()
55827f02ce8SMatthew G. Knepley */
559*2dce792eSToby 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[])
560d71ae5a4SJacob Faibussowitsch {
56127f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
5626528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
56327f02ce8SMatthew G. Knepley   PetscFE            fe;
5646528b96dSMatthew G. Knepley   PetscWeakForm      wf;
5656528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
5666528b96dSMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
56727f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
5680e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
56907218a29SMatthew G. Knepley   PetscTabulation   *Tf, *TfIn, *TfAux = NULL;
57027f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
57127f02ce8SMatthew G. Knepley   const PetscScalar *constants;
57227f02ce8SMatthew G. Knepley   PetscReal         *x;
573665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
57407218a29SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
5756587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
57627f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
5776587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
57827f02ce8SMatthew G. Knepley 
57927f02ce8SMatthew G. Knepley   PetscFunctionBegin;
58027f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
5819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
5829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
5839566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
5849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
58607218a29SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
587429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets
58807218a29SMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x));
5899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
5909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
5919566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
5923ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
5939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
5949566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
5959566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
59627f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
5979566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
59807218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
5999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
60027f02ce8SMatthew G. Knepley   if (dsAux) {
6019566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
6029566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6039566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6049566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6059566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6069566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
60701907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
6089566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
6099566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
61063a3b9bcSJacob 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);
61127f02ce8SMatthew G. Knepley   }
6129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
613665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
614c2b7495fSMatthew G. Knepley   NcS = NcI;
6159566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
6160e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
61763a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
61827f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
61963a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
62027f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
62127f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
6220e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
6230e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
6244e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
62527f02ce8SMatthew G. Knepley 
6266587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6279566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcS));
6289566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
62927f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
6300e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
63127f02ce8SMatthew G. Knepley       PetscReal w;
63227f02ce8SMatthew G. Knepley       PetscInt  c, d;
63327f02ce8SMatthew G. Knepley 
6344e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0]));
6354e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1]));
63607218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
63727f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
6386587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
63963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
64027f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
6419566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
64227f02ce8SMatthew G. Knepley #endif
64327f02ce8SMatthew G. Knepley       }
644a4158a15SMatthew 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]));
64527f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
64607218a29SMatthew G. Knepley       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], &coefficients_t[cOffsetIn], u, u_x, u_t));
64707218a29SMatthew 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));
6486528b96dSMatthew 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]);
64927f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
6509ee2af8cSMatthew 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]);
6519371c9d4SSatish Balay       for (c = 0; c < NcS; ++c)
6529371c9d4SSatish Balay         for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
65327f02ce8SMatthew G. Knepley     }
6549371c9d4SSatish Balay     if (isCohesiveField) {
6553ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
6569371c9d4SSatish Balay     } else {
6573ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
6589371c9d4SSatish Balay     }
65927f02ce8SMatthew G. Knepley     cOffset += totDim;
66007218a29SMatthew G. Knepley     cOffsetIn += totDimIn;
66127f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
66227f02ce8SMatthew G. Knepley   }
6633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66427f02ce8SMatthew G. Knepley }
66527f02ce8SMatthew G. Knepley 
666d71ae5a4SJacob 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[])
667d71ae5a4SJacob Faibussowitsch {
66820cf1dd8SToby Isaac   const PetscInt     debug = 0;
6694bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
6706528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6716528b96dSMatthew G. Knepley   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
6726528b96dSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
67320cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
67420cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
67520cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
67620cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
67720cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
67820cf1dd8SToby Isaac   PetscQuadrature    quad;
679ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
6804bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
68120cf1dd8SToby Isaac   const PetscScalar *constants;
68220cf1dd8SToby Isaac   PetscReal         *x;
683ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
684ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
6856528b96dSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
68620cf1dd8SToby Isaac   PetscInt           dE, Np;
68720cf1dd8SToby Isaac   PetscBool          isAffine;
68820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
68920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
69020cf1dd8SToby Isaac 
69120cf1dd8SToby Isaac   PetscFunctionBegin;
6929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6936528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
6946528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
6959566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
6969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
6979566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
6989566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
6999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
7009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
7019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
7029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
70320cf1dd8SToby Isaac   switch (jtype) {
704d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
705d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
706d71ae5a4SJacob Faibussowitsch     break;
707d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
708d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
709d71ae5a4SJacob Faibussowitsch     break;
710d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
711d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
712d71ae5a4SJacob Faibussowitsch     break;
71320cf1dd8SToby Isaac   }
7143ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
7159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
7169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
7179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
7189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
7199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
7209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
7219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
7224bee2e38SMatthew G. Knepley   if (dsAux) {
7239566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7249566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
7259566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
7269566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
7279566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
7289566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
72963a3b9bcSJacob 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);
73020cf1dd8SToby Isaac   }
73127f02ce8SMatthew G. Knepley   NcI      = T[fieldI]->Nc;
73227f02ce8SMatthew G. Knepley   NcJ      = T[fieldJ]->Nc;
7334bee2e38SMatthew G. Knepley   Np       = cgeom->numPoints;
7344bee2e38SMatthew G. Knepley   dE       = cgeom->dimEmbed;
7354bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
73627f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
7379566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
7389566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7399566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7409566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7419566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
74263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
7434bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
7444bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
7454bee2e38SMatthew G. Knepley 
74627f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
74727f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
7484bee2e38SMatthew G. Knepley     if (isAffine) {
7494bee2e38SMatthew G. Knepley       fegeom.v    = x;
7504bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
7517132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
7527132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
7537132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
7544bee2e38SMatthew G. Knepley     }
75520cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
75620cf1dd8SToby Isaac       PetscReal w;
7574bee2e38SMatthew G. Knepley       PetscInt  c;
75820cf1dd8SToby Isaac 
75920cf1dd8SToby Isaac       if (isAffine) {
7607132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
76120cf1dd8SToby Isaac       } else {
7624bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
7634bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
7644bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
7654bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
76620cf1dd8SToby Isaac       }
7679566063dSJacob 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]));
7684bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
7699566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
7709566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
771ea672e62SMatthew G. Knepley       if (n0) {
7729566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
7736528b96dSMatthew 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);
77420cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
77520cf1dd8SToby Isaac       }
776ea672e62SMatthew G. Knepley       if (n1) {
7779566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7786528b96dSMatthew 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);
7794bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
78020cf1dd8SToby Isaac       }
781ea672e62SMatthew G. Knepley       if (n2) {
7829566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7836528b96dSMatthew 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);
7844bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
78520cf1dd8SToby Isaac       }
786ea672e62SMatthew G. Knepley       if (n3) {
7879566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7886528b96dSMatthew 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);
7894bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
79020cf1dd8SToby Isaac       }
79120cf1dd8SToby Isaac 
7929566063dSJacob 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));
79320cf1dd8SToby Isaac     }
79420cf1dd8SToby Isaac     if (debug > 1) {
79520cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
79620cf1dd8SToby Isaac 
79763a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
798ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
799ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
800ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
801ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
802ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
803ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
80463a3b9bcSJacob 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])));
80520cf1dd8SToby Isaac             }
80620cf1dd8SToby Isaac           }
8079566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
80820cf1dd8SToby Isaac         }
80920cf1dd8SToby Isaac       }
81020cf1dd8SToby Isaac     }
81120cf1dd8SToby Isaac     cOffset += totDim;
81220cf1dd8SToby Isaac     cOffsetAux += totDimAux;
81320cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
81420cf1dd8SToby Isaac   }
8153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81620cf1dd8SToby Isaac }
81720cf1dd8SToby Isaac 
818*2dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
819d71ae5a4SJacob Faibussowitsch {
82020cf1dd8SToby Isaac   const PetscInt     debug = 0;
8214bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
82245480ffeSMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
82345480ffeSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
82420cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
82520cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
82620cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
82720cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
82820cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
82920cf1dd8SToby Isaac   PetscQuadrature    quad;
830ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
8314bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
83220cf1dd8SToby Isaac   const PetscScalar *constants;
83320cf1dd8SToby Isaac   PetscReal         *x;
834ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
835ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
83645480ffeSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
83720cf1dd8SToby Isaac   PetscBool          isAffine;
83820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
83920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
84020cf1dd8SToby Isaac 
84120cf1dd8SToby Isaac   PetscFunctionBegin;
8429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
84345480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
84445480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
8459566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
8469566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
8479566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
8489566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
8499566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
8509566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
8519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
8529566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
8539566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
8549566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
8553ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
8569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
8579566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
8589566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
8599566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
8609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
8614bee2e38SMatthew G. Knepley   if (dsAux) {
8629566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8639566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
8649566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
8659566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
8669566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
8679566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
86820cf1dd8SToby Isaac   }
869ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
87020cf1dd8SToby Isaac   Np       = fgeom->numPoints;
87120cf1dd8SToby Isaac   dE       = fgeom->dimEmbed;
87220cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
87327f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
8749566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
8759566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
8769566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
8779566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8789566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
87963a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
88020cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
8819f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
88220cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
883ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
884ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
885ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
886ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
88727f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
88827f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
88927f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
89027f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
8914bee2e38SMatthew G. Knepley     if (isAffine) {
8924bee2e38SMatthew G. Knepley       fegeom.v    = x;
8934bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
8947132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
8957132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
8967132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
8977132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
8989f209ee4SMatthew G. Knepley 
8997132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
9007132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
9017132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
9024bee2e38SMatthew G. Knepley     }
90320cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
90420cf1dd8SToby Isaac       PetscReal w;
9054bee2e38SMatthew G. Knepley       PetscInt  c;
90620cf1dd8SToby Isaac 
90763a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
90820cf1dd8SToby Isaac       if (isAffine) {
9097132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
91020cf1dd8SToby Isaac       } else {
9113fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
9129f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
9139f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
9144bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
9154bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
9169f209ee4SMatthew G. Knepley 
9179f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
9189f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
9199f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
92020cf1dd8SToby Isaac       }
9214bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
9229566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
9239566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
924ea672e62SMatthew G. Knepley       if (n0) {
9259566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
92645480ffeSMatthew 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);
92720cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
92820cf1dd8SToby Isaac       }
929ea672e62SMatthew G. Knepley       if (n1) {
9309566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
93145480ffeSMatthew 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);
9324bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
93320cf1dd8SToby Isaac       }
934ea672e62SMatthew G. Knepley       if (n2) {
9359566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
93645480ffeSMatthew 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);
9374bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
93820cf1dd8SToby Isaac       }
939ea672e62SMatthew G. Knepley       if (n3) {
9409566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
94145480ffeSMatthew 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);
9424bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
94320cf1dd8SToby Isaac       }
94420cf1dd8SToby Isaac 
9459566063dSJacob 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));
94620cf1dd8SToby Isaac     }
94720cf1dd8SToby Isaac     if (debug > 1) {
94820cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
94920cf1dd8SToby Isaac 
95063a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
951ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
952ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
953ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
954ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
955ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
956ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
95763a3b9bcSJacob 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])));
95820cf1dd8SToby Isaac             }
95920cf1dd8SToby Isaac           }
9609566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
96120cf1dd8SToby Isaac         }
96220cf1dd8SToby Isaac       }
96320cf1dd8SToby Isaac     }
96420cf1dd8SToby Isaac     cOffset += totDim;
96520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
96620cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
96720cf1dd8SToby Isaac   }
9683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
96920cf1dd8SToby Isaac }
97020cf1dd8SToby Isaac 
971*2dce792eSToby 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[])
972d71ae5a4SJacob Faibussowitsch {
97327f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
97427f02ce8SMatthew G. Knepley   PetscFE            feI, feJ;
975148442b3SMatthew G. Knepley   PetscWeakForm      wf;
976148442b3SMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
977148442b3SMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
97827f02ce8SMatthew G. Knepley   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
97927f02ce8SMatthew G. Knepley   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
98027f02ce8SMatthew G. Knepley   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
98127f02ce8SMatthew G. Knepley   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
98227f02ce8SMatthew G. Knepley   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
983665f567fSMatthew G. Knepley   PetscQuadrature    quad;
9840e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
98507218a29SMatthew G. Knepley   PetscTabulation   *T, *TfIn, *TAux = NULL;
98627f02ce8SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
98727f02ce8SMatthew G. Knepley   const PetscScalar *constants;
98827f02ce8SMatthew G. Knepley   PetscReal         *x;
989665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
990665f567fSMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
99145480ffeSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
99207218a29SMatthew G. Knepley   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE;
99327f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
9940502970dSMatthew G. Knepley   PetscInt           qNc, Nq, q;
99527f02ce8SMatthew G. Knepley 
99627f02ce8SMatthew G. Knepley   PetscFunctionBegin;
9979566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
99845480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
99945480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
100027f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
10019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
10029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
10039566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
10049566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
10059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1006429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets
10079566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
10089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
100927f02ce8SMatthew G. Knepley   switch (jtype) {
1010d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
1011d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1012d71ae5a4SJacob Faibussowitsch     break;
1013d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
1014d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1015d71ae5a4SJacob Faibussowitsch     break;
1016d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
1017d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
101827f02ce8SMatthew G. Knepley   }
10193ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
10209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
10219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
10229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
10239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
102407218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
10259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
10269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
10279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
102827f02ce8SMatthew G. Knepley   if (dsAux) {
10299566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
10309566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
10319566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
10329566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
10339566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
10349566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
103501907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
10369566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
10379566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
103863a3b9bcSJacob 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);
103927f02ce8SMatthew G. Knepley   }
10409566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
10419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1042665f567fSMatthew G. Knepley   NcI = T[fieldI]->Nc;
1043665f567fSMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
104427f02ce8SMatthew G. Knepley   NcS = isCohesiveFieldI ? NcI : 2 * NcI;
104527f02ce8SMatthew G. Knepley   NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ;
10460502970dSMatthew G. Knepley   // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though
10470502970dSMatthew G. Knepley   // the coordinates are in dE dimensions
10489566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS * NcT));
10490502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
10500502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
10510502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
10529566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
10530e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
105463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
105527f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
105627f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
10570e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
10580e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
10594e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
106027f02ce8SMatthew G. Knepley 
106107218a29SMatthew G. Knepley     fegeom.v = x; /* Workspace */
106227f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
10630e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
106427f02ce8SMatthew G. Knepley       PetscReal w;
106527f02ce8SMatthew G. Knepley       PetscInt  c;
106627f02ce8SMatthew G. Knepley 
10674e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0]));
10684e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1]));
106907218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
107027f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
107107218a29SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
107263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
107327f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
10749566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
107527f02ce8SMatthew G. Knepley #endif
107627f02ce8SMatthew G. Knepley       }
107763a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
107807218a29SMatthew G. Knepley       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
107907218a29SMatthew 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));
1080ea672e62SMatthew G. Knepley       if (n0) {
10819566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS * NcT));
1082148442b3SMatthew 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);
108327f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT; ++c) g0[c] *= w;
108427f02ce8SMatthew G. Knepley       }
1085ea672e62SMatthew G. Knepley       if (n1) {
10860502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
1087148442b3SMatthew 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);
10880502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w;
108927f02ce8SMatthew G. Knepley       }
1090ea672e62SMatthew G. Knepley       if (n2) {
10910502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
1092148442b3SMatthew 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);
10930502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w;
109427f02ce8SMatthew G. Knepley       }
1095ea672e62SMatthew G. Knepley       if (n3) {
10960502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
1097148442b3SMatthew 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);
10980502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w;
109927f02ce8SMatthew G. Knepley       }
110027f02ce8SMatthew G. Knepley 
11015fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
11025fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
11039566063dSJacob 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));
110427f02ce8SMatthew G. Knepley         } else {
11059566063dSJacob Faibussowitsch           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
11069566063dSJacob Faibussowitsch           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
11075fedec97SMatthew G. Knepley         }
11089371c9d4SSatish Balay       } else
11099371c9d4SSatish Balay         PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
111027f02ce8SMatthew G. Knepley     }
111127f02ce8SMatthew G. Knepley     if (debug > 1) {
11124e913f38SMatthew G. Knepley       const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb));
11134e913f38SMatthew G. Knepley       const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb));
11144e913f38SMatthew G. Knepley       const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb));
11154e913f38SMatthew G. Knepley       const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb));
11164e913f38SMatthew G. Knepley       PetscInt       f, g;
111727f02ce8SMatthew G. Knepley 
11184e913f38SMatthew 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));
11194e913f38SMatthew G. Knepley       for (f = fS; f < fE; ++f) {
11204e913f38SMatthew G. Knepley         const PetscInt i = offsetI + f;
11214e913f38SMatthew G. Knepley         for (g = gS; g < gE; ++g) {
11224e913f38SMatthew G. Knepley           const PetscInt j = offsetJ + g;
11234e913f38SMatthew G. Knepley           PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", f, i, g, j);
11244e913f38SMatthew 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])));
112527f02ce8SMatthew G. Knepley         }
11269566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
112727f02ce8SMatthew G. Knepley       }
112827f02ce8SMatthew G. Knepley     }
112927f02ce8SMatthew G. Knepley     cOffset += totDim;
113027f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
113127f02ce8SMatthew G. Knepley     eOffset += PetscSqr(totDim);
113227f02ce8SMatthew G. Knepley   }
11333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
113427f02ce8SMatthew G. Knepley }
113527f02ce8SMatthew G. Knepley 
1136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1137d71ae5a4SJacob Faibussowitsch {
113820cf1dd8SToby Isaac   PetscFunctionBegin;
113920cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
114020cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
114120cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
114220cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
114320cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1144ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
114520cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1146afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
114720cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
114820cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
114927f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
115020cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
115120cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
115220cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
115327f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
11543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115520cf1dd8SToby Isaac }
115620cf1dd8SToby Isaac 
115720cf1dd8SToby Isaac /*MC
1158dce8aebaSBarry Smith   PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization
115920cf1dd8SToby Isaac 
116020cf1dd8SToby Isaac   Level: intermediate
116120cf1dd8SToby Isaac 
1162dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
116320cf1dd8SToby Isaac M*/
116420cf1dd8SToby Isaac 
1165d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1166d71ae5a4SJacob Faibussowitsch {
116720cf1dd8SToby Isaac   PetscFE_Basic *b;
116820cf1dd8SToby Isaac 
116920cf1dd8SToby Isaac   PetscFunctionBegin;
117020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
11714dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
117220cf1dd8SToby Isaac   fem->data = b;
117320cf1dd8SToby Isaac 
11749566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
11753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117620cf1dd8SToby Isaac }
1177