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; 103b9d4cb8dSJed Brown for (i = 0; i < m; i++) { 104b9d4cb8dSJed Brown PetscBLASInt n_, p_, k_, lda, ldb, ldc; 105b9d4cb8dSJed Brown PetscReal one = 1, zero = 0; 106b9d4cb8dSJed Brown /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n] 107b9d4cb8dSJed Brown * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k) 108b9d4cb8dSJed Brown */ 1099566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(n, &n_)); 1109566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(p, &p_)); 1119566063dSJacob Faibussowitsch PetscCall(PetscBLASIntCast(k, &k_)); 112b9d4cb8dSJed Brown lda = p_; 113b9d4cb8dSJed Brown ldb = n_; 114b9d4cb8dSJed Brown ldc = p_; 115792fecdfSBarry Smith PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc)); 116b9d4cb8dSJed Brown } 1179566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2. * m * n * p * k)); 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119b9d4cb8dSJed Brown } 120b9d4cb8dSJed Brown 121d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PetscFECreateTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 122d71ae5a4SJacob Faibussowitsch { 12320cf1dd8SToby Isaac DM dm; 12420cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 12520cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 12620cf1dd8SToby Isaac PetscInt Nc; /* Field components */ 127ef0bb6c7SMatthew G. Knepley PetscReal *B = K >= 0 ? T->T[0] : NULL; 128ef0bb6c7SMatthew G. Knepley PetscReal *D = K >= 1 ? T->T[1] : NULL; 129ef0bb6c7SMatthew G. Knepley PetscReal *H = K >= 2 ? T->T[2] : NULL; 130ef0bb6c7SMatthew G. Knepley PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL; 13120cf1dd8SToby Isaac 13220cf1dd8SToby Isaac PetscFunctionBegin; 1339566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 1349566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1359566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 1369566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fem, &Nc)); 13720cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 1389566063dSJacob Faibussowitsch if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB)); 1399566063dSJacob Faibussowitsch if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD)); 1409566063dSJacob Faibussowitsch if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH)); 142b9d4cb8dSJed Brown /* Translate from prime to nodal basis */ 14320cf1dd8SToby Isaac if (B) { 144b9d4cb8dSJed Brown /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */ 1459566063dSJacob Faibussowitsch PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B)); 14620cf1dd8SToby Isaac } 14720cf1dd8SToby Isaac if (D) { 148b9d4cb8dSJed Brown /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */ 1499566063dSJacob Faibussowitsch PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D)); 15020cf1dd8SToby Isaac } 15120cf1dd8SToby Isaac if (H) { 152b9d4cb8dSJed Brown /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */ 1539566063dSJacob Faibussowitsch PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H)); 15420cf1dd8SToby Isaac } 1559566063dSJacob Faibussowitsch if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB)); 1569566063dSJacob Faibussowitsch if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD)); 1579566063dSJacob Faibussowitsch if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH)); 1583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15920cf1dd8SToby Isaac } 16020cf1dd8SToby Isaac 161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 162d71ae5a4SJacob Faibussowitsch { 16320cf1dd8SToby Isaac const PetscInt debug = 0; 1644bee2e38SMatthew G. Knepley PetscFE fe; 16520cf1dd8SToby Isaac PetscPointFunc obj_func; 16620cf1dd8SToby Isaac PetscQuadrature quad; 167ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 1684bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x; 16920cf1dd8SToby Isaac const PetscScalar *constants; 17020cf1dd8SToby Isaac PetscReal *x; 171ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 17220cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 17320cf1dd8SToby Isaac PetscBool isAffine; 17420cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 17520cf1dd8SToby Isaac PetscInt qNc, Nq, q; 17620cf1dd8SToby Isaac 17720cf1dd8SToby Isaac PetscFunctionBegin; 1789566063dSJacob Faibussowitsch PetscCall(PetscDSGetObjective(ds, field, &obj_func)); 1793ba16761SJacob Faibussowitsch if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 1809566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 1819566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 1829566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 1839566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1849566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1859566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 1869566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 1879566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 1889566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 1899566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL)); 1909566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 1914bee2e38SMatthew G. Knepley if (dsAux) { 1929566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 1939566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 1949566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 1959566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 1969566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 19863a3b9bcSJacob 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); 19920cf1dd8SToby Isaac } 2009566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 20163a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 20220cf1dd8SToby Isaac Np = cgeom->numPoints; 20320cf1dd8SToby Isaac dE = cgeom->dimEmbed; 20420cf1dd8SToby Isaac isAffine = cgeom->isAffine; 20520cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 2064bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 20720cf1dd8SToby Isaac 20827f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 20927f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 21020cf1dd8SToby Isaac if (isAffine) { 2114bee2e38SMatthew G. Knepley fegeom.v = x; 2124bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 2137132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 2147132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 2157132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 21620cf1dd8SToby Isaac } 2174bee2e38SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 2184bee2e38SMatthew G. Knepley PetscScalar integrand; 2194bee2e38SMatthew G. Knepley PetscReal w; 2204bee2e38SMatthew G. Knepley 2214bee2e38SMatthew G. Knepley if (isAffine) { 2227132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 2234bee2e38SMatthew G. Knepley } else { 2244bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 2254bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 2264bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 2274bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 2284bee2e38SMatthew G. Knepley } 2294bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 23020cf1dd8SToby Isaac if (debug > 1 && q < Np) { 23163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 2327be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 2339566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 23420cf1dd8SToby Isaac #endif 23520cf1dd8SToby Isaac } 23663a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 2379566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 2389566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 2394bee2e38SMatthew 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); 2404bee2e38SMatthew G. Knepley integrand *= w; 24120cf1dd8SToby Isaac integral[e * Nf + field] += integrand; 2429566063dSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[field]))); 24320cf1dd8SToby Isaac } 24420cf1dd8SToby Isaac cOffset += totDim; 24520cf1dd8SToby Isaac cOffsetAux += totDimAux; 24620cf1dd8SToby Isaac } 2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24820cf1dd8SToby Isaac } 24920cf1dd8SToby Isaac 250d71ae5a4SJacob 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[]) 251d71ae5a4SJacob Faibussowitsch { 252afe6d6adSToby Isaac const PetscInt debug = 0; 2534bee2e38SMatthew G. Knepley PetscFE fe; 254afe6d6adSToby Isaac PetscQuadrature quad; 255ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 2564bee2e38SMatthew G. Knepley PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 257afe6d6adSToby Isaac const PetscScalar *constants; 258afe6d6adSToby Isaac PetscReal *x; 259ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 260afe6d6adSToby Isaac PetscBool isAffine, auxOnBd; 261afe6d6adSToby Isaac const PetscReal *quadPoints, *quadWeights; 262afe6d6adSToby Isaac PetscInt qNc, Nq, q, Np, dE; 263afe6d6adSToby Isaac PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 264afe6d6adSToby Isaac 265afe6d6adSToby Isaac PetscFunctionBegin; 2663ba16761SJacob Faibussowitsch if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 2679566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 2689566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 2709566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 2719566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 2729566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 2739566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 2749566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 2759566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 2769566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 2779566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 2784bee2e38SMatthew G. Knepley if (dsAux) { 2799566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 2809566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 2819566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 2829566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 2839566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 285afe6d6adSToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 2869566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 2879566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 28863a3b9bcSJacob 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); 289afe6d6adSToby Isaac } 2909566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 29163a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 292afe6d6adSToby Isaac Np = fgeom->numPoints; 293afe6d6adSToby Isaac dE = fgeom->dimEmbed; 294afe6d6adSToby Isaac isAffine = fgeom->isAffine; 295afe6d6adSToby Isaac for (e = 0; e < Ne; ++e) { 2969f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 297afe6d6adSToby Isaac const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 298ea78f98cSLisandro Dalcin fegeom.n = NULL; 299ea78f98cSLisandro Dalcin fegeom.v = NULL; 300ea78f98cSLisandro Dalcin fegeom.J = NULL; 301ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 30227f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 30327f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 30427f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 30527f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 3064bee2e38SMatthew G. Knepley if (isAffine) { 3074bee2e38SMatthew G. Knepley fegeom.v = x; 3084bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 3097132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 3107132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 3117132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 3127132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 3139f209ee4SMatthew G. Knepley 3147132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 3157132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 3167132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 3174bee2e38SMatthew G. Knepley } 318afe6d6adSToby Isaac for (q = 0; q < Nq; ++q) { 319afe6d6adSToby Isaac PetscScalar integrand; 3204bee2e38SMatthew G. Knepley PetscReal w; 321afe6d6adSToby Isaac 322afe6d6adSToby Isaac if (isAffine) { 3237132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 324afe6d6adSToby Isaac } else { 3253fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 3269f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 3279f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 3284bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 3294bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 3309f209ee4SMatthew G. Knepley 3319f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 3329f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 3339f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 334afe6d6adSToby Isaac } 3354bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 336afe6d6adSToby Isaac if (debug > 1 && q < Np) { 33763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 338afe6d6adSToby Isaac #ifndef PETSC_USE_COMPLEX 3399566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 340afe6d6adSToby Isaac #endif 341afe6d6adSToby Isaac } 34263a3b9bcSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 3439566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 3449566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 3454bee2e38SMatthew 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); 3464bee2e38SMatthew G. Knepley integrand *= w; 347afe6d6adSToby Isaac integral[e * Nf + field] += integrand; 3489566063dSJacob Faibussowitsch if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field]))); 349afe6d6adSToby Isaac } 350afe6d6adSToby Isaac cOffset += totDim; 351afe6d6adSToby Isaac cOffsetAux += totDimAux; 352afe6d6adSToby Isaac } 3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 354afe6d6adSToby Isaac } 355afe6d6adSToby Isaac 356d71ae5a4SJacob 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[]) 357d71ae5a4SJacob Faibussowitsch { 35820cf1dd8SToby Isaac const PetscInt debug = 0; 3596528b96dSMatthew G. Knepley const PetscInt field = key.field; 3604bee2e38SMatthew G. Knepley PetscFE fe; 3616528b96dSMatthew G. Knepley PetscWeakForm wf; 3626528b96dSMatthew G. Knepley PetscInt n0, n1, i; 3636528b96dSMatthew G. Knepley PetscPointFunc *f0_func, *f1_func; 36420cf1dd8SToby Isaac PetscQuadrature quad; 365ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 3664bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 36720cf1dd8SToby Isaac const PetscScalar *constants; 36820cf1dd8SToby Isaac PetscReal *x; 369ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 370ef0bb6c7SMatthew G. Knepley PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 37120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 3726587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 37320cf1dd8SToby Isaac 37420cf1dd8SToby Isaac PetscFunctionBegin; 3759566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 3769566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 3779566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 3789566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 3799566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 3809566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 3819566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 3829566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 3839566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 3849566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 3853ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 3869566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 3879566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 3889566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 3909566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 3914bee2e38SMatthew G. Knepley if (dsAux) { 3929566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 3939566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 3949566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 3959566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 3969566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 3979566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 39863a3b9bcSJacob 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); 39920cf1dd8SToby Isaac } 4009566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 40163a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 40220cf1dd8SToby Isaac dE = cgeom->dimEmbed; 40363a3b9bcSJacob Faibussowitsch PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim); 40420cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 4054bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 40620cf1dd8SToby Isaac 4076587ee25SMatthew G. Knepley fegeom.v = x; /* workspace */ 4089566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc)); 4099566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE)); 41020cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 4114bee2e38SMatthew G. Knepley PetscReal w; 4124bee2e38SMatthew G. Knepley PetscInt c, d; 41320cf1dd8SToby Isaac 4149566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom)); 4154bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 4166587ee25SMatthew G. Knepley if (debug > 1 && q < cgeom->numPoints) { 41763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 4187be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 4199566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 42020cf1dd8SToby Isaac #endif 42120cf1dd8SToby Isaac } 4229566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 4239566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 4246528b96dSMatthew 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]); 425ef0bb6c7SMatthew G. Knepley for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w; 4266528b96dSMatthew 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]); 4279371c9d4SSatish Balay for (c = 0; c < T[field]->Nc; ++c) 4289371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[(q * T[field]->Nc + c) * dim + d] *= w; 429b8025e53SMatthew G. Knepley if (debug) { 43063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " wt %g\n", q, (double)quadWeights[q])); 431b8025e53SMatthew G. Knepley if (debug > 2) { 43263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " field %" PetscInt_FMT ":", field)); 43363a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c]))); 4349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 43563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " resid %" PetscInt_FMT ":", field)); 43663a3b9bcSJacob Faibussowitsch for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c]))); 4379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 438b8025e53SMatthew G. Knepley } 439b8025e53SMatthew G. Knepley } 44020cf1dd8SToby Isaac } 4419566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset])); 44220cf1dd8SToby Isaac cOffset += totDim; 44320cf1dd8SToby Isaac cOffsetAux += totDimAux; 44420cf1dd8SToby Isaac } 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44620cf1dd8SToby Isaac } 44720cf1dd8SToby Isaac 448d71ae5a4SJacob 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[]) 449d71ae5a4SJacob Faibussowitsch { 45020cf1dd8SToby Isaac const PetscInt debug = 0; 45106d8a0d3SMatthew G. Knepley const PetscInt field = key.field; 4524bee2e38SMatthew G. Knepley PetscFE fe; 45306d8a0d3SMatthew G. Knepley PetscInt n0, n1, i; 45406d8a0d3SMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 45520cf1dd8SToby Isaac PetscQuadrature quad; 456ef0bb6c7SMatthew G. Knepley PetscTabulation *Tf, *TfAux = NULL; 4574bee2e38SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 45820cf1dd8SToby Isaac const PetscScalar *constants; 45920cf1dd8SToby Isaac PetscReal *x; 460ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 461ef0bb6c7SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 4626587ee25SMatthew G. Knepley PetscBool auxOnBd = PETSC_FALSE; 46320cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 4646587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 46520cf1dd8SToby Isaac 46620cf1dd8SToby Isaac PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 4689566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 4699566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 4709566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 4719566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 4729566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 4739566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 4749566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 4759566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 4763ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 4779566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 4789566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 4799566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 4809566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 4819566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 4824bee2e38SMatthew G. Knepley if (dsAux) { 4839566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 4849566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 4859566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 4869566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 4879566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 4889566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 4897be5e748SToby Isaac auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 4909566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 4919566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 49263a3b9bcSJacob 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); 49320cf1dd8SToby Isaac } 494ef0bb6c7SMatthew G. Knepley NcI = Tf[field]->Nc; 4959566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 49663a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 49720cf1dd8SToby Isaac dE = fgeom->dimEmbed; 4986587ee25SMatthew G. Knepley /* TODO FIX THIS */ 4996587ee25SMatthew G. Knepley fgeom->dim = dim - 1; 50063a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 50120cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 5029f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 50320cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 5049f209ee4SMatthew G. Knepley 5056587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 5069566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcI)); 5079566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcI * dE)); 50820cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 5094bee2e38SMatthew G. Knepley PetscReal w; 5104bee2e38SMatthew G. Knepley PetscInt c, d; 5114bee2e38SMatthew G. Knepley 5129566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom)); 5139566063dSJacob Faibussowitsch PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom)); 5144bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 51562bd480fSMatthew G. Knepley if (debug > 1) { 5166587ee25SMatthew G. Knepley if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 51763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 5187be5e748SToby Isaac #if !defined(PETSC_USE_COMPLEX) 5199566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 5209566063dSJacob Faibussowitsch PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n)); 52120cf1dd8SToby Isaac #endif 52220cf1dd8SToby Isaac } 52362bd480fSMatthew G. Knepley } 5249566063dSJacob Faibussowitsch PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 5259566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 52606d8a0d3SMatthew 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]); 5274bee2e38SMatthew G. Knepley for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w; 52806d8a0d3SMatthew 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]); 5299371c9d4SSatish Balay for (c = 0; c < NcI; ++c) 5309371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[(q * NcI + c) * dim + d] *= w; 53162bd480fSMatthew G. Knepley if (debug) { 53263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q)); 53362bd480fSMatthew G. Knepley for (c = 0; c < NcI; ++c) { 53463a3b9bcSJacob Faibussowitsch if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c]))); 53562bd480fSMatthew G. Knepley if (n1) { 53663a3b9bcSJacob 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]))); 5379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 53862bd480fSMatthew G. Knepley } 53962bd480fSMatthew G. Knepley } 54062bd480fSMatthew G. Knepley } 54120cf1dd8SToby Isaac } 5429566063dSJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 54320cf1dd8SToby Isaac cOffset += totDim; 54420cf1dd8SToby Isaac cOffsetAux += totDimAux; 54520cf1dd8SToby Isaac } 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54720cf1dd8SToby Isaac } 54820cf1dd8SToby Isaac 54927f02ce8SMatthew G. Knepley /* 55027f02ce8SMatthew 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 55127f02ce8SMatthew G. Knepley all transforms operate in the full space and are square. 55227f02ce8SMatthew G. Knepley 55327f02ce8SMatthew G. Knepley HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 55427f02ce8SMatthew G. Knepley 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 55527f02ce8SMatthew G. Knepley 2) We need to assume that the orientation is 0 for both 55627f02ce8SMatthew 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() 55727f02ce8SMatthew G. Knepley */ 558*07218a29SMatthew 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[]) 559d71ae5a4SJacob Faibussowitsch { 56027f02ce8SMatthew G. Knepley const PetscInt debug = 0; 5616528b96dSMatthew G. Knepley const PetscInt field = key.field; 56227f02ce8SMatthew G. Knepley PetscFE fe; 5636528b96dSMatthew G. Knepley PetscWeakForm wf; 5646528b96dSMatthew G. Knepley PetscInt n0, n1, i; 5656528b96dSMatthew G. Knepley PetscBdPointFunc *f0_func, *f1_func; 56627f02ce8SMatthew G. Knepley PetscQuadrature quad; 567*07218a29SMatthew G. Knepley PetscTabulation *Tf, *TfIn, *TfAux = NULL; 56827f02ce8SMatthew G. Knepley PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 56927f02ce8SMatthew G. Knepley const PetscScalar *constants; 57027f02ce8SMatthew G. Knepley PetscReal *x; 571665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 572*07218a29SMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 5736587ee25SMatthew G. Knepley PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 57427f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 5756587ee25SMatthew G. Knepley PetscInt qdim, qNc, Nq, q, dE; 57627f02ce8SMatthew G. Knepley 57727f02ce8SMatthew G. Knepley PetscFunctionBegin; 57827f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 5799566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 5809566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 5819566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 5829566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 5839566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 584*07218a29SMatthew G. Knepley PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn)); 585*07218a29SMatthew G. Knepley PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, s, &uOff)); 586*07218a29SMatthew G. Knepley PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x)); 5879566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset)); 5889566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 5899566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 5903ba16761SJacob Faibussowitsch if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 5919566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 5929566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 5939566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 59427f02ce8SMatthew G. Knepley /* NOTE This is a bulk tabulation because the DS is a face discretization */ 5959566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &Tf)); 596*07218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 5979566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 59827f02ce8SMatthew G. Knepley if (dsAux) { 5999566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 6009566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 6019566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 6029566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 6039566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 6049566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 60501907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 6069566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 6079566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 60863a3b9bcSJacob 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); 60927f02ce8SMatthew G. Knepley } 6109566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField)); 611665f567fSMatthew G. Knepley NcI = Tf[field]->Nc; 612c2b7495fSMatthew G. Knepley NcS = NcI; 6139566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 61463a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 61527f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 61663a3b9bcSJacob Faibussowitsch PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 61727f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 61827f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 619*07218a29SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][1]}; 62027f02ce8SMatthew G. Knepley 6216587ee25SMatthew G. Knepley fegeom.v = x; /* Workspace */ 6229566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f0, Nq * NcS)); 6239566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(f1, Nq * NcS * dE)); 62427f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 625*07218a29SMatthew G. Knepley const PetscInt qpt[2] = {q, Nq - q - 1}; 62627f02ce8SMatthew G. Knepley PetscReal w; 62727f02ce8SMatthew G. Knepley PetscInt c, d; 62827f02ce8SMatthew G. Knepley 629*07218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 63027f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 6316587ee25SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 63263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 63327f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 6349566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ)); 63527f02ce8SMatthew G. Knepley #endif 63627f02ce8SMatthew G. Knepley } 637a4158a15SMatthew 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])); 63827f02ce8SMatthew G. Knepley /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 639*07218a29SMatthew 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)); 640*07218a29SMatthew 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)); 6416528b96dSMatthew 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]); 64227f02ce8SMatthew G. Knepley for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w; 6439ee2af8cSMatthew 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]); 6449371c9d4SSatish Balay for (c = 0; c < NcS; ++c) 6459371c9d4SSatish Balay for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w; 64627f02ce8SMatthew G. Knepley } 6479371c9d4SSatish Balay if (isCohesiveField) { 6483ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 6499371c9d4SSatish Balay } else { 6503ba16761SJacob Faibussowitsch PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 6519371c9d4SSatish Balay } 65227f02ce8SMatthew G. Knepley cOffset += totDim; 653*07218a29SMatthew G. Knepley cOffsetIn += totDimIn; 65427f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 65527f02ce8SMatthew G. Knepley } 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65727f02ce8SMatthew G. Knepley } 65827f02ce8SMatthew G. Knepley 659d71ae5a4SJacob 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[]) 660d71ae5a4SJacob Faibussowitsch { 66120cf1dd8SToby Isaac const PetscInt debug = 0; 6624bee2e38SMatthew G. Knepley PetscFE feI, feJ; 6636528b96dSMatthew G. Knepley PetscWeakForm wf; 6646528b96dSMatthew G. Knepley PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 6656528b96dSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 66620cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 66720cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 66820cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 66920cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 67020cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 67120cf1dd8SToby Isaac PetscQuadrature quad; 672ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 6734bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 67420cf1dd8SToby Isaac const PetscScalar *constants; 67520cf1dd8SToby Isaac PetscReal *x; 676ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 677ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 6786528b96dSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 67920cf1dd8SToby Isaac PetscInt dE, Np; 68020cf1dd8SToby Isaac PetscBool isAffine; 68120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 68220cf1dd8SToby Isaac PetscInt qNc, Nq, q; 68320cf1dd8SToby Isaac 68420cf1dd8SToby Isaac PetscFunctionBegin; 6859566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 6866528b96dSMatthew G. Knepley fieldI = key.field / Nf; 6876528b96dSMatthew G. Knepley fieldJ = key.field % Nf; 6889566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 6899566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 6909566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 6919566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 6929566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 6939566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 6949566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 6959566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 69620cf1dd8SToby Isaac switch (jtype) { 697d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 698d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 699d71ae5a4SJacob Faibussowitsch break; 700d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 701d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 702d71ae5a4SJacob Faibussowitsch break; 703d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 704d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 705d71ae5a4SJacob Faibussowitsch break; 70620cf1dd8SToby Isaac } 7073ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 7089566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 7099566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 7109566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 7119566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 7129566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 7139566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 7149566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 7154bee2e38SMatthew G. Knepley if (dsAux) { 7169566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 7179566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 7189566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 7199566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 7209566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 7219566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 72263a3b9bcSJacob 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); 72320cf1dd8SToby Isaac } 72427f02ce8SMatthew G. Knepley NcI = T[fieldI]->Nc; 72527f02ce8SMatthew G. Knepley NcJ = T[fieldJ]->Nc; 7264bee2e38SMatthew G. Knepley Np = cgeom->numPoints; 7274bee2e38SMatthew G. Knepley dE = cgeom->dimEmbed; 7284bee2e38SMatthew G. Knepley isAffine = cgeom->isAffine; 72927f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 7309566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 7319566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 7329566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 7339566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 7349566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 73563a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 7364bee2e38SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 7374bee2e38SMatthew G. Knepley PetscFEGeom fegeom; 7384bee2e38SMatthew G. Knepley 73927f02ce8SMatthew G. Knepley fegeom.dim = cgeom->dim; 74027f02ce8SMatthew G. Knepley fegeom.dimEmbed = cgeom->dimEmbed; 7414bee2e38SMatthew G. Knepley if (isAffine) { 7424bee2e38SMatthew G. Knepley fegeom.v = x; 7434bee2e38SMatthew G. Knepley fegeom.xi = cgeom->xi; 7447132c3f7SMatthew G. Knepley fegeom.J = &cgeom->J[e * Np * dE * dE]; 7457132c3f7SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 7467132c3f7SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np]; 7474bee2e38SMatthew G. Knepley } 74820cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 74920cf1dd8SToby Isaac PetscReal w; 7504bee2e38SMatthew G. Knepley PetscInt c; 75120cf1dd8SToby Isaac 75220cf1dd8SToby Isaac if (isAffine) { 7537132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 75420cf1dd8SToby Isaac } else { 7554bee2e38SMatthew G. Knepley fegeom.v = &cgeom->v[(e * Np + q) * dE]; 7564bee2e38SMatthew G. Knepley fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 7574bee2e38SMatthew G. Knepley fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 7584bee2e38SMatthew G. Knepley fegeom.detJ = &cgeom->detJ[e * Np + q]; 75920cf1dd8SToby Isaac } 7609566063dSJacob 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])); 7614bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 7629566063dSJacob Faibussowitsch if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 7639566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 764ea672e62SMatthew G. Knepley if (n0) { 7659566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 7666528b96dSMatthew 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); 76720cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 76820cf1dd8SToby Isaac } 769ea672e62SMatthew G. Knepley if (n1) { 7709566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 7716528b96dSMatthew 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); 7724bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 77320cf1dd8SToby Isaac } 774ea672e62SMatthew G. Knepley if (n2) { 7759566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 7766528b96dSMatthew 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); 7774bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 77820cf1dd8SToby Isaac } 779ea672e62SMatthew G. Knepley if (n3) { 7809566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 7816528b96dSMatthew 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); 7824bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 78320cf1dd8SToby Isaac } 78420cf1dd8SToby Isaac 7859566063dSJacob 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)); 78620cf1dd8SToby Isaac } 78720cf1dd8SToby Isaac if (debug > 1) { 78820cf1dd8SToby Isaac PetscInt fc, f, gc, g; 78920cf1dd8SToby Isaac 79063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 791ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 792ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 793ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 794ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 795ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 796ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 79763a3b9bcSJacob 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]))); 79820cf1dd8SToby Isaac } 79920cf1dd8SToby Isaac } 8009566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 80120cf1dd8SToby Isaac } 80220cf1dd8SToby Isaac } 80320cf1dd8SToby Isaac } 80420cf1dd8SToby Isaac cOffset += totDim; 80520cf1dd8SToby Isaac cOffsetAux += totDimAux; 80620cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 80720cf1dd8SToby Isaac } 8083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80920cf1dd8SToby Isaac } 81020cf1dd8SToby Isaac 811d71ae5a4SJacob 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[]) 812d71ae5a4SJacob Faibussowitsch { 81320cf1dd8SToby Isaac const PetscInt debug = 0; 8144bee2e38SMatthew G. Knepley PetscFE feI, feJ; 81545480ffeSMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 81645480ffeSMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 81720cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 81820cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 81920cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 82020cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 82120cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 82220cf1dd8SToby Isaac PetscQuadrature quad; 823ef0bb6c7SMatthew G. Knepley PetscTabulation *T, *TAux = NULL; 8244bee2e38SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 82520cf1dd8SToby Isaac const PetscScalar *constants; 82620cf1dd8SToby Isaac PetscReal *x; 827ef0bb6c7SMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 828ef0bb6c7SMatthew G. Knepley PetscInt NcI = 0, NcJ = 0; 82945480ffeSMatthew G. Knepley PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 83020cf1dd8SToby Isaac PetscBool isAffine; 83120cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 83220cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 83320cf1dd8SToby Isaac 83420cf1dd8SToby Isaac PetscFunctionBegin; 8359566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 83645480ffeSMatthew G. Knepley fieldI = key.field / Nf; 83745480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 8389566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 8399566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 8409566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 8419566063dSJacob Faibussowitsch PetscCall(PetscFEGetFaceQuadrature(feI, &quad)); 8429566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 8439566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 8449566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 8459566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 8469566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 8479566063dSJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 8483ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 8499566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 8509566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 8519566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 8529566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(ds, &T)); 8539566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 8544bee2e38SMatthew G. Knepley if (dsAux) { 8559566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 8569566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 8579566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 8589566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 8599566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 8609566063dSJacob Faibussowitsch PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 86120cf1dd8SToby Isaac } 862ef0bb6c7SMatthew G. Knepley NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 86320cf1dd8SToby Isaac Np = fgeom->numPoints; 86420cf1dd8SToby Isaac dE = fgeom->dimEmbed; 86520cf1dd8SToby Isaac isAffine = fgeom->isAffine; 86627f02ce8SMatthew G. Knepley /* Initialize here in case the function is not defined */ 8679566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 8689566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 8699566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 8709566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 8719566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 87263a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 87320cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 8749f209ee4SMatthew G. Knepley PetscFEGeom fegeom, cgeom; 87520cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 876ea78f98cSLisandro Dalcin fegeom.n = NULL; 877ea78f98cSLisandro Dalcin fegeom.v = NULL; 878ea78f98cSLisandro Dalcin fegeom.J = NULL; 879ea78f98cSLisandro Dalcin fegeom.detJ = NULL; 88027f02ce8SMatthew G. Knepley fegeom.dim = fgeom->dim; 88127f02ce8SMatthew G. Knepley fegeom.dimEmbed = fgeom->dimEmbed; 88227f02ce8SMatthew G. Knepley cgeom.dim = fgeom->dim; 88327f02ce8SMatthew G. Knepley cgeom.dimEmbed = fgeom->dimEmbed; 8844bee2e38SMatthew G. Knepley if (isAffine) { 8854bee2e38SMatthew G. Knepley fegeom.v = x; 8864bee2e38SMatthew G. Knepley fegeom.xi = fgeom->xi; 8877132c3f7SMatthew G. Knepley fegeom.J = &fgeom->J[e * Np * dE * dE]; 8887132c3f7SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 8897132c3f7SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np]; 8907132c3f7SMatthew G. Knepley fegeom.n = &fgeom->n[e * Np * dE]; 8919f209ee4SMatthew G. Knepley 8927132c3f7SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 8937132c3f7SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 8947132c3f7SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 8954bee2e38SMatthew G. Knepley } 89620cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 89720cf1dd8SToby Isaac PetscReal w; 8984bee2e38SMatthew G. Knepley PetscInt c; 89920cf1dd8SToby Isaac 90063a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 90120cf1dd8SToby Isaac if (isAffine) { 9027132c3f7SMatthew G. Knepley CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 90320cf1dd8SToby Isaac } else { 9043fe841f2SMatthew G. Knepley fegeom.v = &fgeom->v[(e * Np + q) * dE]; 9059f209ee4SMatthew G. Knepley fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 9069f209ee4SMatthew G. Knepley fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 9074bee2e38SMatthew G. Knepley fegeom.detJ = &fgeom->detJ[e * Np + q]; 9084bee2e38SMatthew G. Knepley fegeom.n = &fgeom->n[(e * Np + q) * dE]; 9099f209ee4SMatthew G. Knepley 9109f209ee4SMatthew G. Knepley cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 9119f209ee4SMatthew G. Knepley cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 9129f209ee4SMatthew G. Knepley cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 91320cf1dd8SToby Isaac } 9144bee2e38SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 9159566063dSJacob Faibussowitsch if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 9169566063dSJacob Faibussowitsch if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 917ea672e62SMatthew G. Knepley if (n0) { 9189566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcI * NcJ)); 91945480ffeSMatthew 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); 92020cf1dd8SToby Isaac for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 92120cf1dd8SToby Isaac } 922ea672e62SMatthew G. Knepley if (n1) { 9239566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 92445480ffeSMatthew 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); 9254bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 92620cf1dd8SToby Isaac } 927ea672e62SMatthew G. Knepley if (n2) { 9289566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 92945480ffeSMatthew 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); 9304bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 93120cf1dd8SToby Isaac } 932ea672e62SMatthew G. Knepley if (n3) { 9339566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 93445480ffeSMatthew 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); 9354bee2e38SMatthew G. Knepley for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 93620cf1dd8SToby Isaac } 93720cf1dd8SToby Isaac 9389566063dSJacob 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)); 93920cf1dd8SToby Isaac } 94020cf1dd8SToby Isaac if (debug > 1) { 94120cf1dd8SToby Isaac PetscInt fc, f, gc, g; 94220cf1dd8SToby Isaac 94363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 944ef0bb6c7SMatthew G. Knepley for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 945ef0bb6c7SMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 946ef0bb6c7SMatthew G. Knepley const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 947ef0bb6c7SMatthew G. Knepley for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 948ef0bb6c7SMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 949ef0bb6c7SMatthew G. Knepley const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 95063a3b9bcSJacob 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]))); 95120cf1dd8SToby Isaac } 95220cf1dd8SToby Isaac } 9539566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 95420cf1dd8SToby Isaac } 95520cf1dd8SToby Isaac } 95620cf1dd8SToby Isaac } 95720cf1dd8SToby Isaac cOffset += totDim; 95820cf1dd8SToby Isaac cOffsetAux += totDimAux; 95920cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 96020cf1dd8SToby Isaac } 9613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96220cf1dd8SToby Isaac } 96320cf1dd8SToby Isaac 964*07218a29SMatthew G. Knepley 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[]) 965d71ae5a4SJacob Faibussowitsch { 96627f02ce8SMatthew G. Knepley const PetscInt debug = 0; 96727f02ce8SMatthew G. Knepley PetscFE feI, feJ; 968148442b3SMatthew G. Knepley PetscWeakForm wf; 969148442b3SMatthew G. Knepley PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 970148442b3SMatthew G. Knepley PetscInt n0, n1, n2, n3, i; 97127f02ce8SMatthew G. Knepley PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 97227f02ce8SMatthew G. Knepley PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 97327f02ce8SMatthew G. Knepley PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 97427f02ce8SMatthew G. Knepley PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 97527f02ce8SMatthew G. Knepley PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 976665f567fSMatthew G. Knepley PetscQuadrature quad; 977*07218a29SMatthew G. Knepley PetscTabulation *T, *TfIn, *TAux = NULL; 97827f02ce8SMatthew G. Knepley PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 97927f02ce8SMatthew G. Knepley const PetscScalar *constants; 98027f02ce8SMatthew G. Knepley PetscReal *x; 981665f567fSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 982665f567fSMatthew G. Knepley PetscInt NcI = 0, NcJ = 0, NcS, NcT; 98345480ffeSMatthew G. Knepley PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 984*07218a29SMatthew G. Knepley PetscBool isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE; 98527f02ce8SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 986*07218a29SMatthew G. Knepley PetscInt qNc, Nq, q, dE; 98727f02ce8SMatthew G. Knepley 98827f02ce8SMatthew G. Knepley PetscFunctionBegin; 9899566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 99045480ffeSMatthew G. Knepley fieldI = key.field / Nf; 99145480ffeSMatthew G. Knepley fieldJ = key.field % Nf; 99227f02ce8SMatthew G. Knepley /* Hybrid discretization is posed directly on faces */ 9939566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 9949566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 9959566063dSJacob Faibussowitsch PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 9969566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(feI, &quad)); 9979566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 9989566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsetsCohesive(ds, s, &uOff)); 9999566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x)); 10009566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 100127f02ce8SMatthew G. Knepley switch (jtype) { 1002d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_PRE: 1003d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1004d71ae5a4SJacob Faibussowitsch break; 1005d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN: 1006d71ae5a4SJacob Faibussowitsch PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1007d71ae5a4SJacob Faibussowitsch break; 1008d71ae5a4SJacob Faibussowitsch case PETSCFE_JACOBIAN_DYN: 1009d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 101027f02ce8SMatthew G. Knepley } 10113ba16761SJacob Faibussowitsch if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 10129566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 10139566063dSJacob Faibussowitsch PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 10149566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 10159566063dSJacob Faibussowitsch PetscCall(PetscDSGetTabulation(ds, &T)); 1016*07218a29SMatthew G. Knepley PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 10179566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI)); 10189566063dSJacob Faibussowitsch PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ)); 10199566063dSJacob Faibussowitsch PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 102027f02ce8SMatthew G. Knepley if (dsAux) { 10219566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 10229566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 10239566063dSJacob Faibussowitsch PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 10249566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 10259566063dSJacob Faibussowitsch PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 10269566063dSJacob Faibussowitsch PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 102701907d53SMatthew G. Knepley auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 10289566063dSJacob Faibussowitsch if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 10299566063dSJacob Faibussowitsch else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 103063a3b9bcSJacob 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); 103127f02ce8SMatthew G. Knepley } 10329566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI)); 10339566063dSJacob Faibussowitsch PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ)); 1034665f567fSMatthew G. Knepley NcI = T[fieldI]->Nc; 1035665f567fSMatthew G. Knepley NcJ = T[fieldJ]->Nc; 103627f02ce8SMatthew G. Knepley NcS = isCohesiveFieldI ? NcI : 2 * NcI; 103727f02ce8SMatthew G. Knepley NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ; 103827f02ce8SMatthew G. Knepley dE = fgeom->dimEmbed; 10399566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 10409566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcS * NcT * dE)); 10419566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcS * NcT * dE)); 10429566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcS * NcT * dE * dE)); 10439566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 104463a3b9bcSJacob Faibussowitsch PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 104527f02ce8SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 104627f02ce8SMatthew G. Knepley PetscFEGeom fegeom; 1047*07218a29SMatthew G. Knepley const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][1]}; 104827f02ce8SMatthew G. Knepley 1049*07218a29SMatthew G. Knepley fegeom.v = x; /* Workspace */ 105027f02ce8SMatthew G. Knepley for (q = 0; q < Nq; ++q) { 1051*07218a29SMatthew G. Knepley const PetscInt qpt[2] = {q, Nq - q - 1}; 105227f02ce8SMatthew G. Knepley PetscReal w; 105327f02ce8SMatthew G. Knepley PetscInt c; 105427f02ce8SMatthew G. Knepley 1055*07218a29SMatthew G. Knepley PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 105627f02ce8SMatthew G. Knepley w = fegeom.detJ[0] * quadWeights[q]; 1057*07218a29SMatthew G. Knepley if (debug > 1 && q < fgeom->numPoints) { 105863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 105927f02ce8SMatthew G. Knepley #if !defined(PETSC_USE_COMPLEX) 10609566063dSJacob Faibussowitsch PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 106127f02ce8SMatthew G. Knepley #endif 106227f02ce8SMatthew G. Knepley } 106363a3b9bcSJacob Faibussowitsch if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 1064*07218a29SMatthew 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)); 1065*07218a29SMatthew 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)); 1066ea672e62SMatthew G. Knepley if (n0) { 10679566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g0, NcS * NcT)); 1068148442b3SMatthew 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); 106927f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT; ++c) g0[c] *= w; 107027f02ce8SMatthew G. Knepley } 1071ea672e62SMatthew G. Knepley if (n1) { 10729566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g1, NcS * NcT * dE)); 1073148442b3SMatthew 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); 107427f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT * dE; ++c) g1[c] *= w; 107527f02ce8SMatthew G. Knepley } 1076ea672e62SMatthew G. Knepley if (n2) { 10779566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g2, NcS * NcT * dE)); 1078148442b3SMatthew 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); 107927f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT * dE; ++c) g2[c] *= w; 108027f02ce8SMatthew G. Knepley } 1081ea672e62SMatthew G. Knepley if (n3) { 10829566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(g3, NcS * NcT * dE * dE)); 1083148442b3SMatthew 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); 108427f02ce8SMatthew G. Knepley for (c = 0; c < NcS * NcT * dE * dE; ++c) g3[c] *= w; 108527f02ce8SMatthew G. Knepley } 108627f02ce8SMatthew G. Knepley 10875fedec97SMatthew G. Knepley if (isCohesiveFieldI) { 10885fedec97SMatthew G. Knepley if (isCohesiveFieldJ) { 10899566063dSJacob 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)); 109027f02ce8SMatthew G. Knepley } else { 10919566063dSJacob 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)); 10929566063dSJacob 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)); 10935fedec97SMatthew G. Knepley } 10949371c9d4SSatish Balay } else 10959371c9d4SSatish 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)); 109627f02ce8SMatthew G. Knepley } 109727f02ce8SMatthew G. Knepley if (debug > 1) { 109827f02ce8SMatthew G. Knepley PetscInt fc, f, gc, g; 109927f02ce8SMatthew G. Knepley 110063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 110127f02ce8SMatthew G. Knepley for (fc = 0; fc < NcI; ++fc) { 1102665f567fSMatthew G. Knepley for (f = 0; f < T[fieldI]->Nb; ++f) { 110327f02ce8SMatthew G. Knepley const PetscInt i = offsetI + f * NcI + fc; 110427f02ce8SMatthew G. Knepley for (gc = 0; gc < NcJ; ++gc) { 1105665f567fSMatthew G. Knepley for (g = 0; g < T[fieldJ]->Nb; ++g) { 110627f02ce8SMatthew G. Knepley const PetscInt j = offsetJ + g * NcJ + gc; 110763a3b9bcSJacob 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]))); 110827f02ce8SMatthew G. Knepley } 110927f02ce8SMatthew G. Knepley } 11109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 111127f02ce8SMatthew G. Knepley } 111227f02ce8SMatthew G. Knepley } 111327f02ce8SMatthew G. Knepley } 111427f02ce8SMatthew G. Knepley cOffset += totDim; 111527f02ce8SMatthew G. Knepley cOffsetAux += totDimAux; 111627f02ce8SMatthew G. Knepley eOffset += PetscSqr(totDim); 111727f02ce8SMatthew G. Knepley } 11183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 111927f02ce8SMatthew G. Knepley } 112027f02ce8SMatthew G. Knepley 1121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 1122d71ae5a4SJacob Faibussowitsch { 112320cf1dd8SToby Isaac PetscFunctionBegin; 112420cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 112520cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 112620cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 112720cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 112820cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 1129ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Basic; 113020cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 1131afe6d6adSToby Isaac fem->ops->integratebd = PetscFEIntegrateBd_Basic; 113220cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 113320cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 113427f02ce8SMatthew G. Knepley fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 113520cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */; 113620cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 113720cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 113827f02ce8SMatthew G. Knepley fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 11393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114020cf1dd8SToby Isaac } 114120cf1dd8SToby Isaac 114220cf1dd8SToby Isaac /*MC 1143dce8aebaSBarry Smith PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization 114420cf1dd8SToby Isaac 114520cf1dd8SToby Isaac Level: intermediate 114620cf1dd8SToby Isaac 1147dce8aebaSBarry Smith .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 114820cf1dd8SToby Isaac M*/ 114920cf1dd8SToby Isaac 1150d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 1151d71ae5a4SJacob Faibussowitsch { 115220cf1dd8SToby Isaac PetscFE_Basic *b; 115320cf1dd8SToby Isaac 115420cf1dd8SToby Isaac PetscFunctionBegin; 115520cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 11564dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 115720cf1dd8SToby Isaac fem->data = b; 115820cf1dd8SToby Isaac 11599566063dSJacob Faibussowitsch PetscCall(PetscFEInitialize_Basic(fem)); 11603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116120cf1dd8SToby Isaac } 1162