Actual source code: fecomposite.c

petsc-master 2020-02-18
Report Typos and Errors
  1:  #include <petsc/private/petscfeimpl.h>
  2:  #include <petsc/private/dtimpl.h>
  3:  #include <petsc/private/dmpleximpl.h>
  4:  #include <petscblaslapack.h>

  6: static PetscErrorCode PetscFEDestroy_Composite(PetscFE fem)
  7: {
  8:   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
  9:   PetscErrorCode     ierr;

 12:   CellRefinerRestoreAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);
 13:   PetscFree(cmp->embedding);
 14:   PetscFree(cmp);
 15:   return(0);
 16: }

 18: static PetscErrorCode PetscFESetUp_Composite(PetscFE fem)
 19: {
 20:   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
 21:   DM                 K;
 22:   PetscReal         *subpoint;
 23:   PetscBLASInt      *pivots;
 24:   PetscBLASInt       n, info;
 25:   PetscScalar       *work, *invVscalar;
 26:   PetscInt           dim, pdim, spdim, j, s;
 27:   PetscErrorCode     ierr;

 30:   /* Get affine mapping from reference cell to each subcell */
 31:   PetscDualSpaceGetDM(fem->dualSpace, &K);
 32:   DMGetDimension(K, &dim);
 33:   DMPlexGetCellRefiner_Internal(K, &cmp->cellRefiner);
 34:   CellRefinerGetAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);
 35:   /* Determine dof embedding into subelements */
 36:   PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
 37:   PetscSpaceGetDimension(fem->basisSpace, &spdim);
 38:   PetscMalloc1(cmp->numSubelements*spdim,&cmp->embedding);
 39:   DMGetWorkArray(K, dim, MPIU_REAL, &subpoint);
 40:   for (s = 0; s < cmp->numSubelements; ++s) {
 41:     PetscInt sd = 0;

 43:     for (j = 0; j < pdim; ++j) {
 44:       PetscBool       inside;
 45:       PetscQuadrature f;
 46:       PetscInt        d, e;

 48:       PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);
 49:       /* Apply transform to first point, and check that point is inside subcell */
 50:       for (d = 0; d < dim; ++d) {
 51:         subpoint[d] = -1.0;
 52:         for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(f->points[e] - cmp->v0[s*dim+e]);
 53:       }
 54:       CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);
 55:       if (inside) {cmp->embedding[s*spdim+sd++] = j;}
 56:     }
 57:     if (sd != spdim) SETERRQ3(PetscObjectComm((PetscObject) fem), PETSC_ERR_PLIB, "Subelement %d has %d dual basis vectors != %d", s, sd, spdim);
 58:   }
 59:   DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint);
 60:   /* Construct the change of basis from prime basis to nodal basis for each subelement */
 61:   PetscMalloc1(cmp->numSubelements*spdim*spdim,&fem->invV);
 62:   PetscMalloc2(spdim,&pivots,spdim,&work);
 63: #if defined(PETSC_USE_COMPLEX)
 64:   PetscMalloc1(cmp->numSubelements*spdim*spdim,&invVscalar);
 65: #else
 66:   invVscalar = fem->invV;
 67: #endif
 68:   for (s = 0; s < cmp->numSubelements; ++s) {
 69:     for (j = 0; j < spdim; ++j) {
 70:       PetscReal       *Bf;
 71:       PetscQuadrature  f;
 72:       const PetscReal *points, *weights;
 73:       PetscInt         Nc, Nq, q, k;

 75:       PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s*spdim+j], &f);
 76:       PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);
 77:       PetscMalloc1(f->numPoints*spdim*Nc,&Bf);
 78:       PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);
 79:       for (k = 0; k < spdim; ++k) {
 80:         /* n_j \cdot \phi_k */
 81:         invVscalar[(s*spdim + j)*spdim+k] = 0.0;
 82:         for (q = 0; q < Nq; ++q) {
 83:           invVscalar[(s*spdim + j)*spdim+k] += Bf[q*spdim+k]*weights[q];
 84:         }
 85:       }
 86:       PetscFree(Bf);
 87:     }
 88:     n = spdim;
 89:     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s*spdim*spdim], &n, pivots, &info));
 90:     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s*spdim*spdim], &n, pivots, work, &n, &info));
 91:   }
 92: #if defined(PETSC_USE_COMPLEX)
 93:   for (s = 0; s <cmp->numSubelements*spdim*spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]);
 94:   PetscFree(invVscalar);
 95: #endif
 96:   PetscFree2(pivots,work);
 97:   return(0);
 98: }

