120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscblaslapack.h> 320cf1dd8SToby Isaac 42b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem) 520cf1dd8SToby Isaac { 620cf1dd8SToby Isaac PetscFE_Basic *b = (PetscFE_Basic *) fem->data; 720cf1dd8SToby Isaac PetscErrorCode ierr; 820cf1dd8SToby Isaac 920cf1dd8SToby Isaac PetscFunctionBegin; 1020cf1dd8SToby Isaac ierr = PetscFree(b);CHKERRQ(ierr); 1120cf1dd8SToby Isaac PetscFunctionReturn(0); 1220cf1dd8SToby Isaac } 1320cf1dd8SToby Isaac 142b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v) 1520cf1dd8SToby Isaac { 16d9bac1caSLisandro Dalcin PetscInt dim, Nc; 17d9bac1caSLisandro Dalcin PetscSpace basis = NULL; 18d9bac1caSLisandro Dalcin PetscDualSpace dual = NULL; 19d9bac1caSLisandro Dalcin PetscQuadrature quad = NULL; 2020cf1dd8SToby Isaac PetscErrorCode ierr; 2120cf1dd8SToby Isaac 2220cf1dd8SToby Isaac PetscFunctionBegin; 23d9bac1caSLisandro Dalcin ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 24d9bac1caSLisandro Dalcin ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 2520cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr); 2620cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr); 27d9bac1caSLisandro Dalcin ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 28d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 29d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "Basic Finite Element in %D dimensions with %D components\n",dim,Nc);CHKERRQ(ierr); 30d9bac1caSLisandro Dalcin if (basis) {ierr = PetscSpaceView(basis, v);CHKERRQ(ierr);} 31d9bac1caSLisandro Dalcin if (dual) {ierr = PetscDualSpaceView(dual, v);CHKERRQ(ierr);} 32d9bac1caSLisandro Dalcin if (quad) {ierr = PetscQuadratureView(quad, v);CHKERRQ(ierr);} 33d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 3420cf1dd8SToby Isaac PetscFunctionReturn(0); 3520cf1dd8SToby Isaac } 3620cf1dd8SToby Isaac 372b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v) 3820cf1dd8SToby Isaac { 3920cf1dd8SToby Isaac PetscBool iascii; 4020cf1dd8SToby Isaac PetscErrorCode ierr; 4120cf1dd8SToby Isaac 4220cf1dd8SToby Isaac PetscFunctionBegin; 43d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 44d9bac1caSLisandro Dalcin if (iascii) {ierr = PetscFEView_Basic_Ascii(fe, v);CHKERRQ(ierr);} 4520cf1dd8SToby Isaac PetscFunctionReturn(0); 4620cf1dd8SToby Isaac } 4720cf1dd8SToby Isaac 4820cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */ 4915310ec5SMatthew G. Knepley PetscErrorCode PetscFESetUp_Basic(PetscFE fem) 5020cf1dd8SToby Isaac { 5120cf1dd8SToby Isaac PetscScalar *work, *invVscalar; 5220cf1dd8SToby Isaac PetscBLASInt *pivots; 5320cf1dd8SToby Isaac PetscBLASInt n, info; 5420cf1dd8SToby Isaac PetscInt pdim, j; 5520cf1dd8SToby Isaac PetscErrorCode ierr; 5620cf1dd8SToby Isaac 5720cf1dd8SToby Isaac PetscFunctionBegin; 5820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 5920cf1dd8SToby Isaac ierr = PetscMalloc1(pdim*pdim,&fem->invV);CHKERRQ(ierr); 6020cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 6120cf1dd8SToby Isaac ierr = PetscMalloc1(pdim*pdim,&invVscalar);CHKERRQ(ierr); 6220cf1dd8SToby Isaac #else 6320cf1dd8SToby Isaac invVscalar = fem->invV; 6420cf1dd8SToby Isaac #endif 6520cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 6620cf1dd8SToby Isaac PetscReal *Bf; 6720cf1dd8SToby Isaac PetscQuadrature f; 6820cf1dd8SToby Isaac const PetscReal *points, *weights; 6920cf1dd8SToby Isaac PetscInt Nc, Nq, q, k, c; 7020cf1dd8SToby Isaac 7120cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr); 7220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 7320cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*Nq*pdim,&Bf);CHKERRQ(ierr); 7420cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr); 7520cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 7620cf1dd8SToby Isaac /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */ 7720cf1dd8SToby Isaac invVscalar[j*pdim+k] = 0.0; 7820cf1dd8SToby Isaac 7920cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 8020cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) invVscalar[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c]; 8120cf1dd8SToby Isaac } 8220cf1dd8SToby Isaac } 8320cf1dd8SToby Isaac ierr = PetscFree(Bf);CHKERRQ(ierr); 8420cf1dd8SToby Isaac } 85*ea2bdf6dSBarry Smith 8620cf1dd8SToby Isaac ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr); 8720cf1dd8SToby Isaac n = pdim; 8820cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invVscalar, &n, pivots, &info)); 89*ea2bdf6dSBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 90*ea2bdf6dSBarry Smith PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invVscalar, &n, pivots, work, &pdim, &info)); 91*ea2bdf6dSBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 9220cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 9320cf1dd8SToby Isaac for (j = 0; j < pdim*pdim; j++) fem->invV[j] = PetscRealPart(invVscalar[j]); 9420cf1dd8SToby Isaac ierr = PetscFree(invVscalar);CHKERRQ(ierr); 9520cf1dd8SToby Isaac #endif 9620cf1dd8SToby Isaac ierr = PetscFree2(pivots,work);CHKERRQ(ierr); 9720cf1dd8SToby Isaac PetscFunctionReturn(0); 9820cf1dd8SToby Isaac } 9920cf1dd8SToby Isaac 10020cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim) 10120cf1dd8SToby Isaac { 10220cf1dd8SToby Isaac PetscErrorCode ierr; 10320cf1dd8SToby Isaac 10420cf1dd8SToby Isaac PetscFunctionBegin; 10520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr); 10620cf1dd8SToby Isaac PetscFunctionReturn(0); 10720cf1dd8SToby Isaac } 10820cf1dd8SToby Isaac 10915310ec5SMatthew G. Knepley PetscErrorCode PetscFEGetTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H) 11020cf1dd8SToby Isaac { 11120cf1dd8SToby Isaac DM dm; 11220cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 11320cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 11420cf1dd8SToby Isaac PetscInt Nc; /* Field components */ 11520cf1dd8SToby Isaac PetscReal *tmpB, *tmpD, *tmpH; 11620cf1dd8SToby Isaac PetscInt p, d, j, k, c; 11720cf1dd8SToby Isaac PetscErrorCode ierr; 11820cf1dd8SToby Isaac 11920cf1dd8SToby Isaac PetscFunctionBegin; 12020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 12120cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 12220cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 12320cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 12420cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 12520cf1dd8SToby Isaac if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 12620cf1dd8SToby Isaac if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 12720cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 12820cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);CHKERRQ(ierr); 12920cf1dd8SToby Isaac /* Translate to the nodal basis */ 13020cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 13120cf1dd8SToby Isaac if (B) { 13220cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 13320cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 13420cf1dd8SToby Isaac const PetscInt i = (p*pdim + j)*Nc; 13520cf1dd8SToby Isaac 13620cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 13720cf1dd8SToby Isaac B[i+c] = 0.0; 13820cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 13920cf1dd8SToby Isaac B[i+c] += fem->invV[k*pdim+j] * tmpB[(p*pdim + k)*Nc+c]; 14020cf1dd8SToby Isaac } 14120cf1dd8SToby Isaac } 14220cf1dd8SToby Isaac } 14320cf1dd8SToby Isaac } 14420cf1dd8SToby Isaac if (D) { 14520cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 14620cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 14720cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 14820cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 14920cf1dd8SToby Isaac const PetscInt i = ((p*pdim + j)*Nc + c)*dim + d; 15020cf1dd8SToby Isaac 15120cf1dd8SToby Isaac D[i] = 0.0; 15220cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 15320cf1dd8SToby Isaac D[i] += fem->invV[k*pdim+j] * tmpD[((p*pdim + k)*Nc + c)*dim + d]; 15420cf1dd8SToby Isaac } 15520cf1dd8SToby Isaac } 15620cf1dd8SToby Isaac } 15720cf1dd8SToby Isaac } 15820cf1dd8SToby Isaac } 15920cf1dd8SToby Isaac if (H) { 16020cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 16120cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 16220cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 16320cf1dd8SToby Isaac for (d = 0; d < dim*dim; ++d) { 16420cf1dd8SToby Isaac const PetscInt i = ((p*pdim + j)*Nc + c)*dim*dim + d; 16520cf1dd8SToby Isaac 16620cf1dd8SToby Isaac H[i] = 0.0; 16720cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 16820cf1dd8SToby Isaac H[i] += fem->invV[k*pdim+j] * tmpH[((p*pdim + k)*Nc + c)*dim*dim + d]; 16920cf1dd8SToby Isaac } 17020cf1dd8SToby Isaac } 17120cf1dd8SToby Isaac } 17220cf1dd8SToby Isaac } 17320cf1dd8SToby Isaac } 17420cf1dd8SToby Isaac } 17520cf1dd8SToby Isaac if (B) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 17620cf1dd8SToby Isaac if (D) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 17720cf1dd8SToby Isaac if (H) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 17820cf1dd8SToby Isaac PetscFunctionReturn(0); 17920cf1dd8SToby Isaac } 18020cf1dd8SToby Isaac 1812b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1824bee2e38SMatthew G. Knepley const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 18320cf1dd8SToby Isaac { 18420cf1dd8SToby Isaac const PetscInt debug = 0; 1854bee2e38SMatthew G. Knepley PetscFE fe; 18620cf1dd8SToby Isaac PetscPointFunc obj_func; 18720cf1dd8SToby Isaac PetscQuadrature quad; 1884bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x; 18920cf1dd8SToby Isaac const PetscScalar *constants; 19020cf1dd8SToby Isaac PetscReal *x; 19120cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL; 19220cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 19320cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 19420cf1dd8SToby Isaac PetscBool isAffine; 19520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 19620cf1dd8SToby Isaac PetscInt qNc, Nq, q; 19720cf1dd8SToby Isaac PetscErrorCode ierr; 19820cf1dd8SToby Isaac 19920cf1dd8SToby Isaac PetscFunctionBegin; 2004bee2e38SMatthew G. Knepley ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr); 20120cf1dd8SToby Isaac if (!obj_func) PetscFunctionReturn(0); 2024bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 2034bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 2044bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 2054bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 2064bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 2074bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 2084bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 2094bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 2104bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 2114bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 2124bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 2134bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr); 2144bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 2154bee2e38SMatthew G. Knepley if (dsAux) { 2164bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 2174bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 2184bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 2194bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 2204bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 2214bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 2224bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 2234bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 22420cf1dd8SToby Isaac } 22520cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 2264bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 22720cf1dd8SToby Isaac Np = cgeom->numPoints; 22820cf1dd8SToby Isaac dE = cgeom->dimEmbed; 22920cf1dd8SToby Isaac isAffine = cgeom->isAffine; 23020cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 2314bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 23220cf1dd8SToby Isaac 23320cf1dd8SToby Isaac if (isAffine) { 2344bee2e38SMatthew G. Knepley fegeom.v = x; 2354bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2364bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[e*dE*dE]; 2374bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*dE*dE]; 2384bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e]; 23920cf1dd8SToby Isaac } 2404bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 2414bee2e38SMatthew G. Knepley PetscScalar integrand; 2424bee2e38SMatthew G. Knepley PetscReal w; 2434bee2e38SMatthew G. Knepley 2444bee2e38SMatthew G. Knepley if (isAffine) { 2454bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 2464bee2e38SMatthew G. Knepley } else { 2474bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 2484bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 2494bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 2504bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 2514bee2e38SMatthew G. Knepley } 2524bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 25320cf1dd8SToby Isaac if (debug > 1 && q < Np) { 2544bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 2557be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 2564bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 25720cf1dd8SToby Isaac #endif 25820cf1dd8SToby Isaac } 25920cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 260a8f1f9e5SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, dim, Nf, Nb, Nc, q, B, D, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 261a8f1f9e5SMatthew 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);} 2624bee2e38SMatthew 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); 2634bee2e38SMatthew G. Knepley integrand *= w; 26420cf1dd8SToby Isaac integral[e*Nf+field] += integrand; 26520cf1dd8SToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);} 26620cf1dd8SToby Isaac } 26720cf1dd8SToby Isaac cOffset += totDim; 26820cf1dd8SToby Isaac cOffsetAux += totDimAux; 26920cf1dd8SToby Isaac } 27020cf1dd8SToby Isaac PetscFunctionReturn(0); 27120cf1dd8SToby Isaac } 27220cf1dd8SToby Isaac 2732b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, 274afe6d6adSToby Isaac PetscBdPointFunc obj_func, 2754bee2e38SMatthew G. Knepley PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 276afe6d6adSToby Isaac { 277afe6d6adSToby Isaac const PetscInt debug = 0; 2784bee2e38SMatthew G. Knepley PetscFE fe; 279afe6d6adSToby Isaac PetscQuadrature quad; 2804bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 281afe6d6adSToby Isaac const PetscScalar *constants; 282afe6d6adSToby Isaac PetscReal *x; 283afe6d6adSToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL; 284afe6d6adSToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 285afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 286afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 287afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 288afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 289afe6d6adSToby Isaac PetscErrorCode ierr; 290afe6d6adSToby Isaac 291afe6d6adSToby Isaac PetscFunctionBegin; 292afe6d6adSToby Isaac if (!obj_func) PetscFunctionReturn(0); 2934bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 2944bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 2954bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 2964bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 2974bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 2984bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 2994bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 3004bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 3014bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 3024bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 3034bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 3044bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr); 3054bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 3064bee2e38SMatthew G. Knepley if (dsAux) { 3074bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 3084bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 3094bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 3104bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 3114bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 3124bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 3134bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 3144bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 315afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 3164bee2e38SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 3174bee2e38SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 318afe6d6adSToby Isaac } 319afe6d6adSToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 320afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 321afe6d6adSToby Isaac Np = fgeom->numPoints; 322afe6d6adSToby Isaac dE = fgeom->dimEmbed; 323afe6d6adSToby Isaac isAffine = fgeom->isAffine; 324afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 3259f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 326afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 3274dcf62a8SSatish Balay fegeom.n = 0; 3284dcf62a8SSatish Balay fegeom.v = 0; 3294dcf62a8SSatish Balay fegeom.J = 0; 3304dcf62a8SSatish Balay fegeom.detJ = 0; 3314bee2e38SMatthew G. Knepley if (isAffine) { 3324bee2e38SMatthew G. Knepley fegeom.v = x; 3334bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3344bee2e38SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 3359f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 3364bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 3374bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 3389f209ee4SMatthew G. Knepley 3399f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*dE*dE]; 3409f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*dE*dE]; 3419f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e]; 3424bee2e38SMatthew G. Knepley } 343afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 344afe6d6adSToby Isaac PetscScalar integrand; 3454bee2e38SMatthew G. Knepley PetscReal w; 346afe6d6adSToby Isaac 347afe6d6adSToby Isaac if (isAffine) { 3484bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 349afe6d6adSToby Isaac } else { 3503fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 3519f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 3529f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 3534bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 3544bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 3559f209ee4SMatthew G. Knepley 3569f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 3579f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 3589f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 359afe6d6adSToby Isaac } 3604bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 361afe6d6adSToby Isaac if (debug > 1 && q < Np) { 3624bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 363afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3644bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 365afe6d6adSToby Isaac #endif 366afe6d6adSToby Isaac } 367afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 3689f209ee4SMatthew 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); 3699f209ee4SMatthew 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);} 3704bee2e38SMatthew 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); 3714bee2e38SMatthew G. Knepley integrand *= w; 372afe6d6adSToby Isaac integral[e*Nf+field] += integrand; 373afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);} 374afe6d6adSToby Isaac } 375afe6d6adSToby Isaac cOffset += totDim; 376afe6d6adSToby Isaac cOffsetAux += totDimAux; 377afe6d6adSToby Isaac } 378afe6d6adSToby Isaac PetscFunctionReturn(0); 379afe6d6adSToby Isaac } 380afe6d6adSToby Isaac 3814bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 3824bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 38320cf1dd8SToby Isaac { 38420cf1dd8SToby Isaac const PetscInt debug = 0; 3854bee2e38SMatthew G. Knepley PetscFE fe; 38620cf1dd8SToby Isaac PetscPointFunc f0_func; 38720cf1dd8SToby Isaac PetscPointFunc f1_func; 38820cf1dd8SToby Isaac PetscQuadrature quad; 3894bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 39020cf1dd8SToby Isaac const PetscScalar *constants; 39120cf1dd8SToby Isaac PetscReal *x; 39220cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 39320cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 39420cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 39520cf1dd8SToby Isaac PetscBool isAffine; 39620cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3974bee2e38SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 39820cf1dd8SToby Isaac PetscErrorCode ierr; 39920cf1dd8SToby Isaac 40020cf1dd8SToby Isaac PetscFunctionBegin; 4014bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 4024bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 4034bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 4044bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 4054bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 4064bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 4074bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 4084bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 4094bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 4104bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 4114bee2e38SMatthew G. Knepley ierr = PetscDSGetResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 4124bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 4134bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 4144bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 4154bee2e38SMatthew G. Knepley if (!f0_func && !f1_func) PetscFunctionReturn(0); 4164bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr); 4174bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 4184bee2e38SMatthew G. Knepley if (dsAux) { 4194bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 4204bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 4214bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 4224bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 4234bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 4244bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 4254bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 4264bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 42720cf1dd8SToby Isaac } 42820cf1dd8SToby Isaac NbI = Nb[field]; 42920cf1dd8SToby Isaac NcI = Nc[field]; 43020cf1dd8SToby Isaac BI = B[field]; 43120cf1dd8SToby Isaac DI = D[field]; 43220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 4334bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 43420cf1dd8SToby Isaac Np = cgeom->numPoints; 43520cf1dd8SToby Isaac dE = cgeom->dimEmbed; 43620cf1dd8SToby Isaac isAffine = cgeom->isAffine; 43720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4384bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 43920cf1dd8SToby Isaac 4404bee2e38SMatthew G. Knepley if (isAffine) { 4414bee2e38SMatthew G. Knepley fegeom.v = x; 4424bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 4434bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[e*dE*dE]; 4444bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*dE*dE]; 4454bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e]; 4464bee2e38SMatthew G. Knepley } 447580bdb30SBarry Smith ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr); 448580bdb30SBarry Smith ierr = PetscArrayzero(f1, Nq*NcI*dim);CHKERRQ(ierr); 44920cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4504bee2e38SMatthew G. Knepley PetscReal w; 4514bee2e38SMatthew G. Knepley PetscInt c, d; 45220cf1dd8SToby Isaac 45320cf1dd8SToby Isaac if (isAffine) { 4544bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 45520cf1dd8SToby Isaac } else { 4564bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 4574bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 4584bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 4594bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 46020cf1dd8SToby Isaac } 4614bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 46220cf1dd8SToby Isaac if (debug > 1 && q < Np) { 4634bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 4647be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4654bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 46620cf1dd8SToby Isaac #endif 46720cf1dd8SToby Isaac } 46820cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 469a8f1f9e5SMatthew 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); 470a8f1f9e5SMatthew 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);} 4714bee2e38SMatthew G. Knepley if (f0_func) { 4724bee2e38SMatthew 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]); 4734bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 4744bee2e38SMatthew G. Knepley } 47520cf1dd8SToby Isaac if (f1_func) { 4764bee2e38SMatthew 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]); 4774bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 47820cf1dd8SToby Isaac } 47920cf1dd8SToby Isaac } 480a8f1f9e5SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, dim, Nq, NbI, NcI, BI, DI, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 48120cf1dd8SToby Isaac cOffset += totDim; 48220cf1dd8SToby Isaac cOffsetAux += totDimAux; 48320cf1dd8SToby Isaac } 48420cf1dd8SToby Isaac PetscFunctionReturn(0); 48520cf1dd8SToby Isaac } 48620cf1dd8SToby Isaac 4874bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 4884bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 48920cf1dd8SToby Isaac { 49020cf1dd8SToby Isaac const PetscInt debug = 0; 4914bee2e38SMatthew G. Knepley PetscFE fe; 49220cf1dd8SToby Isaac PetscBdPointFunc f0_func; 49320cf1dd8SToby Isaac PetscBdPointFunc f1_func; 49420cf1dd8SToby Isaac PetscQuadrature quad; 4954bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 49620cf1dd8SToby Isaac const PetscScalar *constants; 49720cf1dd8SToby Isaac PetscReal *x; 49820cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 49920cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 5007be5e748SToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 501b07bd890SMatthew G. Knepley PetscBool isAffine, auxOnBd = PETSC_FALSE; 50220cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 50320cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 50420cf1dd8SToby Isaac PetscErrorCode ierr; 50520cf1dd8SToby Isaac 50620cf1dd8SToby Isaac PetscFunctionBegin; 5074bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 5084bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 5094bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 5104bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 5114bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 5124bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 5134bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 5144bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 5154bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 5164bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 5174bee2e38SMatthew G. Knepley ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 51820cf1dd8SToby Isaac if (!f0_func && !f1_func) PetscFunctionReturn(0); 5194bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 5204bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 5214bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 5224bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr); 5234bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 5244bee2e38SMatthew G. Knepley if (dsAux) { 5254bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 5264bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 5274bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 5284bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 5294bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 5304bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 5314bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 5324bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 5337be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 5344bee2e38SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 5354bee2e38SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr);} 53620cf1dd8SToby Isaac } 53720cf1dd8SToby Isaac NbI = Nb[field]; 53820cf1dd8SToby Isaac NcI = Nc[field]; 53920cf1dd8SToby Isaac BI = B[field]; 54020cf1dd8SToby Isaac DI = D[field]; 54120cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 542afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 54320cf1dd8SToby Isaac Np = fgeom->numPoints; 54420cf1dd8SToby Isaac dE = fgeom->dimEmbed; 54520cf1dd8SToby Isaac isAffine = fgeom->isAffine; 54620cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5479f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 54820cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5494dcf62a8SSatish Balay fegeom.n = 0; 5504dcf62a8SSatish Balay fegeom.v = 0; 5514dcf62a8SSatish Balay fegeom.J = 0; 5524dcf62a8SSatish Balay fegeom.detJ = 0; 5534bee2e38SMatthew G. Knepley if (isAffine) { 5544bee2e38SMatthew G. Knepley fegeom.v = x; 5554bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 5564bee2e38SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 5579f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 5584bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 5594bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 5609f209ee4SMatthew G. Knepley 5619f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*dE*dE]; 5629f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*dE*dE]; 5639f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e]; 5644bee2e38SMatthew G. Knepley } 565580bdb30SBarry Smith ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr); 566580bdb30SBarry Smith ierr = PetscArrayzero(f1, Nq*NcI*dim);CHKERRQ(ierr); 56720cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5684bee2e38SMatthew G. Knepley PetscReal w; 5694bee2e38SMatthew G. Knepley PetscInt c, d; 5704bee2e38SMatthew G. Knepley 57120cf1dd8SToby Isaac if (isAffine) { 5724bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 57320cf1dd8SToby Isaac } else { 5743fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 5759f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 5769f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 5774bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 5784bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 5799f209ee4SMatthew G. Knepley 5809f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 5819f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 5829f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 58320cf1dd8SToby Isaac } 5844bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 58520cf1dd8SToby Isaac if (debug > 1 && q < Np) { 5864bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 5877be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5884bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 58920cf1dd8SToby Isaac #endif 59020cf1dd8SToby Isaac } 59120cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 5929f209ee4SMatthew 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); 5939f209ee4SMatthew 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);} 5944bee2e38SMatthew G. Knepley if (f0_func) { 5954bee2e38SMatthew 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]); 5964bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 5974bee2e38SMatthew G. Knepley } 59820cf1dd8SToby Isaac if (f1_func) { 5994bee2e38SMatthew 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]); 6004bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 60120cf1dd8SToby Isaac } 60220cf1dd8SToby Isaac } 6039f209ee4SMatthew 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); 60420cf1dd8SToby Isaac cOffset += totDim; 60520cf1dd8SToby Isaac cOffsetAux += totDimAux; 60620cf1dd8SToby Isaac } 60720cf1dd8SToby Isaac PetscFunctionReturn(0); 60820cf1dd8SToby Isaac } 60920cf1dd8SToby Isaac 6104bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 6114bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 61220cf1dd8SToby Isaac { 61320cf1dd8SToby Isaac const PetscInt debug = 0; 6144bee2e38SMatthew G. Knepley PetscFE feI, feJ; 6154bee2e38SMatthew G. Knepley PetscPointJac g0_func, g1_func, g2_func, g3_func; 61620cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 61720cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 61820cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 61920cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 62020cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 62120cf1dd8SToby Isaac PetscQuadrature quad; 6224bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 62320cf1dd8SToby Isaac const PetscScalar *constants; 62420cf1dd8SToby Isaac PetscReal *x; 62520cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 62620cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 62720cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 62820cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 62920cf1dd8SToby Isaac PetscInt dE, Np; 63020cf1dd8SToby Isaac PetscBool isAffine; 63120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 63220cf1dd8SToby Isaac PetscInt qNc, Nq, q; 63320cf1dd8SToby Isaac PetscErrorCode ierr; 63420cf1dd8SToby Isaac 63520cf1dd8SToby Isaac PetscFunctionBegin; 6364bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 6374bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 6384bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 6394bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 6404bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 6414bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 6424bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 6434bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 6444bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 6454bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 64620cf1dd8SToby Isaac switch(jtype) { 6474bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 6484bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 6494bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscDSGetJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 65020cf1dd8SToby Isaac } 65120cf1dd8SToby Isaac if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 6524bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 6534bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 6544bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 6554bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &B, &D);CHKERRQ(ierr); 6564bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 6574bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 6584bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 6594bee2e38SMatthew G. Knepley if (dsAux) { 6604bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 6614bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 6624bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 6634bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 6644bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 6654bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 6664bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 6674bee2e38SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 66820cf1dd8SToby Isaac } 66920cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 67020cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 67120cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 67220cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 67320cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 674580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 675580bdb30SBarry Smith ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr); 676580bdb30SBarry Smith ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr); 677580bdb30SBarry Smith ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr); 67820cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 67920cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 6804bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 6814bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 6824bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 6834bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 6844bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 6854bee2e38SMatthew G. Knepley 6864bee2e38SMatthew G. Knepley if (isAffine) { 6874bee2e38SMatthew G. Knepley fegeom.v = x; 6884bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 6894bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[e*dE*dE]; 6904bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*dE*dE]; 6914bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e]; 6924bee2e38SMatthew G. Knepley } 69320cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 69420cf1dd8SToby Isaac const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ]; 69520cf1dd8SToby Isaac const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim]; 69620cf1dd8SToby Isaac PetscReal w; 6974bee2e38SMatthew G. Knepley PetscInt c; 69820cf1dd8SToby Isaac 69920cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 70020cf1dd8SToby Isaac if (isAffine) { 7014bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 70220cf1dd8SToby Isaac } else { 7034bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 7044bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 7054bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 7064bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 70720cf1dd8SToby Isaac } 7084bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 709a8f1f9e5SMatthew 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);} 710a8f1f9e5SMatthew 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);} 71120cf1dd8SToby Isaac if (g0_func) { 712580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 7134bee2e38SMatthew 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); 71420cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 71520cf1dd8SToby Isaac } 71620cf1dd8SToby Isaac if (g1_func) { 717580bdb30SBarry Smith ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr); 7184bee2e38SMatthew 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); 7194bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 72020cf1dd8SToby Isaac } 72120cf1dd8SToby Isaac if (g2_func) { 722580bdb30SBarry Smith ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr); 7234bee2e38SMatthew 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); 7244bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 72520cf1dd8SToby Isaac } 72620cf1dd8SToby Isaac if (g3_func) { 727580bdb30SBarry Smith ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr); 7284bee2e38SMatthew 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); 7294bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 73020cf1dd8SToby Isaac } 73120cf1dd8SToby Isaac 732a8f1f9e5SMatthew 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); 73320cf1dd8SToby Isaac } 73420cf1dd8SToby Isaac if (debug > 1) { 73520cf1dd8SToby Isaac PetscInt fc, f, gc, g; 73620cf1dd8SToby Isaac 73720cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 73820cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 73920cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 74020cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 74120cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 74220cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 74320cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 74420cf1dd8SToby 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); 74520cf1dd8SToby Isaac } 74620cf1dd8SToby Isaac } 74720cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 74820cf1dd8SToby Isaac } 74920cf1dd8SToby Isaac } 75020cf1dd8SToby Isaac } 75120cf1dd8SToby Isaac cOffset += totDim; 75220cf1dd8SToby Isaac cOffsetAux += totDimAux; 75320cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 75420cf1dd8SToby Isaac } 75520cf1dd8SToby Isaac PetscFunctionReturn(0); 75620cf1dd8SToby Isaac } 75720cf1dd8SToby Isaac 7582b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 7594bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 76020cf1dd8SToby Isaac { 76120cf1dd8SToby Isaac const PetscInt debug = 0; 7624bee2e38SMatthew G. Knepley PetscFE feI, feJ; 7634bee2e38SMatthew G. Knepley PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 76420cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 76520cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 76620cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 76720cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 76820cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 76920cf1dd8SToby Isaac PetscQuadrature quad; 7704bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 77120cf1dd8SToby Isaac const PetscScalar *constants; 77220cf1dd8SToby Isaac PetscReal *x; 77320cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 77420cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 77520cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 77620cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 77720cf1dd8SToby Isaac PetscBool isAffine; 77820cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 77920cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 78020cf1dd8SToby Isaac PetscErrorCode ierr; 78120cf1dd8SToby Isaac 78220cf1dd8SToby Isaac PetscFunctionBegin; 7834bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 7844bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 7854bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 7864bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 7874bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 7884bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 7894bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(ds, &Nb);CHKERRQ(ierr); 7904bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(ds, &Nc);CHKERRQ(ierr); 7914bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 7924bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 7934bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 7944bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 7954bee2e38SMatthew G. Knepley ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr); 7964bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 7974bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 7984bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 7994bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &B, &D);CHKERRQ(ierr); 8004bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 8014bee2e38SMatthew G. Knepley if (dsAux) { 8024bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 8034bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 8044bee2e38SMatthew G. Knepley ierr = PetscDSGetDimensions(dsAux, &NbAux);CHKERRQ(ierr); 8054bee2e38SMatthew G. Knepley ierr = PetscDSGetComponents(dsAux, &NcAux);CHKERRQ(ierr); 8064bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 8074bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 8084bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 8094bee2e38SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(dsAux, &BAux, &DAux);CHKERRQ(ierr); 81020cf1dd8SToby Isaac } 81120cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 81220cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 81320cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 81420cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 81520cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 816580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 817580bdb30SBarry Smith ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr); 818580bdb30SBarry Smith ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr); 819580bdb30SBarry Smith ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr); 82020cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 8214bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 82220cf1dd8SToby Isaac Np = fgeom->numPoints; 82320cf1dd8SToby Isaac dE = fgeom->dimEmbed; 82420cf1dd8SToby Isaac isAffine = fgeom->isAffine; 82520cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 8269f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 82720cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 8284dcf62a8SSatish Balay fegeom.n = 0; 8294dcf62a8SSatish Balay fegeom.v = 0; 8304dcf62a8SSatish Balay fegeom.J = 0; 8314dcf62a8SSatish Balay fegeom.detJ = 0; 8324bee2e38SMatthew G. Knepley if (isAffine) { 8334bee2e38SMatthew G. Knepley fegeom.v = x; 8344bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 8354bee2e38SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 8369f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 8374bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 8384bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 8399f209ee4SMatthew G. Knepley 8409f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*dE*dE]; 8419f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*dE*dE]; 8429f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e]; 8434bee2e38SMatthew G. Knepley } 84420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 84520cf1dd8SToby Isaac const PetscReal *BIq = &BI[(face*Nq+q)*NbI*NcI], *BJq = &BJ[(face*Nq+q)*NbJ*NcJ]; 84620cf1dd8SToby Isaac const PetscReal *DIq = &DI[(face*Nq+q)*NbI*NcI*dim], *DJq = &DJ[(face*Nq+q)*NbJ*NcJ*dim]; 84720cf1dd8SToby Isaac PetscReal w; 8484bee2e38SMatthew G. Knepley PetscInt c; 84920cf1dd8SToby Isaac 85020cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 85120cf1dd8SToby Isaac if (isAffine) { 8524bee2e38SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 85320cf1dd8SToby Isaac } else { 8543fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 8559f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 8569f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 8574bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 8584bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 8599f209ee4SMatthew G. Knepley 8609f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 8619f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 8629f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 86320cf1dd8SToby Isaac } 8644bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 8659f209ee4SMatthew 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);} 8669f209ee4SMatthew 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);} 86720cf1dd8SToby Isaac if (g0_func) { 868580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 8694bee2e38SMatthew 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); 87020cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 87120cf1dd8SToby Isaac } 87220cf1dd8SToby Isaac if (g1_func) { 873580bdb30SBarry Smith ierr = PetscArrayzero(g1, NcI*NcJ*dim);CHKERRQ(ierr); 8744bee2e38SMatthew 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); 8754bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 87620cf1dd8SToby Isaac } 87720cf1dd8SToby Isaac if (g2_func) { 878580bdb30SBarry Smith ierr = PetscArrayzero(g2, NcI*NcJ*dim);CHKERRQ(ierr); 8794bee2e38SMatthew 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); 8804bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 88120cf1dd8SToby Isaac } 88220cf1dd8SToby Isaac if (g3_func) { 883580bdb30SBarry Smith ierr = PetscArrayzero(g3, NcI*NcJ*dim*dim);CHKERRQ(ierr); 8844bee2e38SMatthew 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); 8854bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 88620cf1dd8SToby Isaac } 88720cf1dd8SToby Isaac 8889f209ee4SMatthew 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); 88920cf1dd8SToby Isaac } 89020cf1dd8SToby Isaac if (debug > 1) { 89120cf1dd8SToby Isaac PetscInt fc, f, gc, g; 89220cf1dd8SToby Isaac 89320cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 89420cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 89520cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 89620cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 89720cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 89820cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 89920cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 90020cf1dd8SToby 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); 90120cf1dd8SToby Isaac } 90220cf1dd8SToby Isaac } 90320cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 90420cf1dd8SToby Isaac } 90520cf1dd8SToby Isaac } 90620cf1dd8SToby Isaac } 90720cf1dd8SToby Isaac cOffset += totDim; 90820cf1dd8SToby Isaac cOffsetAux += totDimAux; 90920cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 91020cf1dd8SToby Isaac } 91120cf1dd8SToby Isaac PetscFunctionReturn(0); 91220cf1dd8SToby Isaac } 91320cf1dd8SToby Isaac 9142b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 91520cf1dd8SToby Isaac { 91620cf1dd8SToby Isaac PetscFunctionBegin; 91720cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 91820cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 91920cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 92020cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 92120cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 92220cf1dd8SToby Isaac fem->ops->gettabulation = PetscFEGetTabulation_Basic; 92320cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 924afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 92520cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 92620cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 92720cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 92820cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 92920cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 93020cf1dd8SToby Isaac PetscFunctionReturn(0); 93120cf1dd8SToby Isaac } 93220cf1dd8SToby Isaac 93320cf1dd8SToby Isaac /*MC 93420cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 93520cf1dd8SToby Isaac 93620cf1dd8SToby Isaac Level: intermediate 93720cf1dd8SToby Isaac 93820cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 93920cf1dd8SToby Isaac M*/ 94020cf1dd8SToby Isaac 94120cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 94220cf1dd8SToby Isaac { 94320cf1dd8SToby Isaac PetscFE_Basic *b; 94420cf1dd8SToby Isaac PetscErrorCode ierr; 94520cf1dd8SToby Isaac 94620cf1dd8SToby Isaac PetscFunctionBegin; 94720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 94820cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 94920cf1dd8SToby Isaac fem->data = b; 95020cf1dd8SToby Isaac 95120cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 95220cf1dd8SToby Isaac PetscFunctionReturn(0); 95320cf1dd8SToby Isaac } 954