xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
49371c9d4SSatish Balay static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem) {
520cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *)fem->data;
620cf1dd8SToby Isaac 
720cf1dd8SToby Isaac   PetscFunctionBegin;
89566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
920cf1dd8SToby Isaac   PetscFunctionReturn(0);
1020cf1dd8SToby Isaac }
1120cf1dd8SToby Isaac 
129371c9d4SSatish Balay static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v) {
13d9bac1caSLisandro Dalcin   PetscInt        dim, Nc;
14d9bac1caSLisandro Dalcin   PetscSpace      basis = NULL;
15d9bac1caSLisandro Dalcin   PetscDualSpace  dual  = NULL;
16d9bac1caSLisandro Dalcin   PetscQuadrature quad  = NULL;
1720cf1dd8SToby Isaac 
1820cf1dd8SToby Isaac   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
209566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe, &Nc));
219566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &basis));
229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dual));
239566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(v));
2563a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc));
269566063dSJacob Faibussowitsch   if (basis) PetscCall(PetscSpaceView(basis, v));
279566063dSJacob Faibussowitsch   if (dual) PetscCall(PetscDualSpaceView(dual, v));
289566063dSJacob Faibussowitsch   if (quad) PetscCall(PetscQuadratureView(quad, v));
299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(v));
3020cf1dd8SToby Isaac   PetscFunctionReturn(0);
3120cf1dd8SToby Isaac }
3220cf1dd8SToby Isaac 
339371c9d4SSatish Balay static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v) {
3420cf1dd8SToby Isaac   PetscBool iascii;
3520cf1dd8SToby Isaac 
3620cf1dd8SToby Isaac   PetscFunctionBegin;
379566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii));
389566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscFEView_Basic_Ascii(fe, v));
3920cf1dd8SToby Isaac   PetscFunctionReturn(0);
4020cf1dd8SToby Isaac }
4120cf1dd8SToby Isaac 
4220cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
439371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem) {
44b9d4cb8dSJed Brown   PetscReal    *work;
4520cf1dd8SToby Isaac   PetscBLASInt *pivots;
4620cf1dd8SToby Isaac   PetscBLASInt  n, info;
4720cf1dd8SToby Isaac   PetscInt      pdim, j;
4820cf1dd8SToby Isaac 
4920cf1dd8SToby Isaac   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pdim * pdim, &fem->invV));
5220cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
5320cf1dd8SToby Isaac     PetscReal       *Bf;
5420cf1dd8SToby Isaac     PetscQuadrature  f;
5520cf1dd8SToby Isaac     const PetscReal *points, *weights;
5620cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
5720cf1dd8SToby Isaac 
589566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f));
599566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf));
619566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
6220cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
6320cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
64b9d4cb8dSJed Brown       fem->invV[j * pdim + k] = 0.0;
6520cf1dd8SToby Isaac 
6620cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
67b9d4cb8dSJed Brown         for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c];
6820cf1dd8SToby Isaac       }
6920cf1dd8SToby Isaac     }
709566063dSJacob Faibussowitsch     PetscCall(PetscFree(Bf));
7120cf1dd8SToby Isaac   }
72ea2bdf6dSBarry Smith 
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work));
7420cf1dd8SToby Isaac   n = pdim;
75792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
7663a3b9bcSJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info);
77792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
7863a3b9bcSJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)info);
799566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pivots, work));
8020cf1dd8SToby Isaac   PetscFunctionReturn(0);
8120cf1dd8SToby Isaac }
8220cf1dd8SToby Isaac 
839371c9d4SSatish Balay PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim) {
8420cf1dd8SToby Isaac   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim));
8620cf1dd8SToby Isaac   PetscFunctionReturn(0);
8720cf1dd8SToby Isaac }
8820cf1dd8SToby Isaac 
89b9d4cb8dSJed Brown /* Tensor contraction on the middle index,
90b9d4cb8dSJed Brown  *    C[m,n,p] = A[m,k,p] * B[k,n]
91b9d4cb8dSJed Brown  * where all matrices use C-style ordering.
92b9d4cb8dSJed Brown  */
939371c9d4SSatish Balay static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C) {
94b9d4cb8dSJed Brown   PetscInt i;
95b9d4cb8dSJed Brown 
96b9d4cb8dSJed Brown   PetscFunctionBegin;
97b9d4cb8dSJed Brown   for (i = 0; i < m; i++) {
98b9d4cb8dSJed Brown     PetscBLASInt n_, p_, k_, lda, ldb, ldc;
99b9d4cb8dSJed Brown     PetscReal    one = 1, zero = 0;
100b9d4cb8dSJed Brown     /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
101b9d4cb8dSJed Brown      * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
102b9d4cb8dSJed Brown      */
1039566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &n_));
1049566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(p, &p_));
1059566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(k, &k_));
106b9d4cb8dSJed Brown     lda = p_;
107b9d4cb8dSJed Brown     ldb = n_;
108b9d4cb8dSJed Brown     ldc = p_;
109792fecdfSBarry Smith     PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc));
110b9d4cb8dSJed Brown   }
1119566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2. * m * n * p * k));
112b9d4cb8dSJed Brown   PetscFunctionReturn(0);
113b9d4cb8dSJed Brown }
114b9d4cb8dSJed Brown 
1159371c9d4SSatish Balay PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) {
11620cf1dd8SToby Isaac   DM         dm;
11720cf1dd8SToby Isaac   PetscInt   pdim; /* Dimension of FE space P */
11820cf1dd8SToby Isaac   PetscInt   dim;  /* Spatial dimension */
11920cf1dd8SToby Isaac   PetscInt   Nc;   /* Field components */
120ef0bb6c7SMatthew G. Knepley   PetscReal *B    = K >= 0 ? T->T[0] : NULL;
121ef0bb6c7SMatthew G. Knepley   PetscReal *D    = K >= 1 ? T->T[1] : NULL;
122ef0bb6c7SMatthew G. Knepley   PetscReal *H    = K >= 2 ? T->T[2] : NULL;
123ef0bb6c7SMatthew G. Knepley   PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
12420cf1dd8SToby Isaac 
12520cf1dd8SToby Isaac   PetscFunctionBegin;
1269566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
1279566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1289566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
1299566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
13020cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
1319566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1329566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1339566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1349566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
135b9d4cb8dSJed Brown   /* Translate from prime to nodal basis */
13620cf1dd8SToby Isaac   if (B) {
137b9d4cb8dSJed Brown     /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
1389566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B));
13920cf1dd8SToby Isaac   }
14020cf1dd8SToby Isaac   if (D) {
141b9d4cb8dSJed Brown     /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
1429566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D));
14320cf1dd8SToby Isaac   }
14420cf1dd8SToby Isaac   if (H) {
145b9d4cb8dSJed Brown     /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
1469566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H));
14720cf1dd8SToby Isaac   }
1489566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1499566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1509566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
15120cf1dd8SToby Isaac   PetscFunctionReturn(0);
15220cf1dd8SToby Isaac }
15320cf1dd8SToby Isaac 
1549371c9d4SSatish Balay static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) {
15520cf1dd8SToby Isaac   const PetscInt     debug = 0;
1564bee2e38SMatthew G. Knepley   PetscFE            fe;
15720cf1dd8SToby Isaac   PetscPointFunc     obj_func;
15820cf1dd8SToby Isaac   PetscQuadrature    quad;
159ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
1604bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
16120cf1dd8SToby Isaac   const PetscScalar *constants;
16220cf1dd8SToby Isaac   PetscReal         *x;
163ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
16420cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
16520cf1dd8SToby Isaac   PetscBool          isAffine;
16620cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
16720cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
16820cf1dd8SToby Isaac 
16920cf1dd8SToby Isaac   PetscFunctionBegin;
1709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetObjective(ds, field, &obj_func));
17120cf1dd8SToby Isaac   if (!obj_func) PetscFunctionReturn(0);
1729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
1749566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
1759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
1789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
1799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
1809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
1819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
1829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1834bee2e38SMatthew G. Knepley   if (dsAux) {
1849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
1899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
19063a3b9bcSJacob 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);
19120cf1dd8SToby Isaac   }
1929566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
19363a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
19420cf1dd8SToby Isaac   Np       = cgeom->numPoints;
19520cf1dd8SToby Isaac   dE       = cgeom->dimEmbed;
19620cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
19720cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
1984bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
19920cf1dd8SToby Isaac 
20027f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
20127f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
20220cf1dd8SToby Isaac     if (isAffine) {
2034bee2e38SMatthew G. Knepley       fegeom.v    = x;
2044bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2057132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
2067132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
2077132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
20820cf1dd8SToby Isaac     }
2094bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
2104bee2e38SMatthew G. Knepley       PetscScalar integrand;
2114bee2e38SMatthew G. Knepley       PetscReal   w;
2124bee2e38SMatthew G. Knepley 
2134bee2e38SMatthew G. Knepley       if (isAffine) {
2147132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
2154bee2e38SMatthew G. Knepley       } else {
2164bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
2174bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
2184bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
2194bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
2204bee2e38SMatthew G. Knepley       }
2214bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
22220cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
22363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
2247be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2259566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
22620cf1dd8SToby Isaac #endif
22720cf1dd8SToby Isaac       }
22863a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
2299566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
2309566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
2314bee2e38SMatthew 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);
2324bee2e38SMatthew G. Knepley       integrand *= w;
23320cf1dd8SToby Isaac       integral[e * Nf + field] += integrand;
2349566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[field])));
23520cf1dd8SToby Isaac     }
23620cf1dd8SToby Isaac     cOffset += totDim;
23720cf1dd8SToby Isaac     cOffsetAux += totDimAux;
23820cf1dd8SToby Isaac   }
23920cf1dd8SToby Isaac   PetscFunctionReturn(0);
24020cf1dd8SToby Isaac }
24120cf1dd8SToby Isaac 
2429371c9d4SSatish Balay static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFunc obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) {
243afe6d6adSToby Isaac   const PetscInt     debug = 0;
2444bee2e38SMatthew G. Knepley   PetscFE            fe;
245afe6d6adSToby Isaac   PetscQuadrature    quad;
246ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2474bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
248afe6d6adSToby Isaac   const PetscScalar *constants;
249afe6d6adSToby Isaac   PetscReal         *x;
250ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
251afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
252afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
253afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
254afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
255afe6d6adSToby Isaac 
256afe6d6adSToby Isaac   PetscFunctionBegin;
257afe6d6adSToby Isaac   if (!obj_func) PetscFunctionReturn(0);
2589566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
2599566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
2609566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
2619566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2629566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2639566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
2649566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
2659566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
2669566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
2679566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
2689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
2694bee2e38SMatthew G. Knepley   if (dsAux) {
2709566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
2719566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2729566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
2739566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2749566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2759566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
276afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
2779566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
2789566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
27963a3b9bcSJacob 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);
280afe6d6adSToby Isaac   }
2819566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
28263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
283afe6d6adSToby Isaac   Np       = fgeom->numPoints;
284afe6d6adSToby Isaac   dE       = fgeom->dimEmbed;
285afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
286afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
2879f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
288afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
289ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
290ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
291ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
292ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
29327f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
29427f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
29527f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
29627f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
2974bee2e38SMatthew G. Knepley     if (isAffine) {
2984bee2e38SMatthew G. Knepley       fegeom.v    = x;
2994bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3007132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
3017132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
3027132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
3037132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
3049f209ee4SMatthew G. Knepley 
3057132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
3067132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
3077132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
3084bee2e38SMatthew G. Knepley     }
309afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
310afe6d6adSToby Isaac       PetscScalar integrand;
3114bee2e38SMatthew G. Knepley       PetscReal   w;
312afe6d6adSToby Isaac 
313afe6d6adSToby Isaac       if (isAffine) {
3147132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
315afe6d6adSToby Isaac       } else {
3163fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
3179f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
3189f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
3194bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
3204bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
3219f209ee4SMatthew G. Knepley 
3229f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
3239f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
3249f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
325afe6d6adSToby Isaac       }
3264bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
327afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
32863a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
329afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3309566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
331afe6d6adSToby Isaac #endif
332afe6d6adSToby Isaac       }
33363a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
3349566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3359566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3364bee2e38SMatthew 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);
3374bee2e38SMatthew G. Knepley       integrand *= w;
338afe6d6adSToby Isaac       integral[e * Nf + field] += integrand;
3399566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
340afe6d6adSToby Isaac     }
341afe6d6adSToby Isaac     cOffset += totDim;
342afe6d6adSToby Isaac     cOffsetAux += totDimAux;
343afe6d6adSToby Isaac   }
344afe6d6adSToby Isaac   PetscFunctionReturn(0);
345afe6d6adSToby Isaac }
346afe6d6adSToby Isaac 
3479371c9d4SSatish Balay 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[]) {
34820cf1dd8SToby Isaac   const PetscInt     debug = 0;
3496528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3504bee2e38SMatthew G. Knepley   PetscFE            fe;
3516528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3526528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
3536528b96dSMatthew G. Knepley   PetscPointFunc    *f0_func, *f1_func;
35420cf1dd8SToby Isaac   PetscQuadrature    quad;
355ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3564bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
35720cf1dd8SToby Isaac   const PetscScalar *constants;
35820cf1dd8SToby Isaac   PetscReal         *x;
359ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
360ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
36120cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3626587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
36320cf1dd8SToby Isaac 
36420cf1dd8SToby Isaac   PetscFunctionBegin;
3659566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
3669566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
3679566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
3689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
3709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
3739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
3749566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
3756528b96dSMatthew G. Knepley   if (!n0 && !n1) PetscFunctionReturn(0);
3769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
3779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
3789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
3799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
3809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3814bee2e38SMatthew G. Knepley   if (dsAux) {
3829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
3849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
3869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
3879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
38863a3b9bcSJacob 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);
38920cf1dd8SToby Isaac   }
3909566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
39163a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
39220cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
39363a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
39420cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
3954bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
39620cf1dd8SToby Isaac 
3976587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
3989566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
3999566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
40020cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4014bee2e38SMatthew G. Knepley       PetscReal w;
4024bee2e38SMatthew G. Knepley       PetscInt  c, d;
40320cf1dd8SToby Isaac 
4049566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
4054bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
4066587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
40763a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4087be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4099566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
41020cf1dd8SToby Isaac #endif
41120cf1dd8SToby Isaac       }
4129566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
4139566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
4146528b96dSMatthew 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]);
415ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
4166528b96dSMatthew 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]);
4179371c9d4SSatish Balay       for (c = 0; c < T[field]->Nc; ++c)
4189371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w;
419b8025e53SMatthew G. Knepley       if (debug) {
42063a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g\n", q, (double)quadWeights[q]));
421b8025e53SMatthew G. Knepley         if (debug > 2) {
42263a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
42363a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
4249566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
42563a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
42663a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
4279566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
428b8025e53SMatthew G. Knepley         }
429b8025e53SMatthew G. Knepley       }
43020cf1dd8SToby Isaac     }
4319566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
43220cf1dd8SToby Isaac     cOffset += totDim;
43320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
43420cf1dd8SToby Isaac   }
43520cf1dd8SToby Isaac   PetscFunctionReturn(0);
43620cf1dd8SToby Isaac }
43720cf1dd8SToby Isaac 
4389371c9d4SSatish Balay 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[]) {
43920cf1dd8SToby Isaac   const PetscInt     debug = 0;
44006d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
4414bee2e38SMatthew G. Knepley   PetscFE            fe;
44206d8a0d3SMatthew G. Knepley   PetscInt           n0, n1, i;
44306d8a0d3SMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
44420cf1dd8SToby Isaac   PetscQuadrature    quad;
445ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
4464bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
44720cf1dd8SToby Isaac   const PetscScalar *constants;
44820cf1dd8SToby Isaac   PetscReal         *x;
449ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
450ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
4516587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
45220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
4536587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
45420cf1dd8SToby Isaac 
45520cf1dd8SToby Isaac   PetscFunctionBegin;
4569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
4579566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
4589566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
4599566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4619566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4629566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4639566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4649566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
46506d8a0d3SMatthew G. Knepley   if (!n0 && !n1) PetscFunctionReturn(0);
4669566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4679566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
4709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4714bee2e38SMatthew G. Knepley   if (dsAux) {
4729566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
4739566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4749566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4759566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4769566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4779566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4787be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
4799566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
4809566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
48163a3b9bcSJacob 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);
48220cf1dd8SToby Isaac   }
483ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
4849566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
48563a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
48620cf1dd8SToby Isaac   dE         = fgeom->dimEmbed;
4876587ee25SMatthew G. Knepley   /* TODO FIX THIS */
4886587ee25SMatthew G. Knepley   fgeom->dim = dim - 1;
48963a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
49020cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4919f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
49220cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
4939f209ee4SMatthew G. Knepley 
4946587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
4959566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcI));
4969566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
49720cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4984bee2e38SMatthew G. Knepley       PetscReal w;
4994bee2e38SMatthew G. Knepley       PetscInt  c, d;
5004bee2e38SMatthew G. Knepley 
5019566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
5029566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
5034bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
50462bd480fSMatthew G. Knepley       if (debug > 1) {
5056587ee25SMatthew G. Knepley         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
50663a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5077be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5089566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5099566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
51020cf1dd8SToby Isaac #endif
51120cf1dd8SToby Isaac         }
51262bd480fSMatthew G. Knepley       }
5139566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
5149566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
51506d8a0d3SMatthew 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]);
5164bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w;
51706d8a0d3SMatthew 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]);
5189371c9d4SSatish Balay       for (c = 0; c < NcI; ++c)
5199371c9d4SSatish Balay         for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w;
52062bd480fSMatthew G. Knepley       if (debug) {
52163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
52262bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
52363a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c])));
52462bd480fSMatthew G. Knepley           if (n1) {
52563a3b9bcSJacob 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])));
5269566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
52762bd480fSMatthew G. Knepley           }
52862bd480fSMatthew G. Knepley         }
52962bd480fSMatthew G. Knepley       }
53020cf1dd8SToby Isaac     }
5319566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
53220cf1dd8SToby Isaac     cOffset += totDim;
53320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
53420cf1dd8SToby Isaac   }
53520cf1dd8SToby Isaac   PetscFunctionReturn(0);
53620cf1dd8SToby Isaac }
53720cf1dd8SToby Isaac 
53827f02ce8SMatthew G. Knepley /*
53927f02ce8SMatthew 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
54027f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
54127f02ce8SMatthew G. Knepley 
54227f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
54327f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
54427f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
54527f02ce8SMatthew 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()
54627f02ce8SMatthew G. Knepley */
5479371c9d4SSatish Balay static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) {
54827f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
5496528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
55027f02ce8SMatthew G. Knepley   PetscFE            fe;
5516528b96dSMatthew G. Knepley   PetscWeakForm      wf;
5526528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
5536528b96dSMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
55427f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
555665f567fSMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
55627f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
55727f02ce8SMatthew G. Knepley   const PetscScalar *constants;
55827f02ce8SMatthew G. Knepley   PetscReal         *x;
559665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
560665f567fSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
5616587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
56227f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
5636587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
56427f02ce8SMatthew G. Knepley 
56527f02ce8SMatthew G. Knepley   PetscFunctionBegin;
56627f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
5679566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
5699566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
5709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, s, &uOff));
5739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
5749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
5759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
5769566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
5776528b96dSMatthew G. Knepley   if (!n0 && !n1) PetscFunctionReturn(0);
5789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
5799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
5809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
58127f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
5829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
5839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
58427f02ce8SMatthew G. Knepley   if (dsAux) {
5859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
5869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
5879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
5899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
5909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
59101907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
5929566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
5939566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
59463a3b9bcSJacob 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);
59527f02ce8SMatthew G. Knepley   }
5969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
597665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
598c2b7495fSMatthew G. Knepley   NcS = NcI;
5999566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
60063a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
60127f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
60263a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
60327f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
60427f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
60527f02ce8SMatthew G. Knepley     const PetscInt face = fgeom->face[e][0];
60627f02ce8SMatthew G. Knepley 
6076587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6089566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcS));
6099566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
61027f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
61127f02ce8SMatthew G. Knepley       PetscReal w;
61227f02ce8SMatthew G. Knepley       PetscInt  c, d;
61327f02ce8SMatthew G. Knepley 
6149566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
61527f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
6166587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
61763a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
61827f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
6199566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
62027f02ce8SMatthew G. Knepley #endif
62127f02ce8SMatthew G. Knepley       }
622a4158a15SMatthew 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]));
62327f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
6249566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
6259566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, auxOnBd ? q : face * Nq + q, TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
6266528b96dSMatthew 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]);
62727f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
6289ee2af8cSMatthew 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]);
6299371c9d4SSatish Balay       for (c = 0; c < NcS; ++c)
6309371c9d4SSatish Balay         for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
63127f02ce8SMatthew G. Knepley     }
6329371c9d4SSatish Balay     if (isCohesiveField) {
6339371c9d4SSatish Balay       PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]);
6349371c9d4SSatish Balay     } else {
6359371c9d4SSatish Balay       PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]);
6369371c9d4SSatish Balay     }
63727f02ce8SMatthew G. Knepley     cOffset += totDim;
63827f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
63927f02ce8SMatthew G. Knepley   }
64027f02ce8SMatthew G. Knepley   PetscFunctionReturn(0);
64127f02ce8SMatthew G. Knepley }
64227f02ce8SMatthew G. Knepley 
6439371c9d4SSatish Balay 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[]) {
64420cf1dd8SToby Isaac   const PetscInt     debug = 0;
6454bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
6466528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6476528b96dSMatthew G. Knepley   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
6486528b96dSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
64920cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
65020cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
65120cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
65220cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
65320cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
65420cf1dd8SToby Isaac   PetscQuadrature    quad;
655ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
6564bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
65720cf1dd8SToby Isaac   const PetscScalar *constants;
65820cf1dd8SToby Isaac   PetscReal         *x;
659ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
660ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
6616528b96dSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
66220cf1dd8SToby Isaac   PetscInt           dE, Np;
66320cf1dd8SToby Isaac   PetscBool          isAffine;
66420cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
66520cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
66620cf1dd8SToby Isaac 
66720cf1dd8SToby Isaac   PetscFunctionBegin;
6689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6696528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
6706528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
6719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
6729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
6739566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
6749566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
6759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
6769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
6779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
6789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
67920cf1dd8SToby Isaac   switch (jtype) {
6809566063dSJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN: PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); break;
6819566063dSJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE: PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); break;
6829566063dSJacob Faibussowitsch   case PETSCFE_JACOBIAN: PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); break;
68320cf1dd8SToby Isaac   }
6846528b96dSMatthew G. Knepley   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
6859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
6869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
6879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
6889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
6899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
6909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
6919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
6924bee2e38SMatthew G. Knepley   if (dsAux) {
6939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
6989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
69963a3b9bcSJacob 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);
70020cf1dd8SToby Isaac   }
70127f02ce8SMatthew G. Knepley   NcI      = T[fieldI]->Nc;
70227f02ce8SMatthew G. Knepley   NcJ      = T[fieldJ]->Nc;
7034bee2e38SMatthew G. Knepley   Np       = cgeom->numPoints;
7044bee2e38SMatthew G. Knepley   dE       = cgeom->dimEmbed;
7054bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
70627f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
7079566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
7089566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7099566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7109566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7119566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
71263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
7134bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
7144bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
7154bee2e38SMatthew G. Knepley 
71627f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
71727f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
7184bee2e38SMatthew G. Knepley     if (isAffine) {
7194bee2e38SMatthew G. Knepley       fegeom.v    = x;
7204bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
7217132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
7227132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
7237132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
7244bee2e38SMatthew G. Knepley     }
72520cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
72620cf1dd8SToby Isaac       PetscReal w;
7274bee2e38SMatthew G. Knepley       PetscInt  c;
72820cf1dd8SToby Isaac 
72920cf1dd8SToby Isaac       if (isAffine) {
7307132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
73120cf1dd8SToby Isaac       } else {
7324bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
7334bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
7344bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
7354bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
73620cf1dd8SToby Isaac       }
7379566063dSJacob 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]));
7384bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
7399566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
7409566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
741ea672e62SMatthew G. Knepley       if (n0) {
7429566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
7436528b96dSMatthew 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);
74420cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
74520cf1dd8SToby Isaac       }
746ea672e62SMatthew G. Knepley       if (n1) {
7479566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
7486528b96dSMatthew 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);
7494bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
75020cf1dd8SToby Isaac       }
751ea672e62SMatthew G. Knepley       if (n2) {
7529566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
7536528b96dSMatthew 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);
7544bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
75520cf1dd8SToby Isaac       }
756ea672e62SMatthew G. Knepley       if (n3) {
7579566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
7586528b96dSMatthew 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);
7594bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
76020cf1dd8SToby Isaac       }
76120cf1dd8SToby Isaac 
7629566063dSJacob 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));
76320cf1dd8SToby Isaac     }
76420cf1dd8SToby Isaac     if (debug > 1) {
76520cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
76620cf1dd8SToby Isaac 
76763a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
768ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
769ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
770ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
771ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
772ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
773ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
77463a3b9bcSJacob 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])));
77520cf1dd8SToby Isaac             }
77620cf1dd8SToby Isaac           }
7779566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
77820cf1dd8SToby Isaac         }
77920cf1dd8SToby Isaac       }
78020cf1dd8SToby Isaac     }
78120cf1dd8SToby Isaac     cOffset += totDim;
78220cf1dd8SToby Isaac     cOffsetAux += totDimAux;
78320cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
78420cf1dd8SToby Isaac   }
78520cf1dd8SToby Isaac   PetscFunctionReturn(0);
78620cf1dd8SToby Isaac }
78720cf1dd8SToby Isaac 
7889371c9d4SSatish Balay static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) {
78920cf1dd8SToby Isaac   const PetscInt     debug = 0;
7904bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
79145480ffeSMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
79245480ffeSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
79320cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
79420cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
79520cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
79620cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
79720cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
79820cf1dd8SToby Isaac   PetscQuadrature    quad;
799ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
8004bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
80120cf1dd8SToby Isaac   const PetscScalar *constants;
80220cf1dd8SToby Isaac   PetscReal         *x;
803ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
804ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
80545480ffeSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
80620cf1dd8SToby Isaac   PetscBool          isAffine;
80720cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
80820cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
80920cf1dd8SToby Isaac 
81020cf1dd8SToby Isaac   PetscFunctionBegin;
8119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
81245480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
81345480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
8149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
8159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
8169566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
8179566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
8189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
8199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
8209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
8219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
8229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
8239566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
82445480ffeSMatthew G. Knepley   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
8259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
8269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
8279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
8289566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
8299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
8304bee2e38SMatthew G. Knepley   if (dsAux) {
8319566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8329566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
8339566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
8349566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
8359566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
8369566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
83720cf1dd8SToby Isaac   }
838ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
83920cf1dd8SToby Isaac   Np       = fgeom->numPoints;
84020cf1dd8SToby Isaac   dE       = fgeom->dimEmbed;
84120cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
84227f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
8439566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
8449566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
8459566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
8469566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8479566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
84863a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
84920cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
8509f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
85120cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
852ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
853ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
854ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
855ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
85627f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
85727f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
85827f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
85927f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
8604bee2e38SMatthew G. Knepley     if (isAffine) {
8614bee2e38SMatthew G. Knepley       fegeom.v    = x;
8624bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
8637132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
8647132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
8657132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
8667132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
8679f209ee4SMatthew G. Knepley 
8687132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
8697132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
8707132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
8714bee2e38SMatthew G. Knepley     }
87220cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
87320cf1dd8SToby Isaac       PetscReal w;
8744bee2e38SMatthew G. Knepley       PetscInt  c;
87520cf1dd8SToby Isaac 
87663a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
87720cf1dd8SToby Isaac       if (isAffine) {
8787132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
87920cf1dd8SToby Isaac       } else {
8803fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
8819f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
8829f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
8834bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
8844bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
8859f209ee4SMatthew G. Knepley 
8869f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
8879f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
8889f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
88920cf1dd8SToby Isaac       }
8904bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
8919566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
8929566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
893ea672e62SMatthew G. Knepley       if (n0) {
8949566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
89545480ffeSMatthew 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);
89620cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
89720cf1dd8SToby Isaac       }
898ea672e62SMatthew G. Knepley       if (n1) {
8999566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
90045480ffeSMatthew 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);
9014bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
90220cf1dd8SToby Isaac       }
903ea672e62SMatthew G. Knepley       if (n2) {
9049566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
90545480ffeSMatthew 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);
9064bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
90720cf1dd8SToby Isaac       }
908ea672e62SMatthew G. Knepley       if (n3) {
9099566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
91045480ffeSMatthew 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);
9114bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
91220cf1dd8SToby Isaac       }
91320cf1dd8SToby Isaac 
9149566063dSJacob 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));
91520cf1dd8SToby Isaac     }
91620cf1dd8SToby Isaac     if (debug > 1) {
91720cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
91820cf1dd8SToby Isaac 
91963a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
920ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
921ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
922ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
923ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
924ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
925ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
92663a3b9bcSJacob 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])));
92720cf1dd8SToby Isaac             }
92820cf1dd8SToby Isaac           }
9299566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
93020cf1dd8SToby Isaac         }
93120cf1dd8SToby Isaac       }
93220cf1dd8SToby Isaac     }
93320cf1dd8SToby Isaac     cOffset += totDim;
93420cf1dd8SToby Isaac     cOffsetAux += totDimAux;
93520cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
93620cf1dd8SToby Isaac   }
93720cf1dd8SToby Isaac   PetscFunctionReturn(0);
93820cf1dd8SToby Isaac }
93920cf1dd8SToby Isaac 
9409371c9d4SSatish Balay PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) {
94127f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
94227f02ce8SMatthew G. Knepley   PetscFE            feI, feJ;
943148442b3SMatthew G. Knepley   PetscWeakForm      wf;
944148442b3SMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
945148442b3SMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
94627f02ce8SMatthew G. Knepley   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
94727f02ce8SMatthew G. Knepley   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
94827f02ce8SMatthew G. Knepley   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
94927f02ce8SMatthew G. Knepley   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
95027f02ce8SMatthew G. Knepley   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
951665f567fSMatthew G. Knepley   PetscQuadrature    quad;
952665f567fSMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
95327f02ce8SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
95427f02ce8SMatthew G. Knepley   const PetscScalar *constants;
95527f02ce8SMatthew G. Knepley   PetscReal         *x;
956665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
957665f567fSMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
95845480ffeSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
959665f567fSMatthew G. Knepley   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE;
96027f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
96127f02ce8SMatthew G. Knepley   PetscInt           qNc, Nq, q, Np, dE;
96227f02ce8SMatthew G. Knepley 
96327f02ce8SMatthew G. Knepley   PetscFunctionBegin;
9649566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
96545480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
96645480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
96727f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
9689566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
9699566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
9709566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
9719566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
9729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
9739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, s, &uOff));
9749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
9759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
97627f02ce8SMatthew G. Knepley   switch (jtype) {
9779566063dSJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE: PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); break;
9789566063dSJacob Faibussowitsch   case PETSCFE_JACOBIAN: PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); break;
979665f567fSMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
98027f02ce8SMatthew G. Knepley   }
981148442b3SMatthew G. Knepley   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
9829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
9839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
9849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
9859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
9869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
9879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
9889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
98927f02ce8SMatthew G. Knepley   if (dsAux) {
9909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
9919566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
9929566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
9939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
9949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
9959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
99601907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
9979566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
9989566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
99963a3b9bcSJacob 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);
100027f02ce8SMatthew G. Knepley   }
10019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
10029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1003665f567fSMatthew G. Knepley   NcI      = T[fieldI]->Nc;
1004665f567fSMatthew G. Knepley   NcJ      = T[fieldJ]->Nc;
100527f02ce8SMatthew G. Knepley   NcS      = isCohesiveFieldI ? NcI : 2 * NcI;
100627f02ce8SMatthew G. Knepley   NcT      = isCohesiveFieldJ ? NcJ : 2 * NcJ;
100727f02ce8SMatthew G. Knepley   Np       = fgeom->numPoints;
100827f02ce8SMatthew G. Knepley   dE       = fgeom->dimEmbed;
100927f02ce8SMatthew G. Knepley   isAffine = fgeom->isAffine;
10109566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS * NcT));
10119566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcS * NcT * dE));
10129566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcS * NcT * dE));
10139566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcS * NcT * dE * dE));
10149566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
101563a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
101627f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
101727f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
101827f02ce8SMatthew G. Knepley     const PetscInt face = fgeom->face[e][0];
101927f02ce8SMatthew G. Knepley 
102027f02ce8SMatthew G. Knepley     fegeom.dim      = fgeom->dim;
102127f02ce8SMatthew G. Knepley     fegeom.dimEmbed = fgeom->dimEmbed;
102227f02ce8SMatthew G. Knepley     if (isAffine) {
102327f02ce8SMatthew G. Knepley       fegeom.v    = x;
102427f02ce8SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
1025a4158a15SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
1026a4158a15SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
1027a4158a15SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
1028a4158a15SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * dE * Np];
102927f02ce8SMatthew G. Knepley     }
103027f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
103127f02ce8SMatthew G. Knepley       PetscReal w;
103227f02ce8SMatthew G. Knepley       PetscInt  c;
103327f02ce8SMatthew G. Knepley 
103427f02ce8SMatthew G. Knepley       if (isAffine) {
103527f02ce8SMatthew G. Knepley         /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */
1036a4158a15SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
103727f02ce8SMatthew G. Knepley       } else {
103827f02ce8SMatthew G. Knepley         fegeom.v    = &fegeom.v[(e * Np + q) * dE];
103927f02ce8SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
104027f02ce8SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
104127f02ce8SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
104227f02ce8SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
104327f02ce8SMatthew G. Knepley       }
104427f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
104527f02ce8SMatthew G. Knepley       if (debug > 1 && q < Np) {
104663a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
104727f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
10489566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
104927f02ce8SMatthew G. Knepley #endif
105027f02ce8SMatthew G. Knepley       }
105163a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
10529566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
10539566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, auxOnBd ? q : face * Nq + q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1054ea672e62SMatthew G. Knepley       if (n0) {
10559566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS * NcT));
1056148442b3SMatthew 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);
105727f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT; ++c) g0[c] *= w;
105827f02ce8SMatthew G. Knepley       }
1059ea672e62SMatthew G. Knepley       if (n1) {
10609566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcS * NcT * dE));
1061148442b3SMatthew 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);
106227f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT * dE; ++c) g1[c] *= w;
106327f02ce8SMatthew G. Knepley       }
1064ea672e62SMatthew G. Knepley       if (n2) {
10659566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcS * NcT * dE));
1066148442b3SMatthew 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);
106727f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT * dE; ++c) g2[c] *= w;
106827f02ce8SMatthew G. Knepley       }
1069ea672e62SMatthew G. Knepley       if (n3) {
10709566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcS * NcT * dE * dE));
1071148442b3SMatthew 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);
107227f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT * dE * dE; ++c) g3[c] *= w;
107327f02ce8SMatthew G. Knepley       }
107427f02ce8SMatthew G. Knepley 
10755fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
10765fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
10779566063dSJacob 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));
107827f02ce8SMatthew G. Knepley         } else {
10799566063dSJacob 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));
10809566063dSJacob 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));
10815fedec97SMatthew G. Knepley         }
10829371c9d4SSatish Balay       } else
10839371c9d4SSatish 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));
108427f02ce8SMatthew G. Knepley     }
108527f02ce8SMatthew G. Knepley     if (debug > 1) {
108627f02ce8SMatthew G. Knepley       PetscInt fc, f, gc, g;
108727f02ce8SMatthew G. Knepley 
108863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
108927f02ce8SMatthew G. Knepley       for (fc = 0; fc < NcI; ++fc) {
1090665f567fSMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
109127f02ce8SMatthew G. Knepley           const PetscInt i = offsetI + f * NcI + fc;
109227f02ce8SMatthew G. Knepley           for (gc = 0; gc < NcJ; ++gc) {
1093665f567fSMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
109427f02ce8SMatthew G. Knepley               const PetscInt j = offsetJ + g * NcJ + gc;
109563a3b9bcSJacob 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])));
109627f02ce8SMatthew G. Knepley             }
109727f02ce8SMatthew G. Knepley           }
10989566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
109927f02ce8SMatthew G. Knepley         }
110027f02ce8SMatthew G. Knepley       }
110127f02ce8SMatthew G. Knepley     }
110227f02ce8SMatthew G. Knepley     cOffset += totDim;
110327f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
110427f02ce8SMatthew G. Knepley     eOffset += PetscSqr(totDim);
110527f02ce8SMatthew G. Knepley   }
110627f02ce8SMatthew G. Knepley   PetscFunctionReturn(0);
110727f02ce8SMatthew G. Knepley }
110827f02ce8SMatthew G. Knepley 
11099371c9d4SSatish Balay static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) {
111020cf1dd8SToby Isaac   PetscFunctionBegin;
111120cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
111220cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
111320cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
111420cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
111520cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1116ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
111720cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1118afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
111920cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
112020cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
112127f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
112220cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
112320cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
112420cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
112527f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
112620cf1dd8SToby Isaac   PetscFunctionReturn(0);
112720cf1dd8SToby Isaac }
112820cf1dd8SToby Isaac 
112920cf1dd8SToby Isaac /*MC
113020cf1dd8SToby Isaac   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
113120cf1dd8SToby Isaac 
113220cf1dd8SToby Isaac   Level: intermediate
113320cf1dd8SToby Isaac 
1134db781477SPatrick Sanan .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
113520cf1dd8SToby Isaac M*/
113620cf1dd8SToby Isaac 
11379371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) {
113820cf1dd8SToby Isaac   PetscFE_Basic *b;
113920cf1dd8SToby Isaac 
114020cf1dd8SToby Isaac   PetscFunctionBegin;
114120cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1142*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
114320cf1dd8SToby Isaac   fem->data = b;
114420cf1dd8SToby Isaac 
11459566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
114620cf1dd8SToby Isaac   PetscFunctionReturn(0);
114720cf1dd8SToby Isaac }
1148