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 178*4bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 179*4bee2e38SMatthew G. Knepley const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 18020cf1dd8SToby Isaac { 18120cf1dd8SToby Isaac const PetscInt debug = 0; 182*4bee2e38SMatthew G. Knepley PetscFE fe; 18320cf1dd8SToby Isaac PetscPointFunc obj_func; 18420cf1dd8SToby Isaac PetscQuadrature quad; 185*4bee2e38SMatthew 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; 197*4bee2e38SMatthew G. Knepley ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr); 19820cf1dd8SToby Isaac if (!obj_func) PetscFunctionReturn(0); 199*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 200*4bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 201*4bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 202*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 203*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 204*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 205*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 206*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 207*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 208*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 209*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 210*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr); 211*4bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 212*4bee2e38SMatthew G. Knepley if (dsAux) { 213*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 214*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 215*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 216*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 217*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 218*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 219*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 220*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 22120cf1dd8SToby Isaac } 22220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 223*4bee2e38SMatthew 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) { 228*4bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 22920cf1dd8SToby Isaac 23020cf1dd8SToby Isaac if (isAffine) { 231*4bee2e38SMatthew G. Knepley fegeom.v = x; 232*4bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 233*4bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[e*dE*dE]; 234*4bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*dE*dE]; 235*4bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e]; 23620cf1dd8SToby Isaac } 237*4bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 238*4bee2e38SMatthew G. Knepley PetscScalar integrand; 239*4bee2e38SMatthew G. Knepley PetscReal w; 240*4bee2e38SMatthew G. Knepley 241*4bee2e38SMatthew G. Knepley if (isAffine) { 242*4bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 243*4bee2e38SMatthew G. Knepley } else { 244*4bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 245*4bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 246*4bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 247*4bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 248*4bee2e38SMatthew G. Knepley } 249*4bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 25020cf1dd8SToby Isaac if (debug > 1 && q < Np) { 251*4bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 2527be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 253*4bee2e38SMatthew 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);} 257*4bee2e38SMatthew G. Knepley ierr = EvaluateFieldJets(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 258*4bee2e38SMatthew G. Knepley if (dsAux) {ierr = EvaluateFieldJets(dsAux, dim, NfAux, NbAux, NcAux, q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 259*4bee2e38SMatthew 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); 260*4bee2e38SMatthew 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 270*4bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, 271afe6d6adSToby Isaac PetscBdPointFunc obj_func, 272*4bee2e38SMatthew G. Knepley PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 273afe6d6adSToby Isaac { 274afe6d6adSToby Isaac const PetscInt debug = 0; 275*4bee2e38SMatthew G. Knepley PetscFE fe; 276afe6d6adSToby Isaac PetscQuadrature quad; 277*4bee2e38SMatthew 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); 290*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 291*4bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 292*4bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 293*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 294*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 295*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 296*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 297*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 298*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 299*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 300*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 301*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr); 302*4bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 303*4bee2e38SMatthew G. Knepley if (dsAux) { 304*4bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 305*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 306*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 307*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 308*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 309*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 310*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 311*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 312afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 313*4bee2e38SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 314*4bee2e38SMatthew 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) { 322*4bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 323afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 324afe6d6adSToby Isaac 325*4bee2e38SMatthew G. Knepley if (isAffine) { 326*4bee2e38SMatthew G. Knepley fegeom.v = x; 327*4bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 328*4bee2e38SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 329*4bee2e38SMatthew G. Knepley fegeom.invJ = &fgeom->suppInvJ[0][e*dE*dE]; 330*4bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 331*4bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 332*4bee2e38SMatthew G. Knepley } 333afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 334afe6d6adSToby Isaac PetscScalar integrand; 335*4bee2e38SMatthew G. Knepley PetscReal w; 336afe6d6adSToby Isaac 337afe6d6adSToby Isaac if (isAffine) { 338*4bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 339afe6d6adSToby Isaac } else { 340*4bee2e38SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 341*4bee2e38SMatthew G. Knepley fegeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 342*4bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 343*4bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 344afe6d6adSToby Isaac } 345*4bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 346afe6d6adSToby Isaac if (debug > 1 && q < Np) { 347*4bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 348afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 349*4bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 350afe6d6adSToby Isaac #endif 351afe6d6adSToby Isaac } 352afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 353*4bee2e38SMatthew G. Knepley ierr = EvaluateFieldJets(ds, dim, Nf, Nb, Nc, face*Nq+q, B, D, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 354*4bee2e38SMatthew G. Knepley if (dsAux) {ierr = EvaluateFieldJets(dsAux, dimAux, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 355*4bee2e38SMatthew 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); 356*4bee2e38SMatthew G. Knepley integrand *= w; 357afe6d6adSToby Isaac integral[e*Nf+field] += integrand; 358afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);} 359afe6d6adSToby Isaac } 360afe6d6adSToby Isaac cOffset += totDim; 361afe6d6adSToby Isaac cOffsetAux += totDimAux; 362afe6d6adSToby Isaac } 363afe6d6adSToby Isaac PetscFunctionReturn(0); 364afe6d6adSToby Isaac } 365afe6d6adSToby Isaac 366*4bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 367*4bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 36820cf1dd8SToby Isaac { 36920cf1dd8SToby Isaac const PetscInt debug = 0; 370*4bee2e38SMatthew G. Knepley PetscFE fe; 37120cf1dd8SToby Isaac PetscPointFunc f0_func; 37220cf1dd8SToby Isaac PetscPointFunc f1_func; 37320cf1dd8SToby Isaac PetscQuadrature quad; 374*4bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 37520cf1dd8SToby Isaac const PetscScalar *constants; 37620cf1dd8SToby Isaac PetscReal *x; 37720cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 37820cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 37920cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 38020cf1dd8SToby Isaac PetscBool isAffine; 38120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 382*4bee2e38SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 38320cf1dd8SToby Isaac PetscErrorCode ierr; 38420cf1dd8SToby Isaac 38520cf1dd8SToby Isaac PetscFunctionBegin; 386*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 387*4bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 388*4bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 389*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 390*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 391*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 392*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 393*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 394*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 395*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 396*4bee2e38SMatthew G. Knepley ierr = PetscDSGetResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 397*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 398*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 399*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 400*4bee2e38SMatthew G. Knepley if (!f0_func && !f1_func) PetscFunctionReturn(0); 401*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr); 402*4bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 403*4bee2e38SMatthew G. Knepley if (dsAux) { 404*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 405*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 406*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 407*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 408*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 409*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 410*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 411*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 41220cf1dd8SToby Isaac } 41320cf1dd8SToby Isaac NbI = Nb[field]; 41420cf1dd8SToby Isaac NcI = Nc[field]; 41520cf1dd8SToby Isaac BI = B[field]; 41620cf1dd8SToby Isaac DI = D[field]; 41720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 418*4bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 41920cf1dd8SToby Isaac Np = cgeom->numPoints; 42020cf1dd8SToby Isaac dE = cgeom->dimEmbed; 42120cf1dd8SToby Isaac isAffine = cgeom->isAffine; 42220cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 423*4bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 42420cf1dd8SToby Isaac 425*4bee2e38SMatthew G. Knepley if (isAffine) { 426*4bee2e38SMatthew G. Knepley fegeom.v = x; 427*4bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 428*4bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[e*dE*dE]; 429*4bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*dE*dE]; 430*4bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e]; 431*4bee2e38SMatthew G. Knepley } 43220cf1dd8SToby Isaac ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr); 43320cf1dd8SToby Isaac ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 43420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 435*4bee2e38SMatthew G. Knepley PetscReal w; 436*4bee2e38SMatthew G. Knepley PetscInt c, d; 43720cf1dd8SToby Isaac 43820cf1dd8SToby Isaac if (isAffine) { 439*4bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 44020cf1dd8SToby Isaac } else { 441*4bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 442*4bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 443*4bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 444*4bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 44520cf1dd8SToby Isaac } 446*4bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 44720cf1dd8SToby Isaac if (debug > 1 && q < Np) { 448*4bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 4497be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 450*4bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 45120cf1dd8SToby Isaac #endif 45220cf1dd8SToby Isaac } 45320cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 454*4bee2e38SMatthew G. Knepley ierr = EvaluateFieldJets(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 455*4bee2e38SMatthew G. Knepley if (dsAux) {ierr = EvaluateFieldJets(dsAux, dim, NfAux, NbAux, NcAux, q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 456*4bee2e38SMatthew G. Knepley if (f0_func) { 457*4bee2e38SMatthew 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]); 458*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 459*4bee2e38SMatthew G. Knepley } 46020cf1dd8SToby Isaac if (f1_func) { 461*4bee2e38SMatthew 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]); 462*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 46320cf1dd8SToby Isaac } 46420cf1dd8SToby Isaac } 465*4bee2e38SMatthew G. Knepley ierr = UpdateElementVec(fe, dim, Nq, NbI, NcI, BI, DI, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 46620cf1dd8SToby Isaac cOffset += totDim; 46720cf1dd8SToby Isaac cOffsetAux += totDimAux; 46820cf1dd8SToby Isaac } 46920cf1dd8SToby Isaac PetscFunctionReturn(0); 47020cf1dd8SToby Isaac } 47120cf1dd8SToby Isaac 472*4bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 473*4bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 47420cf1dd8SToby Isaac { 47520cf1dd8SToby Isaac const PetscInt debug = 0; 476*4bee2e38SMatthew G. Knepley PetscFE fe; 47720cf1dd8SToby Isaac PetscBdPointFunc f0_func; 47820cf1dd8SToby Isaac PetscBdPointFunc f1_func; 47920cf1dd8SToby Isaac PetscQuadrature quad; 480*4bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 48120cf1dd8SToby Isaac const PetscScalar *constants; 48220cf1dd8SToby Isaac PetscReal *x; 48320cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 48420cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 4857be5e748SToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 486b07bd890SMatthew G. Knepley PetscBool isAffine, auxOnBd = PETSC_FALSE; 48720cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 48820cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 48920cf1dd8SToby Isaac PetscErrorCode ierr; 49020cf1dd8SToby Isaac 49120cf1dd8SToby Isaac PetscFunctionBegin; 492*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 493*4bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 494*4bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 495*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 496*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 497*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 498*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 499*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 500*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 501*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 502*4bee2e38SMatthew G. Knepley ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 50320cf1dd8SToby Isaac if (!f0_func && !f1_func) PetscFunctionReturn(0); 504*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 505*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 506*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 507*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr); 508*4bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 509*4bee2e38SMatthew G. Knepley if (dsAux) { 510*4bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 511*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 512*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 513*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 514*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 515*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 516*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 517*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 5187be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 519*4bee2e38SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 520*4bee2e38SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 52120cf1dd8SToby Isaac } 52220cf1dd8SToby Isaac NbI = Nb[field]; 52320cf1dd8SToby Isaac NcI = Nc[field]; 52420cf1dd8SToby Isaac BI = B[field]; 52520cf1dd8SToby Isaac DI = D[field]; 52620cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 527afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 52820cf1dd8SToby Isaac Np = fgeom->numPoints; 52920cf1dd8SToby Isaac dE = fgeom->dimEmbed; 53020cf1dd8SToby Isaac isAffine = fgeom->isAffine; 53120cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 532*4bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 53320cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 53420cf1dd8SToby Isaac 535*4bee2e38SMatthew G. Knepley if (isAffine) { 536*4bee2e38SMatthew G. Knepley fegeom.v = x; 537*4bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 538*4bee2e38SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 539*4bee2e38SMatthew G. Knepley fegeom.invJ = &fgeom->suppInvJ[0][e*dE*dE]; 540*4bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 541*4bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 542*4bee2e38SMatthew G. Knepley } 54320cf1dd8SToby Isaac ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr); 54420cf1dd8SToby Isaac ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 54520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 546*4bee2e38SMatthew G. Knepley PetscReal w; 547*4bee2e38SMatthew G. Knepley PetscInt c, d; 548*4bee2e38SMatthew G. Knepley 54920cf1dd8SToby Isaac if (isAffine) { 550*4bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 55120cf1dd8SToby Isaac } else { 552*4bee2e38SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 553*4bee2e38SMatthew G. Knepley fegeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 554*4bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 555*4bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 55620cf1dd8SToby Isaac } 557*4bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 55820cf1dd8SToby Isaac if (debug > 1 && q < Np) { 559*4bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 5607be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 561*4bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 56220cf1dd8SToby Isaac #endif 56320cf1dd8SToby Isaac } 56420cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 565*4bee2e38SMatthew G. Knepley ierr = EvaluateFieldJets(ds, dim, Nf, Nb, Nc, face*Nq+q, B, D, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 566*4bee2e38SMatthew G. Knepley if (dsAux) {ierr = EvaluateFieldJets(dsAux, dimAux, NfAux, NbAux, NcAux, auxOnBd ? q : face*Nq+q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 567*4bee2e38SMatthew G. Knepley if (f0_func) { 568*4bee2e38SMatthew 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]); 569*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 570*4bee2e38SMatthew G. Knepley } 57120cf1dd8SToby Isaac if (f1_func) { 572*4bee2e38SMatthew 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]); 573*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 57420cf1dd8SToby Isaac } 57520cf1dd8SToby Isaac } 576*4bee2e38SMatthew G. Knepley ierr = UpdateElementVec(fe, dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 57720cf1dd8SToby Isaac cOffset += totDim; 57820cf1dd8SToby Isaac cOffsetAux += totDimAux; 57920cf1dd8SToby Isaac } 58020cf1dd8SToby Isaac PetscFunctionReturn(0); 58120cf1dd8SToby Isaac } 58220cf1dd8SToby Isaac 583*4bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 584*4bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 58520cf1dd8SToby Isaac { 58620cf1dd8SToby Isaac const PetscInt debug = 0; 587*4bee2e38SMatthew G. Knepley PetscFE feI, feJ; 588*4bee2e38SMatthew G. Knepley PetscPointJac g0_func, g1_func, g2_func, g3_func; 58920cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 59020cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 59120cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 59220cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 59320cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 59420cf1dd8SToby Isaac PetscQuadrature quad; 595*4bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 59620cf1dd8SToby Isaac const PetscScalar *constants; 59720cf1dd8SToby Isaac PetscReal *x; 59820cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 59920cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 60020cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 60120cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 60220cf1dd8SToby Isaac PetscInt dE, Np; 60320cf1dd8SToby Isaac PetscBool isAffine; 60420cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 60520cf1dd8SToby Isaac PetscInt qNc, Nq, q; 60620cf1dd8SToby Isaac PetscErrorCode ierr; 60720cf1dd8SToby Isaac 60820cf1dd8SToby Isaac PetscFunctionBegin; 609*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 610*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 611*4bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 612*4bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 613*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 614*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 615*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 616*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 617*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 618*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 61920cf1dd8SToby Isaac switch(jtype) { 620*4bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 621*4bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 622*4bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscDSGetJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 62320cf1dd8SToby Isaac } 62420cf1dd8SToby Isaac if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 625*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 626*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 627*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 628*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr); 629*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 630*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 631*4bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 632*4bee2e38SMatthew G. Knepley if (dsAux) { 633*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 634*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 635*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 636*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 637*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 638*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 639*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 640*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 64120cf1dd8SToby Isaac } 64220cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 64320cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 64420cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 64520cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 64620cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 64720cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 64820cf1dd8SToby Isaac ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 64920cf1dd8SToby Isaac ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 65020cf1dd8SToby Isaac ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));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); 653*4bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 654*4bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 655*4bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 656*4bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 657*4bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 658*4bee2e38SMatthew G. Knepley 659*4bee2e38SMatthew G. Knepley if (isAffine) { 660*4bee2e38SMatthew G. Knepley fegeom.v = x; 661*4bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 662*4bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[e*dE*dE]; 663*4bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*dE*dE]; 664*4bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e]; 665*4bee2e38SMatthew G. Knepley } 66620cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 66720cf1dd8SToby Isaac const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ]; 66820cf1dd8SToby Isaac const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim]; 66920cf1dd8SToby Isaac PetscReal w; 670*4bee2e38SMatthew G. Knepley PetscInt c; 67120cf1dd8SToby Isaac 67220cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 67320cf1dd8SToby Isaac if (isAffine) { 674*4bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 67520cf1dd8SToby Isaac } else { 676*4bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 677*4bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 678*4bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 679*4bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 68020cf1dd8SToby Isaac } 681*4bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 682*4bee2e38SMatthew G. Knepley if (coefficients) {ierr = EvaluateFieldJets(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 683*4bee2e38SMatthew G. Knepley if (dsAux) {ierr = EvaluateFieldJets(dsAux, dim, NfAux, NbAux, NcAux, q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 68420cf1dd8SToby Isaac if (g0_func) { 68520cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 686*4bee2e38SMatthew 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); 68720cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 68820cf1dd8SToby Isaac } 68920cf1dd8SToby Isaac if (g1_func) { 690*4bee2e38SMatthew G. Knepley ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 691*4bee2e38SMatthew 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); 692*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 69320cf1dd8SToby Isaac } 69420cf1dd8SToby Isaac if (g2_func) { 695*4bee2e38SMatthew G. Knepley ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 696*4bee2e38SMatthew 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); 697*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 69820cf1dd8SToby Isaac } 69920cf1dd8SToby Isaac if (g3_func) { 700*4bee2e38SMatthew G. Knepley ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 701*4bee2e38SMatthew 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); 702*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 70320cf1dd8SToby Isaac } 70420cf1dd8SToby Isaac 705*4bee2e38SMatthew G. Knepley ierr = UpdateElementMat(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); 70620cf1dd8SToby Isaac } 70720cf1dd8SToby Isaac if (debug > 1) { 70820cf1dd8SToby Isaac PetscInt fc, f, gc, g; 70920cf1dd8SToby Isaac 71020cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 71120cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 71220cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 71320cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 71420cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 71520cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 71620cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 71720cf1dd8SToby 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); 71820cf1dd8SToby Isaac } 71920cf1dd8SToby Isaac } 72020cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 72120cf1dd8SToby Isaac } 72220cf1dd8SToby Isaac } 72320cf1dd8SToby Isaac } 72420cf1dd8SToby Isaac cOffset += totDim; 72520cf1dd8SToby Isaac cOffsetAux += totDimAux; 72620cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 72720cf1dd8SToby Isaac } 72820cf1dd8SToby Isaac PetscFunctionReturn(0); 72920cf1dd8SToby Isaac } 73020cf1dd8SToby Isaac 731*4bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 732*4bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 73320cf1dd8SToby Isaac { 73420cf1dd8SToby Isaac const PetscInt debug = 0; 735*4bee2e38SMatthew G. Knepley PetscFE feI, feJ; 736*4bee2e38SMatthew G. Knepley PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 73720cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 73820cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 73920cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 74020cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 74120cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 74220cf1dd8SToby Isaac PetscQuadrature quad; 743*4bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 74420cf1dd8SToby Isaac const PetscScalar *constants; 74520cf1dd8SToby Isaac PetscReal *x; 74620cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 74720cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 74820cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 74920cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 75020cf1dd8SToby Isaac PetscBool isAffine; 75120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 75220cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 75320cf1dd8SToby Isaac PetscErrorCode ierr; 75420cf1dd8SToby Isaac 75520cf1dd8SToby Isaac PetscFunctionBegin; 756*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 757*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 758*4bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 759*4bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 760*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 761*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 762*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 763*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 764*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 765*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 766*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 767*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 768*4bee2e38SMatthew G. Knepley ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr); 769*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 770*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 771*4bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 772*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr); 773*4bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 774*4bee2e38SMatthew G. Knepley if (dsAux) { 775*4bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 776*4bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 777*4bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 778*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 779*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 780*4bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 781*4bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 782*4bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 78320cf1dd8SToby Isaac } 78420cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 78520cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 78620cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 78720cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 78820cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 78920cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 79020cf1dd8SToby Isaac ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 79120cf1dd8SToby Isaac ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 79220cf1dd8SToby Isaac ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 79320cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 794*4bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 79520cf1dd8SToby Isaac Np = fgeom->numPoints; 79620cf1dd8SToby Isaac dE = fgeom->dimEmbed; 79720cf1dd8SToby Isaac isAffine = fgeom->isAffine; 79820cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 799*4bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 80020cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 80120cf1dd8SToby Isaac 802*4bee2e38SMatthew G. Knepley if (isAffine) { 803*4bee2e38SMatthew G. Knepley fegeom.v = x; 804*4bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 805*4bee2e38SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 806*4bee2e38SMatthew G. Knepley fegeom.invJ = &fgeom->suppInvJ[0][e*dE*dE]; 807*4bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 808*4bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 809*4bee2e38SMatthew G. Knepley } 81020cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 81120cf1dd8SToby Isaac const PetscReal *BIq = &BI[(face*Nq+q)*NbI*NcI], *BJq = &BJ[(face*Nq+q)*NbJ*NcJ]; 81220cf1dd8SToby Isaac const PetscReal *DIq = &DI[(face*Nq+q)*NbI*NcI*dim], *DJq = &DJ[(face*Nq+q)*NbJ*NcJ*dim]; 81320cf1dd8SToby Isaac PetscReal w; 814*4bee2e38SMatthew G. Knepley PetscInt c; 81520cf1dd8SToby Isaac 81620cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 81720cf1dd8SToby Isaac if (isAffine) { 818*4bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 81920cf1dd8SToby Isaac } else { 820*4bee2e38SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 821*4bee2e38SMatthew G. Knepley fegeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 822*4bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 823*4bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 82420cf1dd8SToby Isaac } 825*4bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 826*4bee2e38SMatthew G. Knepley if (coefficients) {ierr = EvaluateFieldJets(ds, dim, Nf, Nb, Nc, face*Nq+q, B, D, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 827*4bee2e38SMatthew G. Knepley if (dsAux) {ierr = EvaluateFieldJets(dsAux, dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 82820cf1dd8SToby Isaac if (g0_func) { 82920cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 830*4bee2e38SMatthew 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); 83120cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 83220cf1dd8SToby Isaac } 83320cf1dd8SToby Isaac if (g1_func) { 834*4bee2e38SMatthew G. Knepley ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 835*4bee2e38SMatthew 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); 836*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 83720cf1dd8SToby Isaac } 83820cf1dd8SToby Isaac if (g2_func) { 839*4bee2e38SMatthew G. Knepley ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 840*4bee2e38SMatthew 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); 841*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 84220cf1dd8SToby Isaac } 84320cf1dd8SToby Isaac if (g3_func) { 844*4bee2e38SMatthew G. Knepley ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 845*4bee2e38SMatthew 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); 846*4bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 84720cf1dd8SToby Isaac } 84820cf1dd8SToby Isaac 849*4bee2e38SMatthew G. Knepley ierr = UpdateElementMat(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); 85020cf1dd8SToby Isaac } 85120cf1dd8SToby Isaac if (debug > 1) { 85220cf1dd8SToby Isaac PetscInt fc, f, gc, g; 85320cf1dd8SToby Isaac 85420cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 85520cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 85620cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 85720cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 85820cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 85920cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 86020cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 86120cf1dd8SToby 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); 86220cf1dd8SToby Isaac } 86320cf1dd8SToby Isaac } 86420cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 86520cf1dd8SToby Isaac } 86620cf1dd8SToby Isaac } 86720cf1dd8SToby Isaac } 86820cf1dd8SToby Isaac cOffset += totDim; 86920cf1dd8SToby Isaac cOffsetAux += totDimAux; 87020cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 87120cf1dd8SToby Isaac } 87220cf1dd8SToby Isaac PetscFunctionReturn(0); 87320cf1dd8SToby Isaac } 87420cf1dd8SToby Isaac 87520cf1dd8SToby Isaac PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 87620cf1dd8SToby Isaac { 87720cf1dd8SToby Isaac PetscFunctionBegin; 87820cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 87920cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 88020cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 88120cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 88220cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 88320cf1dd8SToby Isaac fem->ops->gettabulation = PetscFEGetTabulation_Basic; 88420cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 885afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 88620cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 88720cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 88820cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 88920cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 89020cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 89120cf1dd8SToby Isaac PetscFunctionReturn(0); 89220cf1dd8SToby Isaac } 89320cf1dd8SToby Isaac 89420cf1dd8SToby Isaac /*MC 89520cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 89620cf1dd8SToby Isaac 89720cf1dd8SToby Isaac Level: intermediate 89820cf1dd8SToby Isaac 89920cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 90020cf1dd8SToby Isaac M*/ 90120cf1dd8SToby Isaac 90220cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 90320cf1dd8SToby Isaac { 90420cf1dd8SToby Isaac PetscFE_Basic *b; 90520cf1dd8SToby Isaac PetscErrorCode ierr; 90620cf1dd8SToby Isaac 90720cf1dd8SToby Isaac PetscFunctionBegin; 90820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 90920cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 91020cf1dd8SToby Isaac fem->data = b; 91120cf1dd8SToby Isaac 91220cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 91320cf1dd8SToby Isaac PetscFunctionReturn(0); 91420cf1dd8SToby Isaac } 915