120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscblaslapack.h> 320cf1dd8SToby Isaac 42b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem) 520cf1dd8SToby Isaac { 620cf1dd8SToby Isaac PetscFE_Basic *b = (PetscFE_Basic *) fem->data; 720cf1dd8SToby Isaac PetscErrorCode ierr; 820cf1dd8SToby Isaac 920cf1dd8SToby Isaac PetscFunctionBegin; 1020cf1dd8SToby Isaac ierr = PetscFree(b);CHKERRQ(ierr); 1120cf1dd8SToby Isaac PetscFunctionReturn(0); 1220cf1dd8SToby Isaac } 1320cf1dd8SToby Isaac 142b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v) 1520cf1dd8SToby Isaac { 16d9bac1caSLisandro Dalcin PetscInt dim, Nc; 17d9bac1caSLisandro Dalcin PetscSpace basis = NULL; 18d9bac1caSLisandro Dalcin PetscDualSpace dual = NULL; 19d9bac1caSLisandro Dalcin PetscQuadrature quad = NULL; 2020cf1dd8SToby Isaac PetscErrorCode ierr; 2120cf1dd8SToby Isaac 2220cf1dd8SToby Isaac PetscFunctionBegin; 23d9bac1caSLisandro Dalcin ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 24d9bac1caSLisandro Dalcin ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 2520cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr); 2620cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr); 27d9bac1caSLisandro Dalcin ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 28d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 29d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPrintf(v, "Basic Finite Element in %D dimensions with %D components\n",dim,Nc);CHKERRQ(ierr); 30d9bac1caSLisandro Dalcin if (basis) {ierr = PetscSpaceView(basis, v);CHKERRQ(ierr);} 31d9bac1caSLisandro Dalcin if (dual) {ierr = PetscDualSpaceView(dual, v);CHKERRQ(ierr);} 32d9bac1caSLisandro Dalcin if (quad) {ierr = PetscQuadratureView(quad, v);CHKERRQ(ierr);} 33d9bac1caSLisandro Dalcin ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 3420cf1dd8SToby Isaac PetscFunctionReturn(0); 3520cf1dd8SToby Isaac } 3620cf1dd8SToby Isaac 372b99622eSMatthew G. Knepley static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v) 3820cf1dd8SToby Isaac { 3920cf1dd8SToby Isaac PetscBool iascii; 4020cf1dd8SToby Isaac PetscErrorCode ierr; 4120cf1dd8SToby Isaac 4220cf1dd8SToby Isaac PetscFunctionBegin; 43d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 44d9bac1caSLisandro Dalcin if (iascii) {ierr = PetscFEView_Basic_Ascii(fe, v);CHKERRQ(ierr);} 4520cf1dd8SToby Isaac PetscFunctionReturn(0); 4620cf1dd8SToby Isaac } 4720cf1dd8SToby Isaac 4820cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */ 4915310ec5SMatthew G. Knepley PetscErrorCode PetscFESetUp_Basic(PetscFE fem) 5020cf1dd8SToby Isaac { 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 */ 104b9d4cb8dSJed Brown static PetscErrorCode TensorContract_Private(PetscInt m,PetscInt n,PetscInt p,PetscInt k,const PetscReal *A,const PetscReal *B,PetscReal *C) { 105b9d4cb8dSJed Brown PetscErrorCode ierr; 106b9d4cb8dSJed Brown PetscInt i; 107b9d4cb8dSJed Brown 108b9d4cb8dSJed Brown PetscFunctionBegin; 109b9d4cb8dSJed Brown for (i=0; i<m; i++) { 110b9d4cb8dSJed Brown PetscBLASInt n_,p_,k_,lda,ldb,ldc; 111b9d4cb8dSJed Brown PetscReal one = 1, zero = 0; 112b9d4cb8dSJed Brown /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n] 113b9d4cb8dSJed Brown * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k) 114b9d4cb8dSJed Brown */ 115b9d4cb8dSJed Brown ierr = PetscBLASIntCast(n,&n_);CHKERRQ(ierr); 116b9d4cb8dSJed Brown ierr = PetscBLASIntCast(p,&p_);CHKERRQ(ierr); 117b9d4cb8dSJed Brown ierr = PetscBLASIntCast(k,&k_);CHKERRQ(ierr); 118b9d4cb8dSJed Brown lda = p_; 119b9d4cb8dSJed Brown ldb = n_; 120b9d4cb8dSJed Brown ldc = p_; 121b9d4cb8dSJed Brown PetscStackCallBLAS("BLASgemm",BLASREALgemm_("N","T",&p_,&n_,&k_,&one,A+i*k*p,&lda,B,&ldb,&zero,C+i*n*p,&ldc)); 122b9d4cb8dSJed Brown } 1236485a3bbSJed Brown ierr = PetscLogFlops(2.*m*n*p*k);CHKERRQ(ierr); 124b9d4cb8dSJed Brown PetscFunctionReturn(0); 125b9d4cb8dSJed Brown } 126b9d4cb8dSJed Brown 127ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 12820cf1dd8SToby Isaac { 12920cf1dd8SToby Isaac DM dm; 13020cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 13120cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 13220cf1dd8SToby Isaac PetscInt Nc; /* Field components */ 133ef0bb6c7SMatthew G. Knepley PetscReal *B = K >= 0 ? T->T[0] : NULL; 134ef0bb6c7SMatthew G. Knepley PetscReal *D = K >= 1 ? T->T[1] : NULL; 135ef0bb6c7SMatthew G. Knepley PetscReal *H = K >= 2 ? T->T[2] : NULL; 136ef0bb6c7SMatthew G. Knepley PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL; 13720cf1dd8SToby Isaac PetscErrorCode ierr; 13820cf1dd8SToby Isaac 13920cf1dd8SToby Isaac PetscFunctionBegin; 14020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 14120cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 14220cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 14320cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 14420cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 145ef0bb6c7SMatthew G. Knepley if (K >= 0) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 146ef0bb6c7SMatthew G. Knepley if (K >= 1) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 147ef0bb6c7SMatthew G. Knepley if (K >= 2) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 148ef0bb6c7SMatthew G. Knepley ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH);CHKERRQ(ierr); 149b9d4cb8dSJed Brown /* Translate from prime to nodal basis */ 15020cf1dd8SToby Isaac if (B) { 151b9d4cb8dSJed Brown /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */ 152b9d4cb8dSJed Brown ierr = TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B);CHKERRQ(ierr); 15320cf1dd8SToby Isaac } 15420cf1dd8SToby Isaac if (D) { 155b9d4cb8dSJed Brown /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */ 156b9d4cb8dSJed Brown ierr = TensorContract_Private(npoints, pdim, Nc*dim, pdim, tmpD, fem->invV, D);CHKERRQ(ierr); 15720cf1dd8SToby Isaac } 15820cf1dd8SToby Isaac if (H) { 159b9d4cb8dSJed Brown /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */ 160b9d4cb8dSJed Brown ierr = TensorContract_Private(npoints, pdim, Nc*dim*dim, pdim, tmpH, fem->invV, H);CHKERRQ(ierr); 16120cf1dd8SToby Isaac } 162ef0bb6c7SMatthew G. Knepley if (K >= 0) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 163ef0bb6c7SMatthew G. Knepley if (K >= 1) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 164ef0bb6c7SMatthew G. Knepley if (K >= 2) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 16520cf1dd8SToby Isaac PetscFunctionReturn(0); 16620cf1dd8SToby Isaac } 16720cf1dd8SToby Isaac 1682b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1694bee2e38SMatthew G. Knepley const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 17020cf1dd8SToby Isaac { 17120cf1dd8SToby Isaac const PetscInt debug = 0; 1724bee2e38SMatthew G. Knepley PetscFE fe; 17320cf1dd8SToby Isaac PetscPointFunc obj_func; 17420cf1dd8SToby Isaac PetscQuadrature quad; 175ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 1764bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x; 17720cf1dd8SToby Isaac const PetscScalar *constants; 17820cf1dd8SToby Isaac PetscReal *x; 179ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 18020cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 18120cf1dd8SToby Isaac PetscBool isAffine; 18220cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 18320cf1dd8SToby Isaac PetscInt qNc, Nq, q; 18420cf1dd8SToby Isaac PetscErrorCode ierr; 18520cf1dd8SToby Isaac 18620cf1dd8SToby Isaac PetscFunctionBegin; 1874bee2e38SMatthew G. Knepley ierr = PetscDSGetObjective(ds, field, &obj_func);CHKERRQ(ierr); 18820cf1dd8SToby Isaac if (!obj_func) PetscFunctionReturn(0); 1894bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 1904bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 1914bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 1924bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 1934bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 1944bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 1954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 196ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1974bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 1984bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 1994bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 2004bee2e38SMatthew G. Knepley if (dsAux) { 2014bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 2024bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 2034bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 2044bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 205ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 2064bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 20701907d53SMatthew 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); 20820cf1dd8SToby Isaac } 20920cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 2104bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 21120cf1dd8SToby Isaac Np = cgeom->numPoints; 21220cf1dd8SToby Isaac dE = cgeom->dimEmbed; 21320cf1dd8SToby Isaac isAffine = cgeom->isAffine; 21420cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 2154bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 21620cf1dd8SToby Isaac 21727f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 21827f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 21920cf1dd8SToby Isaac if (isAffine) { 2204bee2e38SMatthew G. Knepley fegeom.v = x; 2214bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2227132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 2237132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 2247132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 22520cf1dd8SToby Isaac } 2264bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 2274bee2e38SMatthew G. Knepley PetscScalar integrand; 2284bee2e38SMatthew G. Knepley PetscReal w; 2294bee2e38SMatthew G. Knepley 2304bee2e38SMatthew G. Knepley if (isAffine) { 2317132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 2324bee2e38SMatthew G. Knepley } else { 2334bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 2344bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 2354bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 2364bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 2374bee2e38SMatthew G. Knepley } 2384bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 23920cf1dd8SToby Isaac if (debug > 1 && q < Np) { 2404bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 2417be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 2424bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 24320cf1dd8SToby Isaac #endif 24420cf1dd8SToby Isaac } 24520cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 246ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 247ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 2484bee2e38SMatthew 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); 2494bee2e38SMatthew G. Knepley integrand *= w; 25020cf1dd8SToby Isaac integral[e*Nf+field] += integrand; 25120cf1dd8SToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);} 25220cf1dd8SToby Isaac } 25320cf1dd8SToby Isaac cOffset += totDim; 25420cf1dd8SToby Isaac cOffsetAux += totDimAux; 25520cf1dd8SToby Isaac } 25620cf1dd8SToby Isaac PetscFunctionReturn(0); 25720cf1dd8SToby Isaac } 25820cf1dd8SToby Isaac 2592b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, 260afe6d6adSToby Isaac PetscBdPointFunc obj_func, 2614bee2e38SMatthew G. Knepley PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 262afe6d6adSToby Isaac { 263afe6d6adSToby Isaac const PetscInt debug = 0; 2644bee2e38SMatthew G. Knepley PetscFE fe; 265afe6d6adSToby Isaac PetscQuadrature quad; 266ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 2674bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 268afe6d6adSToby Isaac const PetscScalar *constants; 269afe6d6adSToby Isaac PetscReal *x; 270ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 271afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 272afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 273afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 274afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 275afe6d6adSToby Isaac PetscErrorCode ierr; 276afe6d6adSToby Isaac 277afe6d6adSToby Isaac PetscFunctionBegin; 278afe6d6adSToby Isaac if (!obj_func) PetscFunctionReturn(0); 2794bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 2804bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 2814bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 2824bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 2834bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 2844bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 2854bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 2864bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x);CHKERRQ(ierr); 2874bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 288ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr); 2894bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 2904bee2e38SMatthew G. Knepley if (dsAux) { 2914bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 2924bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 2934bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 2944bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 2954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 2964bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 297afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 298ef0bb6c7SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 299ef0bb6c7SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 30001907d53SMatthew 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); 301afe6d6adSToby Isaac } 302afe6d6adSToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 303afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 304afe6d6adSToby Isaac Np = fgeom->numPoints; 305afe6d6adSToby Isaac dE = fgeom->dimEmbed; 306afe6d6adSToby Isaac isAffine = fgeom->isAffine; 307afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 3089f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 309afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 310*ea78f98cSLisandro Dalcin fegeom.n = NULL; 311*ea78f98cSLisandro Dalcin fegeom.v = NULL; 312*ea78f98cSLisandro Dalcin fegeom.J = NULL; 313*ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 31427f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 31527f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 31627f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 31727f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3184bee2e38SMatthew G. Knepley if (isAffine) { 3194bee2e38SMatthew G. Knepley fegeom.v = x; 3204bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3217132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 3227132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 3237132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 3247132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 3259f209ee4SMatthew G. Knepley 3267132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 3277132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 3287132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 3294bee2e38SMatthew G. Knepley } 330afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 331afe6d6adSToby Isaac PetscScalar integrand; 3324bee2e38SMatthew G. Knepley PetscReal w; 333afe6d6adSToby Isaac 334afe6d6adSToby Isaac if (isAffine) { 3357132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 336afe6d6adSToby Isaac } else { 3373fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 3389f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 3399f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 3404bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 3414bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 3429f209ee4SMatthew G. Knepley 3439f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 3449f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 3459f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 346afe6d6adSToby Isaac } 3474bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 348afe6d6adSToby Isaac if (debug > 1 && q < Np) { 3494bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 350afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3514bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 352afe6d6adSToby Isaac #endif 353afe6d6adSToby Isaac } 354afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 355ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL);CHKERRQ(ierr); 356ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 3574bee2e38SMatthew 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); 3584bee2e38SMatthew G. Knepley integrand *= w; 359afe6d6adSToby Isaac integral[e*Nf+field] += integrand; 360afe6d6adSToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[e*Nf+field]));CHKERRQ(ierr);} 361afe6d6adSToby Isaac } 362afe6d6adSToby Isaac cOffset += totDim; 363afe6d6adSToby Isaac cOffsetAux += totDimAux; 364afe6d6adSToby Isaac } 365afe6d6adSToby Isaac PetscFunctionReturn(0); 366afe6d6adSToby Isaac } 367afe6d6adSToby Isaac 3684bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 3694bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 37020cf1dd8SToby Isaac { 37120cf1dd8SToby Isaac const PetscInt debug = 0; 3724bee2e38SMatthew G. Knepley PetscFE fe; 37320cf1dd8SToby Isaac PetscPointFunc f0_func; 37420cf1dd8SToby Isaac PetscPointFunc f1_func; 37520cf1dd8SToby Isaac PetscQuadrature quad; 376ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 3774bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 37820cf1dd8SToby Isaac const PetscScalar *constants; 37920cf1dd8SToby Isaac PetscReal *x; 380ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 381ef0bb6c7SMatthew G. Knepley PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 38220cf1dd8SToby Isaac PetscBool isAffine; 38320cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3844bee2e38SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 38520cf1dd8SToby Isaac PetscErrorCode ierr; 38620cf1dd8SToby Isaac 38720cf1dd8SToby Isaac PetscFunctionBegin; 3884bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 3894bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 3904bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 3914bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 3924bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 3934bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 3944bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 3954bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 3964bee2e38SMatthew G. Knepley ierr = PetscDSGetResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 3974bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 3984bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 3994bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 4004bee2e38SMatthew G. Knepley if (!f0_func && !f1_func) PetscFunctionReturn(0); 401ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 4024bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 4034bee2e38SMatthew G. Knepley if (dsAux) { 4044bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 4054bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 4064bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 4074bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 4084bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 409ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 41001907d53SMatthew 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); 41120cf1dd8SToby Isaac } 41220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 4134bee2e38SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 41420cf1dd8SToby Isaac Np = cgeom->numPoints; 41520cf1dd8SToby Isaac dE = cgeom->dimEmbed; 41620cf1dd8SToby Isaac isAffine = cgeom->isAffine; 41720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4184bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 41920cf1dd8SToby Isaac 42027f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 42127f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 4224bee2e38SMatthew G. Knepley if (isAffine) { 4234bee2e38SMatthew G. Knepley fegeom.v = x; 4244bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 4257132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 4267132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 4277132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 4284bee2e38SMatthew G. Knepley } 429ef0bb6c7SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr); 43027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dE);CHKERRQ(ierr); 43120cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4324bee2e38SMatthew G. Knepley PetscReal w; 4334bee2e38SMatthew G. Knepley PetscInt c, d; 43420cf1dd8SToby Isaac 43520cf1dd8SToby Isaac if (isAffine) { 4367132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 43720cf1dd8SToby Isaac } else { 4384bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 4394bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 4404bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 4414bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 44220cf1dd8SToby Isaac } 4434bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 44420cf1dd8SToby Isaac if (debug > 1 && q < Np) { 4454bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 4467be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4474bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 44820cf1dd8SToby Isaac #endif 44920cf1dd8SToby Isaac } 45020cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 451ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 452ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 4534bee2e38SMatthew G. Knepley if (f0_func) { 454ef0bb6c7SMatthew G. Knepley f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q*T[field]->Nc]); 455ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w; 4564bee2e38SMatthew G. Knepley } 45720cf1dd8SToby Isaac if (f1_func) { 458ef0bb6c7SMatthew G. Knepley f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q*T[field]->Nc*dim]); 459ef0bb6c7SMatthew 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; 46020cf1dd8SToby Isaac } 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 4694bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, 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; 4734bee2e38SMatthew G. Knepley PetscFE fe; 47420cf1dd8SToby Isaac PetscBdPointFunc f0_func; 47520cf1dd8SToby Isaac PetscBdPointFunc f1_func; 47620cf1dd8SToby Isaac PetscQuadrature quad; 477ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4784bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 47920cf1dd8SToby Isaac const PetscScalar *constants; 48020cf1dd8SToby Isaac PetscReal *x; 481ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 482ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 483b07bd890SMatthew G. Knepley PetscBool isAffine, auxOnBd = PETSC_FALSE; 48420cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 48520cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 48620cf1dd8SToby Isaac PetscErrorCode ierr; 48720cf1dd8SToby Isaac 48820cf1dd8SToby Isaac PetscFunctionBegin; 4894bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 4904bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 4914bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(fe, &quad);CHKERRQ(ierr); 4924bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 4934bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 4944bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 4954bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 4964bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 4974bee2e38SMatthew G. Knepley ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 49820cf1dd8SToby Isaac if (!f0_func && !f1_func) PetscFunctionReturn(0); 4994bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 5004bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 5014bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 502ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &Tf);CHKERRQ(ierr); 5034bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 5044bee2e38SMatthew G. Knepley if (dsAux) { 5054bee2e38SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 5064bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 5074bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 5084bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 5094bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 5104bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 5117be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 512ef0bb6c7SMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 513ef0bb6c7SMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 51401907d53SMatthew 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); 51520cf1dd8SToby Isaac } 516ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 51720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 518afe6d6adSToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 51920cf1dd8SToby Isaac Np = fgeom->numPoints; 52020cf1dd8SToby Isaac dE = fgeom->dimEmbed; 52120cf1dd8SToby Isaac isAffine = fgeom->isAffine; 52220cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5239f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 52420cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 525*ea78f98cSLisandro Dalcin fegeom.n = NULL; 526*ea78f98cSLisandro Dalcin fegeom.v = NULL; 527*ea78f98cSLisandro Dalcin fegeom.J = NULL; 528*ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 52927f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 53027f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 53127f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 53227f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 5334bee2e38SMatthew G. Knepley if (isAffine) { 5344bee2e38SMatthew G. Knepley fegeom.v = x; 5354bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 5367132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 5377132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 5387132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 5397132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 5409f209ee4SMatthew G. Knepley 5417132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 5427132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 5437132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 5444bee2e38SMatthew G. Knepley } 545580bdb30SBarry Smith ierr = PetscArrayzero(f0, Nq*NcI);CHKERRQ(ierr); 54627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcI*dE);CHKERRQ(ierr); 54720cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5484bee2e38SMatthew G. Knepley PetscReal w; 5494bee2e38SMatthew G. Knepley PetscInt c, d; 5504bee2e38SMatthew G. Knepley 55120cf1dd8SToby Isaac if (isAffine) { 5527132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 55320cf1dd8SToby Isaac } else { 5543fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 5559f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 5569f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 5574bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 5584bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 5599f209ee4SMatthew G. Knepley 5609f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 5619f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 5629f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 56320cf1dd8SToby Isaac } 5644bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 56520cf1dd8SToby Isaac if (debug > 1 && q < Np) { 5664bee2e38SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 5677be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5684bee2e38SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 56920cf1dd8SToby Isaac #endif 57020cf1dd8SToby Isaac } 57120cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 572ef0bb6c7SMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 573ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 5744bee2e38SMatthew G. Knepley if (f0_func) { 5754bee2e38SMatthew G. Knepley f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q*NcI]); 5764bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 5774bee2e38SMatthew G. Knepley } 57820cf1dd8SToby Isaac if (f1_func) { 5794bee2e38SMatthew G. Knepley f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q*NcI*dim]); 5804bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 58120cf1dd8SToby Isaac } 58220cf1dd8SToby Isaac } 583ef0bb6c7SMatthew G. Knepley ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, &cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 58420cf1dd8SToby Isaac cOffset += totDim; 58520cf1dd8SToby Isaac cOffsetAux += totDimAux; 58620cf1dd8SToby Isaac } 58720cf1dd8SToby Isaac PetscFunctionReturn(0); 58820cf1dd8SToby Isaac } 58920cf1dd8SToby Isaac 59027f02ce8SMatthew G. Knepley /* 59127f02ce8SMatthew 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 59227f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 59327f02ce8SMatthew G. Knepley 59427f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 59527f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 59627f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 59727f02ce8SMatthew 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() 59827f02ce8SMatthew G. Knepley */ 59927f02ce8SMatthew G. Knepley static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 60027f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 60127f02ce8SMatthew G. Knepley { 60227f02ce8SMatthew G. Knepley const PetscInt debug = 0; 60327f02ce8SMatthew G. Knepley PetscFE fe; 60427f02ce8SMatthew G. Knepley PetscBdPointFunc f0_func; 60527f02ce8SMatthew G. Knepley PetscBdPointFunc f1_func; 60627f02ce8SMatthew G. Knepley PetscQuadrature quad; 607665f567fSMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 60827f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 60927f02ce8SMatthew G. Knepley const PetscScalar *constants; 61027f02ce8SMatthew G. Knepley PetscReal *x; 611665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 612665f567fSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 613665f567fSMatthew G. Knepley PetscBool isCohesiveField, isAffine, auxOnBd = PETSC_FALSE; 61427f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 61527f02ce8SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 61627f02ce8SMatthew G. Knepley PetscErrorCode ierr; 61727f02ce8SMatthew G. Knepley 61827f02ce8SMatthew G. Knepley PetscFunctionBegin; 61927f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 62027f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 62127f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 62227f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 62327f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 62427f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 62527f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 62627f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 62727f02ce8SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 62827f02ce8SMatthew G. Knepley ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 62927f02ce8SMatthew G. Knepley if (!f0_func && !f1_func) PetscFunctionReturn(0); 63027f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 63127f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 63227f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 63327f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 634665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr); 63527f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 63627f02ce8SMatthew G. Knepley if (dsAux) { 63727f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 63827f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 63927f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 64027f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 64127f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 64227f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 64301907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 644665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 645665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 64601907d53SMatthew 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); 64727f02ce8SMatthew G. Knepley } 64827f02ce8SMatthew G. Knepley isCohesiveField = field == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 649665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 65027f02ce8SMatthew G. Knepley NcS = isCohesiveField ? NcI : 2*NcI; 65127f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 65227f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 65327f02ce8SMatthew G. Knepley Np = fgeom->numPoints; 65427f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 65527f02ce8SMatthew G. Knepley isAffine = fgeom->isAffine; 65627f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 65727f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 65827f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 65927f02ce8SMatthew G. Knepley 66027f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 66127f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 66227f02ce8SMatthew G. Knepley if (isAffine) { 66327f02ce8SMatthew G. Knepley fegeom.v = x; 66427f02ce8SMatthew G. Knepley fegeom.xi = fgeom->xi; 66527f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 66627f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 66727f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 66827f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 66927f02ce8SMatthew G. Knepley } 67027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr); 67127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr); 67227f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 67327f02ce8SMatthew G. Knepley PetscReal w; 67427f02ce8SMatthew G. Knepley PetscInt c, d; 67527f02ce8SMatthew G. Knepley 67627f02ce8SMatthew G. Knepley if (isAffine) { 67727f02ce8SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 67827f02ce8SMatthew G. Knepley } else { 67927f02ce8SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 68027f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 68127f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 68227f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 68327f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 68427f02ce8SMatthew G. Knepley } 68527f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 68627f02ce8SMatthew G. Knepley if (debug > 1 && q < Np) { 68727f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 68827f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 68927f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr); 69027f02ce8SMatthew G. Knepley #endif 69127f02ce8SMatthew G. Knepley } 69227f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 69327f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 694665f567fSMatthew G. Knepley ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 695665f567fSMatthew 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);} 69627f02ce8SMatthew G. Knepley if (f0_func) { 69727f02ce8SMatthew G. Knepley f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q*NcS]); 69827f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w; 69927f02ce8SMatthew G. Knepley } 70027f02ce8SMatthew G. Knepley if (f1_func) { 70127f02ce8SMatthew G. Knepley f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q*NcS*dim]); 70227f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w; 70327f02ce8SMatthew G. Knepley } 70427f02ce8SMatthew G. Knepley } 705665f567fSMatthew G. Knepley if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 706665f567fSMatthew G. Knepley else {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 70727f02ce8SMatthew G. Knepley cOffset += totDim; 70827f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 70927f02ce8SMatthew G. Knepley } 71027f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 71127f02ce8SMatthew G. Knepley } 71227f02ce8SMatthew G. Knepley 7134bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 7144bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 71520cf1dd8SToby Isaac { 71620cf1dd8SToby Isaac const PetscInt debug = 0; 7174bee2e38SMatthew G. Knepley PetscFE feI, feJ; 7184bee2e38SMatthew G. Knepley PetscPointJac g0_func, g1_func, g2_func, g3_func; 71920cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 72020cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 72120cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 72220cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 72320cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 72420cf1dd8SToby Isaac PetscQuadrature quad; 725ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 7264bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 72720cf1dd8SToby Isaac const PetscScalar *constants; 72820cf1dd8SToby Isaac PetscReal *x; 729ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 730ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 73120cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 73220cf1dd8SToby Isaac PetscInt dE, Np; 73320cf1dd8SToby Isaac PetscBool isAffine; 73420cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 73520cf1dd8SToby Isaac PetscInt qNc, Nq, q; 73620cf1dd8SToby Isaac PetscErrorCode ierr; 73720cf1dd8SToby Isaac 73820cf1dd8SToby Isaac PetscFunctionBegin; 7394bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 7404bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 7414bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 7424bee2e38SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 7434bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 7444bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 7454bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 7464bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 74720cf1dd8SToby Isaac switch(jtype) { 7484bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 7494bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 7504bee2e38SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscDSGetJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 75120cf1dd8SToby Isaac } 75220cf1dd8SToby Isaac if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 7534bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 7544bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 7554bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 756ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 7574bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 7584bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 7594bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 7604bee2e38SMatthew G. Knepley if (dsAux) { 7614bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 7624bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 7634bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 7644bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 7654bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 766ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 76701907d53SMatthew 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); 76820cf1dd8SToby Isaac } 76927f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 77027f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7714bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7724bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7734bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 77427f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 77527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 77627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 77727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 77827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 77927f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 78027f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 7814bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7824bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7834bee2e38SMatthew G. Knepley 78427f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 78527f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7864bee2e38SMatthew G. Knepley if (isAffine) { 7874bee2e38SMatthew G. Knepley fegeom.v = x; 7884bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 7897132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e*Np*dE*dE]; 7907132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 7917132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np]; 7924bee2e38SMatthew G. Knepley } 79320cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 79420cf1dd8SToby Isaac PetscReal w; 7954bee2e38SMatthew G. Knepley PetscInt c; 79620cf1dd8SToby Isaac 79720cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 79820cf1dd8SToby Isaac if (isAffine) { 7997132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 80020cf1dd8SToby Isaac } else { 8014bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e*Np+q)*dE]; 8024bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 8034bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 8044bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e*Np+q]; 80520cf1dd8SToby Isaac } 8064bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 807ef0bb6c7SMatthew 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);} 808ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 80920cf1dd8SToby Isaac if (g0_func) { 810580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 8114bee2e38SMatthew G. Knepley g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g0); 81220cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 81320cf1dd8SToby Isaac } 81420cf1dd8SToby Isaac if (g1_func) { 81527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 8164bee2e38SMatthew G. Knepley g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g1); 8174bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 81820cf1dd8SToby Isaac } 81920cf1dd8SToby Isaac if (g2_func) { 82027f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 8214bee2e38SMatthew G. Knepley g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g2); 8224bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 82320cf1dd8SToby Isaac } 82420cf1dd8SToby Isaac if (g3_func) { 82527f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 8264bee2e38SMatthew G. Knepley g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g3); 8274bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 82820cf1dd8SToby Isaac } 82920cf1dd8SToby Isaac 830ef0bb6c7SMatthew 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); 83120cf1dd8SToby Isaac } 83220cf1dd8SToby Isaac if (debug > 1) { 83320cf1dd8SToby Isaac PetscInt fc, f, gc, g; 83420cf1dd8SToby Isaac 83520cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 836ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 837ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 838ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 839ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 840ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 841ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 84220cf1dd8SToby 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); 84320cf1dd8SToby Isaac } 84420cf1dd8SToby Isaac } 84520cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 84620cf1dd8SToby Isaac } 84720cf1dd8SToby Isaac } 84820cf1dd8SToby Isaac } 84920cf1dd8SToby Isaac cOffset += totDim; 85020cf1dd8SToby Isaac cOffsetAux += totDimAux; 85120cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 85220cf1dd8SToby Isaac } 85320cf1dd8SToby Isaac PetscFunctionReturn(0); 85420cf1dd8SToby Isaac } 85520cf1dd8SToby Isaac 8562b99622eSMatthew G. Knepley static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 8574bee2e38SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 85820cf1dd8SToby Isaac { 85920cf1dd8SToby Isaac const PetscInt debug = 0; 8604bee2e38SMatthew G. Knepley PetscFE feI, feJ; 8614bee2e38SMatthew G. Knepley PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 86220cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 86320cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 86420cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 86520cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 86620cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 86720cf1dd8SToby Isaac PetscQuadrature quad; 868ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8694bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 87020cf1dd8SToby Isaac const PetscScalar *constants; 87120cf1dd8SToby Isaac PetscReal *x; 872ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 873ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 87420cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 87520cf1dd8SToby Isaac PetscBool isAffine; 87620cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 87720cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 87820cf1dd8SToby Isaac PetscErrorCode ierr; 87920cf1dd8SToby Isaac 88020cf1dd8SToby Isaac PetscFunctionBegin; 8814bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 8824bee2e38SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 8834bee2e38SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 8844bee2e38SMatthew G. Knepley ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 8854bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 8864bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 8874bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 8884bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 8894bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 8904bee2e38SMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 8914bee2e38SMatthew G. Knepley ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr); 89231f073d7SMatthew G. Knepley if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 8934bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 8944bee2e38SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 8954bee2e38SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 896ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr); 8974bee2e38SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 8984bee2e38SMatthew G. Knepley if (dsAux) { 8994bee2e38SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 9004bee2e38SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 9014bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 9024bee2e38SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 9034bee2e38SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 904ef0bb6c7SMatthew G. Knepley ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr); 90520cf1dd8SToby Isaac } 906ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 90720cf1dd8SToby Isaac Np = fgeom->numPoints; 90820cf1dd8SToby Isaac dE = fgeom->dimEmbed; 90920cf1dd8SToby Isaac isAffine = fgeom->isAffine; 91027f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 91131f073d7SMatthew G. Knepley ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 91227f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 91327f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 91427f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 91527f02ce8SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 91627f02ce8SMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 91720cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 9189f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 91920cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 920*ea78f98cSLisandro Dalcin fegeom.n = NULL; 921*ea78f98cSLisandro Dalcin fegeom.v = NULL; 922*ea78f98cSLisandro Dalcin fegeom.J = NULL; 923*ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 92427f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 92527f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 92627f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 92727f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 9284bee2e38SMatthew G. Knepley if (isAffine) { 9294bee2e38SMatthew G. Knepley fegeom.v = x; 9304bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 9317132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e*Np*dE*dE]; 9327132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 9337132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np]; 9347132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e*Np*dE]; 9359f209ee4SMatthew G. Knepley 9367132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 9377132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 9387132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 9394bee2e38SMatthew G. Knepley } 94020cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 94120cf1dd8SToby Isaac PetscReal w; 9424bee2e38SMatthew G. Knepley PetscInt c; 94320cf1dd8SToby Isaac 94420cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 94520cf1dd8SToby Isaac if (isAffine) { 9467132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 94720cf1dd8SToby Isaac } else { 9483fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e*Np+q)*dE]; 9499f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 9509f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 9514bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 9524bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 9539f209ee4SMatthew G. Knepley 9549f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 9559f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 9569f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 95720cf1dd8SToby Isaac } 9584bee2e38SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 959ef0bb6c7SMatthew 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);} 960ef0bb6c7SMatthew G. Knepley if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 96120cf1dd8SToby Isaac if (g0_func) { 962580bdb30SBarry Smith ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 9634bee2e38SMatthew 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); 96420cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 96520cf1dd8SToby Isaac } 96620cf1dd8SToby Isaac if (g1_func) { 96727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 9684bee2e38SMatthew 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); 9694bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 97020cf1dd8SToby Isaac } 97120cf1dd8SToby Isaac if (g2_func) { 97227f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 9734bee2e38SMatthew 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); 9744bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 97520cf1dd8SToby Isaac } 97620cf1dd8SToby Isaac if (g3_func) { 97727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 9784bee2e38SMatthew 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); 9794bee2e38SMatthew G. Knepley for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 98020cf1dd8SToby Isaac } 98120cf1dd8SToby Isaac 982ef0bb6c7SMatthew 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); 98320cf1dd8SToby Isaac } 98420cf1dd8SToby Isaac if (debug > 1) { 98520cf1dd8SToby Isaac PetscInt fc, f, gc, g; 98620cf1dd8SToby Isaac 98720cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 988ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 989ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 990ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 991ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 992ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 993ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 99420cf1dd8SToby 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); 99520cf1dd8SToby Isaac } 99620cf1dd8SToby Isaac } 99720cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 99820cf1dd8SToby Isaac } 99920cf1dd8SToby Isaac } 100020cf1dd8SToby Isaac } 100120cf1dd8SToby Isaac cOffset += totDim; 100220cf1dd8SToby Isaac cOffsetAux += totDimAux; 100320cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 100420cf1dd8SToby Isaac } 100520cf1dd8SToby Isaac PetscFunctionReturn(0); 100620cf1dd8SToby Isaac } 100720cf1dd8SToby Isaac 100827f02ce8SMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 100927f02ce8SMatthew G. Knepley const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 101027f02ce8SMatthew G. Knepley { 101127f02ce8SMatthew G. Knepley const PetscInt debug = 0; 101227f02ce8SMatthew G. Knepley PetscFE feI, feJ; 1013665f567fSMatthew G. Knepley PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 101427f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 101527f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 101627f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 101727f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 101827f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 1019665f567fSMatthew G. Knepley PetscQuadrature quad; 1020665f567fSMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 102127f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 102227f02ce8SMatthew G. Knepley const PetscScalar *constants; 102327f02ce8SMatthew G. Knepley PetscReal *x; 1024665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1025665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 102627f02ce8SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 1027665f567fSMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE; 102827f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 102927f02ce8SMatthew G. Knepley PetscInt qNc, Nq, q, Np, dE; 103027f02ce8SMatthew G. Knepley PetscErrorCode ierr; 103127f02ce8SMatthew G. Knepley 103227f02ce8SMatthew G. Knepley PetscFunctionBegin; 103327f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 103427f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 103527f02ce8SMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 103627f02ce8SMatthew G. Knepley ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 103727f02ce8SMatthew G. Knepley ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 103827f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 103927f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 104027f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 104127f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 104227f02ce8SMatthew G. Knepley switch(jtype) { 104327f02ce8SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetBdJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 104427f02ce8SMatthew G. Knepley case PETSCFE_JACOBIAN: ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 1045665f567fSMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 104627f02ce8SMatthew G. Knepley } 104727f02ce8SMatthew G. Knepley if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 104827f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 104927f02ce8SMatthew G. Knepley ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 105027f02ce8SMatthew G. Knepley ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 1051665f567fSMatthew G. Knepley ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1052665f567fSMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 1053665f567fSMatthew G. Knepley ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 105427f02ce8SMatthew G. Knepley ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 105527f02ce8SMatthew G. Knepley if (dsAux) { 105627f02ce8SMatthew G. Knepley ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 105727f02ce8SMatthew G. Knepley ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 105827f02ce8SMatthew G. Knepley ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 105927f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 106027f02ce8SMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 106127f02ce8SMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 106201907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 1063665f567fSMatthew G. Knepley if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);} 1064665f567fSMatthew G. Knepley else {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);} 106501907d53SMatthew 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); 106627f02ce8SMatthew G. Knepley } 106727f02ce8SMatthew G. Knepley isCohesiveFieldI = fieldI == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 106827f02ce8SMatthew G. Knepley isCohesiveFieldJ = fieldJ == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 1069665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1070665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 107127f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2*NcI; 107227f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2*NcJ; 107327f02ce8SMatthew G. Knepley Np = fgeom->numPoints; 107427f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 107527f02ce8SMatthew G. Knepley isAffine = fgeom->isAffine; 107627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 107727f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 107827f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 107927f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1080665f567fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1081665f567fSMatthew G. Knepley if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 108227f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 108327f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 108427f02ce8SMatthew G. Knepley const PetscInt face = fgeom->face[e][0]; 108527f02ce8SMatthew G. Knepley 108627f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 108727f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 108827f02ce8SMatthew G. Knepley if (isAffine) { 108927f02ce8SMatthew G. Knepley fegeom.v = x; 109027f02ce8SMatthew G. Knepley fegeom.xi = fgeom->xi; 109127f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[e*dE*dE]; 109227f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e*dE*dE]; 109327f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e]; 109427f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[e*dE]; 109527f02ce8SMatthew G. Knepley } 109627f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 109727f02ce8SMatthew G. Knepley PetscReal w; 109827f02ce8SMatthew G. Knepley PetscInt c; 109927f02ce8SMatthew G. Knepley 110027f02ce8SMatthew G. Knepley if (isAffine) { 110127f02ce8SMatthew G. Knepley /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */ 110227f02ce8SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 110327f02ce8SMatthew G. Knepley } else { 110427f02ce8SMatthew G. Knepley fegeom.v = &fegeom.v[(e*Np+q)*dE]; 110527f02ce8SMatthew G. Knepley fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 110627f02ce8SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 110727f02ce8SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e*Np+q]; 110827f02ce8SMatthew G. Knepley fegeom.n = &fgeom->n[(e*Np+q)*dE]; 110927f02ce8SMatthew G. Knepley } 111027f02ce8SMatthew G. Knepley w = fegeom.detJ[0]*quadWeights[q]; 111127f02ce8SMatthew G. Knepley if (debug > 1 && q < Np) { 111227f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 111327f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 111427f02ce8SMatthew G. Knepley ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 111527f02ce8SMatthew G. Knepley #endif 111627f02ce8SMatthew G. Knepley } 111727f02ce8SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 1118665f567fSMatthew 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);} 1119665f567fSMatthew 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);} 112027f02ce8SMatthew G. Knepley if (g0_func) { 112127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 112227f02ce8SMatthew 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); 112327f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT; ++c) g0[c] *= w; 112427f02ce8SMatthew G. Knepley } 112527f02ce8SMatthew G. Knepley if (g1_func) { 112627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 112727f02ce8SMatthew 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); 112827f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w; 112927f02ce8SMatthew G. Knepley } 113027f02ce8SMatthew G. Knepley if (g2_func) { 113127f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 113227f02ce8SMatthew 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); 113327f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w; 113427f02ce8SMatthew G. Knepley } 113527f02ce8SMatthew G. Knepley if (g3_func) { 113627f02ce8SMatthew G. Knepley ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 113727f02ce8SMatthew 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); 113827f02ce8SMatthew G. Knepley for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w; 113927f02ce8SMatthew G. Knepley } 114027f02ce8SMatthew G. Knepley 114127f02ce8SMatthew G. Knepley if (isCohesiveFieldI && isCohesiveFieldJ) { 1142665f567fSMatthew 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); 114327f02ce8SMatthew G. Knepley } else { 1144665f567fSMatthew 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); 114527f02ce8SMatthew G. Knepley } 114627f02ce8SMatthew G. Knepley } 114727f02ce8SMatthew G. Knepley if (debug > 1) { 114827f02ce8SMatthew G. Knepley PetscInt fc, f, gc, g; 114927f02ce8SMatthew G. Knepley 115027f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 115127f02ce8SMatthew G. Knepley for (fc = 0; fc < NcI; ++fc) { 1152665f567fSMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 115327f02ce8SMatthew G. Knepley const PetscInt i = offsetI + f*NcI+fc; 115427f02ce8SMatthew G. Knepley for (gc = 0; gc < NcJ; ++gc) { 1155665f567fSMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 115627f02ce8SMatthew G. Knepley const PetscInt j = offsetJ + g*NcJ+gc; 115727f02ce8SMatthew 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); 115827f02ce8SMatthew G. Knepley } 115927f02ce8SMatthew G. Knepley } 116027f02ce8SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 116127f02ce8SMatthew G. Knepley } 116227f02ce8SMatthew G. Knepley } 116327f02ce8SMatthew G. Knepley } 116427f02ce8SMatthew G. Knepley cOffset += totDim; 116527f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 116627f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 116727f02ce8SMatthew G. Knepley } 116827f02ce8SMatthew G. Knepley PetscFunctionReturn(0); 116927f02ce8SMatthew G. Knepley } 117027f02ce8SMatthew G. Knepley 11712b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 117220cf1dd8SToby Isaac { 117320cf1dd8SToby Isaac PetscFunctionBegin; 117420cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 117520cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 117620cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 117720cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 117820cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1179ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 118020cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1181afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 118220cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 118320cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 118427f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 118520cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 118620cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 118720cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 118827f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 118920cf1dd8SToby Isaac PetscFunctionReturn(0); 119020cf1dd8SToby Isaac } 119120cf1dd8SToby Isaac 119220cf1dd8SToby Isaac /*MC 119320cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 119420cf1dd8SToby Isaac 119520cf1dd8SToby Isaac Level: intermediate 119620cf1dd8SToby Isaac 119720cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 119820cf1dd8SToby Isaac M*/ 119920cf1dd8SToby Isaac 120020cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 120120cf1dd8SToby Isaac { 120220cf1dd8SToby Isaac PetscFE_Basic *b; 120320cf1dd8SToby Isaac PetscErrorCode ierr; 120420cf1dd8SToby Isaac 120520cf1dd8SToby Isaac PetscFunctionBegin; 120620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 120720cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 120820cf1dd8SToby Isaac fem->data = b; 120920cf1dd8SToby Isaac 121020cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 121120cf1dd8SToby Isaac PetscFunctionReturn(0); 121220cf1dd8SToby Isaac } 1213