xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
42b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
520cf1dd8SToby Isaac {
620cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *) fem->data;
720cf1dd8SToby Isaac 
820cf1dd8SToby Isaac   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
1020cf1dd8SToby Isaac   PetscFunctionReturn(0);
1120cf1dd8SToby Isaac }
1220cf1dd8SToby Isaac 
132b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
1420cf1dd8SToby Isaac {
15d9bac1caSLisandro Dalcin   PetscInt          dim, Nc;
16d9bac1caSLisandro Dalcin   PetscSpace        basis = NULL;
17d9bac1caSLisandro Dalcin   PetscDualSpace    dual = NULL;
18d9bac1caSLisandro Dalcin   PetscQuadrature   quad = NULL;
1920cf1dd8SToby Isaac 
2020cf1dd8SToby Isaac   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe, &Nc));
239566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &basis));
249566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dual));
259566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(v));
2763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n",dim,Nc));
289566063dSJacob Faibussowitsch   if (basis) PetscCall(PetscSpaceView(basis, v));
299566063dSJacob Faibussowitsch   if (dual)  PetscCall(PetscDualSpaceView(dual, v));
309566063dSJacob Faibussowitsch   if (quad)  PetscCall(PetscQuadratureView(quad, v));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(v));
3220cf1dd8SToby Isaac   PetscFunctionReturn(0);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac 
352b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
3620cf1dd8SToby Isaac {
3720cf1dd8SToby Isaac   PetscBool      iascii;
3820cf1dd8SToby Isaac 
3920cf1dd8SToby Isaac   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii));
419566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscFEView_Basic_Ascii(fe, v));
4220cf1dd8SToby Isaac   PetscFunctionReturn(0);
4320cf1dd8SToby Isaac }
4420cf1dd8SToby Isaac 
4520cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
46526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
4720cf1dd8SToby Isaac {
48b9d4cb8dSJed Brown   PetscReal     *work;
4920cf1dd8SToby Isaac   PetscBLASInt  *pivots;
5020cf1dd8SToby Isaac   PetscBLASInt   n, info;
5120cf1dd8SToby Isaac   PetscInt       pdim, j;
5220cf1dd8SToby Isaac 
5320cf1dd8SToby Isaac   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pdim*pdim,&fem->invV));
5620cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
5720cf1dd8SToby Isaac     PetscReal       *Bf;
5820cf1dd8SToby Isaac     PetscQuadrature  f;
5920cf1dd8SToby Isaac     const PetscReal *points, *weights;
6020cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
6120cf1dd8SToby Isaac 
629566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f));
639566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc*Nq*pdim,&Bf));
659566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
6620cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
6720cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
68b9d4cb8dSJed Brown       fem->invV[j*pdim+k] = 0.0;
6920cf1dd8SToby Isaac 
7020cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
71b9d4cb8dSJed Brown         for (c = 0; c < Nc; ++c) fem->invV[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c];
7220cf1dd8SToby Isaac       }
7320cf1dd8SToby Isaac     }
749566063dSJacob Faibussowitsch     PetscCall(PetscFree(Bf));
7520cf1dd8SToby Isaac   }
76ea2bdf6dSBarry Smith 
779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pdim,&pivots,pdim,&work));
7820cf1dd8SToby Isaac   n = pdim;
79*792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
8063a3b9bcSJacob Faibussowitsch   PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %" PetscInt_FMT,(PetscInt)info);
81*792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
8263a3b9bcSJacob Faibussowitsch   PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %" PetscInt_FMT,(PetscInt)info);
839566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pivots,work));
8420cf1dd8SToby Isaac   PetscFunctionReturn(0);
8520cf1dd8SToby Isaac }
8620cf1dd8SToby Isaac 
8720cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
8820cf1dd8SToby Isaac {
8920cf1dd8SToby Isaac   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim));
9120cf1dd8SToby Isaac   PetscFunctionReturn(0);
9220cf1dd8SToby Isaac }
9320cf1dd8SToby Isaac 
94b9d4cb8dSJed Brown /* Tensor contraction on the middle index,
95b9d4cb8dSJed Brown  *    C[m,n,p] = A[m,k,p] * B[k,n]
96b9d4cb8dSJed Brown  * where all matrices use C-style ordering.
97b9d4cb8dSJed Brown  */
98bdb10af2SPierre Jolivet static PetscErrorCode TensorContract_Private(PetscInt m,PetscInt n,PetscInt p,PetscInt k,const PetscReal *A,const PetscReal *B,PetscReal *C)
99bdb10af2SPierre Jolivet {
100b9d4cb8dSJed Brown   PetscInt i;
101b9d4cb8dSJed Brown 
102b9d4cb8dSJed Brown   PetscFunctionBegin;
103b9d4cb8dSJed Brown   for (i=0; i<m; i++) {
104b9d4cb8dSJed Brown     PetscBLASInt n_,p_,k_,lda,ldb,ldc;
105b9d4cb8dSJed Brown     PetscReal one = 1, zero = 0;
106b9d4cb8dSJed Brown     /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
107b9d4cb8dSJed Brown      * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
108b9d4cb8dSJed Brown      */
1099566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n,&n_));
1109566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(p,&p_));
1119566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(k,&k_));
112b9d4cb8dSJed Brown     lda = p_;
113b9d4cb8dSJed Brown     ldb = n_;
114b9d4cb8dSJed Brown     ldc = p_;
115*792fecdfSBarry Smith     PetscCallBLAS("BLASgemm",BLASREALgemm_("N","T",&p_,&n_,&k_,&one,A+i*k*p,&lda,B,&ldb,&zero,C+i*n*p,&ldc));
116b9d4cb8dSJed Brown   }
1179566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2.*m*n*p*k));
118b9d4cb8dSJed Brown   PetscFunctionReturn(0);
119b9d4cb8dSJed Brown }
120b9d4cb8dSJed Brown 
121526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
12220cf1dd8SToby Isaac {
12320cf1dd8SToby Isaac   DM               dm;
12420cf1dd8SToby Isaac   PetscInt         pdim; /* Dimension of FE space P */
12520cf1dd8SToby Isaac   PetscInt         dim;  /* Spatial dimension */
12620cf1dd8SToby Isaac   PetscInt         Nc;   /* Field components */
127ef0bb6c7SMatthew G. Knepley   PetscReal       *B = K >= 0 ? T->T[0] : NULL;
128ef0bb6c7SMatthew G. Knepley   PetscReal       *D = K >= 1 ? T->T[1] : NULL;
129ef0bb6c7SMatthew G. Knepley   PetscReal       *H = K >= 2 ? T->T[2] : NULL;
130ef0bb6c7SMatthew G. Knepley   PetscReal       *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
13120cf1dd8SToby Isaac 
13220cf1dd8SToby Isaac   PetscFunctionBegin;
1339566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
1349566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1359566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
13720cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
1389566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB));
1399566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD));
1409566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH));
1419566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
142b9d4cb8dSJed Brown   /* Translate from prime to nodal basis */
14320cf1dd8SToby Isaac   if (B) {
144b9d4cb8dSJed Brown     /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
1459566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B));
14620cf1dd8SToby Isaac   }
14720cf1dd8SToby Isaac   if (D) {
148b9d4cb8dSJed Brown     /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
1499566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc*dim, pdim, tmpD, fem->invV, D));
15020cf1dd8SToby Isaac   }
15120cf1dd8SToby Isaac   if (H) {
152b9d4cb8dSJed Brown     /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
1539566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc*dim*dim, pdim, tmpH, fem->invV, H));
15420cf1dd8SToby Isaac   }
1559566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB));
1569566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD));
1579566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH));
15820cf1dd8SToby Isaac   PetscFunctionReturn(0);
15920cf1dd8SToby Isaac }
16020cf1dd8SToby Isaac 
1612b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1624bee2e38SMatthew G. Knepley                                              const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
16320cf1dd8SToby Isaac {
16420cf1dd8SToby Isaac   const PetscInt     debug = 0;
1654bee2e38SMatthew G. Knepley   PetscFE            fe;
16620cf1dd8SToby Isaac   PetscPointFunc     obj_func;
16720cf1dd8SToby Isaac   PetscQuadrature    quad;
168ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
1694bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
17020cf1dd8SToby Isaac   const PetscScalar *constants;
17120cf1dd8SToby Isaac   PetscReal         *x;
172ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
17320cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
17420cf1dd8SToby Isaac   PetscBool          isAffine;
17520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
17620cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
17720cf1dd8SToby Isaac 
17820cf1dd8SToby Isaac   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetObjective(ds, field, &obj_func));
18020cf1dd8SToby Isaac   if (!obj_func) PetscFunctionReturn(0);
1819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *) &fe));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
1839566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
1849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
1879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
1889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
1899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
1909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1924bee2e38SMatthew G. Knepley   if (dsAux) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1949566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
1989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
19963a3b9bcSJacob Faibussowitsch     PetscCheck(T[0]->Np == TAux[0]->Np,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np);
20020cf1dd8SToby Isaac   }
2019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
20263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
20320cf1dd8SToby Isaac   Np = cgeom->numPoints;
20420cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
20520cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
20620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2074bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
20820cf1dd8SToby Isaac 
20927f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
21027f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
21120cf1dd8SToby Isaac     if (isAffine) {
2124bee2e38SMatthew G. Knepley       fegeom.v    = x;
2134bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2147132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*Np*dE*dE];
2157132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*Np*dE*dE];
2167132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e*Np];
21720cf1dd8SToby Isaac     }
2184bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
2194bee2e38SMatthew G. Knepley       PetscScalar integrand;
2204bee2e38SMatthew G. Knepley       PetscReal   w;
2214bee2e38SMatthew G. Knepley 
2224bee2e38SMatthew G. Knepley       if (isAffine) {
2237132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
2244bee2e38SMatthew G. Knepley       } else {
2254bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
2264bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
2274bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
2284bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
2294bee2e38SMatthew G. Knepley       }
2304bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
23120cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
23263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
2337be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2349566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
23520cf1dd8SToby Isaac #endif
23620cf1dd8SToby Isaac       }
23763a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
2389566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
2399566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
2404bee2e38SMatthew G. Knepley       obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, numConstants, constants, &integrand);
2414bee2e38SMatthew G. Knepley       integrand *= w;
24220cf1dd8SToby Isaac       integral[e*Nf+field] += integrand;
2439566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field])));
24420cf1dd8SToby Isaac     }
24520cf1dd8SToby Isaac     cOffset    += totDim;
24620cf1dd8SToby Isaac     cOffsetAux += totDimAux;
24720cf1dd8SToby Isaac   }
24820cf1dd8SToby Isaac   PetscFunctionReturn(0);
24920cf1dd8SToby Isaac }
25020cf1dd8SToby Isaac 
2512b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field,
252afe6d6adSToby Isaac                                                PetscBdPointFunc obj_func,
2534bee2e38SMatthew G. Knepley                                                PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
254afe6d6adSToby Isaac {
255afe6d6adSToby Isaac   const PetscInt     debug = 0;
2564bee2e38SMatthew G. Knepley   PetscFE            fe;
257afe6d6adSToby Isaac   PetscQuadrature    quad;
258ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2594bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
260afe6d6adSToby Isaac   const PetscScalar *constants;
261afe6d6adSToby Isaac   PetscReal         *x;
262ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
263afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
264afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
265afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
266afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
267afe6d6adSToby Isaac 
268afe6d6adSToby Isaac   PetscFunctionBegin;
269afe6d6adSToby Isaac   if (!obj_func) PetscFunctionReturn(0);
2709566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *) &fe));
2719566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
2729566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
2739566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
2769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
2779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
2789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
2799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
2809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
2814bee2e38SMatthew G. Knepley   if (dsAux) {
2829566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
2839566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
2859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
288afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
2899566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
2909566063dSJacob Faibussowitsch     else         PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
29163a3b9bcSJacob 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);
292afe6d6adSToby Isaac   }
2939566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
29463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
295afe6d6adSToby Isaac   Np = fgeom->numPoints;
296afe6d6adSToby Isaac   dE = fgeom->dimEmbed;
297afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
298afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
2999f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
300afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
301ea78f98cSLisandro Dalcin     fegeom.n = NULL;
302ea78f98cSLisandro Dalcin     fegeom.v = NULL;
303ea78f98cSLisandro Dalcin     fegeom.J = NULL;
304ea78f98cSLisandro Dalcin     fegeom.detJ = NULL;
30527f02ce8SMatthew G. Knepley     fegeom.dim      = fgeom->dim;
30627f02ce8SMatthew G. Knepley     fegeom.dimEmbed = fgeom->dimEmbed;
30727f02ce8SMatthew G. Knepley     cgeom.dim       = fgeom->dim;
30827f02ce8SMatthew G. Knepley     cgeom.dimEmbed  = fgeom->dimEmbed;
3094bee2e38SMatthew G. Knepley     if (isAffine) {
3104bee2e38SMatthew G. Knepley       fegeom.v    = x;
3114bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3127132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*Np*dE*dE];
3137132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
3147132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e*Np];
3157132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*Np*dE];
3169f209ee4SMatthew G. Knepley 
3177132c3f7SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*Np*dE*dE];
3187132c3f7SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*Np*dE*dE];
3197132c3f7SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e*Np];
3204bee2e38SMatthew G. Knepley     }
321afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
322afe6d6adSToby Isaac       PetscScalar integrand;
3234bee2e38SMatthew G. Knepley       PetscReal   w;
324afe6d6adSToby Isaac 
325afe6d6adSToby Isaac       if (isAffine) {
3267132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
327afe6d6adSToby Isaac       } else {
3283fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
3299f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
3309f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
3314bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
3324bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
3339f209ee4SMatthew G. Knepley 
3349f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
3359f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
3369f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
337afe6d6adSToby Isaac       }
3384bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
339afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
34063a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
341afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3429566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
343afe6d6adSToby Isaac #endif
344afe6d6adSToby Isaac       }
34563a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
3469566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3479566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3484bee2e38SMatthew 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);
3494bee2e38SMatthew G. Knepley       integrand *= w;
350afe6d6adSToby Isaac       integral[e*Nf+field] += integrand;
3519566063dSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field])));
352afe6d6adSToby Isaac     }
353afe6d6adSToby Isaac     cOffset    += totDim;
354afe6d6adSToby Isaac     cOffsetAux += totDimAux;
355afe6d6adSToby Isaac   }
356afe6d6adSToby Isaac   PetscFunctionReturn(0);
357afe6d6adSToby Isaac }
358afe6d6adSToby Isaac 
35906ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
3604bee2e38SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
36120cf1dd8SToby Isaac {
36220cf1dd8SToby Isaac   const PetscInt     debug = 0;
3636528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3644bee2e38SMatthew G. Knepley   PetscFE            fe;
3656528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3666528b96dSMatthew G. Knepley   PetscInt           n0,       n1, i;
3676528b96dSMatthew G. Knepley   PetscPointFunc    *f0_func, *f1_func;
36820cf1dd8SToby Isaac   PetscQuadrature    quad;
369ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3704bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
37120cf1dd8SToby Isaac   const PetscScalar *constants;
37220cf1dd8SToby Isaac   PetscReal         *x;
373ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
374ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
37520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3766587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
37720cf1dd8SToby Isaac 
37820cf1dd8SToby Isaac   PetscFunctionBegin;
3799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *) &fe));
3809566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
3819566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
3829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
3839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
3849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
3859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
3869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
3879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
3889566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
3896528b96dSMatthew G. Knepley   if (!n0 && !n1) PetscFunctionReturn(0);
3909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
3919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
3929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
3939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
3949566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
3954bee2e38SMatthew G. Knepley   if (dsAux) {
3969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
3979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
3989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
3999566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4009566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4019566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
40263a3b9bcSJacob 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);
40320cf1dd8SToby Isaac   }
4049566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
40563a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
40620cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
40763a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
40820cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4094bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
41020cf1dd8SToby Isaac 
4116587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
4129566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq*T[field]->Nc));
4139566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq*T[field]->Nc*dE));
41420cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4154bee2e38SMatthew G. Knepley       PetscReal w;
4164bee2e38SMatthew G. Knepley       PetscInt  c, d;
41720cf1dd8SToby Isaac 
4189566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q*cgeom->dim], &fegeom));
4194bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
4206587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
42163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4227be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4239566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
42420cf1dd8SToby Isaac #endif
42520cf1dd8SToby Isaac       }
4269566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
4279566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
4286528b96dSMatthew 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]);
429ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w;
4306528b96dSMatthew 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]);
431ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) for (d = 0; d < dim; ++d) f1[(q*T[field]->Nc+c)*dim+d] *= w;
432b8025e53SMatthew G. Knepley       if (debug) {
43363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g\n", q, (double)quadWeights[q]));
434b8025e53SMatthew G. Knepley         if (debug > 2) {
43563a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
43663a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field]+c])));
4379566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
43863a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
43963a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q*T[field]->Nc+c])));
4409566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
441b8025e53SMatthew G. Knepley         }
442b8025e53SMatthew G. Knepley       }
44320cf1dd8SToby Isaac     }
4449566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset+fOffset]));
44520cf1dd8SToby Isaac     cOffset    += totDim;
44620cf1dd8SToby Isaac     cOffsetAux += totDimAux;
44720cf1dd8SToby Isaac   }
44820cf1dd8SToby Isaac   PetscFunctionReturn(0);
44920cf1dd8SToby Isaac }
45020cf1dd8SToby Isaac 
45106ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
4524bee2e38SMatthew G. Knepley                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
45320cf1dd8SToby Isaac {
45420cf1dd8SToby Isaac   const PetscInt     debug = 0;
45506d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
4564bee2e38SMatthew G. Knepley   PetscFE            fe;
45706d8a0d3SMatthew G. Knepley   PetscInt           n0,       n1, i;
45806d8a0d3SMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
45920cf1dd8SToby Isaac   PetscQuadrature    quad;
460ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
4614bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
46220cf1dd8SToby Isaac   const PetscScalar *constants;
46320cf1dd8SToby Isaac   PetscReal         *x;
464ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
465ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
4666587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
46720cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
4686587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
46920cf1dd8SToby Isaac 
47020cf1dd8SToby Isaac   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *) &fe));
4729566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
4739566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
4749566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4759566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4799566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
48006d8a0d3SMatthew G. Knepley   if (!n0 && !n1) PetscFunctionReturn(0);
4819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
4859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4864bee2e38SMatthew G. Knepley   if (dsAux) {
4879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
4889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4919566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4929566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4937be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
4949566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
4959566063dSJacob Faibussowitsch     else         PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
49663a3b9bcSJacob 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);
49720cf1dd8SToby Isaac   }
498ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
4999566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
50063a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
50120cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
5026587ee25SMatthew G. Knepley   /* TODO FIX THIS */
5036587ee25SMatthew G. Knepley   fgeom->dim = dim-1;
50463a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
50520cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5069f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
50720cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5089f209ee4SMatthew G. Knepley 
5096587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
5109566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq*NcI));
5119566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq*NcI*dE));
51220cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5134bee2e38SMatthew G. Knepley       PetscReal w;
5144bee2e38SMatthew G. Knepley       PetscInt  c, d;
5154bee2e38SMatthew G. Knepley 
5169566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom));
5179566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
5184bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
51962bd480fSMatthew G. Knepley       if (debug > 1) {
5206587ee25SMatthew G. Knepley         if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) {
52163a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5227be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5239566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5249566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
52520cf1dd8SToby Isaac #endif
52620cf1dd8SToby Isaac         }
52762bd480fSMatthew G. Knepley       }
5289566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
5299566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
53006d8a0d3SMatthew 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]);
5314bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w;
53206d8a0d3SMatthew 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]);
5334bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w;
53462bd480fSMatthew G. Knepley       if (debug) {
53563a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
53662bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
53763a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q*NcI+c])));
53862bd480fSMatthew G. Knepley           if (n1) {
53963a3b9bcSJacob 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])));
5409566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
54162bd480fSMatthew G. Knepley           }
54262bd480fSMatthew G. Knepley         }
54362bd480fSMatthew G. Knepley       }
54420cf1dd8SToby Isaac     }
5459566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset]));
54620cf1dd8SToby Isaac     cOffset    += totDim;
54720cf1dd8SToby Isaac     cOffsetAux += totDimAux;
54820cf1dd8SToby Isaac   }
54920cf1dd8SToby Isaac   PetscFunctionReturn(0);
55020cf1dd8SToby Isaac }
55120cf1dd8SToby Isaac 
55227f02ce8SMatthew G. Knepley /*
55327f02ce8SMatthew 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
55427f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
55527f02ce8SMatthew G. Knepley 
55627f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
55727f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
55827f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
55927f02ce8SMatthew 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()
56027f02ce8SMatthew G. Knepley */
561c2b7495fSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom,
56227f02ce8SMatthew G. Knepley                                                            const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
56327f02ce8SMatthew G. Knepley {
56427f02ce8SMatthew G. Knepley   const PetscInt     debug = 0;
5656528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
56627f02ce8SMatthew G. Knepley   PetscFE            fe;
5676528b96dSMatthew G. Knepley   PetscWeakForm      wf;
5686528b96dSMatthew G. Knepley   PetscInt           n0,      n1, i;
5696528b96dSMatthew G. Knepley   PetscBdPointFunc  *f0_func, *f1_func;
57027f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
571665f567fSMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
57227f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
57327f02ce8SMatthew G. Knepley   const PetscScalar *constants;
57427f02ce8SMatthew G. Knepley   PetscReal         *x;
575665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
576665f567fSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
5776587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
57827f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
5796587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
58027f02ce8SMatthew G. Knepley 
58127f02ce8SMatthew G. Knepley   PetscFunctionBegin;
58227f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
5839566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *) &fe));
5849566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
5859566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
5869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, s, &uOff));
5899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
5909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
5919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
5929566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
5936528b96dSMatthew G. Knepley   if (!n0 && !n1) PetscFunctionReturn(0);
5949566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
5959566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
5969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
59727f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
5989566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
5999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
60027f02ce8SMatthew G. Knepley   if (dsAux) {
6019566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
6029566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6039566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6049566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6059566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6069566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
60701907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
6089566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
6099566063dSJacob Faibussowitsch     else         PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
61063a3b9bcSJacob Faibussowitsch     PetscCheck(Tf[0]->Np == TfAux[0]->Np,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np);
61127f02ce8SMatthew G. Knepley   }
6129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
613665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
614c2b7495fSMatthew G. Knepley   NcS = NcI;
6159566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
61663a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
61727f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
61863a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
61927f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
62027f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
62127f02ce8SMatthew G. Knepley     const PetscInt face = fgeom->face[e][0];
62227f02ce8SMatthew G. Knepley 
6236587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6249566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq*NcS));
6259566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq*NcS*dE));
62627f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
62727f02ce8SMatthew G. Knepley       PetscReal w;
62827f02ce8SMatthew G. Knepley       PetscInt  c, d;
62927f02ce8SMatthew G. Knepley 
6309566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom));
63127f02ce8SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
6326587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
63363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
63427f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
6359566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
63627f02ce8SMatthew G. Knepley #endif
63727f02ce8SMatthew G. Knepley       }
638a4158a15SMatthew 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]));
63927f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
6409566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
6419566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
6426528b96dSMatthew 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]);
64327f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w;
6449ee2af8cSMatthew 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]);
6459ee2af8cSMatthew G. Knepley       for (c = 0; c < NcS; ++c) for (d = 0; d < dE; ++d) f1[(q*NcS+c)*dE+d] *= w;
64627f02ce8SMatthew G. Knepley     }
6475fedec97SMatthew G. Knepley     if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset]);}
6485fedec97SMatthew G. Knepley     else                 {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset+fOffset]);}
64927f02ce8SMatthew G. Knepley     cOffset    += totDim;
65027f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
65127f02ce8SMatthew G. Knepley   }
65227f02ce8SMatthew G. Knepley   PetscFunctionReturn(0);
65327f02ce8SMatthew G. Knepley }
65427f02ce8SMatthew G. Knepley 
65506ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
6564bee2e38SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
65720cf1dd8SToby Isaac {
65820cf1dd8SToby Isaac   const PetscInt     debug      = 0;
6594bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
6606528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6616528b96dSMatthew G. Knepley   PetscPointJac     *g0_func, *g1_func, *g2_func, *g3_func;
6626528b96dSMatthew G. Knepley   PetscInt           n0, n1, n2, n3, i;
66320cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
66420cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
66520cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
66620cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
66720cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
66820cf1dd8SToby Isaac   PetscQuadrature    quad;
669ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
6704bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
67120cf1dd8SToby Isaac   const PetscScalar *constants;
67220cf1dd8SToby Isaac   PetscReal         *x;
673ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
674ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
6756528b96dSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
67620cf1dd8SToby Isaac   PetscInt           dE, Np;
67720cf1dd8SToby Isaac   PetscBool          isAffine;
67820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
67920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
68020cf1dd8SToby Isaac 
68120cf1dd8SToby Isaac   PetscFunctionBegin;
6829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6836528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
6846528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
6859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI));
6869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ));
6879566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
6889566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
6899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
6909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
6919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
6929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
69320cf1dd8SToby Isaac   switch(jtype) {
6949566063dSJacob 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;
6959566063dSJacob 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;
6969566063dSJacob 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;
69720cf1dd8SToby Isaac   }
6986528b96dSMatthew G. Knepley   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
6999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
7009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
7019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
7029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
7039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
7049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
7059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
7064bee2e38SMatthew G. Knepley   if (dsAux) {
7079566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
7089566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
7099566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
7109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
7119566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
7129566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
71363a3b9bcSJacob 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);
71420cf1dd8SToby Isaac   }
71527f02ce8SMatthew G. Knepley   NcI = T[fieldI]->Nc;
71627f02ce8SMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
7174bee2e38SMatthew G. Knepley   Np  = cgeom->numPoints;
7184bee2e38SMatthew G. Knepley   dE  = cgeom->dimEmbed;
7194bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
72027f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
7219566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI*NcJ));
7229566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI*NcJ*dE));
7239566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI*NcJ*dE));
7249566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI*NcJ*dE*dE));
7259566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
72663a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
7274bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
7284bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
7294bee2e38SMatthew G. Knepley 
73027f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
73127f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
7324bee2e38SMatthew G. Knepley     if (isAffine) {
7334bee2e38SMatthew G. Knepley       fegeom.v    = x;
7344bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
7357132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*Np*dE*dE];
7367132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*Np*dE*dE];
7377132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e*Np];
7384bee2e38SMatthew G. Knepley     }
73920cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
74020cf1dd8SToby Isaac       PetscReal w;
7414bee2e38SMatthew G. Knepley       PetscInt  c;
74220cf1dd8SToby Isaac 
74320cf1dd8SToby Isaac       if (isAffine) {
7447132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
74520cf1dd8SToby Isaac       } else {
7464bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
7474bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
7484bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
7494bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
75020cf1dd8SToby Isaac       }
7519566063dSJacob 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]));
7524bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
7539566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
7549566063dSJacob Faibussowitsch       if (dsAux)        PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
755ea672e62SMatthew G. Knepley       if (n0) {
7569566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI*NcJ));
7576528b96dSMatthew 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);
75820cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
75920cf1dd8SToby Isaac       }
760ea672e62SMatthew G. Knepley       if (n1) {
7619566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI*NcJ*dE));
7626528b96dSMatthew 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);
7634bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
76420cf1dd8SToby Isaac       }
765ea672e62SMatthew G. Knepley       if (n2) {
7669566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI*NcJ*dE));
7676528b96dSMatthew 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);
7684bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
76920cf1dd8SToby Isaac       }
770ea672e62SMatthew G. Knepley       if (n3) {
7719566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI*NcJ*dE*dE));
7726528b96dSMatthew 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);
7734bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
77420cf1dd8SToby Isaac       }
77520cf1dd8SToby Isaac 
7769566063dSJacob 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));
77720cf1dd8SToby Isaac     }
77820cf1dd8SToby Isaac     if (debug > 1) {
77920cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
78020cf1dd8SToby Isaac 
78163a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
782ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
783ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
784ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
785ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
786ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
787ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
78863a3b9bcSJacob 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])));
78920cf1dd8SToby Isaac             }
79020cf1dd8SToby Isaac           }
7919566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
79220cf1dd8SToby Isaac         }
79320cf1dd8SToby Isaac       }
79420cf1dd8SToby Isaac     }
79520cf1dd8SToby Isaac     cOffset    += totDim;
79620cf1dd8SToby Isaac     cOffsetAux += totDimAux;
79720cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
79820cf1dd8SToby Isaac   }
79920cf1dd8SToby Isaac   PetscFunctionReturn(0);
80020cf1dd8SToby Isaac }
80120cf1dd8SToby Isaac 
80206ad1575SMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
8034bee2e38SMatthew G. Knepley                                                        const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
80420cf1dd8SToby Isaac {
80520cf1dd8SToby Isaac   const PetscInt     debug      = 0;
8064bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
80745480ffeSMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
80845480ffeSMatthew G. Knepley   PetscInt           n0,       n1,       n2,       n3, i;
80920cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
81020cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
81120cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
81220cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
81320cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
81420cf1dd8SToby Isaac   PetscQuadrature    quad;
815ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
8164bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
81720cf1dd8SToby Isaac   const PetscScalar *constants;
81820cf1dd8SToby Isaac   PetscReal         *x;
819ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
820ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
82145480ffeSMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
82220cf1dd8SToby Isaac   PetscBool          isAffine;
82320cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
82420cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
82520cf1dd8SToby Isaac 
82620cf1dd8SToby Isaac   PetscFunctionBegin;
8279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
82845480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
82945480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
8309566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI));
8319566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ));
8329566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
8339566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
8349566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
8359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
8369566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
8379566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
8389566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
8399566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
84045480ffeSMatthew G. Knepley   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
8419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
8429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
8439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
8449566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
8459566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
8464bee2e38SMatthew G. Knepley   if (dsAux) {
8479566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8489566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
8499566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
8509566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
8519566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
8529566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
85320cf1dd8SToby Isaac   }
854ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
85520cf1dd8SToby Isaac   Np = fgeom->numPoints;
85620cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
85720cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
85827f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
8599566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI*NcJ));
8609566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI*NcJ*dE));
8619566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI*NcJ*dE));
8629566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI*NcJ*dE*dE));
8639566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
86463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
86520cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
8669f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
86720cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
868ea78f98cSLisandro Dalcin     fegeom.n = NULL;
869ea78f98cSLisandro Dalcin     fegeom.v = NULL;
870ea78f98cSLisandro Dalcin     fegeom.J = NULL;
871ea78f98cSLisandro Dalcin     fegeom.detJ = NULL;
87227f02ce8SMatthew G. Knepley     fegeom.dim      = fgeom->dim;
87327f02ce8SMatthew G. Knepley     fegeom.dimEmbed = fgeom->dimEmbed;
87427f02ce8SMatthew G. Knepley     cgeom.dim       = fgeom->dim;
87527f02ce8SMatthew G. Knepley     cgeom.dimEmbed  = fgeom->dimEmbed;
8764bee2e38SMatthew G. Knepley     if (isAffine) {
8774bee2e38SMatthew G. Knepley       fegeom.v    = x;
8784bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
8797132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*Np*dE*dE];
8807132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
8817132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e*Np];
8827132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*Np*dE];
8839f209ee4SMatthew G. Knepley 
8847132c3f7SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*Np*dE*dE];
8857132c3f7SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*Np*dE*dE];
8867132c3f7SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e*Np];
8874bee2e38SMatthew G. Knepley     }
88820cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
88920cf1dd8SToby Isaac       PetscReal w;
8904bee2e38SMatthew G. Knepley       PetscInt  c;
89120cf1dd8SToby Isaac 
89263a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
89320cf1dd8SToby Isaac       if (isAffine) {
8947132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
89520cf1dd8SToby Isaac       } else {
8963fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
8979f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
8989f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
8994bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
9004bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
9019f209ee4SMatthew G. Knepley 
9029f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
9039f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
9049f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
90520cf1dd8SToby Isaac       }
9064bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
9079566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
9089566063dSJacob Faibussowitsch       if (dsAux)        PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
909ea672e62SMatthew G. Knepley       if (n0) {
9109566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI*NcJ));
91145480ffeSMatthew 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);
91220cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
91320cf1dd8SToby Isaac       }
914ea672e62SMatthew G. Knepley       if (n1) {
9159566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI*NcJ*dE));
91645480ffeSMatthew 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);
9174bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
91820cf1dd8SToby Isaac       }
919ea672e62SMatthew G. Knepley       if (n2) {
9209566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI*NcJ*dE));
92145480ffeSMatthew 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);
9224bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
92320cf1dd8SToby Isaac       }
924ea672e62SMatthew G. Knepley       if (n3) {
9259566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI*NcJ*dE*dE));
92645480ffeSMatthew 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);
9274bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
92820cf1dd8SToby Isaac       }
92920cf1dd8SToby Isaac 
9309566063dSJacob 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));
93120cf1dd8SToby Isaac     }
93220cf1dd8SToby Isaac     if (debug > 1) {
93320cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
93420cf1dd8SToby Isaac 
93563a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
936ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
937ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
938ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
939ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
940ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
941ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
94263a3b9bcSJacob 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])));
94320cf1dd8SToby Isaac             }
94420cf1dd8SToby Isaac           }
9459566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
94620cf1dd8SToby Isaac         }
94720cf1dd8SToby Isaac       }
94820cf1dd8SToby Isaac     }
94920cf1dd8SToby Isaac     cOffset    += totDim;
95020cf1dd8SToby Isaac     cOffsetAux += totDimAux;
95120cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
95220cf1dd8SToby Isaac   }
95320cf1dd8SToby Isaac   PetscFunctionReturn(0);
95420cf1dd8SToby Isaac }
95520cf1dd8SToby Isaac 
9565fedec97SMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom,
95727f02ce8SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
95827f02ce8SMatthew G. Knepley {
95927f02ce8SMatthew G. Knepley   const PetscInt     debug      = 0;
96027f02ce8SMatthew G. Knepley   PetscFE            feI, feJ;
961148442b3SMatthew G. Knepley   PetscWeakForm      wf;
962148442b3SMatthew G. Knepley   PetscBdPointJac   *g0_func, *g1_func, *g2_func, *g3_func;
963148442b3SMatthew G. Knepley   PetscInt           n0,       n1,       n2,       n3, i;
96427f02ce8SMatthew G. Knepley   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
96527f02ce8SMatthew G. Knepley   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
96627f02ce8SMatthew G. Knepley   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
96727f02ce8SMatthew G. Knepley   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
96827f02ce8SMatthew G. Knepley   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
969665f567fSMatthew G. Knepley   PetscQuadrature    quad;
970665f567fSMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
97127f02ce8SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
97227f02ce8SMatthew G. Knepley   const PetscScalar *constants;
97327f02ce8SMatthew G. Knepley   PetscReal         *x;
974665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
975665f567fSMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0, NcS, NcT;
97645480ffeSMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
977665f567fSMatthew G. Knepley   PetscBool          isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE;
97827f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
97927f02ce8SMatthew G. Knepley   PetscInt           qNc, Nq, q, Np, dE;
98027f02ce8SMatthew G. Knepley 
98127f02ce8SMatthew G. Knepley   PetscFunctionBegin;
9829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
98345480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
98445480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
98527f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
9869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI));
9879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ));
9889566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
9899566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
9909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
9919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, s, &uOff));
9929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
9939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
99427f02ce8SMatthew G. Knepley   switch(jtype) {
9959566063dSJacob 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;
9969566063dSJacob 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;
997665f567fSMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
99827f02ce8SMatthew G. Knepley   }
999148442b3SMatthew G. Knepley   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0);
10009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
10019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
10029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
10039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
10049566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
10059566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
10069566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
100727f02ce8SMatthew G. Knepley   if (dsAux) {
10089566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
10099566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
10109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
10119566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
10129566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
10139566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
101401907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
10159566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
10169566063dSJacob Faibussowitsch     else         PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
101763a3b9bcSJacob 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);
101827f02ce8SMatthew G. Knepley   }
10199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
10209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
1021665f567fSMatthew G. Knepley   NcI = T[fieldI]->Nc;
1022665f567fSMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
102327f02ce8SMatthew G. Knepley   NcS = isCohesiveFieldI ? NcI : 2*NcI;
102427f02ce8SMatthew G. Knepley   NcT = isCohesiveFieldJ ? NcJ : 2*NcJ;
102527f02ce8SMatthew G. Knepley   Np = fgeom->numPoints;
102627f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
102727f02ce8SMatthew G. Knepley   isAffine = fgeom->isAffine;
10289566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS*NcT));
10299566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcS*NcT*dE));
10309566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcS*NcT*dE));
10319566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcS*NcT*dE*dE));
10329566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
103363a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
103427f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
103527f02ce8SMatthew G. Knepley     PetscFEGeom    fegeom;
103627f02ce8SMatthew G. Knepley     const PetscInt face = fgeom->face[e][0];
103727f02ce8SMatthew G. Knepley 
103827f02ce8SMatthew G. Knepley     fegeom.dim      = fgeom->dim;
103927f02ce8SMatthew G. Knepley     fegeom.dimEmbed = fgeom->dimEmbed;
104027f02ce8SMatthew G. Knepley     if (isAffine) {
104127f02ce8SMatthew G. Knepley       fegeom.v    = x;
104227f02ce8SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
1043a4158a15SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*Np*dE*dE];
1044a4158a15SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*Np*dE*dE];
1045a4158a15SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e*Np];
1046a4158a15SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE*Np];
104727f02ce8SMatthew G. Knepley     }
104827f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
104927f02ce8SMatthew G. Knepley       PetscReal w;
105027f02ce8SMatthew G. Knepley       PetscInt  c;
105127f02ce8SMatthew G. Knepley 
105227f02ce8SMatthew G. Knepley       if (isAffine) {
105327f02ce8SMatthew G. Knepley         /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */
1054a4158a15SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x);
105527f02ce8SMatthew G. Knepley       } else {
105627f02ce8SMatthew G. Knepley         fegeom.v    = &fegeom.v[(e*Np+q)*dE];
105727f02ce8SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
105827f02ce8SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
105927f02ce8SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
106027f02ce8SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
106127f02ce8SMatthew G. Knepley       }
106227f02ce8SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
106327f02ce8SMatthew G. Knepley       if (debug > 1 && q < Np) {
106463a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
106527f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
10669566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
106727f02ce8SMatthew G. Knepley #endif
106827f02ce8SMatthew G. Knepley       }
106963a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
10709566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
10719566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1072ea672e62SMatthew G. Knepley       if (n0) {
10739566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS*NcT));
1074148442b3SMatthew 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);
107527f02ce8SMatthew G. Knepley         for (c = 0; c < NcS*NcT; ++c) g0[c] *= w;
107627f02ce8SMatthew G. Knepley       }
1077ea672e62SMatthew G. Knepley       if (n1) {
10789566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcS*NcT*dE));
1079148442b3SMatthew 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);
108027f02ce8SMatthew G. Knepley         for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w;
108127f02ce8SMatthew G. Knepley       }
1082ea672e62SMatthew G. Knepley       if (n2) {
10839566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcS*NcT*dE));
1084148442b3SMatthew 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);
108527f02ce8SMatthew G. Knepley         for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w;
108627f02ce8SMatthew G. Knepley       }
1087ea672e62SMatthew G. Knepley       if (n3) {
10889566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcS*NcT*dE*dE));
1089148442b3SMatthew 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);
109027f02ce8SMatthew G. Knepley         for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w;
109127f02ce8SMatthew G. Knepley       }
109227f02ce8SMatthew G. Knepley 
10935fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
10945fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
10959566063dSJacob 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));
109627f02ce8SMatthew G. Knepley         } else {
10979566063dSJacob 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));
10989566063dSJacob 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));
10995fedec97SMatthew G. Knepley         }
11001baa6e33SBarry Smith       } else 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));
110127f02ce8SMatthew G. Knepley     }
110227f02ce8SMatthew G. Knepley     if (debug > 1) {
110327f02ce8SMatthew G. Knepley       PetscInt fc, f, gc, g;
110427f02ce8SMatthew G. Knepley 
110563a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
110627f02ce8SMatthew G. Knepley       for (fc = 0; fc < NcI; ++fc) {
1107665f567fSMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
110827f02ce8SMatthew G. Knepley           const PetscInt i = offsetI + f*NcI+fc;
110927f02ce8SMatthew G. Knepley           for (gc = 0; gc < NcJ; ++gc) {
1110665f567fSMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
111127f02ce8SMatthew G. Knepley               const PetscInt j = offsetJ + g*NcJ+gc;
111263a3b9bcSJacob 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])));
111327f02ce8SMatthew G. Knepley             }
111427f02ce8SMatthew G. Knepley           }
11159566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
111627f02ce8SMatthew G. Knepley         }
111727f02ce8SMatthew G. Knepley       }
111827f02ce8SMatthew G. Knepley     }
111927f02ce8SMatthew G. Knepley     cOffset    += totDim;
112027f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
112127f02ce8SMatthew G. Knepley     eOffset    += PetscSqr(totDim);
112227f02ce8SMatthew G. Knepley   }
112327f02ce8SMatthew G. Knepley   PetscFunctionReturn(0);
112427f02ce8SMatthew G. Knepley }
112527f02ce8SMatthew G. Knepley 
11262b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
112720cf1dd8SToby Isaac {
112820cf1dd8SToby Isaac   PetscFunctionBegin;
112920cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
113020cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
113120cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
113220cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
113320cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1134ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
113520cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1136afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
113720cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
113820cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
113927f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
114020cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
114120cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
114220cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
114327f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
114420cf1dd8SToby Isaac   PetscFunctionReturn(0);
114520cf1dd8SToby Isaac }
114620cf1dd8SToby Isaac 
114720cf1dd8SToby Isaac /*MC
114820cf1dd8SToby Isaac   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
114920cf1dd8SToby Isaac 
115020cf1dd8SToby Isaac   Level: intermediate
115120cf1dd8SToby Isaac 
1152db781477SPatrick Sanan .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
115320cf1dd8SToby Isaac M*/
115420cf1dd8SToby Isaac 
115520cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
115620cf1dd8SToby Isaac {
115720cf1dd8SToby Isaac   PetscFE_Basic *b;
115820cf1dd8SToby Isaac 
115920cf1dd8SToby Isaac   PetscFunctionBegin;
116020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
11619566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(fem,&b));
116220cf1dd8SToby Isaac   fem->data = b;
116320cf1dd8SToby Isaac 
11649566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
116520cf1dd8SToby Isaac   PetscFunctionReturn(0);
116620cf1dd8SToby Isaac }
1167