xref: /petsc/src/dm/dt/fe/impls/basic/febasic.c (revision 4dcf62a8fdbe4bbf3ef5c1ae035a1ab8dc216862)
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 
14d9bac1caSLisandro Dalcin 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 
37d9bac1caSLisandro Dalcin 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 */
4920cf1dd8SToby Isaac 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   }
8520cf1dd8SToby Isaac   ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr);
8620cf1dd8SToby Isaac   n = pdim;
8720cf1dd8SToby Isaac   PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invVscalar, &n, pivots, &info));
8820cf1dd8SToby Isaac   PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invVscalar, &n, pivots, work, &n, &info));
8920cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX)
9020cf1dd8SToby Isaac   for (j = 0; j < pdim*pdim; j++) fem->invV[j] = PetscRealPart(invVscalar[j]);
9120cf1dd8SToby Isaac   ierr = PetscFree(invVscalar);CHKERRQ(ierr);
9220cf1dd8SToby Isaac #endif
9320cf1dd8SToby Isaac   ierr = PetscFree2(pivots,work);CHKERRQ(ierr);
9420cf1dd8SToby Isaac   PetscFunctionReturn(0);
9520cf1dd8SToby Isaac }
9620cf1dd8SToby Isaac 
9720cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
9820cf1dd8SToby Isaac {
9920cf1dd8SToby Isaac   PetscErrorCode ierr;
10020cf1dd8SToby Isaac 
10120cf1dd8SToby Isaac   PetscFunctionBegin;
10220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr);
10320cf1dd8SToby Isaac   PetscFunctionReturn(0);
10420cf1dd8SToby Isaac }
10520cf1dd8SToby Isaac 
10620cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H)
10720cf1dd8SToby Isaac {
10820cf1dd8SToby Isaac   DM               dm;
10920cf1dd8SToby Isaac   PetscInt         pdim; /* Dimension of FE space P */
11020cf1dd8SToby Isaac   PetscInt         dim;  /* Spatial dimension */
11120cf1dd8SToby Isaac   PetscInt         Nc;   /* Field components */
11220cf1dd8SToby Isaac   PetscReal       *tmpB, *tmpD, *tmpH;
11320cf1dd8SToby Isaac   PetscInt         p, d, j, k, c;
11420cf1dd8SToby Isaac   PetscErrorCode   ierr;
11520cf1dd8SToby Isaac 
11620cf1dd8SToby Isaac   PetscFunctionBegin;
11720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
11820cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
11920cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
12020cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr);
12120cf1dd8SToby Isaac   /* Evaluate the prime basis functions at all points */
12220cf1dd8SToby Isaac   if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
12320cf1dd8SToby Isaac   if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
12420cf1dd8SToby Isaac   if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
12520cf1dd8SToby Isaac   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);CHKERRQ(ierr);
12620cf1dd8SToby Isaac   /* Translate to the nodal basis */
12720cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
12820cf1dd8SToby Isaac     if (B) {
12920cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
13020cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
13120cf1dd8SToby Isaac         const PetscInt i = (p*pdim + j)*Nc;
13220cf1dd8SToby Isaac 
13320cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
13420cf1dd8SToby Isaac           B[i+c] = 0.0;
13520cf1dd8SToby Isaac           for (k = 0; k < pdim; ++k) {
13620cf1dd8SToby Isaac             B[i+c] += fem->invV[k*pdim+j] * tmpB[(p*pdim + k)*Nc+c];
13720cf1dd8SToby Isaac           }
13820cf1dd8SToby Isaac         }
13920cf1dd8SToby Isaac       }
14020cf1dd8SToby Isaac     }
14120cf1dd8SToby Isaac     if (D) {
14220cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
14320cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
14420cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
14520cf1dd8SToby Isaac           for (d = 0; d < dim; ++d) {
14620cf1dd8SToby Isaac             const PetscInt i = ((p*pdim + j)*Nc + c)*dim + d;
14720cf1dd8SToby Isaac 
14820cf1dd8SToby Isaac             D[i] = 0.0;
14920cf1dd8SToby Isaac             for (k = 0; k < pdim; ++k) {
15020cf1dd8SToby Isaac               D[i] += fem->invV[k*pdim+j] * tmpD[((p*pdim + k)*Nc + c)*dim + d];
15120cf1dd8SToby Isaac             }
15220cf1dd8SToby Isaac           }
15320cf1dd8SToby Isaac         }
15420cf1dd8SToby Isaac       }
15520cf1dd8SToby Isaac     }
15620cf1dd8SToby Isaac     if (H) {
15720cf1dd8SToby Isaac       /* Multiply by V^{-1} (pdim x pdim) */
15820cf1dd8SToby Isaac       for (j = 0; j < pdim; ++j) {
15920cf1dd8SToby Isaac         for (c = 0; c < Nc; ++c) {
16020cf1dd8SToby Isaac           for (d = 0; d < dim*dim; ++d) {
16120cf1dd8SToby Isaac             const PetscInt i = ((p*pdim + j)*Nc + c)*dim*dim + d;
16220cf1dd8SToby Isaac 
16320cf1dd8SToby Isaac             H[i] = 0.0;
16420cf1dd8SToby Isaac             for (k = 0; k < pdim; ++k) {
16520cf1dd8SToby Isaac               H[i] += fem->invV[k*pdim+j] * tmpH[((p*pdim + k)*Nc + c)*dim*dim + d];
16620cf1dd8SToby Isaac             }
16720cf1dd8SToby Isaac           }
16820cf1dd8SToby Isaac         }
16920cf1dd8SToby Isaac       }
17020cf1dd8SToby Isaac     }
17120cf1dd8SToby Isaac   }
17220cf1dd8SToby Isaac   if (B) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);}
17320cf1dd8SToby Isaac   if (D) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);}
17420cf1dd8SToby Isaac   if (H) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);}
17520cf1dd8SToby Isaac   PetscFunctionReturn(0);
17620cf1dd8SToby Isaac }
17720cf1dd8SToby Isaac 
1784bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
1794bee2e38SMatthew G. Knepley                                       const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
18020cf1dd8SToby Isaac {
18120cf1dd8SToby Isaac   const PetscInt     debug = 0;
1824bee2e38SMatthew G. Knepley   PetscFE            fe;
18320cf1dd8SToby Isaac   PetscPointFunc     obj_func;
18420cf1dd8SToby Isaac   PetscQuadrature    quad;
1854bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x;
18620cf1dd8SToby Isaac   const PetscScalar *constants;
18720cf1dd8SToby Isaac   PetscReal         *x;
18820cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL;
18920cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
19020cf1dd8SToby Isaac   PetscInt           dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
19120cf1dd8SToby Isaac   PetscBool          isAffine;
19220cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
19320cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
19420cf1dd8SToby Isaac   PetscErrorCode     ierr;
19520cf1dd8SToby Isaac 
19620cf1dd8SToby Isaac   PetscFunctionBegin;
1974bee2e38SMatthew G. Knepley   ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr);
19820cf1dd8SToby Isaac   if (!obj_func) PetscFunctionReturn(0);
1994bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
2004bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
2014bee2e38SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
2024bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
2034bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
2044bee2e38SMatthew G. Knepley   ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr);
2054bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
2064bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
2074bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
2084bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr);
2094bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
2104bee2e38SMatthew G. Knepley   ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr);
2114bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
2124bee2e38SMatthew G. Knepley   if (dsAux) {
2134bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
2144bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
2154bee2e38SMatthew G. Knepley     ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr);
2164bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr);
2174bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
2184bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
2194bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
2204bee2e38SMatthew G. Knepley     ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);
22120cf1dd8SToby Isaac   }
22220cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
2234bee2e38SMatthew G. Knepley   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
22420cf1dd8SToby Isaac   Np = cgeom->numPoints;
22520cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
22620cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
22720cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
2284bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
22920cf1dd8SToby Isaac 
23020cf1dd8SToby Isaac     if (isAffine) {
2314bee2e38SMatthew G. Knepley       fegeom.v    = x;
2324bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
2334bee2e38SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*dE*dE];
2344bee2e38SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*dE*dE];
2354bee2e38SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e];
23620cf1dd8SToby Isaac     }
2374bee2e38SMatthew G. Knepley     for (q = 0; q < Nq; ++q) {
2384bee2e38SMatthew G. Knepley       PetscScalar integrand;
2394bee2e38SMatthew G. Knepley       PetscReal   w;
2404bee2e38SMatthew G. Knepley 
2414bee2e38SMatthew G. Knepley       if (isAffine) {
2424bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
2434bee2e38SMatthew G. Knepley       } else {
2444bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
2454bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
2464bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
2474bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
2484bee2e38SMatthew G. Knepley       }
2494bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
25020cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
2514bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
2527be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
2534bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
25420cf1dd8SToby Isaac #endif
25520cf1dd8SToby Isaac       }
25620cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
257a8f1f9e5SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr);
258a8f1f9e5SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dim, NfAux, NbAux, NcAux, q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
2594bee2e38SMatthew 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);
2604bee2e38SMatthew G. Knepley       integrand *= w;
26120cf1dd8SToby Isaac       integral[e*Nf+field] += integrand;
26220cf1dd8SToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);}
26320cf1dd8SToby Isaac     }
26420cf1dd8SToby Isaac     cOffset    += totDim;
26520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
26620cf1dd8SToby Isaac   }
26720cf1dd8SToby Isaac   PetscFunctionReturn(0);
26820cf1dd8SToby Isaac }
26920cf1dd8SToby Isaac 
2704bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field,
271afe6d6adSToby Isaac                                         PetscBdPointFunc obj_func,
2724bee2e38SMatthew G. Knepley                                         PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
273afe6d6adSToby Isaac {
274afe6d6adSToby Isaac   const PetscInt     debug = 0;
2754bee2e38SMatthew G. Knepley   PetscFE            fe;
276afe6d6adSToby Isaac   PetscQuadrature    quad;
2774bee2e38SMatthew G. Knepley   PetscScalar       *u, *u_x, *a, *a_x, *basisReal, *basisDerReal;
278afe6d6adSToby Isaac   const PetscScalar *constants;
279afe6d6adSToby Isaac   PetscReal         *x;
280afe6d6adSToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL;
281afe6d6adSToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
282afe6d6adSToby Isaac   PetscBool          isAffine, auxOnBd;
283afe6d6adSToby Isaac   const PetscReal   *quadPoints, *quadWeights;
284afe6d6adSToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
285afe6d6adSToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
286afe6d6adSToby Isaac   PetscErrorCode     ierr;
287afe6d6adSToby Isaac 
288afe6d6adSToby Isaac   PetscFunctionBegin;
289afe6d6adSToby Isaac   if (!obj_func) PetscFunctionReturn(0);
2904bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
2914bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
2924bee2e38SMatthew G. Knepley   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
2934bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
2944bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
2954bee2e38SMatthew G. Knepley   ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr);
2964bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);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);
3014bee2e38SMatthew G. Knepley   ierr = PetscDSGetFaceTabulation(ds, &B, &D);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 = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr);
3084bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr);
3094bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
3104bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
3114bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
312afe6d6adSToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
3134bee2e38SMatthew G. Knepley     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);}
3144bee2e38SMatthew G. Knepley     else         {ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);}
315afe6d6adSToby Isaac   }
316afe6d6adSToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
317afe6d6adSToby Isaac   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
318afe6d6adSToby Isaac   Np = fgeom->numPoints;
319afe6d6adSToby Isaac   dE = fgeom->dimEmbed;
320afe6d6adSToby Isaac   isAffine = fgeom->isAffine;
321afe6d6adSToby Isaac   for (e = 0; e < Ne; ++e) {
3229f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
323afe6d6adSToby Isaac     const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */
324*4dcf62a8SSatish Balay     fegeom.n = 0;
325*4dcf62a8SSatish Balay     fegeom.v = 0;
326*4dcf62a8SSatish Balay     fegeom.J = 0;
327*4dcf62a8SSatish Balay     fegeom.detJ = 0;
3284bee2e38SMatthew G. Knepley     if (isAffine) {
3294bee2e38SMatthew G. Knepley       fegeom.v    = x;
3304bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
3314bee2e38SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*dE*dE];
3329f209ee4SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*dE*dE];
3334bee2e38SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e];
3344bee2e38SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE];
3359f209ee4SMatthew G. Knepley 
3369f209ee4SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*dE*dE];
3379f209ee4SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*dE*dE];
3389f209ee4SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e];
3394bee2e38SMatthew G. Knepley     }
340afe6d6adSToby Isaac     for (q = 0; q < Nq; ++q) {
341afe6d6adSToby Isaac       PetscScalar integrand;
3424bee2e38SMatthew G. Knepley       PetscReal   w;
343afe6d6adSToby Isaac 
344afe6d6adSToby Isaac       if (isAffine) {
3454bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
346afe6d6adSToby Isaac       } else {
3473fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
3489f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
3499f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
3504bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
3514bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
3529f209ee4SMatthew G. Knepley 
3539f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
3549f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
3559f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
356afe6d6adSToby Isaac       }
3574bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
358afe6d6adSToby Isaac       if (debug > 1 && q < Np) {
3594bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
360afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX
3614bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
362afe6d6adSToby Isaac #endif
363afe6d6adSToby Isaac       }
364afe6d6adSToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
3659f209ee4SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, face*Nq+q, B, D, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr);
3669f209ee4SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dimAux, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
3674bee2e38SMatthew 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);
3684bee2e38SMatthew G. Knepley       integrand *= w;
369afe6d6adSToby Isaac       integral[e*Nf+field] += integrand;
370afe6d6adSToby Isaac       if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, "    int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);}
371afe6d6adSToby Isaac     }
372afe6d6adSToby Isaac     cOffset    += totDim;
373afe6d6adSToby Isaac     cOffsetAux += totDimAux;
374afe6d6adSToby Isaac   }
375afe6d6adSToby Isaac   PetscFunctionReturn(0);
376afe6d6adSToby Isaac }
377afe6d6adSToby Isaac 
3784bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
3794bee2e38SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
38020cf1dd8SToby Isaac {
38120cf1dd8SToby Isaac   const PetscInt     debug = 0;
3824bee2e38SMatthew G. Knepley   PetscFE            fe;
38320cf1dd8SToby Isaac   PetscPointFunc     f0_func;
38420cf1dd8SToby Isaac   PetscPointFunc     f1_func;
38520cf1dd8SToby Isaac   PetscQuadrature    quad;
3864bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
38720cf1dd8SToby Isaac   const PetscScalar *constants;
38820cf1dd8SToby Isaac   PetscReal         *x;
38920cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
39020cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
39120cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
39220cf1dd8SToby Isaac   PetscBool          isAffine;
39320cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
3944bee2e38SMatthew G. Knepley   PetscInt           qNc, Nq, q, Np, dE;
39520cf1dd8SToby Isaac   PetscErrorCode     ierr;
39620cf1dd8SToby Isaac 
39720cf1dd8SToby Isaac   PetscFunctionBegin;
3984bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
3994bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
4004bee2e38SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
4014bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
4024bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
4034bee2e38SMatthew G. Knepley   ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr);
4044bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
4054bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
4064bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
4074bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
4084bee2e38SMatthew G. Knepley   ierr = PetscDSGetResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr);
4094bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
4104bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
4114bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
4124bee2e38SMatthew G. Knepley   if (!f0_func && !f1_func) PetscFunctionReturn(0);
4134bee2e38SMatthew G. Knepley   ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr);
4144bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
4154bee2e38SMatthew G. Knepley   if (dsAux) {
4164bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
4174bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
4184bee2e38SMatthew G. Knepley     ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr);
4194bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr);
4204bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
4214bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
4224bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
4234bee2e38SMatthew G. Knepley     ierr = PetscDSGetTabulation(dsAux, &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);
4304bee2e38SMatthew G. Knepley   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
43120cf1dd8SToby Isaac   Np = cgeom->numPoints;
43220cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
43320cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
43420cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
4354bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
43620cf1dd8SToby Isaac 
4374bee2e38SMatthew G. Knepley     if (isAffine) {
4384bee2e38SMatthew G. Knepley       fegeom.v    = x;
4394bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
4404bee2e38SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*dE*dE];
4414bee2e38SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*dE*dE];
4424bee2e38SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e];
4434bee2e38SMatthew G. Knepley     }
44420cf1dd8SToby Isaac     ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr);
44520cf1dd8SToby Isaac     ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr);
44620cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
4474bee2e38SMatthew G. Knepley       PetscReal w;
4484bee2e38SMatthew G. Knepley       PetscInt  c, d;
44920cf1dd8SToby Isaac 
45020cf1dd8SToby Isaac       if (isAffine) {
4514bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
45220cf1dd8SToby Isaac       } else {
4534bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
4544bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
4554bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
4564bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
45720cf1dd8SToby Isaac       }
4584bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
45920cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
4604bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
4617be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
4624bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
46320cf1dd8SToby Isaac #endif
46420cf1dd8SToby Isaac       }
46520cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
466a8f1f9e5SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
467a8f1f9e5SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dim, NfAux, NbAux, NcAux, q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
4684bee2e38SMatthew G. Knepley       if (f0_func) {
4694bee2e38SMatthew 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*NcI]);
4704bee2e38SMatthew G. Knepley         for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w;
4714bee2e38SMatthew G. Knepley       }
47220cf1dd8SToby Isaac       if (f1_func) {
4734bee2e38SMatthew 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*NcI*dim]);
4744bee2e38SMatthew G. Knepley         for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w;
47520cf1dd8SToby Isaac       }
47620cf1dd8SToby Isaac     }
477a8f1f9e5SMatthew G. Knepley     ierr = PetscFEUpdateElementVec_Internal(fe, dim, Nq, NbI, NcI, BI, DI, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
47820cf1dd8SToby Isaac     cOffset    += totDim;
47920cf1dd8SToby Isaac     cOffsetAux += totDimAux;
48020cf1dd8SToby Isaac   }
48120cf1dd8SToby Isaac   PetscFunctionReturn(0);
48220cf1dd8SToby Isaac }
48320cf1dd8SToby Isaac 
4844bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
4854bee2e38SMatthew G. Knepley                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
48620cf1dd8SToby Isaac {
48720cf1dd8SToby Isaac   const PetscInt     debug = 0;
4884bee2e38SMatthew G. Knepley   PetscFE            fe;
48920cf1dd8SToby Isaac   PetscBdPointFunc   f0_func;
49020cf1dd8SToby Isaac   PetscBdPointFunc   f1_func;
49120cf1dd8SToby Isaac   PetscQuadrature    quad;
4924bee2e38SMatthew G. Knepley   PetscScalar       *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal;
49320cf1dd8SToby Isaac   const PetscScalar *constants;
49420cf1dd8SToby Isaac   PetscReal         *x;
49520cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
49620cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4977be5e748SToby Isaac   PetscInt           dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
498b07bd890SMatthew G. Knepley   PetscBool          isAffine, auxOnBd = PETSC_FALSE;
49920cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
50020cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
50120cf1dd8SToby Isaac   PetscErrorCode     ierr;
50220cf1dd8SToby Isaac 
50320cf1dd8SToby Isaac   PetscFunctionBegin;
5044bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr);
5054bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr);
5064bee2e38SMatthew G. Knepley   ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr);
5074bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
5084bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
5094bee2e38SMatthew G. Knepley   ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr);
5104bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
5114bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
5124bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
5134bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr);
5144bee2e38SMatthew G. Knepley   ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr);
51520cf1dd8SToby Isaac   if (!f0_func && !f1_func) PetscFunctionReturn(0);
5164bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
5174bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr);
5184bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
5194bee2e38SMatthew G. Knepley   ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr);
5204bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
5214bee2e38SMatthew G. Knepley   if (dsAux) {
5224bee2e38SMatthew G. Knepley     ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr);
5234bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
5244bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
5254bee2e38SMatthew G. Knepley     ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr);
5264bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr);
5274bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
5284bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
5294bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
5307be5e748SToby Isaac     auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE;
5314bee2e38SMatthew G. Knepley     if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);}
5324bee2e38SMatthew G. Knepley     else         {ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);}
53320cf1dd8SToby Isaac   }
53420cf1dd8SToby Isaac   NbI = Nb[field];
53520cf1dd8SToby Isaac   NcI = Nc[field];
53620cf1dd8SToby Isaac   BI  = B[field];
53720cf1dd8SToby Isaac   DI  = D[field];
53820cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
539afe6d6adSToby Isaac   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
54020cf1dd8SToby Isaac   Np = fgeom->numPoints;
54120cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
54220cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
54320cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
5449f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
54520cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
546*4dcf62a8SSatish Balay     fegeom.n = 0;
547*4dcf62a8SSatish Balay     fegeom.v = 0;
548*4dcf62a8SSatish Balay     fegeom.J = 0;
549*4dcf62a8SSatish Balay     fegeom.detJ = 0;
5504bee2e38SMatthew G. Knepley     if (isAffine) {
5514bee2e38SMatthew G. Knepley       fegeom.v    = x;
5524bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
5534bee2e38SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*dE*dE];
5549f209ee4SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*dE*dE];
5554bee2e38SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e];
5564bee2e38SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE];
5579f209ee4SMatthew G. Knepley 
5589f209ee4SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*dE*dE];
5599f209ee4SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*dE*dE];
5609f209ee4SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e];
5614bee2e38SMatthew G. Knepley     }
56220cf1dd8SToby Isaac     ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr);
56320cf1dd8SToby Isaac     ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr);
56420cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
5654bee2e38SMatthew G. Knepley       PetscReal w;
5664bee2e38SMatthew G. Knepley       PetscInt  c, d;
5674bee2e38SMatthew G. Knepley 
56820cf1dd8SToby Isaac       if (isAffine) {
5694bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
57020cf1dd8SToby Isaac       } else {
5713fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
5729f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
5739f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
5744bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
5754bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
5769f209ee4SMatthew G. Knepley 
5779f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
5789f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
5799f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
58020cf1dd8SToby Isaac       }
5814bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
58220cf1dd8SToby Isaac       if (debug > 1 && q < Np) {
5834bee2e38SMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr);
5847be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX)
5854bee2e38SMatthew G. Knepley         ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr);
58620cf1dd8SToby Isaac #endif
58720cf1dd8SToby Isaac       }
58820cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
5899f209ee4SMatthew G. Knepley       ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, face*Nq+q, B, D, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);
5909f209ee4SMatthew G. Knepley       if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dimAux, NfAux, NbAux, NcAux, auxOnBd ? q : face*Nq+q, BAux, DAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
5914bee2e38SMatthew G. Knepley       if (f0_func) {
5924bee2e38SMatthew 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]);
5934bee2e38SMatthew G. Knepley         for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w;
5944bee2e38SMatthew G. Knepley       }
59520cf1dd8SToby Isaac       if (f1_func) {
5964bee2e38SMatthew 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]);
5974bee2e38SMatthew G. Knepley         for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w;
59820cf1dd8SToby Isaac       }
59920cf1dd8SToby Isaac     }
6009f209ee4SMatthew G. Knepley     ierr = PetscFEUpdateElementVec_Internal(fe, dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], basisReal, basisDerReal, &cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr);
60120cf1dd8SToby Isaac     cOffset    += totDim;
60220cf1dd8SToby Isaac     cOffsetAux += totDimAux;
60320cf1dd8SToby Isaac   }
60420cf1dd8SToby Isaac   PetscFunctionReturn(0);
60520cf1dd8SToby Isaac }
60620cf1dd8SToby Isaac 
6074bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
6084bee2e38SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
60920cf1dd8SToby Isaac {
61020cf1dd8SToby Isaac   const PetscInt     debug      = 0;
6114bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
6124bee2e38SMatthew G. Knepley   PetscPointJac      g0_func, g1_func, g2_func, g3_func;
61320cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
61420cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
61520cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
61620cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
61720cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
61820cf1dd8SToby Isaac   PetscQuadrature    quad;
6194bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
62020cf1dd8SToby Isaac   const PetscScalar *constants;
62120cf1dd8SToby Isaac   PetscReal         *x;
62220cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
62320cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
62420cf1dd8SToby Isaac   PetscInt           NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
62520cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
62620cf1dd8SToby Isaac   PetscInt           dE, Np;
62720cf1dd8SToby Isaac   PetscBool          isAffine;
62820cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
62920cf1dd8SToby Isaac   PetscInt           qNc, Nq, q;
63020cf1dd8SToby Isaac   PetscErrorCode     ierr;
63120cf1dd8SToby Isaac 
63220cf1dd8SToby Isaac   PetscFunctionBegin;
6334bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
6344bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
6354bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
6364bee2e38SMatthew G. Knepley   ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr);
6374bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
6384bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
6394bee2e38SMatthew G. Knepley   ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr);
6404bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
6414bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
6424bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
64320cf1dd8SToby Isaac   switch(jtype) {
6444bee2e38SMatthew G. Knepley   case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
6454bee2e38SMatthew G. Knepley   case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
6464bee2e38SMatthew G. Knepley   case PETSCFE_JACOBIAN:     ierr = PetscDSGetJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break;
64720cf1dd8SToby Isaac   }
64820cf1dd8SToby Isaac   if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0);
6494bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
6504bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
6514bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
6524bee2e38SMatthew G. Knepley   ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr);
6534bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
6544bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
6554bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
6564bee2e38SMatthew G. Knepley   if (dsAux) {
6574bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
6584bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
6594bee2e38SMatthew G. Knepley     ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr);
6604bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr);
6614bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
6624bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
6634bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
6644bee2e38SMatthew G. Knepley     ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);
66520cf1dd8SToby Isaac   }
66620cf1dd8SToby Isaac   NbI = Nb[fieldI], NbJ = Nb[fieldJ];
66720cf1dd8SToby Isaac   NcI = Nc[fieldI], NcJ = Nc[fieldJ];
66820cf1dd8SToby Isaac   BI  = B[fieldI],  BJ  = B[fieldJ];
66920cf1dd8SToby Isaac   DI  = D[fieldI],  DJ  = D[fieldJ];
67020cf1dd8SToby Isaac   /* Initialize here in case the function is not defined */
67120cf1dd8SToby Isaac   ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
67220cf1dd8SToby Isaac   ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
67320cf1dd8SToby Isaac   ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
67420cf1dd8SToby Isaac   ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
67520cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
67620cf1dd8SToby Isaac   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
6774bee2e38SMatthew G. Knepley   Np = cgeom->numPoints;
6784bee2e38SMatthew G. Knepley   dE = cgeom->dimEmbed;
6794bee2e38SMatthew G. Knepley   isAffine = cgeom->isAffine;
6804bee2e38SMatthew G. Knepley   for (e = 0; e < Ne; ++e) {
6814bee2e38SMatthew G. Knepley     PetscFEGeom fegeom;
6824bee2e38SMatthew G. Knepley 
6834bee2e38SMatthew G. Knepley     if (isAffine) {
6844bee2e38SMatthew G. Knepley       fegeom.v    = x;
6854bee2e38SMatthew G. Knepley       fegeom.xi   = cgeom->xi;
6864bee2e38SMatthew G. Knepley       fegeom.J    = &cgeom->J[e*dE*dE];
6874bee2e38SMatthew G. Knepley       fegeom.invJ = &cgeom->invJ[e*dE*dE];
6884bee2e38SMatthew G. Knepley       fegeom.detJ = &cgeom->detJ[e];
6894bee2e38SMatthew G. Knepley     }
69020cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
69120cf1dd8SToby Isaac       const PetscReal *BIq = &BI[q*NbI*NcI],     *BJq = &BJ[q*NbJ*NcJ];
69220cf1dd8SToby Isaac       const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim];
69320cf1dd8SToby Isaac       PetscReal        w;
6944bee2e38SMatthew G. Knepley       PetscInt         c;
69520cf1dd8SToby Isaac 
69620cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
69720cf1dd8SToby Isaac       if (isAffine) {
6984bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x);
69920cf1dd8SToby Isaac       } else {
7004bee2e38SMatthew G. Knepley         fegeom.v    = &cgeom->v[(e*Np+q)*dE];
7014bee2e38SMatthew G. Knepley         fegeom.J    = &cgeom->J[(e*Np+q)*dE*dE];
7024bee2e38SMatthew G. Knepley         fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE];
7034bee2e38SMatthew G. Knepley         fegeom.detJ = &cgeom->detJ[e*Np+q];
70420cf1dd8SToby Isaac       }
7054bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
706a8f1f9e5SMatthew G. Knepley       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
707a8f1f9e5SMatthew G. Knepley       if (dsAux)      {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dim, NfAux, NbAux, NcAux, q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
70820cf1dd8SToby Isaac       if (g0_func) {
70920cf1dd8SToby Isaac         ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
7104bee2e38SMatthew 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);
71120cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
71220cf1dd8SToby Isaac       }
71320cf1dd8SToby Isaac       if (g1_func) {
7144bee2e38SMatthew G. Knepley         ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
7154bee2e38SMatthew 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);
7164bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
71720cf1dd8SToby Isaac       }
71820cf1dd8SToby Isaac       if (g2_func) {
7194bee2e38SMatthew G. Knepley         ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
7204bee2e38SMatthew 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);
7214bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
72220cf1dd8SToby Isaac       }
72320cf1dd8SToby Isaac       if (g3_func) {
7244bee2e38SMatthew G. Knepley         ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
7254bee2e38SMatthew 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);
7264bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
72720cf1dd8SToby Isaac       }
72820cf1dd8SToby Isaac 
729a8f1f9e5SMatthew G. Knepley       ierr = PetscFEUpdateElementMat_Internal(feI, feJ, dim, NbI, NcI, BIq, DIq, basisReal, basisDerReal, NbJ, NcJ, BJq, DJq, testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr);
73020cf1dd8SToby Isaac     }
73120cf1dd8SToby Isaac     if (debug > 1) {
73220cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
73320cf1dd8SToby Isaac 
73420cf1dd8SToby Isaac       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
73520cf1dd8SToby Isaac       for (fc = 0; fc < NcI; ++fc) {
73620cf1dd8SToby Isaac         for (f = 0; f < NbI; ++f) {
73720cf1dd8SToby Isaac           const PetscInt i = offsetI + f*NcI+fc;
73820cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
73920cf1dd8SToby Isaac             for (g = 0; g < NbJ; ++g) {
74020cf1dd8SToby Isaac               const PetscInt j = offsetJ + g*NcJ+gc;
74120cf1dd8SToby 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);
74220cf1dd8SToby Isaac             }
74320cf1dd8SToby Isaac           }
74420cf1dd8SToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
74520cf1dd8SToby Isaac         }
74620cf1dd8SToby Isaac       }
74720cf1dd8SToby Isaac     }
74820cf1dd8SToby Isaac     cOffset    += totDim;
74920cf1dd8SToby Isaac     cOffsetAux += totDimAux;
75020cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
75120cf1dd8SToby Isaac   }
75220cf1dd8SToby Isaac   PetscFunctionReturn(0);
75320cf1dd8SToby Isaac }
75420cf1dd8SToby Isaac 
7554bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
7564bee2e38SMatthew G. Knepley                                                 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
75720cf1dd8SToby Isaac {
75820cf1dd8SToby Isaac   const PetscInt     debug      = 0;
7594bee2e38SMatthew G. Knepley   PetscFE            feI, feJ;
7604bee2e38SMatthew G. Knepley   PetscBdPointJac    g0_func, g1_func, g2_func, g3_func;
76120cf1dd8SToby Isaac   PetscInt           cOffset    = 0; /* Offset into coefficients[] for element e */
76220cf1dd8SToby Isaac   PetscInt           cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
76320cf1dd8SToby Isaac   PetscInt           eOffset    = 0; /* Offset into elemMat[] for element e */
76420cf1dd8SToby Isaac   PetscInt           offsetI    = 0; /* Offset into an element vector for fieldI */
76520cf1dd8SToby Isaac   PetscInt           offsetJ    = 0; /* Offset into an element vector for fieldJ */
76620cf1dd8SToby Isaac   PetscQuadrature    quad;
7674bee2e38SMatthew G. Knepley   PetscScalar       *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal;
76820cf1dd8SToby Isaac   const PetscScalar *constants;
76920cf1dd8SToby Isaac   PetscReal         *x;
77020cf1dd8SToby Isaac   PetscReal        **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
77120cf1dd8SToby Isaac   PetscInt          *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
77220cf1dd8SToby Isaac   PetscInt           NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
77320cf1dd8SToby Isaac   PetscInt           dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
77420cf1dd8SToby Isaac   PetscBool          isAffine;
77520cf1dd8SToby Isaac   const PetscReal   *quadPoints, *quadWeights;
77620cf1dd8SToby Isaac   PetscInt           qNc, Nq, q, Np, dE;
77720cf1dd8SToby Isaac   PetscErrorCode     ierr;
77820cf1dd8SToby Isaac 
77920cf1dd8SToby Isaac   PetscFunctionBegin;
7804bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr);
7814bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr);
7824bee2e38SMatthew G. Knepley   ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr);
7834bee2e38SMatthew G. Knepley   ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr);
7844bee2e38SMatthew G. Knepley   ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr);
7854bee2e38SMatthew G. Knepley   ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr);
7864bee2e38SMatthew G. Knepley   ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr);
7874bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr);
7884bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr);
7894bee2e38SMatthew G. Knepley   ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr);
7904bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr);
7914bee2e38SMatthew G. Knepley   ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr);
7924bee2e38SMatthew G. Knepley   ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);
7934bee2e38SMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr);
7944bee2e38SMatthew G. Knepley   ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr);
7954bee2e38SMatthew G. Knepley   ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
7964bee2e38SMatthew G. Knepley   ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr);
7974bee2e38SMatthew G. Knepley   ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr);
7984bee2e38SMatthew G. Knepley   if (dsAux) {
7994bee2e38SMatthew G. Knepley     ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr);
8004bee2e38SMatthew G. Knepley     ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr);
8014bee2e38SMatthew G. Knepley     ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr);
8024bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr);
8034bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr);
8044bee2e38SMatthew G. Knepley     ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr);
8054bee2e38SMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr);
8064bee2e38SMatthew G. Knepley     ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);
80720cf1dd8SToby Isaac   }
80820cf1dd8SToby Isaac   NbI = Nb[fieldI], NbJ = Nb[fieldJ];
80920cf1dd8SToby Isaac   NcI = Nc[fieldI], NcJ = Nc[fieldJ];
81020cf1dd8SToby Isaac   BI  = B[fieldI],  BJ  = B[fieldJ];
81120cf1dd8SToby Isaac   DI  = D[fieldI],  DJ  = D[fieldJ];
81220cf1dd8SToby Isaac   /* Initialize here in case the function is not defined */
81320cf1dd8SToby Isaac   ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
81420cf1dd8SToby Isaac   ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
81520cf1dd8SToby Isaac   ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
81620cf1dd8SToby Isaac   ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
81720cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
8184bee2e38SMatthew G. Knepley   if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
81920cf1dd8SToby Isaac   Np = fgeom->numPoints;
82020cf1dd8SToby Isaac   dE = fgeom->dimEmbed;
82120cf1dd8SToby Isaac   isAffine = fgeom->isAffine;
82220cf1dd8SToby Isaac   for (e = 0; e < Ne; ++e) {
8239f209ee4SMatthew G. Knepley     PetscFEGeom    fegeom, cgeom;
82420cf1dd8SToby Isaac     const PetscInt face = fgeom->face[e][0];
825*4dcf62a8SSatish Balay     fegeom.n = 0;
826*4dcf62a8SSatish Balay     fegeom.v = 0;
827*4dcf62a8SSatish Balay     fegeom.J = 0;
828*4dcf62a8SSatish Balay     fegeom.detJ = 0;
8294bee2e38SMatthew G. Knepley     if (isAffine) {
8304bee2e38SMatthew G. Knepley       fegeom.v    = x;
8314bee2e38SMatthew G. Knepley       fegeom.xi   = fgeom->xi;
8324bee2e38SMatthew G. Knepley       fegeom.J    = &fgeom->J[e*dE*dE];
8339f209ee4SMatthew G. Knepley       fegeom.invJ = &fgeom->invJ[e*dE*dE];
8344bee2e38SMatthew G. Knepley       fegeom.detJ = &fgeom->detJ[e];
8354bee2e38SMatthew G. Knepley       fegeom.n    = &fgeom->n[e*dE];
8369f209ee4SMatthew G. Knepley 
8379f209ee4SMatthew G. Knepley       cgeom.J     = &fgeom->suppJ[0][e*dE*dE];
8389f209ee4SMatthew G. Knepley       cgeom.invJ  = &fgeom->suppInvJ[0][e*dE*dE];
8399f209ee4SMatthew G. Knepley       cgeom.detJ  = &fgeom->suppDetJ[0][e];
8404bee2e38SMatthew G. Knepley     }
84120cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
84220cf1dd8SToby Isaac       const PetscReal *BIq = &BI[(face*Nq+q)*NbI*NcI],     *BJq = &BJ[(face*Nq+q)*NbJ*NcJ];
84320cf1dd8SToby Isaac       const PetscReal *DIq = &DI[(face*Nq+q)*NbI*NcI*dim], *DJq = &DJ[(face*Nq+q)*NbJ*NcJ*dim];
84420cf1dd8SToby Isaac       PetscReal  w;
8454bee2e38SMatthew G. Knepley       PetscInt   c;
84620cf1dd8SToby Isaac 
84720cf1dd8SToby Isaac       if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  quad point %d\n", q);CHKERRQ(ierr);}
84820cf1dd8SToby Isaac       if (isAffine) {
8494bee2e38SMatthew G. Knepley         CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x);
85020cf1dd8SToby Isaac       } else {
8513fe841f2SMatthew G. Knepley         fegeom.v    = &fgeom->v[(e*Np+q)*dE];
8529f209ee4SMatthew G. Knepley         fegeom.J    = &fgeom->J[(e*Np+q)*dE*dE];
8539f209ee4SMatthew G. Knepley         fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE];
8544bee2e38SMatthew G. Knepley         fegeom.detJ = &fgeom->detJ[e*Np+q];
8554bee2e38SMatthew G. Knepley         fegeom.n    = &fgeom->n[(e*Np+q)*dE];
8569f209ee4SMatthew G. Knepley 
8579f209ee4SMatthew G. Knepley         cgeom.J     = &fgeom->suppJ[0][(e*Np+q)*dE*dE];
8589f209ee4SMatthew G. Knepley         cgeom.invJ  = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE];
8599f209ee4SMatthew G. Knepley         cgeom.detJ  = &fgeom->suppDetJ[0][e*Np+q];
86020cf1dd8SToby Isaac       }
8614bee2e38SMatthew G. Knepley       w = fegeom.detJ[0]*quadWeights[q];
8629f209ee4SMatthew G. Knepley       if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, face*Nq+q, B, D, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);}
8639f209ee4SMatthew G. Knepley       if (dsAux)      {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);}
86420cf1dd8SToby Isaac       if (g0_func) {
86520cf1dd8SToby Isaac         ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr);
8664bee2e38SMatthew 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);
86720cf1dd8SToby Isaac         for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
86820cf1dd8SToby Isaac       }
86920cf1dd8SToby Isaac       if (g1_func) {
8704bee2e38SMatthew G. Knepley         ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
8714bee2e38SMatthew 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);
8724bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w;
87320cf1dd8SToby Isaac       }
87420cf1dd8SToby Isaac       if (g2_func) {
8754bee2e38SMatthew G. Knepley         ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr);
8764bee2e38SMatthew 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);
8774bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w;
87820cf1dd8SToby Isaac       }
87920cf1dd8SToby Isaac       if (g3_func) {
8804bee2e38SMatthew G. Knepley         ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr);
8814bee2e38SMatthew 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);
8824bee2e38SMatthew G. Knepley         for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w;
88320cf1dd8SToby Isaac       }
88420cf1dd8SToby Isaac 
8859f209ee4SMatthew G. Knepley       ierr = PetscFEUpdateElementMat_Internal(feI, feJ, dim, NbI, NcI, BIq, DIq, basisReal, basisDerReal, NbJ, NcJ, BJq, DJq, testReal, testDerReal, &cgeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr);
88620cf1dd8SToby Isaac     }
88720cf1dd8SToby Isaac     if (debug > 1) {
88820cf1dd8SToby Isaac       PetscInt fc, f, gc, g;
88920cf1dd8SToby Isaac 
89020cf1dd8SToby Isaac       ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr);
89120cf1dd8SToby Isaac       for (fc = 0; fc < NcI; ++fc) {
89220cf1dd8SToby Isaac         for (f = 0; f < NbI; ++f) {
89320cf1dd8SToby Isaac           const PetscInt i = offsetI + f*NcI+fc;
89420cf1dd8SToby Isaac           for (gc = 0; gc < NcJ; ++gc) {
89520cf1dd8SToby Isaac             for (g = 0; g < NbJ; ++g) {
89620cf1dd8SToby Isaac               const PetscInt j = offsetJ + g*NcJ+gc;
89720cf1dd8SToby 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);
89820cf1dd8SToby Isaac             }
89920cf1dd8SToby Isaac           }
90020cf1dd8SToby Isaac           ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr);
90120cf1dd8SToby Isaac         }
90220cf1dd8SToby Isaac       }
90320cf1dd8SToby Isaac     }
90420cf1dd8SToby Isaac     cOffset    += totDim;
90520cf1dd8SToby Isaac     cOffsetAux += totDimAux;
90620cf1dd8SToby Isaac     eOffset    += PetscSqr(totDim);
90720cf1dd8SToby Isaac   }
90820cf1dd8SToby Isaac   PetscFunctionReturn(0);
90920cf1dd8SToby Isaac }
91020cf1dd8SToby Isaac 
91120cf1dd8SToby Isaac PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
91220cf1dd8SToby Isaac {
91320cf1dd8SToby Isaac   PetscFunctionBegin;
91420cf1dd8SToby Isaac   fem->ops->setfromoptions          = NULL;
91520cf1dd8SToby Isaac   fem->ops->setup                   = PetscFESetUp_Basic;
91620cf1dd8SToby Isaac   fem->ops->view                    = PetscFEView_Basic;
91720cf1dd8SToby Isaac   fem->ops->destroy                 = PetscFEDestroy_Basic;
91820cf1dd8SToby Isaac   fem->ops->getdimension            = PetscFEGetDimension_Basic;
91920cf1dd8SToby Isaac   fem->ops->gettabulation           = PetscFEGetTabulation_Basic;
92020cf1dd8SToby Isaac   fem->ops->integrate               = PetscFEIntegrate_Basic;
921afe6d6adSToby Isaac   fem->ops->integratebd             = PetscFEIntegrateBd_Basic;
92220cf1dd8SToby Isaac   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
92320cf1dd8SToby Isaac   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
92420cf1dd8SToby Isaac   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
92520cf1dd8SToby Isaac   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
92620cf1dd8SToby Isaac   fem->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
92720cf1dd8SToby Isaac   PetscFunctionReturn(0);
92820cf1dd8SToby Isaac }
92920cf1dd8SToby Isaac 
93020cf1dd8SToby Isaac /*MC
93120cf1dd8SToby Isaac   PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
93220cf1dd8SToby Isaac 
93320cf1dd8SToby Isaac   Level: intermediate
93420cf1dd8SToby Isaac 
93520cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
93620cf1dd8SToby Isaac M*/
93720cf1dd8SToby Isaac 
93820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
93920cf1dd8SToby Isaac {
94020cf1dd8SToby Isaac   PetscFE_Basic *b;
94120cf1dd8SToby Isaac   PetscErrorCode ierr;
94220cf1dd8SToby Isaac 
94320cf1dd8SToby Isaac   PetscFunctionBegin;
94420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
94520cf1dd8SToby Isaac   ierr      = PetscNewLog(fem,&b);CHKERRQ(ierr);
94620cf1dd8SToby Isaac   fem->data = b;
94720cf1dd8SToby Isaac 
94820cf1dd8SToby Isaac   ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr);
94920cf1dd8SToby Isaac   PetscFunctionReturn(0);
95020cf1dd8SToby Isaac }
951