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 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 163d71ae5a4SJacob Faibussowitsch { 16420cf1dd8SToby Isaac const PetscInt debug = 0; 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) { 2194bee2e38SMatthew G. Knepley PetscScalar integrand; 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 251d71ae5a4SJacob Faibussowitsch static 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 { 253afe6d6adSToby Isaac const PetscInt debug = 0; 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; 302ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 30327f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 30427f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 30527f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 30627f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3074bee2e38SMatthew G. Knepley if (isAffine) { 3084bee2e38SMatthew G. Knepley fegeom.v = x; 3094bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3107132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 3117132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 3127132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 3137132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 3149f209ee4SMatthew G. Knepley 3157132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 3167132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 3177132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 3184bee2e38SMatthew G. Knepley } 319afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 320afe6d6adSToby Isaac PetscScalar integrand; 3214bee2e38SMatthew G. Knepley PetscReal w; 322afe6d6adSToby Isaac 323afe6d6adSToby Isaac if (isAffine) { 3247132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 325afe6d6adSToby Isaac } else { 3263fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 3279f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 3289f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 3294bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 3304bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 3319f209ee4SMatthew G. Knepley 3329f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 3339f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 3349f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 335afe6d6adSToby Isaac } 3364bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 337afe6d6adSToby Isaac if (debug > 1 && q < Np) { 33863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 339afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3409566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 341afe6d6adSToby Isaac #endif 342afe6d6adSToby Isaac } 34363a3b9bcSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 3449566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 3459566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 3464bee2e38SMatthew 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); 3474bee2e38SMatthew G. Knepley integrand *= w; 348afe6d6adSToby Isaac integral[e * Nf + field] += integrand; 3499566063dSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field]))); 350afe6d6adSToby Isaac } 351afe6d6adSToby Isaac cOffset += totDim; 352afe6d6adSToby Isaac cOffsetAux += totDimAux; 353afe6d6adSToby Isaac } 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355afe6d6adSToby Isaac } 356afe6d6adSToby Isaac 357d71ae5a4SJacob 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[]) 358d71ae5a4SJacob Faibussowitsch { 35920cf1dd8SToby Isaac const PetscInt debug = 0; 3606528b96dSMatthew G. Knepley const PetscInt field = key.field; 3614bee2e38SMatthew G. Knepley PetscFE fe; 3626528b96dSMatthew G. Knepley PetscWeakForm wf; 3636528b96dSMatthew G. Knepley PetscInt n0, n1, i; 3646528b96dSMatthew G. Knepley PetscPointFunc *f0_func, *f1_func; 36520cf1dd8SToby Isaac PetscQuadrature quad; 366ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 3674bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 36820cf1dd8SToby Isaac const PetscScalar *constants; 36920cf1dd8SToby Isaac PetscReal *x; 370ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 371ef0bb6c7SMatthew G. Knepley PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 37220cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3736587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 37420cf1dd8SToby Isaac 37520cf1dd8SToby Isaac PetscFunctionBegin; 3769566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 3779566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 3799566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 3809566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 3819566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 3829566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 3839566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 3849566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 3859566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 3863ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 3879566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 3889566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 3909566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 3919566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 3924bee2e38SMatthew G. Knepley if (dsAux) { 3939566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 3949566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 3959566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 3969566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 3979566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 3989566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 39963a3b9bcSJacob 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); 40020cf1dd8SToby Isaac } 4019566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 40263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 40320cf1dd8SToby Isaac dE = cgeom->dimEmbed; 40463a3b9bcSJacob Faibussowitsch PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim); 40520cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4064bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 40720cf1dd8SToby Isaac 4086587ee25SMatthew G. Knepley fegeom.v = x; /* workspace */ 4099566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc)); 4109566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE)); 41120cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4124bee2e38SMatthew G. Knepley PetscReal w; 4134bee2e38SMatthew G. Knepley PetscInt c, d; 41420cf1dd8SToby Isaac 4159566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom)); 4164bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 4176587ee25SMatthew G. Knepley if (debug > 1 && q < cgeom->numPoints) { 41863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 4197be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4209566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 42120cf1dd8SToby Isaac #endif 42220cf1dd8SToby Isaac } 4239566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 4249566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 4256528b96dSMatthew 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]); 426ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w; 4276528b96dSMatthew 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]); 4289371c9d4SSatish Balay for (c = 0; c < T[field]->Nc; ++c) 4299371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w; 430b8025e53SMatthew G. Knepley if (debug) { 431*e8e188d2SZach Atkins // LCOV_EXCL_START 432*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q])); 433*e8e188d2SZach Atkins for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c])); 434*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 435b8025e53SMatthew G. Knepley if (debug > 2) { 43663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " field %" PetscInt_FMT ":", field)); 43763a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c]))); 4389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 439*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " field der %" PetscInt_FMT ":", field)); 440*e8e188d2SZach Atkins for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c]))); 441*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 44263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " resid %" PetscInt_FMT ":", field)); 44363a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c]))); 4449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 445*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, " res der %" PetscInt_FMT ":", field)); 446*e8e188d2SZach Atkins for (c = 0; c < T[field]->Nc; ++c) { 447*e8e188d2SZach Atkins for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dim + d]))); 448b8025e53SMatthew G. Knepley } 449*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 450*e8e188d2SZach Atkins } 451*e8e188d2SZach Atkins // LCOV_EXCL_STOP 452b8025e53SMatthew G. Knepley } 45320cf1dd8SToby Isaac } 4549566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset])); 45520cf1dd8SToby Isaac cOffset += totDim; 45620cf1dd8SToby Isaac cOffsetAux += totDimAux; 45720cf1dd8SToby Isaac } 4583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45920cf1dd8SToby Isaac } 46020cf1dd8SToby Isaac 461d71ae5a4SJacob 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[]) 462d71ae5a4SJacob Faibussowitsch { 46320cf1dd8SToby Isaac const PetscInt debug = 0; 46406d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 4654bee2e38SMatthew G. Knepley PetscFE fe; 46606d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 46706d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 46820cf1dd8SToby Isaac PetscQuadrature quad; 469ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4704bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 47120cf1dd8SToby Isaac const PetscScalar *constants; 47220cf1dd8SToby Isaac PetscReal *x; 473ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 474ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 4756587ee25SMatthew G. Knepley PetscBool auxOnBd = PETSC_FALSE; 47620cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 4776587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 47820cf1dd8SToby Isaac 47920cf1dd8SToby Isaac PetscFunctionBegin; 4809566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 4819566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 4829566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 4839566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 4849566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 4859566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 4869566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 4879566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 4889566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 4893ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 4909566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 4919566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 4929566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 4939566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 4949566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 4954bee2e38SMatthew G. Knepley if (dsAux) { 4969566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 4979566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 4989566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 4999566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 5009566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 5019566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 5027be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 5039566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 5049566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 50563a3b9bcSJacob 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); 50620cf1dd8SToby Isaac } 507ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 5089566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 50963a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 51020cf1dd8SToby Isaac dE = fgeom->dimEmbed; 5116587ee25SMatthew G. Knepley /* TODO FIX THIS */ 5126587ee25SMatthew G. Knepley fgeom->dim = dim - 1; 51363a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 51420cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5159f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 51620cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5179f209ee4SMatthew G. Knepley 5186587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 5199566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcI)); 5209566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcI * dE)); 52120cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5224bee2e38SMatthew G. Knepley PetscReal w; 5234bee2e38SMatthew G. Knepley PetscInt c, d; 5244bee2e38SMatthew G. Knepley 5259566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom)); 5269566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom)); 5274bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 52862bd480fSMatthew G. Knepley if (debug > 1) { 5296587ee25SMatthew G. Knepley if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 53063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 5317be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5329566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 5339566063dSJacob Faibussowitsch PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n)); 53420cf1dd8SToby Isaac #endif 53520cf1dd8SToby Isaac } 53662bd480fSMatthew G. Knepley } 5379566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 5389566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 53906d8a0d3SMatthew 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]); 5404bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w; 54106d8a0d3SMatthew 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]); 5429371c9d4SSatish Balay for (c = 0; c < NcI; ++c) 5439371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w; 54462bd480fSMatthew G. Knepley if (debug) { 54563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q)); 54662bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 54763a3b9bcSJacob Faibussowitsch if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c]))); 54862bd480fSMatthew G. Knepley if (n1) { 54963a3b9bcSJacob 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]))); 5509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 55162bd480fSMatthew G. Knepley } 55262bd480fSMatthew G. Knepley } 55362bd480fSMatthew G. Knepley } 55420cf1dd8SToby Isaac } 5559566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 55620cf1dd8SToby Isaac cOffset += totDim; 55720cf1dd8SToby Isaac cOffsetAux += totDimAux; 55820cf1dd8SToby Isaac } 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56020cf1dd8SToby Isaac } 56120cf1dd8SToby Isaac 56227f02ce8SMatthew G. Knepley /* 56327f02ce8SMatthew 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 56427f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 56527f02ce8SMatthew G. Knepley 56627f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 56727f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 56827f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 56927f02ce8SMatthew 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() 57027f02ce8SMatthew G. Knepley */ 57107218a29SMatthew G. Knepley static 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[]) 572d71ae5a4SJacob Faibussowitsch { 57327f02ce8SMatthew G. Knepley const PetscInt debug = 0; 5746528b96dSMatthew G. Knepley const PetscInt field = key.field; 57527f02ce8SMatthew G. Knepley PetscFE fe; 5766528b96dSMatthew G. Knepley PetscWeakForm wf; 5776528b96dSMatthew G. Knepley PetscInt n0, n1, i; 5786528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 57927f02ce8SMatthew G. Knepley PetscQuadrature quad; 5800e18dc48SMatthew G. Knepley DMPolytopeType ct; 58107218a29SMatthew G. Knepley PetscTabulation *Tf, *TfIn, *TfAux = NULL; 58227f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 58327f02ce8SMatthew G. Knepley const PetscScalar *constants; 58427f02ce8SMatthew G. Knepley PetscReal *x; 585665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 58607218a29SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 5876587ee25SMatthew G. Knepley PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 58827f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 5896587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 59027f02ce8SMatthew G. Knepley 59127f02ce8SMatthew G. Knepley PetscFunctionBegin; 59227f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 5939566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 5949566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 5959566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 5969566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 5979566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 59807218a29SMatthew G. Knepley PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn)); 599429ebbe4SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets 60007218a29SMatthew G. Knepley PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x)); 6019566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset)); 6029566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 6039566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 6043ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 6059566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 6069566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 6079566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 60827f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 6099566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &Tf)); 61007218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 6119566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 61227f02ce8SMatthew G. Knepley if (dsAux) { 6139566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 6149566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 6159566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 6169566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 6179566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 6189566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 61901907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 6209566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 6219566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 62263a3b9bcSJacob 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); 62327f02ce8SMatthew G. Knepley } 6249566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField)); 625665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 626c2b7495fSMatthew G. Knepley NcS = NcI; 6279566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 6280e18dc48SMatthew G. Knepley PetscCall(PetscQuadratureGetCellType(quad, &ct)); 62963a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 63027f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 63163a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 63227f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 63327f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 6340e18dc48SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 6350e18dc48SMatthew G. Knepley const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 6364e913f38SMatthew G. Knepley const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 63727f02ce8SMatthew G. Knepley 6386587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 6399566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcS)); 6409566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcS * dE)); 64127f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 6420e18dc48SMatthew G. Knepley PetscInt qpt[2]; 64327f02ce8SMatthew G. Knepley PetscReal w; 64427f02ce8SMatthew G. Knepley PetscInt c, d; 64527f02ce8SMatthew G. Knepley 6464e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0])); 6474e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1])); 64807218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 64927f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 6506587ee25SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 65163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 65227f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 6539566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ)); 65427f02ce8SMatthew G. Knepley #endif 65527f02ce8SMatthew G. Knepley } 656a4158a15SMatthew 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])); 65727f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 65807218a29SMatthew G. Knepley PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], &coefficients_t[cOffsetIn], u, u_x, u_t)); 65907218a29SMatthew 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)); 6606528b96dSMatthew 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]); 66127f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w; 6629ee2af8cSMatthew 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]); 6639371c9d4SSatish Balay for (c = 0; c < NcS; ++c) 6649371c9d4SSatish Balay for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w; 66527f02ce8SMatthew G. Knepley } 6669371c9d4SSatish Balay if (isCohesiveField) { 6673ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 6689371c9d4SSatish Balay } else { 6693ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 6709371c9d4SSatish Balay } 67127f02ce8SMatthew G. Knepley cOffset += totDim; 67207218a29SMatthew G. Knepley cOffsetIn += totDimIn; 67327f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 67427f02ce8SMatthew G. Knepley } 6753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67627f02ce8SMatthew G. Knepley } 67727f02ce8SMatthew G. Knepley 678d71ae5a4SJacob 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[]) 679d71ae5a4SJacob Faibussowitsch { 68020cf1dd8SToby Isaac const PetscInt debug = 0; 6814bee2e38SMatthew G. Knepley PetscFE feI, feJ; 6826528b96dSMatthew G. Knepley PetscWeakForm wf; 6836528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 6846528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 68520cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 68620cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 68720cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 68820cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 68920cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 69020cf1dd8SToby Isaac PetscQuadrature quad; 691ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 6924bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 69320cf1dd8SToby Isaac const PetscScalar *constants; 69420cf1dd8SToby Isaac PetscReal *x; 695ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 696ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 6976528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 69820cf1dd8SToby Isaac PetscInt dE, Np; 69920cf1dd8SToby Isaac PetscBool isAffine; 70020cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 70120cf1dd8SToby Isaac PetscInt qNc, Nq, q; 70220cf1dd8SToby Isaac 70320cf1dd8SToby Isaac PetscFunctionBegin; 7049566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 7056528b96dSMatthew G. Knepley fieldI = key.field / Nf; 7066528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 7079566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 7089566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 7099566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 7109566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 7119566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 7129566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 7139566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 7149566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 71520cf1dd8SToby Isaac switch (jtype) { 716d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 717d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 718d71ae5a4SJacob Faibussowitsch break; 719d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 720d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 721d71ae5a4SJacob Faibussowitsch break; 722d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 723d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 724d71ae5a4SJacob Faibussowitsch break; 72520cf1dd8SToby Isaac } 7263ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 7279566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 7289566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 7299566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 7309566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 7319566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 7329566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 7339566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 7344bee2e38SMatthew G. Knepley if (dsAux) { 7359566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 7369566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 7379566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 7389566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 7399566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 7409566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 74163a3b9bcSJacob 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); 74220cf1dd8SToby Isaac } 74327f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 74427f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7454bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7464bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7474bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 74827f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 7499566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 7509566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 7519566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 7529566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 7539566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 75463a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 7554bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7564bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7574bee2e38SMatthew G. Knepley 75827f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 75927f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7604bee2e38SMatthew G. Knepley if (isAffine) { 7614bee2e38SMatthew G. Knepley fegeom.v = x; 7624bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 7637132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 7647132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 7657132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 7664bee2e38SMatthew G. Knepley } 76720cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 76820cf1dd8SToby Isaac PetscReal w; 7694bee2e38SMatthew G. Knepley PetscInt c; 77020cf1dd8SToby Isaac 77120cf1dd8SToby Isaac if (isAffine) { 7727132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 77320cf1dd8SToby Isaac } else { 7744bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 7754bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 7764bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 7774bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 77820cf1dd8SToby Isaac } 7799566063dSJacob 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])); 7804bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 7819566063dSJacob Faibussowitsch if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 7829566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 783ea672e62SMatthew G. Knepley if (n0) { 7849566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 7856528b96dSMatthew 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); 78620cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 78720cf1dd8SToby Isaac } 788ea672e62SMatthew G. Knepley if (n1) { 7899566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 7906528b96dSMatthew 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); 7914bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 79220cf1dd8SToby Isaac } 793ea672e62SMatthew G. Knepley if (n2) { 7949566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 7956528b96dSMatthew 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); 7964bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 79720cf1dd8SToby Isaac } 798ea672e62SMatthew G. Knepley if (n3) { 7999566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 8006528b96dSMatthew 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); 8014bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 80220cf1dd8SToby Isaac } 80320cf1dd8SToby Isaac 8049566063dSJacob 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)); 80520cf1dd8SToby Isaac } 80620cf1dd8SToby Isaac if (debug > 1) { 80720cf1dd8SToby Isaac PetscInt fc, f, gc, g; 80820cf1dd8SToby Isaac 80963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 810ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 811ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 812ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 813ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 814ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 815ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 81663a3b9bcSJacob 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]))); 81720cf1dd8SToby Isaac } 81820cf1dd8SToby Isaac } 8199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 82020cf1dd8SToby Isaac } 82120cf1dd8SToby Isaac } 82220cf1dd8SToby Isaac } 82320cf1dd8SToby Isaac cOffset += totDim; 82420cf1dd8SToby Isaac cOffsetAux += totDimAux; 82520cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 82620cf1dd8SToby Isaac } 8273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82820cf1dd8SToby Isaac } 82920cf1dd8SToby Isaac 830d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEIntegrateBdJacobian_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, PetscReal u_tshift, PetscScalar elemMat[]) 831d71ae5a4SJacob Faibussowitsch { 83220cf1dd8SToby Isaac const PetscInt debug = 0; 8334bee2e38SMatthew G. Knepley PetscFE feI, feJ; 83445480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 83545480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 83620cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 83720cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 83820cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 83920cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 84020cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 84120cf1dd8SToby Isaac PetscQuadrature quad; 842ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8434bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 84420cf1dd8SToby Isaac const PetscScalar *constants; 84520cf1dd8SToby Isaac PetscReal *x; 846ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 847ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 84845480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 84920cf1dd8SToby Isaac PetscBool isAffine; 85020cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 85120cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 85220cf1dd8SToby Isaac 85320cf1dd8SToby Isaac PetscFunctionBegin; 8549566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 85545480ffeSMatthew G. Knepley fieldI = key.field / Nf; 85645480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8579566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 8589566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 8599566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 8609566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(feI, &quad)); 8619566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 8629566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 8639566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 8649566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 8659566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 8669566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 8673ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 8689566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 8699566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 8709566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 8719566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &T)); 8729566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 8734bee2e38SMatthew G. Knepley if (dsAux) { 8749566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 8759566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 8769566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 8779566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 8789566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 8799566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 88020cf1dd8SToby Isaac } 881ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 88220cf1dd8SToby Isaac Np = fgeom->numPoints; 88320cf1dd8SToby Isaac dE = fgeom->dimEmbed; 88420cf1dd8SToby Isaac isAffine = fgeom->isAffine; 88527f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 8869566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 8879566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 8889566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 8899566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 8909566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 89163a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 89220cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 8939f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 89420cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 895ea78f98cSLisandro Dalcin fegeom.n = NULL; 896ea78f98cSLisandro Dalcin fegeom.v = NULL; 897ea78f98cSLisandro Dalcin fegeom.J = NULL; 898ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 89927f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 90027f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 90127f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 90227f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 9034bee2e38SMatthew G. Knepley if (isAffine) { 9044bee2e38SMatthew G. Knepley fegeom.v = x; 9054bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 9067132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 9077132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 9087132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 9097132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 9109f209ee4SMatthew G. Knepley 9117132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 9127132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 9137132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 9144bee2e38SMatthew G. Knepley } 91520cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 91620cf1dd8SToby Isaac PetscReal w; 9174bee2e38SMatthew G. Knepley PetscInt c; 91820cf1dd8SToby Isaac 91963a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 92020cf1dd8SToby Isaac if (isAffine) { 9217132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 92220cf1dd8SToby Isaac } else { 9233fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 9249f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 9259f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 9264bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 9274bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 9289f209ee4SMatthew G. Knepley 9299f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 9309f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 9319f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 93220cf1dd8SToby Isaac } 9334bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 9349566063dSJacob Faibussowitsch if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 9359566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 936ea672e62SMatthew G. Knepley if (n0) { 9379566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 93845480ffeSMatthew 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); 93920cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 94020cf1dd8SToby Isaac } 941ea672e62SMatthew G. Knepley if (n1) { 9429566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 94345480ffeSMatthew 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); 9444bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 94520cf1dd8SToby Isaac } 946ea672e62SMatthew G. Knepley if (n2) { 9479566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 94845480ffeSMatthew 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); 9494bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 95020cf1dd8SToby Isaac } 951ea672e62SMatthew G. Knepley if (n3) { 9529566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 95345480ffeSMatthew 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); 9544bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 95520cf1dd8SToby Isaac } 95620cf1dd8SToby Isaac 9579566063dSJacob 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)); 95820cf1dd8SToby Isaac } 95920cf1dd8SToby Isaac if (debug > 1) { 96020cf1dd8SToby Isaac PetscInt fc, f, gc, g; 96120cf1dd8SToby Isaac 96263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 963ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 964ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 965ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 966ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 967ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 968ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 96963a3b9bcSJacob 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]))); 97020cf1dd8SToby Isaac } 97120cf1dd8SToby Isaac } 9729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 97320cf1dd8SToby Isaac } 97420cf1dd8SToby Isaac } 97520cf1dd8SToby Isaac } 97620cf1dd8SToby Isaac cOffset += totDim; 97720cf1dd8SToby Isaac cOffsetAux += totDimAux; 97820cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 97920cf1dd8SToby Isaac } 9803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98120cf1dd8SToby Isaac } 98220cf1dd8SToby Isaac 98366976f2fSJacob Faibussowitsch static 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[]) 984d71ae5a4SJacob Faibussowitsch { 98527f02ce8SMatthew G. Knepley const PetscInt debug = 0; 98627f02ce8SMatthew G. Knepley PetscFE feI, feJ; 987148442b3SMatthew G. Knepley PetscWeakForm wf; 988148442b3SMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 989148442b3SMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 99027f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 99127f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 99227f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 99327f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 99427f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 995665f567fSMatthew G. Knepley PetscQuadrature quad; 9960e18dc48SMatthew G. Knepley DMPolytopeType ct; 99707218a29SMatthew G. Knepley PetscTabulation *T, *TfIn, *TAux = NULL; 99827f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 99927f02ce8SMatthew G. Knepley const PetscScalar *constants; 100027f02ce8SMatthew G. Knepley PetscReal *x; 1001665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1002665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 100345480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 100407218a29SMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE; 100527f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 10060502970dSMatthew G. Knepley PetscInt qNc, Nq, q; 100727f02ce8SMatthew G. Knepley 100827f02ce8SMatthew G. Knepley PetscFunctionBegin; 10099566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 101045480ffeSMatthew G. Knepley fieldI = key.field / Nf; 101145480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 101227f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 10139566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 10149566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 10159566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 10169566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 10179566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1018429ebbe4SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets 10199566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x)); 10209566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 102127f02ce8SMatthew G. Knepley switch (jtype) { 1022d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 1023d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1024d71ae5a4SJacob Faibussowitsch break; 1025d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 1026d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1027d71ae5a4SJacob Faibussowitsch break; 1028d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 1029d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 103027f02ce8SMatthew G. Knepley } 10313ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 10329566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 10339566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 10349566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 10359566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 103607218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 10379566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI)); 10389566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ)); 10399566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 104027f02ce8SMatthew G. Knepley if (dsAux) { 10419566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 10429566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 10439566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 10449566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 10459566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 10469566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 104701907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 10489566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 10499566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 105063a3b9bcSJacob 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); 105127f02ce8SMatthew G. Knepley } 10529566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI)); 10539566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ)); 1054665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1055665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 105627f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2 * NcI; 105727f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ; 10580502970dSMatthew G. Knepley // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though 10590502970dSMatthew G. Knepley // the coordinates are in dE dimensions 10609566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 10610502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 10620502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 10630502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 10649566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 10650e18dc48SMatthew G. Knepley PetscCall(PetscQuadratureGetCellType(quad, &ct)); 106663a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 106727f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 106827f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 10690e18dc48SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 10700e18dc48SMatthew G. Knepley const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 10714e913f38SMatthew G. Knepley const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 107227f02ce8SMatthew G. Knepley 107307218a29SMatthew G. Knepley fegeom.v = x; /* Workspace */ 107427f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 10750e18dc48SMatthew G. Knepley PetscInt qpt[2]; 107627f02ce8SMatthew G. Knepley PetscReal w; 107727f02ce8SMatthew G. Knepley PetscInt c; 107827f02ce8SMatthew G. Knepley 10794e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0])); 10804e913f38SMatthew G. Knepley PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1])); 108107218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 108227f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 108307218a29SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 108463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 108527f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 10869566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 108727f02ce8SMatthew G. Knepley #endif 108827f02ce8SMatthew G. Knepley } 108963a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 109007218a29SMatthew G. Knepley if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 109107218a29SMatthew 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)); 1092ea672e62SMatthew G. Knepley if (n0) { 10939566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 1094148442b3SMatthew 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); 109527f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT; ++c) g0[c] *= w; 109627f02ce8SMatthew G. Knepley } 1097ea672e62SMatthew G. Knepley if (n1) { 10980502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 1099148442b3SMatthew 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); 11000502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w; 110127f02ce8SMatthew G. Knepley } 1102ea672e62SMatthew G. Knepley if (n2) { 11030502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 1104148442b3SMatthew 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); 11050502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w; 110627f02ce8SMatthew G. Knepley } 1107ea672e62SMatthew G. Knepley if (n3) { 11080502970dSMatthew G. Knepley PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 1109148442b3SMatthew 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); 11100502970dSMatthew G. Knepley for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w; 111127f02ce8SMatthew G. Knepley } 111227f02ce8SMatthew G. Knepley 11135fedec97SMatthew G. Knepley if (isCohesiveFieldI) { 11145fedec97SMatthew G. Knepley if (isCohesiveFieldJ) { 11159566063dSJacob 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)); 111627f02ce8SMatthew G. Knepley } else { 11179566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 11189566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 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)); 11195fedec97SMatthew G. Knepley } 11209371c9d4SSatish Balay } else 11219371c9d4SSatish Balay PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 112227f02ce8SMatthew G. Knepley } 112327f02ce8SMatthew G. Knepley if (debug > 1) { 11244e913f38SMatthew G. Knepley const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb)); 11254e913f38SMatthew G. Knepley const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb)); 11264e913f38SMatthew G. Knepley const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb)); 11274e913f38SMatthew G. Knepley const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb)); 11284e913f38SMatthew G. Knepley PetscInt f, g; 112927f02ce8SMatthew G. Knepley 11304e913f38SMatthew 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)); 11314e913f38SMatthew G. Knepley for (f = fS; f < fE; ++f) { 11324e913f38SMatthew G. Knepley const PetscInt i = offsetI + f; 11334e913f38SMatthew G. Knepley for (g = gS; g < gE; ++g) { 11344e913f38SMatthew G. Knepley const PetscInt j = offsetJ + g; 11354e913f38SMatthew G. Knepley PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", f, i, g, j); 11364e913f38SMatthew 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]))); 113727f02ce8SMatthew G. Knepley } 11389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 113927f02ce8SMatthew G. Knepley } 114027f02ce8SMatthew G. Knepley } 114127f02ce8SMatthew G. Knepley cOffset += totDim; 114227f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 114327f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 114427f02ce8SMatthew G. Knepley } 11453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114627f02ce8SMatthew G. Knepley } 114727f02ce8SMatthew G. Knepley 1148d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 1149d71ae5a4SJacob Faibussowitsch { 115020cf1dd8SToby Isaac PetscFunctionBegin; 115120cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 115220cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 115320cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 115420cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 115520cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1156ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 115720cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1158afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 115920cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 116020cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 116127f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 116220cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */; 116320cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 116420cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 116527f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 11663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116720cf1dd8SToby Isaac } 116820cf1dd8SToby Isaac 116920cf1dd8SToby Isaac /*MC 1170dce8aebaSBarry Smith PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization 117120cf1dd8SToby Isaac 117220cf1dd8SToby Isaac Level: intermediate 117320cf1dd8SToby Isaac 1174dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 117520cf1dd8SToby Isaac M*/ 117620cf1dd8SToby Isaac 1177d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 1178d71ae5a4SJacob Faibussowitsch { 117920cf1dd8SToby Isaac PetscFE_Basic *b; 118020cf1dd8SToby Isaac 118120cf1dd8SToby Isaac PetscFunctionBegin; 118220cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 11834dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 118420cf1dd8SToby Isaac fem->data = b; 118520cf1dd8SToby Isaac 11869566063dSJacob Faibussowitsch PetscCall(PetscFEInitialize_Basic(fem)); 11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118820cf1dd8SToby Isaac } 1189