120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/ 320cf1dd8SToby Isaac #include <petscblaslapack.h> 4012bc364SMatthew G. Knepley #include <petscdmplextransform.h> 520cf1dd8SToby Isaac 69371c9d4SSatish Balay static PetscErrorCode PetscFEDestroy_Composite(PetscFE fem) { 720cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data; 820cf1dd8SToby Isaac 920cf1dd8SToby Isaac PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(PetscFree(cmp->embedding)); 119566063dSJacob Faibussowitsch PetscCall(PetscFree(cmp)); 1220cf1dd8SToby Isaac PetscFunctionReturn(0); 1320cf1dd8SToby Isaac } 1420cf1dd8SToby Isaac 159371c9d4SSatish Balay static PetscErrorCode PetscFESetUp_Composite(PetscFE fem) { 1620cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data; 1720cf1dd8SToby Isaac DM K; 18412e9a14SMatthew G. Knepley DMPolytopeType ct; 19012bc364SMatthew G. Knepley DMPlexTransform tr; 2020cf1dd8SToby Isaac PetscReal *subpoint; 2120cf1dd8SToby Isaac PetscBLASInt *pivots; 2220cf1dd8SToby Isaac PetscBLASInt n, info; 2320cf1dd8SToby Isaac PetscScalar *work, *invVscalar; 2420cf1dd8SToby Isaac PetscInt dim, pdim, spdim, j, s; 25b4457527SToby Isaac PetscSection section; 2620cf1dd8SToby Isaac 2720cf1dd8SToby Isaac PetscFunctionBegin; 2820cf1dd8SToby Isaac /* Get affine mapping from reference cell to each subcell */ 299566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &K)); 309566063dSJacob Faibussowitsch PetscCall(DMGetDimension(K, &dim)); 319566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(K, 0, &ct)); 329566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &tr)); 339566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR)); 349566063dSJacob Faibussowitsch PetscCall(DMPlexRefineRegularGetAffineTransforms(tr, ct, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac)); 359566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 3620cf1dd8SToby Isaac /* Determine dof embedding into subelements */ 379566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 389566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(fem->basisSpace, &spdim)); 399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cmp->numSubelements * spdim, &cmp->embedding)); 409566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(K, dim, MPIU_REAL, &subpoint)); 419566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSection(fem->dualSpace, §ion)); 4220cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 4320cf1dd8SToby Isaac PetscInt sd = 0; 44b4457527SToby Isaac PetscInt closureSize; 45b4457527SToby Isaac PetscInt *closure = NULL; 4620cf1dd8SToby Isaac 479566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure)); 48b4457527SToby Isaac for (j = 0; j < closureSize; j++) { 49b4457527SToby Isaac PetscInt point = closure[2 * j]; 50b4457527SToby Isaac PetscInt dof, off, k; 5120cf1dd8SToby Isaac 529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, point, &dof)); 539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, point, &off)); 54b4457527SToby Isaac for (k = 0; k < dof; k++) cmp->embedding[s * spdim + sd++] = off + k; 5520cf1dd8SToby Isaac } 569566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure)); 5763a3b9bcSJacob Faibussowitsch PetscCheck(sd == spdim, PetscObjectComm((PetscObject)fem), PETSC_ERR_PLIB, "Subelement %" PetscInt_FMT " has %" PetscInt_FMT " dual basis vectors != %" PetscInt_FMT, s, sd, spdim); 5820cf1dd8SToby Isaac } 599566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint)); 6020cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis for each subelement */ 619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cmp->numSubelements * spdim * spdim, &fem->invV)); 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(spdim, &pivots, spdim, &work)); 6320cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cmp->numSubelements * spdim * spdim, &invVscalar)); 6520cf1dd8SToby Isaac #else 6620cf1dd8SToby Isaac invVscalar = fem->invV; 6720cf1dd8SToby Isaac #endif 6820cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 6920cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 7020cf1dd8SToby Isaac PetscReal *Bf; 7120cf1dd8SToby Isaac PetscQuadrature f; 7220cf1dd8SToby Isaac const PetscReal *points, *weights; 7320cf1dd8SToby Isaac PetscInt Nc, Nq, q, k; 7420cf1dd8SToby Isaac 759566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s * spdim + j], &f)); 769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights)); 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(f->numPoints * spdim * Nc, &Bf)); 789566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL)); 7920cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 8020cf1dd8SToby Isaac /* n_j \cdot \phi_k */ 8120cf1dd8SToby Isaac invVscalar[(s * spdim + j) * spdim + k] = 0.0; 82ad540459SPierre Jolivet for (q = 0; q < Nq; ++q) invVscalar[(s * spdim + j) * spdim + k] += Bf[q * spdim + k] * weights[q]; 8320cf1dd8SToby Isaac } 849566063dSJacob Faibussowitsch PetscCall(PetscFree(Bf)); 8520cf1dd8SToby Isaac } 8620cf1dd8SToby Isaac n = spdim; 87792fecdfSBarry Smith PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s * spdim * spdim], &n, pivots, &info)); 88792fecdfSBarry Smith PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s * spdim * spdim], &n, pivots, work, &n, &info)); 8920cf1dd8SToby Isaac } 9020cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 9120cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements * spdim * spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]); 929566063dSJacob Faibussowitsch PetscCall(PetscFree(invVscalar)); 9320cf1dd8SToby Isaac #endif 949566063dSJacob Faibussowitsch PetscCall(PetscFree2(pivots, work)); 9520cf1dd8SToby Isaac PetscFunctionReturn(0); 9620cf1dd8SToby Isaac } 9720cf1dd8SToby Isaac 989371c9d4SSatish Balay static PetscErrorCode PetscFECreateTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) { 9920cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data; 10020cf1dd8SToby Isaac DM dm; 101412e9a14SMatthew G. Knepley DMPolytopeType ct; 10220cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 10320cf1dd8SToby Isaac PetscInt spdim; /* Dimension of subelement FE space P */ 10420cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 10520cf1dd8SToby Isaac PetscInt comp; /* Field components */ 10620cf1dd8SToby Isaac PetscInt *subpoints; 107ef0bb6c7SMatthew G. Knepley PetscReal *B = K >= 0 ? T->T[0] : NULL; 108ef0bb6c7SMatthew G. Knepley PetscReal *D = K >= 1 ? T->T[1] : NULL; 109ef0bb6c7SMatthew G. Knepley PetscReal *H = K >= 2 ? T->T[2] : NULL; 110ef0bb6c7SMatthew G. Knepley PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL, *subpoint; 11120cf1dd8SToby Isaac PetscInt p, s, d, e, j, k; 11220cf1dd8SToby Isaac 11320cf1dd8SToby Isaac PetscFunctionBegin; 1149566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 1159566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1169566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, 0, &ct)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(fem->basisSpace, &spdim)); 1189566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 1199566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fem, &comp)); 12020cf1dd8SToby Isaac /* Divide points into subelements */ 1219566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints)); 1229566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint)); 12320cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 12420cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 12520cf1dd8SToby Isaac PetscBool inside; 12620cf1dd8SToby Isaac 12720cf1dd8SToby Isaac /* Apply transform, and check that point is inside cell */ 12820cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 12920cf1dd8SToby Isaac subpoint[d] = -1.0; 13020cf1dd8SToby Isaac for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s * dim + d) * dim + e] * (points[p * dim + e] - cmp->v0[s * dim + e]); 13120cf1dd8SToby Isaac } 1329566063dSJacob Faibussowitsch PetscCall(DMPolytopeInCellTest(ct, subpoint, &inside)); 1339371c9d4SSatish Balay if (inside) { 1349371c9d4SSatish Balay subpoints[p] = s; 1359371c9d4SSatish Balay break; 1369371c9d4SSatish Balay } 13720cf1dd8SToby Isaac } 13863a3b9bcSJacob Faibussowitsch PetscCheck(s < cmp->numSubelements, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " was not found in any subelement", p); 13920cf1dd8SToby Isaac } 1409566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint)); 14120cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 1429566063dSJacob Faibussowitsch if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * spdim, MPIU_REAL, &tmpB)); 1439566063dSJacob Faibussowitsch if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * spdim * dim, MPIU_REAL, &tmpD)); 1449566063dSJacob Faibussowitsch if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * spdim * dim * dim, MPIU_REAL, &tmpH)); 1459566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH)); 14620cf1dd8SToby Isaac /* Translate to the nodal basis */ 1479566063dSJacob Faibussowitsch if (K >= 0) PetscCall(PetscArrayzero(B, npoints * pdim * comp)); 1489566063dSJacob Faibussowitsch if (K >= 1) PetscCall(PetscArrayzero(D, npoints * pdim * comp * dim)); 1499566063dSJacob Faibussowitsch if (K >= 2) PetscCall(PetscArrayzero(H, npoints * pdim * comp * dim * dim)); 15020cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 15120cf1dd8SToby Isaac const PetscInt s = subpoints[p]; 15220cf1dd8SToby Isaac 15320cf1dd8SToby Isaac if (B) { 15420cf1dd8SToby Isaac /* Multiply by V^{-1} (spdim x spdim) */ 15520cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 15620cf1dd8SToby Isaac const PetscInt i = (p * pdim + cmp->embedding[s * spdim + j]) * comp; 15720cf1dd8SToby Isaac 15820cf1dd8SToby Isaac B[i] = 0.0; 159ad540459SPierre Jolivet for (k = 0; k < spdim; ++k) B[i] += fem->invV[(s * spdim + k) * spdim + j] * tmpB[p * spdim + k]; 16020cf1dd8SToby Isaac } 16120cf1dd8SToby Isaac } 16220cf1dd8SToby Isaac if (D) { 16320cf1dd8SToby Isaac /* Multiply by V^{-1} (spdim x spdim) */ 16420cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 16520cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 16620cf1dd8SToby Isaac const PetscInt i = ((p * pdim + cmp->embedding[s * spdim + j]) * comp + 0) * dim + d; 16720cf1dd8SToby Isaac 16820cf1dd8SToby Isaac D[i] = 0.0; 169ad540459SPierre Jolivet for (k = 0; k < spdim; ++k) D[i] += fem->invV[(s * spdim + k) * spdim + j] * tmpD[(p * spdim + k) * dim + d]; 17020cf1dd8SToby Isaac } 17120cf1dd8SToby Isaac } 17220cf1dd8SToby Isaac } 17320cf1dd8SToby Isaac if (H) { 17420cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 17520cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 17620cf1dd8SToby Isaac for (d = 0; d < dim * dim; ++d) { 17720cf1dd8SToby Isaac const PetscInt i = ((p * pdim + cmp->embedding[s * spdim + j]) * comp + 0) * dim * dim + d; 17820cf1dd8SToby Isaac 17920cf1dd8SToby Isaac H[i] = 0.0; 180ad540459SPierre Jolivet for (k = 0; k < spdim; ++k) H[i] += fem->invV[(s * spdim + k) * spdim + j] * tmpH[(p * spdim + k) * dim * dim + d]; 18120cf1dd8SToby Isaac } 18220cf1dd8SToby Isaac } 18320cf1dd8SToby Isaac } 18420cf1dd8SToby Isaac } 1859566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints)); 1869566063dSJacob Faibussowitsch if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * spdim, MPIU_REAL, &tmpB)); 1879566063dSJacob Faibussowitsch if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * spdim * dim, MPIU_REAL, &tmpD)); 1889566063dSJacob Faibussowitsch if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * spdim * dim * dim, MPIU_REAL, &tmpH)); 18920cf1dd8SToby Isaac PetscFunctionReturn(0); 19020cf1dd8SToby Isaac } 19120cf1dd8SToby Isaac 1929371c9d4SSatish Balay static PetscErrorCode PetscFEInitialize_Composite(PetscFE fem) { 19320cf1dd8SToby Isaac PetscFunctionBegin; 19420cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 19520cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Composite; 19620cf1dd8SToby Isaac fem->ops->view = NULL; 19720cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Composite; 19820cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 199ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Composite; 20020cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 20120cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 20220cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */; 20320cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 20420cf1dd8SToby Isaac PetscFunctionReturn(0); 20520cf1dd8SToby Isaac } 20620cf1dd8SToby Isaac 20720cf1dd8SToby Isaac /*MC 20820cf1dd8SToby Isaac PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element 20920cf1dd8SToby Isaac 21020cf1dd8SToby Isaac Level: intermediate 21120cf1dd8SToby Isaac 212db781477SPatrick Sanan .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 21320cf1dd8SToby Isaac M*/ 2149371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem) { 21520cf1dd8SToby Isaac PetscFE_Composite *cmp; 21620cf1dd8SToby Isaac 21720cf1dd8SToby Isaac PetscFunctionBegin; 21820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 219*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cmp)); 22020cf1dd8SToby Isaac fem->data = cmp; 22120cf1dd8SToby Isaac 22220cf1dd8SToby Isaac cmp->numSubelements = -1; 22320cf1dd8SToby Isaac cmp->v0 = NULL; 22420cf1dd8SToby Isaac cmp->jac = NULL; 22520cf1dd8SToby Isaac 2269566063dSJacob Faibussowitsch PetscCall(PetscFEInitialize_Composite(fem)); 22720cf1dd8SToby Isaac PetscFunctionReturn(0); 22820cf1dd8SToby Isaac } 22920cf1dd8SToby Isaac 23020cf1dd8SToby Isaac /*@C 23120cf1dd8SToby Isaac PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement 23220cf1dd8SToby Isaac 23320cf1dd8SToby Isaac Not collective 23420cf1dd8SToby Isaac 23520cf1dd8SToby Isaac Input Parameter: 23620cf1dd8SToby Isaac . fem - The PetscFE object 23720cf1dd8SToby Isaac 23820cf1dd8SToby Isaac Output Parameters: 23920cf1dd8SToby Isaac + blockSize - The number of elements in a block 24020cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 24120cf1dd8SToby Isaac . batchSize - The number of elements in a batch 24220cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 24320cf1dd8SToby Isaac 24420cf1dd8SToby Isaac Level: intermediate 24520cf1dd8SToby Isaac 246db781477SPatrick Sanan .seealso: `PetscFECreate()` 24720cf1dd8SToby Isaac @*/ 2489371c9d4SSatish Balay PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[]) { 24920cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *)fem->data; 25020cf1dd8SToby Isaac 25120cf1dd8SToby Isaac PetscFunctionBegin; 25220cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 2539371c9d4SSatish Balay if (numSubelements) { 2549371c9d4SSatish Balay PetscValidIntPointer(numSubelements, 2); 2559371c9d4SSatish Balay *numSubelements = cmp->numSubelements; 2569371c9d4SSatish Balay } 2579371c9d4SSatish Balay if (v0) { 2589371c9d4SSatish Balay PetscValidPointer(v0, 3); 2599371c9d4SSatish Balay *v0 = cmp->v0; 2609371c9d4SSatish Balay } 2619371c9d4SSatish Balay if (jac) { 2629371c9d4SSatish Balay PetscValidPointer(jac, 4); 2639371c9d4SSatish Balay *jac = cmp->jac; 2649371c9d4SSatish Balay } 2659371c9d4SSatish Balay if (invjac) { 2669371c9d4SSatish Balay PetscValidPointer(invjac, 5); 2679371c9d4SSatish Balay *invjac = cmp->invjac; 2689371c9d4SSatish Balay } 26920cf1dd8SToby Isaac PetscFunctionReturn(0); 27020cf1dd8SToby Isaac } 271