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 17820cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrate_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 17920cf1dd8SToby Isaac const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 18020cf1dd8SToby Isaac { 18120cf1dd8SToby Isaac const PetscInt debug = 0; 18220cf1dd8SToby Isaac PetscPointFunc obj_func; 18320cf1dd8SToby Isaac PetscQuadrature quad; 18420cf1dd8SToby Isaac PetscScalar *u, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 18520cf1dd8SToby Isaac const PetscScalar *constants; 18620cf1dd8SToby Isaac PetscReal *x; 18720cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL; 18820cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 18920cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 19020cf1dd8SToby Isaac PetscBool isAffine; 19120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 19220cf1dd8SToby Isaac PetscInt qNc, Nq, q; 19320cf1dd8SToby Isaac PetscErrorCode ierr; 19420cf1dd8SToby Isaac 19520cf1dd8SToby Isaac PetscFunctionBegin; 19620cf1dd8SToby Isaac ierr = PetscDSGetObjective(prob, field, &obj_func);CHKERRQ(ierr); 19720cf1dd8SToby Isaac if (!obj_func) PetscFunctionReturn(0); 19820cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 19920cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 20020cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 20120cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 20220cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 20320cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 20420cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 20520cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 20620cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 20720cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 20820cf1dd8SToby Isaac ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr); 20920cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 21020cf1dd8SToby Isaac if (probAux) { 21120cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 21220cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 21320cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 21420cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 21520cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 21620cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 21720cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 21820cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 21920cf1dd8SToby Isaac ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 22020cf1dd8SToby Isaac } 22120cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 22220cf1dd8SToby Isaac Np = cgeom->numPoints; 22320cf1dd8SToby Isaac dE = cgeom->dimEmbed; 22420cf1dd8SToby Isaac isAffine = cgeom->isAffine; 22520cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 22620cf1dd8SToby Isaac const PetscReal *v0 = &cgeom->v[e*Np*dE]; 22720cf1dd8SToby Isaac const PetscReal *J = &cgeom->J[e*Np*dE*dE]; 22820cf1dd8SToby Isaac 22920cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 23020cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 23120cf1dd8SToby Isaac PetscScalar integrand; 23220cf1dd8SToby Isaac const PetscReal *v; 23320cf1dd8SToby Isaac const PetscReal *invJ; 23420cf1dd8SToby Isaac PetscReal detJ; 23520cf1dd8SToby Isaac 23620cf1dd8SToby Isaac if (isAffine) { 23720cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim, cgeom->xi, v0, J, &quadPoints[q*dim], x); 23820cf1dd8SToby Isaac v = x; 23920cf1dd8SToby Isaac invJ = &cgeom->invJ[e*dE*dE]; 24020cf1dd8SToby Isaac detJ = cgeom->detJ[e]; 24120cf1dd8SToby Isaac } else { 24220cf1dd8SToby Isaac v = &v0[q*dE]; 24320cf1dd8SToby Isaac invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 24420cf1dd8SToby Isaac detJ = cgeom->detJ[e*Np + q]; 24520cf1dd8SToby Isaac } 24620cf1dd8SToby Isaac if (debug > 1 && q < Np) { 24720cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 2487be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 24920cf1dd8SToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 25020cf1dd8SToby Isaac #endif 25120cf1dd8SToby Isaac } 25220cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 25320cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], NULL, u, u_x, NULL); 25420cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 25520cf1dd8SToby Isaac obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, v, numConstants, constants, &integrand); 25620cf1dd8SToby Isaac integrand *= detJ*quadWeights[q]; 25720cf1dd8SToby Isaac integral[e*Nf+field] += integrand; 25820cf1dd8SToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);} 25920cf1dd8SToby Isaac } 26020cf1dd8SToby Isaac cOffset += totDim; 26120cf1dd8SToby Isaac cOffsetAux += totDimAux; 26220cf1dd8SToby Isaac } 26320cf1dd8SToby Isaac PetscFunctionReturn(0); 26420cf1dd8SToby Isaac } 26520cf1dd8SToby Isaac 266afe6d6adSToby Isaac PetscErrorCode PetscFEIntegrateBd_Basic(PetscFE fem, PetscDS prob, PetscInt field, 267afe6d6adSToby Isaac PetscBdPointFunc obj_func, 268afe6d6adSToby Isaac PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 269afe6d6adSToby Isaac { 270afe6d6adSToby Isaac const PetscInt debug = 0; 271afe6d6adSToby Isaac PetscQuadrature quad; 272afe6d6adSToby Isaac PetscScalar *u, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 273afe6d6adSToby Isaac const PetscScalar *constants; 274afe6d6adSToby Isaac PetscReal *x; 275afe6d6adSToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL; 276afe6d6adSToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 277afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 278afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 279afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 280afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 281afe6d6adSToby Isaac PetscErrorCode ierr; 282afe6d6adSToby Isaac 283afe6d6adSToby Isaac PetscFunctionBegin; 284afe6d6adSToby Isaac if (!obj_func) PetscFunctionReturn(0); 285afe6d6adSToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 286afe6d6adSToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr); 287afe6d6adSToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 288afe6d6adSToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 289afe6d6adSToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 290afe6d6adSToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 291afe6d6adSToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 292afe6d6adSToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 293afe6d6adSToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 294afe6d6adSToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 295afe6d6adSToby Isaac ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr); 296afe6d6adSToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 297afe6d6adSToby Isaac if (probAux) { 298afe6d6adSToby Isaac ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 299afe6d6adSToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 300afe6d6adSToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 301afe6d6adSToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 302afe6d6adSToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 303afe6d6adSToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 304afe6d6adSToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 305afe6d6adSToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 306afe6d6adSToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 307afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 308afe6d6adSToby Isaac if (auxOnBd) {ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);} 309afe6d6adSToby Isaac else {ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);} 310afe6d6adSToby Isaac } 311afe6d6adSToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 312afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 313afe6d6adSToby Isaac Np = fgeom->numPoints; 314afe6d6adSToby Isaac dE = fgeom->dimEmbed; 315afe6d6adSToby Isaac isAffine = fgeom->isAffine; 316afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 317afe6d6adSToby Isaac const PetscReal *v0 = &fgeom->v[e*Np*dE]; 318afe6d6adSToby Isaac const PetscReal *J = &fgeom->J[e*Np*dE*dE]; 319afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 320afe6d6adSToby Isaac 321afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 322afe6d6adSToby Isaac const PetscReal *v; 323afe6d6adSToby Isaac const PetscReal *invJ; 324afe6d6adSToby Isaac const PetscReal *n; 325afe6d6adSToby Isaac PetscReal detJ; 326afe6d6adSToby Isaac PetscScalar integrand; 327afe6d6adSToby Isaac 328afe6d6adSToby Isaac if (isAffine) { 329afe6d6adSToby Isaac CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x); 330afe6d6adSToby Isaac v = x; 331afe6d6adSToby Isaac invJ = &fgeom->suppInvJ[0][e*dE*dE]; 332afe6d6adSToby Isaac detJ = fgeom->detJ[e]; 333afe6d6adSToby Isaac n = &fgeom->n[e*dE]; 334afe6d6adSToby Isaac } else { 335afe6d6adSToby Isaac v = &v0[q*dE]; 336afe6d6adSToby Isaac invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 337afe6d6adSToby Isaac detJ = fgeom->detJ[e*Np + q]; 338afe6d6adSToby Isaac n = &fgeom->n[(e*Np+q)*dE]; 339afe6d6adSToby Isaac } 340afe6d6adSToby Isaac if (debug > 1 && q < Np) { 341afe6d6adSToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 342afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 343afe6d6adSToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 344afe6d6adSToby Isaac #endif 345afe6d6adSToby Isaac } 346afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 347afe6d6adSToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], NULL, u, u_x, NULL); 348afe6d6adSToby Isaac if (probAux) EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 349afe6d6adSToby Isaac obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, v, n, numConstants, constants, &integrand); 350afe6d6adSToby Isaac integrand *= detJ*quadWeights[q]; 351afe6d6adSToby Isaac integral[e*Nf+field] += integrand; 352afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);} 353afe6d6adSToby Isaac } 354afe6d6adSToby Isaac cOffset += totDim; 355afe6d6adSToby Isaac cOffsetAux += totDimAux; 356afe6d6adSToby Isaac } 357afe6d6adSToby Isaac PetscFunctionReturn(0); 358afe6d6adSToby Isaac } 359afe6d6adSToby Isaac 36020cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 36120cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 36220cf1dd8SToby Isaac { 36320cf1dd8SToby Isaac const PetscInt debug = 0; 36420cf1dd8SToby Isaac PetscPointFunc f0_func; 36520cf1dd8SToby Isaac PetscPointFunc f1_func; 36620cf1dd8SToby Isaac PetscQuadrature quad; 36720cf1dd8SToby Isaac PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 36820cf1dd8SToby Isaac const PetscScalar *constants; 36920cf1dd8SToby Isaac PetscReal *x; 37020cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 37120cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 37220cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 37320cf1dd8SToby Isaac PetscInt dE, Np; 37420cf1dd8SToby Isaac PetscBool isAffine; 37520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 37620cf1dd8SToby Isaac PetscInt qNc, Nq, q; 37720cf1dd8SToby Isaac PetscErrorCode ierr; 37820cf1dd8SToby Isaac 37920cf1dd8SToby Isaac PetscFunctionBegin; 38020cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 38120cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 38220cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 38320cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 38420cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 38520cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 38620cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 38720cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 38820cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, field, &fOffset);CHKERRQ(ierr); 38920cf1dd8SToby Isaac ierr = PetscDSGetResidual(prob, field, &f0_func, &f1_func);CHKERRQ(ierr); 39020cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 39120cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 39220cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 39320cf1dd8SToby Isaac ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr); 39420cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 39520cf1dd8SToby Isaac if (probAux) { 39620cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 39720cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 39820cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 39920cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 40020cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 40120cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 40220cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 40320cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 40420cf1dd8SToby Isaac ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 40520cf1dd8SToby Isaac } 40620cf1dd8SToby Isaac NbI = Nb[field]; 40720cf1dd8SToby Isaac NcI = Nc[field]; 40820cf1dd8SToby Isaac BI = B[field]; 40920cf1dd8SToby Isaac DI = D[field]; 41020cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 41120cf1dd8SToby Isaac Np = cgeom->numPoints; 41220cf1dd8SToby Isaac dE = cgeom->dimEmbed; 41320cf1dd8SToby Isaac isAffine = cgeom->isAffine; 41420cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 41520cf1dd8SToby Isaac const PetscReal *v0 = &cgeom->v[e*Np*dE]; 41620cf1dd8SToby Isaac const PetscReal *J = &cgeom->J[e*Np*dE*dE]; 41720cf1dd8SToby Isaac 41820cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 41920cf1dd8SToby Isaac ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr); 42020cf1dd8SToby Isaac ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 42120cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 42220cf1dd8SToby Isaac const PetscReal *v; 42320cf1dd8SToby Isaac const PetscReal *invJ; 42420cf1dd8SToby Isaac PetscReal detJ; 42520cf1dd8SToby Isaac 42620cf1dd8SToby Isaac if (isAffine) { 42720cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim, cgeom->xi, v0, J, &quadPoints[q*dim], x); 42820cf1dd8SToby Isaac v = x; 42920cf1dd8SToby Isaac invJ = &cgeom->invJ[e*dE*dE]; 43020cf1dd8SToby Isaac detJ = cgeom->detJ[e]; 43120cf1dd8SToby Isaac } else { 43220cf1dd8SToby Isaac v = &v0[q*dE]; 43320cf1dd8SToby Isaac invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 43420cf1dd8SToby Isaac detJ = cgeom->detJ[e*Np + q]; 43520cf1dd8SToby Isaac } 43620cf1dd8SToby Isaac if (debug > 1 && q < Np) { 43720cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 4387be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 43920cf1dd8SToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 44020cf1dd8SToby Isaac #endif 44120cf1dd8SToby Isaac } 44220cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 44320cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 44420cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 44520cf1dd8SToby Isaac if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, numConstants, constants, &f0[q*NcI]); 44620cf1dd8SToby Isaac if (f1_func) { 44720cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 44820cf1dd8SToby Isaac f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, numConstants, constants, refSpaceDer); 44920cf1dd8SToby Isaac } 45020cf1dd8SToby Isaac TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL); 45120cf1dd8SToby Isaac } 45220cf1dd8SToby Isaac UpdateElementVec(dim, Nq, NbI, NcI, BI, DI, f0, f1, &elemVec[cOffset+fOffset]); 45320cf1dd8SToby Isaac cOffset += totDim; 45420cf1dd8SToby Isaac cOffsetAux += totDimAux; 45520cf1dd8SToby Isaac } 45620cf1dd8SToby Isaac PetscFunctionReturn(0); 45720cf1dd8SToby Isaac } 45820cf1dd8SToby Isaac 45920cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 46020cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 46120cf1dd8SToby Isaac { 46220cf1dd8SToby Isaac const PetscInt debug = 0; 46320cf1dd8SToby Isaac PetscBdPointFunc f0_func; 46420cf1dd8SToby Isaac PetscBdPointFunc f1_func; 46520cf1dd8SToby Isaac PetscQuadrature quad; 46620cf1dd8SToby Isaac PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 46720cf1dd8SToby Isaac const PetscScalar *constants; 46820cf1dd8SToby Isaac PetscReal *x; 46920cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 47020cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 4717be5e748SToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 472*b07bd890SMatthew G. Knepley PetscBool isAffine, auxOnBd = PETSC_FALSE; 47320cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 47420cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 47520cf1dd8SToby Isaac PetscErrorCode ierr; 47620cf1dd8SToby Isaac 47720cf1dd8SToby Isaac PetscFunctionBegin; 47820cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 47920cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr); 48020cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 48120cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 48220cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 48320cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 48420cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 48520cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 48620cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, field, &fOffset);CHKERRQ(ierr); 48720cf1dd8SToby Isaac ierr = PetscDSGetBdResidual(prob, field, &f0_func, &f1_func);CHKERRQ(ierr); 48820cf1dd8SToby Isaac if (!f0_func && !f1_func) PetscFunctionReturn(0); 48920cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 49020cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 49120cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 49220cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr); 49320cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 49420cf1dd8SToby Isaac if (probAux) { 4957be5e748SToby Isaac ierr = PetscDSGetSpatialDimension(probAux, &dimAux);CHKERRQ(ierr); 49620cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 49720cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 49820cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 49920cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 50020cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 50120cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 50220cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 50320cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 5047be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 5057be5e748SToby Isaac if (auxOnBd) {ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);} 5067be5e748SToby Isaac else {ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr);} 50720cf1dd8SToby Isaac } 50820cf1dd8SToby Isaac NbI = Nb[field]; 50920cf1dd8SToby Isaac NcI = Nc[field]; 51020cf1dd8SToby Isaac BI = B[field]; 51120cf1dd8SToby Isaac DI = D[field]; 51220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 513afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 51420cf1dd8SToby Isaac Np = fgeom->numPoints; 51520cf1dd8SToby Isaac dE = fgeom->dimEmbed; 51620cf1dd8SToby Isaac isAffine = fgeom->isAffine; 51720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 51820cf1dd8SToby Isaac const PetscReal *v0 = &fgeom->v[e*Np*dE]; 51920cf1dd8SToby Isaac const PetscReal *J = &fgeom->J[e*Np*dE*dE]; 52020cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 52120cf1dd8SToby Isaac 52220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 52320cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 52420cf1dd8SToby Isaac ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr); 52520cf1dd8SToby Isaac ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 52620cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 52720cf1dd8SToby Isaac const PetscReal *v; 52820cf1dd8SToby Isaac const PetscReal *invJ; 52920cf1dd8SToby Isaac const PetscReal *n; 53020cf1dd8SToby Isaac PetscReal detJ; 53120cf1dd8SToby Isaac if (isAffine) { 53220cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x); 53320cf1dd8SToby Isaac v = x; 53420cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][e*dE*dE]; 53520cf1dd8SToby Isaac detJ = fgeom->detJ[e]; 53620cf1dd8SToby Isaac n = &fgeom->n[e*dE]; 53720cf1dd8SToby Isaac } else { 53820cf1dd8SToby Isaac v = &v0[q*dE]; 53920cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 54020cf1dd8SToby Isaac detJ = fgeom->detJ[e*Np + q]; 54120cf1dd8SToby Isaac n = &fgeom->n[(e*Np+q)*dE]; 54220cf1dd8SToby Isaac } 54320cf1dd8SToby Isaac if (debug > 1 && q < Np) { 54420cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 5457be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 54620cf1dd8SToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 54720cf1dd8SToby Isaac #endif 54820cf1dd8SToby Isaac } 54920cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 55020cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 55147a6c934SMatthew G. Knepley if (probAux) EvaluateFieldJets(dimAux, NfAux, NbAux, NcAux, auxOnBd ? q : face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 55220cf1dd8SToby Isaac if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, n, numConstants, constants, &f0[q*NcI]); 55320cf1dd8SToby Isaac if (f1_func) { 55420cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 55520cf1dd8SToby Isaac f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, n, numConstants, constants, refSpaceDer); 55620cf1dd8SToby Isaac } 55720cf1dd8SToby Isaac TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL); 55820cf1dd8SToby Isaac } 55920cf1dd8SToby Isaac UpdateElementVec(dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], f0, f1, &elemVec[cOffset+fOffset]); 56020cf1dd8SToby Isaac cOffset += totDim; 56120cf1dd8SToby Isaac cOffsetAux += totDimAux; 56220cf1dd8SToby Isaac } 56320cf1dd8SToby Isaac PetscFunctionReturn(0); 56420cf1dd8SToby Isaac } 56520cf1dd8SToby Isaac 56620cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *geom, 56720cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 56820cf1dd8SToby Isaac { 56920cf1dd8SToby Isaac const PetscInt debug = 0; 57020cf1dd8SToby Isaac PetscPointJac g0_func; 57120cf1dd8SToby Isaac PetscPointJac g1_func; 57220cf1dd8SToby Isaac PetscPointJac g2_func; 57320cf1dd8SToby Isaac PetscPointJac g3_func; 57420cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 57520cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 57620cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 57720cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 57820cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 57920cf1dd8SToby Isaac PetscQuadrature quad; 58020cf1dd8SToby Isaac PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 58120cf1dd8SToby Isaac const PetscScalar *constants; 58220cf1dd8SToby Isaac PetscReal *x; 58320cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 58420cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 58520cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 58620cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 58720cf1dd8SToby Isaac PetscInt dE, Np; 58820cf1dd8SToby Isaac PetscBool isAffine; 58920cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 59020cf1dd8SToby Isaac PetscInt qNc, Nq, q; 59120cf1dd8SToby Isaac PetscErrorCode ierr; 59220cf1dd8SToby Isaac 59320cf1dd8SToby Isaac PetscFunctionBegin; 59420cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 59520cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 59620cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 59720cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 59820cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 59920cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 60020cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 60120cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 60220cf1dd8SToby Isaac switch(jtype) { 60320cf1dd8SToby Isaac case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 60420cf1dd8SToby Isaac case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 60520cf1dd8SToby Isaac case PETSCFE_JACOBIAN: ierr = PetscDSGetJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 60620cf1dd8SToby Isaac } 60720cf1dd8SToby Isaac if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 60820cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 60920cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 61020cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 61120cf1dd8SToby Isaac ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr); 61220cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr); 61320cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);CHKERRQ(ierr); 61420cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 61520cf1dd8SToby Isaac if (probAux) { 61620cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 61720cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 61820cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 61920cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 62020cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 62120cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 62220cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 62320cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 62420cf1dd8SToby Isaac ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 62520cf1dd8SToby Isaac } 62620cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 62720cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 62820cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 62920cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 63020cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 63120cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 63220cf1dd8SToby Isaac ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 63320cf1dd8SToby Isaac ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 63420cf1dd8SToby Isaac ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 63520cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 63620cf1dd8SToby Isaac Np = geom->numPoints; 63720cf1dd8SToby Isaac dE = geom->dimEmbed; 63820cf1dd8SToby Isaac isAffine = geom->isAffine; 63920cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 64020cf1dd8SToby Isaac const PetscReal *v0 = &geom->v[e*Np*dE]; 64120cf1dd8SToby Isaac const PetscReal *J = &geom->J[e*Np*dE*dE]; 64220cf1dd8SToby Isaac 64320cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 64420cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 64520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 64620cf1dd8SToby Isaac const PetscReal *v; 64720cf1dd8SToby Isaac const PetscReal *invJ; 64820cf1dd8SToby Isaac PetscReal detJ; 64920cf1dd8SToby Isaac const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ]; 65020cf1dd8SToby Isaac const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim]; 65120cf1dd8SToby Isaac PetscReal w; 65220cf1dd8SToby Isaac PetscInt f, g, fc, gc, c; 65320cf1dd8SToby Isaac 65420cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 65520cf1dd8SToby Isaac if (isAffine) { 65620cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim, geom->xi, v0, J, &quadPoints[q*dim], x); 65720cf1dd8SToby Isaac v = x; 65820cf1dd8SToby Isaac invJ = &geom->invJ[e*dE*dE]; 65920cf1dd8SToby Isaac detJ = geom->detJ[e]; 66020cf1dd8SToby Isaac } else { 66120cf1dd8SToby Isaac v = &v0[q*dE]; 66220cf1dd8SToby Isaac invJ = &geom->invJ[(e*Np+q)*dE*dE]; 66320cf1dd8SToby Isaac detJ = geom->detJ[e*Np + q]; 66420cf1dd8SToby Isaac } 66520cf1dd8SToby Isaac w = detJ*quadWeights[q]; 6667489efa5SMatthew G. Knepley if (coefficients) EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 66720cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 66820cf1dd8SToby Isaac if (g0_func) { 66920cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 67020cf1dd8SToby Isaac g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, g0); 67120cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 67220cf1dd8SToby Isaac } 67320cf1dd8SToby Isaac if (g1_func) { 67420cf1dd8SToby Isaac PetscInt d, d2; 67520cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 67620cf1dd8SToby Isaac g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer); 67720cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 67820cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 67920cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 68020cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] = 0.0; 68120cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 68220cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] *= w; 68320cf1dd8SToby Isaac } 68420cf1dd8SToby Isaac } 68520cf1dd8SToby Isaac } 68620cf1dd8SToby Isaac } 68720cf1dd8SToby Isaac if (g2_func) { 68820cf1dd8SToby Isaac PetscInt d, d2; 68920cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 69020cf1dd8SToby Isaac g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer); 69120cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 69220cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 69320cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 69420cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] = 0.0; 69520cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 69620cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] *= w; 69720cf1dd8SToby Isaac } 69820cf1dd8SToby Isaac } 69920cf1dd8SToby Isaac } 70020cf1dd8SToby Isaac } 70120cf1dd8SToby Isaac if (g3_func) { 70220cf1dd8SToby Isaac PetscInt d, d2, dp, d3; 70320cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 70420cf1dd8SToby Isaac g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer); 70520cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 70620cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 70720cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 70820cf1dd8SToby Isaac for (dp = 0; dp < dim; ++dp) { 70920cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0; 71020cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 71120cf1dd8SToby Isaac for (d3 = 0; d3 < dim; ++d3) { 71220cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3]; 71320cf1dd8SToby Isaac } 71420cf1dd8SToby Isaac } 71520cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w; 71620cf1dd8SToby Isaac } 71720cf1dd8SToby Isaac } 71820cf1dd8SToby Isaac } 71920cf1dd8SToby Isaac } 72020cf1dd8SToby Isaac } 72120cf1dd8SToby Isaac 72220cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 72320cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 72420cf1dd8SToby Isaac const PetscInt fidx = f*NcI+fc; /* Test function basis index */ 72520cf1dd8SToby Isaac const PetscInt i = offsetI+f; /* Element matrix row */ 72620cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 72720cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 72820cf1dd8SToby Isaac const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */ 72920cf1dd8SToby Isaac const PetscInt j = offsetJ+g; /* Element matrix column */ 73020cf1dd8SToby Isaac const PetscInt fOff = eOffset+i*totDim+j; 73120cf1dd8SToby Isaac PetscInt d, d2; 73220cf1dd8SToby Isaac 73320cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx]; 73420cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 73520cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d]; 73620cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx]; 73720cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 73820cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2]; 73920cf1dd8SToby Isaac } 74020cf1dd8SToby Isaac } 74120cf1dd8SToby Isaac } 74220cf1dd8SToby Isaac } 74320cf1dd8SToby Isaac } 74420cf1dd8SToby Isaac } 74520cf1dd8SToby Isaac } 74620cf1dd8SToby Isaac if (debug > 1) { 74720cf1dd8SToby Isaac PetscInt fc, f, gc, g; 74820cf1dd8SToby Isaac 74920cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 75020cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 75120cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 75220cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 75320cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 75420cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 75520cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 75620cf1dd8SToby 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); 75720cf1dd8SToby Isaac } 75820cf1dd8SToby Isaac } 75920cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 76020cf1dd8SToby Isaac } 76120cf1dd8SToby Isaac } 76220cf1dd8SToby Isaac } 76320cf1dd8SToby Isaac cOffset += totDim; 76420cf1dd8SToby Isaac cOffsetAux += totDimAux; 76520cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 76620cf1dd8SToby Isaac } 76720cf1dd8SToby Isaac PetscFunctionReturn(0); 76820cf1dd8SToby Isaac } 76920cf1dd8SToby Isaac 77020cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 77120cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 77220cf1dd8SToby Isaac { 77320cf1dd8SToby Isaac const PetscInt debug = 0; 77420cf1dd8SToby Isaac PetscBdPointJac g0_func; 77520cf1dd8SToby Isaac PetscBdPointJac g1_func; 77620cf1dd8SToby Isaac PetscBdPointJac g2_func; 77720cf1dd8SToby Isaac PetscBdPointJac g3_func; 77820cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 77920cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 78020cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 78120cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 78220cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 78320cf1dd8SToby Isaac PetscQuadrature quad; 78420cf1dd8SToby Isaac PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 78520cf1dd8SToby Isaac const PetscScalar *constants; 78620cf1dd8SToby Isaac PetscReal *x; 78720cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 78820cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 78920cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 79020cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 79120cf1dd8SToby Isaac PetscBool isAffine; 79220cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 79320cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 79420cf1dd8SToby Isaac PetscErrorCode ierr; 79520cf1dd8SToby Isaac 79620cf1dd8SToby Isaac PetscFunctionBegin; 79720cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 79820cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr); 79920cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 80020cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 80120cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 80220cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 80320cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 80420cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 80520cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr); 80620cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);CHKERRQ(ierr); 80720cf1dd8SToby Isaac ierr = PetscDSGetBdJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr); 80820cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 80920cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 81020cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 81120cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr); 81220cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 81320cf1dd8SToby Isaac if (probAux) { 81420cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 81520cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 81620cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 81720cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 81820cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 81920cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 82020cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 82120cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 82220cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 82320cf1dd8SToby Isaac } 82420cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 82520cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 82620cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 82720cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 82820cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 82920cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 83020cf1dd8SToby Isaac ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 83120cf1dd8SToby Isaac ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 83220cf1dd8SToby Isaac ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 83320cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 83420cf1dd8SToby Isaac Np = fgeom->numPoints; 83520cf1dd8SToby Isaac dE = fgeom->dimEmbed; 83620cf1dd8SToby Isaac isAffine = fgeom->isAffine; 83720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 83820cf1dd8SToby Isaac const PetscReal *v0 = &fgeom->v[e*Np*dE]; 83920cf1dd8SToby Isaac const PetscReal *J = &fgeom->J[e*Np*dE*dE]; 84020cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 84120cf1dd8SToby Isaac 84220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 84320cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 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; 84820cf1dd8SToby Isaac PetscInt f, g, fc, gc, c; 84920cf1dd8SToby Isaac const PetscReal *v; 85020cf1dd8SToby Isaac const PetscReal *invJ; 85120cf1dd8SToby Isaac const PetscReal *n; 85220cf1dd8SToby Isaac PetscReal detJ; 85320cf1dd8SToby Isaac 85420cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 85520cf1dd8SToby Isaac if (isAffine) { 85620cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x); 85720cf1dd8SToby Isaac v = x; 85820cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][e*dE*dE]; 85920cf1dd8SToby Isaac detJ = fgeom->detJ[e]; 86020cf1dd8SToby Isaac n = &fgeom->n[e*dE]; 86120cf1dd8SToby Isaac } else { 86220cf1dd8SToby Isaac v = &v0[q*dE]; 86320cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 86420cf1dd8SToby Isaac detJ = fgeom->detJ[e*Np + q]; 86520cf1dd8SToby Isaac n = &fgeom->n[(e*Np+q)*dE]; 86620cf1dd8SToby Isaac } 86720cf1dd8SToby Isaac w = detJ*quadWeights[q]; 86820cf1dd8SToby Isaac 8697489efa5SMatthew G. Knepley if (coefficients) EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 87020cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 87120cf1dd8SToby Isaac if (g0_func) { 87220cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 87320cf1dd8SToby Isaac g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, g0); 87420cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 87520cf1dd8SToby Isaac } 87620cf1dd8SToby Isaac if (g1_func) { 87720cf1dd8SToby Isaac PetscInt d, d2; 87820cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 87920cf1dd8SToby Isaac g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer); 88020cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 88120cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 88220cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 88320cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] = 0.0; 88420cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 88520cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] *= w; 88620cf1dd8SToby Isaac } 88720cf1dd8SToby Isaac } 88820cf1dd8SToby Isaac } 88920cf1dd8SToby Isaac } 89020cf1dd8SToby Isaac if (g2_func) { 89120cf1dd8SToby Isaac PetscInt d, d2; 89220cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 89320cf1dd8SToby Isaac g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer); 89420cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 89520cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 89620cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 89720cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] = 0.0; 89820cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 89920cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] *= w; 90020cf1dd8SToby Isaac } 90120cf1dd8SToby Isaac } 90220cf1dd8SToby Isaac } 90320cf1dd8SToby Isaac } 90420cf1dd8SToby Isaac if (g3_func) { 90520cf1dd8SToby Isaac PetscInt d, d2, dp, d3; 90620cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 90720cf1dd8SToby Isaac g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer); 90820cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 90920cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 91020cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 91120cf1dd8SToby Isaac for (dp = 0; dp < dim; ++dp) { 91220cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0; 91320cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 91420cf1dd8SToby Isaac for (d3 = 0; d3 < dim; ++d3) { 91520cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3]; 91620cf1dd8SToby Isaac } 91720cf1dd8SToby Isaac } 91820cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w; 91920cf1dd8SToby Isaac } 92020cf1dd8SToby Isaac } 92120cf1dd8SToby Isaac } 92220cf1dd8SToby Isaac } 92320cf1dd8SToby Isaac } 92420cf1dd8SToby Isaac 92520cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 92620cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 92720cf1dd8SToby Isaac const PetscInt fidx = f*NcI+fc; /* Test function basis index */ 92820cf1dd8SToby Isaac const PetscInt i = offsetI+f; /* Element matrix row */ 92920cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 93020cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 93120cf1dd8SToby Isaac const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */ 93220cf1dd8SToby Isaac const PetscInt j = offsetJ+g; /* Element matrix column */ 93320cf1dd8SToby Isaac const PetscInt fOff = eOffset+i*totDim+j; 93420cf1dd8SToby Isaac PetscInt d, d2; 93520cf1dd8SToby Isaac 93620cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx]; 93720cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 93820cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d]; 93920cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx]; 94020cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 94120cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2]; 94220cf1dd8SToby Isaac } 94320cf1dd8SToby Isaac } 94420cf1dd8SToby Isaac } 94520cf1dd8SToby Isaac } 94620cf1dd8SToby Isaac } 94720cf1dd8SToby Isaac } 94820cf1dd8SToby Isaac } 94920cf1dd8SToby Isaac if (debug > 1) { 95020cf1dd8SToby Isaac PetscInt fc, f, gc, g; 95120cf1dd8SToby Isaac 95220cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 95320cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 95420cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 95520cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 95620cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 95720cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 95820cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 95920cf1dd8SToby 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); 96020cf1dd8SToby Isaac } 96120cf1dd8SToby Isaac } 96220cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 96320cf1dd8SToby Isaac } 96420cf1dd8SToby Isaac } 96520cf1dd8SToby Isaac } 96620cf1dd8SToby Isaac cOffset += totDim; 96720cf1dd8SToby Isaac cOffsetAux += totDimAux; 96820cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 96920cf1dd8SToby Isaac } 97020cf1dd8SToby Isaac PetscFunctionReturn(0); 97120cf1dd8SToby Isaac } 97220cf1dd8SToby Isaac 97320cf1dd8SToby Isaac PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 97420cf1dd8SToby Isaac { 97520cf1dd8SToby Isaac PetscFunctionBegin; 97620cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 97720cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 97820cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 97920cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 98020cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 98120cf1dd8SToby Isaac fem->ops->gettabulation = PetscFEGetTabulation_Basic; 98220cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 983afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 98420cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 98520cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 98620cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 98720cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 98820cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 98920cf1dd8SToby Isaac PetscFunctionReturn(0); 99020cf1dd8SToby Isaac } 99120cf1dd8SToby Isaac 99220cf1dd8SToby Isaac /*MC 99320cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 99420cf1dd8SToby Isaac 99520cf1dd8SToby Isaac Level: intermediate 99620cf1dd8SToby Isaac 99720cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 99820cf1dd8SToby Isaac M*/ 99920cf1dd8SToby Isaac 100020cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 100120cf1dd8SToby Isaac { 100220cf1dd8SToby Isaac PetscFE_Basic *b; 100320cf1dd8SToby Isaac PetscErrorCode ierr; 100420cf1dd8SToby Isaac 100520cf1dd8SToby Isaac PetscFunctionBegin; 100620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 100720cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 100820cf1dd8SToby Isaac fem->data = b; 100920cf1dd8SToby Isaac 101020cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 101120cf1dd8SToby Isaac PetscFunctionReturn(0); 101220cf1dd8SToby Isaac } 101320cf1dd8SToby Isaac 101420cf1dd8SToby Isaac 1015