xref: /petsc/src/dm/impls/swarm/swarmpic_plex.c (revision 84467f862f2de26368b63ea79dccd665bcda30cd)
1 #include <petscdm.h>
2 #include <petscdmplex.h>
3 #include <petscdmswarm.h>
4 #include "../src/dm/impls/swarm/data_bucket.h"
5 
6 PetscBool  SwarmProjcite       = PETSC_FALSE;
7 const char SwarmProjCitation[] = "@article{PusztayKnepleyAdams2022,\n"
8                                  "title   = {Conservative Projection Between FEM and Particle Bases},\n"
9                                  "author  = {Joseph V. Pusztay and Matthew G. Knepley and Mark F. Adams},\n"
10                                  "journal = {SIAM Journal on Scientific Computing},\n"
11                                  "volume  = {44},\n"
12                                  "number  = {4},\n"
13                                  "pages   = {C310--C319},\n"
14                                  "doi     = {10.1137/21M145407},\n"
15                                  "year    = {2022}\n}\n";
16 
17 PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM, DM, PetscInt, PetscReal *xi);
18 
19 static PetscErrorCode private_PetscFECreateDefault_scalar_pk1(DM dm, PetscInt dim, PetscBool isSimplex, PetscInt qorder, PetscFE *fem)
20 {
21   const PetscInt  Nc = 1;
22   PetscQuadrature q, fq;
23   DM              K;
24   PetscSpace      P;
25   PetscDualSpace  Q;
26   PetscInt        order, quadPointsPerEdge;
27   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
28 
29   PetscFunctionBegin;
30   /* Create space */
31   PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)dm), &P));
32   /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); */
33   PetscCall(PetscSpacePolynomialSetTensor(P, tensor));
34   /* PetscCall(PetscSpaceSetFromOptions(P)); */
35   PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
36   PetscCall(PetscSpaceSetDegree(P, 1, PETSC_DETERMINE));
37   PetscCall(PetscSpaceSetNumComponents(P, Nc));
38   PetscCall(PetscSpaceSetNumVariables(P, dim));
39   PetscCall(PetscSpaceSetUp(P));
40   PetscCall(PetscSpaceGetDegree(P, &order, NULL));
41   PetscCall(PetscSpacePolynomialGetTensor(P, &tensor));
42   /* Create dual space */
43   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)dm), &Q));
44   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
45   /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); */
46   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
47   PetscCall(PetscDualSpaceSetDM(Q, K));
48   PetscCall(DMDestroy(&K));
49   PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
50   PetscCall(PetscDualSpaceSetOrder(Q, order));
51   PetscCall(PetscDualSpaceLagrangeSetTensor(Q, tensor));
52   /* PetscCall(PetscDualSpaceSetFromOptions(Q)); */
53   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
54   PetscCall(PetscDualSpaceSetUp(Q));
55   /* Create element */
56   PetscCall(PetscFECreate(PetscObjectComm((PetscObject)dm), fem));
57   /* PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); */
58   /* PetscCall(PetscFESetFromOptions(*fem)); */
59   PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
60   PetscCall(PetscFESetBasisSpace(*fem, P));
61   PetscCall(PetscFESetDualSpace(*fem, Q));
62   PetscCall(PetscFESetNumComponents(*fem, Nc));
63   PetscCall(PetscFESetUp(*fem));
64   PetscCall(PetscSpaceDestroy(&P));
65   PetscCall(PetscDualSpaceDestroy(&Q));
66   /* Create quadrature (with specified order if given) */
67   qorder            = qorder >= 0 ? qorder : order;
68   quadPointsPerEdge = PetscMax(qorder + 1, 1);
69   if (isSimplex) {
70     PetscCall(PetscDTStroudConicalQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q));
71     PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq));
72   } else {
73     PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q));
74     PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, &fq));
75   }
76   PetscCall(PetscFESetQuadrature(*fem, q));
77   PetscCall(PetscFESetFaceQuadrature(*fem, fq));
78   PetscCall(PetscQuadratureDestroy(&q));
79   PetscCall(PetscQuadratureDestroy(&fq));
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 
83 static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(DM dm, DM dmc, PetscInt nsub)
84 {
85   PetscInt         dim, nfaces, nbasis;
86   PetscInt         q, npoints_q, e, nel, pcnt, ps, pe, d, k, r;
87   PetscTabulation  T;
88   Vec              coorlocal;
89   PetscSection     coordSection;
90   PetscScalar     *elcoor = NULL;
91   PetscReal       *swarm_coor;
92   PetscInt        *swarm_cellid;
93   const PetscReal *xiq;
94   PetscQuadrature  quadrature;
95   PetscFE          fe, feRef;
96   PetscBool        is_simplex;
97   const char      *coordname;
98 
99   PetscFunctionBegin;
100   PetscCall(DMGetDimension(dmc, &dim));
101   is_simplex = PETSC_FALSE;
102   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
103   PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
104   if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
105 
106   PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
107 
108   for (r = 0; r < nsub; r++) {
109     PetscCall(PetscFERefine(fe, &feRef));
110     PetscCall(PetscFECopyQuadrature(feRef, fe));
111     PetscCall(PetscFEDestroy(&feRef));
112   }
113 
114   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
115   PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &xiq, NULL));
116   PetscCall(PetscFEGetDimension(fe, &nbasis));
117   PetscCall(PetscFEGetCellTabulation(fe, 1, &T));
118 
119   /* 0->cell, 1->edge, 2->vert */
120   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
121   nel = pe - ps;
122 
123   PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
124   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
125   PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
126   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
127 
128   PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
129   PetscCall(DMGetCoordinateSection(dmc, &coordSection));
130 
131   pcnt = 0;
132   for (e = 0; e < nel; e++) {
133     PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
134 
135     for (q = 0; q < npoints_q; q++) {
136       for (d = 0; d < dim; d++) {
137         swarm_coor[dim * pcnt + d] = 0.0;
138         for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][q * nbasis + k] * PetscRealPart(elcoor[dim * k + d]);
139       }
140       swarm_cellid[pcnt] = e;
141       pcnt++;
142     }
143     PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
144   }
145   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
146   PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
147 
148   PetscCall(PetscFEDestroy(&fe));
149   PetscFunctionReturn(PETSC_SUCCESS);
150 }
151 
152 static PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(DM dm, DM dmc, PetscInt npoints)
153 {
154   PetscInt     dim;
155   PetscInt     ii, jj, q, npoints_q, e, nel, npe, pcnt, ps, pe, d, k, nfaces;
156   PetscReal   *xi, ds, ds2;
157   PetscReal  **basis;
158   Vec          coorlocal;
159   PetscSection coordSection;
160   PetscScalar *elcoor = NULL;
161   PetscReal   *swarm_coor;
162   PetscInt    *swarm_cellid;
163   PetscBool    is_simplex;
164   const char  *coordname;
165 
166   PetscFunctionBegin;
167   PetscCall(DMGetDimension(dmc, &dim));
168   PetscCheck(dim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only 2D is supported");
169   is_simplex = PETSC_FALSE;
170   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
171   PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
172   if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
173   PetscCheck(is_simplex, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only the simplex is supported");
174 
175   PetscCall(PetscMalloc1(dim * npoints * npoints, &xi));
176   pcnt = 0;
177   ds   = 1.0 / (PetscReal)(npoints - 1);
178   ds2  = 1.0 / (PetscReal)npoints;
179   for (jj = 0; jj < npoints; jj++) {
180     for (ii = 0; ii < npoints - jj; ii++) {
181       xi[dim * pcnt + 0] = ii * ds;
182       xi[dim * pcnt + 1] = jj * ds;
183 
184       xi[dim * pcnt + 0] *= (1.0 - 1.2 * ds2);
185       xi[dim * pcnt + 1] *= (1.0 - 1.2 * ds2);
186 
187       xi[dim * pcnt + 0] += 0.35 * ds2;
188       xi[dim * pcnt + 1] += 0.35 * ds2;
189       pcnt++;
190     }
191   }
192   npoints_q = pcnt;
193 
194   npe = 3; /* nodes per element (triangle) */
195   PetscCall(PetscMalloc1(npoints_q, &basis));
196   for (q = 0; q < npoints_q; q++) {
197     PetscCall(PetscMalloc1(npe, &basis[q]));
198 
199     basis[q][0] = 1.0 - xi[dim * q + 0] - xi[dim * q + 1];
200     basis[q][1] = xi[dim * q + 0];
201     basis[q][2] = xi[dim * q + 1];
202   }
203 
204   /* 0->cell, 1->edge, 2->vert */
205   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
206   nel = pe - ps;
207 
208   PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1));
209   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
210   PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
211   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
212 
213   PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
214   PetscCall(DMGetCoordinateSection(dmc, &coordSection));
215 
216   pcnt = 0;
217   for (e = 0; e < nel; e++) {
218     PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor));
219 
220     for (q = 0; q < npoints_q; q++) {
221       for (d = 0; d < dim; d++) {
222         swarm_coor[dim * pcnt + d] = 0.0;
223         for (k = 0; k < npe; k++) swarm_coor[dim * pcnt + d] += basis[q][k] * PetscRealPart(elcoor[dim * k + d]);
224       }
225       swarm_cellid[pcnt] = e;
226       pcnt++;
227     }
228     PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, e, NULL, &elcoor));
229   }
230   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
231   PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
232 
233   PetscCall(PetscFree(xi));
234   for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q]));
235   PetscCall(PetscFree(basis));
236   PetscFunctionReturn(PETSC_SUCCESS);
237 }
238 
239 PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_PLEX(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param)
240 {
241   PetscInt dim;
242 
243   PetscFunctionBegin;
244   PetscCall(DMGetDimension(celldm, &dim));
245   switch (layout) {
246   case DMSWARMPIC_LAYOUT_REGULAR:
247     PetscCheck(dim != 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No 3D support for REGULAR+PLEX");
248     PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX2D_Regular(dm, celldm, layout_param));
249     break;
250   case DMSWARMPIC_LAYOUT_GAUSS: {
251     PetscQuadrature  quad, facequad;
252     const PetscReal *xi;
253     DMPolytopeType   ct;
254     PetscInt         cStart, Nq;
255 
256     PetscCall(DMPlexGetHeightStratum(celldm, 0, &cStart, NULL));
257     PetscCall(DMPlexGetCellType(celldm, cStart, &ct));
258     PetscCall(PetscDTCreateDefaultQuadrature(ct, layout_param, &quad, &facequad));
259     PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &xi, NULL));
260     PetscCall(private_DMSwarmSetPointCoordinatesCellwise_PLEX(dm, celldm, Nq, (PetscReal *)xi));
261     PetscCall(PetscQuadratureDestroy(&quad));
262     PetscCall(PetscQuadratureDestroy(&facequad));
263   } break;
264   case DMSWARMPIC_LAYOUT_SUBDIVISION:
265     PetscCall(private_DMSwarmInsertPointsUsingCellDM_PLEX_SubDivide(dm, celldm, layout_param));
266     break;
267   }
268   PetscFunctionReturn(PETSC_SUCCESS);
269 }
270 
271 PetscErrorCode private_DMSwarmSetPointCoordinatesCellwise_PLEX(DM dm, DM dmc, PetscInt npoints, PetscReal xi[])
272 {
273   PetscBool       is_simplex, is_tensorcell;
274   PetscInt        dim, nfaces, ps, pe, p, d, nbasis, pcnt, e, k, nel;
275   PetscFE         fe;
276   PetscQuadrature quadrature;
277   PetscTabulation T;
278   PetscReal      *xiq;
279   Vec             coorlocal;
280   PetscSection    coordSection;
281   PetscScalar    *elcoor = NULL;
282   PetscReal      *swarm_coor;
283   PetscInt       *swarm_cellid;
284   const char     *coordname;
285 
286   PetscFunctionBegin;
287   PetscCall(DMGetDimension(dmc, &dim));
288 
289   is_simplex    = PETSC_FALSE;
290   is_tensorcell = PETSC_FALSE;
291   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
292   PetscCall(DMPlexGetConeSize(dmc, ps, &nfaces));
293 
294   if (nfaces == (dim + 1)) is_simplex = PETSC_TRUE;
295 
296   switch (dim) {
297   case 2:
298     if (nfaces == 4) is_tensorcell = PETSC_TRUE;
299     break;
300   case 3:
301     if (nfaces == 6) is_tensorcell = PETSC_TRUE;
302     break;
303   default:
304     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for 2D, 3D");
305   }
306 
307   /* check points provided fail inside the reference cell */
308   if (is_simplex) {
309     for (p = 0; p < npoints; p++) {
310       PetscReal sum;
311       for (d = 0; d < dim; d++) PetscCheck(xi[dim * p + d] >= -1.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain");
312       sum = 0.0;
313       for (d = 0; d < dim; d++) sum += xi[dim * p + d];
314       PetscCheck(sum <= 0.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the simplex domain");
315     }
316   } else if (is_tensorcell) {
317     for (p = 0; p < npoints; p++) {
318       for (d = 0; d < dim; d++) PetscCheck(PetscAbsReal(xi[dim * p + d]) <= 1.0, PetscObjectComm((PetscObject)dm), PETSC_ERR_USER, "Points do not fail inside the tensor domain [-1,1]^d");
319     }
320   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only support for d-simplex and d-tensorcell");
321 
322   PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)dm), &quadrature));
323   PetscCall(PetscMalloc1(npoints * dim, &xiq));
324   PetscCall(PetscArraycpy(xiq, xi, npoints * dim));
325   PetscCall(PetscQuadratureSetData(quadrature, dim, 1, npoints, (const PetscReal *)xiq, NULL));
326   PetscCall(private_PetscFECreateDefault_scalar_pk1(dmc, dim, is_simplex, 0, &fe));
327   PetscCall(PetscFESetQuadrature(fe, quadrature));
328   PetscCall(PetscFEGetDimension(fe, &nbasis));
329   PetscCall(PetscFEGetCellTabulation(fe, 1, &T));
330 
331   /* for each cell, interpolate coordaintes and insert the interpolated points coordinates into swarm */
332   /* 0->cell, 1->edge, 2->vert */
333   PetscCall(DMPlexGetHeightStratum(dmc, 0, &ps, &pe));
334   nel = pe - ps;
335 
336   PetscCall(DMSwarmSetLocalSizes(dm, npoints * nel, -1));
337   PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
338   PetscCall(DMSwarmGetField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
339   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
340 
341   PetscCall(DMGetCoordinatesLocal(dmc, &coorlocal));
342   PetscCall(DMGetCoordinateSection(dmc, &coordSection));
343 
344   pcnt = 0;
345   for (e = 0; e < nel; e++) {
346     PetscCall(DMPlexVecGetClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
347 
348     for (p = 0; p < npoints; p++) {
349       for (d = 0; d < dim; d++) {
350         swarm_coor[dim * pcnt + d] = 0.0;
351         for (k = 0; k < nbasis; k++) swarm_coor[dim * pcnt + d] += T->T[0][p * nbasis + k] * PetscRealPart(elcoor[dim * k + d]);
352       }
353       swarm_cellid[pcnt] = e;
354       pcnt++;
355     }
356     PetscCall(DMPlexVecRestoreClosure(dmc, coordSection, coorlocal, ps + e, NULL, &elcoor));
357   }
358   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid));
359   PetscCall(DMSwarmRestoreField(dm, coordname, NULL, NULL, (void **)&swarm_coor));
360 
361   PetscCall(PetscQuadratureDestroy(&quadrature));
362   PetscCall(PetscFEDestroy(&fe));
363   PetscFunctionReturn(PETSC_SUCCESS);
364 }
365