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 <petsc/private/dmpleximpl.h> /* For CellRefiner */ 420cf1dd8SToby Isaac #include <petscblaslapack.h> 520cf1dd8SToby Isaac 6*2b99622eSMatthew G. Knepley static PetscErrorCode PetscFEDestroy_Composite(PetscFE fem) 720cf1dd8SToby Isaac { 820cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 920cf1dd8SToby Isaac PetscErrorCode ierr; 1020cf1dd8SToby Isaac 1120cf1dd8SToby Isaac PetscFunctionBegin; 1220cf1dd8SToby Isaac ierr = CellRefinerRestoreAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);CHKERRQ(ierr); 1320cf1dd8SToby Isaac ierr = PetscFree(cmp->embedding);CHKERRQ(ierr); 1420cf1dd8SToby Isaac ierr = PetscFree(cmp);CHKERRQ(ierr); 1520cf1dd8SToby Isaac PetscFunctionReturn(0); 1620cf1dd8SToby Isaac } 1720cf1dd8SToby Isaac 18*2b99622eSMatthew G. Knepley static PetscErrorCode PetscFESetUp_Composite(PetscFE fem) 1920cf1dd8SToby Isaac { 2020cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 2120cf1dd8SToby Isaac DM K; 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; 2720cf1dd8SToby Isaac PetscErrorCode ierr; 2820cf1dd8SToby Isaac 2920cf1dd8SToby Isaac PetscFunctionBegin; 3020cf1dd8SToby Isaac /* Get affine mapping from reference cell to each subcell */ 3120cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &K);CHKERRQ(ierr); 3220cf1dd8SToby Isaac ierr = DMGetDimension(K, &dim);CHKERRQ(ierr); 3320cf1dd8SToby Isaac ierr = DMPlexGetCellRefiner_Internal(K, &cmp->cellRefiner);CHKERRQ(ierr); 3420cf1dd8SToby Isaac ierr = CellRefinerGetAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);CHKERRQ(ierr); 3520cf1dd8SToby Isaac /* Determine dof embedding into subelements */ 3620cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 3720cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(fem->basisSpace, &spdim);CHKERRQ(ierr); 3820cf1dd8SToby Isaac ierr = PetscMalloc1(cmp->numSubelements*spdim,&cmp->embedding);CHKERRQ(ierr); 3920cf1dd8SToby Isaac ierr = DMGetWorkArray(K, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr); 4020cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 4120cf1dd8SToby Isaac PetscInt sd = 0; 4220cf1dd8SToby Isaac 4320cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 4420cf1dd8SToby Isaac PetscBool inside; 4520cf1dd8SToby Isaac PetscQuadrature f; 4620cf1dd8SToby Isaac PetscInt d, e; 4720cf1dd8SToby Isaac 4820cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr); 4920cf1dd8SToby Isaac /* Apply transform to first point, and check that point is inside subcell */ 5020cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 5120cf1dd8SToby Isaac subpoint[d] = -1.0; 5220cf1dd8SToby Isaac for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(f->points[e] - cmp->v0[s*dim+e]); 5320cf1dd8SToby Isaac } 5420cf1dd8SToby Isaac ierr = CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);CHKERRQ(ierr); 5520cf1dd8SToby Isaac if (inside) {cmp->embedding[s*spdim+sd++] = j;} 5620cf1dd8SToby Isaac } 5720cf1dd8SToby Isaac if (sd != spdim) SETERRQ3(PetscObjectComm((PetscObject) fem), PETSC_ERR_PLIB, "Subelement %d has %d dual basis vectors != %d", s, sd, spdim); 5820cf1dd8SToby Isaac } 5920cf1dd8SToby Isaac ierr = DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr); 6020cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis for each subelement */ 6120cf1dd8SToby Isaac ierr = PetscMalloc1(cmp->numSubelements*spdim*spdim,&fem->invV);CHKERRQ(ierr); 6220cf1dd8SToby Isaac ierr = PetscMalloc2(spdim,&pivots,spdim,&work);CHKERRQ(ierr); 6320cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 6420cf1dd8SToby Isaac ierr = PetscMalloc1(cmp->numSubelements*spdim*spdim,&invVscalar);CHKERRQ(ierr); 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 7520cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s*spdim+j], &f);CHKERRQ(ierr); 7620cf1dd8SToby Isaac ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 7720cf1dd8SToby Isaac ierr = PetscMalloc1(f->numPoints*spdim*Nc,&Bf);CHKERRQ(ierr); 7820cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr); 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; 8220cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 8320cf1dd8SToby Isaac invVscalar[(s*spdim + j)*spdim+k] += Bf[q*spdim+k]*weights[q]; 8420cf1dd8SToby Isaac } 8520cf1dd8SToby Isaac } 8620cf1dd8SToby Isaac ierr = PetscFree(Bf);CHKERRQ(ierr); 8720cf1dd8SToby Isaac } 8820cf1dd8SToby Isaac n = spdim; 8920cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s*spdim*spdim], &n, pivots, &info)); 9020cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s*spdim*spdim], &n, pivots, work, &n, &info)); 9120cf1dd8SToby Isaac } 9220cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 9320cf1dd8SToby Isaac for (s = 0; s <cmp->numSubelements*spdim*spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]); 9420cf1dd8SToby Isaac ierr = PetscFree(invVscalar);CHKERRQ(ierr); 9520cf1dd8SToby Isaac #endif 9620cf1dd8SToby Isaac ierr = PetscFree2(pivots,work);CHKERRQ(ierr); 9720cf1dd8SToby Isaac PetscFunctionReturn(0); 9820cf1dd8SToby Isaac } 9920cf1dd8SToby Isaac 100*2b99622eSMatthew G. Knepley static PetscErrorCode PetscFEGetTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H) 10120cf1dd8SToby Isaac { 10220cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 10320cf1dd8SToby Isaac DM dm; 10420cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 10520cf1dd8SToby Isaac PetscInt spdim; /* Dimension of subelement FE space P */ 10620cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 10720cf1dd8SToby Isaac PetscInt comp; /* Field components */ 10820cf1dd8SToby Isaac PetscInt *subpoints; 10920cf1dd8SToby Isaac PetscReal *tmpB, *tmpD, *tmpH, *subpoint; 11020cf1dd8SToby Isaac PetscInt p, s, d, e, j, k; 11120cf1dd8SToby Isaac PetscErrorCode ierr; 11220cf1dd8SToby Isaac 11320cf1dd8SToby Isaac PetscFunctionBegin; 11420cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 11520cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 11620cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(fem->basisSpace, &spdim);CHKERRQ(ierr); 11720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 11820cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr); 11920cf1dd8SToby Isaac /* Divide points into subelements */ 12020cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints);CHKERRQ(ierr); 12120cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr); 12220cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 12320cf1dd8SToby Isaac for (s = 0; s < cmp->numSubelements; ++s) { 12420cf1dd8SToby Isaac PetscBool inside; 12520cf1dd8SToby Isaac 12620cf1dd8SToby Isaac /* Apply transform, and check that point is inside cell */ 12720cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 12820cf1dd8SToby Isaac subpoint[d] = -1.0; 12920cf1dd8SToby 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]); 13020cf1dd8SToby Isaac } 13120cf1dd8SToby Isaac ierr = CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);CHKERRQ(ierr); 13220cf1dd8SToby Isaac if (inside) {subpoints[p] = s; break;} 13320cf1dd8SToby Isaac } 13420cf1dd8SToby Isaac if (s >= cmp->numSubelements) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d was not found in any subelement", p); 13520cf1dd8SToby Isaac } 13620cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint);CHKERRQ(ierr); 13720cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 13820cf1dd8SToby Isaac if (B) {ierr = DMGetWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 13920cf1dd8SToby Isaac if (D) {ierr = DMGetWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 14020cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 14120cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);CHKERRQ(ierr); 14220cf1dd8SToby Isaac /* Translate to the nodal basis */ 143580bdb30SBarry Smith if (B) {ierr = PetscArrayzero(B, npoints*pdim*comp);CHKERRQ(ierr);} 144580bdb30SBarry Smith if (D) {ierr = PetscArrayzero(D, npoints*pdim*comp*dim);CHKERRQ(ierr);} 145580bdb30SBarry Smith if (H) {ierr = PetscArrayzero(H, npoints*pdim*comp*dim*dim);CHKERRQ(ierr);} 14620cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 14720cf1dd8SToby Isaac const PetscInt s = subpoints[p]; 14820cf1dd8SToby Isaac 14920cf1dd8SToby Isaac if (B) { 15020cf1dd8SToby Isaac /* Multiply by V^{-1} (spdim x spdim) */ 15120cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 15220cf1dd8SToby Isaac const PetscInt i = (p*pdim + cmp->embedding[s*spdim+j])*comp; 15320cf1dd8SToby Isaac 15420cf1dd8SToby Isaac B[i] = 0.0; 15520cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 15620cf1dd8SToby Isaac B[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpB[p*spdim + k]; 15720cf1dd8SToby Isaac } 15820cf1dd8SToby Isaac } 15920cf1dd8SToby Isaac } 16020cf1dd8SToby Isaac if (D) { 16120cf1dd8SToby Isaac /* Multiply by V^{-1} (spdim x spdim) */ 16220cf1dd8SToby Isaac for (j = 0; j < spdim; ++j) { 16320cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 16420cf1dd8SToby Isaac const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim + d; 16520cf1dd8SToby Isaac 16620cf1dd8SToby Isaac D[i] = 0.0; 16720cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 16820cf1dd8SToby Isaac D[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpD[(p*spdim + k)*dim + d]; 16920cf1dd8SToby Isaac } 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; 18020cf1dd8SToby Isaac for (k = 0; k < spdim; ++k) { 18120cf1dd8SToby Isaac H[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpH[(p*spdim + k)*dim*dim + d]; 18220cf1dd8SToby Isaac } 18320cf1dd8SToby Isaac } 18420cf1dd8SToby Isaac } 18520cf1dd8SToby Isaac } 18620cf1dd8SToby Isaac } 18720cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints);CHKERRQ(ierr); 18820cf1dd8SToby Isaac if (B) {ierr = DMRestoreWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 18920cf1dd8SToby Isaac if (D) {ierr = DMRestoreWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 19020cf1dd8SToby Isaac if (H) {ierr = DMRestoreWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 19120cf1dd8SToby Isaac PetscFunctionReturn(0); 19220cf1dd8SToby Isaac } 19320cf1dd8SToby Isaac 194*2b99622eSMatthew G. Knepley static PetscErrorCode PetscFEInitialize_Composite(PetscFE fem) 19520cf1dd8SToby Isaac { 19620cf1dd8SToby Isaac PetscFunctionBegin; 19720cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 19820cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Composite; 19920cf1dd8SToby Isaac fem->ops->view = NULL; 20020cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Composite; 20120cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 20220cf1dd8SToby Isaac fem->ops->gettabulation = PetscFEGetTabulation_Composite; 20320cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 20420cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 20520cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 20620cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 20720cf1dd8SToby Isaac PetscFunctionReturn(0); 20820cf1dd8SToby Isaac } 20920cf1dd8SToby Isaac 21020cf1dd8SToby Isaac /*MC 21120cf1dd8SToby Isaac PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element 21220cf1dd8SToby Isaac 21320cf1dd8SToby Isaac Level: intermediate 21420cf1dd8SToby Isaac 21520cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 21620cf1dd8SToby Isaac M*/ 21720cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem) 21820cf1dd8SToby Isaac { 21920cf1dd8SToby Isaac PetscFE_Composite *cmp; 22020cf1dd8SToby Isaac PetscErrorCode ierr; 22120cf1dd8SToby Isaac 22220cf1dd8SToby Isaac PetscFunctionBegin; 22320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 22420cf1dd8SToby Isaac ierr = PetscNewLog(fem, &cmp);CHKERRQ(ierr); 22520cf1dd8SToby Isaac fem->data = cmp; 22620cf1dd8SToby Isaac 22720cf1dd8SToby Isaac cmp->cellRefiner = REFINER_NOOP; 22820cf1dd8SToby Isaac cmp->numSubelements = -1; 22920cf1dd8SToby Isaac cmp->v0 = NULL; 23020cf1dd8SToby Isaac cmp->jac = NULL; 23120cf1dd8SToby Isaac 23220cf1dd8SToby Isaac ierr = PetscFEInitialize_Composite(fem);CHKERRQ(ierr); 23320cf1dd8SToby Isaac PetscFunctionReturn(0); 23420cf1dd8SToby Isaac } 23520cf1dd8SToby Isaac 23620cf1dd8SToby Isaac /*@C 23720cf1dd8SToby Isaac PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement 23820cf1dd8SToby Isaac 23920cf1dd8SToby Isaac Not collective 24020cf1dd8SToby Isaac 24120cf1dd8SToby Isaac Input Parameter: 24220cf1dd8SToby Isaac . fem - The PetscFE object 24320cf1dd8SToby Isaac 24420cf1dd8SToby Isaac Output Parameters: 24520cf1dd8SToby Isaac + blockSize - The number of elements in a block 24620cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 24720cf1dd8SToby Isaac . batchSize - The number of elements in a batch 24820cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 24920cf1dd8SToby Isaac 25020cf1dd8SToby Isaac Level: intermediate 25120cf1dd8SToby Isaac 25220cf1dd8SToby Isaac .seealso: PetscFECreate() 25320cf1dd8SToby Isaac @*/ 25420cf1dd8SToby Isaac PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[]) 25520cf1dd8SToby Isaac { 25620cf1dd8SToby Isaac PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data; 25720cf1dd8SToby Isaac 25820cf1dd8SToby Isaac PetscFunctionBegin; 25920cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 26020cf1dd8SToby Isaac if (numSubelements) {PetscValidPointer(numSubelements, 2); *numSubelements = cmp->numSubelements;} 26120cf1dd8SToby Isaac if (v0) {PetscValidPointer(v0, 3); *v0 = cmp->v0;} 26220cf1dd8SToby Isaac if (jac) {PetscValidPointer(jac, 4); *jac = cmp->jac;} 26320cf1dd8SToby Isaac if (invjac) {PetscValidPointer(invjac, 5); *invjac = cmp->invjac;} 26420cf1dd8SToby Isaac PetscFunctionReturn(0); 26520cf1dd8SToby Isaac } 266