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 62b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Composite(PetscFE fem) 720cf1dd8SToby Isaac { 820cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 920cf1dd8SToby Isaac 1020cf1dd8SToby Isaac PetscFunctionBegin; 115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cmp->embedding)); 125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cmp)); 1320cf1dd8SToby Isaac PetscFunctionReturn(0); 1420cf1dd8SToby Isaac } 1520cf1dd8SToby Isaac 162b99622eSMatthew G. Knepley static PetscErrorCode PetscFESetUp_Composite(PetscFE fem) 1720cf1dd8SToby Isaac { 1820cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 1920cf1dd8SToby Isaac DM K; 20412e9a14SMatthew G. Knepley DMPolytopeType ct; 21012bc364SMatthew G. Knepley DMPlexTransform tr; 2220cf1dd8SToby Isaac PetscReal *subpoint; 2320cf1dd8SToby Isaac PetscBLASInt *pivots; 2420cf1dd8SToby Isaac PetscBLASInt n, info; 2520cf1dd8SToby Isaac PetscScalar *work, *invVscalar; 2620cf1dd8SToby Isaac PetscInt dim, pdim, spdim, j, s; 27b4457527SToby Isaac PetscSection section; 2820cf1dd8SToby Isaac 2920cf1dd8SToby Isaac PetscFunctionBegin; 3020cf1dd8SToby Isaac /* Get affine mapping from reference cell to each subcell */ 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDM(fem->dualSpace, &K)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(K, &dim)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(K, 0, &ct)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformCreate(PETSC_COMM_SELF, &tr)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRefineRegularGetAffineTransforms(tr, ct, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTransformDestroy(&tr)); 3820cf1dd8SToby Isaac /* Determine dof embedding into subelements */ 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSpaceGetDimension(fem->basisSpace, &spdim)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(cmp->numSubelements*spdim,&cmp->embedding)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetWorkArray(K, dim, MPIU_REAL, &subpoint)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetSection(fem->dualSpace, §ion)); 4420cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 4520cf1dd8SToby Isaac PetscInt sd = 0; 46b4457527SToby Isaac PetscInt closureSize; 47b4457527SToby Isaac PetscInt *closure = NULL; 4820cf1dd8SToby Isaac 495f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure)); 50b4457527SToby Isaac for (j = 0; j < closureSize; j++) { 51b4457527SToby Isaac PetscInt point = closure[2*j]; 52b4457527SToby Isaac PetscInt dof, off, k; 5320cf1dd8SToby Isaac 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(section, point, &dof)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(section, point, &off)); 56b4457527SToby Isaac for (k = 0; k < dof; k++) cmp->embedding[s*spdim+sd++] = off + k; 5720cf1dd8SToby Isaac } 585f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(K, s, PETSC_TRUE, &closureSize, &closure)); 592c71b3e2SJacob Faibussowitsch PetscCheckFalse(sd != spdim,PetscObjectComm((PetscObject) fem), PETSC_ERR_PLIB, "Subelement %d has %d dual basis vectors != %d", s, sd, spdim); 6020cf1dd8SToby Isaac } 615f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint)); 6220cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis for each subelement */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(cmp->numSubelements*spdim*spdim,&fem->invV)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(spdim,&pivots,spdim,&work)); 6520cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(cmp->numSubelements*spdim*spdim,&invVscalar)); 6720cf1dd8SToby Isaac #else 6820cf1dd8SToby Isaac invVscalar = fem->invV; 6920cf1dd8SToby Isaac #endif 7020cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 7120cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 7220cf1dd8SToby Isaac PetscReal *Bf; 7320cf1dd8SToby Isaac PetscQuadrature f; 7420cf1dd8SToby Isaac const PetscReal *points, *weights; 7520cf1dd8SToby Isaac PetscInt Nc, Nq, q, k; 7620cf1dd8SToby Isaac 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s*spdim+j], &f)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(f->numPoints*spdim*Nc,&Bf)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL)); 8120cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 8220cf1dd8SToby Isaac /* n_j \cdot \phi_k */ 8320cf1dd8SToby Isaac invVscalar[(s*spdim + j)*spdim+k] = 0.0; 8420cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 8520cf1dd8SToby Isaac invVscalar[(s*spdim + j)*spdim+k] += Bf[q*spdim+k]*weights[q]; 8620cf1dd8SToby Isaac } 8720cf1dd8SToby Isaac } 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(Bf)); 8920cf1dd8SToby Isaac } 9020cf1dd8SToby Isaac n = spdim; 9120cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s*spdim*spdim], &n, pivots, &info)); 9220cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s*spdim*spdim], &n, pivots, work, &n, &info)); 9320cf1dd8SToby Isaac } 9420cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 9520cf1dd8SToby Isaac for (s = 0; s <cmp->numSubelements*spdim*spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(invVscalar)); 9720cf1dd8SToby Isaac #endif 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(pivots,work)); 9920cf1dd8SToby Isaac PetscFunctionReturn(0); 10020cf1dd8SToby Isaac } 10120cf1dd8SToby Isaac 102ef0bb6c7SMatthew G. Knepley static PetscErrorCode PetscFECreateTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 10320cf1dd8SToby Isaac { 10420cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 10520cf1dd8SToby Isaac DM dm; 106412e9a14SMatthew G. Knepley DMPolytopeType ct; 10720cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 10820cf1dd8SToby Isaac PetscInt spdim; /* Dimension of subelement FE space P */ 10920cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 11020cf1dd8SToby Isaac PetscInt comp; /* Field components */ 11120cf1dd8SToby Isaac PetscInt *subpoints; 112ef0bb6c7SMatthew G. Knepley PetscReal *B = K >= 0 ? T->T[0] : NULL; 113ef0bb6c7SMatthew G. Knepley PetscReal *D = K >= 1 ? T->T[1] : NULL; 114ef0bb6c7SMatthew G. Knepley PetscReal *H = K >= 2 ? T->T[2] : NULL; 115ef0bb6c7SMatthew G. Knepley PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL, *subpoint; 11620cf1dd8SToby Isaac PetscInt p, s, d, e, j, k; 11720cf1dd8SToby Isaac 11820cf1dd8SToby Isaac PetscFunctionBegin; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, 0, &ct)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSpaceGetDimension(fem->basisSpace, &spdim)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetNumComponents(fem, &comp)); 12520cf1dd8SToby Isaac /* Divide points into subelements */ 1265f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint)); 12820cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 12920cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 13020cf1dd8SToby Isaac PetscBool inside; 13120cf1dd8SToby Isaac 13220cf1dd8SToby Isaac /* Apply transform, and check that point is inside cell */ 13320cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 13420cf1dd8SToby Isaac subpoint[d] = -1.0; 13520cf1dd8SToby 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]); 13620cf1dd8SToby Isaac } 1375f80ce2aSJacob Faibussowitsch CHKERRQ(DMPolytopeInCellTest(ct, subpoint, &inside)); 13820cf1dd8SToby Isaac if (inside) {subpoints[p] = s; break;} 13920cf1dd8SToby Isaac } 1402c71b3e2SJacob Faibussowitsch PetscCheckFalse(s >= cmp->numSubelements,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d was not found in any subelement", p); 14120cf1dd8SToby Isaac } 1425f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint)); 14320cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 1445f80ce2aSJacob Faibussowitsch if (K >= 0) CHKERRQ(DMGetWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB)); 1455f80ce2aSJacob Faibussowitsch if (K >= 1) CHKERRQ(DMGetWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD)); 1465f80ce2aSJacob Faibussowitsch if (K >= 2) CHKERRQ(DMGetWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH)); 14820cf1dd8SToby Isaac /* Translate to the nodal basis */ 1495f80ce2aSJacob Faibussowitsch if (K >= 0) CHKERRQ(PetscArrayzero(B, npoints*pdim*comp)); 1505f80ce2aSJacob Faibussowitsch if (K >= 1) CHKERRQ(PetscArrayzero(D, npoints*pdim*comp*dim)); 1515f80ce2aSJacob Faibussowitsch if (K >= 2) CHKERRQ(PetscArrayzero(H, npoints*pdim*comp*dim*dim)); 15220cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 15320cf1dd8SToby Isaac const PetscInt s = subpoints[p]; 15420cf1dd8SToby Isaac 15520cf1dd8SToby Isaac if (B) { 15620cf1dd8SToby Isaac /* Multiply by V^{-1} (spdim x spdim) */ 15720cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 15820cf1dd8SToby Isaac const PetscInt i = (p*pdim + cmp->embedding[s*spdim+j])*comp; 15920cf1dd8SToby Isaac 16020cf1dd8SToby Isaac B[i] = 0.0; 16120cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 16220cf1dd8SToby Isaac B[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpB[p*spdim + k]; 16320cf1dd8SToby Isaac } 16420cf1dd8SToby Isaac } 16520cf1dd8SToby Isaac } 16620cf1dd8SToby Isaac if (D) { 16720cf1dd8SToby Isaac /* Multiply by V^{-1} (spdim x spdim) */ 16820cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 16920cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 17020cf1dd8SToby Isaac const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim + d; 17120cf1dd8SToby Isaac 17220cf1dd8SToby Isaac D[i] = 0.0; 17320cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 17420cf1dd8SToby Isaac D[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpD[(p*spdim + k)*dim + d]; 17520cf1dd8SToby Isaac } 17620cf1dd8SToby Isaac } 17720cf1dd8SToby Isaac } 17820cf1dd8SToby Isaac } 17920cf1dd8SToby Isaac if (H) { 18020cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 18120cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 18220cf1dd8SToby Isaac for (d = 0; d < dim*dim; ++d) { 18320cf1dd8SToby Isaac const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim*dim + d; 18420cf1dd8SToby Isaac 18520cf1dd8SToby Isaac H[i] = 0.0; 18620cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 18720cf1dd8SToby Isaac H[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpH[(p*spdim + k)*dim*dim + d]; 18820cf1dd8SToby Isaac } 18920cf1dd8SToby Isaac } 19020cf1dd8SToby Isaac } 19120cf1dd8SToby Isaac } 19220cf1dd8SToby Isaac } 1935f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints)); 1945f80ce2aSJacob Faibussowitsch if (K >= 0) CHKERRQ(DMRestoreWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB)); 1955f80ce2aSJacob Faibussowitsch if (K >= 1) CHKERRQ(DMRestoreWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD)); 1965f80ce2aSJacob Faibussowitsch if (K >= 2) CHKERRQ(DMRestoreWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH)); 19720cf1dd8SToby Isaac PetscFunctionReturn(0); 19820cf1dd8SToby Isaac } 19920cf1dd8SToby Isaac 2002b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Composite(PetscFE fem) 20120cf1dd8SToby Isaac { 20220cf1dd8SToby Isaac PetscFunctionBegin; 20320cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 20420cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Composite; 20520cf1dd8SToby Isaac fem->ops->view = NULL; 20620cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Composite; 20720cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 208ef0bb6c7SMatthew G. Knepley fem->ops->createtabulation = PetscFECreateTabulation_Composite; 20920cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 21020cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 21120cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 21220cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 21320cf1dd8SToby Isaac PetscFunctionReturn(0); 21420cf1dd8SToby Isaac } 21520cf1dd8SToby Isaac 21620cf1dd8SToby Isaac /*MC 21720cf1dd8SToby Isaac PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element 21820cf1dd8SToby Isaac 21920cf1dd8SToby Isaac Level: intermediate 22020cf1dd8SToby Isaac 22120cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 22220cf1dd8SToby Isaac M*/ 22320cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem) 22420cf1dd8SToby Isaac { 22520cf1dd8SToby Isaac PetscFE_Composite *cmp; 22620cf1dd8SToby Isaac 22720cf1dd8SToby Isaac PetscFunctionBegin; 22820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(fem, &cmp)); 23020cf1dd8SToby Isaac fem->data = cmp; 23120cf1dd8SToby Isaac 23220cf1dd8SToby Isaac cmp->numSubelements = -1; 23320cf1dd8SToby Isaac cmp->v0 = NULL; 23420cf1dd8SToby Isaac cmp->jac = NULL; 23520cf1dd8SToby Isaac 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEInitialize_Composite(fem)); 23720cf1dd8SToby Isaac PetscFunctionReturn(0); 23820cf1dd8SToby Isaac } 23920cf1dd8SToby Isaac 24020cf1dd8SToby Isaac /*@C 24120cf1dd8SToby Isaac PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement 24220cf1dd8SToby Isaac 24320cf1dd8SToby Isaac Not collective 24420cf1dd8SToby Isaac 24520cf1dd8SToby Isaac Input Parameter: 24620cf1dd8SToby Isaac . fem - The PetscFE object 24720cf1dd8SToby Isaac 24820cf1dd8SToby Isaac Output Parameters: 24920cf1dd8SToby Isaac + blockSize - The number of elements in a block 25020cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 25120cf1dd8SToby Isaac . batchSize - The number of elements in a batch 25220cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 25320cf1dd8SToby Isaac 25420cf1dd8SToby Isaac Level: intermediate 25520cf1dd8SToby Isaac 25620cf1dd8SToby Isaac .seealso: PetscFECreate() 25720cf1dd8SToby Isaac @*/ 25820cf1dd8SToby Isaac PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[]) 25920cf1dd8SToby Isaac { 26020cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 26120cf1dd8SToby Isaac 26220cf1dd8SToby Isaac PetscFunctionBegin; 26320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 264*dadcf809SJacob Faibussowitsch if (numSubelements) {PetscValidIntPointer(numSubelements, 2); *numSubelements = cmp->numSubelements;} 26520cf1dd8SToby Isaac if (v0) {PetscValidPointer(v0, 3); *v0 = cmp->v0;} 26620cf1dd8SToby Isaac if (jac) {PetscValidPointer(jac, 4); *jac = cmp->jac;} 26720cf1dd8SToby Isaac if (invjac) {PetscValidPointer(invjac, 5); *invjac = cmp->invjac;} 26820cf1dd8SToby Isaac PetscFunctionReturn(0); 26920cf1dd8SToby Isaac } 270