120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscblaslapack.h> 320cf1dd8SToby Isaac 4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem) 5d71ae5a4SJacob Faibussowitsch { 620cf1dd8SToby Isaac PetscFE_Basic *b = (PetscFE_Basic *)fem->data; 720cf1dd8SToby Isaac 820cf1dd8SToby Isaac PetscFunctionBegin; 99566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1120cf1dd8SToby Isaac } 1220cf1dd8SToby Isaac 13d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v) 14d71ae5a4SJacob Faibussowitsch { 15d9bac1caSLisandro Dalcin PetscInt dim, Nc; 16d9bac1caSLisandro Dalcin PetscSpace basis = NULL; 17d9bac1caSLisandro Dalcin PetscDualSpace dual = NULL; 18d9bac1caSLisandro Dalcin PetscQuadrature quad = NULL; 1920cf1dd8SToby Isaac 2020cf1dd8SToby Isaac PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 229566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fe, &Nc)); 239566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &basis)); 249566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &dual)); 259566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(v)); 2763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc)); 289566063dSJacob Faibussowitsch if (basis) PetscCall(PetscSpaceView(basis, v)); 299566063dSJacob Faibussowitsch if (dual) PetscCall(PetscDualSpaceView(dual, v)); 309566063dSJacob Faibussowitsch if (quad) PetscCall(PetscQuadratureView(quad, v)); 319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(v)); 323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3320cf1dd8SToby Isaac } 3420cf1dd8SToby Isaac 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v) 36d71ae5a4SJacob Faibussowitsch { 3720cf1dd8SToby Isaac PetscBool iascii; 3820cf1dd8SToby Isaac 3920cf1dd8SToby Isaac PetscFunctionBegin; 409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 419566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscFEView_Basic_Ascii(fe, v)); 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4320cf1dd8SToby Isaac } 4420cf1dd8SToby Isaac 4520cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */ 46d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem) 47d71ae5a4SJacob Faibussowitsch { 48b9d4cb8dSJed Brown PetscReal *work; 4920cf1dd8SToby Isaac PetscBLASInt *pivots; 5020cf1dd8SToby Isaac PetscBLASInt n, info; 5120cf1dd8SToby Isaac PetscInt pdim, j; 5220cf1dd8SToby Isaac 5320cf1dd8SToby Isaac PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pdim * pdim, &fem->invV)); 5620cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 5720cf1dd8SToby Isaac PetscReal *Bf; 5820cf1dd8SToby Isaac PetscQuadrature f; 5920cf1dd8SToby Isaac const PetscReal *points, *weights; 6020cf1dd8SToby Isaac PetscInt Nc, Nq, q, k, c; 6120cf1dd8SToby Isaac 629566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f)); 639566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights)); 649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf)); 659566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL)); 6620cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 6720cf1dd8SToby Isaac /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */ 68b9d4cb8dSJed Brown fem->invV[j * pdim + k] = 0.0; 6920cf1dd8SToby Isaac 7020cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 71b9d4cb8dSJed Brown for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c]; 7220cf1dd8SToby Isaac } 7320cf1dd8SToby Isaac } 749566063dSJacob Faibussowitsch PetscCall(PetscFree(Bf)); 7520cf1dd8SToby Isaac } 76ea2bdf6dSBarry Smith 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work)); 786497c311SBarry Smith PetscCall(PetscBLASIntCast(pdim, &n)); 79792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info)); 80835f2295SStefano Zampini PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info); 81792fecdfSBarry Smith PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info)); 82835f2295SStefano Zampini PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info); 839566063dSJacob Faibussowitsch PetscCall(PetscFree2(pivots, work)); 843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8520cf1dd8SToby Isaac } 8620cf1dd8SToby Isaac 87d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim) 88d71ae5a4SJacob Faibussowitsch { 8920cf1dd8SToby Isaac PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim)); 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9220cf1dd8SToby Isaac } 9320cf1dd8SToby Isaac 94b9d4cb8dSJed Brown /* Tensor contraction on the middle index, 95b9d4cb8dSJed Brown * C[m,n,p] = A[m,k,p] * B[k,n] 96b9d4cb8dSJed Brown * where all matrices use C-style ordering. 97b9d4cb8dSJed Brown */ 98d71ae5a4SJacob Faibussowitsch static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C) 99d71ae5a4SJacob Faibussowitsch { 100b9d4cb8dSJed Brown PetscInt i; 101b9d4cb8dSJed Brown 102b9d4cb8dSJed Brown PetscFunctionBegin; 103aa9788aaSMatthew G. Knepley PetscCheck(n && p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Empty tensor is not allowed %" PetscInt_FMT " %" PetscInt_FMT, n, p); 104b9d4cb8dSJed Brown for (i = 0; i < m; i++) { 105b9d4cb8dSJed Brown PetscBLASInt n_, p_, k_, lda, ldb, ldc; 106b9d4cb8dSJed Brown PetscReal one = 1, zero = 0; 107b9d4cb8dSJed Brown /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n] 108b9d4cb8dSJed Brown * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k) 109b9d4cb8dSJed Brown */ 1109566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &n_)); 1119566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(p, &p_)); 1129566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(k, &k_)); 113b9d4cb8dSJed Brown lda = p_; 114b9d4cb8dSJed Brown ldb = n_; 115b9d4cb8dSJed Brown ldc = p_; 116792fecdfSBarry Smith PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc)); 117b9d4cb8dSJed Brown } 1189566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2. * m * n * p * k)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120b9d4cb8dSJed Brown } 121b9d4cb8dSJed Brown 122ce78bad3SBarry Smith PETSC_INTERN PetscErrorCode PetscFEComputeTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 123d71ae5a4SJacob Faibussowitsch { 12420cf1dd8SToby Isaac DM dm; 12520cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 12620cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 12720cf1dd8SToby Isaac PetscInt Nc; /* Field components */ 128ef0bb6c7SMatthew G. Knepley PetscReal *B = K >= 0 ? T->T[0] : NULL; 129ef0bb6c7SMatthew G. Knepley PetscReal *D = K >= 1 ? T->T[1] : NULL; 130ef0bb6c7SMatthew G. Knepley PetscReal *H = K >= 2 ? T->T[2] : NULL; 131ef0bb6c7SMatthew G. Knepley PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL; 13220cf1dd8SToby Isaac 13320cf1dd8SToby Isaac PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 1359566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1369566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fem, &Nc)); 13820cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 1399566063dSJacob Faibussowitsch if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB)); 1409566063dSJacob Faibussowitsch if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD)); 1419566063dSJacob Faibussowitsch if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH)); 143b9d4cb8dSJed Brown /* Translate from prime to nodal basis */ 14420cf1dd8SToby Isaac if (B) { 145b9d4cb8dSJed Brown /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */ 1469566063dSJacob Faibussowitsch PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B)); 14720cf1dd8SToby Isaac } 148aa9788aaSMatthew G. Knepley if (D && dim) { 149b9d4cb8dSJed Brown /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */ 1509566063dSJacob Faibussowitsch PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D)); 15120cf1dd8SToby Isaac } 152aa9788aaSMatthew G. Knepley if (H && dim) { 153b9d4cb8dSJed Brown /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */ 1549566063dSJacob Faibussowitsch PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H)); 15520cf1dd8SToby Isaac } 1569566063dSJacob Faibussowitsch if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB)); 1579566063dSJacob Faibussowitsch if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD)); 1589566063dSJacob Faibussowitsch if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH)); 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16020cf1dd8SToby Isaac } 16120cf1dd8SToby Isaac 1622dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 163d71ae5a4SJacob Faibussowitsch { 164b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 1654bee2e38SMatthew G. Knepley PetscFE fe; 16620cf1dd8SToby Isaac PetscPointFunc obj_func; 16720cf1dd8SToby Isaac PetscQuadrature quad; 168ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 1694bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x; 17020cf1dd8SToby Isaac const PetscScalar *constants; 17187510d7dSMatthew G. Knepley PetscReal *x, cellScale; 172ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 17320cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 17420cf1dd8SToby Isaac PetscBool isAffine; 17520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 17620cf1dd8SToby Isaac PetscInt qNc, Nq, q; 17720cf1dd8SToby Isaac 17820cf1dd8SToby Isaac PetscFunctionBegin; 1799566063dSJacob Faibussowitsch PetscCall(PetscDSGetObjective(ds, field, &obj_func)); 1803ba16761SJacob Faibussowitsch if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 1819566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 1829566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 18387510d7dSMatthew G. Knepley cellScale = (PetscReal)PetscPowInt(2, dim); 1849566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 1859566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1869566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1879566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 1889566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 1899566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 1909566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 1919566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL)); 19287510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 1939566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 1944bee2e38SMatthew G. Knepley if (dsAux) { 1959566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 1969566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 20163a3b9bcSJacob Faibussowitsch PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 20220cf1dd8SToby Isaac } 2039566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 20463a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 20520cf1dd8SToby Isaac Np = cgeom->numPoints; 20620cf1dd8SToby Isaac dE = cgeom->dimEmbed; 20720cf1dd8SToby Isaac isAffine = cgeom->isAffine; 20820cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 2094bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 21020cf1dd8SToby Isaac 21127f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 21227f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 213*a7be0fb8SMatthew G. Knepley fegeom.xi = NULL; 21420cf1dd8SToby Isaac if (isAffine) { 2154bee2e38SMatthew G. Knepley fegeom.v = x; 2164bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2177132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 2187132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 2197132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 22020cf1dd8SToby Isaac } 2214bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 222d627b919SMatthew G. Knepley PetscScalar integrand = 0.; 2234bee2e38SMatthew G. Knepley PetscReal w; 2244bee2e38SMatthew G. Knepley 2254bee2e38SMatthew G. Knepley if (isAffine) { 2267132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 2274bee2e38SMatthew G. Knepley } else { 2284bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 2294bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 2304bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 2314bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 2324bee2e38SMatthew G. Knepley } 23387510d7dSMatthew G. Knepley PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 2344bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 23520cf1dd8SToby Isaac if (debug > 1 && q < Np) { 23663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 2377be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 2389566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 23920cf1dd8SToby Isaac #endif 24020cf1dd8SToby Isaac } 24163a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 2429566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 2439566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 2444bee2e38SMatthew 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); 2454bee2e38SMatthew G. Knepley integrand *= w; 24620cf1dd8SToby Isaac integral[e * Nf + field] += integrand; 24720cf1dd8SToby Isaac } 248699b48e5SStefano Zampini if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Element Field %" PetscInt_FMT " integral: %g\n", Nf, (double)PetscRealPart(integral[e * Nf + field]))); 24920cf1dd8SToby Isaac cOffset += totDim; 25020cf1dd8SToby Isaac cOffsetAux += totDimAux; 25120cf1dd8SToby Isaac } 2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25320cf1dd8SToby Isaac } 25420cf1dd8SToby Isaac 2552dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFunc obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 256d71ae5a4SJacob Faibussowitsch { 257b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 2584bee2e38SMatthew G. Knepley PetscFE fe; 259afe6d6adSToby Isaac PetscQuadrature quad; 260ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 2614bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 262afe6d6adSToby Isaac const PetscScalar *constants; 26387510d7dSMatthew G. Knepley PetscReal *x, cellScale; 264ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 265afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 266afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 267afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 268afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 269afe6d6adSToby Isaac 270afe6d6adSToby Isaac PetscFunctionBegin; 2713ba16761SJacob Faibussowitsch if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 2729566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 2739566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 27487510d7dSMatthew G. Knepley cellScale = (PetscReal)PetscPowInt(2, dim); 2759566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 2789566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 2799566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 2809566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 2819566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 2829566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 28387510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 2854bee2e38SMatthew G. Knepley if (dsAux) { 2869566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 2879566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 2889566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 2899566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 2909566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 2919566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 292afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 2939566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 2949566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 29563a3b9bcSJacob Faibussowitsch PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 296afe6d6adSToby Isaac } 2979566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 29863a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 29979ab67a3SMatthew G. Knepley if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Field: %" PetscInt_FMT " Nface: %" PetscInt_FMT " Nq: %" PetscInt_FMT "\n", field, Ne, Nq)); 300afe6d6adSToby Isaac Np = fgeom->numPoints; 301afe6d6adSToby Isaac dE = fgeom->dimEmbed; 302afe6d6adSToby Isaac isAffine = fgeom->isAffine; 303afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 3049f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 305afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 306ea78f98cSLisandro Dalcin fegeom.n = NULL; 307ea78f98cSLisandro Dalcin fegeom.v = NULL; 308*a7be0fb8SMatthew G. Knepley fegeom.xi = NULL; 309ea78f98cSLisandro Dalcin fegeom.J = NULL; 310b2deab97SMatthew G. Knepley fegeom.invJ = NULL; 311ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 31227f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 31327f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 31427f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 31527f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3164bee2e38SMatthew G. Knepley if (isAffine) { 3174bee2e38SMatthew G. Knepley fegeom.v = x; 3184bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3197132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 3207132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 3217132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 3227132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 3239f209ee4SMatthew G. Knepley 3247132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 3257132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 3267132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 3274bee2e38SMatthew G. Knepley } 328afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 32979ab67a3SMatthew G. Knepley PetscScalar integrand = 0.; 3304bee2e38SMatthew G. Knepley PetscReal w; 331afe6d6adSToby Isaac 332afe6d6adSToby Isaac if (isAffine) { 3337132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 334afe6d6adSToby Isaac } else { 3353fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 3369f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 3379f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 3384bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 3394bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 3409f209ee4SMatthew G. Knepley 3419f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 3429f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 3439f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 344afe6d6adSToby Isaac } 34587510d7dSMatthew G. Knepley PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 3464bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 347afe6d6adSToby Isaac if (debug > 1 && q < Np) { 34863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 349afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3509566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 351afe6d6adSToby Isaac #endif 352afe6d6adSToby Isaac } 35363a3b9bcSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 35479ab67a3SMatthew G. Knepley if (debug > 3) { 35579ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, " x_q (")); 35679ab67a3SMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) { 35779ab67a3SMatthew G. Knepley if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 35879ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.v[d])); 35979ab67a3SMatthew G. Knepley } 36079ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 36179ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, " n_q (")); 36279ab67a3SMatthew G. Knepley for (PetscInt d = 0; d < dE; ++d) { 36379ab67a3SMatthew G. Knepley if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 36479ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.n[d])); 36579ab67a3SMatthew G. Knepley } 36679ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 36779ab67a3SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) { 36879ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, " u_%" PetscInt_FMT " (", f)); 36979ab67a3SMatthew G. Knepley for (PetscInt c = 0; c < uOff[f + 1] - uOff[f]; ++c) { 37079ab67a3SMatthew G. Knepley if (c) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 37179ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(u[uOff[f] + c]))); 37279ab67a3SMatthew G. Knepley } 37379ab67a3SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 37479ab67a3SMatthew G. Knepley } 37579ab67a3SMatthew G. Knepley } 3769566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 3779566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 3784bee2e38SMatthew 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); 3794bee2e38SMatthew G. Knepley integrand *= w; 380afe6d6adSToby Isaac integral[e * Nf + field] += integrand; 38179ab67a3SMatthew G. Knepley if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g tot: %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field]))); 382afe6d6adSToby Isaac } 383afe6d6adSToby Isaac cOffset += totDim; 384afe6d6adSToby Isaac cOffsetAux += totDimAux; 385afe6d6adSToby Isaac } 3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 387afe6d6adSToby Isaac } 388afe6d6adSToby Isaac 389d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 390d71ae5a4SJacob Faibussowitsch { 391b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 3926528b96dSMatthew G. Knepley const PetscInt field = key.field; 3934bee2e38SMatthew G. Knepley PetscFE fe; 3946528b96dSMatthew G. Knepley PetscWeakForm wf; 3956528b96dSMatthew G. Knepley PetscInt n0, n1, i; 3966528b96dSMatthew G. Knepley PetscPointFunc *f0_func, *f1_func; 39720cf1dd8SToby Isaac PetscQuadrature quad; 398ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 3994bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 40020cf1dd8SToby Isaac const PetscScalar *constants; 40187510d7dSMatthew G. Knepley PetscReal *x, cellScale; 402ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 403ef0bb6c7SMatthew G. Knepley PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 40420cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 4056587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 40620cf1dd8SToby Isaac 40720cf1dd8SToby Isaac PetscFunctionBegin; 4089566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 41087510d7dSMatthew G. Knepley cellScale = (PetscReal)PetscPowInt(2, dim); 4119566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 4129566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 4139566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 4149566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 4159566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 4169566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 4179566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 4189566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 4193ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 4209566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 4219566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 4229566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 4239566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 42487510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 4259566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 4264bee2e38SMatthew G. Knepley if (dsAux) { 4279566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 4289566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 4299566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 4309566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 4319566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 4329566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 43363a3b9bcSJacob Faibussowitsch PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 43420cf1dd8SToby Isaac } 4359566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 43663a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 43720cf1dd8SToby Isaac dE = cgeom->dimEmbed; 43863a3b9bcSJacob Faibussowitsch PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim); 43920cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4404bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 44120cf1dd8SToby Isaac 4426587ee25SMatthew G. Knepley fegeom.v = x; /* workspace */ 4439566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc)); 4449566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE)); 44520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4464bee2e38SMatthew G. Knepley PetscReal w; 4474bee2e38SMatthew G. Knepley PetscInt c, d; 44820cf1dd8SToby Isaac 4499566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom)); 45087510d7dSMatthew G. Knepley PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 4514bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 4526587ee25SMatthew G. Knepley if (debug > 1 && q < cgeom->numPoints) { 45363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 4547be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 455ac9d17c7SMatthew G. Knepley PetscCall(DMPrintCellMatrix(e, "invJ", dE, dE, fegeom.invJ)); 45620cf1dd8SToby Isaac #endif 45720cf1dd8SToby Isaac } 45816cd844bSPierre Jolivet PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 4599566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 460ac9d17c7SMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dE, 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]); 461ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w; 462ac9d17c7SMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dE, 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 * dE]); 4639371c9d4SSatish Balay for (c = 0; c < T[field]->Nc; ++c) 464ac9d17c7SMatthew G. Knepley for (d = 0; d < dE; ++d) f1[(q * T[field]->Nc + c) * dE + d] *= w; 465b8025e53SMatthew G. Knepley if (debug) { 466e8e188d2SZach Atkins // LCOV_EXCL_START 467e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q])); 468e8e188d2SZach Atkins for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c])); 469e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 470b8025e53SMatthew G. Knepley if (debug > 2) { 47163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " field %" PetscInt_FMT ":", field)); 47263a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c]))); 4739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 474e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " field der %" PetscInt_FMT ":", field)); 475e8e188d2SZach Atkins for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c]))); 476e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 47763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " resid %" PetscInt_FMT ":", field)); 47863a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c]))); 4799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 480e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " res der %" PetscInt_FMT ":", field)); 481e8e188d2SZach Atkins for (c = 0; c < T[field]->Nc; ++c) { 482ac9d17c7SMatthew G. Knepley for (d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dE + d]))); 483b8025e53SMatthew G. Knepley } 484e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 485e8e188d2SZach Atkins } 486e8e188d2SZach Atkins // LCOV_EXCL_STOP 487b8025e53SMatthew G. Knepley } 48820cf1dd8SToby Isaac } 4899566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset])); 49020cf1dd8SToby Isaac cOffset += totDim; 49120cf1dd8SToby Isaac cOffsetAux += totDimAux; 49220cf1dd8SToby Isaac } 4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49420cf1dd8SToby Isaac } 49520cf1dd8SToby Isaac 496d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 497d71ae5a4SJacob Faibussowitsch { 498b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 49906d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 5004bee2e38SMatthew G. Knepley PetscFE fe; 50106d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 50206d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 50320cf1dd8SToby Isaac PetscQuadrature quad; 504ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 5054bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 50620cf1dd8SToby Isaac const PetscScalar *constants; 50787510d7dSMatthew G. Knepley PetscReal *x, cellScale; 508ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 509ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 5106587ee25SMatthew G. Knepley PetscBool auxOnBd = PETSC_FALSE; 51120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 5126587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 51320cf1dd8SToby Isaac 51420cf1dd8SToby Isaac PetscFunctionBegin; 5159566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 5169566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 51787510d7dSMatthew G. Knepley cellScale = (PetscReal)PetscPowInt(2, dim); 5189566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 5199566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 5209566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 5219566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 5229566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 5239566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 5249566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 5253ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 5269566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 5279566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 5289566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 5299566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 53087510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 5319566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 5324bee2e38SMatthew G. Knepley if (dsAux) { 5339566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 5349566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 5359566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 5369566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 5379566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 5389566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 5397be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 5409566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 5419566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 54263a3b9bcSJacob Faibussowitsch PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 54320cf1dd8SToby Isaac } 544ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 5459566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 54663a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 54720cf1dd8SToby Isaac dE = fgeom->dimEmbed; 5486587ee25SMatthew G. Knepley /* TODO FIX THIS */ 5496587ee25SMatthew G. Knepley fgeom->dim = dim - 1; 55063a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 55120cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5529f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 55320cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5549f209ee4SMatthew G. Knepley 5556587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 5569566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcI)); 5579566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcI * dE)); 55820cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5594bee2e38SMatthew G. Knepley PetscReal w; 5604bee2e38SMatthew G. Knepley PetscInt c, d; 5614bee2e38SMatthew G. Knepley 5629566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom)); 5639566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom)); 56487510d7dSMatthew G. Knepley PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 5654bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 56662bd480fSMatthew G. Knepley if (debug > 1) { 5676587ee25SMatthew G. Knepley if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 56863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 5697be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5709566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 5719566063dSJacob Faibussowitsch PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n)); 57220cf1dd8SToby Isaac #endif 57320cf1dd8SToby Isaac } 57462bd480fSMatthew G. Knepley } 5758e3a54c0SPierre Jolivet PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 5769566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 577ac9d17c7SMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q * NcI]); 5784bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w; 579ac9d17c7SMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dE, 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 * dE]); 5809371c9d4SSatish Balay for (c = 0; c < NcI; ++c) 581ac9d17c7SMatthew G. Knepley for (d = 0; d < dE; ++d) f1[(q * NcI + c) * dE + d] *= w; 58262bd480fSMatthew G. Knepley if (debug) { 58363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q)); 58462bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 58563a3b9bcSJacob Faibussowitsch if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c]))); 58662bd480fSMatthew G. Knepley if (n1) { 58763a3b9bcSJacob Faibussowitsch for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f1[%" PetscInt_FMT ",%" PetscInt_FMT "] %g", c, d, (double)PetscRealPart(f1[(q * NcI + c) * dim + d]))); 5889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 58962bd480fSMatthew G. Knepley } 59062bd480fSMatthew G. Knepley } 59162bd480fSMatthew G. Knepley } 59220cf1dd8SToby Isaac } 5939566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 59420cf1dd8SToby Isaac cOffset += totDim; 59520cf1dd8SToby Isaac cOffsetAux += totDimAux; 59620cf1dd8SToby Isaac } 5973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59820cf1dd8SToby Isaac } 59920cf1dd8SToby Isaac 60027f02ce8SMatthew G. Knepley /* 60127f02ce8SMatthew 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 60227f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 60327f02ce8SMatthew G. Knepley 60427f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 60527f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 60627f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 60727f02ce8SMatthew 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() 60827f02ce8SMatthew G. Knepley */ 6092dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 610d71ae5a4SJacob Faibussowitsch { 611b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 6126528b96dSMatthew G. Knepley const PetscInt field = key.field; 61327f02ce8SMatthew G. Knepley PetscFE fe; 6146528b96dSMatthew G. Knepley PetscWeakForm wf; 6156528b96dSMatthew G. Knepley PetscInt n0, n1, i; 6166528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 61727f02ce8SMatthew G. Knepley PetscQuadrature quad; 6180e18dc48SMatthew G. Knepley DMPolytopeType ct; 61907218a29SMatthew G. Knepley PetscTabulation *Tf, *TfIn, *TfAux = NULL; 62027f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 62127f02ce8SMatthew G. Knepley const PetscScalar *constants; 62227f02ce8SMatthew G. Knepley PetscReal *x; 623665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 62407218a29SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 6256587ee25SMatthew G. Knepley PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 62627f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 6276587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 62827f02ce8SMatthew G. Knepley 62927f02ce8SMatthew G. Knepley PetscFunctionBegin; 63027f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 6319566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 6329566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 6339566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 6349566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 6359566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 63607218a29SMatthew G. Knepley PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn)); 637429ebbe4SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets 63807218a29SMatthew G. Knepley PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x)); 6399566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset)); 6409566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 6419566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 6423ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 6439566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 6449566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 6459566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 64627f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 6479566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &Tf)); 64807218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 64987510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 6509566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 65127f02ce8SMatthew G. Knepley if (dsAux) { 6529566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 6539566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 6549566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 6559566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 6569566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 6579566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 65801907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 6599566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 6609566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 66163a3b9bcSJacob Faibussowitsch PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 66227f02ce8SMatthew G. Knepley } 6639566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField)); 664665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 665c2b7495fSMatthew G. Knepley NcS = NcI; 6660abb75b6SMatthew G. Knepley if (!isCohesiveField && s == 2) { 6670abb75b6SMatthew G. Knepley // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 6680abb75b6SMatthew G. Knepley NcS *= 2; 6690abb75b6SMatthew G. Knepley } 6709566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 6710e18dc48SMatthew G. Knepley PetscCall(PetscQuadratureGetCellType(quad, &ct)); 67263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 67327f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 67463a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 67527f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 67627f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 6770e18dc48SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 6780e18dc48SMatthew G. Knepley const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 6794e913f38SMatthew G. Knepley const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 68027f02ce8SMatthew G. Knepley 6816587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 6829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcS)); 6839566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcS * dE)); 68427f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 6850e18dc48SMatthew G. Knepley PetscInt qpt[2]; 68627f02ce8SMatthew G. Knepley PetscReal w; 68727f02ce8SMatthew G. Knepley PetscInt c, d; 68827f02ce8SMatthew G. Knepley 6894e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0])); 6904e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1])); 69107218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 69227f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 6936587ee25SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 69463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 69527f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 6969566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ)); 69727f02ce8SMatthew G. Knepley #endif 69827f02ce8SMatthew G. Knepley } 699a4158a15SMatthew G. Knepley if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0])); 70027f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 7018e3a54c0SPierre Jolivet PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], PetscSafePointerPlusOffset(coefficients_t, cOffsetIn), u, u_x, u_t)); 70207218a29SMatthew G. Knepley if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 703ac9d17c7SMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dE, 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]); 70427f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w; 705ac9d17c7SMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dE, 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 * dE]); 7069371c9d4SSatish Balay for (c = 0; c < NcS; ++c) 7079371c9d4SSatish Balay for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w; 70827f02ce8SMatthew G. Knepley } 7099371c9d4SSatish Balay if (isCohesiveField) { 7103ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 7119371c9d4SSatish Balay } else { 7123ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 7139371c9d4SSatish Balay } 71427f02ce8SMatthew G. Knepley cOffset += totDim; 71507218a29SMatthew G. Knepley cOffsetIn += totDimIn; 71627f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 71727f02ce8SMatthew G. Knepley } 7183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71927f02ce8SMatthew G. Knepley } 72027f02ce8SMatthew G. Knepley 721d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 722d71ae5a4SJacob Faibussowitsch { 723b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 7244bee2e38SMatthew G. Knepley PetscFE feI, feJ; 7256528b96dSMatthew G. Knepley PetscWeakForm wf; 7266528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 7276528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 72820cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 72920cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 73020cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 73120cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 73220cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 73320cf1dd8SToby Isaac PetscQuadrature quad; 734ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 7351a768569SStefano Zampini PetscScalar *g0 = NULL, *g1 = NULL, *g2 = NULL, *g3 = NULL, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 73620cf1dd8SToby Isaac const PetscScalar *constants; 73787510d7dSMatthew G. Knepley PetscReal *x, cellScale; 738ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 739ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 7406528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 74120cf1dd8SToby Isaac PetscInt dE, Np; 74220cf1dd8SToby Isaac PetscBool isAffine; 74320cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 74420cf1dd8SToby Isaac PetscInt qNc, Nq, q; 74520cf1dd8SToby Isaac 74620cf1dd8SToby Isaac PetscFunctionBegin; 7479566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 7486528b96dSMatthew G. Knepley fieldI = key.field / Nf; 7496528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 7509566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 7519566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 7529566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 75387510d7dSMatthew G. Knepley cellScale = (PetscReal)PetscPowInt(2, dim); 7549566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 7559566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 7569566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 7579566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 7589566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 75920cf1dd8SToby Isaac switch (jtype) { 760d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 761d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 762d71ae5a4SJacob Faibussowitsch break; 763d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 764d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 765d71ae5a4SJacob Faibussowitsch break; 766d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 767d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 768d71ae5a4SJacob Faibussowitsch break; 76920cf1dd8SToby Isaac } 7703ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 7719566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 7729566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 7731a768569SStefano Zampini PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, n0 ? &g0 : NULL, n1 ? &g1 : NULL, n2 ? &g2 : NULL, n3 ? &g3 : NULL)); 7741a768569SStefano Zampini 7759566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 7769566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 7779566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 778a102dd69SStefano Zampini PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ)); 7799566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 7804bee2e38SMatthew G. Knepley if (dsAux) { 7819566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 7829566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 7839566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 7849566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 7859566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 7869566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 78763a3b9bcSJacob Faibussowitsch PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 78820cf1dd8SToby Isaac } 78927f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 79027f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7914bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7924bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7934bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 7949566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 79563a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 7961a768569SStefano Zampini 7974bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7984bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7994bee2e38SMatthew G. Knepley 80027f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 80127f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 802*a7be0fb8SMatthew G. Knepley fegeom.xi = NULL; 8034bee2e38SMatthew G. Knepley if (isAffine) { 8044bee2e38SMatthew G. Knepley fegeom.v = x; 8054bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 8067132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 8077132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 8087132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 8094bee2e38SMatthew G. Knepley } 81020cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 81120cf1dd8SToby Isaac PetscReal w; 8124bee2e38SMatthew G. Knepley PetscInt c; 81320cf1dd8SToby Isaac 81420cf1dd8SToby Isaac if (isAffine) { 8157132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 81620cf1dd8SToby Isaac } else { 8174bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 8184bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 8194bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 8204bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 82120cf1dd8SToby Isaac } 82287510d7dSMatthew G. Knepley PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 8239566063dSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0])); 8244bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 82516cd844bSPierre Jolivet if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 8269566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 827ea672e62SMatthew G. Knepley if (n0) { 8289566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 829ac9d17c7SMatthew G. Knepley for (i = 0; i < n0; ++i) g0_func[i](dE, 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); 83020cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 83120cf1dd8SToby Isaac } 832ea672e62SMatthew G. Knepley if (n1) { 8339566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 834ac9d17c7SMatthew G. Knepley for (i = 0; i < n1; ++i) g1_func[i](dE, 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); 8359b529ac9SStefano Zampini for (c = 0; c < NcI * NcJ * dE; ++c) g1[c] *= w; 83620cf1dd8SToby Isaac } 837ea672e62SMatthew G. Knepley if (n2) { 8389566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 839ac9d17c7SMatthew G. Knepley for (i = 0; i < n2; ++i) g2_func[i](dE, 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); 8409b529ac9SStefano Zampini for (c = 0; c < NcI * NcJ * dE; ++c) g2[c] *= w; 84120cf1dd8SToby Isaac } 842ea672e62SMatthew G. Knepley if (n3) { 8439566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 844ac9d17c7SMatthew G. Knepley for (i = 0; i < n3; ++i) g3_func[i](dE, 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); 8459b529ac9SStefano Zampini for (c = 0; c < NcI * NcJ * dE * dE; ++c) g3[c] *= w; 84620cf1dd8SToby Isaac } 84720cf1dd8SToby Isaac 8481a768569SStefano Zampini PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset)); 84920cf1dd8SToby Isaac } 85020cf1dd8SToby Isaac if (debug > 1) { 851699b48e5SStefano Zampini PetscInt f, g; 85220cf1dd8SToby Isaac 85363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 854ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 855699b48e5SStefano Zampini const PetscInt i = offsetI + f; 856ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 857699b48e5SStefano Zampini const PetscInt j = offsetJ + g; 858699b48e5SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ", %" PetscInt_FMT "]: %g\n", f, g, (double)PetscRealPart(elemMat[eOffset + i * totDim + j]))); 85920cf1dd8SToby Isaac } 8609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 86120cf1dd8SToby Isaac } 86220cf1dd8SToby Isaac } 86320cf1dd8SToby Isaac cOffset += totDim; 86420cf1dd8SToby Isaac cOffsetAux += totDimAux; 86520cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 86620cf1dd8SToby Isaac } 8673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86820cf1dd8SToby Isaac } 86920cf1dd8SToby Isaac 870e3d591f2SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 871d71ae5a4SJacob Faibussowitsch { 872b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 8734bee2e38SMatthew G. Knepley PetscFE feI, feJ; 87445480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 87545480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 87620cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 87720cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 87820cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 87920cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 88020cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 88120cf1dd8SToby Isaac PetscQuadrature quad; 882ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8834bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 88420cf1dd8SToby Isaac const PetscScalar *constants; 88587510d7dSMatthew G. Knepley PetscReal *x, cellScale; 886ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 887ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 88845480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 88920cf1dd8SToby Isaac PetscBool isAffine; 89020cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 89120cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 89220cf1dd8SToby Isaac 89320cf1dd8SToby Isaac PetscFunctionBegin; 8949566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 89545480ffeSMatthew G. Knepley fieldI = key.field / Nf; 89645480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8979566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 8989566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 8999566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 90087510d7dSMatthew G. Knepley cellScale = (PetscReal)PetscPowInt(2, dim); 9019566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(feI, &quad)); 9029566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 9039566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 9049566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 9059566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 9069566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 907e3d591f2SMatthew G. Knepley switch (jtype) { 908e3d591f2SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: 909e3d591f2SMatthew G. Knepley PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 910e3d591f2SMatthew G. Knepley break; 911e3d591f2SMatthew G. Knepley case PETSCFE_JACOBIAN: 9129566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 913e3d591f2SMatthew G. Knepley break; 914e3d591f2SMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: 915e3d591f2SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()"); 916e3d591f2SMatthew G. Knepley } 9173ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 9189566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 9199566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 9209566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 9219566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &T)); 92287510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ)); 9239566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 9244bee2e38SMatthew G. Knepley if (dsAux) { 9259566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 9269566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 9279566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 9289566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 9299566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 9309566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 93120cf1dd8SToby Isaac } 932ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 93320cf1dd8SToby Isaac Np = fgeom->numPoints; 93420cf1dd8SToby Isaac dE = fgeom->dimEmbed; 93520cf1dd8SToby Isaac isAffine = fgeom->isAffine; 93627f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 9379566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 9389566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 9399566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 9409566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 9419566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 94263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 94320cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 9449f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 94520cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 946ea78f98cSLisandro Dalcin fegeom.n = NULL; 947ea78f98cSLisandro Dalcin fegeom.v = NULL; 948*a7be0fb8SMatthew G. Knepley fegeom.xi = NULL; 949ea78f98cSLisandro Dalcin fegeom.J = NULL; 950ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 95127f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 95227f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 95327f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 95427f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 9554bee2e38SMatthew G. Knepley if (isAffine) { 9564bee2e38SMatthew G. Knepley fegeom.v = x; 9574bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 9587132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 9597132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 9607132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 9617132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 9629f209ee4SMatthew G. Knepley 9637132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 9647132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 9657132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 9664bee2e38SMatthew G. Knepley } 96720cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 96820cf1dd8SToby Isaac PetscReal w; 9694bee2e38SMatthew G. Knepley PetscInt c; 97020cf1dd8SToby Isaac 97163a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 97220cf1dd8SToby Isaac if (isAffine) { 9737132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 97420cf1dd8SToby Isaac } else { 9753fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 9769f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 9779f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 9784bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 9794bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 9809f209ee4SMatthew G. Knepley 9819f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 9829f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 9839f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 98420cf1dd8SToby Isaac } 98587510d7dSMatthew G. Knepley PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 9864bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 9879566063dSJacob Faibussowitsch if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 9889566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 989ea672e62SMatthew G. Knepley if (n0) { 9909566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 991ac9d17c7SMatthew G. Knepley for (i = 0; i < n0; ++i) g0_func[i](dE, 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); 99220cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 99320cf1dd8SToby Isaac } 994ea672e62SMatthew G. Knepley if (n1) { 9959566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 996ac9d17c7SMatthew G. Knepley for (i = 0; i < n1; ++i) g1_func[i](dE, 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); 9974bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 99820cf1dd8SToby Isaac } 999ea672e62SMatthew G. Knepley if (n2) { 10009566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 1001ac9d17c7SMatthew G. Knepley for (i = 0; i < n2; ++i) g2_func[i](dE, 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); 10024bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 100320cf1dd8SToby Isaac } 1004ea672e62SMatthew G. Knepley if (n3) { 10059566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 1006ac9d17c7SMatthew G. Knepley for (i = 0; i < n3; ++i) g3_func[i](dE, 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); 10074bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 100820cf1dd8SToby Isaac } 100920cf1dd8SToby Isaac 10101a768569SStefano Zampini PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset)); 101120cf1dd8SToby Isaac } 101220cf1dd8SToby Isaac if (debug > 1) { 101320cf1dd8SToby Isaac PetscInt fc, f, gc, g; 101420cf1dd8SToby Isaac 101563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 1016ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 1017ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 1018ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 1019ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 1020ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 1021ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 102263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f, fc, g, gc, (double)PetscRealPart(elemMat[eOffset + i * totDim + j]))); 102320cf1dd8SToby Isaac } 102420cf1dd8SToby Isaac } 10259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 102620cf1dd8SToby Isaac } 102720cf1dd8SToby Isaac } 102820cf1dd8SToby Isaac } 102920cf1dd8SToby Isaac cOffset += totDim; 103020cf1dd8SToby Isaac cOffsetAux += totDimAux; 103120cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 103220cf1dd8SToby Isaac } 10333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103420cf1dd8SToby Isaac } 103520cf1dd8SToby Isaac 10362dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1037d71ae5a4SJacob Faibussowitsch { 1038b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 103927f02ce8SMatthew G. Knepley PetscFE feI, feJ; 1040148442b3SMatthew G. Knepley PetscWeakForm wf; 1041148442b3SMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 1042148442b3SMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 104327f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 104427f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 104527f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 104627f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 104727f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 1048665f567fSMatthew G. Knepley PetscQuadrature quad; 10490e18dc48SMatthew G. Knepley DMPolytopeType ct; 105007218a29SMatthew G. Knepley PetscTabulation *T, *TfIn, *TAux = NULL; 105127f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 105227f02ce8SMatthew G. Knepley const PetscScalar *constants; 105327f02ce8SMatthew G. Knepley PetscReal *x; 1054665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1055665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 105645480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 105707218a29SMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE; 105827f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 10590502970dSMatthew G. Knepley PetscInt qNc, Nq, q; 106027f02ce8SMatthew G. Knepley 106127f02ce8SMatthew G. Knepley PetscFunctionBegin; 10629566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 106345480ffeSMatthew G. Knepley fieldI = key.field / Nf; 106445480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 106527f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 10669566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 10679566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 10689566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 10699566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 10709566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1071429ebbe4SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets 10729566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x)); 10739566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 107427f02ce8SMatthew G. Knepley switch (jtype) { 1075d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 1076d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1077d71ae5a4SJacob Faibussowitsch break; 1078d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 1079d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1080d71ae5a4SJacob Faibussowitsch break; 1081d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 1082d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 108327f02ce8SMatthew G. Knepley } 10843ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 10859566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 10869566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 10879566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 10889566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 108907218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 10909566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI)); 10919566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ)); 109287510d7dSMatthew G. Knepley PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ)); 10939566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 109427f02ce8SMatthew G. Knepley if (dsAux) { 10959566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 10969566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 10979566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 10989566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 10999566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 11009566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 110101907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 11029566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 11039566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 110463a3b9bcSJacob Faibussowitsch PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 110527f02ce8SMatthew G. Knepley } 11069566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI)); 11079566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ)); 1108665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1109665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 111027f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2 * NcI; 111127f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ; 11120abb75b6SMatthew G. Knepley if (!isCohesiveFieldI && s == 2) { 11130abb75b6SMatthew G. Knepley // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 11140abb75b6SMatthew G. Knepley NcS *= 2; 11150abb75b6SMatthew G. Knepley } 11160abb75b6SMatthew G. Knepley if (!isCohesiveFieldJ && s == 2) { 11170abb75b6SMatthew G. Knepley // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 11180abb75b6SMatthew G. Knepley NcT *= 2; 11190abb75b6SMatthew G. Knepley } 11200502970dSMatthew G. Knepley // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though 11210502970dSMatthew G. Knepley // the coordinates are in dE dimensions 11229566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 11230502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 11240502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 11250502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 11269566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 11270e18dc48SMatthew G. Knepley PetscCall(PetscQuadratureGetCellType(quad, &ct)); 112863a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 112927f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 113027f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 11310e18dc48SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 11320e18dc48SMatthew G. Knepley const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 11334e913f38SMatthew G. Knepley const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 113427f02ce8SMatthew G. Knepley 113507218a29SMatthew G. Knepley fegeom.v = x; /* Workspace */ 113627f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 11370e18dc48SMatthew G. Knepley PetscInt qpt[2]; 113827f02ce8SMatthew G. Knepley PetscReal w; 113927f02ce8SMatthew G. Knepley PetscInt c; 114027f02ce8SMatthew G. Knepley 11414e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0])); 11424e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1])); 114307218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 114427f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 114507218a29SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 114663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 114727f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 11489566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 114927f02ce8SMatthew G. Knepley #endif 115027f02ce8SMatthew G. Knepley } 115163a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 11528e3a54c0SPierre Jolivet if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 115307218a29SMatthew G. Knepley if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 1154ea672e62SMatthew G. Knepley if (n0) { 11559566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 1156148442b3SMatthew 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); 115727f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT; ++c) g0[c] *= w; 115827f02ce8SMatthew G. Knepley } 1159ea672e62SMatthew G. Knepley if (n1) { 11600502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 1161148442b3SMatthew 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); 11620502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w; 116327f02ce8SMatthew G. Knepley } 1164ea672e62SMatthew G. Knepley if (n2) { 11650502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 1166148442b3SMatthew 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); 11670502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w; 116827f02ce8SMatthew G. Knepley } 1169ea672e62SMatthew G. Knepley if (n3) { 11700502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 1171148442b3SMatthew 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); 11720502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w; 117327f02ce8SMatthew G. Knepley } 117427f02ce8SMatthew G. Knepley 11755fedec97SMatthew G. Knepley if (isCohesiveFieldI) { 11765fedec97SMatthew G. Knepley if (isCohesiveFieldJ) { 11771a768569SStefano Zampini PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset)); 117827f02ce8SMatthew G. Knepley } else { 11790abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 11800abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat)); 11810abb75b6SMatthew G. Knepley } 11820abb75b6SMatthew G. Knepley } else { 11830abb75b6SMatthew G. Knepley if (s == 2) { 11840abb75b6SMatthew G. Knepley if (isCohesiveFieldJ) { 11850abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 11860abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat)); 11870abb75b6SMatthew G. Knepley } else { 11880abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 11890abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat)); 11900abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 2], &g1[NcI * NcJ * dim * 2], &g2[NcI * NcJ * dim * 2], &g3[NcI * NcJ * dim * dim * 2], eOffset, totDim, offsetI, offsetJ, elemMat)); 11910abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 3], &g1[NcI * NcJ * dim * 3], &g2[NcI * NcJ * dim * 3], &g3[NcI * NcJ * dim * dim * 3], eOffset, totDim, offsetI, offsetJ, elemMat)); 11925fedec97SMatthew G. Knepley } 11939371c9d4SSatish Balay } else 11940abb75b6SMatthew G. Knepley PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 11950abb75b6SMatthew G. Knepley } 119627f02ce8SMatthew G. Knepley } 119727f02ce8SMatthew G. Knepley if (debug > 1) { 11984e913f38SMatthew G. Knepley const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb)); 11994e913f38SMatthew G. Knepley const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb)); 12004e913f38SMatthew G. Knepley const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb)); 12014e913f38SMatthew G. Knepley const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb)); 12024e913f38SMatthew G. Knepley PetscInt f, g; 120327f02ce8SMatthew G. Knepley 12044e913f38SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT " s %s totDim %" PetscInt_FMT " offsets (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", fieldI, fieldJ, s ? (s > 1 ? "Coh" : "Pos") : "Neg", totDim, eOffset, offsetI, offsetJ)); 12054e913f38SMatthew G. Knepley for (f = fS; f < fE; ++f) { 12064e913f38SMatthew G. Knepley const PetscInt i = offsetI + f; 12074e913f38SMatthew G. Knepley for (g = gS; g < gE; ++g) { 12084e913f38SMatthew G. Knepley const PetscInt j = offsetJ + g; 1209e978a55eSPierre Jolivet PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, f, i, g, j); 12104e913f38SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f / NcI, f % NcI, g / NcJ, g % NcJ, (double)PetscRealPart(elemMat[eOffset + i * totDim + j]))); 121127f02ce8SMatthew G. Knepley } 12129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 121327f02ce8SMatthew G. Knepley } 121427f02ce8SMatthew G. Knepley } 121527f02ce8SMatthew G. Knepley cOffset += totDim; 121627f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 121727f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 121827f02ce8SMatthew G. Knepley } 12193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122027f02ce8SMatthew G. Knepley } 122127f02ce8SMatthew G. Knepley 1222d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 1223d71ae5a4SJacob Faibussowitsch { 122420cf1dd8SToby Isaac PetscFunctionBegin; 122520cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 122620cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 122720cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 122820cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 122920cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1230ce78bad3SBarry Smith fem->ops->computetabulation = PetscFEComputeTabulation_Basic; 123120cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1232afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 123320cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 123420cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 123527f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 123620cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */; 123720cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 123820cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 123927f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 12403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124120cf1dd8SToby Isaac } 124220cf1dd8SToby Isaac 124320cf1dd8SToby Isaac /*MC 1244dce8aebaSBarry Smith PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization 124520cf1dd8SToby Isaac 124620cf1dd8SToby Isaac Level: intermediate 124720cf1dd8SToby Isaac 1248dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 124920cf1dd8SToby Isaac M*/ 125020cf1dd8SToby Isaac 1251d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 1252d71ae5a4SJacob Faibussowitsch { 125320cf1dd8SToby Isaac PetscFE_Basic *b; 125420cf1dd8SToby Isaac 125520cf1dd8SToby Isaac PetscFunctionBegin; 125620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 12574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 125820cf1dd8SToby Isaac fem->data = b; 125920cf1dd8SToby Isaac 12609566063dSJacob Faibussowitsch PetscCall(PetscFEInitialize_Basic(fem)); 12613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126220cf1dd8SToby Isaac } 1263