xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 9d47de495d3c23378050c1b4a410c12a375cb6c6)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
PetscFEDestroy_Basic(PetscFE fem)4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
5d71ae5a4SJacob Faibussowitsch {
620cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *)fem->data;
720cf1dd8SToby Isaac 
820cf1dd8SToby Isaac   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1120cf1dd8SToby Isaac }
1220cf1dd8SToby Isaac 
PetscFEView_Basic_Ascii(PetscFE fe,PetscViewer v)13d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
14d71ae5a4SJacob Faibussowitsch {
15d9bac1caSLisandro Dalcin   PetscInt        dim, Nc;
16d9bac1caSLisandro Dalcin   PetscSpace      basis = NULL;
17d9bac1caSLisandro Dalcin   PetscDualSpace  dual  = NULL;
18d9bac1caSLisandro Dalcin   PetscQuadrature quad  = NULL;
1920cf1dd8SToby Isaac 
2020cf1dd8SToby Isaac   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe, &Nc));
239566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &basis));
249566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dual));
259566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(v));
2763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc));
289566063dSJacob Faibussowitsch   if (basis) PetscCall(PetscSpaceView(basis, v));
299566063dSJacob Faibussowitsch   if (dual) PetscCall(PetscDualSpaceView(dual, v));
309566063dSJacob Faibussowitsch   if (quad) PetscCall(PetscQuadratureView(quad, v));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(v));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac 
PetscFEView_Basic(PetscFE fe,PetscViewer v)35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
36d71ae5a4SJacob Faibussowitsch {
379f196a02SMartin Diehl   PetscBool isascii;
3820cf1dd8SToby Isaac 
3920cf1dd8SToby Isaac   PetscFunctionBegin;
409f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
419f196a02SMartin Diehl   if (isascii) PetscCall(PetscFEView_Basic_Ascii(fe, v));
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4320cf1dd8SToby Isaac }
4420cf1dd8SToby Isaac 
4520cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
PetscFESetUp_Basic(PetscFE fem)46d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
47d71ae5a4SJacob Faibussowitsch {
48b9d4cb8dSJed Brown   PetscReal    *work;
4920cf1dd8SToby Isaac   PetscBLASInt *pivots;
5020cf1dd8SToby Isaac   PetscBLASInt  n, info;
5120cf1dd8SToby Isaac   PetscInt      pdim, j;
5220cf1dd8SToby Isaac 
5320cf1dd8SToby Isaac   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pdim * pdim, &fem->invV));
5620cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
5720cf1dd8SToby Isaac     PetscReal       *Bf;
5820cf1dd8SToby Isaac     PetscQuadrature  f;
5920cf1dd8SToby Isaac     const PetscReal *points, *weights;
6020cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
6120cf1dd8SToby Isaac 
629566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f));
639566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights));
649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf));
659566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL));
6620cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
6720cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
68b9d4cb8dSJed Brown       fem->invV[j * pdim + k] = 0.0;
6920cf1dd8SToby Isaac 
7020cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
71b9d4cb8dSJed Brown         for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c];
7220cf1dd8SToby Isaac       }
7320cf1dd8SToby Isaac     }
749566063dSJacob Faibussowitsch     PetscCall(PetscFree(Bf));
7520cf1dd8SToby Isaac   }
76ea2bdf6dSBarry Smith 
779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work));
786497c311SBarry Smith   PetscCall(PetscBLASIntCast(pdim, &n));
79792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info));
80835f2295SStefano Zampini   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info);
81792fecdfSBarry Smith   PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info));
82835f2295SStefano Zampini   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info);
839566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pivots, work));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8520cf1dd8SToby Isaac }
8620cf1dd8SToby Isaac 
PetscFEGetDimension_Basic(PetscFE fem,PetscInt * dim)87d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
88d71ae5a4SJacob Faibussowitsch {
8920cf1dd8SToby Isaac   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9220cf1dd8SToby Isaac }
9320cf1dd8SToby Isaac 
94b9d4cb8dSJed Brown /* Tensor contraction on the middle index,
95b9d4cb8dSJed Brown  *    C[m,n,p] = A[m,k,p] * B[k,n]
96b9d4cb8dSJed Brown  * where all matrices use C-style ordering.
97b9d4cb8dSJed Brown  */
TensorContract_Private(PetscInt m,PetscInt n,PetscInt p,PetscInt k,const PetscReal * A,const PetscReal * B,PetscReal * C)98d71ae5a4SJacob Faibussowitsch static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C)
99d71ae5a4SJacob Faibussowitsch {
100b9d4cb8dSJed Brown   PetscInt i;
101b9d4cb8dSJed Brown 
102b9d4cb8dSJed Brown   PetscFunctionBegin;
103aa9788aaSMatthew G. Knepley   PetscCheck(n && p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Empty tensor is not allowed %" PetscInt_FMT " %" PetscInt_FMT, n, p);
104b9d4cb8dSJed Brown   for (i = 0; i < m; i++) {
105b9d4cb8dSJed Brown     PetscBLASInt n_, p_, k_, lda, ldb, ldc;
106b9d4cb8dSJed Brown     PetscReal    one = 1, zero = 0;
107b9d4cb8dSJed Brown     /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n]
108b9d4cb8dSJed Brown      * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k)
109b9d4cb8dSJed Brown      */
1109566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(n, &n_));
1119566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(p, &p_));
1129566063dSJacob Faibussowitsch     PetscCall(PetscBLASIntCast(k, &k_));
113b9d4cb8dSJed Brown     lda = p_;
114b9d4cb8dSJed Brown     ldb = n_;
115b9d4cb8dSJed Brown     ldc = p_;
116792fecdfSBarry Smith     PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc));
117b9d4cb8dSJed Brown   }
1189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(2. * m * n * p * k));
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
120b9d4cb8dSJed Brown }
121b9d4cb8dSJed Brown 
PetscFEComputeTabulation_Basic(PetscFE fem,PetscInt npoints,const PetscReal points[],PetscInt K,PetscTabulation T)122ce78bad3SBarry Smith PETSC_INTERN PetscErrorCode PetscFEComputeTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
123d71ae5a4SJacob Faibussowitsch {
12420cf1dd8SToby Isaac   DM         dm;
12520cf1dd8SToby Isaac   PetscInt   pdim; /* Dimension of FE space P */
12620cf1dd8SToby Isaac   PetscInt   dim;  /* Spatial dimension */
12720cf1dd8SToby Isaac   PetscInt   Nc;   /* Field components */
128ef0bb6c7SMatthew G. Knepley   PetscReal *B    = K >= 0 ? T->T[0] : NULL;
129ef0bb6c7SMatthew G. Knepley   PetscReal *D    = K >= 1 ? T->T[1] : NULL;
130ef0bb6c7SMatthew G. Knepley   PetscReal *H    = K >= 2 ? T->T[2] : NULL;
131ef0bb6c7SMatthew G. Knepley   PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
13220cf1dd8SToby Isaac 
13320cf1dd8SToby Isaac   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
1359566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1369566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
13820cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
1399566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1409566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1419566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1429566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH));
143b9d4cb8dSJed Brown   /* Translate from prime to nodal basis */
14420cf1dd8SToby Isaac   if (B) {
145b9d4cb8dSJed Brown     /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */
1469566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B));
14720cf1dd8SToby Isaac   }
148aa9788aaSMatthew G. Knepley   if (D && dim) {
149b9d4cb8dSJed Brown     /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */
1509566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D));
15120cf1dd8SToby Isaac   }
152aa9788aaSMatthew G. Knepley   if (H && dim) {
153b9d4cb8dSJed Brown     /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */
1549566063dSJacob Faibussowitsch     PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H));
15520cf1dd8SToby Isaac   }
1569566063dSJacob Faibussowitsch   if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB));
1579566063dSJacob Faibussowitsch   if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD));
1589566063dSJacob Faibussowitsch   if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH));
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16020cf1dd8SToby Isaac }
16120cf1dd8SToby Isaac 
PetscFEIntegrate_Basic(PetscDS ds,PetscInt field,PetscInt Ne,PetscFEGeom * cgeom,const PetscScalar coefficients[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscScalar integral[])1622dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
163d71ae5a4SJacob Faibussowitsch {
164b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
1654bee2e38SMatthew G. Knepley   PetscFE            fe;
1662192575eSBarry Smith   PetscPointFn      *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;
17187510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
172ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
17320cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
17420cf1dd8SToby Isaac   PetscBool          isAffine;
17520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
17620cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
17720cf1dd8SToby Isaac 
17820cf1dd8SToby Isaac   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetObjective(ds, field, &obj_func));
1803ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
1819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
18387510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
1849566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
1859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
1869566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1879566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
1889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
1899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
1909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
1919566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL));
19287510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
1939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
1944bee2e38SMatthew G. Knepley   if (dsAux) {
1959566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
1969566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
1979566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
1989566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
1999566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
2009566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
20163a3b9bcSJacob 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);
20220cf1dd8SToby Isaac   }
2039566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
20463a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
20520cf1dd8SToby Isaac   Np       = cgeom->numPoints;
20620cf1dd8SToby Isaac   dE       = cgeom->dimEmbed;
20720cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
20820cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2094bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
21020cf1dd8SToby Isaac 
21127f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
21227f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
213a7be0fb8SMatthew G. Knepley     fegeom.xi       = NULL;
21420cf1dd8SToby Isaac     if (isAffine) {
2154bee2e38SMatthew G. Knepley       fegeom.v    = x;
2164bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2177132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
2187132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
2197132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
22025c09e31SMark Adams     } else fegeom.xi = NULL;
2214bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
222d627b919SMatthew G. Knepley       PetscScalar integrand = 0.;
2234bee2e38SMatthew G. Knepley       PetscReal   w;
2244bee2e38SMatthew G. Knepley 
2254bee2e38SMatthew G. Knepley       if (isAffine) {
2267132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
2274bee2e38SMatthew G. Knepley       } else {
2284bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
2294bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
2304bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
2314bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
2324bee2e38SMatthew G. Knepley       }
23387510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
2344bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
23520cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
23663a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
2377be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2389566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
23920cf1dd8SToby Isaac #endif
24020cf1dd8SToby Isaac       }
24163a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
2429566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL));
2439566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
2444bee2e38SMatthew 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);
2454bee2e38SMatthew G. Knepley       integrand *= w;
24620cf1dd8SToby Isaac       integral[e * Nf + field] += integrand;
24720cf1dd8SToby Isaac     }
248699b48e5SStefano Zampini     if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    Element Field %" PetscInt_FMT " integral: %g\n", Nf, (double)PetscRealPart(integral[e * Nf + field])));
24920cf1dd8SToby Isaac     cOffset += totDim;
25020cf1dd8SToby Isaac     cOffsetAux += totDimAux;
25120cf1dd8SToby Isaac   }
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25320cf1dd8SToby Isaac }
25420cf1dd8SToby Isaac 
PetscFEIntegrateBd_Basic(PetscDS ds,PetscInt field,PetscBdPointFn * obj_func,PetscInt Ne,PetscFEGeom * fgeom,const PetscScalar coefficients[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscScalar integral[])2552192575eSBarry Smith PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFn *obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
256d71ae5a4SJacob Faibussowitsch {
257b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
2584bee2e38SMatthew G. Knepley   PetscFE            fe;
259afe6d6adSToby Isaac   PetscQuadrature    quad;
260ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2614bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
262afe6d6adSToby Isaac   const PetscScalar *constants;
26387510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
264ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
265afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
266afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
267afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
268afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
269afe6d6adSToby Isaac 
270afe6d6adSToby Isaac   PetscFunctionBegin;
2713ba16761SJacob Faibussowitsch   if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS);
2729566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
2739566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
27487510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
2759566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
2769566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
2779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
2789566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
2799566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
2809566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x));
2819566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
2829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
28387510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
2849566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
2854bee2e38SMatthew G. Knepley   if (dsAux) {
2869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
2879566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
2889566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
2899566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
2909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
2919566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
292afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
2939566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
2949566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
29563a3b9bcSJacob 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);
296afe6d6adSToby Isaac   }
2979566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
29863a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
29979ab67a3SMatthew G. Knepley   if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Field: %" PetscInt_FMT " Nface: %" PetscInt_FMT " Nq: %" PetscInt_FMT "\n", field, Ne, Nq));
300afe6d6adSToby Isaac   Np       = fgeom->numPoints;
301afe6d6adSToby Isaac   dE       = fgeom->dimEmbed;
302afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
303afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
3049f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
305afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
306ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
307ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
308a7be0fb8SMatthew G. Knepley     fegeom.xi           = NULL;
309ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
310b2deab97SMatthew G. Knepley     fegeom.invJ         = NULL;
311ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
31227f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
31327f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
31427f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
31527f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
3164bee2e38SMatthew G. Knepley     if (isAffine) {
3174bee2e38SMatthew G. Knepley       fegeom.v    = x;
3184bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3197132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
3207132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
3217132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
3227132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
3239f209ee4SMatthew G. Knepley 
3247132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
3257132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
3267132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
32725c09e31SMark Adams     } else fegeom.xi = NULL;
328afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
32979ab67a3SMatthew G. Knepley       PetscScalar integrand = 0.;
3304bee2e38SMatthew G. Knepley       PetscReal   w;
331afe6d6adSToby Isaac 
332afe6d6adSToby Isaac       if (isAffine) {
3337132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
334afe6d6adSToby Isaac       } else {
3353fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
3369f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
3379f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
3384bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
3394bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
3409f209ee4SMatthew G. Knepley 
3419f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
3429f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
3439f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
344afe6d6adSToby Isaac       }
34587510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
3464bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
347afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
34863a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
349*beceaeb6SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3509566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
351afe6d6adSToby Isaac #endif
352afe6d6adSToby Isaac       }
35363a3b9bcSJacob Faibussowitsch       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
35479ab67a3SMatthew G. Knepley       if (debug > 3) {
35579ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    x_q ("));
35679ab67a3SMatthew G. Knepley         for (PetscInt d = 0; d < dE; ++d) {
35779ab67a3SMatthew G. Knepley           if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
35879ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.v[d]));
35979ab67a3SMatthew G. Knepley         }
36079ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
36179ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "    n_q ("));
36279ab67a3SMatthew G. Knepley         for (PetscInt d = 0; d < dE; ++d) {
36379ab67a3SMatthew G. Knepley           if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
36479ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.n[d]));
36579ab67a3SMatthew G. Knepley         }
36679ab67a3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
36779ab67a3SMatthew G. Knepley         for (PetscInt f = 0; f < Nf; ++f) {
36879ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "    u_%" PetscInt_FMT " (", f));
36979ab67a3SMatthew G. Knepley           for (PetscInt c = 0; c < uOff[f + 1] - uOff[f]; ++c) {
37079ab67a3SMatthew G. Knepley             if (c) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
37179ab67a3SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(u[uOff[f] + c])));
37279ab67a3SMatthew G. Knepley           }
37379ab67a3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
37479ab67a3SMatthew G. Knepley         }
37579ab67a3SMatthew G. Knepley       }
3769566063dSJacob Faibussowitsch       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL));
3779566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
3784bee2e38SMatthew 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);
3794bee2e38SMatthew G. Knepley       integrand *= w;
380afe6d6adSToby Isaac       integral[e * Nf + field] += integrand;
38179ab67a3SMatthew G. Knepley       if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "    int: %g tot: %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field])));
382afe6d6adSToby Isaac     }
383afe6d6adSToby Isaac     cOffset += totDim;
384afe6d6adSToby Isaac     cOffsetAux += totDimAux;
385afe6d6adSToby Isaac   }
3863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
387afe6d6adSToby Isaac }
388afe6d6adSToby Isaac 
PetscFEIntegrateResidual_Basic(PetscDS ds,PetscFormKey key,PetscInt Ne,PetscFEGeom * cgeom,const PetscScalar coefficients[],const PetscScalar coefficients_t[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscReal t,PetscScalar elemVec[])389d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
390d71ae5a4SJacob Faibussowitsch {
391b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
3926528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
3934bee2e38SMatthew G. Knepley   PetscFE            fe;
3946528b96dSMatthew G. Knepley   PetscWeakForm      wf;
3956528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
3962192575eSBarry Smith   PetscPointFn     **f0_func, **f1_func;
39720cf1dd8SToby Isaac   PetscQuadrature    quad;
398ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3994bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
40020cf1dd8SToby Isaac   const PetscScalar *constants;
40187510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
402ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
403ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
40420cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
4056587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
40620cf1dd8SToby Isaac 
40720cf1dd8SToby Isaac   PetscFunctionBegin;
4089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
4099566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
41087510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
4119566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
4129566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
4139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
4149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
4159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
4169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
4179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
4189566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
4193ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
4209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
4219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
4229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
4239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
42487510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
4259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
4264bee2e38SMatthew G. Knepley   if (dsAux) {
4279566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
4289566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
4299566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
4309566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
4319566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
4329566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
43363a3b9bcSJacob 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);
43420cf1dd8SToby Isaac   }
4359566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
43663a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
43720cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
43863a3b9bcSJacob Faibussowitsch   PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim);
43920cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4404bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
44120cf1dd8SToby Isaac 
4426587ee25SMatthew G. Knepley     fegeom.v = x; /* workspace */
4439566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc));
4449566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE));
44520cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4464bee2e38SMatthew G. Knepley       PetscReal w;
4474bee2e38SMatthew G. Knepley       PetscInt  c, d;
44820cf1dd8SToby Isaac 
4499566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom));
45087510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
4514bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
4526587ee25SMatthew G. Knepley       if (debug > 1 && q < cgeom->numPoints) {
45363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
4547be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
455ac9d17c7SMatthew G. Knepley         PetscCall(DMPrintCellMatrix(e, "invJ", dE, dE, fegeom.invJ));
45620cf1dd8SToby Isaac #endif
45720cf1dd8SToby Isaac       }
45816cd844bSPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
4599566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
460ac9d17c7SMatthew G. Knepley       for (i = 0; i < n0; ++i) f0_func[i](dE, 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]);
461ef0bb6c7SMatthew G. Knepley       for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w;
462ac9d17c7SMatthew G. Knepley       for (i = 0; i < n1; ++i) f1_func[i](dE, 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 * dE]);
4639371c9d4SSatish Balay       for (c = 0; c < T[field]->Nc; ++c)
464ac9d17c7SMatthew G. Knepley         for (d = 0; d < dE; ++d) f1[(q * T[field]->Nc + c) * dE + d] *= w;
465b8025e53SMatthew G. Knepley       if (debug) {
466e8e188d2SZach Atkins         // LCOV_EXCL_START
467e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q]));
468e8e188d2SZach Atkins         for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c]));
469e8e188d2SZach Atkins         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
470b8025e53SMatthew G. Knepley         if (debug > 2) {
47163a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field %" PetscInt_FMT ":", field));
47263a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c])));
4739566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
474e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  field der %" PetscInt_FMT ":", field));
475e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c])));
476e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
47763a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  resid %" PetscInt_FMT ":", field));
47863a3b9bcSJacob Faibussowitsch           for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c])));
4799566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
480e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  res der %" PetscInt_FMT ":", field));
481e8e188d2SZach Atkins           for (c = 0; c < T[field]->Nc; ++c) {
482ac9d17c7SMatthew G. Knepley             for (d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dE + d])));
483b8025e53SMatthew G. Knepley           }
484e8e188d2SZach Atkins           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
485e8e188d2SZach Atkins         }
486e8e188d2SZach Atkins         // LCOV_EXCL_STOP
487b8025e53SMatthew G. Knepley       }
48820cf1dd8SToby Isaac     }
4899566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset]));
49020cf1dd8SToby Isaac     cOffset += totDim;
49120cf1dd8SToby Isaac     cOffsetAux += totDimAux;
49220cf1dd8SToby Isaac   }
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49420cf1dd8SToby Isaac }
49520cf1dd8SToby Isaac 
PetscFEIntegrateBdResidual_Basic(PetscDS ds,PetscWeakForm wf,PetscFormKey key,PetscInt Ne,PetscFEGeom * fgeom,const PetscScalar coefficients[],const PetscScalar coefficients_t[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscReal t,PetscScalar elemVec[])496d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
497d71ae5a4SJacob Faibussowitsch {
498b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
49906d8a0d3SMatthew G. Knepley   const PetscInt     field = key.field;
5004bee2e38SMatthew G. Knepley   PetscFE            fe;
50106d8a0d3SMatthew G. Knepley   PetscInt           n0, n1, i;
5022192575eSBarry Smith   PetscBdPointFn   **f0_func, **f1_func;
50320cf1dd8SToby Isaac   PetscQuadrature    quad;
504ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
5054bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
50620cf1dd8SToby Isaac   const PetscScalar *constants;
50787510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
508ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
509ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
5106587ee25SMatthew G. Knepley   PetscBool          auxOnBd = PETSC_FALSE;
51120cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
5126587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
51320cf1dd8SToby Isaac 
51420cf1dd8SToby Isaac   PetscFunctionBegin;
5159566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
51787510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
5189566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &quad));
5199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
5209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
5219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
5229566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
5239566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset));
5249566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
5253ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
5269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
5279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
5289566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
5299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &Tf));
53087510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
5319566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
5324bee2e38SMatthew G. Knepley   if (dsAux) {
5339566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
5349566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
5359566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
5369566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
5379566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
5389566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
5397be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
5409566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
5419566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
54263a3b9bcSJacob 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);
54320cf1dd8SToby Isaac   }
544ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
5459566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
54663a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
54720cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
5486587ee25SMatthew G. Knepley   /* TODO FIX THIS */
5496587ee25SMatthew G. Knepley   fgeom->dim = dim - 1;
55063a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
55120cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5529f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
55320cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5549f209ee4SMatthew G. Knepley 
5556587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
5569566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcI));
5579566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcI * dE));
55820cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5594bee2e38SMatthew G. Knepley       PetscReal w;
5604bee2e38SMatthew G. Knepley       PetscInt  c, d;
5614bee2e38SMatthew G. Knepley 
5629566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom));
5639566063dSJacob Faibussowitsch       PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom));
56487510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
5654bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
56662bd480fSMatthew G. Knepley       if (debug > 1) {
567b6555650SPierre Jolivet         if ((fgeom->isAffine && q == 0) || !fgeom->isAffine) {
56863a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
5697be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5709566063dSJacob Faibussowitsch           PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
5719566063dSJacob Faibussowitsch           PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n));
57220cf1dd8SToby Isaac #endif
57320cf1dd8SToby Isaac         }
57462bd480fSMatthew G. Knepley       }
5758e3a54c0SPierre Jolivet       PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
5769566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
577ac9d17c7SMatthew G. Knepley       for (i = 0; i < n0; ++i) f0_func[i](dE, 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]);
5784bee2e38SMatthew G. Knepley       for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w;
579ac9d17c7SMatthew G. Knepley       for (i = 0; i < n1; ++i) f1_func[i](dE, 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 * dE]);
5809371c9d4SSatish Balay       for (c = 0; c < NcI; ++c)
581ac9d17c7SMatthew G. Knepley         for (d = 0; d < dE; ++d) f1[(q * NcI + c) * dE + d] *= w;
58262bd480fSMatthew G. Knepley       if (debug) {
58363a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q));
58462bd480fSMatthew G. Knepley         for (c = 0; c < NcI; ++c) {
58563a3b9bcSJacob Faibussowitsch           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c])));
58662bd480fSMatthew G. Knepley           if (n1) {
58763a3b9bcSJacob 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])));
5889566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
58962bd480fSMatthew G. Knepley           }
59062bd480fSMatthew G. Knepley         }
59162bd480fSMatthew G. Knepley       }
59220cf1dd8SToby Isaac     }
5939566063dSJacob Faibussowitsch     PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
59420cf1dd8SToby Isaac     cOffset += totDim;
59520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
59620cf1dd8SToby Isaac   }
5973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59820cf1dd8SToby Isaac }
59920cf1dd8SToby Isaac 
60027f02ce8SMatthew G. Knepley /*
60127f02ce8SMatthew 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
60227f02ce8SMatthew G. Knepley               all transforms operate in the full space and are square.
60327f02ce8SMatthew G. Knepley 
60427f02ce8SMatthew G. Knepley   HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square.
60527f02ce8SMatthew G. Knepley     1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces
60627f02ce8SMatthew G. Knepley     2) We need to assume that the orientation is 0 for both
60727f02ce8SMatthew 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()
60827f02ce8SMatthew G. Knepley */
PetscFEIntegrateHybridResidual_Basic(PetscDS ds,PetscDS dsIn,PetscFormKey key,PetscInt s,PetscInt Ne,PetscFEGeom * fgeom,PetscFEGeom * nbrgeom,const PetscScalar coefficients[],const PetscScalar coefficients_t[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscReal t,PetscScalar elemVec[])609989fa639SBrad Aagaard PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *nbrgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
610d71ae5a4SJacob Faibussowitsch {
611b2deab97SMatthew G. Knepley   const PetscInt     debug = ds->printIntegrate;
6126528b96dSMatthew G. Knepley   const PetscInt     field = key.field;
61327f02ce8SMatthew G. Knepley   PetscFE            fe;
6146528b96dSMatthew G. Knepley   PetscWeakForm      wf;
6156528b96dSMatthew G. Knepley   PetscInt           n0, n1, i;
6162192575eSBarry Smith   PetscBdPointFn   **f0_func, **f1_func;
61727f02ce8SMatthew G. Knepley   PetscQuadrature    quad;
6180e18dc48SMatthew G. Knepley   DMPolytopeType     ct;
61907218a29SMatthew G. Knepley   PetscTabulation   *Tf, *TfIn, *TfAux = NULL;
62027f02ce8SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
62127f02ce8SMatthew G. Knepley   const PetscScalar *constants;
62227f02ce8SMatthew G. Knepley   PetscReal         *x;
623665f567fSMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
62407218a29SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS;
6256587ee25SMatthew G. Knepley   PetscBool          isCohesiveField, auxOnBd = PETSC_FALSE;
62627f02ce8SMatthew G. Knepley   const PetscReal   *quadPoints, *quadWeights;
6276587ee25SMatthew G. Knepley   PetscInt           qdim, qNc, Nq, q, dE;
62827f02ce8SMatthew G. Knepley 
62927f02ce8SMatthew G. Knepley   PetscFunctionBegin;
63027f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
6319566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe));
6329566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
6339566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quad));
6349566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
6359566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
63607218a29SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn));
637429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets
63807218a29SMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x));
6399566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset));
6409566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
6419566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func));
6423ba16761SJacob Faibussowitsch   if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS);
6439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
6449566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL));
6459566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL));
64627f02ce8SMatthew G. Knepley   /* NOTE This is a bulk tabulation because the DS is a face discretization */
6479566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &Tf));
64807218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
64987510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE));
6509566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
65127f02ce8SMatthew G. Knepley   if (dsAux) {
6529566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
6539566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
6549566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
6559566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
6569566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
6579566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
65801907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
6599566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux));
6609566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux));
66163a3b9bcSJacob 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);
66227f02ce8SMatthew G. Knepley   }
6639566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField));
664665f567fSMatthew G. Knepley   NcI = Tf[field]->Nc;
665c2b7495fSMatthew G. Knepley   NcS = NcI;
6660abb75b6SMatthew G. Knepley   if (!isCohesiveField && s == 2) {
6670abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
6680abb75b6SMatthew G. Knepley     NcS *= 2;
6690abb75b6SMatthew G. Knepley   }
6709566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights));
6710e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
67263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
67327f02ce8SMatthew G. Knepley   dE = fgeom->dimEmbed;
67463a3b9bcSJacob Faibussowitsch   PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim);
67527f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
676989fa639SBrad Aagaard     // In order for the face information to be correct, the support of endcap faces _must_ be correctly oriented
677989fa639SBrad Aagaard     PetscFEGeom    fegeom, fegeomN[2];
6780e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
6790e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
6804e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
68127f02ce8SMatthew G. Knepley 
6826587ee25SMatthew G. Knepley     fegeom.v = x; /* Workspace */
6839566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f0, Nq * NcS));
6849566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(f1, Nq * NcS * dE));
6850847595cSMatthew G. Knepley     if (debug > 2) {
6860847595cSMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Negative %s face: %" PetscInt_FMT " (%" PetscInt_FMT ") (%" PetscInt_FMT ") perm %" PetscInt_FMT "\n", DMPolytopeTypes[ct], face[0], ornt[0], cornt[0], DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0])));
6870847595cSMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Positive %s face: %" PetscInt_FMT " (%" PetscInt_FMT ") (%" PetscInt_FMT ") perm %" PetscInt_FMT "\n", DMPolytopeTypes[ct], face[1], ornt[1], cornt[1], DMPolytopeTypeComposeOrientationInv(ct, cornt[1], ornt[1])));
6880847595cSMatthew G. Knepley     }
68927f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
6900e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
69127f02ce8SMatthew G. Knepley       PetscReal w;
69227f02ce8SMatthew G. Knepley       PetscInt  c, d;
69327f02ce8SMatthew G. Knepley 
6944e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0]));
6950847595cSMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[1], ornt[1]), field, q, &qpt[1]));
69607218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
697989fa639SBrad Aagaard       PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2, q, NULL, &fegeomN[0]));
698989fa639SBrad Aagaard       PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2 + 1, q, NULL, &fegeomN[1]));
69927f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
7006587ee25SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
70163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
70227f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
7039566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ));
70427f02ce8SMatthew G. Knepley #endif
70527f02ce8SMatthew G. Knepley       }
706a4158a15SMatthew 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]));
70727f02ce8SMatthew G. Knepley       /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */
708989fa639SBrad Aagaard       PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, fegeomN, &coefficients[cOffsetIn], PetscSafePointerPlusOffset(coefficients_t, cOffsetIn), u, u_x, u_t));
70907218a29SMatthew G. Knepley       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
710ac9d17c7SMatthew G. Knepley       for (i = 0; i < n0; ++i) f0_func[i](dE, 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]);
71127f02ce8SMatthew G. Knepley       for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w;
712ac9d17c7SMatthew G. Knepley       for (i = 0; i < n1; ++i) f1_func[i](dE, 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]);
7139371c9d4SSatish Balay       for (c = 0; c < NcS; ++c)
7149371c9d4SSatish Balay         for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w;
71588409cc3SMatthew G. Knepley       if (debug) {
71688409cc3SMatthew G. Knepley         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  elem %" PetscInt_FMT " quad point %" PetscInt_FMT " field %" PetscInt_FMT " side %" PetscInt_FMT "\n", e, q, field, s));
71788409cc3SMatthew G. Knepley         for (PetscInt f = 0; f < Nf; ++f) {
71888409cc3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Field %" PetscInt_FMT ":", f));
71988409cc3SMatthew G. Knepley           for (PetscInt c = uOff[f]; c < uOff[f + 1]; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  %g", (double)PetscRealPart(u[c])));
72088409cc3SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
72188409cc3SMatthew G. Knepley         }
72288409cc3SMatthew G. Knepley         for (c = 0; c < NcS; ++c) {
72388409cc3SMatthew G. Knepley           if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcS + c])));
72488409cc3SMatthew G. Knepley           if (n1) {
72588409cc3SMatthew G. Knepley             for (d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  f1[%" PetscInt_FMT ",%" PetscInt_FMT "] %g", c, d, (double)PetscRealPart(f1[(q * NcS + c) * dE + d])));
72688409cc3SMatthew G. Knepley             PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
72788409cc3SMatthew G. Knepley           }
72888409cc3SMatthew G. Knepley         }
72988409cc3SMatthew G. Knepley       }
73027f02ce8SMatthew G. Knepley     }
7319371c9d4SSatish Balay     if (isCohesiveField) {
7323ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
7339371c9d4SSatish Balay     } else {
7343ba16761SJacob Faibussowitsch       PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset]));
7359371c9d4SSatish Balay     }
73627f02ce8SMatthew G. Knepley     cOffset += totDim;
73707218a29SMatthew G. Knepley     cOffsetIn += totDimIn;
73827f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
73927f02ce8SMatthew G. Knepley   }
7403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74127f02ce8SMatthew G. Knepley }
74227f02ce8SMatthew G. Knepley 
PetscFEIntegrateJacobian_Basic(PetscDS rds,PetscDS cds,PetscFEJacobianType jtype,PetscFormKey key,PetscInt Ne,PetscFEGeom * cgeom,const PetscScalar coefficients[],const PetscScalar coefficients_t[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscReal t,PetscReal u_tshift,PetscScalar elemMat[])7434561e6c9SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS rds, PetscDS cds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
744d71ae5a4SJacob Faibussowitsch {
7454561e6c9SMatthew G. Knepley   const PetscInt     debug = rds->printIntegrate;
7464bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
7476528b96dSMatthew G. Knepley   PetscWeakForm      wf;
7482192575eSBarry Smith   PetscPointJacFn  **g0_func, **g1_func, **g2_func, **g3_func;
7494561e6c9SMatthew G. Knepley   PetscInt           n0, n1, n2, n3;
75020cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
75120cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
75220cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
75320cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
75420cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
75520cf1dd8SToby Isaac   PetscQuadrature    quad;
7564561e6c9SMatthew G. Knepley   PetscTabulation   *rT, *cT, *TAux = NULL;
7571a768569SStefano Zampini   PetscScalar       *g0 = NULL, *g1 = NULL, *g2 = NULL, *g3 = NULL, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
75820cf1dd8SToby Isaac   const PetscScalar *constants;
75987510d7dSMatthew G. Knepley   PetscReal         *x, cellScale;
760ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
761ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
7624561e6c9SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, rtotDim, ctotDim, totDimAux = 0;
76320cf1dd8SToby Isaac   PetscInt           dE, Np;
76420cf1dd8SToby Isaac   PetscBool          isAffine;
76520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
7664561e6c9SMatthew G. Knepley   PetscInt           qNc, Nq;
76720cf1dd8SToby Isaac 
76820cf1dd8SToby Isaac   PetscFunctionBegin;
7694561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetNumFields(rds, &Nf));
7706528b96dSMatthew G. Knepley   fieldI = key.field / Nf;
7716528b96dSMatthew G. Knepley   fieldJ = key.field % Nf;
7724561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(rds, fieldI, (PetscObject *)&feI));
7734561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetDiscretization(cds, fieldJ, (PetscObject *)&feJ));
7749566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
77587510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
7769566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
7774561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(rds, &rtotDim));
7784561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetTotalDimension(cds, &ctotDim));
7794561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsets(rds, &uOff));
7804561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetComponentDerivativeOffsets(rds, &uOff_x));
7814561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetWeakForm(rds, &wf));
78220cf1dd8SToby Isaac   switch (jtype) {
783d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
784d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
785d71ae5a4SJacob Faibussowitsch     break;
786d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
787d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
788d71ae5a4SJacob Faibussowitsch     break;
789d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
790d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
791d71ae5a4SJacob Faibussowitsch     break;
79220cf1dd8SToby Isaac   }
7933ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
7944561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetEvaluationArrays(rds, &u, coefficients_t ? &u_t : NULL, &u_x));
7954561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetWorkspace(rds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
7964561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetWeakFormArrays(rds, NULL, NULL, n0 ? &g0 : NULL, n1 ? &g1 : NULL, n2 ? &g2 : NULL, n3 ? &g3 : NULL));
7971a768569SStefano Zampini 
7984561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetTabulation(rds, &rT));
7994561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetTabulation(cds, &cT));
8004561e6c9SMatthew G. Knepley   PetscCheck(rT[0]->Np == cT[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of row tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of col tabulation points", rT[0]->Np, cT[0]->Np);
8014561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetFieldOffset(rds, fieldI, &offsetI));
8024561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetFieldOffset(cds, fieldJ, &offsetJ));
8034561e6c9SMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(rds, fieldI, fieldJ));
8044561e6c9SMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(cds, fieldI, fieldJ));
8054561e6c9SMatthew G. Knepley   PetscCall(PetscDSGetConstants(rds, &numConstants, &constants));
8064bee2e38SMatthew G. Knepley   if (dsAux) {
8079566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
8089566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
8099566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
8109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
8119566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
8129566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTabulation(dsAux, &TAux));
8134561e6c9SMatthew G. Knepley     PetscCheck(rT[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", rT[0]->Np, TAux[0]->Np);
81420cf1dd8SToby Isaac   }
8154561e6c9SMatthew G. Knepley   NcI      = rT[fieldI]->Nc;
8164561e6c9SMatthew G. Knepley   NcJ      = cT[fieldJ]->Nc;
8174bee2e38SMatthew G. Knepley   Np       = cgeom->numPoints;
8184bee2e38SMatthew G. Knepley   dE       = cgeom->dimEmbed;
8194bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
8209566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
82163a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
8221a768569SStefano Zampini 
8234561e6c9SMatthew G. Knepley   for (PetscInt e = 0; e < Ne; ++e) {
8244bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
8254bee2e38SMatthew G. Knepley 
82627f02ce8SMatthew G. Knepley     fegeom.dim      = cgeom->dim;
82727f02ce8SMatthew G. Knepley     fegeom.dimEmbed = cgeom->dimEmbed;
828a7be0fb8SMatthew G. Knepley     fegeom.xi       = NULL;
8294bee2e38SMatthew G. Knepley     if (isAffine) {
8304bee2e38SMatthew G. Knepley       fegeom.v    = x;
8314bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
8327132c3f7SMatthew G. Knepley       fegeom.J    = &cgeom->J[e * Np * dE * dE];
8337132c3f7SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e * Np * dE * dE];
8347132c3f7SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e * Np];
83525c09e31SMark Adams     } else fegeom.xi = NULL;
8364561e6c9SMatthew G. Knepley     for (PetscInt q = 0; q < Nq; ++q) {
83720cf1dd8SToby Isaac       PetscReal w;
83820cf1dd8SToby Isaac 
83920cf1dd8SToby Isaac       if (isAffine) {
8407132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x);
84120cf1dd8SToby Isaac       } else {
8424bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e * Np + q) * dE];
8434bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e * Np + q) * dE * dE];
8444bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE];
8454bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e * Np + q];
84620cf1dd8SToby Isaac       }
8474561e6c9SMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(rds, fegeom.detJ[0] * cellScale));
8489566063dSJacob 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]));
8494bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
8504561e6c9SMatthew G. Knepley       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(rds, Nf, 0, q, rT, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
8519566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
852ea672e62SMatthew G. Knepley       if (n0) {
8539566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
8544561e6c9SMatthew G. Knepley         for (PetscInt i = 0; i < n0; ++i) g0_func[i](dE, 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);
8554561e6c9SMatthew G. Knepley         for (PetscInt c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
85620cf1dd8SToby Isaac       }
857ea672e62SMatthew G. Knepley       if (n1) {
8589566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
8594561e6c9SMatthew G. Knepley         for (PetscInt i = 0; i < n1; ++i) g1_func[i](dE, 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);
8604561e6c9SMatthew G. Knepley         for (PetscInt c = 0; c < NcI * NcJ * dE; ++c) g1[c] *= w;
86120cf1dd8SToby Isaac       }
862ea672e62SMatthew G. Knepley       if (n2) {
8639566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
8644561e6c9SMatthew G. Knepley         for (PetscInt i = 0; i < n2; ++i) g2_func[i](dE, 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);
8654561e6c9SMatthew G. Knepley         for (PetscInt c = 0; c < NcI * NcJ * dE; ++c) g2[c] *= w;
86620cf1dd8SToby Isaac       }
867ea672e62SMatthew G. Knepley       if (n3) {
8689566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
8694561e6c9SMatthew G. Knepley         for (PetscInt i = 0; i < n3; ++i) g3_func[i](dE, 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);
8704561e6c9SMatthew G. Knepley         for (PetscInt c = 0; c < NcI * NcJ * dE * dE; ++c) g3[c] *= w;
87120cf1dd8SToby Isaac       }
87220cf1dd8SToby Isaac 
8734561e6c9SMatthew G. Knepley       PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, rT[fieldI], basisReal, basisDerReal, cT[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, ctotDim, offsetI, offsetJ, elemMat + eOffset));
87420cf1dd8SToby Isaac     }
87520cf1dd8SToby Isaac     if (debug > 1) {
87663a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
8774561e6c9SMatthew G. Knepley       for (PetscInt f = 0; f < rT[fieldI]->Nb; ++f) {
878699b48e5SStefano Zampini         const PetscInt i = offsetI + f;
8794561e6c9SMatthew G. Knepley         for (PetscInt g = 0; g < cT[fieldJ]->Nb; ++g) {
880699b48e5SStefano Zampini           const PetscInt j = offsetJ + g;
8814561e6c9SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "    elemMat[%" PetscInt_FMT ", %" PetscInt_FMT "]: %g\n", f, g, (double)PetscRealPart(elemMat[eOffset + i * ctotDim + j])));
88220cf1dd8SToby Isaac         }
8839566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
88420cf1dd8SToby Isaac       }
88520cf1dd8SToby Isaac     }
8864561e6c9SMatthew G. Knepley     cOffset += rtotDim;
88720cf1dd8SToby Isaac     cOffsetAux += totDimAux;
8884561e6c9SMatthew G. Knepley     eOffset += rtotDim * ctotDim;
88920cf1dd8SToby Isaac   }
8903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89120cf1dd8SToby Isaac }
89220cf1dd8SToby Isaac 
PetscFEIntegrateBdJacobian_Basic(PetscDS ds,PetscWeakForm wf,PetscFEJacobianType jtype,PetscFormKey key,PetscInt Ne,PetscFEGeom * fgeom,const PetscScalar coefficients[],const PetscScalar coefficients_t[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscReal t,PetscReal u_tshift,PetscScalar elemMat[])893e3d591f2SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
894d71ae5a4SJacob Faibussowitsch {
895b2deab97SMatthew G. Knepley   const PetscInt      debug = ds->printIntegrate;
8964bee2e38SMatthew G. Knepley   PetscFE             feI, feJ;
8972192575eSBarry Smith   PetscBdPointJacFn **g0_func, **g1_func, **g2_func, **g3_func;
89845480ffeSMatthew G. Knepley   PetscInt            n0, n1, n2, n3, i;
89920cf1dd8SToby Isaac   PetscInt            cOffset    = 0; /* Offset into coefficients[] for element e */
90020cf1dd8SToby Isaac   PetscInt            cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
90120cf1dd8SToby Isaac   PetscInt            eOffset    = 0; /* Offset into elemMat[] for element e */
90220cf1dd8SToby Isaac   PetscInt            offsetI    = 0; /* Offset into an element vector for fieldI */
90320cf1dd8SToby Isaac   PetscInt            offsetJ    = 0; /* Offset into an element vector for fieldJ */
90420cf1dd8SToby Isaac   PetscQuadrature     quad;
905ef0bb6c7SMatthew G. Knepley   PetscTabulation    *T, *TAux = NULL;
9064bee2e38SMatthew G. Knepley   PetscScalar        *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
90720cf1dd8SToby Isaac   const PetscScalar  *constants;
90887510d7dSMatthew G. Knepley   PetscReal          *x, cellScale;
909ef0bb6c7SMatthew G. Knepley   PetscInt           *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
910ef0bb6c7SMatthew G. Knepley   PetscInt            NcI = 0, NcJ = 0;
91145480ffeSMatthew G. Knepley   PetscInt            dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
91220cf1dd8SToby Isaac   PetscBool           isAffine;
91320cf1dd8SToby Isaac   const PetscReal    *quadPoints, *quadWeights;
91420cf1dd8SToby Isaac   PetscInt            qNc, Nq, q, Np, dE;
91520cf1dd8SToby Isaac 
91620cf1dd8SToby Isaac   PetscFunctionBegin;
9179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
91845480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
91945480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
9209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
9219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
9229566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
92387510d7dSMatthew G. Knepley   cellScale = (PetscReal)PetscPowInt(2, dim);
9249566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(feI, &quad));
9259566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
9269566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentOffsets(ds, &uOff));
9279566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x));
9289566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI));
9299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ));
930e3d591f2SMatthew G. Knepley   switch (jtype) {
931e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN_PRE:
932e3d591f2SMatthew G. Knepley     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
933e3d591f2SMatthew G. Knepley     break;
934e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN:
9359566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
936e3d591f2SMatthew G. Knepley     break;
937e3d591f2SMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN:
938e3d591f2SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()");
939e3d591f2SMatthew G. Knepley   }
9403ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
9419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
9429566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
9439566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
9449566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFaceTabulation(ds, &T));
94587510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
9469566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
9474bee2e38SMatthew G. Knepley   if (dsAux) {
9489566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
9499566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
9509566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
9519566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
9529566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
9539566063dSJacob Faibussowitsch     PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
95420cf1dd8SToby Isaac   }
955ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
95620cf1dd8SToby Isaac   Np       = fgeom->numPoints;
95720cf1dd8SToby Isaac   dE       = fgeom->dimEmbed;
95820cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
95927f02ce8SMatthew G. Knepley   /* Initialize here in case the function is not defined */
9609566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcI * NcJ));
9619566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
9629566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
9639566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
9649566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
96563a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
96620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
9679f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
96820cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
969ea78f98cSLisandro Dalcin     fegeom.n            = NULL;
970ea78f98cSLisandro Dalcin     fegeom.v            = NULL;
971a7be0fb8SMatthew G. Knepley     fegeom.xi           = NULL;
972ea78f98cSLisandro Dalcin     fegeom.J            = NULL;
973ea78f98cSLisandro Dalcin     fegeom.detJ         = NULL;
97427f02ce8SMatthew G. Knepley     fegeom.dim          = fgeom->dim;
97527f02ce8SMatthew G. Knepley     fegeom.dimEmbed     = fgeom->dimEmbed;
97627f02ce8SMatthew G. Knepley     cgeom.dim           = fgeom->dim;
97727f02ce8SMatthew G. Knepley     cgeom.dimEmbed      = fgeom->dimEmbed;
9784bee2e38SMatthew G. Knepley     if (isAffine) {
9794bee2e38SMatthew G. Knepley       fegeom.v    = x;
9804bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
9817132c3f7SMatthew G. Knepley       fegeom.J    = &fgeom->J[e * Np * dE * dE];
9827132c3f7SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e * Np * dE * dE];
9837132c3f7SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e * Np];
9847132c3f7SMatthew G. Knepley       fegeom.n    = &fgeom->n[e * Np * dE];
9859f209ee4SMatthew G. Knepley 
9867132c3f7SMatthew G. Knepley       cgeom.J    = &fgeom->suppJ[0][e * Np * dE * dE];
9877132c3f7SMatthew G. Knepley       cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE];
9887132c3f7SMatthew G. Knepley       cgeom.detJ = &fgeom->suppDetJ[0][e * Np];
98925c09e31SMark Adams     } else fegeom.xi = NULL;
99020cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
99120cf1dd8SToby Isaac       PetscReal w;
9924bee2e38SMatthew G. Knepley       PetscInt  c;
99320cf1dd8SToby Isaac 
99463a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
99520cf1dd8SToby Isaac       if (isAffine) {
9967132c3f7SMatthew G. Knepley         CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x);
99720cf1dd8SToby Isaac       } else {
9983fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e * Np + q) * dE];
9999f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e * Np + q) * dE * dE];
10009f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE];
10014bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e * Np + q];
10024bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e * Np + q) * dE];
10039f209ee4SMatthew G. Knepley 
10049f209ee4SMatthew G. Knepley         cgeom.J    = &fgeom->suppJ[0][(e * Np + q) * dE * dE];
10059f209ee4SMatthew G. Knepley         cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE];
10069f209ee4SMatthew G. Knepley         cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q];
100720cf1dd8SToby Isaac       }
100887510d7dSMatthew G. Knepley       PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale));
10094bee2e38SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
10109566063dSJacob Faibussowitsch       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t));
10119566063dSJacob Faibussowitsch       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1012ea672e62SMatthew G. Knepley       if (n0) {
10139566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcI * NcJ));
1014ac9d17c7SMatthew G. Knepley         for (i = 0; i < n0; ++i) g0_func[i](dE, 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);
101520cf1dd8SToby Isaac         for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w;
101620cf1dd8SToby Isaac       }
1017ea672e62SMatthew G. Knepley       if (n1) {
10189566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g1, NcI * NcJ * dE));
1019ac9d17c7SMatthew G. Knepley         for (i = 0; i < n1; ++i) g1_func[i](dE, 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);
10204bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w;
102120cf1dd8SToby Isaac       }
1022ea672e62SMatthew G. Knepley       if (n2) {
10239566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g2, NcI * NcJ * dE));
1024ac9d17c7SMatthew G. Knepley         for (i = 0; i < n2; ++i) g2_func[i](dE, 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);
10254bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w;
102620cf1dd8SToby Isaac       }
1027ea672e62SMatthew G. Knepley       if (n3) {
10289566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE));
1029ac9d17c7SMatthew G. Knepley         for (i = 0; i < n3; ++i) g3_func[i](dE, 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);
10304bee2e38SMatthew G. Knepley         for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w;
103120cf1dd8SToby Isaac       }
103220cf1dd8SToby Isaac 
10331a768569SStefano Zampini       PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset));
103420cf1dd8SToby Isaac     }
103520cf1dd8SToby Isaac     if (debug > 1) {
103620cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
103720cf1dd8SToby Isaac 
103863a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ));
1039ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
1040ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
1041ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f * T[fieldI]->Nc + fc;
1042ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
1043ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
1044ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc;
104563a3b9bcSJacob 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])));
104620cf1dd8SToby Isaac             }
104720cf1dd8SToby Isaac           }
10489566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
104920cf1dd8SToby Isaac         }
105020cf1dd8SToby Isaac       }
105120cf1dd8SToby Isaac     }
105220cf1dd8SToby Isaac     cOffset += totDim;
105320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
105420cf1dd8SToby Isaac     eOffset += PetscSqr(totDim);
105520cf1dd8SToby Isaac   }
10563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
105720cf1dd8SToby Isaac }
105820cf1dd8SToby Isaac 
PetscFEIntegrateHybridJacobian_Basic(PetscDS ds,PetscDS dsIn,PetscFEJacobianType jtype,PetscFormKey key,PetscInt s,PetscInt Ne,PetscFEGeom * fgeom,PetscFEGeom * nbrgeom,const PetscScalar coefficients[],const PetscScalar coefficients_t[],PetscDS dsAux,const PetscScalar coefficientsAux[],PetscReal t,PetscReal u_tshift,PetscScalar elemMat[])1059989fa639SBrad Aagaard PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *nbrgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
1060d71ae5a4SJacob Faibussowitsch {
1061b2deab97SMatthew G. Knepley   const PetscInt      debug = ds->printIntegrate;
106227f02ce8SMatthew G. Knepley   PetscFE             feI, feJ;
1063148442b3SMatthew G. Knepley   PetscWeakForm       wf;
10642192575eSBarry Smith   PetscBdPointJacFn **g0_func, **g1_func, **g2_func, **g3_func;
1065148442b3SMatthew G. Knepley   PetscInt            n0, n1, n2, n3, i;
106627f02ce8SMatthew G. Knepley   PetscInt            cOffset    = 0; /* Offset into coefficients[] for element e */
106727f02ce8SMatthew G. Knepley   PetscInt            cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
106827f02ce8SMatthew G. Knepley   PetscInt            eOffset    = 0; /* Offset into elemMat[] for element e */
106927f02ce8SMatthew G. Knepley   PetscInt            offsetI    = 0; /* Offset into an element vector for fieldI */
107027f02ce8SMatthew G. Knepley   PetscInt            offsetJ    = 0; /* Offset into an element vector for fieldJ */
1071665f567fSMatthew G. Knepley   PetscQuadrature     quad;
10720e18dc48SMatthew G. Knepley   DMPolytopeType      ct;
107307218a29SMatthew G. Knepley   PetscTabulation    *T, *TfIn, *TAux = NULL;
107427f02ce8SMatthew G. Knepley   PetscScalar        *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
107527f02ce8SMatthew G. Knepley   const PetscScalar  *constants;
107627f02ce8SMatthew G. Knepley   PetscReal          *x;
1077665f567fSMatthew G. Knepley   PetscInt           *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
1078665f567fSMatthew G. Knepley   PetscInt            NcI = 0, NcJ = 0, NcS, NcT;
107945480ffeSMatthew G. Knepley   PetscInt            dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e;
108007218a29SMatthew G. Knepley   PetscBool           isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE;
108127f02ce8SMatthew G. Knepley   const PetscReal    *quadPoints, *quadWeights;
10820563a546SMatthew G. Knepley   PetscInt            qNc, Nq, q, dE;
108327f02ce8SMatthew G. Knepley 
108427f02ce8SMatthew G. Knepley   PetscFunctionBegin;
10859566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
108645480ffeSMatthew G. Knepley   fieldI = key.field / Nf;
108745480ffeSMatthew G. Knepley   fieldJ = key.field % Nf;
108827f02ce8SMatthew G. Knepley   /* Hybrid discretization is posed directly on faces */
10899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI));
10909566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ));
10919566063dSJacob Faibussowitsch   PetscCall(PetscFEGetSpatialDimension(feI, &dim));
10929566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(feI, &quad));
10939566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(ds, &totDim));
1094429ebbe4SMatthew G. Knepley   PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets
10959566063dSJacob Faibussowitsch   PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x));
10969566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
109727f02ce8SMatthew G. Knepley   switch (jtype) {
1098d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_PRE:
1099d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1100d71ae5a4SJacob Faibussowitsch     break;
1101d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN:
1102d71ae5a4SJacob Faibussowitsch     PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func));
1103d71ae5a4SJacob Faibussowitsch     break;
1104d71ae5a4SJacob Faibussowitsch   case PETSCFE_JACOBIAN_DYN:
1105d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)");
110627f02ce8SMatthew G. Knepley   }
11073ba16761SJacob Faibussowitsch   if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS);
11089566063dSJacob Faibussowitsch   PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x));
11099566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal));
11109566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3));
11119566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTabulation(ds, &T));
111207218a29SMatthew G. Knepley   PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn));
11139566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI));
11149566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ));
111587510d7dSMatthew G. Knepley   PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ));
11169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(ds, &numConstants, &constants));
111727f02ce8SMatthew G. Knepley   if (dsAux) {
11189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux));
11199566063dSJacob Faibussowitsch     PetscCall(PetscDSGetNumFields(dsAux, &NfAux));
11209566063dSJacob Faibussowitsch     PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux));
11219566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff));
11229566063dSJacob Faibussowitsch     PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x));
11239566063dSJacob Faibussowitsch     PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x));
112401907d53SMatthew G. Knepley     auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE;
11259566063dSJacob Faibussowitsch     if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux));
11269566063dSJacob Faibussowitsch     else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux));
112763a3b9bcSJacob 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);
112827f02ce8SMatthew G. Knepley   }
11299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI));
11309566063dSJacob Faibussowitsch   PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ));
11310563a546SMatthew G. Knepley   dE  = fgeom->dimEmbed;
1132665f567fSMatthew G. Knepley   NcI = T[fieldI]->Nc;
1133665f567fSMatthew G. Knepley   NcJ = T[fieldJ]->Nc;
113427f02ce8SMatthew G. Knepley   NcS = isCohesiveFieldI ? NcI : 2 * NcI;
113527f02ce8SMatthew G. Knepley   NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ;
11360abb75b6SMatthew G. Knepley   if (!isCohesiveFieldI && s == 2) {
11370abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
11380abb75b6SMatthew G. Knepley     NcS *= 2;
11390abb75b6SMatthew G. Knepley   }
11400abb75b6SMatthew G. Knepley   if (!isCohesiveFieldJ && s == 2) {
11410abb75b6SMatthew G. Knepley     // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides
11420abb75b6SMatthew G. Knepley     NcT *= 2;
11430abb75b6SMatthew G. Knepley   }
11440502970dSMatthew G. Knepley   // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though
11450502970dSMatthew G. Knepley   // the coordinates are in dE dimensions
11469566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(g0, NcS * NcT));
11470502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
11480502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
11490502970dSMatthew G. Knepley   PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
11509566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
11510e18dc48SMatthew G. Knepley   PetscCall(PetscQuadratureGetCellType(quad, &ct));
115263a3b9bcSJacob Faibussowitsch   PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc);
115327f02ce8SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
1154989fa639SBrad Aagaard     PetscFEGeom    fegeom, fegeomN[2];
11550e18dc48SMatthew G. Knepley     const PetscInt face[2]  = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]};
11560e18dc48SMatthew G. Knepley     const PetscInt ornt[2]  = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]};
11574e913f38SMatthew G. Knepley     const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]};
115827f02ce8SMatthew G. Knepley 
115907218a29SMatthew G. Knepley     fegeom.v = x; /* Workspace */
116027f02ce8SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
11610e18dc48SMatthew G. Knepley       PetscInt  qpt[2];
116227f02ce8SMatthew G. Knepley       PetscReal w;
116327f02ce8SMatthew G. Knepley       PetscInt  c;
116427f02ce8SMatthew G. Knepley 
11654e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0]));
11664e913f38SMatthew G. Knepley       PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1]));
116707218a29SMatthew G. Knepley       PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom));
1168989fa639SBrad Aagaard       PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2, q, NULL, &fegeomN[0]));
1169989fa639SBrad Aagaard       PetscCall(PetscFEGeomGetPoint(nbrgeom, e * 2 + 1, q, NULL, &fegeomN[1]));
117027f02ce8SMatthew G. Knepley       w = fegeom.detJ[0] * quadWeights[q];
117107218a29SMatthew G. Knepley       if (debug > 1 && q < fgeom->numPoints) {
117263a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", (double)fegeom.detJ[0]));
117327f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX)
11749566063dSJacob Faibussowitsch         PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ));
117527f02ce8SMatthew G. Knepley #endif
117627f02ce8SMatthew G. Knepley       }
117763a3b9bcSJacob Faibussowitsch       if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, "  quad point %" PetscInt_FMT "\n", q));
1178989fa639SBrad Aagaard       if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(dsIn, Nf, 0, q, T, face, qpt, TfIn, &fegeom, fegeomN, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t));
117907218a29SMatthew G. Knepley       if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL));
1180ea672e62SMatthew G. Knepley       if (n0) {
11819566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(g0, NcS * NcT));
11820563a546SMatthew G. Knepley         for (i = 0; i < n0; ++i) g0_func[i](dE, 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);
118327f02ce8SMatthew G. Knepley         for (c = 0; c < NcS * NcT; ++c) g0[c] *= w;
118427f02ce8SMatthew G. Knepley       }
1185ea672e62SMatthew G. Knepley       if (n1) {
11860502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g1, NcS * NcT * dim));
11870563a546SMatthew G. Knepley         for (i = 0; i < n1; ++i) g1_func[i](dE, 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);
11880502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w;
118927f02ce8SMatthew G. Knepley       }
1190ea672e62SMatthew G. Knepley       if (n2) {
11910502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g2, NcS * NcT * dim));
11920563a546SMatthew G. Knepley         for (i = 0; i < n2; ++i) g2_func[i](dE, 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);
11930502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w;
119427f02ce8SMatthew G. Knepley       }
1195ea672e62SMatthew G. Knepley       if (n3) {
11960502970dSMatthew G. Knepley         PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim));
11970563a546SMatthew G. Knepley         for (i = 0; i < n3; ++i) g3_func[i](dE, 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);
11980502970dSMatthew G. Knepley         for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w;
119927f02ce8SMatthew G. Knepley       }
120027f02ce8SMatthew G. Knepley 
12015fedec97SMatthew G. Knepley       if (isCohesiveFieldI) {
12025fedec97SMatthew G. Knepley         if (isCohesiveFieldJ) {
1203989fa639SBrad Aagaard           //PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset));
1204989fa639SBrad Aagaard           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
120527f02ce8SMatthew G. Knepley         } else {
12060abb75b6SMatthew G. Knepley           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
12070abb75b6SMatthew G. Knepley           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
12080abb75b6SMatthew G. Knepley         }
12090abb75b6SMatthew G. Knepley       } else {
12100abb75b6SMatthew G. Knepley         if (s == 2) {
12110abb75b6SMatthew G. Knepley           if (isCohesiveFieldJ) {
12120abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
12130abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
12140abb75b6SMatthew G. Knepley           } else {
12150abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
12160abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat));
12170abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 2], &g1[NcI * NcJ * dim * 2], &g2[NcI * NcJ * dim * 2], &g3[NcI * NcJ * dim * dim * 2], eOffset, totDim, offsetI, offsetJ, elemMat));
12180abb75b6SMatthew G. Knepley             PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 3], &g1[NcI * NcJ * dim * 3], &g2[NcI * NcJ * dim * 3], &g3[NcI * NcJ * dim * dim * 3], eOffset, totDim, offsetI, offsetJ, elemMat));
12195fedec97SMatthew G. Knepley           }
12209371c9d4SSatish Balay         } else
12210abb75b6SMatthew G. Knepley           PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat));
12220abb75b6SMatthew G. Knepley       }
122327f02ce8SMatthew G. Knepley     }
122427f02ce8SMatthew G. Knepley     if (debug > 1) {
12254e913f38SMatthew G. Knepley       const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb));
12264e913f38SMatthew G. Knepley       const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb));
12274e913f38SMatthew G. Knepley       const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb));
12284e913f38SMatthew G. Knepley       const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb));
12294e913f38SMatthew G. Knepley       PetscInt       f, g;
123027f02ce8SMatthew G. Knepley 
12314e913f38SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT " s %s totDim %" PetscInt_FMT " offsets (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", fieldI, fieldJ, s ? (s > 1 ? "Coh" : "Pos") : "Neg", totDim, eOffset, offsetI, offsetJ));
12324e913f38SMatthew G. Knepley       for (f = fS; f < fE; ++f) {
12334e913f38SMatthew G. Knepley         const PetscInt i = offsetI + f;
12344e913f38SMatthew G. Knepley         for (g = gS; g < gE; ++g) {
12354e913f38SMatthew G. Knepley           const PetscInt j = offsetJ + g;
1236e978a55eSPierre Jolivet           PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, f, i, g, j);
12374e913f38SMatthew G. Knepley           PetscCall(PetscPrintf(PETSC_COMM_SELF, "    elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f / NcI, f % NcI, g / NcJ, g % NcJ, (double)PetscRealPart(elemMat[eOffset + i * totDim + j])));
123827f02ce8SMatthew G. Knepley         }
12399566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
124027f02ce8SMatthew G. Knepley       }
124127f02ce8SMatthew G. Knepley     }
124227f02ce8SMatthew G. Knepley     cOffset += totDim;
124327f02ce8SMatthew G. Knepley     cOffsetAux += totDimAux;
124427f02ce8SMatthew G. Knepley     eOffset += PetscSqr(totDim);
124527f02ce8SMatthew G. Knepley   }
12463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124727f02ce8SMatthew G. Knepley }
124827f02ce8SMatthew G. Knepley 
PetscFEInitialize_Basic(PetscFE fem)1249d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
1250d71ae5a4SJacob Faibussowitsch {
125120cf1dd8SToby Isaac   PetscFunctionBegin;
125220cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
125320cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
125420cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
125520cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
125620cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
1257ce78bad3SBarry Smith   fem->ops->computetabulation       = PetscFEComputeTabulation_Basic;
125820cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
1259afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
126020cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
126120cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
126227f02ce8SMatthew G. Knepley   fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
126320cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */;
126420cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
126520cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
126627f02ce8SMatthew G. Knepley   fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
12673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126820cf1dd8SToby Isaac }
126920cf1dd8SToby Isaac 
127020cf1dd8SToby Isaac /*MC
1271dce8aebaSBarry Smith   PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization
127220cf1dd8SToby Isaac 
127320cf1dd8SToby Isaac   Level: intermediate
127420cf1dd8SToby Isaac 
1275dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
127620cf1dd8SToby Isaac M*/
127720cf1dd8SToby Isaac 
PetscFECreate_Basic(PetscFE fem)1278d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
1279d71ae5a4SJacob Faibussowitsch {
128020cf1dd8SToby Isaac   PetscFE_Basic *b;
128120cf1dd8SToby Isaac 
128220cf1dd8SToby Isaac   PetscFunctionBegin;
128320cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
12844dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
128520cf1dd8SToby Isaac   fem->data = b;
128620cf1dd8SToby Isaac 
12879566063dSJacob Faibussowitsch   PetscCall(PetscFEInitialize_Basic(fem));
12883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
128920cf1dd8SToby Isaac }
1290