xref: /petsc/src/dm/tutorials/ex20.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 
2 static char help[] = "DMSwarm-PIC demonstator of inserting points into a cell DM \n\
3 Options: \n\
4 -mode {0,1} : 0 ==> DMDA, 1 ==> DMPLEX cell DM \n\
5 -dim {2,3}  : spatial dimension\n";
6 
7 #include <petsc.h>
8 #include <petscdm.h>
9 #include <petscdmda.h>
10 #include <petscdmplex.h>
11 #include <petscdmswarm.h>
12 
13 PetscErrorCode pic_insert_DMDA(PetscInt dim) {
14   DM        celldm = NULL, swarm;
15   PetscInt  dof, stencil_width;
16   PetscReal min[3], max[3];
17   PetscInt  ndir[3];
18 
19   PetscFunctionBegin;
20   /* Create the background cell DM */
21   dof           = 1;
22   stencil_width = 1;
23   if (dim == 2) { PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 25, 13, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, &celldm)); }
24   if (dim == 3) { PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 25, 13, 19, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, NULL, &celldm)); }
25 
26   PetscCall(DMDASetElementType(celldm, DMDA_ELEMENT_Q1));
27   PetscCall(DMSetFromOptions(celldm));
28   PetscCall(DMSetUp(celldm));
29 
30   PetscCall(DMDASetUniformCoordinates(celldm, 0.0, 2.0, 0.0, 1.0, 0.0, 1.5));
31 
32   /* Create the DMSwarm */
33   PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
34   PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm"));
35   PetscCall(DMSetType(swarm, DMSWARM));
36   PetscCall(DMSetDimension(swarm, dim));
37 
38   /* Configure swarm to be of type PIC */
39   PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
40   PetscCall(DMSwarmSetCellDM(swarm, celldm));
41 
42   /* Register two scalar fields within the DMSwarm */
43   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE));
44   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE));
45   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
46 
47   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
48   PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
49 
50   /* Insert swarm coordinates cell-wise */
51   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_REGULAR, 3));
52   min[0]  = 0.5;
53   max[0]  = 0.7;
54   min[1]  = 0.5;
55   max[1]  = 0.8;
56   min[2]  = 0.5;
57   max[2]  = 0.9;
58   ndir[0] = ndir[1] = ndir[2] = 30;
59   PetscCall(DMSwarmSetPointsUniformCoordinates(swarm, min, max, ndir, ADD_VALUES));
60 
61   /* This should be dispatched from a regular DMView() */
62   PetscCall(DMSwarmViewXDMF(swarm, "ex20.xmf"));
63   PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD));
64   PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD));
65 
66   {
67     PetscInt    npoints, *list;
68     PetscMPIInt rank;
69 
70     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
71     PetscCall(DMSwarmSortGetAccess(swarm));
72     PetscCall(DMSwarmSortGetNumberOfPointsPerCell(swarm, 0, &npoints));
73     PetscCall(DMSwarmSortGetPointsPerCell(swarm, rank, &npoints, &list));
74     PetscCall(PetscFree(list));
75     PetscCall(DMSwarmSortRestoreAccess(swarm));
76   }
77   PetscCall(DMSwarmMigrate(swarm, PETSC_FALSE));
78   PetscCall(DMDestroy(&celldm));
79   PetscCall(DMDestroy(&swarm));
80   PetscFunctionReturn(0);
81 }
82 
83 PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim) {
84   DM          celldm = NULL, swarm, distributedMesh = NULL;
85   const char *fieldnames[] = {"viscosity"};
86 
87   PetscFunctionBegin;
88   /* Create the background cell DM */
89   if (dim == 2) {
90     PetscInt   cells_per_dim[2], nx[2];
91     PetscInt   n_tricells;
92     PetscInt   n_trivert;
93     PetscInt  *tricells;
94     PetscReal *trivert, dx, dy;
95     PetscInt   ii, jj, cnt;
96 
97     cells_per_dim[0] = 4;
98     cells_per_dim[1] = 4;
99     n_tricells       = cells_per_dim[0] * cells_per_dim[1] * 2;
100     nx[0]            = cells_per_dim[0] + 1;
101     nx[1]            = cells_per_dim[1] + 1;
102     n_trivert        = nx[0] * nx[1];
103 
104     PetscCall(PetscMalloc1(n_tricells * 3, &tricells));
105     PetscCall(PetscMalloc1(nx[0] * nx[1] * 2, &trivert));
106 
107     /* verts */
108     cnt = 0;
109     dx  = 2.0 / ((PetscReal)cells_per_dim[0]);
110     dy  = 1.0 / ((PetscReal)cells_per_dim[1]);
111     for (jj = 0; jj < nx[1]; jj++) {
112       for (ii = 0; ii < nx[0]; ii++) {
113         trivert[2 * cnt + 0] = 0.0 + ii * dx;
114         trivert[2 * cnt + 1] = 0.0 + jj * dy;
115         cnt++;
116       }
117     }
118 
119     /* connectivity */
120     cnt = 0;
121     for (jj = 0; jj < cells_per_dim[1]; jj++) {
122       for (ii = 0; ii < cells_per_dim[0]; ii++) {
123         PetscInt idx, idx0, idx1, idx2, idx3;
124 
125         idx  = (ii) + (jj)*nx[0];
126         idx0 = idx;
127         idx1 = idx0 + 1;
128         idx2 = idx1 + nx[0];
129         idx3 = idx0 + nx[0];
130 
131         tricells[3 * cnt + 0] = idx0;
132         tricells[3 * cnt + 1] = idx1;
133         tricells[3 * cnt + 2] = idx2;
134         cnt++;
135 
136         tricells[3 * cnt + 0] = idx0;
137         tricells[3 * cnt + 1] = idx2;
138         tricells[3 * cnt + 2] = idx3;
139         cnt++;
140       }
141     }
142     PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, n_tricells, n_trivert, 3, PETSC_TRUE, tricells, dim, trivert, &celldm));
143     PetscCall(PetscFree(trivert));
144     PetscCall(PetscFree(tricells));
145   }
146   PetscCheck(dim != 3, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Only 2D PLEX example supported");
147 
148   /* Distribute mesh over processes */
149   PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh));
150   if (distributedMesh) {
151     PetscCall(DMDestroy(&celldm));
152     celldm = distributedMesh;
153   }
154   PetscCall(PetscObjectSetName((PetscObject)celldm, "Cells"));
155   PetscCall(DMSetFromOptions(celldm));
156   {
157     PetscInt     numComp[] = {1};
158     PetscInt     numDof[]  = {1, 0, 0}; /* vert, edge, cell */
159     PetscInt     numBC     = 0;
160     PetscSection section;
161 
162     PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &section));
163     PetscCall(DMSetLocalSection(celldm, section));
164     PetscCall(PetscSectionDestroy(&section));
165   }
166   PetscCall(DMSetUp(celldm));
167   {
168     PetscViewer viewer;
169 
170     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
171     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK));
172     PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
173     PetscCall(PetscViewerFileSetName(viewer, "ex20plex.vtk"));
174     PetscCall(DMView(celldm, viewer));
175     PetscCall(PetscViewerDestroy(&viewer));
176   }
177 
178   /* Create the DMSwarm */
179   PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
180   PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm"));
181   PetscCall(DMSetType(swarm, DMSWARM));
182   PetscCall(DMSetDimension(swarm, dim));
183 
184   PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
185   PetscCall(DMSwarmSetCellDM(swarm, celldm));
186 
187   /* Register two scalar fields within the DMSwarm */
188   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE));
189   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE));
190   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
191 
192   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
193   PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
194 
195   /* Insert swarm coordinates cell-wise */
196   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_SUBDIVISION, 2));
197   PetscCall(DMSwarmViewFieldsXDMF(swarm, "ex20.xmf", 1, fieldnames));
198   PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD));
199   PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD));
200   PetscCall(DMDestroy(&celldm));
201   PetscCall(DMDestroy(&swarm));
202   PetscFunctionReturn(0);
203 }
204 
205 PetscErrorCode pic_insert_DMPLEX(PetscBool is_simplex, PetscInt dim) {
206   DM          celldm, swarm, distributedMesh = NULL;
207   const char *fieldnames[] = {"viscosity", "DMSwarm_rank"};
208 
209   PetscFunctionBegin;
210 
211   /* Create the background cell DM */
212   {
213     PetscInt faces[3] = {4, 2, 4};
214     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm));
215   }
216 
217   /* Distribute mesh over processes */
218   PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh));
219   if (distributedMesh) {
220     PetscCall(DMDestroy(&celldm));
221     celldm = distributedMesh;
222   }
223   PetscCall(PetscObjectSetName((PetscObject)celldm, "Cells"));
224   PetscCall(DMSetFromOptions(celldm));
225   {
226     PetscInt     numComp[] = {1};
227     PetscInt     numDof[]  = {1, 0, 0}; /* vert, edge, cell */
228     PetscInt     numBC     = 0;
229     PetscSection section;
230 
231     PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &section));
232     PetscCall(DMSetLocalSection(celldm, section));
233     PetscCall(PetscSectionDestroy(&section));
234   }
235   PetscCall(DMSetUp(celldm));
236   {
237     PetscViewer viewer;
238 
239     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
240     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK));
241     PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
242     PetscCall(PetscViewerFileSetName(viewer, "ex20plex.vtk"));
243     PetscCall(DMView(celldm, viewer));
244     PetscCall(PetscViewerDestroy(&viewer));
245   }
246 
247   PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
248   PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm"));
249   PetscCall(DMSetType(swarm, DMSWARM));
250   PetscCall(DMSetDimension(swarm, dim));
251 
252   PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
253   PetscCall(DMSwarmSetCellDM(swarm, celldm));
254 
255   /* Register two scalar fields within the DMSwarm */
256   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE));
257   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE));
258   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
259 
260   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
261   PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
262 
263   /* Insert swarm coordinates cell-wise */
264   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_GAUSS, 3));
265   PetscCall(DMSwarmViewFieldsXDMF(swarm, "ex20.xmf", 2, fieldnames));
266   PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD));
267   PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD));
268   PetscCall(DMDestroy(&celldm));
269   PetscCall(DMDestroy(&swarm));
270   PetscFunctionReturn(0);
271 }
272 
273 int main(int argc, char **args) {
274   PetscInt mode = 0;
275   PetscInt dim  = 2;
276 
277   PetscFunctionBeginUser;
278   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
279   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mode", &mode, NULL));
280   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
281   switch (mode) {
282   case 0: PetscCall(pic_insert_DMDA(dim)); break;
283   case 1:
284     /* tri / tet */
285     PetscCall(pic_insert_DMPLEX(PETSC_TRUE, dim));
286     break;
287   case 2:
288     /* quad / hex */
289     PetscCall(pic_insert_DMPLEX(PETSC_FALSE, dim));
290     break;
291   default: PetscCall(pic_insert_DMDA(dim)); break;
292   }
293   PetscCall(PetscFinalize());
294   return 0;
295 }
296 
297 /*TEST
298 
299    test:
300       args:
301       requires: !complex double
302       filter: grep -v atomic
303       filter_output: grep -v atomic
304 
305    test:
306       suffix: 2
307       requires: triangle double !complex
308       args: -mode 1
309       filter: grep -v atomic
310       filter_output: grep -v atomic
311 
312 TEST*/
313