xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision ef0bb6c736604ce380bf8bea4ebd4a7bda431d97)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscblaslapack.h>
320cf1dd8SToby Isaac 
42b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
520cf1dd8SToby Isaac {
620cf1dd8SToby Isaac   PetscFE_Basic *b = (PetscFE_Basic *) fem->data;
720cf1dd8SToby Isaac   PetscErrorCode ierr;
820cf1dd8SToby Isaac 
920cf1dd8SToby Isaac   PetscFunctionBegin;
1020cf1dd8SToby Isaac   ierr = PetscFree(b);CHKERRQ(ierr);
1120cf1dd8SToby Isaac   PetscFunctionReturn(0);
1220cf1dd8SToby Isaac }
1320cf1dd8SToby Isaac 
142b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v)
1520cf1dd8SToby Isaac {
16d9bac1caSLisandro Dalcin   PetscInt          dim, Nc;
17d9bac1caSLisandro Dalcin   PetscSpace        basis = NULL;
18d9bac1caSLisandro Dalcin   PetscDualSpace    dual = NULL;
19d9bac1caSLisandro Dalcin   PetscQuadrature   quad = NULL;
2020cf1dd8SToby Isaac   PetscErrorCode    ierr;
2120cf1dd8SToby Isaac 
2220cf1dd8SToby Isaac   PetscFunctionBegin;
23d9bac1caSLisandro Dalcin   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
24d9bac1caSLisandro Dalcin   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
2520cf1dd8SToby Isaac   ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr);
2620cf1dd8SToby Isaac   ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr);
27d9bac1caSLisandro Dalcin   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
28d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
29d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPrintf(v, "Basic Finite Element in %D dimensions with %D components\n",dim,Nc);CHKERRQ(ierr);
30d9bac1caSLisandro Dalcin   if (basis) {ierr = PetscSpaceView(basis, v);CHKERRQ(ierr);}
31d9bac1caSLisandro Dalcin   if (dual)  {ierr = PetscDualSpaceView(dual, v);CHKERRQ(ierr);}
32d9bac1caSLisandro Dalcin   if (quad)  {ierr = PetscQuadratureView(quad, v);CHKERRQ(ierr);}
33d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
3420cf1dd8SToby Isaac   PetscFunctionReturn(0);
3520cf1dd8SToby Isaac }
3620cf1dd8SToby Isaac 
372b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v)
3820cf1dd8SToby Isaac {
3920cf1dd8SToby Isaac   PetscBool      iascii;
4020cf1dd8SToby Isaac   PetscErrorCode ierr;
4120cf1dd8SToby Isaac 
4220cf1dd8SToby Isaac   PetscFunctionBegin;
43d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
44d9bac1caSLisandro Dalcin   if (iascii) {ierr = PetscFEView_Basic_Ascii(fe, v);CHKERRQ(ierr);}
4520cf1dd8SToby Isaac   PetscFunctionReturn(0);
4620cf1dd8SToby Isaac }
4720cf1dd8SToby Isaac 
4820cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */
4915310ec5SMatthew G. Knepley PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
5020cf1dd8SToby Isaac {
5120cf1dd8SToby Isaac   PetscScalar   *work, *invVscalar;
5220cf1dd8SToby Isaac   PetscBLASInt  *pivots;
5320cf1dd8SToby Isaac   PetscBLASInt   n, info;
5420cf1dd8SToby Isaac   PetscInt       pdim, j;
5520cf1dd8SToby Isaac   PetscErrorCode ierr;
5620cf1dd8SToby Isaac 
5720cf1dd8SToby Isaac   PetscFunctionBegin;
5820cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
5920cf1dd8SToby Isaac   ierr = PetscMalloc1(pdim*pdim,&fem->invV);CHKERRQ(ierr);
6020cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
6120cf1dd8SToby Isaac   ierr = PetscMalloc1(pdim*pdim,&invVscalar);CHKERRQ(ierr);
6220cf1dd8SToby Isaac #else
6320cf1dd8SToby Isaac   invVscalar = fem->invV;
6420cf1dd8SToby Isaac #endif
6520cf1dd8SToby Isaac   for (j = 0; j < pdim; ++j) {
6620cf1dd8SToby Isaac     PetscReal       *Bf;
6720cf1dd8SToby Isaac     PetscQuadrature  f;
6820cf1dd8SToby Isaac     const PetscReal *points, *weights;
6920cf1dd8SToby Isaac     PetscInt         Nc, Nq, q, k, c;
7020cf1dd8SToby Isaac 
7120cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr);
7220cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
7320cf1dd8SToby Isaac     ierr = PetscMalloc1(Nc*Nq*pdim,&Bf);CHKERRQ(ierr);
7420cf1dd8SToby Isaac     ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr);
7520cf1dd8SToby Isaac     for (k = 0; k < pdim; ++k) {
7620cf1dd8SToby Isaac       /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
7720cf1dd8SToby Isaac       invVscalar[j*pdim+k] = 0.0;
7820cf1dd8SToby Isaac 
7920cf1dd8SToby Isaac       for (q = 0; q < Nq; ++q) {
8020cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) invVscalar[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c];
8120cf1dd8SToby Isaac       }
8220cf1dd8SToby Isaac     }
8320cf1dd8SToby Isaac     ierr = PetscFree(Bf);CHKERRQ(ierr);
8420cf1dd8SToby Isaac   }
85ea2bdf6dSBarry Smith 
8620cf1dd8SToby Isaac   ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr);
8720cf1dd8SToby Isaac   n = pdim;
8820cf1dd8SToby Isaac   PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invVscalar, &n, pivots, &info));
89ea2bdf6dSBarry Smith   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
90449a1c5fSBarry Smith   PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invVscalar, &n, pivots, work, &n, &info));
91ea2bdf6dSBarry Smith   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
9220cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
9320cf1dd8SToby Isaac   for (j = 0; j < pdim*pdim; j++) fem->invV[j] = PetscRealPart(invVscalar[j]);
9420cf1dd8SToby Isaac   ierr = PetscFree(invVscalar);CHKERRQ(ierr);
9520cf1dd8SToby Isaac #endif
9620cf1dd8SToby Isaac   ierr = PetscFree2(pivots,work);CHKERRQ(ierr);
9720cf1dd8SToby Isaac   PetscFunctionReturn(0);
9820cf1dd8SToby Isaac }
9920cf1dd8SToby Isaac 
10020cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
10120cf1dd8SToby Isaac {
10220cf1dd8SToby Isaac   PetscErrorCode ierr;
10320cf1dd8SToby Isaac 
10420cf1dd8SToby Isaac   PetscFunctionBegin;
10520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr);
10620cf1dd8SToby Isaac   PetscFunctionReturn(0);
10720cf1dd8SToby Isaac }
10820cf1dd8SToby Isaac 
109*ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
11020cf1dd8SToby Isaac {
11120cf1dd8SToby Isaac   DM               dm;
11220cf1dd8SToby Isaac   PetscInt         pdim; /* Dimension of FE space P */
11320cf1dd8SToby Isaac   PetscInt         dim;  /* Spatial dimension */
11420cf1dd8SToby Isaac   PetscInt         Nc;   /* Field components */
115*ef0bb6c7SMatthew G. Knepley   PetscReal       *B = K >= 0 ? T->T[0] : NULL;
116*ef0bb6c7SMatthew G. Knepley   PetscReal       *D = K >= 1 ? T->T[1] : NULL;
117*ef0bb6c7SMatthew G. Knepley   PetscReal       *H = K >= 2 ? T->T[2] : NULL;
118*ef0bb6c7SMatthew G. Knepley   PetscReal       *tmpB = NULL, *tmpD = NULL, *tmpH = NULL;
11920cf1dd8SToby Isaac   PetscInt         p, d, j, k, c;
12020cf1dd8SToby Isaac   PetscErrorCode   ierr;
12120cf1dd8SToby Isaac 
12220cf1dd8SToby Isaac   PetscFunctionBegin;
12320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
12420cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
12520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
12620cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr);
12720cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
128*ef0bb6c7SMatthew G. Knepley   if (K >= 0) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
129*ef0bb6c7SMatthew G. Knepley   if (K >= 1) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
130*ef0bb6c7SMatthew G. Knepley   if (K >= 2) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
131*ef0bb6c7SMatthew G. Knepley   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH);CHKERRQ(ierr);
13220cf1dd8SToby Isaac   /* Translate to the nodal basis */
13320cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
13420cf1dd8SToby Isaac     if (B) {
13520cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
13620cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
13720cf1dd8SToby Isaac         const PetscInt i = (p*pdim + j)*Nc;
13820cf1dd8SToby Isaac 
13920cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
14020cf1dd8SToby Isaac           B[i+c] = 0.0;
14120cf1dd8SToby Isaac           for (k = 0; k < pdim; ++k) {
14220cf1dd8SToby Isaac             B[i+c] += fem->invV[k*pdim+j] * tmpB[(p*pdim + k)*Nc+c];
14320cf1dd8SToby Isaac           }
14420cf1dd8SToby Isaac         }
14520cf1dd8SToby Isaac       }
14620cf1dd8SToby Isaac     }
14720cf1dd8SToby Isaac     if (D) {
14820cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
14920cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
15020cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
15120cf1dd8SToby Isaac           for (d = 0; d < dim; ++d) {
15220cf1dd8SToby Isaac             const PetscInt i = ((p*pdim + j)*Nc + c)*dim + d;
15320cf1dd8SToby Isaac 
15420cf1dd8SToby Isaac             D[i] = 0.0;
15520cf1dd8SToby Isaac             for (k = 0; k < pdim; ++k) {
15620cf1dd8SToby Isaac               D[i] += fem->invV[k*pdim+j] * tmpD[((p*pdim + k)*Nc + c)*dim + d];
15720cf1dd8SToby Isaac             }
15820cf1dd8SToby Isaac           }
15920cf1dd8SToby Isaac         }
16020cf1dd8SToby Isaac       }
16120cf1dd8SToby Isaac     }
16220cf1dd8SToby Isaac     if (H) {
16320cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
16420cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
16520cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
16620cf1dd8SToby Isaac           for (d = 0; d < dim*dim; ++d) {
16720cf1dd8SToby Isaac             const PetscInt i = ((p*pdim + j)*Nc + c)*dim*dim + d;
16820cf1dd8SToby Isaac 
16920cf1dd8SToby Isaac             H[i] = 0.0;
17020cf1dd8SToby Isaac             for (k = 0; k < pdim; ++k) {
17120cf1dd8SToby Isaac               H[i] += fem->invV[k*pdim+j] * tmpH[((p*pdim + k)*Nc + c)*dim*dim + d];
17220cf1dd8SToby Isaac             }
17320cf1dd8SToby Isaac           }
17420cf1dd8SToby Isaac         }
17520cf1dd8SToby Isaac       }
17620cf1dd8SToby Isaac     }
17720cf1dd8SToby Isaac   }
178*ef0bb6c7SMatthew G. Knepley   if (K >= 0) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
179*ef0bb6c7SMatthew G. Knepley   if (K >= 1) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
180*ef0bb6c7SMatthew G. Knepley   if (K >= 2) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
18120cf1dd8SToby Isaac   PetscFunctionReturn(0);
18220cf1dd8SToby Isaac }
18320cf1dd8SToby Isaac 
1842b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1854bee2e38SMatthew G. Knepley                                              const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
18620cf1dd8SToby Isaac {
18720cf1dd8SToby Isaac   const PetscInt     debug = 0;
1884bee2e38SMatthew G. Knepley   PetscFE            fe;
18920cf1dd8SToby Isaac   PetscPointFunc     obj_func;
19020cf1dd8SToby Isaac   PetscQuadrature    quad;
191*ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
1924bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
19320cf1dd8SToby Isaac   const PetscScalar *constants;
19420cf1dd8SToby Isaac   PetscReal         *x;
195*ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
19620cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
19720cf1dd8SToby Isaac   PetscBool          isAffine;
19820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
19920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
20020cf1dd8SToby Isaac   PetscErrorCode     ierr;
20120cf1dd8SToby Isaac 
20220cf1dd8SToby Isaac   PetscFunctionBegin;
2034bee2e38SMatthew G. Knepley   ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr);
20420cf1dd8SToby Isaac   if (!obj_func) PetscFunctionReturn(0);
2054bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
2064bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
2074bee2e38SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2084bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
2094bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
2104bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
2114bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
212*ef0bb6c7SMatthew G. Knepley   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
2134bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr);
2144bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
2154bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
2164bee2e38SMatthew G. Knepley   if (dsAux) {
2174bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
2184bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
2194bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
2204bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
221*ef0bb6c7SMatthew G. Knepley     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
2224bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
22320cf1dd8SToby Isaac   }
22420cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
2254bee2e38SMatthew G. Knepley   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
22620cf1dd8SToby Isaac   Np = cgeom->numPoints;
22720cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
22820cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
22920cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2304bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
23120cf1dd8SToby Isaac 
23220cf1dd8SToby Isaac     if (isAffine) {
2334bee2e38SMatthew G. Knepley       fegeom.v    = x;
2344bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2354bee2e38SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*dE*dE];
2364bee2e38SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*dE*dE];
2374bee2e38SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e];
23820cf1dd8SToby Isaac     }
2394bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
2404bee2e38SMatthew G. Knepley       PetscScalar integrand;
2414bee2e38SMatthew G. Knepley       PetscReal   w;
2424bee2e38SMatthew G. Knepley 
2434bee2e38SMatthew G. Knepley       if (isAffine) {
2444bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
2454bee2e38SMatthew G. Knepley       } else {
2464bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
2474bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
2484bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
2494bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
2504bee2e38SMatthew G. Knepley       }
2514bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
25220cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
2534bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
2547be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2554bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
25620cf1dd8SToby Isaac #endif
25720cf1dd8SToby Isaac       }
25820cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
259*ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr);
260*ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
2614bee2e38SMatthew 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);
2624bee2e38SMatthew G. Knepley       integrand *= w;
26320cf1dd8SToby Isaac       integral[e*Nf+field] += integrand;
26420cf1dd8SToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);}
26520cf1dd8SToby Isaac     }
26620cf1dd8SToby Isaac     cOffset    += totDim;
26720cf1dd8SToby Isaac     cOffsetAux += totDimAux;
26820cf1dd8SToby Isaac   }
26920cf1dd8SToby Isaac   PetscFunctionReturn(0);
27020cf1dd8SToby Isaac }
27120cf1dd8SToby Isaac 
2722b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field,
273afe6d6adSToby Isaac                                                PetscBdPointFunc obj_func,
2744bee2e38SMatthew G. Knepley                                                PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
275afe6d6adSToby Isaac {
276afe6d6adSToby Isaac   const PetscInt     debug = 0;
2774bee2e38SMatthew G. Knepley   PetscFE            fe;
278afe6d6adSToby Isaac   PetscQuadrature    quad;
279*ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
2804bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
281afe6d6adSToby Isaac   const PetscScalar *constants;
282afe6d6adSToby Isaac   PetscReal         *x;
283*ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
284afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
285afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
286afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
287afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
288afe6d6adSToby Isaac   PetscErrorCode     ierr;
289afe6d6adSToby Isaac 
290afe6d6adSToby Isaac   PetscFunctionBegin;
291afe6d6adSToby Isaac   if (!obj_func) PetscFunctionReturn(0);
2924bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
2934bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
2944bee2e38SMatthew G. Knepley   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
2954bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
2964bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
2974bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
2984bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
2994bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr);
3004bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
301*ef0bb6c7SMatthew G. Knepley   ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr);
3024bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
3034bee2e38SMatthew G. Knepley   if (dsAux) {
3044bee2e38SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
3054bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
3064bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
3074bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
3084bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
3094bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
310afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
311*ef0bb6c7SMatthew G. Knepley     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
312*ef0bb6c7SMatthew G. Knepley     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
313afe6d6adSToby Isaac   }
314afe6d6adSToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
315afe6d6adSToby Isaac   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
316afe6d6adSToby Isaac   Np = fgeom->numPoints;
317afe6d6adSToby Isaac   dE = fgeom->dimEmbed;
318afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
319afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
3209f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
321afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
3224dcf62a8SSatish Balay     fegeom.n = 0;
3234dcf62a8SSatish Balay     fegeom.v = 0;
3244dcf62a8SSatish Balay     fegeom.J = 0;
3254dcf62a8SSatish Balay     fegeom.detJ = 0;
3264bee2e38SMatthew G. Knepley     if (isAffine) {
3274bee2e38SMatthew G. Knepley       fegeom.v    = x;
3284bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3294bee2e38SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*dE*dE];
3309f209ee4SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*dE*dE];
3314bee2e38SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e];
3324bee2e38SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE];
3339f209ee4SMatthew G. Knepley 
3349f209ee4SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*dE*dE];
3359f209ee4SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*dE*dE];
3369f209ee4SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e];
3374bee2e38SMatthew G. Knepley     }
338afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
339afe6d6adSToby Isaac       PetscScalar integrand;
3404bee2e38SMatthew G. Knepley       PetscReal   w;
341afe6d6adSToby Isaac 
342afe6d6adSToby Isaac       if (isAffine) {
3434bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
344afe6d6adSToby Isaac       } else {
3453fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
3469f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
3479f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
3484bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
3494bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
3509f209ee4SMatthew G. Knepley 
3519f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
3529f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
3539f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
354afe6d6adSToby Isaac       }
3554bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
356afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
3574bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
358afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3594bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
360afe6d6adSToby Isaac #endif
361afe6d6adSToby Isaac       }
362afe6d6adSToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
363*ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr);
364*ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
3654bee2e38SMatthew 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);
3664bee2e38SMatthew G. Knepley       integrand *= w;
367afe6d6adSToby Isaac       integral[e*Nf+field] += integrand;
368afe6d6adSToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);}
369afe6d6adSToby Isaac     }
370afe6d6adSToby Isaac     cOffset    += totDim;
371afe6d6adSToby Isaac     cOffsetAux += totDimAux;
372afe6d6adSToby Isaac   }
373afe6d6adSToby Isaac   PetscFunctionReturn(0);
374afe6d6adSToby Isaac }
375afe6d6adSToby Isaac 
3764bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
3774bee2e38SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
37820cf1dd8SToby Isaac {
37920cf1dd8SToby Isaac   const PetscInt     debug = 0;
3804bee2e38SMatthew G. Knepley   PetscFE            fe;
38120cf1dd8SToby Isaac   PetscPointFunc     f0_func;
38220cf1dd8SToby Isaac   PetscPointFunc     f1_func;
38320cf1dd8SToby Isaac   PetscQuadrature    quad;
384*ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
3854bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
38620cf1dd8SToby Isaac   const PetscScalar *constants;
38720cf1dd8SToby Isaac   PetscReal         *x;
388*ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
389*ef0bb6c7SMatthew G. Knepley   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e;
39020cf1dd8SToby Isaac   PetscBool          isAffine;
39120cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3924bee2e38SMatthew G. Knepley   PetscInt           qNc, Nq, q, Np, dE;
39320cf1dd8SToby Isaac   PetscErrorCode     ierr;
39420cf1dd8SToby Isaac 
39520cf1dd8SToby Isaac   PetscFunctionBegin;
3964bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
3974bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
3984bee2e38SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
3994bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
4004bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
4014bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
4024bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
4034bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
4044bee2e38SMatthew G. Knepley   ierr = PetscDSGetResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr);
4054bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
4064bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
4074bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
4084bee2e38SMatthew G. Knepley   if (!f0_func && !f1_func) PetscFunctionReturn(0);
409*ef0bb6c7SMatthew G. Knepley   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
4104bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
4114bee2e38SMatthew G. Knepley   if (dsAux) {
4124bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
4134bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
4144bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
4154bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
4164bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
417*ef0bb6c7SMatthew G. Knepley     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
41820cf1dd8SToby Isaac   }
41920cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
4204bee2e38SMatthew G. Knepley   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
42120cf1dd8SToby Isaac   Np = cgeom->numPoints;
42220cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
42320cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
42420cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4254bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
42620cf1dd8SToby Isaac 
4274bee2e38SMatthew G. Knepley     if (isAffine) {
4284bee2e38SMatthew G. Knepley       fegeom.v    = x;
4294bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
4304bee2e38SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*dE*dE];
4314bee2e38SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*dE*dE];
4324bee2e38SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e];
4334bee2e38SMatthew G. Knepley     }
434*ef0bb6c7SMatthew G. Knepley     ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr);
435*ef0bb6c7SMatthew G. Knepley     ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dim);CHKERRQ(ierr);
43620cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4374bee2e38SMatthew G. Knepley       PetscReal w;
4384bee2e38SMatthew G. Knepley       PetscInt  c, d;
43920cf1dd8SToby Isaac 
44020cf1dd8SToby Isaac       if (isAffine) {
4414bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
44220cf1dd8SToby Isaac       } else {
4434bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
4444bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
4454bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
4464bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
44720cf1dd8SToby Isaac       }
4484bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
44920cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
4504bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
4517be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4524bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
45320cf1dd8SToby Isaac #endif
45420cf1dd8SToby Isaac       }
45520cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
456*ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
457*ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
4584bee2e38SMatthew G. Knepley       if (f0_func) {
459*ef0bb6c7SMatthew G. Knepley         f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q*T[field]->Nc]);
460*ef0bb6c7SMatthew G. Knepley         for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w;
4614bee2e38SMatthew G. Knepley       }
46220cf1dd8SToby Isaac       if (f1_func) {
463*ef0bb6c7SMatthew G. Knepley         f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q*T[field]->Nc*dim]);
464*ef0bb6c7SMatthew G. Knepley         for (c = 0; c < T[field]->Nc; ++c) for (d = 0; d < dim; ++d) f1[(q*T[field]->Nc+c)*dim+d] *= w;
46520cf1dd8SToby Isaac       }
46620cf1dd8SToby Isaac     }
467*ef0bb6c7SMatthew G. Knepley     ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
46820cf1dd8SToby Isaac     cOffset    += totDim;
46920cf1dd8SToby Isaac     cOffsetAux += totDimAux;
47020cf1dd8SToby Isaac   }
47120cf1dd8SToby Isaac   PetscFunctionReturn(0);
47220cf1dd8SToby Isaac }
47320cf1dd8SToby Isaac 
4744bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
4754bee2e38SMatthew G. Knepley                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
47620cf1dd8SToby Isaac {
47720cf1dd8SToby Isaac   const PetscInt     debug = 0;
4784bee2e38SMatthew G. Knepley   PetscFE            fe;
47920cf1dd8SToby Isaac   PetscBdPointFunc   f0_func;
48020cf1dd8SToby Isaac   PetscBdPointFunc   f1_func;
48120cf1dd8SToby Isaac   PetscQuadrature    quad;
482*ef0bb6c7SMatthew G. Knepley   PetscTabulation   *Tf, *TfAux = NULL;
4834bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
48420cf1dd8SToby Isaac   const PetscScalar *constants;
48520cf1dd8SToby Isaac   PetscReal         *x;
486*ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
487*ef0bb6c7SMatthew G. Knepley   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI;
488b07bd890SMatthew G. Knepley   PetscBool          isAffine, auxOnBd = PETSC_FALSE;
48920cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
49020cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
49120cf1dd8SToby Isaac   PetscErrorCode     ierr;
49220cf1dd8SToby Isaac 
49320cf1dd8SToby Isaac   PetscFunctionBegin;
4944bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
4954bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
4964bee2e38SMatthew G. Knepley   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
4974bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
4984bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
4994bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
5004bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
5014bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
5024bee2e38SMatthew G. Knepley   ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr);
50320cf1dd8SToby Isaac   if (!f0_func && !f1_func) PetscFunctionReturn(0);
5044bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
5054bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
5064bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
507*ef0bb6c7SMatthew G. Knepley   ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr);
5084bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
5094bee2e38SMatthew G. Knepley   if (dsAux) {
5104bee2e38SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
5114bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
5124bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
5134bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
5144bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
5154bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
5167be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
517*ef0bb6c7SMatthew G. Knepley     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
518*ef0bb6c7SMatthew G. Knepley     else         {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);}
51920cf1dd8SToby Isaac   }
520*ef0bb6c7SMatthew G. Knepley   NcI = Tf[field]->Nc;
52120cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
522afe6d6adSToby Isaac   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
52320cf1dd8SToby Isaac   Np = fgeom->numPoints;
52420cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
52520cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
52620cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5279f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
52820cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
5294dcf62a8SSatish Balay     fegeom.n = 0;
5304dcf62a8SSatish Balay     fegeom.v = 0;
5314dcf62a8SSatish Balay     fegeom.J = 0;
5324dcf62a8SSatish Balay     fegeom.detJ = 0;
5334bee2e38SMatthew G. Knepley     if (isAffine) {
5344bee2e38SMatthew G. Knepley       fegeom.v    = x;
5354bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
5364bee2e38SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*dE*dE];
5379f209ee4SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*dE*dE];
5384bee2e38SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e];
5394bee2e38SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE];
5409f209ee4SMatthew G. Knepley 
5419f209ee4SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*dE*dE];
5429f209ee4SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*dE*dE];
5439f209ee4SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e];
5444bee2e38SMatthew G. Knepley     }
545580bdb30SBarry Smith     ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr);
546580bdb30SBarry Smith     ierr = PetscArrayzero(f1, Nq*NcI*dim);CHKERRQ(ierr);
54720cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5484bee2e38SMatthew G. Knepley       PetscReal w;
5494bee2e38SMatthew G. Knepley       PetscInt  c, d;
5504bee2e38SMatthew G. Knepley 
55120cf1dd8SToby Isaac       if (isAffine) {
5524bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
55320cf1dd8SToby Isaac       } else {
5543fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
5559f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
5569f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
5574bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
5584bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
5599f209ee4SMatthew G. Knepley 
5609f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
5619f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
5629f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
56320cf1dd8SToby Isaac       }
5644bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
56520cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
5664bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
5677be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5684bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
56920cf1dd8SToby Isaac #endif
57020cf1dd8SToby Isaac       }
57120cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
572*ef0bb6c7SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
573*ef0bb6c7SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
5744bee2e38SMatthew G. Knepley       if (f0_func) {
5754bee2e38SMatthew G. Knepley         f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q*NcI]);
5764bee2e38SMatthew G. Knepley         for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w;
5774bee2e38SMatthew G. Knepley       }
57820cf1dd8SToby Isaac       if (f1_func) {
5794bee2e38SMatthew G. Knepley         f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q*NcI*dim]);
5804bee2e38SMatthew G. Knepley         for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w;
58120cf1dd8SToby Isaac       }
58220cf1dd8SToby Isaac     }
583*ef0bb6c7SMatthew G. Knepley     ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, &cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
58420cf1dd8SToby Isaac     cOffset    += totDim;
58520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
58620cf1dd8SToby Isaac   }
58720cf1dd8SToby Isaac   PetscFunctionReturn(0);
58820cf1dd8SToby Isaac }
58920cf1dd8SToby Isaac 
5904bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
5914bee2e38SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
59220cf1dd8SToby Isaac {
59320cf1dd8SToby Isaac   const PetscInt     debug      = 0;
5944bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
5954bee2e38SMatthew G. Knepley   PetscPointJac      g0_func, g1_func, g2_func, g3_func;
59620cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
59720cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
59820cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
59920cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
60020cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
60120cf1dd8SToby Isaac   PetscQuadrature    quad;
602*ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
6034bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
60420cf1dd8SToby Isaac   const PetscScalar *constants;
60520cf1dd8SToby Isaac   PetscReal         *x;
606*ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
607*ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
60820cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
60920cf1dd8SToby Isaac   PetscInt           dE, Np;
61020cf1dd8SToby Isaac   PetscBool          isAffine;
61120cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
61220cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
61320cf1dd8SToby Isaac   PetscErrorCode     ierr;
61420cf1dd8SToby Isaac 
61520cf1dd8SToby Isaac   PetscFunctionBegin;
6164bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
6174bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
6184bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
6194bee2e38SMatthew G. Knepley   ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr);
6204bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
6214bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
6224bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
6234bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
62420cf1dd8SToby Isaac   switch(jtype) {
6254bee2e38SMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
6264bee2e38SMatthew G. Knepley   case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
6274bee2e38SMatthew G. Knepley   case PETSCFE_JACOBIAN:     ierr = PetscDSGetJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
62820cf1dd8SToby Isaac   }
62920cf1dd8SToby Isaac   if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0);
6304bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
6314bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
6324bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
633*ef0bb6c7SMatthew G. Knepley   ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr);
6344bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
6354bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
6364bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
6374bee2e38SMatthew G. Knepley   if (dsAux) {
6384bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
6394bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
6404bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
6414bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
6424bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
643*ef0bb6c7SMatthew G. Knepley     ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);
64420cf1dd8SToby Isaac   }
645*ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
64620cf1dd8SToby Isaac   /* Initialize here in case the function is not defined */
647580bdb30SBarry Smith   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
648580bdb30SBarry Smith   ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr);
649580bdb30SBarry Smith   ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr);
650580bdb30SBarry Smith   ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr);
65120cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
65220cf1dd8SToby Isaac   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
6534bee2e38SMatthew G. Knepley   Np = cgeom->numPoints;
6544bee2e38SMatthew G. Knepley   dE = cgeom->dimEmbed;
6554bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
6564bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
6574bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
6584bee2e38SMatthew G. Knepley 
6594bee2e38SMatthew G. Knepley     if (isAffine) {
6604bee2e38SMatthew G. Knepley       fegeom.v    = x;
6614bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
6624bee2e38SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*dE*dE];
6634bee2e38SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*dE*dE];
6644bee2e38SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e];
6654bee2e38SMatthew G. Knepley     }
66620cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
66720cf1dd8SToby Isaac       PetscReal w;
6684bee2e38SMatthew G. Knepley       PetscInt  c;
66920cf1dd8SToby Isaac 
67020cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
67120cf1dd8SToby Isaac       if (isAffine) {
6724bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
67320cf1dd8SToby Isaac       } else {
6744bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
6754bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
6764bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
6774bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
67820cf1dd8SToby Isaac       }
6794bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
680*ef0bb6c7SMatthew G. Knepley       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
681*ef0bb6c7SMatthew G. Knepley       if (dsAux)      {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
68220cf1dd8SToby Isaac       if (g0_func) {
683580bdb30SBarry Smith         ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
6844bee2e38SMatthew G. Knepley         g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g0);
68520cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
68620cf1dd8SToby Isaac       }
68720cf1dd8SToby Isaac       if (g1_func) {
688580bdb30SBarry Smith         ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr);
6894bee2e38SMatthew G. Knepley         g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g1);
6904bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
69120cf1dd8SToby Isaac       }
69220cf1dd8SToby Isaac       if (g2_func) {
693580bdb30SBarry Smith         ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr);
6944bee2e38SMatthew G. Knepley         g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g2);
6954bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
69620cf1dd8SToby Isaac       }
69720cf1dd8SToby Isaac       if (g3_func) {
698580bdb30SBarry Smith         ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr);
6994bee2e38SMatthew G. Knepley         g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g3);
7004bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
70120cf1dd8SToby Isaac       }
70220cf1dd8SToby Isaac 
703*ef0bb6c7SMatthew G. Knepley       ierr = PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr);
70420cf1dd8SToby Isaac     }
70520cf1dd8SToby Isaac     if (debug > 1) {
70620cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
70720cf1dd8SToby Isaac 
70820cf1dd8SToby Isaac       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
709*ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
710*ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
711*ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
712*ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
713*ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
714*ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
71520cf1dd8SToby 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);
71620cf1dd8SToby Isaac             }
71720cf1dd8SToby Isaac           }
71820cf1dd8SToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
71920cf1dd8SToby Isaac         }
72020cf1dd8SToby Isaac       }
72120cf1dd8SToby Isaac     }
72220cf1dd8SToby Isaac     cOffset    += totDim;
72320cf1dd8SToby Isaac     cOffsetAux += totDimAux;
72420cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
72520cf1dd8SToby Isaac   }
72620cf1dd8SToby Isaac   PetscFunctionReturn(0);
72720cf1dd8SToby Isaac }
72820cf1dd8SToby Isaac 
7292b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
7304bee2e38SMatthew G. Knepley                                                        const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
73120cf1dd8SToby Isaac {
73220cf1dd8SToby Isaac   const PetscInt     debug      = 0;
7334bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
7344bee2e38SMatthew G. Knepley   PetscBdPointJac    g0_func, g1_func, g2_func, g3_func;
73520cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
73620cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
73720cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
73820cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
73920cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
74020cf1dd8SToby Isaac   PetscQuadrature    quad;
741*ef0bb6c7SMatthew G. Knepley   PetscTabulation   *T, *TAux = NULL;
7424bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
74320cf1dd8SToby Isaac   const PetscScalar *constants;
74420cf1dd8SToby Isaac   PetscReal         *x;
745*ef0bb6c7SMatthew G. Knepley   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL;
746*ef0bb6c7SMatthew G. Knepley   PetscInt           NcI = 0, NcJ = 0;
74720cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
74820cf1dd8SToby Isaac   PetscBool          isAffine;
74920cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
75020cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
75120cf1dd8SToby Isaac   PetscErrorCode     ierr;
75220cf1dd8SToby Isaac 
75320cf1dd8SToby Isaac   PetscFunctionBegin;
7544bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
7554bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
7564bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
7574bee2e38SMatthew G. Knepley   ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr);
7584bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
7594bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
7604bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
7614bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
7624bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
7634bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
7644bee2e38SMatthew G. Knepley   ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);
76531f073d7SMatthew G. Knepley   if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0);
7664bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
7674bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
7684bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
769*ef0bb6c7SMatthew G. Knepley   ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr);
7704bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
7714bee2e38SMatthew G. Knepley   if (dsAux) {
7724bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
7734bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
7744bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
7754bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
7764bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
777*ef0bb6c7SMatthew G. Knepley     ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);
77820cf1dd8SToby Isaac   }
779*ef0bb6c7SMatthew G. Knepley   NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc;
78020cf1dd8SToby Isaac   /* Initialize here in case the function is not defined */
781580bdb30SBarry Smith   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
782580bdb30SBarry Smith   ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr);
783580bdb30SBarry Smith   ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr);
784580bdb30SBarry Smith   ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr);
78520cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
7864bee2e38SMatthew G. Knepley   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
78720cf1dd8SToby Isaac   Np = fgeom->numPoints;
78820cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
78920cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
79031f073d7SMatthew G. Knepley   ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
79131f073d7SMatthew G. Knepley   ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr);
79231f073d7SMatthew G. Knepley   ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr);
79331f073d7SMatthew G. Knepley   ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr);
79420cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
7959f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
79620cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
7974dcf62a8SSatish Balay     fegeom.n = 0;
7984dcf62a8SSatish Balay     fegeom.v = 0;
7994dcf62a8SSatish Balay     fegeom.J = 0;
8004dcf62a8SSatish Balay     fegeom.detJ = 0;
8014bee2e38SMatthew G. Knepley     if (isAffine) {
8024bee2e38SMatthew G. Knepley       fegeom.v    = x;
8034bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
8044bee2e38SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*dE*dE];
8059f209ee4SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*dE*dE];
8064bee2e38SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e];
8074bee2e38SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE];
8089f209ee4SMatthew G. Knepley 
8099f209ee4SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*dE*dE];
8109f209ee4SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*dE*dE];
8119f209ee4SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e];
8124bee2e38SMatthew G. Knepley     }
81320cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
81420cf1dd8SToby Isaac       PetscReal w;
8154bee2e38SMatthew G. Knepley       PetscInt  c;
81620cf1dd8SToby Isaac 
81720cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
81820cf1dd8SToby Isaac       if (isAffine) {
8194bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
82020cf1dd8SToby Isaac       } else {
8213fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
8229f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
8239f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
8244bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
8254bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
8269f209ee4SMatthew G. Knepley 
8279f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
8289f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
8299f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
83020cf1dd8SToby Isaac       }
8314bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
832*ef0bb6c7SMatthew G. Knepley       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
833*ef0bb6c7SMatthew G. Knepley       if (dsAux)      {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
83420cf1dd8SToby Isaac       if (g0_func) {
835580bdb30SBarry Smith         ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr);
8364bee2e38SMatthew G. Knepley         g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0);
83720cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
83820cf1dd8SToby Isaac       }
83920cf1dd8SToby Isaac       if (g1_func) {
840580bdb30SBarry Smith         ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr);
8414bee2e38SMatthew G. Knepley         g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1);
8424bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
84320cf1dd8SToby Isaac       }
84420cf1dd8SToby Isaac       if (g2_func) {
845580bdb30SBarry Smith         ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr);
8464bee2e38SMatthew G. Knepley         g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2);
8474bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
84820cf1dd8SToby Isaac       }
84920cf1dd8SToby Isaac       if (g3_func) {
850580bdb30SBarry Smith         ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr);
8514bee2e38SMatthew G. Knepley         g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3);
8524bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
85320cf1dd8SToby Isaac       }
85420cf1dd8SToby Isaac 
855*ef0bb6c7SMatthew G. Knepley       ierr = PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr);
85620cf1dd8SToby Isaac     }
85720cf1dd8SToby Isaac     if (debug > 1) {
85820cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
85920cf1dd8SToby Isaac 
86020cf1dd8SToby Isaac       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
861*ef0bb6c7SMatthew G. Knepley       for (fc = 0; fc < T[fieldI]->Nc; ++fc) {
862*ef0bb6c7SMatthew G. Knepley         for (f = 0; f < T[fieldI]->Nb; ++f) {
863*ef0bb6c7SMatthew G. Knepley           const PetscInt i = offsetI + f*T[fieldI]->Nc+fc;
864*ef0bb6c7SMatthew G. Knepley           for (gc = 0; gc < T[fieldJ]->Nc; ++gc) {
865*ef0bb6c7SMatthew G. Knepley             for (g = 0; g < T[fieldJ]->Nb; ++g) {
866*ef0bb6c7SMatthew G. Knepley               const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc;
86720cf1dd8SToby 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);
86820cf1dd8SToby Isaac             }
86920cf1dd8SToby Isaac           }
87020cf1dd8SToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
87120cf1dd8SToby Isaac         }
87220cf1dd8SToby Isaac       }
87320cf1dd8SToby Isaac     }
87420cf1dd8SToby Isaac     cOffset    += totDim;
87520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
87620cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
87720cf1dd8SToby Isaac   }
87820cf1dd8SToby Isaac   PetscFunctionReturn(0);
87920cf1dd8SToby Isaac }
88020cf1dd8SToby Isaac 
8812b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
88220cf1dd8SToby Isaac {
88320cf1dd8SToby Isaac   PetscFunctionBegin;
88420cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
88520cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
88620cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
88720cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
88820cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
889*ef0bb6c7SMatthew G. Knepley   fem->ops->createtabulation        = PetscFECreateTabulation_Basic;
89020cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
891afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
89220cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
89320cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
89420cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
89520cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
89620cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
89720cf1dd8SToby Isaac   PetscFunctionReturn(0);
89820cf1dd8SToby Isaac }
89920cf1dd8SToby Isaac 
90020cf1dd8SToby Isaac /*MC
90120cf1dd8SToby Isaac   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
90220cf1dd8SToby Isaac 
90320cf1dd8SToby Isaac   Level: intermediate
90420cf1dd8SToby Isaac 
90520cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
90620cf1dd8SToby Isaac M*/
90720cf1dd8SToby Isaac 
90820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
90920cf1dd8SToby Isaac {
91020cf1dd8SToby Isaac   PetscFE_Basic *b;
91120cf1dd8SToby Isaac   PetscErrorCode ierr;
91220cf1dd8SToby Isaac 
91320cf1dd8SToby Isaac   PetscFunctionBegin;
91420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
91520cf1dd8SToby Isaac   ierr      = PetscNewLog(fem,&b);CHKERRQ(ierr);
91620cf1dd8SToby Isaac   fem->data = b;
91720cf1dd8SToby Isaac 
91820cf1dd8SToby Isaac   ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr);
91920cf1dd8SToby Isaac   PetscFunctionReturn(0);
92020cf1dd8SToby Isaac }
921