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 3696528b96dSMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscHashFormKey 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 PetscBool isAffine; 38620cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3874bee2e38SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 38820cf1dd8SToby Isaac PetscErrorCode ierr; 38920cf1dd8SToby Isaac 39020cf1dd8SToby Isaac PetscFunctionBegin; 3914bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 3924bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 3934bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 3944bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 3954bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 3964bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 3974bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 3984bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 3996528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 4006528b96dSMatthew G. Knepley ierr = PetscWeakFormGetResidual(wf, key.label, key.value, key.field, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 4016528b96dSMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 4024bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 4034bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 4044bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 405ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 4064bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 4074bee2e38SMatthew G. Knepley if (dsAux) { 4084bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 4094bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 4104bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 4114bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 4124bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 413ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 41401907d53SMatthew 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); 41520cf1dd8SToby Isaac } 41620cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 4174bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 41820cf1dd8SToby Isaac Np = cgeom->numPoints; 41920cf1dd8SToby Isaac dE = cgeom->dimEmbed; 42020cf1dd8SToby Isaac isAffine = cgeom->isAffine; 42120cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4224bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 42320cf1dd8SToby Isaac 42427f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 42527f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 4264bee2e38SMatthew G. Knepley if (isAffine) { 4274bee2e38SMatthew G. Knepley fegeom.v = x; 4284bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 4297132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 4307132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 4317132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 4324bee2e38SMatthew G. Knepley } 433ef0bb6c7SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr); 43427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dE);CHKERRQ(ierr); 43520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4364bee2e38SMatthew G. Knepley PetscReal w; 4374bee2e38SMatthew G. Knepley PetscInt c, d; 43820cf1dd8SToby Isaac 43920cf1dd8SToby Isaac if (isAffine) { 4407132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 44120cf1dd8SToby Isaac } else { 4424bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 4434bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 4444bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 4454bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 44620cf1dd8SToby Isaac } 4474bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 44820cf1dd8SToby Isaac if (debug > 1 && q < Np) { 4494bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 4507be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4514bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 45220cf1dd8SToby Isaac #endif 45320cf1dd8SToby Isaac } 45420cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 455ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 456ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 4576528b96dSMatthew 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]); 458ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w; 4596528b96dSMatthew 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]); 460ef0bb6c7SMatthew 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; 46120cf1dd8SToby Isaac } 462ef0bb6c7SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 46320cf1dd8SToby Isaac cOffset += totDim; 46420cf1dd8SToby Isaac cOffsetAux += totDimAux; 46520cf1dd8SToby Isaac } 46620cf1dd8SToby Isaac PetscFunctionReturn(0); 46720cf1dd8SToby Isaac } 46820cf1dd8SToby Isaac 46945480ffeSMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 4704bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 47120cf1dd8SToby Isaac { 47220cf1dd8SToby Isaac const PetscInt debug = 0; 47306d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 4744bee2e38SMatthew G. Knepley PetscFE fe; 47506d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 47606d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 47720cf1dd8SToby Isaac PetscQuadrature quad; 478ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4794bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 48020cf1dd8SToby Isaac const PetscScalar *constants; 48120cf1dd8SToby Isaac PetscReal *x; 482ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 483ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 484b07bd890SMatthew G. Knepley PetscBool isAffine, auxOnBd = PETSC_FALSE; 48520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 48620cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 48720cf1dd8SToby Isaac PetscErrorCode ierr; 48820cf1dd8SToby Isaac 48920cf1dd8SToby Isaac PetscFunctionBegin; 4904bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 4914bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 4924bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 4934bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 4944bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 4954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 4964bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 4974bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 49806d8a0d3SMatthew G. Knepley ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 49906d8a0d3SMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 5004bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 5014bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 5024bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 503ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr); 5044bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 5054bee2e38SMatthew G. Knepley if (dsAux) { 5064bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 5074bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 5084bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 5094bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 5104bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 5114bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 5127be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 513ef0bb6c7SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 514ef0bb6c7SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 51501907d53SMatthew 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); 51620cf1dd8SToby Isaac } 517ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 51820cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 519afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 52020cf1dd8SToby Isaac Np = fgeom->numPoints; 52120cf1dd8SToby Isaac dE = fgeom->dimEmbed; 52220cf1dd8SToby Isaac isAffine = fgeom->isAffine; 52320cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5249f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 52520cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 526ea78f98cSLisandro Dalcin fegeom.n = NULL; 527ea78f98cSLisandro Dalcin fegeom.v = NULL; 528ea78f98cSLisandro Dalcin fegeom.J = NULL; 529ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 53027f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 53127f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 53227f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 53327f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 5344bee2e38SMatthew G. Knepley if (isAffine) { 5354bee2e38SMatthew G. Knepley fegeom.v = x; 5364bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 5377132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 5387132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 5397132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 5407132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 5419f209ee4SMatthew G. Knepley 5427132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 5437132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 5447132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 5454bee2e38SMatthew G. Knepley } 546580bdb30SBarry Smith ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr); 54727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcI*dE);CHKERRQ(ierr); 54820cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5494bee2e38SMatthew G. Knepley PetscReal w; 5504bee2e38SMatthew G. Knepley PetscInt c, d; 5514bee2e38SMatthew G. Knepley 55220cf1dd8SToby Isaac if (isAffine) { 5537132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 55420cf1dd8SToby Isaac } else { 5553fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 5569f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 5579f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 5584bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 5594bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 5609f209ee4SMatthew G. Knepley 5619f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 5629f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 5639f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 56420cf1dd8SToby Isaac } 5654bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 56662bd480fSMatthew G. Knepley if (debug > 1) { 56762bd480fSMatthew G. Knepley if ((isAffine && q == 0) || (!isAffine)) { 5684bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 5697be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5704bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 57162bd480fSMatthew G. Knepley ierr = DMPrintCellVector(e, "n", dim, fegeom.n);CHKERRQ(ierr); 57220cf1dd8SToby Isaac #endif 57320cf1dd8SToby Isaac } 57462bd480fSMatthew G. Knepley } 575ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 576ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 57706d8a0d3SMatthew 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]); 5784bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 57906d8a0d3SMatthew 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]); 5804bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 58162bd480fSMatthew G. Knepley if (debug) { 58262bd480fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " elem %D quad point %d\n", e, q);CHKERRQ(ierr); 58362bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 58462bd480fSMatthew G. Knepley if (n0) {ierr = PetscPrintf(PETSC_COMM_SELF, " f0[%D] %g\n", c, f0[q*NcI+c]);CHKERRQ(ierr);} 58562bd480fSMatthew G. Knepley if (n1) { 58662bd480fSMatthew 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);} 58762bd480fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 58862bd480fSMatthew G. Knepley } 58962bd480fSMatthew G. Knepley } 59062bd480fSMatthew G. Knepley } 59120cf1dd8SToby Isaac } 592ef0bb6c7SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, &cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 59320cf1dd8SToby Isaac cOffset += totDim; 59420cf1dd8SToby Isaac cOffsetAux += totDimAux; 59520cf1dd8SToby Isaac } 59620cf1dd8SToby Isaac PetscFunctionReturn(0); 59720cf1dd8SToby Isaac } 59820cf1dd8SToby Isaac 59927f02ce8SMatthew G. Knepley /* 60027f02ce8SMatthew 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 60127f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 60227f02ce8SMatthew G. Knepley 60327f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 60427f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 60527f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 60627f02ce8SMatthew 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() 60727f02ce8SMatthew G. Knepley */ 6086528b96dSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 60927f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 61027f02ce8SMatthew G. Knepley { 61127f02ce8SMatthew G. Knepley const PetscInt debug = 0; 6126528b96dSMatthew G. Knepley const PetscInt field = key.field; 61327f02ce8SMatthew G. Knepley PetscFE fe; 6146528b96dSMatthew G. Knepley PetscWeakForm wf; 6156528b96dSMatthew G. Knepley PetscInt n0, n1, i; 6166528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 61727f02ce8SMatthew G. Knepley PetscQuadrature quad; 618665f567fSMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 61927f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 62027f02ce8SMatthew G. Knepley const PetscScalar *constants; 62127f02ce8SMatthew G. Knepley PetscReal *x; 622665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 623665f567fSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 624665f567fSMatthew G. Knepley PetscBool isCohesiveField, isAffine, auxOnBd = PETSC_FALSE; 62527f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 62627f02ce8SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 62727f02ce8SMatthew G. Knepley PetscErrorCode ierr; 62827f02ce8SMatthew G. Knepley 62927f02ce8SMatthew G. Knepley PetscFunctionBegin; 63027f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 63127f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 63227f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 63327f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 63427f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 63527f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 63627f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 63727f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 63827f02ce8SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 6396528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 6406528b96dSMatthew G. Knepley ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 6416528b96dSMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 64227f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 64327f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 64427f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 64527f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 646665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr); 64727f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 64827f02ce8SMatthew G. Knepley if (dsAux) { 64927f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 65027f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 65127f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 65227f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 65327f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 65427f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 65501907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 656665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 657665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 65801907d53SMatthew 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); 65927f02ce8SMatthew G. Knepley } 66027f02ce8SMatthew G. Knepley isCohesiveField = field == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 661665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 66227f02ce8SMatthew G. Knepley NcS = isCohesiveField ? NcI : 2*NcI; 66327f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 66427f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 66527f02ce8SMatthew G. Knepley Np = fgeom->numPoints; 66627f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 66727f02ce8SMatthew G. Knepley isAffine = fgeom->isAffine; 66827f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 66927f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 67027f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 67127f02ce8SMatthew G. Knepley 67227f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 67327f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 67427f02ce8SMatthew G. Knepley if (isAffine) { 67527f02ce8SMatthew G. Knepley fegeom.v = x; 67627f02ce8SMatthew G. Knepley fegeom.xi = fgeom->xi; 67727f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 67827f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 67927f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 68027f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 68127f02ce8SMatthew G. Knepley } 68227f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr); 68327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr); 68427f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 68527f02ce8SMatthew G. Knepley PetscReal w; 68627f02ce8SMatthew G. Knepley PetscInt c, d; 68727f02ce8SMatthew G. Knepley 68827f02ce8SMatthew G. Knepley if (isAffine) { 68927f02ce8SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 69027f02ce8SMatthew G. Knepley } else { 69127f02ce8SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 69227f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 69327f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 69427f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 69527f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 69627f02ce8SMatthew G. Knepley } 69727f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 69827f02ce8SMatthew G. Knepley if (debug > 1 && q < Np) { 69927f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 70027f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 70127f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr); 70227f02ce8SMatthew G. Knepley #endif 70327f02ce8SMatthew G. Knepley } 70427f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 70527f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 706665f567fSMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 707665f567fSMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 7086528b96dSMatthew 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]); 70927f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w; 7106528b96dSMatthew 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]); 71127f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w; 71227f02ce8SMatthew G. Knepley } 713665f567fSMatthew G. Knepley if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 714665f567fSMatthew G. Knepley else {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 71527f02ce8SMatthew G. Knepley cOffset += totDim; 71627f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 71727f02ce8SMatthew G. Knepley } 71827f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 71927f02ce8SMatthew G. Knepley } 72027f02ce8SMatthew G. Knepley 7216528b96dSMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *cgeom, 7224bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 72320cf1dd8SToby Isaac { 72420cf1dd8SToby Isaac const PetscInt debug = 0; 7254bee2e38SMatthew G. Knepley PetscFE feI, feJ; 7266528b96dSMatthew G. Knepley PetscWeakForm wf; 7276528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 7286528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 72920cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 73020cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 73120cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 73220cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 73320cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 73420cf1dd8SToby Isaac PetscQuadrature quad; 735ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 7364bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 73720cf1dd8SToby Isaac const PetscScalar *constants; 73820cf1dd8SToby Isaac PetscReal *x; 739ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 740ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 7416528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 74220cf1dd8SToby Isaac PetscInt dE, Np; 74320cf1dd8SToby Isaac PetscBool isAffine; 74420cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 74520cf1dd8SToby Isaac PetscInt qNc, Nq, q; 74620cf1dd8SToby Isaac PetscErrorCode ierr; 74720cf1dd8SToby Isaac 74820cf1dd8SToby Isaac PetscFunctionBegin; 7496528b96dSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 7506528b96dSMatthew G. Knepley fieldI = key.field / Nf; 7516528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 7524bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 7534bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 7544bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 7554bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 7564bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 7574bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 7584bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 7596528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 76020cf1dd8SToby Isaac switch(jtype) { 7616528b96dSMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: ierr = PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 7626528b96dSMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 7636528b96dSMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr);break; 76420cf1dd8SToby Isaac } 7656528b96dSMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 7664bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 7674bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 7684bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 769ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 7704bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 7714bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 7724bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 7734bee2e38SMatthew G. Knepley if (dsAux) { 7744bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 7754bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 7764bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 7774bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 7784bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 779ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 78001907d53SMatthew 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); 78120cf1dd8SToby Isaac } 78227f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 78327f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7844bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7854bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7864bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 78727f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 78827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 78927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 79027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 79127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 79227f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 79327f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 7944bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7954bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7964bee2e38SMatthew G. Knepley 79727f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 79827f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7994bee2e38SMatthew G. Knepley if (isAffine) { 8004bee2e38SMatthew G. Knepley fegeom.v = x; 8014bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 8027132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 8037132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 8047132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 8054bee2e38SMatthew G. Knepley } 80620cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 80720cf1dd8SToby Isaac PetscReal w; 8084bee2e38SMatthew G. Knepley PetscInt c; 80920cf1dd8SToby Isaac 81020cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 81120cf1dd8SToby Isaac if (isAffine) { 8127132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 81320cf1dd8SToby Isaac } else { 8144bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 8154bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 8164bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 8174bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 81820cf1dd8SToby Isaac } 8194bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 820ef0bb6c7SMatthew 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);} 821ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 822*ea672e62SMatthew G. Knepley if (n0) { 823580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 8246528b96dSMatthew 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); 82520cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 82620cf1dd8SToby Isaac } 827*ea672e62SMatthew G. Knepley if (n1) { 82827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 8296528b96dSMatthew 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); 8304bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 83120cf1dd8SToby Isaac } 832*ea672e62SMatthew G. Knepley if (n2) { 83327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 8346528b96dSMatthew 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); 8354bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 83620cf1dd8SToby Isaac } 837*ea672e62SMatthew G. Knepley if (n3) { 83827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 8396528b96dSMatthew 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); 8404bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 84120cf1dd8SToby Isaac } 84220cf1dd8SToby Isaac 843ef0bb6c7SMatthew 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); 84420cf1dd8SToby Isaac } 84520cf1dd8SToby Isaac if (debug > 1) { 84620cf1dd8SToby Isaac PetscInt fc, f, gc, g; 84720cf1dd8SToby Isaac 84820cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 849ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 850ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 851ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 852ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 853ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 854ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 85520cf1dd8SToby 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); 85620cf1dd8SToby Isaac } 85720cf1dd8SToby Isaac } 85820cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 85920cf1dd8SToby Isaac } 86020cf1dd8SToby Isaac } 86120cf1dd8SToby Isaac } 86220cf1dd8SToby Isaac cOffset += totDim; 86320cf1dd8SToby Isaac cOffsetAux += totDimAux; 86420cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 86520cf1dd8SToby Isaac } 86620cf1dd8SToby Isaac PetscFunctionReturn(0); 86720cf1dd8SToby Isaac } 86820cf1dd8SToby Isaac 86945480ffeSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 8704bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 87120cf1dd8SToby Isaac { 87220cf1dd8SToby Isaac const PetscInt debug = 0; 8734bee2e38SMatthew G. Knepley PetscFE feI, feJ; 87445480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 87545480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 87620cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 87720cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 87820cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 87920cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 88020cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 88120cf1dd8SToby Isaac PetscQuadrature quad; 882ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8834bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 88420cf1dd8SToby Isaac const PetscScalar *constants; 88520cf1dd8SToby Isaac PetscReal *x; 886ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 887ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 88845480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 88920cf1dd8SToby Isaac PetscBool isAffine; 89020cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 89120cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 89220cf1dd8SToby Isaac PetscErrorCode ierr; 89320cf1dd8SToby Isaac 89420cf1dd8SToby Isaac PetscFunctionBegin; 89545480ffeSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 89645480ffeSMatthew G. Knepley fieldI = key.field / Nf; 89745480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8984bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 8994bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 9004bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 9014bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 9024bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 9034bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 9044bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 9054bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 9064bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 90745480ffeSMatthew G. Knepley ierr = PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func);CHKERRQ(ierr); 90845480ffeSMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 9094bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 9104bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 9114bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 912ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr); 9134bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 9144bee2e38SMatthew G. Knepley if (dsAux) { 9154bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 9164bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 9174bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 9184bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 9194bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 920ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr); 92120cf1dd8SToby Isaac } 922ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 92320cf1dd8SToby Isaac Np = fgeom->numPoints; 92420cf1dd8SToby Isaac dE = fgeom->dimEmbed; 92520cf1dd8SToby Isaac isAffine = fgeom->isAffine; 92627f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 92731f073d7SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 92827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 92927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 93027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 93127f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 93227f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 93320cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 9349f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 93520cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 936ea78f98cSLisandro Dalcin fegeom.n = NULL; 937ea78f98cSLisandro Dalcin fegeom.v = NULL; 938ea78f98cSLisandro Dalcin fegeom.J = NULL; 939ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 94027f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 94127f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 94227f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 94327f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 9444bee2e38SMatthew G. Knepley if (isAffine) { 9454bee2e38SMatthew G. Knepley fegeom.v = x; 9464bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 9477132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 9487132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 9497132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 9507132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 9519f209ee4SMatthew G. Knepley 9527132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 9537132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 9547132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 9554bee2e38SMatthew G. Knepley } 95620cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 95720cf1dd8SToby Isaac PetscReal w; 9584bee2e38SMatthew G. Knepley PetscInt c; 95920cf1dd8SToby Isaac 96020cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 96120cf1dd8SToby Isaac if (isAffine) { 9627132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 96320cf1dd8SToby Isaac } else { 9643fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 9659f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 9669f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 9674bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 9684bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 9699f209ee4SMatthew G. Knepley 9709f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 9719f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 9729f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 97320cf1dd8SToby Isaac } 9744bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 975ef0bb6c7SMatthew 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);} 976ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 977*ea672e62SMatthew G. Knepley if (n0) { 978580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 97945480ffeSMatthew 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); 98020cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 98120cf1dd8SToby Isaac } 982*ea672e62SMatthew G. Knepley if (n1) { 98327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 98445480ffeSMatthew 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); 9854bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 98620cf1dd8SToby Isaac } 987*ea672e62SMatthew G. Knepley if (n2) { 98827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 98945480ffeSMatthew 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); 9904bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 99120cf1dd8SToby Isaac } 992*ea672e62SMatthew G. Knepley if (n3) { 99327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 99445480ffeSMatthew 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); 9954bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 99620cf1dd8SToby Isaac } 99720cf1dd8SToby Isaac 998ef0bb6c7SMatthew 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); 99920cf1dd8SToby Isaac } 100020cf1dd8SToby Isaac if (debug > 1) { 100120cf1dd8SToby Isaac PetscInt fc, f, gc, g; 100220cf1dd8SToby Isaac 100320cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 1004ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 1005ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 1006ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 1007ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 1008ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 1009ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 101020cf1dd8SToby 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); 101120cf1dd8SToby Isaac } 101220cf1dd8SToby Isaac } 101320cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 101420cf1dd8SToby Isaac } 101520cf1dd8SToby Isaac } 101620cf1dd8SToby Isaac } 101720cf1dd8SToby Isaac cOffset += totDim; 101820cf1dd8SToby Isaac cOffsetAux += totDimAux; 101920cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 102020cf1dd8SToby Isaac } 102120cf1dd8SToby Isaac PetscFunctionReturn(0); 102220cf1dd8SToby Isaac } 102320cf1dd8SToby Isaac 102445480ffeSMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscHashFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 102527f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 102627f02ce8SMatthew G. Knepley { 102727f02ce8SMatthew G. Knepley const PetscInt debug = 0; 102827f02ce8SMatthew G. Knepley PetscFE feI, feJ; 1029665f567fSMatthew G. Knepley PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 103027f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 103127f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 103227f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 103327f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 103427f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 1035665f567fSMatthew G. Knepley PetscQuadrature quad; 1036665f567fSMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 103727f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 103827f02ce8SMatthew G. Knepley const PetscScalar *constants; 103927f02ce8SMatthew G. Knepley PetscReal *x; 1040665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1041665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 104245480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 1043665f567fSMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE; 104427f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 104527f02ce8SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 104627f02ce8SMatthew G. Knepley PetscErrorCode ierr; 104727f02ce8SMatthew G. Knepley 104827f02ce8SMatthew G. Knepley PetscFunctionBegin; 104945480ffeSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 105045480ffeSMatthew G. Knepley fieldI = key.field / Nf; 105145480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 105227f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 105327f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 105427f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 105527f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 105627f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 105727f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 105827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 105927f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 106027f02ce8SMatthew G. Knepley switch(jtype) { 106127f02ce8SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetBdJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 106227f02ce8SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 1063665f567fSMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 106427f02ce8SMatthew G. Knepley } 106527f02ce8SMatthew G. Knepley if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 106627f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 106727f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 106827f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 1069665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1070665f567fSMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 1071665f567fSMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 107227f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 107327f02ce8SMatthew G. Knepley if (dsAux) { 107427f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 107527f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 107627f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 107727f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 107827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 107927f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 108001907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 1081665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);} 1082665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);} 108301907d53SMatthew 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); 108427f02ce8SMatthew G. Knepley } 108527f02ce8SMatthew G. Knepley isCohesiveFieldI = fieldI == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 108627f02ce8SMatthew G. Knepley isCohesiveFieldJ = fieldJ == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 1087665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1088665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 108927f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2*NcI; 109027f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2*NcJ; 109127f02ce8SMatthew G. Knepley Np = fgeom->numPoints; 109227f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 109327f02ce8SMatthew G. Knepley isAffine = fgeom->isAffine; 109427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 109527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 109627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 109727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1098665f567fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1099665f567fSMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 110027f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 110127f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 110227f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 110327f02ce8SMatthew G. Knepley 110427f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 110527f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 110627f02ce8SMatthew G. Knepley if (isAffine) { 110727f02ce8SMatthew G. Knepley fegeom.v = x; 110827f02ce8SMatthew G. Knepley fegeom.xi = fgeom->xi; 110927f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 111027f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 111127f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 111227f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 111327f02ce8SMatthew G. Knepley } 111427f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 111527f02ce8SMatthew G. Knepley PetscReal w; 111627f02ce8SMatthew G. Knepley PetscInt c; 111727f02ce8SMatthew G. Knepley 111827f02ce8SMatthew G. Knepley if (isAffine) { 111927f02ce8SMatthew G. Knepley /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */ 112027f02ce8SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 112127f02ce8SMatthew G. Knepley } else { 112227f02ce8SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 112327f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 112427f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 112527f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 112627f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 112727f02ce8SMatthew G. Knepley } 112827f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 112927f02ce8SMatthew G. Knepley if (debug > 1 && q < Np) { 113027f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 113127f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 113227f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 113327f02ce8SMatthew G. Knepley #endif 113427f02ce8SMatthew G. Knepley } 113527f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 1136665f567fSMatthew 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);} 1137665f567fSMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(dsAux, NfAux, 0, auxOnBd ? q : face*Nq+q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 1138*ea672e62SMatthew G. Knepley if (n0) { 113927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 114027f02ce8SMatthew G. Knepley g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0); 114127f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT; ++c) g0[c] *= w; 114227f02ce8SMatthew G. Knepley } 1143*ea672e62SMatthew G. Knepley if (n1) { 114427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 114527f02ce8SMatthew G. Knepley g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1); 114627f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w; 114727f02ce8SMatthew G. Knepley } 1148*ea672e62SMatthew G. Knepley if (n2) { 114927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 115027f02ce8SMatthew G. Knepley g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2); 115127f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w; 115227f02ce8SMatthew G. Knepley } 1153*ea672e62SMatthew G. Knepley if (n3) { 115427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 115527f02ce8SMatthew G. Knepley g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3); 115627f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w; 115727f02ce8SMatthew G. Knepley } 115827f02ce8SMatthew G. Knepley 115927f02ce8SMatthew G. Knepley if (isCohesiveFieldI && isCohesiveFieldJ) { 1160665f567fSMatthew 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*2, offsetJ*2, elemMat);CHKERRQ(ierr); 116127f02ce8SMatthew G. Knepley } else { 1162665f567fSMatthew G. Knepley ierr = PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI*2, offsetJ*2, elemMat);CHKERRQ(ierr); 116327f02ce8SMatthew G. Knepley } 116427f02ce8SMatthew G. Knepley } 116527f02ce8SMatthew G. Knepley if (debug > 1) { 116627f02ce8SMatthew G. Knepley PetscInt fc, f, gc, g; 116727f02ce8SMatthew G. Knepley 116827f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 116927f02ce8SMatthew G. Knepley for (fc = 0; fc < NcI; ++fc) { 1170665f567fSMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 117127f02ce8SMatthew G. Knepley const PetscInt i = offsetI + f*NcI+fc; 117227f02ce8SMatthew G. Knepley for (gc = 0; gc < NcJ; ++gc) { 1173665f567fSMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 117427f02ce8SMatthew G. Knepley const PetscInt j = offsetJ + g*NcJ+gc; 117527f02ce8SMatthew 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); 117627f02ce8SMatthew G. Knepley } 117727f02ce8SMatthew G. Knepley } 117827f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 117927f02ce8SMatthew G. Knepley } 118027f02ce8SMatthew G. Knepley } 118127f02ce8SMatthew G. Knepley } 118227f02ce8SMatthew G. Knepley cOffset += totDim; 118327f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 118427f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 118527f02ce8SMatthew G. Knepley } 118627f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 118727f02ce8SMatthew G. Knepley } 118827f02ce8SMatthew G. Knepley 11892b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 119020cf1dd8SToby Isaac { 119120cf1dd8SToby Isaac PetscFunctionBegin; 119220cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 119320cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 119420cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 119520cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 119620cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1197ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 119820cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1199afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 120020cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 120120cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 120227f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 120320cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 120420cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 120520cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 120627f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 120720cf1dd8SToby Isaac PetscFunctionReturn(0); 120820cf1dd8SToby Isaac } 120920cf1dd8SToby Isaac 121020cf1dd8SToby Isaac /*MC 121120cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 121220cf1dd8SToby Isaac 121320cf1dd8SToby Isaac Level: intermediate 121420cf1dd8SToby Isaac 121520cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 121620cf1dd8SToby Isaac M*/ 121720cf1dd8SToby Isaac 121820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 121920cf1dd8SToby Isaac { 122020cf1dd8SToby Isaac PetscFE_Basic *b; 122120cf1dd8SToby Isaac PetscErrorCode ierr; 122220cf1dd8SToby Isaac 122320cf1dd8SToby Isaac PetscFunctionBegin; 122420cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 122520cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 122620cf1dd8SToby Isaac fem->data = b; 122720cf1dd8SToby Isaac 122820cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 122920cf1dd8SToby Isaac PetscFunctionReturn(0); 123020cf1dd8SToby Isaac } 1231