100: static PetscErrorCode PetscFECreateTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
101: {
102:   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
103:   DM                 dm;
104:   PetscInt           pdim;  /* Dimension of FE space P */
105:   PetscInt           spdim; /* Dimension of subelement FE space P */
106:   PetscInt           dim;   /* Spatial dimension */
107:   PetscInt           comp;  /* Field components */
108:   PetscInt          *subpoints;
109:   PetscReal         *B = K >= 0 ? T->T[0] : NULL;
110:   PetscReal         *D = K >= 1 ? T->T[1] : NULL;
111:   PetscReal         *H = K >= 2 ? T->T[2] : NULL;
112:   PetscReal         *tmpB = NULL, *tmpD = NULL, *tmpH = NULL, *subpoint;
113:   PetscInt           p, s, d, e, j, k;
114:   PetscErrorCode     ierr;

117:   PetscDualSpaceGetDM(fem->dualSpace, &dm);
118:   DMGetDimension(dm, &dim);
119:   PetscSpaceGetDimension(fem->basisSpace, &spdim);
120:   PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
121:   PetscFEGetNumComponents(fem, &comp);
122:   /* Divide points into subelements */
123:   DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints);
124:   DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint);
125:   for (p = 0; p < npoints; ++p) {
126:     for (s = 0; s < cmp->numSubelements; ++s) {
127:       PetscBool inside;

129:       /* Apply transform, and check that point is inside cell */
130:       for (d = 0; d < dim; ++d) {
131:         subpoint[d] = -1.0;
132:         for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(points[p*dim+e] - cmp->v0[s*dim+e]);
133:       }
134:       CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);
135:       if (inside) {subpoints[p] = s; break;}
136:     }
137:     if (s >= cmp->numSubelements) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d was not found in any subelement", p);
138:   }
139:   DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint);
140:   /* Evaluate the prime basis functions at all points */
141:   if (K >= 0) {DMGetWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);}
142:   if (K >= 1) {DMGetWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);}
143:   if (K >= 2) {DMGetWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);}
144:   PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH);
145:   /* Translate to the nodal basis */
146:   if (K >= 0) {PetscArrayzero(B, npoints*pdim*comp);}
147:   if (K >= 1) {PetscArrayzero(D, npoints*pdim*comp*dim);}
148:   if (K >= 2) {PetscArrayzero(H, npoints*pdim*comp*dim*dim);}
149:   for (p = 0; p < npoints; ++p) {
150:     const PetscInt s = subpoints[p];

152:     if (B) {
153:       /* Multiply by V^{-1} (spdim x spdim) */
154:       for (j = 0; j < spdim; ++j) {
155:         const PetscInt i = (p*pdim + cmp->embedding[s*spdim+j])*comp;

157:         B[i] = 0.0;
158:         for (k = 0; k < spdim; ++k) {
159:           B[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpB[p*spdim + k];
160:         }
161:       }
162:     }
163:     if (D) {
164:       /* Multiply by V^{-1} (spdim x spdim) */
165:       for (j = 0; j < spdim; ++j) {
166:         for (d = 0; d < dim; ++d) {
167:           const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim + d;

169:           D[i] = 0.0;
170:           for (k = 0; k < spdim; ++k) {
171:             D[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpD[(p*spdim + k)*dim + d];
172:           }
173:         }
174:       }
175:     }
176:     if (H) {
177:       /* Multiply by V^{-1} (pdim x pdim) */
178:       for (j = 0; j < spdim; ++j) {
179:         for (d = 0; d < dim*dim; ++d) {
180:           const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim*dim + d;

182:           H[i] = 0.0;
183:           for (k = 0; k < spdim; ++k) {
184:             H[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpH[(p*spdim + k)*dim*dim + d];
185:           }
186:         }
187:       }
188:     }
189:   }
190:   DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints);
191:   if (K >= 0) {DMRestoreWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);}
192:   if (K >= 1) {DMRestoreWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);}
193:   if (K >= 2) {DMRestoreWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);}
194:   return(0);
195: }

197: static PetscErrorCode PetscFEInitialize_Composite(PetscFE fem)
198: {
200:   fem->ops->setfromoptions          = NULL;
201:   fem->ops->setup                   = PetscFESetUp_Composite;
202:   fem->ops->view                    = NULL;
203:   fem->ops->destroy                 = PetscFEDestroy_Composite;
204:   fem->ops->getdimension            = PetscFEGetDimension_Basic;
205:   fem->ops->createtabulation        = PetscFECreateTabulation_Composite;
206:   fem->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
207:   fem->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
208:   fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
209:   fem->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
210:   return(0);
211: }

213: /*MC
214:   PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element

216:   Level: intermediate

218: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
219: M*/
220: PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem)
221: {
222:   PetscFE_Composite *cmp;
223:   PetscErrorCode     ierr;

227:   PetscNewLog(fem, &cmp);
228:   fem->data = cmp;

230:   cmp->cellRefiner    = REFINER_NOOP;
231:   cmp->numSubelements = -1;
232:   cmp->v0             = NULL;
233:   cmp->jac            = NULL;

235:   PetscFEInitialize_Composite(fem);
236:   return(0);
237: }

239: /*@C
240:   PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement

242:   Not collective

244:   Input Parameter:
245: . fem - The PetscFE object

247:   Output Parameters:
248: + blockSize - The number of elements in a block
249: . numBlocks - The number of blocks in a batch
250: . batchSize - The number of elements in a batch
251: - numBatches - The number of batches in a chunk

253:   Level: intermediate

255: .seealso: PetscFECreate()
256: @*/
257: PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[])
258: {
259:   PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;

267:   return(0);
268: }