xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 7be5e7486ad090eef6dbc1e09e6fd40be450d1fc)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
420cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
520cf1dd8SToby Isaac {
620cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *) fem->data;
720cf1dd8SToby Isaac   PetscErrorCode ierr;
820cf1dd8SToby Isaac 
920cf1dd8SToby Isaac   PetscFunctionBegin;
1020cf1dd8SToby Isaac   ierr = PetscFree(b);CHKERRQ(ierr);
1120cf1dd8SToby Isaac   PetscFunctionReturn(0);
1220cf1dd8SToby Isaac }
1320cf1dd8SToby Isaac 
1420cf1dd8SToby Isaac PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer viewer)
1520cf1dd8SToby Isaac {
1620cf1dd8SToby Isaac   PetscSpace        basis;
1720cf1dd8SToby Isaac   PetscDualSpace    dual;
1820cf1dd8SToby Isaac   PetscQuadrature   q = NULL;
1920cf1dd8SToby Isaac   PetscInt          dim, Nc, Nq;
2020cf1dd8SToby Isaac   PetscViewerFormat format;
2120cf1dd8SToby Isaac   PetscErrorCode    ierr;
2220cf1dd8SToby Isaac 
2320cf1dd8SToby Isaac   PetscFunctionBegin;
2420cf1dd8SToby Isaac   ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr);
2520cf1dd8SToby Isaac   ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr);
2620cf1dd8SToby Isaac   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
2720cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
28*7be5e748SToby Isaac   if (q) {ierr = PetscQuadratureGetData(q, &dim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
2920cf1dd8SToby Isaac   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
3020cf1dd8SToby Isaac   ierr = PetscViewerASCIIPrintf(viewer, "Basic Finite Element:\n");CHKERRQ(ierr);
3120cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "  dimension:       %d\n", dim);CHKERRQ(ierr);
3220cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "  components:      %d\n", Nc);CHKERRQ(ierr);
3320cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "  num quad points: %d\n", Nq);CHKERRQ(ierr);
3420cf1dd8SToby Isaac   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
3520cf1dd8SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
36*7be5e748SToby Isaac     if (q) {ierr = PetscQuadratureView(q, viewer);CHKERRQ(ierr);}
3720cf1dd8SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3820cf1dd8SToby Isaac   }
3920cf1dd8SToby Isaac   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
4020cf1dd8SToby Isaac   ierr = PetscSpaceView(basis, viewer);CHKERRQ(ierr);
4120cf1dd8SToby Isaac   ierr = PetscDualSpaceView(dual, viewer);CHKERRQ(ierr);
4220cf1dd8SToby Isaac   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
4320cf1dd8SToby Isaac   PetscFunctionReturn(0);
4420cf1dd8SToby Isaac }
4520cf1dd8SToby Isaac 
4620cf1dd8SToby Isaac PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer viewer)
4720cf1dd8SToby Isaac {
4820cf1dd8SToby Isaac   PetscBool      iascii;
4920cf1dd8SToby Isaac   PetscErrorCode ierr;
5020cf1dd8SToby Isaac 
5120cf1dd8SToby Isaac   PetscFunctionBegin;
5220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
5320cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
5420cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
5520cf1dd8SToby Isaac   if (iascii) {ierr = PetscFEView_Basic_Ascii(fe, viewer);CHKERRQ(ierr);}
5620cf1dd8SToby Isaac   PetscFunctionReturn(0);
5720cf1dd8SToby Isaac }
5820cf1dd8SToby Isaac 
5920cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
6020cf1dd8SToby Isaac PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
6120cf1dd8SToby Isaac {
6220cf1dd8SToby Isaac   PetscScalar   *work, *invVscalar;
6320cf1dd8SToby Isaac   PetscBLASInt  *pivots;
6420cf1dd8SToby Isaac   PetscBLASInt   n, info;
6520cf1dd8SToby Isaac   PetscInt       pdim, j;
6620cf1dd8SToby Isaac   PetscErrorCode ierr;
6720cf1dd8SToby Isaac 
6820cf1dd8SToby Isaac   PetscFunctionBegin;
6920cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
7020cf1dd8SToby Isaac   ierr = PetscMalloc1(pdim*pdim,&fem->invV);CHKERRQ(ierr);
7120cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
7220cf1dd8SToby Isaac   ierr = PetscMalloc1(pdim*pdim,&invVscalar);CHKERRQ(ierr);
7320cf1dd8SToby Isaac #else
7420cf1dd8SToby Isaac   invVscalar = fem->invV;
7520cf1dd8SToby Isaac #endif
7620cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
7720cf1dd8SToby Isaac     PetscReal       *Bf;
7820cf1dd8SToby Isaac     PetscQuadrature  f;
7920cf1dd8SToby Isaac     const PetscReal *points, *weights;
8020cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
8120cf1dd8SToby Isaac 
8220cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr);
8320cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
8420cf1dd8SToby Isaac     ierr = PetscMalloc1(Nc*Nq*pdim,&Bf);CHKERRQ(ierr);
8520cf1dd8SToby Isaac     ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr);
8620cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
8720cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
8820cf1dd8SToby Isaac       invVscalar[j*pdim+k] = 0.0;
8920cf1dd8SToby Isaac 
9020cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
9120cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) invVscalar[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c];
9220cf1dd8SToby Isaac       }
9320cf1dd8SToby Isaac     }
9420cf1dd8SToby Isaac     ierr = PetscFree(Bf);CHKERRQ(ierr);
9520cf1dd8SToby Isaac   }
9620cf1dd8SToby Isaac   ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr);
9720cf1dd8SToby Isaac   n = pdim;
9820cf1dd8SToby Isaac   PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invVscalar, &n, pivots, &info));
9920cf1dd8SToby Isaac   PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invVscalar, &n, pivots, work, &n, &info));
10020cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
10120cf1dd8SToby Isaac   for (j = 0; j < pdim*pdim; j++) fem->invV[j] = PetscRealPart(invVscalar[j]);
10220cf1dd8SToby Isaac   ierr = PetscFree(invVscalar);CHKERRQ(ierr);
10320cf1dd8SToby Isaac #endif
10420cf1dd8SToby Isaac   ierr = PetscFree2(pivots,work);CHKERRQ(ierr);
10520cf1dd8SToby Isaac   PetscFunctionReturn(0);
10620cf1dd8SToby Isaac }
10720cf1dd8SToby Isaac 
10820cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
10920cf1dd8SToby Isaac {
11020cf1dd8SToby Isaac   PetscErrorCode ierr;
11120cf1dd8SToby Isaac 
11220cf1dd8SToby Isaac   PetscFunctionBegin;
11320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr);
11420cf1dd8SToby Isaac   PetscFunctionReturn(0);
11520cf1dd8SToby Isaac }
11620cf1dd8SToby Isaac 
11720cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H)
11820cf1dd8SToby Isaac {
11920cf1dd8SToby Isaac   DM               dm;
12020cf1dd8SToby Isaac   PetscInt         pdim; /* Dimension of FE space P */
12120cf1dd8SToby Isaac   PetscInt         dim;  /* Spatial dimension */
12220cf1dd8SToby Isaac   PetscInt         Nc;   /* Field components */
12320cf1dd8SToby Isaac   PetscReal       *tmpB, *tmpD, *tmpH;
12420cf1dd8SToby Isaac   PetscInt         p, d, j, k, c;
12520cf1dd8SToby Isaac   PetscErrorCode   ierr;
12620cf1dd8SToby Isaac 
12720cf1dd8SToby Isaac   PetscFunctionBegin;
12820cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
12920cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
13020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
13120cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr);
13220cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
13320cf1dd8SToby Isaac   if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
13420cf1dd8SToby Isaac   if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
13520cf1dd8SToby Isaac   if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
13620cf1dd8SToby Isaac   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);CHKERRQ(ierr);
13720cf1dd8SToby Isaac   /* Translate to the nodal basis */
13820cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
13920cf1dd8SToby Isaac     if (B) {
14020cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
14120cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
14220cf1dd8SToby Isaac         const PetscInt i = (p*pdim + j)*Nc;
14320cf1dd8SToby Isaac 
14420cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
14520cf1dd8SToby Isaac           B[i+c] = 0.0;
14620cf1dd8SToby Isaac           for (k = 0; k < pdim; ++k) {
14720cf1dd8SToby Isaac             B[i+c] += fem->invV[k*pdim+j] * tmpB[(p*pdim + k)*Nc+c];
14820cf1dd8SToby Isaac           }
14920cf1dd8SToby Isaac         }
15020cf1dd8SToby Isaac       }
15120cf1dd8SToby Isaac     }
15220cf1dd8SToby Isaac     if (D) {
15320cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
15420cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
15520cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
15620cf1dd8SToby Isaac           for (d = 0; d < dim; ++d) {
15720cf1dd8SToby Isaac             const PetscInt i = ((p*pdim + j)*Nc + c)*dim + d;
15820cf1dd8SToby Isaac 
15920cf1dd8SToby Isaac             D[i] = 0.0;
16020cf1dd8SToby Isaac             for (k = 0; k < pdim; ++k) {
16120cf1dd8SToby Isaac               D[i] += fem->invV[k*pdim+j] * tmpD[((p*pdim + k)*Nc + c)*dim + d];
16220cf1dd8SToby Isaac             }
16320cf1dd8SToby Isaac           }
16420cf1dd8SToby Isaac         }
16520cf1dd8SToby Isaac       }
16620cf1dd8SToby Isaac     }
16720cf1dd8SToby Isaac     if (H) {
16820cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
16920cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
17020cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
17120cf1dd8SToby Isaac           for (d = 0; d < dim*dim; ++d) {
17220cf1dd8SToby Isaac             const PetscInt i = ((p*pdim + j)*Nc + c)*dim*dim + d;
17320cf1dd8SToby Isaac 
17420cf1dd8SToby Isaac             H[i] = 0.0;
17520cf1dd8SToby Isaac             for (k = 0; k < pdim; ++k) {
17620cf1dd8SToby Isaac               H[i] += fem->invV[k*pdim+j] * tmpH[((p*pdim + k)*Nc + c)*dim*dim + d];
17720cf1dd8SToby Isaac             }
17820cf1dd8SToby Isaac           }
17920cf1dd8SToby Isaac         }
18020cf1dd8SToby Isaac       }
18120cf1dd8SToby Isaac     }
18220cf1dd8SToby Isaac   }
18320cf1dd8SToby Isaac   if (B) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
18420cf1dd8SToby Isaac   if (D) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
18520cf1dd8SToby Isaac   if (H) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
18620cf1dd8SToby Isaac   PetscFunctionReturn(0);
18720cf1dd8SToby Isaac }
18820cf1dd8SToby Isaac 
18920cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrate_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
19020cf1dd8SToby Isaac                                       const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
19120cf1dd8SToby Isaac {
19220cf1dd8SToby Isaac   const PetscInt     debug = 0;
19320cf1dd8SToby Isaac   PetscPointFunc     obj_func;
19420cf1dd8SToby Isaac   PetscQuadrature    quad;
19520cf1dd8SToby Isaac   PetscScalar       *u, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
19620cf1dd8SToby Isaac   const PetscScalar *constants;
19720cf1dd8SToby Isaac   PetscReal         *x;
19820cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL;
19920cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
20020cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
20120cf1dd8SToby Isaac   PetscBool          isAffine;
20220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
20320cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
20420cf1dd8SToby Isaac   PetscErrorCode     ierr;
20520cf1dd8SToby Isaac 
20620cf1dd8SToby Isaac   PetscFunctionBegin;
20720cf1dd8SToby Isaac   ierr = PetscDSGetObjective(prob, field, &obj_func);CHKERRQ(ierr);
20820cf1dd8SToby Isaac   if (!obj_func) PetscFunctionReturn(0);
20920cf1dd8SToby Isaac   ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr);
21020cf1dd8SToby Isaac   ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr);
21120cf1dd8SToby Isaac   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
21220cf1dd8SToby Isaac   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
21320cf1dd8SToby Isaac   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
21420cf1dd8SToby Isaac   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
21520cf1dd8SToby Isaac   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
21620cf1dd8SToby Isaac   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
21720cf1dd8SToby Isaac   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
21820cf1dd8SToby Isaac   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
21920cf1dd8SToby Isaac   ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr);
22020cf1dd8SToby Isaac   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
22120cf1dd8SToby Isaac   if (probAux) {
22220cf1dd8SToby Isaac     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
22320cf1dd8SToby Isaac     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
22420cf1dd8SToby Isaac     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
22520cf1dd8SToby Isaac     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
22620cf1dd8SToby Isaac     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
22720cf1dd8SToby Isaac     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
22820cf1dd8SToby Isaac     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
22920cf1dd8SToby Isaac     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
23020cf1dd8SToby Isaac     ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);
23120cf1dd8SToby Isaac   }
23220cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
23320cf1dd8SToby Isaac   Np = cgeom->numPoints;
23420cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
23520cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
23620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
23720cf1dd8SToby Isaac     const PetscReal *v0   = &cgeom->v[e*Np*dE];
23820cf1dd8SToby Isaac     const PetscReal *J    = &cgeom->J[e*Np*dE*dE];
23920cf1dd8SToby Isaac 
24020cf1dd8SToby Isaac     if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
24120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
24220cf1dd8SToby Isaac       PetscScalar integrand;
24320cf1dd8SToby Isaac       const PetscReal *v;
24420cf1dd8SToby Isaac       const PetscReal *invJ;
24520cf1dd8SToby Isaac       PetscReal detJ;
24620cf1dd8SToby Isaac 
24720cf1dd8SToby Isaac       if (isAffine) {
24820cf1dd8SToby Isaac         CoordinatesRefToReal(dE, dim, cgeom->xi, v0, J, &quadPoints[q*dim], x);
24920cf1dd8SToby Isaac         v = x;
25020cf1dd8SToby Isaac         invJ = &cgeom->invJ[e*dE*dE];
25120cf1dd8SToby Isaac         detJ = cgeom->detJ[e];
25220cf1dd8SToby Isaac       } else {
25320cf1dd8SToby Isaac         v = &v0[q*dE];
25420cf1dd8SToby Isaac         invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
25520cf1dd8SToby Isaac         detJ = cgeom->detJ[e*Np + q];
25620cf1dd8SToby Isaac       }
25720cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
25820cf1dd8SToby Isaac         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", detJ);CHKERRQ(ierr);
259*7be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
26020cf1dd8SToby Isaac         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr);
26120cf1dd8SToby Isaac #endif
26220cf1dd8SToby Isaac       }
26320cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
26420cf1dd8SToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], NULL, u, u_x, NULL);
26520cf1dd8SToby Isaac       if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
26620cf1dd8SToby Isaac       obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, v, numConstants, constants, &integrand);
26720cf1dd8SToby Isaac       integrand *= detJ*quadWeights[q];
26820cf1dd8SToby Isaac       integral[e*Nf+field] += integrand;
26920cf1dd8SToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);}
27020cf1dd8SToby Isaac     }
27120cf1dd8SToby Isaac     cOffset    += totDim;
27220cf1dd8SToby Isaac     cOffsetAux += totDimAux;
27320cf1dd8SToby Isaac   }
27420cf1dd8SToby Isaac   PetscFunctionReturn(0);
27520cf1dd8SToby Isaac }
27620cf1dd8SToby Isaac 
27720cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
27820cf1dd8SToby Isaac                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
27920cf1dd8SToby Isaac {
28020cf1dd8SToby Isaac   const PetscInt     debug = 0;
28120cf1dd8SToby Isaac   PetscPointFunc     f0_func;
28220cf1dd8SToby Isaac   PetscPointFunc     f1_func;
28320cf1dd8SToby Isaac   PetscQuadrature    quad;
28420cf1dd8SToby Isaac   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
28520cf1dd8SToby Isaac   const PetscScalar *constants;
28620cf1dd8SToby Isaac   PetscReal         *x;
28720cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
28820cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
28920cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
29020cf1dd8SToby Isaac   PetscInt           dE, Np;
29120cf1dd8SToby Isaac   PetscBool          isAffine;
29220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
29320cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
29420cf1dd8SToby Isaac   PetscErrorCode     ierr;
29520cf1dd8SToby Isaac 
29620cf1dd8SToby Isaac   PetscFunctionBegin;
29720cf1dd8SToby Isaac   ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr);
29820cf1dd8SToby Isaac   ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr);
29920cf1dd8SToby Isaac   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
30020cf1dd8SToby Isaac   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
30120cf1dd8SToby Isaac   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
30220cf1dd8SToby Isaac   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
30320cf1dd8SToby Isaac   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
30420cf1dd8SToby Isaac   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
30520cf1dd8SToby Isaac   ierr = PetscDSGetFieldOffset(prob, field, &fOffset);CHKERRQ(ierr);
30620cf1dd8SToby Isaac   ierr = PetscDSGetResidual(prob, field, &f0_func, &f1_func);CHKERRQ(ierr);
30720cf1dd8SToby Isaac   ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
30820cf1dd8SToby Isaac   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
30920cf1dd8SToby Isaac   ierr = PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
31020cf1dd8SToby Isaac   ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr);
31120cf1dd8SToby Isaac   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
31220cf1dd8SToby Isaac   if (probAux) {
31320cf1dd8SToby Isaac     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
31420cf1dd8SToby Isaac     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
31520cf1dd8SToby Isaac     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
31620cf1dd8SToby Isaac     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
31720cf1dd8SToby Isaac     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
31820cf1dd8SToby Isaac     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
31920cf1dd8SToby Isaac     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
32020cf1dd8SToby Isaac     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
32120cf1dd8SToby Isaac     ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);
32220cf1dd8SToby Isaac   }
32320cf1dd8SToby Isaac   NbI = Nb[field];
32420cf1dd8SToby Isaac   NcI = Nc[field];
32520cf1dd8SToby Isaac   BI  = B[field];
32620cf1dd8SToby Isaac   DI  = D[field];
32720cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
32820cf1dd8SToby Isaac   Np = cgeom->numPoints;
32920cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
33020cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
33120cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
33220cf1dd8SToby Isaac     const PetscReal *v0   = &cgeom->v[e*Np*dE];
33320cf1dd8SToby Isaac     const PetscReal *J    = &cgeom->J[e*Np*dE*dE];
33420cf1dd8SToby Isaac 
33520cf1dd8SToby Isaac     if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
33620cf1dd8SToby Isaac     ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr);
33720cf1dd8SToby Isaac     ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr);
33820cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
33920cf1dd8SToby Isaac       const PetscReal *v;
34020cf1dd8SToby Isaac       const PetscReal *invJ;
34120cf1dd8SToby Isaac       PetscReal detJ;
34220cf1dd8SToby Isaac 
34320cf1dd8SToby Isaac       if (isAffine) {
34420cf1dd8SToby Isaac         CoordinatesRefToReal(dE, dim, cgeom->xi, v0, J, &quadPoints[q*dim], x);
34520cf1dd8SToby Isaac         v = x;
34620cf1dd8SToby Isaac         invJ = &cgeom->invJ[e*dE*dE];
34720cf1dd8SToby Isaac         detJ = cgeom->detJ[e];
34820cf1dd8SToby Isaac       } else {
34920cf1dd8SToby Isaac         v = &v0[q*dE];
35020cf1dd8SToby Isaac         invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
35120cf1dd8SToby Isaac         detJ = cgeom->detJ[e*Np + q];
35220cf1dd8SToby Isaac       }
35320cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
35420cf1dd8SToby Isaac         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", detJ);CHKERRQ(ierr);
355*7be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
35620cf1dd8SToby Isaac         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr);
35720cf1dd8SToby Isaac #endif
35820cf1dd8SToby Isaac       }
35920cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
36020cf1dd8SToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
36120cf1dd8SToby Isaac       if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
36220cf1dd8SToby Isaac       if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, numConstants, constants, &f0[q*NcI]);
36320cf1dd8SToby Isaac       if (f1_func) {
36420cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr);
36520cf1dd8SToby Isaac         f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, numConstants, constants, refSpaceDer);
36620cf1dd8SToby Isaac       }
36720cf1dd8SToby Isaac       TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL);
36820cf1dd8SToby Isaac     }
36920cf1dd8SToby Isaac     UpdateElementVec(dim, Nq, NbI, NcI, BI, DI, f0, f1, &elemVec[cOffset+fOffset]);
37020cf1dd8SToby Isaac     cOffset    += totDim;
37120cf1dd8SToby Isaac     cOffsetAux += totDimAux;
37220cf1dd8SToby Isaac   }
37320cf1dd8SToby Isaac   PetscFunctionReturn(0);
37420cf1dd8SToby Isaac }
37520cf1dd8SToby Isaac 
37620cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
37720cf1dd8SToby Isaac                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
37820cf1dd8SToby Isaac {
37920cf1dd8SToby Isaac   const PetscInt     debug = 0;
38020cf1dd8SToby Isaac   PetscBdPointFunc   f0_func;
38120cf1dd8SToby Isaac   PetscBdPointFunc   f1_func;
38220cf1dd8SToby Isaac   PetscQuadrature    quad;
38320cf1dd8SToby Isaac   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
38420cf1dd8SToby Isaac   const PetscScalar *constants;
38520cf1dd8SToby Isaac   PetscReal         *x;
38620cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
38720cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
388*7be5e748SToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
389*7be5e748SToby Isaac   PetscBool          isAffine, auxOnBd;
39020cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
39120cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
39220cf1dd8SToby Isaac   PetscErrorCode     ierr;
39320cf1dd8SToby Isaac 
39420cf1dd8SToby Isaac   PetscFunctionBegin;
39520cf1dd8SToby Isaac   ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr);
39620cf1dd8SToby Isaac   ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr);
39720cf1dd8SToby Isaac   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
39820cf1dd8SToby Isaac   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
39920cf1dd8SToby Isaac   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
40020cf1dd8SToby Isaac   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
40120cf1dd8SToby Isaac   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
40220cf1dd8SToby Isaac   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
40320cf1dd8SToby Isaac   ierr = PetscDSGetFieldOffset(prob, field, &fOffset);CHKERRQ(ierr);
40420cf1dd8SToby Isaac   ierr = PetscDSGetBdResidual(prob, field, &f0_func, &f1_func);CHKERRQ(ierr);
40520cf1dd8SToby Isaac   if (!f0_func && !f1_func) PetscFunctionReturn(0);
40620cf1dd8SToby Isaac   ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
40720cf1dd8SToby Isaac   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
40820cf1dd8SToby Isaac   ierr = PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
40920cf1dd8SToby Isaac   ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr);
41020cf1dd8SToby Isaac   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
41120cf1dd8SToby Isaac   if (probAux) {
412*7be5e748SToby Isaac     ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr);
41320cf1dd8SToby Isaac     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
41420cf1dd8SToby Isaac     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
41520cf1dd8SToby Isaac     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
41620cf1dd8SToby Isaac     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
41720cf1dd8SToby Isaac     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
41820cf1dd8SToby Isaac     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
41920cf1dd8SToby Isaac     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
42020cf1dd8SToby Isaac     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
421*7be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
422*7be5e748SToby Isaac     if (auxOnBd) {ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);}
423*7be5e748SToby Isaac     else         {ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);}
42420cf1dd8SToby Isaac   }
42520cf1dd8SToby Isaac   NbI = Nb[field];
42620cf1dd8SToby Isaac   NcI = Nc[field];
42720cf1dd8SToby Isaac   BI  = B[field];
42820cf1dd8SToby Isaac   DI  = D[field];
42920cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
43020cf1dd8SToby Isaac   Np = fgeom->numPoints;
43120cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
43220cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
43320cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
43420cf1dd8SToby Isaac     const PetscReal *v0   = &fgeom->v[e*Np*dE];
43520cf1dd8SToby Isaac     const PetscReal *J    = &fgeom->J[e*Np*dE*dE];
43620cf1dd8SToby Isaac     const PetscInt   face = fgeom->face[e][0];
43720cf1dd8SToby Isaac 
43820cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
43920cf1dd8SToby Isaac     if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
44020cf1dd8SToby Isaac     ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr);
44120cf1dd8SToby Isaac     ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr);
44220cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
44320cf1dd8SToby Isaac       const PetscReal *v;
44420cf1dd8SToby Isaac       const PetscReal *invJ;
44520cf1dd8SToby Isaac       const PetscReal *n;
44620cf1dd8SToby Isaac       PetscReal detJ;
44720cf1dd8SToby Isaac       if (isAffine) {
44820cf1dd8SToby Isaac         CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x);
44920cf1dd8SToby Isaac         v = x;
45020cf1dd8SToby Isaac         invJ = &fgeom->suppInvJ[0][e*dE*dE];
45120cf1dd8SToby Isaac         detJ = fgeom->detJ[e];
45220cf1dd8SToby Isaac         n    = &fgeom->n[e*dE];
45320cf1dd8SToby Isaac       } else {
45420cf1dd8SToby Isaac         v = &v0[q*dE];
45520cf1dd8SToby Isaac         invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
45620cf1dd8SToby Isaac         detJ = fgeom->detJ[e*Np + q];
45720cf1dd8SToby Isaac         n    = &fgeom->n[(e*Np+q)*dE];
45820cf1dd8SToby Isaac       }
45920cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
46020cf1dd8SToby Isaac         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", detJ);CHKERRQ(ierr);
461*7be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
46220cf1dd8SToby Isaac         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr);
46320cf1dd8SToby Isaac #endif
46420cf1dd8SToby Isaac       }
46520cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
46620cf1dd8SToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
467*7be5e748SToby Isaac       if (probAux) EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
46820cf1dd8SToby Isaac       if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, n, numConstants, constants, &f0[q*NcI]);
46920cf1dd8SToby Isaac       if (f1_func) {
47020cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr);
47120cf1dd8SToby Isaac         f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, n, numConstants, constants, refSpaceDer);
47220cf1dd8SToby Isaac       }
47320cf1dd8SToby Isaac       TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL);
47420cf1dd8SToby Isaac     }
47520cf1dd8SToby Isaac     UpdateElementVec(dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], f0, f1, &elemVec[cOffset+fOffset]);
47620cf1dd8SToby Isaac     cOffset    += totDim;
47720cf1dd8SToby Isaac     cOffsetAux += totDimAux;
47820cf1dd8SToby Isaac   }
47920cf1dd8SToby Isaac   PetscFunctionReturn(0);
48020cf1dd8SToby Isaac }
48120cf1dd8SToby Isaac 
48220cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *geom,
48320cf1dd8SToby Isaac                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
48420cf1dd8SToby Isaac {
48520cf1dd8SToby Isaac   const PetscInt     debug      = 0;
48620cf1dd8SToby Isaac   PetscPointJac      g0_func;
48720cf1dd8SToby Isaac   PetscPointJac      g1_func;
48820cf1dd8SToby Isaac   PetscPointJac      g2_func;
48920cf1dd8SToby Isaac   PetscPointJac      g3_func;
49020cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
49120cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
49220cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
49320cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
49420cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
49520cf1dd8SToby Isaac   PetscQuadrature    quad;
49620cf1dd8SToby Isaac   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
49720cf1dd8SToby Isaac   const PetscScalar *constants;
49820cf1dd8SToby Isaac   PetscReal         *x;
49920cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
50020cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
50120cf1dd8SToby Isaac   PetscInt           NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
50220cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
50320cf1dd8SToby Isaac   PetscInt           dE, Np;
50420cf1dd8SToby Isaac   PetscBool          isAffine;
50520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
50620cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
50720cf1dd8SToby Isaac   PetscErrorCode     ierr;
50820cf1dd8SToby Isaac 
50920cf1dd8SToby Isaac   PetscFunctionBegin;
51020cf1dd8SToby Isaac   ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr);
51120cf1dd8SToby Isaac   ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr);
51220cf1dd8SToby Isaac   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
51320cf1dd8SToby Isaac   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
51420cf1dd8SToby Isaac   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
51520cf1dd8SToby Isaac   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
51620cf1dd8SToby Isaac   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
51720cf1dd8SToby Isaac   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
51820cf1dd8SToby Isaac   switch(jtype) {
51920cf1dd8SToby Isaac   case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
52020cf1dd8SToby Isaac   case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
52120cf1dd8SToby Isaac   case PETSCFE_JACOBIAN:     ierr = PetscDSGetJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
52220cf1dd8SToby Isaac   }
52320cf1dd8SToby Isaac   if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0);
52420cf1dd8SToby Isaac   ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
52520cf1dd8SToby Isaac   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
52620cf1dd8SToby Isaac   ierr = PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
52720cf1dd8SToby Isaac   ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr);
52820cf1dd8SToby Isaac   ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr);
52920cf1dd8SToby Isaac   ierr = PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);CHKERRQ(ierr);
53020cf1dd8SToby Isaac   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
53120cf1dd8SToby Isaac   if (probAux) {
53220cf1dd8SToby Isaac     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
53320cf1dd8SToby Isaac     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
53420cf1dd8SToby Isaac     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
53520cf1dd8SToby Isaac     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
53620cf1dd8SToby Isaac     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
53720cf1dd8SToby Isaac     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
53820cf1dd8SToby Isaac     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
53920cf1dd8SToby Isaac     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
54020cf1dd8SToby Isaac     ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);
54120cf1dd8SToby Isaac   }
54220cf1dd8SToby Isaac   NbI = Nb[fieldI], NbJ = Nb[fieldJ];
54320cf1dd8SToby Isaac   NcI = Nc[fieldI], NcJ = Nc[fieldJ];
54420cf1dd8SToby Isaac   BI  = B[fieldI],  BJ  = B[fieldJ];
54520cf1dd8SToby Isaac   DI  = D[fieldI],  DJ  = D[fieldJ];
54620cf1dd8SToby Isaac   /* Initialize here in case the function is not defined */
54720cf1dd8SToby Isaac   ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
54820cf1dd8SToby Isaac   ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
54920cf1dd8SToby Isaac   ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
55020cf1dd8SToby Isaac   ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
55120cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
55220cf1dd8SToby Isaac   Np = geom->numPoints;
55320cf1dd8SToby Isaac   dE = geom->dimEmbed;
55420cf1dd8SToby Isaac   isAffine = geom->isAffine;
55520cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
55620cf1dd8SToby Isaac     const PetscReal *v0   = &geom->v[e*Np*dE];
55720cf1dd8SToby Isaac     const PetscReal *J    = &geom->J[e*Np*dE*dE];
55820cf1dd8SToby Isaac 
55920cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
56020cf1dd8SToby Isaac     if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
56120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
56220cf1dd8SToby Isaac       const PetscReal *v;
56320cf1dd8SToby Isaac       const PetscReal *invJ;
56420cf1dd8SToby Isaac       PetscReal detJ;
56520cf1dd8SToby Isaac       const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ];
56620cf1dd8SToby Isaac       const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim];
56720cf1dd8SToby Isaac       PetscReal  w;
56820cf1dd8SToby Isaac       PetscInt f, g, fc, gc, c;
56920cf1dd8SToby Isaac 
57020cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
57120cf1dd8SToby Isaac       if (isAffine) {
57220cf1dd8SToby Isaac         CoordinatesRefToReal(dE, dim, geom->xi, v0, J, &quadPoints[q*dim], x);
57320cf1dd8SToby Isaac         v = x;
57420cf1dd8SToby Isaac         invJ = &geom->invJ[e*dE*dE];
57520cf1dd8SToby Isaac         detJ = geom->detJ[e];
57620cf1dd8SToby Isaac       } else {
57720cf1dd8SToby Isaac         v = &v0[q*dE];
57820cf1dd8SToby Isaac         invJ = &geom->invJ[(e*Np+q)*dE*dE];
57920cf1dd8SToby Isaac         detJ = geom->detJ[e*Np + q];
58020cf1dd8SToby Isaac       }
58120cf1dd8SToby Isaac       w = detJ*quadWeights[q];
58220cf1dd8SToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
58320cf1dd8SToby Isaac       if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
58420cf1dd8SToby Isaac       if (g0_func) {
58520cf1dd8SToby Isaac         ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
58620cf1dd8SToby Isaac         g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, g0);
58720cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
58820cf1dd8SToby Isaac       }
58920cf1dd8SToby Isaac       if (g1_func) {
59020cf1dd8SToby Isaac         PetscInt d, d2;
59120cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
59220cf1dd8SToby Isaac         g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer);
59320cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
59420cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
59520cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
59620cf1dd8SToby Isaac               g1[(fc*NcJ+gc)*dim+d] = 0.0;
59720cf1dd8SToby Isaac               for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
59820cf1dd8SToby Isaac               g1[(fc*NcJ+gc)*dim+d] *= w;
59920cf1dd8SToby Isaac             }
60020cf1dd8SToby Isaac           }
60120cf1dd8SToby Isaac         }
60220cf1dd8SToby Isaac       }
60320cf1dd8SToby Isaac       if (g2_func) {
60420cf1dd8SToby Isaac         PetscInt d, d2;
60520cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
60620cf1dd8SToby Isaac         g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer);
60720cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
60820cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
60920cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
61020cf1dd8SToby Isaac               g2[(fc*NcJ+gc)*dim+d] = 0.0;
61120cf1dd8SToby Isaac               for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
61220cf1dd8SToby Isaac               g2[(fc*NcJ+gc)*dim+d] *= w;
61320cf1dd8SToby Isaac             }
61420cf1dd8SToby Isaac           }
61520cf1dd8SToby Isaac         }
61620cf1dd8SToby Isaac       }
61720cf1dd8SToby Isaac       if (g3_func) {
61820cf1dd8SToby Isaac         PetscInt d, d2, dp, d3;
61920cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
62020cf1dd8SToby Isaac         g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer);
62120cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
62220cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
62320cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
62420cf1dd8SToby Isaac               for (dp = 0; dp < dim; ++dp) {
62520cf1dd8SToby Isaac                 g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0;
62620cf1dd8SToby Isaac                 for (d2 = 0; d2 < dim; ++d2) {
62720cf1dd8SToby Isaac                   for (d3 = 0; d3 < dim; ++d3) {
62820cf1dd8SToby Isaac                     g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3];
62920cf1dd8SToby Isaac                   }
63020cf1dd8SToby Isaac                 }
63120cf1dd8SToby Isaac                 g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w;
63220cf1dd8SToby Isaac               }
63320cf1dd8SToby Isaac             }
63420cf1dd8SToby Isaac           }
63520cf1dd8SToby Isaac         }
63620cf1dd8SToby Isaac       }
63720cf1dd8SToby Isaac 
63820cf1dd8SToby Isaac       for (f = 0; f < NbI; ++f) {
63920cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
64020cf1dd8SToby Isaac           const PetscInt fidx = f*NcI+fc; /* Test function basis index */
64120cf1dd8SToby Isaac           const PetscInt i    = offsetI+f; /* Element matrix row */
64220cf1dd8SToby Isaac           for (g = 0; g < NbJ; ++g) {
64320cf1dd8SToby Isaac             for (gc = 0; gc < NcJ; ++gc) {
64420cf1dd8SToby Isaac               const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
64520cf1dd8SToby Isaac               const PetscInt j    = offsetJ+g; /* Element matrix column */
64620cf1dd8SToby Isaac               const PetscInt fOff = eOffset+i*totDim+j;
64720cf1dd8SToby Isaac               PetscInt       d, d2;
64820cf1dd8SToby Isaac 
64920cf1dd8SToby Isaac               elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx];
65020cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
65120cf1dd8SToby Isaac                 elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d];
65220cf1dd8SToby Isaac                 elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx];
65320cf1dd8SToby Isaac                 for (d2 = 0; d2 < dim; ++d2) {
65420cf1dd8SToby Isaac                   elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2];
65520cf1dd8SToby Isaac                 }
65620cf1dd8SToby Isaac               }
65720cf1dd8SToby Isaac             }
65820cf1dd8SToby Isaac           }
65920cf1dd8SToby Isaac         }
66020cf1dd8SToby Isaac       }
66120cf1dd8SToby Isaac     }
66220cf1dd8SToby Isaac     if (debug > 1) {
66320cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
66420cf1dd8SToby Isaac 
66520cf1dd8SToby Isaac       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
66620cf1dd8SToby Isaac       for (fc = 0; fc < NcI; ++fc) {
66720cf1dd8SToby Isaac         for (f = 0; f < NbI; ++f) {
66820cf1dd8SToby Isaac           const PetscInt i = offsetI + f*NcI+fc;
66920cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
67020cf1dd8SToby Isaac             for (g = 0; g < NbJ; ++g) {
67120cf1dd8SToby Isaac               const PetscInt j = offsetJ + g*NcJ+gc;
67220cf1dd8SToby Isaac               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
67320cf1dd8SToby Isaac             }
67420cf1dd8SToby Isaac           }
67520cf1dd8SToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
67620cf1dd8SToby Isaac         }
67720cf1dd8SToby Isaac       }
67820cf1dd8SToby Isaac     }
67920cf1dd8SToby Isaac     cOffset    += totDim;
68020cf1dd8SToby Isaac     cOffsetAux += totDimAux;
68120cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
68220cf1dd8SToby Isaac   }
68320cf1dd8SToby Isaac   PetscFunctionReturn(0);
68420cf1dd8SToby Isaac }
68520cf1dd8SToby Isaac 
68620cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
68720cf1dd8SToby Isaac                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
68820cf1dd8SToby Isaac {
68920cf1dd8SToby Isaac   const PetscInt     debug      = 0;
69020cf1dd8SToby Isaac   PetscBdPointJac    g0_func;
69120cf1dd8SToby Isaac   PetscBdPointJac    g1_func;
69220cf1dd8SToby Isaac   PetscBdPointJac    g2_func;
69320cf1dd8SToby Isaac   PetscBdPointJac    g3_func;
69420cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
69520cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
69620cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
69720cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
69820cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
69920cf1dd8SToby Isaac   PetscQuadrature    quad;
70020cf1dd8SToby Isaac   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
70120cf1dd8SToby Isaac   const PetscScalar *constants;
70220cf1dd8SToby Isaac   PetscReal         *x;
70320cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
70420cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
70520cf1dd8SToby Isaac   PetscInt           NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
70620cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
70720cf1dd8SToby Isaac   PetscBool          isAffine;
70820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
70920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
71020cf1dd8SToby Isaac   PetscErrorCode     ierr;
71120cf1dd8SToby Isaac 
71220cf1dd8SToby Isaac   PetscFunctionBegin;
71320cf1dd8SToby Isaac   ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr);
71420cf1dd8SToby Isaac   ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr);
71520cf1dd8SToby Isaac   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
71620cf1dd8SToby Isaac   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
71720cf1dd8SToby Isaac   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
71820cf1dd8SToby Isaac   ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr);
71920cf1dd8SToby Isaac   ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr);
72020cf1dd8SToby Isaac   ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr);
72120cf1dd8SToby Isaac   ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr);
72220cf1dd8SToby Isaac   ierr = PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);CHKERRQ(ierr);
72320cf1dd8SToby Isaac   ierr = PetscDSGetBdJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);
72420cf1dd8SToby Isaac   ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
72520cf1dd8SToby Isaac   ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr);
72620cf1dd8SToby Isaac   ierr = PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
72720cf1dd8SToby Isaac   ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr);
72820cf1dd8SToby Isaac   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
72920cf1dd8SToby Isaac   if (probAux) {
73020cf1dd8SToby Isaac     ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr);
73120cf1dd8SToby Isaac     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
73220cf1dd8SToby Isaac     ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr);
73320cf1dd8SToby Isaac     ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr);
73420cf1dd8SToby Isaac     ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr);
73520cf1dd8SToby Isaac     ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr);
73620cf1dd8SToby Isaac     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
73720cf1dd8SToby Isaac     ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr);
73820cf1dd8SToby Isaac     ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);
73920cf1dd8SToby Isaac   }
74020cf1dd8SToby Isaac   NbI = Nb[fieldI], NbJ = Nb[fieldJ];
74120cf1dd8SToby Isaac   NcI = Nc[fieldI], NcJ = Nc[fieldJ];
74220cf1dd8SToby Isaac   BI  = B[fieldI],  BJ  = B[fieldJ];
74320cf1dd8SToby Isaac   DI  = D[fieldI],  DJ  = D[fieldJ];
74420cf1dd8SToby Isaac   /* Initialize here in case the function is not defined */
74520cf1dd8SToby Isaac   ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
74620cf1dd8SToby Isaac   ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
74720cf1dd8SToby Isaac   ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
74820cf1dd8SToby Isaac   ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
74920cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
75020cf1dd8SToby Isaac   Np = fgeom->numPoints;
75120cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
75220cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
75320cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
75420cf1dd8SToby Isaac     const PetscReal *v0   = &fgeom->v[e*Np*dE];
75520cf1dd8SToby Isaac     const PetscReal *J    = &fgeom->J[e*Np*dE*dE];
75620cf1dd8SToby Isaac     const PetscInt   face = fgeom->face[e][0];
75720cf1dd8SToby Isaac 
75820cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
75920cf1dd8SToby Isaac     if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
76020cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
76120cf1dd8SToby Isaac       const PetscReal *BIq = &BI[(face*Nq+q)*NbI*NcI], *BJq = &BJ[(face*Nq+q)*NbJ*NcJ];
76220cf1dd8SToby Isaac       const PetscReal *DIq = &DI[(face*Nq+q)*NbI*NcI*dim], *DJq = &DJ[(face*Nq+q)*NbJ*NcJ*dim];
76320cf1dd8SToby Isaac       PetscReal  w;
76420cf1dd8SToby Isaac       PetscInt f, g, fc, gc, c;
76520cf1dd8SToby Isaac       const PetscReal *v;
76620cf1dd8SToby Isaac       const PetscReal *invJ;
76720cf1dd8SToby Isaac       const PetscReal *n;
76820cf1dd8SToby Isaac       PetscReal detJ;
76920cf1dd8SToby Isaac 
77020cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
77120cf1dd8SToby Isaac       if (isAffine) {
77220cf1dd8SToby Isaac         CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x);
77320cf1dd8SToby Isaac         v = x;
77420cf1dd8SToby Isaac         invJ = &fgeom->suppInvJ[0][e*dE*dE];
77520cf1dd8SToby Isaac         detJ = fgeom->detJ[e];
77620cf1dd8SToby Isaac         n    = &fgeom->n[e*dE];
77720cf1dd8SToby Isaac       } else {
77820cf1dd8SToby Isaac         v = &v0[q*dE];
77920cf1dd8SToby Isaac         invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
78020cf1dd8SToby Isaac         detJ = fgeom->detJ[e*Np + q];
78120cf1dd8SToby Isaac         n    = &fgeom->n[(e*Np+q)*dE];
78220cf1dd8SToby Isaac       }
78320cf1dd8SToby Isaac       w = detJ*quadWeights[q];
78420cf1dd8SToby Isaac 
78520cf1dd8SToby Isaac       EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
78620cf1dd8SToby Isaac       if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
78720cf1dd8SToby Isaac       if (g0_func) {
78820cf1dd8SToby Isaac         ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
78920cf1dd8SToby Isaac         g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, g0);
79020cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
79120cf1dd8SToby Isaac       }
79220cf1dd8SToby Isaac       if (g1_func) {
79320cf1dd8SToby Isaac         PetscInt d, d2;
79420cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
79520cf1dd8SToby Isaac         g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer);
79620cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
79720cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
79820cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
79920cf1dd8SToby Isaac               g1[(fc*NcJ+gc)*dim+d] = 0.0;
80020cf1dd8SToby Isaac               for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
80120cf1dd8SToby Isaac               g1[(fc*NcJ+gc)*dim+d] *= w;
80220cf1dd8SToby Isaac             }
80320cf1dd8SToby Isaac           }
80420cf1dd8SToby Isaac         }
80520cf1dd8SToby Isaac       }
80620cf1dd8SToby Isaac       if (g2_func) {
80720cf1dd8SToby Isaac         PetscInt d, d2;
80820cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
80920cf1dd8SToby Isaac         g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer);
81020cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
81120cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
81220cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
81320cf1dd8SToby Isaac               g2[(fc*NcJ+gc)*dim+d] = 0.0;
81420cf1dd8SToby Isaac               for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
81520cf1dd8SToby Isaac               g2[(fc*NcJ+gc)*dim+d] *= w;
81620cf1dd8SToby Isaac             }
81720cf1dd8SToby Isaac           }
81820cf1dd8SToby Isaac         }
81920cf1dd8SToby Isaac       }
82020cf1dd8SToby Isaac       if (g3_func) {
82120cf1dd8SToby Isaac         PetscInt d, d2, dp, d3;
82220cf1dd8SToby Isaac         ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
82320cf1dd8SToby Isaac         g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer);
82420cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
82520cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
82620cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
82720cf1dd8SToby Isaac               for (dp = 0; dp < dim; ++dp) {
82820cf1dd8SToby Isaac                 g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0;
82920cf1dd8SToby Isaac                 for (d2 = 0; d2 < dim; ++d2) {
83020cf1dd8SToby Isaac                   for (d3 = 0; d3 < dim; ++d3) {
83120cf1dd8SToby Isaac                     g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3];
83220cf1dd8SToby Isaac                   }
83320cf1dd8SToby Isaac                 }
83420cf1dd8SToby Isaac                 g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w;
83520cf1dd8SToby Isaac               }
83620cf1dd8SToby Isaac             }
83720cf1dd8SToby Isaac           }
83820cf1dd8SToby Isaac         }
83920cf1dd8SToby Isaac       }
84020cf1dd8SToby Isaac 
84120cf1dd8SToby Isaac       for (f = 0; f < NbI; ++f) {
84220cf1dd8SToby Isaac         for (fc = 0; fc < NcI; ++fc) {
84320cf1dd8SToby Isaac           const PetscInt fidx = f*NcI+fc; /* Test function basis index */
84420cf1dd8SToby Isaac           const PetscInt i    = offsetI+f; /* Element matrix row */
84520cf1dd8SToby Isaac           for (g = 0; g < NbJ; ++g) {
84620cf1dd8SToby Isaac             for (gc = 0; gc < NcJ; ++gc) {
84720cf1dd8SToby Isaac               const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
84820cf1dd8SToby Isaac               const PetscInt j    = offsetJ+g; /* Element matrix column */
84920cf1dd8SToby Isaac               const PetscInt fOff = eOffset+i*totDim+j;
85020cf1dd8SToby Isaac               PetscInt       d, d2;
85120cf1dd8SToby Isaac 
85220cf1dd8SToby Isaac               elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx];
85320cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
85420cf1dd8SToby Isaac                 elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d];
85520cf1dd8SToby Isaac                 elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx];
85620cf1dd8SToby Isaac                 for (d2 = 0; d2 < dim; ++d2) {
85720cf1dd8SToby Isaac                   elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2];
85820cf1dd8SToby Isaac                 }
85920cf1dd8SToby Isaac               }
86020cf1dd8SToby Isaac             }
86120cf1dd8SToby Isaac           }
86220cf1dd8SToby Isaac         }
86320cf1dd8SToby Isaac       }
86420cf1dd8SToby Isaac     }
86520cf1dd8SToby Isaac     if (debug > 1) {
86620cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
86720cf1dd8SToby Isaac 
86820cf1dd8SToby Isaac       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
86920cf1dd8SToby Isaac       for (fc = 0; fc < NcI; ++fc) {
87020cf1dd8SToby Isaac         for (f = 0; f < NbI; ++f) {
87120cf1dd8SToby Isaac           const PetscInt i = offsetI + f*NcI+fc;
87220cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
87320cf1dd8SToby Isaac             for (g = 0; g < NbJ; ++g) {
87420cf1dd8SToby Isaac               const PetscInt j = offsetJ + g*NcJ+gc;
87520cf1dd8SToby Isaac               ierr = PetscPrintf(PETSC_COMM_SELF, "    elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr);
87620cf1dd8SToby Isaac             }
87720cf1dd8SToby Isaac           }
87820cf1dd8SToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
87920cf1dd8SToby Isaac         }
88020cf1dd8SToby Isaac       }
88120cf1dd8SToby Isaac     }
88220cf1dd8SToby Isaac     cOffset    += totDim;
88320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
88420cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
88520cf1dd8SToby Isaac   }
88620cf1dd8SToby Isaac   PetscFunctionReturn(0);
88720cf1dd8SToby Isaac }
88820cf1dd8SToby Isaac 
88920cf1dd8SToby Isaac PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
89020cf1dd8SToby Isaac {
89120cf1dd8SToby Isaac   PetscFunctionBegin;
89220cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
89320cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
89420cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
89520cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
89620cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
89720cf1dd8SToby Isaac   fem->ops->gettabulation           = PetscFEGetTabulation_Basic;
89820cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
89920cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
90020cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
90120cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
90220cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
90320cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
90420cf1dd8SToby Isaac   PetscFunctionReturn(0);
90520cf1dd8SToby Isaac }
90620cf1dd8SToby Isaac 
90720cf1dd8SToby Isaac /*MC
90820cf1dd8SToby Isaac   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
90920cf1dd8SToby Isaac 
91020cf1dd8SToby Isaac   Level: intermediate
91120cf1dd8SToby Isaac 
91220cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
91320cf1dd8SToby Isaac M*/
91420cf1dd8SToby Isaac 
91520cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
91620cf1dd8SToby Isaac {
91720cf1dd8SToby Isaac   PetscFE_Basic *b;
91820cf1dd8SToby Isaac   PetscErrorCode ierr;
91920cf1dd8SToby Isaac 
92020cf1dd8SToby Isaac   PetscFunctionBegin;
92120cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
92220cf1dd8SToby Isaac   ierr      = PetscNewLog(fem,&b);CHKERRQ(ierr);
92320cf1dd8SToby Isaac   fem->data = b;
92420cf1dd8SToby Isaac 
92520cf1dd8SToby Isaac   ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr);
92620cf1dd8SToby Isaac   PetscFunctionReturn(0);
92720cf1dd8SToby Isaac }
92820cf1dd8SToby Isaac 
92920cf1dd8SToby Isaac 
930