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)); 7820cf1dd8SToby Isaac n = pdim; 79792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info)); 8063a3b9bcSJacob Faibussowitsch PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscInt_FMT, (PetscInt)info); 81792fecdfSBarry Smith PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info)); 8263a3b9bcSJacob Faibussowitsch PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscInt_FMT, (PetscInt)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 122d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFECreateTabulation_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 { 164*b2deab97SMatthew 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; 17120cf1dd8SToby Isaac PetscReal *x; 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)); 1839566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 1849566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1859566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1869566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 1879566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 1889566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 1899566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 1909566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL)); 1919566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 1924bee2e38SMatthew G. Knepley if (dsAux) { 1939566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 1949566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 1959566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 1969566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 19963a3b9bcSJacob 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); 20020cf1dd8SToby Isaac } 2019566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 20263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 20320cf1dd8SToby Isaac Np = cgeom->numPoints; 20420cf1dd8SToby Isaac dE = cgeom->dimEmbed; 20520cf1dd8SToby Isaac isAffine = cgeom->isAffine; 20620cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 2074bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 20820cf1dd8SToby Isaac 20927f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 21027f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 21120cf1dd8SToby Isaac if (isAffine) { 2124bee2e38SMatthew G. Knepley fegeom.v = x; 2134bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2147132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 2157132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 2167132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 21720cf1dd8SToby Isaac } 2184bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 219d627b919SMatthew G. Knepley PetscScalar integrand = 0.; 2204bee2e38SMatthew G. Knepley PetscReal w; 2214bee2e38SMatthew G. Knepley 2224bee2e38SMatthew G. Knepley if (isAffine) { 2237132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 2244bee2e38SMatthew G. Knepley } else { 2254bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 2264bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 2274bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 2284bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 2294bee2e38SMatthew G. Knepley } 2304bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 23120cf1dd8SToby Isaac if (debug > 1 && q < Np) { 23263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 2337be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 2349566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 23520cf1dd8SToby Isaac #endif 23620cf1dd8SToby Isaac } 23763a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 2389566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 2399566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 2404bee2e38SMatthew 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); 2414bee2e38SMatthew G. Knepley integrand *= w; 24220cf1dd8SToby Isaac integral[e * Nf + field] += integrand; 2439566063dSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[field]))); 24420cf1dd8SToby Isaac } 24520cf1dd8SToby Isaac cOffset += totDim; 24620cf1dd8SToby Isaac cOffsetAux += totDimAux; 24720cf1dd8SToby Isaac } 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24920cf1dd8SToby Isaac } 25020cf1dd8SToby Isaac 2512dce792eSToby 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[]) 252d71ae5a4SJacob Faibussowitsch { 253*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 2544bee2e38SMatthew G. Knepley PetscFE fe; 255afe6d6adSToby Isaac PetscQuadrature quad; 256ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 2574bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 258afe6d6adSToby Isaac const PetscScalar *constants; 259afe6d6adSToby Isaac PetscReal *x; 260ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 261afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 262afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 263afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 264afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 265afe6d6adSToby Isaac 266afe6d6adSToby Isaac PetscFunctionBegin; 2673ba16761SJacob Faibussowitsch if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 2689566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 2709566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 2719566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2729566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 2739566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 2749566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 2759566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 2789566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 2794bee2e38SMatthew G. Knepley if (dsAux) { 2809566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 2819566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 2829566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 2839566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 2859566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 286afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 2879566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 2889566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 28963a3b9bcSJacob 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); 290afe6d6adSToby Isaac } 2919566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 29263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 293afe6d6adSToby Isaac Np = fgeom->numPoints; 294afe6d6adSToby Isaac dE = fgeom->dimEmbed; 295afe6d6adSToby Isaac isAffine = fgeom->isAffine; 296afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 2979f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 298afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 299ea78f98cSLisandro Dalcin fegeom.n = NULL; 300ea78f98cSLisandro Dalcin fegeom.v = NULL; 301ea78f98cSLisandro Dalcin fegeom.J = NULL; 302*b2deab97SMatthew G. Knepley fegeom.invJ = NULL; 303ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 30427f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 30527f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 30627f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 30727f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3084bee2e38SMatthew G. Knepley if (isAffine) { 3094bee2e38SMatthew G. Knepley fegeom.v = x; 3104bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3117132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 3127132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 3137132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 3147132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 3159f209ee4SMatthew G. Knepley 3167132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 3177132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 3187132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 3194bee2e38SMatthew G. Knepley } 320afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 321afe6d6adSToby Isaac PetscScalar integrand; 3224bee2e38SMatthew G. Knepley PetscReal w; 323afe6d6adSToby Isaac 324afe6d6adSToby Isaac if (isAffine) { 3257132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 326afe6d6adSToby Isaac } else { 3273fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 3289f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 3299f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 3304bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 3314bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 3329f209ee4SMatthew G. Knepley 3339f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 3349f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 3359f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 336afe6d6adSToby Isaac } 3374bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 338afe6d6adSToby Isaac if (debug > 1 && q < Np) { 33963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 340afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3419566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 342afe6d6adSToby Isaac #endif 343afe6d6adSToby Isaac } 34463a3b9bcSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 3459566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 3469566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 3474bee2e38SMatthew 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); 3484bee2e38SMatthew G. Knepley integrand *= w; 349afe6d6adSToby Isaac integral[e * Nf + field] += integrand; 3509566063dSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field]))); 351afe6d6adSToby Isaac } 352afe6d6adSToby Isaac cOffset += totDim; 353afe6d6adSToby Isaac cOffsetAux += totDimAux; 354afe6d6adSToby Isaac } 3553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 356afe6d6adSToby Isaac } 357afe6d6adSToby Isaac 358d71ae5a4SJacob 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[]) 359d71ae5a4SJacob Faibussowitsch { 360*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 3616528b96dSMatthew G. Knepley const PetscInt field = key.field; 3624bee2e38SMatthew G. Knepley PetscFE fe; 3636528b96dSMatthew G. Knepley PetscWeakForm wf; 3646528b96dSMatthew G. Knepley PetscInt n0, n1, i; 3656528b96dSMatthew G. Knepley PetscPointFunc *f0_func, *f1_func; 36620cf1dd8SToby Isaac PetscQuadrature quad; 367ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 3684bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 36920cf1dd8SToby Isaac const PetscScalar *constants; 37020cf1dd8SToby Isaac PetscReal *x; 371ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 372ef0bb6c7SMatthew G. Knepley PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 37320cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3746587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 37520cf1dd8SToby Isaac 37620cf1dd8SToby Isaac PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 3799566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 3809566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 3819566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 3829566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 3839566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 3849566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 3859566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 3869566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 3873ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 3889566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 3909566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 3919566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 3929566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 3934bee2e38SMatthew G. Knepley if (dsAux) { 3949566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 3959566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 3969566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 3979566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 3989566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 3999566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 40063a3b9bcSJacob 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); 40120cf1dd8SToby Isaac } 4029566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 40363a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 40420cf1dd8SToby Isaac dE = cgeom->dimEmbed; 40563a3b9bcSJacob Faibussowitsch PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim); 40620cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4074bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 40820cf1dd8SToby Isaac 4096587ee25SMatthew G. Knepley fegeom.v = x; /* workspace */ 4109566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc)); 4119566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE)); 41220cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4134bee2e38SMatthew G. Knepley PetscReal w; 4144bee2e38SMatthew G. Knepley PetscInt c, d; 41520cf1dd8SToby Isaac 4169566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom)); 4174bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 4186587ee25SMatthew G. Knepley if (debug > 1 && q < cgeom->numPoints) { 41963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 4207be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4219566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 42220cf1dd8SToby Isaac #endif 42320cf1dd8SToby Isaac } 42416cd844bSPierre Jolivet PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 4259566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 4266528b96dSMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q * T[field]->Nc]); 427ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w; 4286528b96dSMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q * T[field]->Nc * dim]); 4299371c9d4SSatish Balay for (c = 0; c < T[field]->Nc; ++c) 4309371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w; 431b8025e53SMatthew G. Knepley if (debug) { 432e8e188d2SZach Atkins // LCOV_EXCL_START 433e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q])); 434e8e188d2SZach Atkins for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c])); 435e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 436b8025e53SMatthew G. Knepley if (debug > 2) { 43763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " field %" PetscInt_FMT ":", field)); 43863a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c]))); 4399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 440e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " field der %" PetscInt_FMT ":", field)); 441e8e188d2SZach Atkins for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c]))); 442e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 44363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " resid %" PetscInt_FMT ":", field)); 44463a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c]))); 4459566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 446e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " res der %" PetscInt_FMT ":", field)); 447e8e188d2SZach Atkins for (c = 0; c < T[field]->Nc; ++c) { 448e8e188d2SZach Atkins for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dim + d]))); 449b8025e53SMatthew G. Knepley } 450e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 451e8e188d2SZach Atkins } 452e8e188d2SZach Atkins // LCOV_EXCL_STOP 453b8025e53SMatthew G. Knepley } 45420cf1dd8SToby Isaac } 4559566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset])); 45620cf1dd8SToby Isaac cOffset += totDim; 45720cf1dd8SToby Isaac cOffsetAux += totDimAux; 45820cf1dd8SToby Isaac } 4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46020cf1dd8SToby Isaac } 46120cf1dd8SToby Isaac 462d71ae5a4SJacob 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[]) 463d71ae5a4SJacob Faibussowitsch { 464*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 46506d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 4664bee2e38SMatthew G. Knepley PetscFE fe; 46706d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 46806d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 46920cf1dd8SToby Isaac PetscQuadrature quad; 470ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4714bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 47220cf1dd8SToby Isaac const PetscScalar *constants; 47320cf1dd8SToby Isaac PetscReal *x; 474ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 475ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 4766587ee25SMatthew G. Knepley PetscBool auxOnBd = PETSC_FALSE; 47720cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 4786587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 47920cf1dd8SToby Isaac 48020cf1dd8SToby Isaac PetscFunctionBegin; 4819566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 4829566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 4839566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 4849566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 4859566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 4869566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 4879566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 4889566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 4899566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 4903ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 4919566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 4929566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 4939566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 4949566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 4959566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 4964bee2e38SMatthew G. Knepley if (dsAux) { 4979566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 4989566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 4999566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 5009566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 5019566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 5029566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 5037be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 5049566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 5059566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 50663a3b9bcSJacob 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); 50720cf1dd8SToby Isaac } 508ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 5099566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 51063a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 51120cf1dd8SToby Isaac dE = fgeom->dimEmbed; 5126587ee25SMatthew G. Knepley /* TODO FIX THIS */ 5136587ee25SMatthew G. Knepley fgeom->dim = dim - 1; 51463a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 51520cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5169f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 51720cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5189f209ee4SMatthew G. Knepley 5196587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 5209566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcI)); 5219566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcI * dE)); 52220cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5234bee2e38SMatthew G. Knepley PetscReal w; 5244bee2e38SMatthew G. Knepley PetscInt c, d; 5254bee2e38SMatthew G. Knepley 5269566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom)); 5279566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom)); 5284bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 52962bd480fSMatthew G. Knepley if (debug > 1) { 5306587ee25SMatthew G. Knepley if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 53163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 5327be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5339566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 5349566063dSJacob Faibussowitsch PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n)); 53520cf1dd8SToby Isaac #endif 53620cf1dd8SToby Isaac } 53762bd480fSMatthew G. Knepley } 5388e3a54c0SPierre Jolivet PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 5399566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 54006d8a0d3SMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q * NcI]); 5414bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w; 54206d8a0d3SMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q * NcI * dim]); 5439371c9d4SSatish Balay for (c = 0; c < NcI; ++c) 5449371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w; 54562bd480fSMatthew G. Knepley if (debug) { 54663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q)); 54762bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 54863a3b9bcSJacob Faibussowitsch if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c]))); 54962bd480fSMatthew G. Knepley if (n1) { 55063a3b9bcSJacob 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]))); 5519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 55262bd480fSMatthew G. Knepley } 55362bd480fSMatthew G. Knepley } 55462bd480fSMatthew G. Knepley } 55520cf1dd8SToby Isaac } 5569566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 55720cf1dd8SToby Isaac cOffset += totDim; 55820cf1dd8SToby Isaac cOffsetAux += totDimAux; 55920cf1dd8SToby Isaac } 5603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56120cf1dd8SToby Isaac } 56220cf1dd8SToby Isaac 56327f02ce8SMatthew G. Knepley /* 56427f02ce8SMatthew 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 56527f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 56627f02ce8SMatthew G. Knepley 56727f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 56827f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 56927f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 57027f02ce8SMatthew 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() 57127f02ce8SMatthew G. Knepley */ 5722dce792eSToby 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[]) 573d71ae5a4SJacob Faibussowitsch { 574*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 5756528b96dSMatthew G. Knepley const PetscInt field = key.field; 57627f02ce8SMatthew G. Knepley PetscFE fe; 5776528b96dSMatthew G. Knepley PetscWeakForm wf; 5786528b96dSMatthew G. Knepley PetscInt n0, n1, i; 5796528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 58027f02ce8SMatthew G. Knepley PetscQuadrature quad; 5810e18dc48SMatthew G. Knepley DMPolytopeType ct; 58207218a29SMatthew G. Knepley PetscTabulation *Tf, *TfIn, *TfAux = NULL; 58327f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 58427f02ce8SMatthew G. Knepley const PetscScalar *constants; 58527f02ce8SMatthew G. Knepley PetscReal *x; 586665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 58707218a29SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 5886587ee25SMatthew G. Knepley PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 58927f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 5906587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 59127f02ce8SMatthew G. Knepley 59227f02ce8SMatthew G. Knepley PetscFunctionBegin; 59327f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 5949566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 5959566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 5969566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 5979566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 5989566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 59907218a29SMatthew G. Knepley PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn)); 600429ebbe4SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets 60107218a29SMatthew G. Knepley PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x)); 6029566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset)); 6039566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 6049566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 6053ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 6069566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 6079566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 6089566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 60927f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 6109566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &Tf)); 61107218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 6129566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 61327f02ce8SMatthew G. Knepley if (dsAux) { 6149566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 6159566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 6169566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 6179566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 6189566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 6199566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 62001907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 6219566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 6229566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 62363a3b9bcSJacob 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); 62427f02ce8SMatthew G. Knepley } 6259566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField)); 626665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 627c2b7495fSMatthew G. Knepley NcS = NcI; 6280abb75b6SMatthew G. Knepley if (!isCohesiveField && s == 2) { 6290abb75b6SMatthew G. Knepley // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 6300abb75b6SMatthew G. Knepley NcS *= 2; 6310abb75b6SMatthew G. Knepley } 6329566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 6330e18dc48SMatthew G. Knepley PetscCall(PetscQuadratureGetCellType(quad, &ct)); 63463a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 63527f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 63663a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 63727f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 63827f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 6390e18dc48SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 6400e18dc48SMatthew G. Knepley const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 6414e913f38SMatthew G. Knepley const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 64227f02ce8SMatthew G. Knepley 6436587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 6449566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcS)); 6459566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcS * dE)); 64627f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 6470e18dc48SMatthew G. Knepley PetscInt qpt[2]; 64827f02ce8SMatthew G. Knepley PetscReal w; 64927f02ce8SMatthew G. Knepley PetscInt c, d; 65027f02ce8SMatthew G. Knepley 6514e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0])); 6524e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1])); 65307218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 65427f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 6556587ee25SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 65663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 65727f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 6589566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ)); 65927f02ce8SMatthew G. Knepley #endif 66027f02ce8SMatthew G. Knepley } 661a4158a15SMatthew 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])); 66227f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 6638e3a54c0SPierre 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)); 66407218a29SMatthew 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)); 6656528b96dSMatthew G. Knepley for (i = 0; i < n0; ++i) f0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q * NcS]); 66627f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w; 6679ee2af8cSMatthew G. Knepley for (i = 0; i < n1; ++i) f1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q * NcS * dE]); 6689371c9d4SSatish Balay for (c = 0; c < NcS; ++c) 6699371c9d4SSatish Balay for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w; 67027f02ce8SMatthew G. Knepley } 6719371c9d4SSatish Balay if (isCohesiveField) { 6723ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 6739371c9d4SSatish Balay } else { 6743ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 6759371c9d4SSatish Balay } 67627f02ce8SMatthew G. Knepley cOffset += totDim; 67707218a29SMatthew G. Knepley cOffsetIn += totDimIn; 67827f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 67927f02ce8SMatthew G. Knepley } 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68127f02ce8SMatthew G. Knepley } 68227f02ce8SMatthew G. Knepley 683d71ae5a4SJacob 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[]) 684d71ae5a4SJacob Faibussowitsch { 685*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 6864bee2e38SMatthew G. Knepley PetscFE feI, feJ; 6876528b96dSMatthew G. Knepley PetscWeakForm wf; 6886528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 6896528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 69020cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 69120cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 69220cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 69320cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 69420cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 69520cf1dd8SToby Isaac PetscQuadrature quad; 696ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 6974bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 69820cf1dd8SToby Isaac const PetscScalar *constants; 69920cf1dd8SToby Isaac PetscReal *x; 700ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 701ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 7026528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 70320cf1dd8SToby Isaac PetscInt dE, Np; 70420cf1dd8SToby Isaac PetscBool isAffine; 70520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 70620cf1dd8SToby Isaac PetscInt qNc, Nq, q; 70720cf1dd8SToby Isaac 70820cf1dd8SToby Isaac PetscFunctionBegin; 7099566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 7106528b96dSMatthew G. Knepley fieldI = key.field / Nf; 7116528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 7129566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 7139566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 7149566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 7159566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 7169566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 7179566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 7189566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 7199566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 72020cf1dd8SToby Isaac switch (jtype) { 721d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 722d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 723d71ae5a4SJacob Faibussowitsch break; 724d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 725d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 726d71ae5a4SJacob Faibussowitsch break; 727d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 728d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 729d71ae5a4SJacob Faibussowitsch break; 73020cf1dd8SToby Isaac } 7313ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 7329566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 7339566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 7349566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 7359566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 7369566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 7379566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 7389566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 7394bee2e38SMatthew G. Knepley if (dsAux) { 7409566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 7419566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 7429566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 7439566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 7449566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 7459566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 74663a3b9bcSJacob 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); 74720cf1dd8SToby Isaac } 74827f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 74927f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7504bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7514bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7524bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 75327f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 7549566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 7559566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 7569566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 7579566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 7589566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 75963a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 7604bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7614bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7624bee2e38SMatthew G. Knepley 76327f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 76427f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7654bee2e38SMatthew G. Knepley if (isAffine) { 7664bee2e38SMatthew G. Knepley fegeom.v = x; 7674bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 7687132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 7697132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 7707132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 7714bee2e38SMatthew G. Knepley } 77220cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 77320cf1dd8SToby Isaac PetscReal w; 7744bee2e38SMatthew G. Knepley PetscInt c; 77520cf1dd8SToby Isaac 77620cf1dd8SToby Isaac if (isAffine) { 7777132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 77820cf1dd8SToby Isaac } else { 7794bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 7804bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 7814bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 7824bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 78320cf1dd8SToby Isaac } 7849566063dSJacob 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])); 7854bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 78616cd844bSPierre Jolivet if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 7879566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 788ea672e62SMatthew G. Knepley if (n0) { 7899566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 7906528b96dSMatthew G. Knepley for (i = 0; i < n0; ++i) g0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g0); 79120cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 79220cf1dd8SToby Isaac } 793ea672e62SMatthew G. Knepley if (n1) { 7949566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 7956528b96dSMatthew G. Knepley for (i = 0; i < n1; ++i) g1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g1); 7964bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 79720cf1dd8SToby Isaac } 798ea672e62SMatthew G. Knepley if (n2) { 7999566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 8006528b96dSMatthew G. Knepley for (i = 0; i < n2; ++i) g2_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g2); 8014bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 80220cf1dd8SToby Isaac } 803ea672e62SMatthew G. Knepley if (n3) { 8049566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 8056528b96dSMatthew G. Knepley for (i = 0; i < n3; ++i) g3_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g3); 8064bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 80720cf1dd8SToby Isaac } 80820cf1dd8SToby Isaac 8099566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 81020cf1dd8SToby Isaac } 81120cf1dd8SToby Isaac if (debug > 1) { 81220cf1dd8SToby Isaac PetscInt fc, f, gc, g; 81320cf1dd8SToby Isaac 81463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 815ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 816ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 817ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 818ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 819ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 820ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 82163a3b9bcSJacob 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]))); 82220cf1dd8SToby Isaac } 82320cf1dd8SToby Isaac } 8249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 82520cf1dd8SToby Isaac } 82620cf1dd8SToby Isaac } 82720cf1dd8SToby Isaac } 82820cf1dd8SToby Isaac cOffset += totDim; 82920cf1dd8SToby Isaac cOffsetAux += totDimAux; 83020cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 83120cf1dd8SToby Isaac } 8323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 83320cf1dd8SToby Isaac } 83420cf1dd8SToby Isaac 835e3d591f2SMatthew 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[]) 836d71ae5a4SJacob Faibussowitsch { 837*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 8384bee2e38SMatthew G. Knepley PetscFE feI, feJ; 83945480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 84045480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 84120cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 84220cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 84320cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 84420cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 84520cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 84620cf1dd8SToby Isaac PetscQuadrature quad; 847ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8484bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 84920cf1dd8SToby Isaac const PetscScalar *constants; 85020cf1dd8SToby Isaac PetscReal *x; 851ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 852ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 85345480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 85420cf1dd8SToby Isaac PetscBool isAffine; 85520cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 85620cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 85720cf1dd8SToby Isaac 85820cf1dd8SToby Isaac PetscFunctionBegin; 8599566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 86045480ffeSMatthew G. Knepley fieldI = key.field / Nf; 86145480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8629566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 8639566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 8649566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 8659566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(feI, &quad)); 8669566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 8679566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 8689566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 8699566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 8709566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 871e3d591f2SMatthew G. Knepley switch (jtype) { 872e3d591f2SMatthew G. Knepley case PETSCFE_JACOBIAN_PRE: 873e3d591f2SMatthew 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)); 874e3d591f2SMatthew G. Knepley break; 875e3d591f2SMatthew G. Knepley case PETSCFE_JACOBIAN: 8769566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 877e3d591f2SMatthew G. Knepley break; 878e3d591f2SMatthew G. Knepley case PETSCFE_JACOBIAN_DYN: 879e3d591f2SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()"); 880e3d591f2SMatthew G. Knepley } 8813ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 8829566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 8839566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 8849566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 8859566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &T)); 8869566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 8874bee2e38SMatthew G. Knepley if (dsAux) { 8889566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 8899566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 8909566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 8919566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 8929566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 8939566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 89420cf1dd8SToby Isaac } 895ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 89620cf1dd8SToby Isaac Np = fgeom->numPoints; 89720cf1dd8SToby Isaac dE = fgeom->dimEmbed; 89820cf1dd8SToby Isaac isAffine = fgeom->isAffine; 89927f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 9009566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 9019566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 9029566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 9039566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 9049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 90563a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 90620cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 9079f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 90820cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 909ea78f98cSLisandro Dalcin fegeom.n = NULL; 910ea78f98cSLisandro Dalcin fegeom.v = NULL; 911ea78f98cSLisandro Dalcin fegeom.J = NULL; 912ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 91327f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 91427f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 91527f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 91627f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 9174bee2e38SMatthew G. Knepley if (isAffine) { 9184bee2e38SMatthew G. Knepley fegeom.v = x; 9194bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 9207132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 9217132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 9227132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 9237132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 9249f209ee4SMatthew G. Knepley 9257132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 9267132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 9277132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 9284bee2e38SMatthew G. Knepley } 92920cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 93020cf1dd8SToby Isaac PetscReal w; 9314bee2e38SMatthew G. Knepley PetscInt c; 93220cf1dd8SToby Isaac 93363a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 93420cf1dd8SToby Isaac if (isAffine) { 9357132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 93620cf1dd8SToby Isaac } else { 9373fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 9389f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 9399f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 9404bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 9414bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 9429f209ee4SMatthew G. Knepley 9439f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 9449f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 9459f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 94620cf1dd8SToby Isaac } 9474bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 9489566063dSJacob Faibussowitsch if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 9499566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 950ea672e62SMatthew G. Knepley if (n0) { 9519566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 95245480ffeSMatthew 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); 95320cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 95420cf1dd8SToby Isaac } 955ea672e62SMatthew G. Knepley if (n1) { 9569566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 95745480ffeSMatthew 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); 9584bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 95920cf1dd8SToby Isaac } 960ea672e62SMatthew G. Knepley if (n2) { 9619566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 96245480ffeSMatthew 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); 9634bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 96420cf1dd8SToby Isaac } 965ea672e62SMatthew G. Knepley if (n3) { 9669566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 96745480ffeSMatthew 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); 9684bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 96920cf1dd8SToby Isaac } 97020cf1dd8SToby Isaac 9719566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 97220cf1dd8SToby Isaac } 97320cf1dd8SToby Isaac if (debug > 1) { 97420cf1dd8SToby Isaac PetscInt fc, f, gc, g; 97520cf1dd8SToby Isaac 97663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 977ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 978ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 979ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 980ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 981ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 982ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 98363a3b9bcSJacob 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]))); 98420cf1dd8SToby Isaac } 98520cf1dd8SToby Isaac } 9869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 98720cf1dd8SToby Isaac } 98820cf1dd8SToby Isaac } 98920cf1dd8SToby Isaac } 99020cf1dd8SToby Isaac cOffset += totDim; 99120cf1dd8SToby Isaac cOffsetAux += totDimAux; 99220cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 99320cf1dd8SToby Isaac } 9943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 99520cf1dd8SToby Isaac } 99620cf1dd8SToby Isaac 9972dce792eSToby 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[]) 998d71ae5a4SJacob Faibussowitsch { 999*b2deab97SMatthew G. Knepley const PetscInt debug = ds->printIntegrate; 100027f02ce8SMatthew G. Knepley PetscFE feI, feJ; 1001148442b3SMatthew G. Knepley PetscWeakForm wf; 1002148442b3SMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 1003148442b3SMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 100427f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 100527f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 100627f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 100727f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 100827f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 1009665f567fSMatthew G. Knepley PetscQuadrature quad; 10100e18dc48SMatthew G. Knepley DMPolytopeType ct; 101107218a29SMatthew G. Knepley PetscTabulation *T, *TfIn, *TAux = NULL; 101227f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 101327f02ce8SMatthew G. Knepley const PetscScalar *constants; 101427f02ce8SMatthew G. Knepley PetscReal *x; 1015665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1016665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 101745480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 101807218a29SMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE; 101927f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 10200502970dSMatthew G. Knepley PetscInt qNc, Nq, q; 102127f02ce8SMatthew G. Knepley 102227f02ce8SMatthew G. Knepley PetscFunctionBegin; 10239566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 102445480ffeSMatthew G. Knepley fieldI = key.field / Nf; 102545480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 102627f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 10279566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 10289566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 10299566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 10309566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 10319566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1032429ebbe4SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets 10339566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x)); 10349566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 103527f02ce8SMatthew G. Knepley switch (jtype) { 1036d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 1037d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1038d71ae5a4SJacob Faibussowitsch break; 1039d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 1040d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1041d71ae5a4SJacob Faibussowitsch break; 1042d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 1043d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 104427f02ce8SMatthew G. Knepley } 10453ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 10469566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 10479566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 10489566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 10499566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 105007218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 10519566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI)); 10529566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ)); 10539566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 105427f02ce8SMatthew G. Knepley if (dsAux) { 10559566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 10569566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 10579566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 10589566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 10599566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 10609566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 106101907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 10629566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 10639566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 106463a3b9bcSJacob 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); 106527f02ce8SMatthew G. Knepley } 10669566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI)); 10679566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ)); 1068665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1069665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 107027f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2 * NcI; 107127f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ; 10720abb75b6SMatthew G. Knepley if (!isCohesiveFieldI && s == 2) { 10730abb75b6SMatthew G. Knepley // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 10740abb75b6SMatthew G. Knepley NcS *= 2; 10750abb75b6SMatthew G. Knepley } 10760abb75b6SMatthew G. Knepley if (!isCohesiveFieldJ && s == 2) { 10770abb75b6SMatthew G. Knepley // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 10780abb75b6SMatthew G. Knepley NcT *= 2; 10790abb75b6SMatthew G. Knepley } 10800502970dSMatthew G. Knepley // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though 10810502970dSMatthew G. Knepley // the coordinates are in dE dimensions 10829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 10830502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 10840502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 10850502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 10869566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 10870e18dc48SMatthew G. Knepley PetscCall(PetscQuadratureGetCellType(quad, &ct)); 108863a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 108927f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 109027f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 10910e18dc48SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 10920e18dc48SMatthew G. Knepley const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 10934e913f38SMatthew G. Knepley const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 109427f02ce8SMatthew G. Knepley 109507218a29SMatthew G. Knepley fegeom.v = x; /* Workspace */ 109627f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 10970e18dc48SMatthew G. Knepley PetscInt qpt[2]; 109827f02ce8SMatthew G. Knepley PetscReal w; 109927f02ce8SMatthew G. Knepley PetscInt c; 110027f02ce8SMatthew G. Knepley 11014e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0])); 11024e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1])); 110307218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 110427f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 110507218a29SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 110663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 110727f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 11089566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 110927f02ce8SMatthew G. Knepley #endif 111027f02ce8SMatthew G. Knepley } 111163a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 11128e3a54c0SPierre 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)); 111307218a29SMatthew 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)); 1114ea672e62SMatthew G. Knepley if (n0) { 11159566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 1116148442b3SMatthew 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); 111727f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT; ++c) g0[c] *= w; 111827f02ce8SMatthew G. Knepley } 1119ea672e62SMatthew G. Knepley if (n1) { 11200502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 1121148442b3SMatthew 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); 11220502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w; 112327f02ce8SMatthew G. Knepley } 1124ea672e62SMatthew G. Knepley if (n2) { 11250502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 1126148442b3SMatthew 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); 11270502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w; 112827f02ce8SMatthew G. Knepley } 1129ea672e62SMatthew G. Knepley if (n3) { 11300502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 1131148442b3SMatthew 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); 11320502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w; 113327f02ce8SMatthew G. Knepley } 113427f02ce8SMatthew G. Knepley 11355fedec97SMatthew G. Knepley if (isCohesiveFieldI) { 11365fedec97SMatthew G. Knepley if (isCohesiveFieldJ) { 11379566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 113827f02ce8SMatthew G. Knepley } else { 11390abb75b6SMatthew 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)); 11400abb75b6SMatthew 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)); 11410abb75b6SMatthew G. Knepley } 11420abb75b6SMatthew G. Knepley } else { 11430abb75b6SMatthew G. Knepley if (s == 2) { 11440abb75b6SMatthew G. Knepley if (isCohesiveFieldJ) { 11450abb75b6SMatthew 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)); 11460abb75b6SMatthew 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)); 11470abb75b6SMatthew G. Knepley } else { 11480abb75b6SMatthew 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)); 11490abb75b6SMatthew 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)); 11500abb75b6SMatthew 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)); 11510abb75b6SMatthew 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)); 11525fedec97SMatthew G. Knepley } 11539371c9d4SSatish Balay } else 11540abb75b6SMatthew 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)); 11550abb75b6SMatthew G. Knepley } 115627f02ce8SMatthew G. Knepley } 115727f02ce8SMatthew G. Knepley if (debug > 1) { 11584e913f38SMatthew G. Knepley const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb)); 11594e913f38SMatthew G. Knepley const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb)); 11604e913f38SMatthew G. Knepley const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb)); 11614e913f38SMatthew G. Knepley const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb)); 11624e913f38SMatthew G. Knepley PetscInt f, g; 116327f02ce8SMatthew G. Knepley 11644e913f38SMatthew 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)); 11654e913f38SMatthew G. Knepley for (f = fS; f < fE; ++f) { 11664e913f38SMatthew G. Knepley const PetscInt i = offsetI + f; 11674e913f38SMatthew G. Knepley for (g = gS; g < gE; ++g) { 11684e913f38SMatthew G. Knepley const PetscInt j = offsetJ + g; 1169e978a55eSPierre 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); 11704e913f38SMatthew 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]))); 117127f02ce8SMatthew G. Knepley } 11729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 117327f02ce8SMatthew G. Knepley } 117427f02ce8SMatthew G. Knepley } 117527f02ce8SMatthew G. Knepley cOffset += totDim; 117627f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 117727f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 117827f02ce8SMatthew G. Knepley } 11793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118027f02ce8SMatthew G. Knepley } 118127f02ce8SMatthew G. Knepley 1182d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 1183d71ae5a4SJacob Faibussowitsch { 118420cf1dd8SToby Isaac PetscFunctionBegin; 118520cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 118620cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 118720cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 118820cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 118920cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1190ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 119120cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1192afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 119320cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 119420cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 119527f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 119620cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */; 119720cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 119820cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 119927f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 12003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120120cf1dd8SToby Isaac } 120220cf1dd8SToby Isaac 120320cf1dd8SToby Isaac /*MC 1204dce8aebaSBarry Smith PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization 120520cf1dd8SToby Isaac 120620cf1dd8SToby Isaac Level: intermediate 120720cf1dd8SToby Isaac 1208dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 120920cf1dd8SToby Isaac M*/ 121020cf1dd8SToby Isaac 1211d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 1212d71ae5a4SJacob Faibussowitsch { 121320cf1dd8SToby Isaac PetscFE_Basic *b; 121420cf1dd8SToby Isaac 121520cf1dd8SToby Isaac PetscFunctionBegin; 121620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 12174dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 121820cf1dd8SToby Isaac fem->data = b; 121920cf1dd8SToby Isaac 12209566063dSJacob Faibussowitsch PetscCall(PetscFEInitialize_Basic(fem)); 12213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122220cf1dd8SToby Isaac } 1223