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; 386*6587ee25SMatthew 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 } 415*6587ee25SMatthew 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; 418*6587ee25SMatthew 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 422*6587ee25SMatthew 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 429*6587ee25SMatthew G. Knepley ierr = PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q*cgeom->dim], &fegeom);CHKERRQ(ierr); 4304bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 431*6587ee25SMatthew 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 } 437*6587ee25SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d wt %g\n", q, quadWeights[q]);CHKERRQ(ierr);} 438ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 439ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 4406528b96dSMatthew 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]); 441ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w; 4426528b96dSMatthew 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]); 443ef0bb6c7SMatthew 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; 44420cf1dd8SToby Isaac } 445*6587ee25SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 44620cf1dd8SToby Isaac cOffset += totDim; 44720cf1dd8SToby Isaac cOffsetAux += totDimAux; 44820cf1dd8SToby Isaac } 44920cf1dd8SToby Isaac PetscFunctionReturn(0); 45020cf1dd8SToby Isaac } 45120cf1dd8SToby Isaac 45206ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 4534bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 45420cf1dd8SToby Isaac { 45520cf1dd8SToby Isaac const PetscInt debug = 0; 45606d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 4574bee2e38SMatthew G. Knepley PetscFE fe; 45806d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 45906d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 46020cf1dd8SToby Isaac PetscQuadrature quad; 461ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4624bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 46320cf1dd8SToby Isaac const PetscScalar *constants; 46420cf1dd8SToby Isaac PetscReal *x; 465ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 466ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 467*6587ee25SMatthew G. Knepley PetscBool auxOnBd = PETSC_FALSE; 46820cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 469*6587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 47020cf1dd8SToby Isaac PetscErrorCode ierr; 47120cf1dd8SToby Isaac 47220cf1dd8SToby Isaac PetscFunctionBegin; 4734bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 4744bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 4754bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 4764bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 4774bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 4784bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 4794bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 4804bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 48106ad1575SMatthew G. Knepley ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 48206d8a0d3SMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 4834bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 4844bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 4854bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 486ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr); 4874bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 4884bee2e38SMatthew G. Knepley if (dsAux) { 4894bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 4904bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 4914bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 4924bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 4934bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 4944bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 4957be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 496ef0bb6c7SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 497ef0bb6c7SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 49801907d53SMatthew 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); 49920cf1dd8SToby Isaac } 500ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 501*6587ee25SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 502afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 50320cf1dd8SToby Isaac dE = fgeom->dimEmbed; 504*6587ee25SMatthew G. Knepley /* TODO FIX THIS */ 505*6587ee25SMatthew G. Knepley fgeom->dim = dim-1; 506*6587ee25SMatthew G. Knepley if (fgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", fgeom->dim, qdim); 50720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5089f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 50920cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5109f209ee4SMatthew G. Knepley 511*6587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 512580bdb30SBarry Smith ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr); 51327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcI*dE);CHKERRQ(ierr); 51420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5154bee2e38SMatthew G. Knepley PetscReal w; 5164bee2e38SMatthew G. Knepley PetscInt c, d; 5174bee2e38SMatthew G. Knepley 518*6587ee25SMatthew G. Knepley ierr = PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom);CHKERRQ(ierr); 519*6587ee25SMatthew G. Knepley ierr = PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom);CHKERRQ(ierr); 5204bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 52162bd480fSMatthew G. Knepley if (debug > 1) { 522*6587ee25SMatthew G. Knepley if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 5234bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 5247be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5254bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 52662bd480fSMatthew G. Knepley ierr = DMPrintCellVector(e, "n", dim, fegeom.n);CHKERRQ(ierr); 52720cf1dd8SToby Isaac #endif 52820cf1dd8SToby Isaac } 52962bd480fSMatthew G. Knepley } 530ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 531ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 53206d8a0d3SMatthew 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]); 5334bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 53406d8a0d3SMatthew 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]); 5354bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 53662bd480fSMatthew G. Knepley if (debug) { 53762bd480fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " elem %D quad point %d\n", e, q);CHKERRQ(ierr); 53862bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 53962bd480fSMatthew G. Knepley if (n0) {ierr = PetscPrintf(PETSC_COMM_SELF, " f0[%D] %g\n", c, f0[q*NcI+c]);CHKERRQ(ierr);} 54062bd480fSMatthew G. Knepley if (n1) { 54162bd480fSMatthew 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);} 54262bd480fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 54362bd480fSMatthew G. Knepley } 54462bd480fSMatthew G. Knepley } 54562bd480fSMatthew G. Knepley } 54620cf1dd8SToby Isaac } 547*6587ee25SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 54820cf1dd8SToby Isaac cOffset += totDim; 54920cf1dd8SToby Isaac cOffsetAux += totDimAux; 55020cf1dd8SToby Isaac } 55120cf1dd8SToby Isaac PetscFunctionReturn(0); 55220cf1dd8SToby Isaac } 55320cf1dd8SToby Isaac 55427f02ce8SMatthew G. Knepley /* 55527f02ce8SMatthew 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 55627f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 55727f02ce8SMatthew G. Knepley 55827f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 55927f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 56027f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 56127f02ce8SMatthew 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() 56227f02ce8SMatthew G. Knepley */ 56306ad1575SMatthew G. Knepley static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 56427f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 56527f02ce8SMatthew G. Knepley { 56627f02ce8SMatthew G. Knepley const PetscInt debug = 0; 5676528b96dSMatthew G. Knepley const PetscInt field = key.field; 56827f02ce8SMatthew G. Knepley PetscFE fe; 5696528b96dSMatthew G. Knepley PetscWeakForm wf; 5706528b96dSMatthew G. Knepley PetscInt n0, n1, i; 5716528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 57227f02ce8SMatthew G. Knepley PetscQuadrature quad; 573665f567fSMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 57427f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 57527f02ce8SMatthew G. Knepley const PetscScalar *constants; 57627f02ce8SMatthew G. Knepley PetscReal *x; 577665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 578665f567fSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 579*6587ee25SMatthew G. Knepley PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 58027f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 581*6587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 58227f02ce8SMatthew G. Knepley PetscErrorCode ierr; 58327f02ce8SMatthew G. Knepley 58427f02ce8SMatthew G. Knepley PetscFunctionBegin; 58527f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 58627f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 58727f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 58827f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 58927f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 59027f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 59127f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 59227f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 59327f02ce8SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 5946528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 59506ad1575SMatthew G. Knepley ierr = PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func);CHKERRQ(ierr); 5966528b96dSMatthew G. Knepley if (!n0 && !n1) PetscFunctionReturn(0); 59727f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 59827f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 59927f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 60027f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 601665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr); 60227f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 60327f02ce8SMatthew G. Knepley if (dsAux) { 60427f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 60527f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 60627f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 60727f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 60827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 60927f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 61001907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 611665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 612665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 61301907d53SMatthew 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); 61427f02ce8SMatthew G. Knepley } 61527f02ce8SMatthew G. Knepley isCohesiveField = field == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 616665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 61727f02ce8SMatthew G. Knepley NcS = isCohesiveField ? NcI : 2*NcI; 618*6587ee25SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 61927f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 62027f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 621*6587ee25SMatthew G. Knepley if (fgeom->dim != qdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %D != %D quadrature dim", fgeom->dim, qdim); 62227f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 62327f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 62427f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 62527f02ce8SMatthew G. Knepley 626*6587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 62727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr); 62827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr); 62927f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 63027f02ce8SMatthew G. Knepley PetscReal w; 63127f02ce8SMatthew G. Knepley PetscInt c, d; 63227f02ce8SMatthew G. Knepley 633*6587ee25SMatthew G. Knepley ierr = PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q*fgeom->dim], &fegeom);CHKERRQ(ierr); 63427f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 635*6587ee25SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 63627f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 63727f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 63827f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr); 63927f02ce8SMatthew G. Knepley #endif 64027f02ce8SMatthew G. Knepley } 64127f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 64227f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 643665f567fSMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 644665f567fSMatthew 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);} 6456528b96dSMatthew 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]); 64627f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w; 6476528b96dSMatthew 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]); 64827f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w; 64927f02ce8SMatthew G. Knepley } 650*6587ee25SMatthew G. Knepley if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 651*6587ee25SMatthew G. Knepley else {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 65227f02ce8SMatthew G. Knepley cOffset += totDim; 65327f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 65427f02ce8SMatthew G. Knepley } 65527f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 65627f02ce8SMatthew G. Knepley } 65727f02ce8SMatthew G. Knepley 65806ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, 6594bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 66020cf1dd8SToby Isaac { 66120cf1dd8SToby Isaac const PetscInt debug = 0; 6624bee2e38SMatthew G. Knepley PetscFE feI, feJ; 6636528b96dSMatthew G. Knepley PetscWeakForm wf; 6646528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 6656528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 66620cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 66720cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 66820cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 66920cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 67020cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 67120cf1dd8SToby Isaac PetscQuadrature quad; 672ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 6734bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 67420cf1dd8SToby Isaac const PetscScalar *constants; 67520cf1dd8SToby Isaac PetscReal *x; 676ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 677ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 6786528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 67920cf1dd8SToby Isaac PetscInt dE, Np; 68020cf1dd8SToby Isaac PetscBool isAffine; 68120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 68220cf1dd8SToby Isaac PetscInt qNc, Nq, q; 68320cf1dd8SToby Isaac PetscErrorCode ierr; 68420cf1dd8SToby Isaac 68520cf1dd8SToby Isaac PetscFunctionBegin; 6866528b96dSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 6876528b96dSMatthew G. Knepley fieldI = key.field / Nf; 6886528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 6894bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 6904bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 6914bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 6924bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 6934bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 6944bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 6954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 6966528b96dSMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 69720cf1dd8SToby Isaac switch(jtype) { 69806ad1575SMatthew 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; 69906ad1575SMatthew 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; 70006ad1575SMatthew 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; 70120cf1dd8SToby Isaac } 7026528b96dSMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 7034bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 7044bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 7054bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 706ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 7074bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 7084bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 7094bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 7104bee2e38SMatthew G. Knepley if (dsAux) { 7114bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 7124bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 7134bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 7144bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 7154bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 716ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 71701907d53SMatthew 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); 71820cf1dd8SToby Isaac } 71927f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 72027f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7214bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7224bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7234bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 72427f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 72527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 72627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 72727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 72827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 72927f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 73027f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 7314bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7324bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7334bee2e38SMatthew G. Knepley 73427f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 73527f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7364bee2e38SMatthew G. Knepley if (isAffine) { 7374bee2e38SMatthew G. Knepley fegeom.v = x; 7384bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 7397132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 7407132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 7417132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 7424bee2e38SMatthew G. Knepley } 74320cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 74420cf1dd8SToby Isaac PetscReal w; 7454bee2e38SMatthew G. Knepley PetscInt c; 74620cf1dd8SToby Isaac 74720cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 74820cf1dd8SToby Isaac if (isAffine) { 7497132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 75020cf1dd8SToby Isaac } else { 7514bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 7524bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 7534bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 7544bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 75520cf1dd8SToby Isaac } 7564bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 757ef0bb6c7SMatthew 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);} 758ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 759ea672e62SMatthew G. Knepley if (n0) { 760580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 7616528b96dSMatthew 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); 76220cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 76320cf1dd8SToby Isaac } 764ea672e62SMatthew G. Knepley if (n1) { 76527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 7666528b96dSMatthew 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); 7674bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 76820cf1dd8SToby Isaac } 769ea672e62SMatthew G. Knepley if (n2) { 77027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 7716528b96dSMatthew 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); 7724bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 77320cf1dd8SToby Isaac } 774ea672e62SMatthew G. Knepley if (n3) { 77527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 7766528b96dSMatthew 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); 7774bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 77820cf1dd8SToby Isaac } 77920cf1dd8SToby Isaac 780ef0bb6c7SMatthew 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); 78120cf1dd8SToby Isaac } 78220cf1dd8SToby Isaac if (debug > 1) { 78320cf1dd8SToby Isaac PetscInt fc, f, gc, g; 78420cf1dd8SToby Isaac 78520cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 786ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 787ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 788ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 789ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 790ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 791ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 79220cf1dd8SToby 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); 79320cf1dd8SToby Isaac } 79420cf1dd8SToby Isaac } 79520cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 79620cf1dd8SToby Isaac } 79720cf1dd8SToby Isaac } 79820cf1dd8SToby Isaac } 79920cf1dd8SToby Isaac cOffset += totDim; 80020cf1dd8SToby Isaac cOffsetAux += totDimAux; 80120cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 80220cf1dd8SToby Isaac } 80320cf1dd8SToby Isaac PetscFunctionReturn(0); 80420cf1dd8SToby Isaac } 80520cf1dd8SToby Isaac 80606ad1575SMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 8074bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 80820cf1dd8SToby Isaac { 80920cf1dd8SToby Isaac const PetscInt debug = 0; 8104bee2e38SMatthew G. Knepley PetscFE feI, feJ; 81145480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 81245480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 81320cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 81420cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 81520cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 81620cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 81720cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 81820cf1dd8SToby Isaac PetscQuadrature quad; 819ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8204bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 82120cf1dd8SToby Isaac const PetscScalar *constants; 82220cf1dd8SToby Isaac PetscReal *x; 823ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 824ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 82545480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 82620cf1dd8SToby Isaac PetscBool isAffine; 82720cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 82820cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 82920cf1dd8SToby Isaac PetscErrorCode ierr; 83020cf1dd8SToby Isaac 83120cf1dd8SToby Isaac PetscFunctionBegin; 83245480ffeSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 83345480ffeSMatthew G. Knepley fieldI = key.field / Nf; 83445480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8354bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 8364bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 8374bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 8384bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 8394bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 8404bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 8414bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 8424bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 8434bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 84406ad1575SMatthew 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); 84545480ffeSMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 8464bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 8474bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 8484bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 849ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr); 8504bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 8514bee2e38SMatthew G. Knepley if (dsAux) { 8524bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 8534bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 8544bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 8554bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 8564bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 857ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr); 85820cf1dd8SToby Isaac } 859ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 86020cf1dd8SToby Isaac Np = fgeom->numPoints; 86120cf1dd8SToby Isaac dE = fgeom->dimEmbed; 86220cf1dd8SToby Isaac isAffine = fgeom->isAffine; 86327f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 86431f073d7SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 86527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 86627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 86727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 86827f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 86927f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 87020cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 8719f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 87220cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 873ea78f98cSLisandro Dalcin fegeom.n = NULL; 874ea78f98cSLisandro Dalcin fegeom.v = NULL; 875ea78f98cSLisandro Dalcin fegeom.J = NULL; 876ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 87727f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 87827f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 87927f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 88027f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 8814bee2e38SMatthew G. Knepley if (isAffine) { 8824bee2e38SMatthew G. Knepley fegeom.v = x; 8834bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 8847132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 8857132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 8867132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 8877132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 8889f209ee4SMatthew G. Knepley 8897132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 8907132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 8917132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 8924bee2e38SMatthew G. Knepley } 89320cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 89420cf1dd8SToby Isaac PetscReal w; 8954bee2e38SMatthew G. Knepley PetscInt c; 89620cf1dd8SToby Isaac 89720cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 89820cf1dd8SToby Isaac if (isAffine) { 8997132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 90020cf1dd8SToby Isaac } else { 9013fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 9029f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 9039f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 9044bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 9054bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 9069f209ee4SMatthew G. Knepley 9079f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 9089f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 9099f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 91020cf1dd8SToby Isaac } 9114bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 912ef0bb6c7SMatthew 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);} 913ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 914ea672e62SMatthew G. Knepley if (n0) { 915580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 91645480ffeSMatthew 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); 91720cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 91820cf1dd8SToby Isaac } 919ea672e62SMatthew G. Knepley if (n1) { 92027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 92145480ffeSMatthew 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); 9224bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 92320cf1dd8SToby Isaac } 924ea672e62SMatthew G. Knepley if (n2) { 92527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 92645480ffeSMatthew 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); 9274bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 92820cf1dd8SToby Isaac } 929ea672e62SMatthew G. Knepley if (n3) { 93027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 93145480ffeSMatthew 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); 9324bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 93320cf1dd8SToby Isaac } 93420cf1dd8SToby Isaac 935ef0bb6c7SMatthew 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); 93620cf1dd8SToby Isaac } 93720cf1dd8SToby Isaac if (debug > 1) { 93820cf1dd8SToby Isaac PetscInt fc, f, gc, g; 93920cf1dd8SToby Isaac 94020cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 941ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 942ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 943ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 944ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 945ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 946ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 94720cf1dd8SToby 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); 94820cf1dd8SToby Isaac } 94920cf1dd8SToby Isaac } 95020cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 95120cf1dd8SToby Isaac } 95220cf1dd8SToby Isaac } 95320cf1dd8SToby Isaac } 95420cf1dd8SToby Isaac cOffset += totDim; 95520cf1dd8SToby Isaac cOffsetAux += totDimAux; 95620cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 95720cf1dd8SToby Isaac } 95820cf1dd8SToby Isaac PetscFunctionReturn(0); 95920cf1dd8SToby Isaac } 96020cf1dd8SToby Isaac 96106ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, 96227f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 96327f02ce8SMatthew G. Knepley { 96427f02ce8SMatthew G. Knepley const PetscInt debug = 0; 96527f02ce8SMatthew G. Knepley PetscFE feI, feJ; 966148442b3SMatthew G. Knepley PetscWeakForm wf; 967148442b3SMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 968148442b3SMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 96927f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 97027f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 97127f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 97227f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 97327f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 974665f567fSMatthew G. Knepley PetscQuadrature quad; 975665f567fSMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 97627f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 97727f02ce8SMatthew G. Knepley const PetscScalar *constants; 97827f02ce8SMatthew G. Knepley PetscReal *x; 979665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 980665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 98145480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 982665f567fSMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE; 98327f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 98427f02ce8SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 98527f02ce8SMatthew G. Knepley PetscErrorCode ierr; 98627f02ce8SMatthew G. Knepley 98727f02ce8SMatthew G. Knepley PetscFunctionBegin; 98845480ffeSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 98945480ffeSMatthew G. Knepley fieldI = key.field / Nf; 99045480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 99127f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 99227f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 99327f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 99427f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 99527f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 99627f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 99727f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 99827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 999148442b3SMatthew G. Knepley ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr); 100027f02ce8SMatthew G. Knepley switch(jtype) { 100106ad1575SMatthew 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; 100206ad1575SMatthew 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; 1003665f567fSMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 100427f02ce8SMatthew G. Knepley } 1005148442b3SMatthew G. Knepley if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(0); 100627f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 100727f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 100827f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 1009665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1010665f567fSMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 1011665f567fSMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 101227f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 101327f02ce8SMatthew G. Knepley if (dsAux) { 101427f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 101527f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 101627f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 101727f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 101827f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 101927f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 102001907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 1021665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);} 1022665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);} 102301907d53SMatthew 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); 102427f02ce8SMatthew G. Knepley } 102527f02ce8SMatthew G. Knepley isCohesiveFieldI = fieldI == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 102627f02ce8SMatthew G. Knepley isCohesiveFieldJ = fieldJ == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 1027665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1028665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 102927f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2*NcI; 103027f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2*NcJ; 103127f02ce8SMatthew G. Knepley Np = fgeom->numPoints; 103227f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 103327f02ce8SMatthew G. Knepley isAffine = fgeom->isAffine; 103427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 103527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 103627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 103727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1038665f567fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1039665f567fSMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 104027f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 104127f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 104227f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 104327f02ce8SMatthew G. Knepley 104427f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 104527f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 104627f02ce8SMatthew G. Knepley if (isAffine) { 104727f02ce8SMatthew G. Knepley fegeom.v = x; 104827f02ce8SMatthew G. Knepley fegeom.xi = fgeom->xi; 104927f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 105027f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 105127f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 105227f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 105327f02ce8SMatthew G. Knepley } 105427f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 105527f02ce8SMatthew G. Knepley PetscReal w; 105627f02ce8SMatthew G. Knepley PetscInt c; 105727f02ce8SMatthew G. Knepley 105827f02ce8SMatthew G. Knepley if (isAffine) { 105927f02ce8SMatthew G. Knepley /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */ 106027f02ce8SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 106127f02ce8SMatthew G. Knepley } else { 106227f02ce8SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 106327f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 106427f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 106527f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 106627f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 106727f02ce8SMatthew G. Knepley } 106827f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 106927f02ce8SMatthew G. Knepley if (debug > 1 && q < Np) { 107027f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 107127f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 107227f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 107327f02ce8SMatthew G. Knepley #endif 107427f02ce8SMatthew G. Knepley } 107527f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 1076665f567fSMatthew 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);} 1077665f567fSMatthew 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);} 1078ea672e62SMatthew G. Knepley if (n0) { 107927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 1080148442b3SMatthew 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); 108127f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT; ++c) g0[c] *= w; 108227f02ce8SMatthew G. Knepley } 1083ea672e62SMatthew G. Knepley if (n1) { 108427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 1085148442b3SMatthew 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); 108627f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w; 108727f02ce8SMatthew G. Knepley } 1088ea672e62SMatthew G. Knepley if (n2) { 108927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 1090148442b3SMatthew 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); 109127f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w; 109227f02ce8SMatthew G. Knepley } 1093ea672e62SMatthew G. Knepley if (n3) { 109427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1095148442b3SMatthew 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); 109627f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w; 109727f02ce8SMatthew G. Knepley } 109827f02ce8SMatthew G. Knepley 109927f02ce8SMatthew G. Knepley if (isCohesiveFieldI && isCohesiveFieldJ) { 1100665f567fSMatthew 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); 110127f02ce8SMatthew G. Knepley } else { 1102665f567fSMatthew 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); 110327f02ce8SMatthew G. Knepley } 110427f02ce8SMatthew G. Knepley } 110527f02ce8SMatthew G. Knepley if (debug > 1) { 110627f02ce8SMatthew G. Knepley PetscInt fc, f, gc, g; 110727f02ce8SMatthew G. Knepley 110827f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 110927f02ce8SMatthew G. Knepley for (fc = 0; fc < NcI; ++fc) { 1110665f567fSMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 111127f02ce8SMatthew G. Knepley const PetscInt i = offsetI + f*NcI+fc; 111227f02ce8SMatthew G. Knepley for (gc = 0; gc < NcJ; ++gc) { 1113665f567fSMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 111427f02ce8SMatthew G. Knepley const PetscInt j = offsetJ + g*NcJ+gc; 111527f02ce8SMatthew 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); 111627f02ce8SMatthew G. Knepley } 111727f02ce8SMatthew G. Knepley } 111827f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 111927f02ce8SMatthew G. Knepley } 112027f02ce8SMatthew G. Knepley } 112127f02ce8SMatthew G. Knepley } 112227f02ce8SMatthew G. Knepley cOffset += totDim; 112327f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 112427f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 112527f02ce8SMatthew G. Knepley } 112627f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 112727f02ce8SMatthew G. Knepley } 112827f02ce8SMatthew G. Knepley 11292b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 113020cf1dd8SToby Isaac { 113120cf1dd8SToby Isaac PetscFunctionBegin; 113220cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 113320cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 113420cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 113520cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 113620cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1137ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 113820cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1139afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 114020cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 114120cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 114227f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 114320cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 114420cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 114520cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 114627f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 114720cf1dd8SToby Isaac PetscFunctionReturn(0); 114820cf1dd8SToby Isaac } 114920cf1dd8SToby Isaac 115020cf1dd8SToby Isaac /*MC 115120cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 115220cf1dd8SToby Isaac 115320cf1dd8SToby Isaac Level: intermediate 115420cf1dd8SToby Isaac 115520cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 115620cf1dd8SToby Isaac M*/ 115720cf1dd8SToby Isaac 115820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 115920cf1dd8SToby Isaac { 116020cf1dd8SToby Isaac PetscFE_Basic *b; 116120cf1dd8SToby Isaac PetscErrorCode ierr; 116220cf1dd8SToby Isaac 116320cf1dd8SToby Isaac PetscFunctionBegin; 116420cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 116520cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 116620cf1dd8SToby Isaac fem->data = b; 116720cf1dd8SToby Isaac 116820cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 116920cf1dd8SToby Isaac PetscFunctionReturn(0); 117020cf1dd8SToby Isaac } 1171