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 */ 49526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem) 5020cf1dd8SToby Isaac { 51b9d4cb8dSJed Brown PetscReal *work; 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 for (j = 0; j < pdim; ++j) { 6120cf1dd8SToby Isaac PetscReal *Bf; 6220cf1dd8SToby Isaac PetscQuadrature f; 6320cf1dd8SToby Isaac const PetscReal *points, *weights; 6420cf1dd8SToby Isaac PetscInt Nc, Nq, q, k, c; 6520cf1dd8SToby Isaac 6620cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr); 6720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 6820cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*Nq*pdim,&Bf);CHKERRQ(ierr); 6920cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr); 7020cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 7120cf1dd8SToby Isaac /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */ 72b9d4cb8dSJed Brown fem->invV[j*pdim+k] = 0.0; 7320cf1dd8SToby Isaac 7420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 75b9d4cb8dSJed Brown for (c = 0; c < Nc; ++c) fem->invV[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c]; 7620cf1dd8SToby Isaac } 7720cf1dd8SToby Isaac } 7820cf1dd8SToby Isaac ierr = PetscFree(Bf);CHKERRQ(ierr); 7920cf1dd8SToby Isaac } 80ea2bdf6dSBarry Smith 8120cf1dd8SToby Isaac ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr); 8220cf1dd8SToby Isaac n = pdim; 83b9d4cb8dSJed Brown PetscStackCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info)); 84ea2bdf6dSBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info); 85b9d4cb8dSJed Brown PetscStackCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info)); 86ea2bdf6dSBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info); 8720cf1dd8SToby Isaac ierr = PetscFree2(pivots,work);CHKERRQ(ierr); 8820cf1dd8SToby Isaac PetscFunctionReturn(0); 8920cf1dd8SToby Isaac } 9020cf1dd8SToby Isaac 9120cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim) 9220cf1dd8SToby Isaac { 9320cf1dd8SToby Isaac PetscErrorCode ierr; 9420cf1dd8SToby Isaac 9520cf1dd8SToby Isaac PetscFunctionBegin; 9620cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr); 9720cf1dd8SToby Isaac PetscFunctionReturn(0); 9820cf1dd8SToby Isaac } 9920cf1dd8SToby Isaac 100b9d4cb8dSJed Brown /* Tensor contraction on the middle index, 101b9d4cb8dSJed Brown * C[m,n,p] = A[m,k,p] * B[k,n] 102b9d4cb8dSJed Brown * where all matrices use C-style ordering. 103b9d4cb8dSJed Brown */ 104bdb10af2SPierre Jolivet static PetscErrorCode TensorContract_Private(PetscInt m,PetscInt n,PetscInt p,PetscInt k,const PetscReal *A,const PetscReal *B,PetscReal *C) 105bdb10af2SPierre Jolivet { 106b9d4cb8dSJed Brown PetscErrorCode ierr; 107b9d4cb8dSJed Brown PetscInt i; 108b9d4cb8dSJed Brown 109b9d4cb8dSJed Brown PetscFunctionBegin; 110b9d4cb8dSJed Brown for (i=0; i<m; i++) { 111b9d4cb8dSJed Brown PetscBLASInt n_,p_,k_,lda,ldb,ldc; 112b9d4cb8dSJed Brown PetscReal one = 1, zero = 0; 113b9d4cb8dSJed Brown /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n] 114b9d4cb8dSJed Brown * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k) 115b9d4cb8dSJed Brown */ 116b9d4cb8dSJed Brown ierr = PetscBLASIntCast(n,&n_);CHKERRQ(ierr); 117b9d4cb8dSJed Brown ierr = PetscBLASIntCast(p,&p_);CHKERRQ(ierr); 118b9d4cb8dSJed Brown ierr = PetscBLASIntCast(k,&k_);CHKERRQ(ierr); 119b9d4cb8dSJed Brown lda = p_; 120b9d4cb8dSJed Brown ldb = n_; 121b9d4cb8dSJed Brown ldc = p_; 122b9d4cb8dSJed Brown PetscStackCallBLAS("BLASgemm",BLASREALgemm_("N","T",&p_,&n_,&k_,&one,A+i*k*p,&lda,B,&ldb,&zero,C+i*n*p,&ldc)); 123b9d4cb8dSJed Brown } 1246485a3bbSJed Brown ierr = PetscLogFlops(2.*m*n*p*k);CHKERRQ(ierr); 125b9d4cb8dSJed Brown PetscFunctionReturn(0); 126b9d4cb8dSJed Brown } 127b9d4cb8dSJed Brown 128526996e7SStefano Zampini PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 12920cf1dd8SToby Isaac { 13020cf1dd8SToby Isaac DM dm; 13120cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 13220cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 13320cf1dd8SToby Isaac PetscInt Nc; /* Field components */ 134ef0bb6c7SMatthew G. Knepley PetscReal *B = K >= 0 ? T->T[0] : NULL; 135ef0bb6c7SMatthew G. Knepley PetscReal *D = K >= 1 ? T->T[1] : NULL; 136ef0bb6c7SMatthew G. Knepley PetscReal *H = K >= 2 ? T->T[2] : NULL; 137ef0bb6c7SMatthew G. Knepley PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL; 13820cf1dd8SToby Isaac PetscErrorCode ierr; 13920cf1dd8SToby Isaac 14020cf1dd8SToby Isaac PetscFunctionBegin; 14120cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 14220cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 14320cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 14420cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 14520cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 146ef0bb6c7SMatthew G. Knepley if (K >= 0) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 147ef0bb6c7SMatthew G. Knepley if (K >= 1) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 148ef0bb6c7SMatthew G. Knepley if (K >= 2) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 149ef0bb6c7SMatthew G. Knepley ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH);CHKERRQ(ierr); 150b9d4cb8dSJed Brown /* Translate from prime to nodal basis */ 15120cf1dd8SToby Isaac if (B) { 152b9d4cb8dSJed Brown /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */ 153b9d4cb8dSJed Brown ierr = TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B);CHKERRQ(ierr); 15420cf1dd8SToby Isaac } 15520cf1dd8SToby Isaac if (D) { 156b9d4cb8dSJed Brown /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */ 157b9d4cb8dSJed Brown ierr = TensorContract_Private(npoints, pdim, Nc*dim, pdim, tmpD, fem->invV, D);CHKERRQ(ierr); 15820cf1dd8SToby Isaac } 15920cf1dd8SToby Isaac if (H) { 160b9d4cb8dSJed Brown /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */ 161b9d4cb8dSJed Brown ierr = TensorContract_Private(npoints, pdim, Nc*dim*dim, pdim, tmpH, fem->invV, H);CHKERRQ(ierr); 16220cf1dd8SToby Isaac } 163ef0bb6c7SMatthew G. Knepley if (K >= 0) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 164ef0bb6c7SMatthew G. Knepley if (K >= 1) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 165ef0bb6c7SMatthew G. Knepley if (K >= 2) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 16620cf1dd8SToby Isaac PetscFunctionReturn(0); 16720cf1dd8SToby Isaac } 16820cf1dd8SToby Isaac 1692b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1704bee2e38SMatthew G. Knepley const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 17120cf1dd8SToby Isaac { 17220cf1dd8SToby Isaac const PetscInt debug = 0; 1734bee2e38SMatthew G. Knepley PetscFE fe; 17420cf1dd8SToby Isaac PetscPointFunc obj_func; 17520cf1dd8SToby Isaac PetscQuadrature quad; 176ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 1774bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x; 17820cf1dd8SToby Isaac const PetscScalar *constants; 17920cf1dd8SToby Isaac PetscReal *x; 180ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 18120cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 18220cf1dd8SToby Isaac PetscBool isAffine; 18320cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 18420cf1dd8SToby Isaac PetscInt qNc, Nq, q; 18520cf1dd8SToby Isaac PetscErrorCode ierr; 18620cf1dd8SToby Isaac 18720cf1dd8SToby Isaac PetscFunctionBegin; 1884bee2e38SMatthew G. Knepley ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr); 18920cf1dd8SToby Isaac if (!obj_func) PetscFunctionReturn(0); 1904bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 1914bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 1924bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1934bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 1944bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 1954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 1964bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 197ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1984bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 1994bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 2004bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 2014bee2e38SMatthew G. Knepley if (dsAux) { 2024bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 2034bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 2044bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 2054bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 206ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 2074bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 20801907d53SMatthew G. Knepley if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 20920cf1dd8SToby Isaac } 21020cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 2114bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 21220cf1dd8SToby Isaac Np = cgeom->numPoints; 21320cf1dd8SToby Isaac dE = cgeom->dimEmbed; 21420cf1dd8SToby Isaac isAffine = cgeom->isAffine; 21520cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 2164bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 21720cf1dd8SToby Isaac 21827f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 21927f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 22020cf1dd8SToby Isaac if (isAffine) { 2214bee2e38SMatthew G. Knepley fegeom.v = x; 2224bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2237132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 2247132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 2257132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 22620cf1dd8SToby Isaac } 2274bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 2284bee2e38SMatthew G. Knepley PetscScalar integrand; 2294bee2e38SMatthew G. Knepley PetscReal w; 2304bee2e38SMatthew G. Knepley 2314bee2e38SMatthew G. Knepley if (isAffine) { 2327132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 2334bee2e38SMatthew G. Knepley } else { 2344bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 2354bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 2364bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 2374bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 2384bee2e38SMatthew G. Knepley } 2394bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 24020cf1dd8SToby Isaac if (debug > 1 && q < Np) { 2414bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 2427be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 2434bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 24420cf1dd8SToby Isaac #endif 24520cf1dd8SToby Isaac } 24620cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 247ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 248ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 2494bee2e38SMatthew 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); 2504bee2e38SMatthew G. Knepley integrand *= w; 25120cf1dd8SToby Isaac integral[e*Nf+field] += integrand; 25220cf1dd8SToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);} 25320cf1dd8SToby Isaac } 25420cf1dd8SToby Isaac cOffset += totDim; 25520cf1dd8SToby Isaac cOffsetAux += totDimAux; 25620cf1dd8SToby Isaac } 25720cf1dd8SToby Isaac PetscFunctionReturn(0); 25820cf1dd8SToby Isaac } 25920cf1dd8SToby Isaac 2602b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, 261afe6d6adSToby Isaac PetscBdPointFunc obj_func, 2624bee2e38SMatthew G. Knepley PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 263afe6d6adSToby Isaac { 264afe6d6adSToby Isaac const PetscInt debug = 0; 2654bee2e38SMatthew G. Knepley PetscFE fe; 266afe6d6adSToby Isaac PetscQuadrature quad; 267ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 2684bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 269afe6d6adSToby Isaac const PetscScalar *constants; 270afe6d6adSToby Isaac PetscReal *x; 271ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 272afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 273afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 274afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 275afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 276afe6d6adSToby Isaac PetscErrorCode ierr; 277afe6d6adSToby Isaac 278afe6d6adSToby Isaac PetscFunctionBegin; 279afe6d6adSToby Isaac if (!obj_func) PetscFunctionReturn(0); 2804bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 2814bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 2824bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 2834bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 2844bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 2854bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 2864bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 2874bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 2884bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 289ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr); 2904bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 2914bee2e38SMatthew G. Knepley if (dsAux) { 2924bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 2934bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 2944bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 2954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 2964bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 2974bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 298afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 299ef0bb6c7SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 300ef0bb6c7SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 30101907d53SMatthew G. Knepley if (Tf[0]->Np != TfAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 302afe6d6adSToby Isaac } 303afe6d6adSToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 304afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 305afe6d6adSToby Isaac Np = fgeom->numPoints; 306afe6d6adSToby Isaac dE = fgeom->dimEmbed; 307afe6d6adSToby Isaac isAffine = fgeom->isAffine; 308afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 3099f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 310afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 311ea78f98cSLisandro Dalcin fegeom.n = NULL; 312ea78f98cSLisandro Dalcin fegeom.v = NULL; 313ea78f98cSLisandro Dalcin fegeom.J = NULL; 314ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 31527f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 31627f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 31727f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 31827f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3194bee2e38SMatthew G. Knepley if (isAffine) { 3204bee2e38SMatthew G. Knepley fegeom.v = x; 3214bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3227132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 3237132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 3247132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 3257132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 3269f209ee4SMatthew G. Knepley 3277132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 3287132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 3297132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 3304bee2e38SMatthew G. Knepley } 331afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 332afe6d6adSToby Isaac PetscScalar integrand; 3334bee2e38SMatthew G. Knepley PetscReal w; 334afe6d6adSToby Isaac 335afe6d6adSToby Isaac if (isAffine) { 3367132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 337afe6d6adSToby Isaac } else { 3383fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 3399f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 3409f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 3414bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 3424bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 3439f209ee4SMatthew G. Knepley 3449f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 3459f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 3469f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 347afe6d6adSToby Isaac } 3484bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 349afe6d6adSToby Isaac if (debug > 1 && q < Np) { 3504bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 351afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3524bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 353afe6d6adSToby Isaac #endif 354afe6d6adSToby Isaac } 355afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 356ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 357ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 3584bee2e38SMatthew 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); 3594bee2e38SMatthew G. Knepley integrand *= w; 360afe6d6adSToby Isaac integral[e*Nf+field] += integrand; 361afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);} 362afe6d6adSToby Isaac } 363afe6d6adSToby Isaac cOffset += totDim; 364afe6d6adSToby Isaac cOffsetAux += totDimAux; 365afe6d6adSToby Isaac } 366afe6d6adSToby Isaac PetscFunctionReturn(0); 367afe6d6adSToby Isaac } 368afe6d6adSToby Isaac 36906ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, 3704bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 37120cf1dd8SToby Isaac { 37220cf1dd8SToby Isaac const PetscInt debug = 0; 3736528b96dSMatthew G. Knepley const PetscInt field = key.field; 3744bee2e38SMatthew G. Knepley PetscFE fe; 3756528b96dSMatthew G. Knepley PetscWeakForm wf; 3766528b96dSMatthew G. Knepley PetscInt n0, n1, i; 3776528b96dSMatthew G. Knepley PetscPointFunc *f0_func, *f1_func; 37820cf1dd8SToby Isaac PetscQuadrature quad; 379ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 3804bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 38120cf1dd8SToby Isaac const PetscScalar *constants; 38220cf1dd8SToby Isaac PetscReal *x; 383ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 384ef0bb6c7SMatthew G. Knepley PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 38520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3866587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 38720cf1dd8SToby Isaac PetscErrorCode ierr; 38820cf1dd8SToby Isaac 38920cf1dd8SToby Isaac PetscFunctionBegin; 3904bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 3914bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 3924bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 3934bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 3944bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 3954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 3964bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 3974bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 3986528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 39906ad1575SMatthew G. Knepley ierr = PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 4006528b96dSMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 4014bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 4024bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 4034bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 404ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 4054bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 4064bee2e38SMatthew G. Knepley if (dsAux) { 4074bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 4084bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 4094bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 4104bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 4114bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 412ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 41301907d53SMatthew G. Knepley if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 41420cf1dd8SToby Isaac } 4156587ee25SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 4164bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 41720cf1dd8SToby Isaac dE = cgeom->dimEmbed; 4186587ee25SMatthew G. Knepley if (cgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", cgeom->dim, qdim); 41920cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4204bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 42120cf1dd8SToby Isaac 4226587ee25SMatthew G. Knepley fegeom.v = x; /* workspace */ 423ef0bb6c7SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr); 42427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dE);CHKERRQ(ierr); 42520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4264bee2e38SMatthew G. Knepley PetscReal w; 4274bee2e38SMatthew G. Knepley PetscInt c, d; 42820cf1dd8SToby Isaac 4296587ee25SMatthew G. Knepley ierr = PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q*cgeom->dim], &fegeom);CHKERRQ(ierr); 4304bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 4316587ee25SMatthew G. Knepley if (debug > 1 && q < cgeom->numPoints) { 4324bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 4337be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4344bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 43520cf1dd8SToby Isaac #endif 43620cf1dd8SToby Isaac } 437ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 438ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 4396528b96dSMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q*T[field]->Nc]); 440ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w; 4416528b96dSMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q*T[field]->Nc*dim]); 442ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) for (d = 0; d < dim; ++d) f1[(q*T[field]->Nc+c)*dim+d] *= w; 443b8025e53SMatthew G. Knepley if (debug) { 444b8025e53SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d wt %g\n", q, quadWeights[q]);CHKERRQ(ierr); 445b8025e53SMatthew G. Knepley if (debug > 2) { 446b8025e53SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " field %d:", field);CHKERRQ(ierr); 447b8025e53SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %g", u[uOff[field]+c]);CHKERRQ(ierr);} 448b8025e53SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 449b8025e53SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " resid %d:", field);CHKERRQ(ierr); 450b8025e53SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %g", f0[q*T[field]->Nc+c]);CHKERRQ(ierr);} 451b8025e53SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 452b8025e53SMatthew G. Knepley } 453b8025e53SMatthew G. Knepley } 45420cf1dd8SToby Isaac } 4556587ee25SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 45620cf1dd8SToby Isaac cOffset += totDim; 45720cf1dd8SToby Isaac cOffsetAux += totDimAux; 45820cf1dd8SToby Isaac } 45920cf1dd8SToby Isaac PetscFunctionReturn(0); 46020cf1dd8SToby Isaac } 46120cf1dd8SToby Isaac 46206ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 4634bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 46420cf1dd8SToby Isaac { 46520cf1dd8SToby Isaac const PetscInt debug = 0; 46606d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 4674bee2e38SMatthew G. Knepley PetscFE fe; 46806d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 46906d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 47020cf1dd8SToby Isaac PetscQuadrature quad; 471ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4724bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 47320cf1dd8SToby Isaac const PetscScalar *constants; 47420cf1dd8SToby Isaac PetscReal *x; 475ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 476ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 4776587ee25SMatthew G. Knepley PetscBool auxOnBd = PETSC_FALSE; 47820cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 4796587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 48020cf1dd8SToby Isaac PetscErrorCode ierr; 48120cf1dd8SToby Isaac 48220cf1dd8SToby Isaac PetscFunctionBegin; 4834bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 4844bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 4854bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 4864bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 4874bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 4884bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 4894bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 4904bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 49106ad1575SMatthew G. Knepley ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 49206d8a0d3SMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 4934bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 4944bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 4954bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 496ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr); 4974bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 4984bee2e38SMatthew G. Knepley if (dsAux) { 4994bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 5004bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 5014bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 5024bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 5034bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 5044bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 5057be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 506ef0bb6c7SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 507ef0bb6c7SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 50801907d53SMatthew G. Knepley if (Tf[0]->Np != TfAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 50920cf1dd8SToby Isaac } 510ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 5116587ee25SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 512afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 51320cf1dd8SToby Isaac dE = fgeom->dimEmbed; 5146587ee25SMatthew G. Knepley /* TODO FIX THIS */ 5156587ee25SMatthew G. Knepley fgeom->dim = dim-1; 5166587ee25SMatthew G. Knepley if (fgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", fgeom->dim, qdim); 51720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5189f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 51920cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5209f209ee4SMatthew G. Knepley 5216587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 522580bdb30SBarry Smith ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr); 52327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcI*dE);CHKERRQ(ierr); 52420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5254bee2e38SMatthew G. Knepley PetscReal w; 5264bee2e38SMatthew G. Knepley PetscInt c, d; 5274bee2e38SMatthew G. Knepley 5286587ee25SMatthew G. Knepley ierr = PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom);CHKERRQ(ierr); 5296587ee25SMatthew G. Knepley ierr = PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom);CHKERRQ(ierr); 5304bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 53162bd480fSMatthew G. Knepley if (debug > 1) { 5326587ee25SMatthew G. Knepley if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 5334bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 5347be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5354bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 53662bd480fSMatthew G. Knepley ierr = DMPrintCellVector(e, "n", dim, fegeom.n);CHKERRQ(ierr); 53720cf1dd8SToby Isaac #endif 53820cf1dd8SToby Isaac } 53962bd480fSMatthew G. Knepley } 540ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 541ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 54206d8a0d3SMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](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]); 5434bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 54406d8a0d3SMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](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]); 5454bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 54662bd480fSMatthew G. Knepley if (debug) { 54762bd480fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " elem %D quad point %d\n", e, q);CHKERRQ(ierr); 54862bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 54962bd480fSMatthew G. Knepley if (n0) {ierr = PetscPrintf(PETSC_COMM_SELF, " f0[%D] %g\n", c, f0[q*NcI+c]);CHKERRQ(ierr);} 55062bd480fSMatthew G. Knepley if (n1) { 55162bd480fSMatthew G. Knepley for (d = 0; d < dim; ++d) {ierr = PetscPrintf(PETSC_COMM_SELF, " f1[%D,%D] %g", c, d, f1[(q*NcI + c)*dim + d]);CHKERRQ(ierr);} 55262bd480fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 55362bd480fSMatthew G. Knepley } 55462bd480fSMatthew G. Knepley } 55562bd480fSMatthew G. Knepley } 55620cf1dd8SToby Isaac } 5576587ee25SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 55820cf1dd8SToby Isaac cOffset += totDim; 55920cf1dd8SToby Isaac cOffsetAux += totDimAux; 56020cf1dd8SToby Isaac } 56120cf1dd8SToby Isaac PetscFunctionReturn(0); 56220cf1dd8SToby Isaac } 56320cf1dd8SToby Isaac 56427f02ce8SMatthew G. Knepley /* 56527f02ce8SMatthew G. Knepley BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but 56627f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 56727f02ce8SMatthew G. Knepley 56827f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 56927f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 57027f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 57127f02ce8SMatthew G. Knepley 3) TODO We need to use a non-square Jacobian for the derivative maps, meaning the embedding dimension has to go to EvaluateFieldJets() and UpdateElementVec() 57227f02ce8SMatthew G. Knepley */ 573c2b7495fSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, 57427f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 57527f02ce8SMatthew G. Knepley { 57627f02ce8SMatthew G. Knepley const PetscInt debug = 0; 5776528b96dSMatthew G. Knepley const PetscInt field = key.field; 57827f02ce8SMatthew G. Knepley PetscFE fe; 5796528b96dSMatthew G. Knepley PetscWeakForm wf; 5806528b96dSMatthew G. Knepley PetscInt n0, n1, i; 5816528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 58227f02ce8SMatthew G. Knepley PetscQuadrature quad; 583665f567fSMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 58427f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 58527f02ce8SMatthew G. Knepley const PetscScalar *constants; 58627f02ce8SMatthew G. Knepley PetscReal *x; 587665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 588665f567fSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 5896587ee25SMatthew G. Knepley PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 59027f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 5916587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 59227f02ce8SMatthew G. Knepley PetscErrorCode ierr; 59327f02ce8SMatthew G. Knepley 59427f02ce8SMatthew G. Knepley PetscFunctionBegin; 595c2b7495fSMatthew G. Knepley if (s < 0 || s > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The side %D must be in [0, 1]", s); 59627f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 59727f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 59827f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 59927f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 60027f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 60127f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 60227f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 60327f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 604*5fedec97SMatthew G. Knepley ierr = PetscDSGetFieldOffsetCohesive(ds, field, &fOffset);CHKERRQ(ierr); 6056528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 60606ad1575SMatthew G. Knepley ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 6076528b96dSMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 60827f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 60927f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 61027f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 61127f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 612665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr); 61327f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 61427f02ce8SMatthew G. Knepley if (dsAux) { 61527f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 61627f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 61727f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 61827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 61927f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 62027f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 62101907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 622665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 623665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 62401907d53SMatthew G. Knepley if (Tf[0]->Np != TfAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 62527f02ce8SMatthew G. Knepley } 626*5fedec97SMatthew G. Knepley ierr = PetscDSGetCohesive(ds, field, &isCohesiveField);CHKERRQ(ierr); 627665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 628c2b7495fSMatthew G. Knepley NcS = NcI; 6296587ee25SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 63027f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 63127f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 6326587ee25SMatthew G. Knepley if (fgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", fgeom->dim, qdim); 63327f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 63427f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 63527f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 63627f02ce8SMatthew G. Knepley 6376587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 63827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr); 63927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr); 64027f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 64127f02ce8SMatthew G. Knepley PetscReal w; 64227f02ce8SMatthew G. Knepley PetscInt c, d; 64327f02ce8SMatthew G. Knepley 6446587ee25SMatthew G. Knepley ierr = PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom);CHKERRQ(ierr); 64527f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 6466587ee25SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 64727f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 64827f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 64927f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr); 65027f02ce8SMatthew G. Knepley #endif 65127f02ce8SMatthew G. Knepley } 65227f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 65327f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 654665f567fSMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 655*5fedec97SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 6566528b96dSMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](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*NcS]); 65727f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w; 6586528b96dSMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](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*NcS*dim]); 65927f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w; 66027f02ce8SMatthew G. Knepley } 661*5fedec97SMatthew G. Knepley if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset]);} 662*5fedec97SMatthew G. Knepley else {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset+fOffset]);} 66327f02ce8SMatthew G. Knepley cOffset += totDim; 66427f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 66527f02ce8SMatthew G. Knepley } 66627f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 66727f02ce8SMatthew G. Knepley } 66827f02ce8SMatthew G. Knepley 66906ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, 6704bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 67120cf1dd8SToby Isaac { 67220cf1dd8SToby Isaac const PetscInt debug = 0; 6734bee2e38SMatthew G. Knepley PetscFE feI, feJ; 6746528b96dSMatthew G. Knepley PetscWeakForm wf; 6756528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 6766528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 67720cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 67820cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 67920cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 68020cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 68120cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 68220cf1dd8SToby Isaac PetscQuadrature quad; 683ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 6844bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 68520cf1dd8SToby Isaac const PetscScalar *constants; 68620cf1dd8SToby Isaac PetscReal *x; 687ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 688ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 6896528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 69020cf1dd8SToby Isaac PetscInt dE, Np; 69120cf1dd8SToby Isaac PetscBool isAffine; 69220cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 69320cf1dd8SToby Isaac PetscInt qNc, Nq, q; 69420cf1dd8SToby Isaac PetscErrorCode ierr; 69520cf1dd8SToby Isaac 69620cf1dd8SToby Isaac PetscFunctionBegin; 6976528b96dSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 6986528b96dSMatthew G. Knepley fieldI = key.field / Nf; 6996528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 7004bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 7014bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 7024bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 7034bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 7044bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 7054bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 7064bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 7076528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 70820cf1dd8SToby Isaac switch(jtype) { 70906ad1575SMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: ierr = PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 71006ad1575SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 71106ad1575SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 71220cf1dd8SToby Isaac } 7136528b96dSMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 7144bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 7154bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 7164bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 717ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 7184bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 7194bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 7204bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 7214bee2e38SMatthew G. Knepley if (dsAux) { 7224bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 7234bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 7244bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 7254bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 7264bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 727ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 72801907d53SMatthew G. Knepley if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 72920cf1dd8SToby Isaac } 73027f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 73127f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7324bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7334bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7344bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 73527f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 73627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 73727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 73827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 73927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 74027f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 74127f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 7424bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7434bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7444bee2e38SMatthew G. Knepley 74527f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 74627f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7474bee2e38SMatthew G. Knepley if (isAffine) { 7484bee2e38SMatthew G. Knepley fegeom.v = x; 7494bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 7507132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 7517132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 7527132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 7534bee2e38SMatthew G. Knepley } 75420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 75520cf1dd8SToby Isaac PetscReal w; 7564bee2e38SMatthew G. Knepley PetscInt c; 75720cf1dd8SToby Isaac 75820cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 75920cf1dd8SToby Isaac if (isAffine) { 7607132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 76120cf1dd8SToby Isaac } else { 7624bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 7634bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 7644bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 7654bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 76620cf1dd8SToby Isaac } 7674bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 768ef0bb6c7SMatthew G. Knepley if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 769ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 770ea672e62SMatthew G. Knepley if (n0) { 771580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 7726528b96dSMatthew G. Knepley for (i = 0; i < n0; ++i) g0_func[i](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); 77320cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 77420cf1dd8SToby Isaac } 775ea672e62SMatthew G. Knepley if (n1) { 77627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 7776528b96dSMatthew G. Knepley for (i = 0; i < n1; ++i) g1_func[i](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); 7784bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 77920cf1dd8SToby Isaac } 780ea672e62SMatthew G. Knepley if (n2) { 78127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 7826528b96dSMatthew G. Knepley for (i = 0; i < n2; ++i) g2_func[i](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); 7834bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 78420cf1dd8SToby Isaac } 785ea672e62SMatthew G. Knepley if (n3) { 78627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 7876528b96dSMatthew G. Knepley for (i = 0; i < n3; ++i) g3_func[i](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); 7884bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 78920cf1dd8SToby Isaac } 79020cf1dd8SToby Isaac 791ef0bb6c7SMatthew G. Knepley ierr = PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr); 79220cf1dd8SToby Isaac } 79320cf1dd8SToby Isaac if (debug > 1) { 79420cf1dd8SToby Isaac PetscInt fc, f, gc, g; 79520cf1dd8SToby Isaac 79620cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 797ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 798ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 799ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 800ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 801ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 802ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 80320cf1dd8SToby 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); 80420cf1dd8SToby Isaac } 80520cf1dd8SToby Isaac } 80620cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 80720cf1dd8SToby Isaac } 80820cf1dd8SToby Isaac } 80920cf1dd8SToby Isaac } 81020cf1dd8SToby Isaac cOffset += totDim; 81120cf1dd8SToby Isaac cOffsetAux += totDimAux; 81220cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 81320cf1dd8SToby Isaac } 81420cf1dd8SToby Isaac PetscFunctionReturn(0); 81520cf1dd8SToby Isaac } 81620cf1dd8SToby Isaac 81706ad1575SMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 8184bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 81920cf1dd8SToby Isaac { 82020cf1dd8SToby Isaac const PetscInt debug = 0; 8214bee2e38SMatthew G. Knepley PetscFE feI, feJ; 82245480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 82345480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 82420cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 82520cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 82620cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 82720cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 82820cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 82920cf1dd8SToby Isaac PetscQuadrature quad; 830ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8314bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 83220cf1dd8SToby Isaac const PetscScalar *constants; 83320cf1dd8SToby Isaac PetscReal *x; 834ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 835ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 83645480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 83720cf1dd8SToby Isaac PetscBool isAffine; 83820cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 83920cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 84020cf1dd8SToby Isaac PetscErrorCode ierr; 84120cf1dd8SToby Isaac 84220cf1dd8SToby Isaac PetscFunctionBegin; 84345480ffeSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 84445480ffeSMatthew G. Knepley fieldI = key.field / Nf; 84545480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8464bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 8474bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 8484bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 8494bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 8504bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 8514bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 8524bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 8534bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 8544bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 85506ad1575SMatthew G. Knepley ierr = PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr); 85645480ffeSMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 8574bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 8584bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 8594bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 860ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr); 8614bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 8624bee2e38SMatthew G. Knepley if (dsAux) { 8634bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 8644bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 8654bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 8664bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 8674bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 868ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr); 86920cf1dd8SToby Isaac } 870ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 87120cf1dd8SToby Isaac Np = fgeom->numPoints; 87220cf1dd8SToby Isaac dE = fgeom->dimEmbed; 87320cf1dd8SToby Isaac isAffine = fgeom->isAffine; 87427f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 87531f073d7SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 87627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 87727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 87827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 87927f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 88027f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 88120cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 8829f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 88320cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 884ea78f98cSLisandro Dalcin fegeom.n = NULL; 885ea78f98cSLisandro Dalcin fegeom.v = NULL; 886ea78f98cSLisandro Dalcin fegeom.J = NULL; 887ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 88827f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 88927f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 89027f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 89127f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 8924bee2e38SMatthew G. Knepley if (isAffine) { 8934bee2e38SMatthew G. Knepley fegeom.v = x; 8944bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 8957132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 8967132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 8977132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 8987132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 8999f209ee4SMatthew G. Knepley 9007132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 9017132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 9027132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 9034bee2e38SMatthew G. Knepley } 90420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 90520cf1dd8SToby Isaac PetscReal w; 9064bee2e38SMatthew G. Knepley PetscInt c; 90720cf1dd8SToby Isaac 90820cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 90920cf1dd8SToby Isaac if (isAffine) { 9107132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 91120cf1dd8SToby Isaac } else { 9123fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 9139f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 9149f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 9154bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 9164bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 9179f209ee4SMatthew G. Knepley 9189f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 9199f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 9209f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 92120cf1dd8SToby Isaac } 9224bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 923ef0bb6c7SMatthew G. Knepley if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 924ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 925ea672e62SMatthew G. Knepley if (n0) { 926580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 92745480ffeSMatthew G. Knepley for (i = 0; i < n0; ++i) g0_func[i](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); 92820cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 92920cf1dd8SToby Isaac } 930ea672e62SMatthew G. Knepley if (n1) { 93127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 93245480ffeSMatthew G. Knepley for (i = 0; i < n1; ++i) g1_func[i](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); 9334bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 93420cf1dd8SToby Isaac } 935ea672e62SMatthew G. Knepley if (n2) { 93627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 93745480ffeSMatthew G. Knepley for (i = 0; i < n2; ++i) g2_func[i](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); 9384bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 93920cf1dd8SToby Isaac } 940ea672e62SMatthew G. Knepley if (n3) { 94127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 94245480ffeSMatthew G. Knepley for (i = 0; i < n3; ++i) g3_func[i](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); 9434bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 94420cf1dd8SToby Isaac } 94520cf1dd8SToby Isaac 946ef0bb6c7SMatthew G. Knepley ierr = PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr); 94720cf1dd8SToby Isaac } 94820cf1dd8SToby Isaac if (debug > 1) { 94920cf1dd8SToby Isaac PetscInt fc, f, gc, g; 95020cf1dd8SToby Isaac 95120cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 952ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 953ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 954ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 955ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 956ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 957ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 95820cf1dd8SToby 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); 95920cf1dd8SToby Isaac } 96020cf1dd8SToby Isaac } 96120cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 96220cf1dd8SToby Isaac } 96320cf1dd8SToby Isaac } 96420cf1dd8SToby Isaac } 96520cf1dd8SToby Isaac cOffset += totDim; 96620cf1dd8SToby Isaac cOffsetAux += totDimAux; 96720cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 96820cf1dd8SToby Isaac } 96920cf1dd8SToby Isaac PetscFunctionReturn(0); 97020cf1dd8SToby Isaac } 97120cf1dd8SToby Isaac 972*5fedec97SMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, 97327f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 97427f02ce8SMatthew G. Knepley { 97527f02ce8SMatthew G. Knepley const PetscInt debug = 0; 97627f02ce8SMatthew G. Knepley PetscFE feI, feJ; 977148442b3SMatthew G. Knepley PetscWeakForm wf; 978148442b3SMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 979148442b3SMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 98027f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 98127f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 98227f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 98327f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 98427f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 985665f567fSMatthew G. Knepley PetscQuadrature quad; 986665f567fSMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 98727f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 98827f02ce8SMatthew G. Knepley const PetscScalar *constants; 98927f02ce8SMatthew G. Knepley PetscReal *x; 990665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 991665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 99245480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 993665f567fSMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE; 99427f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 99527f02ce8SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 99627f02ce8SMatthew G. Knepley PetscErrorCode ierr; 99727f02ce8SMatthew G. Knepley 99827f02ce8SMatthew G. Knepley PetscFunctionBegin; 99945480ffeSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 100045480ffeSMatthew G. Knepley fieldI = key.field / Nf; 100145480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 100227f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 100327f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 100427f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 100527f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 100627f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 100727f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 100827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 100927f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 1010148442b3SMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 101127f02ce8SMatthew G. Knepley switch(jtype) { 101206ad1575SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 101306ad1575SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 1014665f567fSMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 101527f02ce8SMatthew G. Knepley } 1016148442b3SMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 101727f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 101827f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 101927f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 1020665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1021*5fedec97SMatthew G. Knepley ierr = PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI);CHKERRQ(ierr); 1022*5fedec97SMatthew G. Knepley ierr = PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 102327f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 102427f02ce8SMatthew G. Knepley if (dsAux) { 102527f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 102627f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 102727f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 102827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 102927f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 103027f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 103101907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 1032665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);} 1033665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);} 103401907d53SMatthew G. Knepley if (T[0]->Np != TAux[0]->Np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %D != %D number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 103527f02ce8SMatthew G. Knepley } 1036*5fedec97SMatthew G. Knepley ierr = PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI);CHKERRQ(ierr); 1037*5fedec97SMatthew G. Knepley ierr = PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ);CHKERRQ(ierr); 1038665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1039665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 104027f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2*NcI; 104127f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2*NcJ; 104227f02ce8SMatthew G. Knepley Np = fgeom->numPoints; 104327f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 104427f02ce8SMatthew G. Knepley isAffine = fgeom->isAffine; 104527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 104627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 104727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 104827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1049665f567fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1050665f567fSMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 105127f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 105227f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 105327f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 105427f02ce8SMatthew G. Knepley 105527f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 105627f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 105727f02ce8SMatthew G. Knepley if (isAffine) { 105827f02ce8SMatthew G. Knepley fegeom.v = x; 105927f02ce8SMatthew G. Knepley fegeom.xi = fgeom->xi; 106027f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 106127f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 106227f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 106327f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 106427f02ce8SMatthew G. Knepley } 106527f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 106627f02ce8SMatthew G. Knepley PetscReal w; 106727f02ce8SMatthew G. Knepley PetscInt c; 106827f02ce8SMatthew G. Knepley 106927f02ce8SMatthew G. Knepley if (isAffine) { 107027f02ce8SMatthew G. Knepley /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */ 107127f02ce8SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 107227f02ce8SMatthew G. Knepley } else { 107327f02ce8SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 107427f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 107527f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 107627f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 107727f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 107827f02ce8SMatthew G. Knepley } 107927f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 108027f02ce8SMatthew G. Knepley if (debug > 1 && q < Np) { 108127f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 108227f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 108327f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 108427f02ce8SMatthew G. Knepley #endif 108527f02ce8SMatthew G. Knepley } 108627f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 1087665f567fSMatthew G. Knepley if (coefficients) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 1088*5fedec97SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 1089ea672e62SMatthew G. Knepley if (n0) { 109027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 1091148442b3SMatthew G. Knepley for (i = 0; i < n0; ++i) g0_func[i](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); 109227f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT; ++c) g0[c] *= w; 109327f02ce8SMatthew G. Knepley } 1094ea672e62SMatthew G. Knepley if (n1) { 109527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 1096148442b3SMatthew G. Knepley for (i = 0; i < n1; ++i) g1_func[i](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); 109727f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w; 109827f02ce8SMatthew G. Knepley } 1099ea672e62SMatthew G. Knepley if (n2) { 110027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 1101148442b3SMatthew G. Knepley for (i = 0; i < n2; ++i) g2_func[i](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); 110227f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w; 110327f02ce8SMatthew G. Knepley } 1104ea672e62SMatthew G. Knepley if (n3) { 110527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1106148442b3SMatthew G. Knepley for (i = 0; i < n3; ++i) g3_func[i](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); 110727f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w; 110827f02ce8SMatthew G. Knepley } 110927f02ce8SMatthew G. Knepley 1110*5fedec97SMatthew G. Knepley if (isCohesiveFieldI) { 1111*5fedec97SMatthew G. Knepley if (isCohesiveFieldJ) { 1112*5fedec97SMatthew G. Knepley ierr = PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr); 111327f02ce8SMatthew G. Knepley } else { 1114*5fedec97SMatthew G. Knepley ierr = PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr); 1115*5fedec97SMatthew G. Knepley ierr = PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI*NcJ], &g1[NcI*NcJ*dim], &g2[NcI*NcJ*dim], &g3[NcI*NcJ*dim*dim], eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr); 1116*5fedec97SMatthew G. Knepley } 1117*5fedec97SMatthew G. Knepley } else { 1118*5fedec97SMatthew G. Knepley ierr = PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat);CHKERRQ(ierr); 111927f02ce8SMatthew G. Knepley } 112027f02ce8SMatthew G. Knepley } 112127f02ce8SMatthew G. Knepley if (debug > 1) { 112227f02ce8SMatthew G. Knepley PetscInt fc, f, gc, g; 112327f02ce8SMatthew G. Knepley 112427f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 112527f02ce8SMatthew G. Knepley for (fc = 0; fc < NcI; ++fc) { 1126665f567fSMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 112727f02ce8SMatthew G. Knepley const PetscInt i = offsetI + f*NcI+fc; 112827f02ce8SMatthew G. Knepley for (gc = 0; gc < NcJ; ++gc) { 1129665f567fSMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 113027f02ce8SMatthew G. Knepley const PetscInt j = offsetJ + g*NcJ+gc; 113127f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr); 113227f02ce8SMatthew G. Knepley } 113327f02ce8SMatthew G. Knepley } 113427f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 113527f02ce8SMatthew G. Knepley } 113627f02ce8SMatthew G. Knepley } 113727f02ce8SMatthew G. Knepley } 113827f02ce8SMatthew G. Knepley cOffset += totDim; 113927f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 114027f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 114127f02ce8SMatthew G. Knepley } 114227f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 114327f02ce8SMatthew G. Knepley } 114427f02ce8SMatthew G. Knepley 11452b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 114620cf1dd8SToby Isaac { 114720cf1dd8SToby Isaac PetscFunctionBegin; 114820cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 114920cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 115020cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 115120cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 115220cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1153ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 115420cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1155afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 115620cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 115720cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 115827f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 115920cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 116020cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 116120cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 116227f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 116320cf1dd8SToby Isaac PetscFunctionReturn(0); 116420cf1dd8SToby Isaac } 116520cf1dd8SToby Isaac 116620cf1dd8SToby Isaac /*MC 116720cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 116820cf1dd8SToby Isaac 116920cf1dd8SToby Isaac Level: intermediate 117020cf1dd8SToby Isaac 117120cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 117220cf1dd8SToby Isaac M*/ 117320cf1dd8SToby Isaac 117420cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 117520cf1dd8SToby Isaac { 117620cf1dd8SToby Isaac PetscFE_Basic *b; 117720cf1dd8SToby Isaac PetscErrorCode ierr; 117820cf1dd8SToby Isaac 117920cf1dd8SToby Isaac PetscFunctionBegin; 118020cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 118120cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 118220cf1dd8SToby Isaac fem->data = b; 118320cf1dd8SToby Isaac 118420cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 118520cf1dd8SToby Isaac PetscFunctionReturn(0); 118620cf1dd8SToby Isaac } 1187