xref: /petsc/src/dm/tutorials/swarm_ex1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests DMSwarm\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown #include <petscdmswarm.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM, PetscErrorCode (*)(DM, void *, PetscInt *, PetscInt **), size_t, void *, PetscInt *);
9c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM, PetscInt *);
10c4762a1bSJed Brown 
11*9371c9d4SSatish Balay PetscErrorCode ex1_1(void) {
12c4762a1bSJed Brown   DM          dms;
13c4762a1bSJed Brown   Vec         x;
14c4762a1bSJed Brown   PetscMPIInt rank, size;
15c4762a1bSJed Brown   PetscInt    p;
16c4762a1bSJed Brown 
17362febeeSStefano Zampini   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2008401ef6SPierre Jolivet   PetscCheck(!(size > 1) || !(size != 4), PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must be run wuth 4 MPI ranks");
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
239566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
246a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
279566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
289566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 1, PETSC_REAL));
299566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
309566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
319566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   {
34c4762a1bSJed Brown     PetscReal *array;
359566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", NULL, NULL, (void **)&array));
36*9371c9d4SSatish Balay     for (p = 0; p < 5 + rank; p++) { array[p] = 11.1 + p * 0.1 + rank * 100.0; }
379566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", NULL, NULL, (void **)&array));
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   {
41c4762a1bSJed Brown     PetscReal *array;
429566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", NULL, NULL, (void **)&array));
43*9371c9d4SSatish Balay     for (p = 0; p < 5 + rank; p++) { array[p] = 2.0e-2 + p * 0.001 + rank * 1.0; }
449566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", NULL, NULL, (void **)&array));
45c4762a1bSJed Brown   }
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
489566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(dms, "strain"));
519566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dms, &x));
529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   {
55c4762a1bSJed Brown     PetscInt *rankval;
56c4762a1bSJed Brown     PetscInt  npoints[2], npoints_orig[2];
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
599566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
61c4762a1bSJed Brown     if ((rank == 0) && (size > 1)) {
62c4762a1bSJed Brown       rankval[0] = 1;
63c4762a1bSJed Brown       rankval[3] = 1;
64c4762a1bSJed Brown     }
65*9371c9d4SSatish Balay     if (rank == 3) { rankval[2] = 1; }
669566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
679566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
689566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
699566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
7163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ") after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1], npoints[0], npoints[1]));
729566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown   {
769566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
779566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
789566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown   {
819566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "strain", &x));
829566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
839566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "strain", &x));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
87c4762a1bSJed Brown   PetscFunctionReturn(0);
88c4762a1bSJed Brown }
89c4762a1bSJed Brown 
90*9371c9d4SSatish Balay PetscErrorCode ex1_2(void) {
91c4762a1bSJed Brown   DM          dms;
92c4762a1bSJed Brown   Vec         x;
93c4762a1bSJed Brown   PetscMPIInt rank, size;
94c4762a1bSJed Brown   PetscInt    p;
95c4762a1bSJed Brown 
96362febeeSStefano Zampini   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
999566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1009566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1016a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
1029566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
1059566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 1, PETSC_REAL));
1069566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
1079566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
1089566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
109c4762a1bSJed Brown   {
110c4762a1bSJed Brown     PetscReal *array;
1119566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", NULL, NULL, (void **)&array));
112*9371c9d4SSatish Balay     for (p = 0; p < 5 + rank; p++) { array[p] = 11.1 + p * 0.1 + rank * 100.0; }
1139566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", NULL, NULL, (void **)&array));
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown   {
116c4762a1bSJed Brown     PetscReal *array;
1179566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", NULL, NULL, (void **)&array));
118*9371c9d4SSatish Balay     for (p = 0; p < 5 + rank; p++) { array[p] = 2.0e-2 + p * 0.001 + rank * 1.0; }
1199566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", NULL, NULL, (void **)&array));
120c4762a1bSJed Brown   }
121c4762a1bSJed Brown   {
122c4762a1bSJed Brown     PetscInt *rankval;
123c4762a1bSJed Brown     PetscInt  npoints[2], npoints_orig[2];
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
1269566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
1279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
12863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
1299566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
133c4762a1bSJed Brown 
134*9371c9d4SSatish Balay     if (rank == 1) { rankval[0] = -1; }
135*9371c9d4SSatish Balay     if (rank == 2) { rankval[1] = -1; }
136c4762a1bSJed Brown     if (rank == 3) {
137c4762a1bSJed Brown       rankval[3] = -1;
138c4762a1bSJed Brown       rankval[4] = -1;
139c4762a1bSJed Brown     }
1409566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
1419566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewCreate(dms));
1429566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
1439566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
1449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
14563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
1469566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
1509566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1519566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewDestroy(dms));
1549566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
1559566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
1569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
15763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after_v(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
1589566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
1629566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1639566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
164c4762a1bSJed Brown   }
1659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown /*
170c4762a1bSJed Brown  splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp"
171c4762a1bSJed Brown */
172*9371c9d4SSatish Balay PetscErrorCode ex1_3(void) {
173c4762a1bSJed Brown   DM          dms;
174c4762a1bSJed Brown   PetscMPIInt rank, size;
175c4762a1bSJed Brown   PetscInt    is, js, ni, nj, overlap;
176c4762a1bSJed Brown   DM          dmcell;
177c4762a1bSJed Brown 
178362febeeSStefano Zampini   PetscFunctionBegin;
1799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
181c4762a1bSJed Brown   overlap = 2;
1829566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 13, 13, PETSC_DECIDE, PETSC_DECIDE, 1, overlap, NULL, NULL, &dmcell));
1839566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
1849566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
1859566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
1869566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell, &is, &js, NULL, &ni, &nj, NULL));
1879566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1889566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1896a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
1909566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(dms, dmcell));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* load in data types */
1939566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
1949566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
1959566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coorx", 1, PETSC_REAL));
1969566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coory", 1, PETSC_REAL));
1979566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
1989566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, ni * nj * 4, 4));
1999566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* set values within the swarm */
202c4762a1bSJed Brown   {
203c4762a1bSJed Brown     PetscReal   *array_x, *array_y;
204c4762a1bSJed Brown     PetscInt     npoints, i, j, cnt;
205c4762a1bSJed Brown     DMDACoor2d **LA_coor;
206c4762a1bSJed Brown     Vec          coor;
207c4762a1bSJed Brown     DM           dmcellcdm;
208c4762a1bSJed Brown 
2099566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell, &dmcellcdm));
2109566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell, &coor));
2119566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm, coor, &LA_coor));
2129566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
2139566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
2149566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
215c4762a1bSJed Brown     cnt = 0;
216c4762a1bSJed Brown     for (j = js; j < js + nj; j++) {
217c4762a1bSJed Brown       for (i = is; i < is + ni; i++) {
218c4762a1bSJed Brown         PetscReal xp, yp;
219c4762a1bSJed Brown         xp                   = PetscRealPart(LA_coor[j][i].x);
220c4762a1bSJed Brown         yp                   = PetscRealPart(LA_coor[j][i].y);
221*9371c9d4SSatish Balay         array_x[4 * cnt + 0] = xp - 0.05;
222*9371c9d4SSatish Balay         if (array_x[4 * cnt + 0] < -1.0) { array_x[4 * cnt + 0] = -1.0 + 1.0e-12; }
223*9371c9d4SSatish Balay         array_x[4 * cnt + 1] = xp + 0.05;
224*9371c9d4SSatish Balay         if (array_x[4 * cnt + 1] > 1.0) { array_x[4 * cnt + 1] = 1.0 - 1.0e-12; }
225*9371c9d4SSatish Balay         array_x[4 * cnt + 2] = xp - 0.05;
226*9371c9d4SSatish Balay         if (array_x[4 * cnt + 2] < -1.0) { array_x[4 * cnt + 2] = -1.0 + 1.0e-12; }
227*9371c9d4SSatish Balay         array_x[4 * cnt + 3] = xp + 0.05;
228*9371c9d4SSatish Balay         if (array_x[4 * cnt + 3] > 1.0) { array_x[4 * cnt + 3] = 1.0 - 1.0e-12; }
229c4762a1bSJed Brown 
230*9371c9d4SSatish Balay         array_y[4 * cnt + 0] = yp - 0.05;
231*9371c9d4SSatish Balay         if (array_y[4 * cnt + 0] < -1.0) { array_y[4 * cnt + 0] = -1.0 + 1.0e-12; }
232*9371c9d4SSatish Balay         array_y[4 * cnt + 1] = yp - 0.05;
233*9371c9d4SSatish Balay         if (array_y[4 * cnt + 1] < -1.0) { array_y[4 * cnt + 1] = -1.0 + 1.0e-12; }
234*9371c9d4SSatish Balay         array_y[4 * cnt + 2] = yp + 0.05;
235*9371c9d4SSatish Balay         if (array_y[4 * cnt + 2] > 1.0) { array_y[4 * cnt + 2] = 1.0 - 1.0e-12; }
236*9371c9d4SSatish Balay         array_y[4 * cnt + 3] = yp + 0.05;
237*9371c9d4SSatish Balay         if (array_y[4 * cnt + 3] > 1.0) { array_y[4 * cnt + 3] = 1.0 - 1.0e-12; }
238c4762a1bSJed Brown         cnt++;
239c4762a1bSJed Brown       }
240c4762a1bSJed Brown     }
2419566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
2429566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
2439566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm, coor, &LA_coor));
244c4762a1bSJed Brown   }
245c4762a1bSJed Brown   {
246c4762a1bSJed Brown     PetscInt npoints[2], npoints_orig[2], ng;
247c4762a1bSJed Brown 
2489566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
2499566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
25163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
2549566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_DMDABoundingBox(dms, &ng));
255c4762a1bSJed Brown 
2569566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
2579566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
2589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
25963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
2609566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
262c4762a1bSJed Brown   }
263c4762a1bSJed Brown   {
264c4762a1bSJed Brown     PetscReal *array_x, *array_y;
265c4762a1bSJed Brown     PetscInt   npoints, p;
266c4762a1bSJed Brown     FILE      *fp = NULL;
267c4762a1bSJed Brown     char       name[PETSC_MAX_PATH_LEN];
268c4762a1bSJed Brown 
2699566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "c-rank%d.gp", rank));
270c4762a1bSJed Brown     fp = fopen(name, "w");
27128b400f6SJacob Faibussowitsch     PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
2729566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
2739566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
2749566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
275*9371c9d4SSatish Balay     for (p = 0; p < npoints; p++) { fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array_x[p], array_y[p], (double)rank); }
2769566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
2779566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
278c4762a1bSJed Brown     fclose(fp);
279c4762a1bSJed Brown   }
2809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
2819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
282c4762a1bSJed Brown   PetscFunctionReturn(0);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown typedef struct {
286c4762a1bSJed Brown   PetscReal cx[2];
287c4762a1bSJed Brown   PetscReal radius;
288c4762a1bSJed Brown } CollectZoneCtx;
289c4762a1bSJed Brown 
290*9371c9d4SSatish Balay PetscErrorCode collect_zone(DM dm, void *ctx, PetscInt *nfound, PetscInt **foundlist) {
291c4762a1bSJed Brown   CollectZoneCtx *zone = (CollectZoneCtx *)ctx;
292c4762a1bSJed Brown   PetscInt        p, npoints;
293c4762a1bSJed Brown   PetscReal      *array_x, *array_y, r2;
294c4762a1bSJed Brown   PetscInt        p2collect, *plist;
295c4762a1bSJed Brown   PetscMPIInt     rank;
296c4762a1bSJed Brown 
297362febeeSStefano Zampini   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
2999566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &npoints));
3009566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
3019566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   r2        = zone->radius * zone->radius;
304c4762a1bSJed Brown   p2collect = 0;
305c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
306c4762a1bSJed Brown     PetscReal sep2;
307c4762a1bSJed Brown 
308c4762a1bSJed Brown     sep2 = (array_x[p] - zone->cx[0]) * (array_x[p] - zone->cx[0]);
309c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1]) * (array_y[p] - zone->cx[1]);
310*9371c9d4SSatish Balay     if (sep2 < r2) { p2collect++; }
311c4762a1bSJed Brown   }
312c4762a1bSJed Brown 
3139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(p2collect + 1, &plist));
314c4762a1bSJed Brown   p2collect = 0;
315c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
316c4762a1bSJed Brown     PetscReal sep2;
317c4762a1bSJed Brown 
318c4762a1bSJed Brown     sep2 = (array_x[p] - zone->cx[0]) * (array_x[p] - zone->cx[0]);
319c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1]) * (array_y[p] - zone->cx[1]);
320c4762a1bSJed Brown     if (sep2 < r2) {
321c4762a1bSJed Brown       plist[p2collect] = p;
322c4762a1bSJed Brown       p2collect++;
323c4762a1bSJed Brown     }
324c4762a1bSJed Brown   }
3259566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
3269566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   *nfound    = p2collect;
329c4762a1bSJed Brown   *foundlist = plist;
330c4762a1bSJed Brown   PetscFunctionReturn(0);
331c4762a1bSJed Brown }
332c4762a1bSJed Brown 
333*9371c9d4SSatish Balay PetscErrorCode ex1_4(void) {
334c4762a1bSJed Brown   DM              dms;
335c4762a1bSJed Brown   PetscMPIInt     rank, size;
336c4762a1bSJed Brown   PetscInt        is, js, ni, nj, overlap, nn;
337c4762a1bSJed Brown   DM              dmcell;
338c4762a1bSJed Brown   CollectZoneCtx *zone;
339c4762a1bSJed Brown   PetscReal       dx;
340c4762a1bSJed Brown 
341362febeeSStefano Zampini   PetscFunctionBegin;
3429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
344c4762a1bSJed Brown   nn      = 101;
345c4762a1bSJed Brown   dx      = 2.0 / (PetscReal)(nn - 1);
346c4762a1bSJed Brown   overlap = 0;
3479566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, nn, nn, PETSC_DECIDE, PETSC_DECIDE, 1, overlap, NULL, NULL, &dmcell));
3489566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
3499566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
3509566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
3519566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell, &is, &js, NULL, &ni, &nj, NULL));
3529566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
3539566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
3546a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   /* load in data types */
3579566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
3589566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
3599566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coorx", 1, PETSC_REAL));
3609566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coory", 1, PETSC_REAL));
3619566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
3629566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, ni * nj * 4, 4));
3639566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   /* set values within the swarm */
366c4762a1bSJed Brown   {
367c4762a1bSJed Brown     PetscReal   *array_x, *array_y;
368c4762a1bSJed Brown     PetscInt     npoints, i, j, cnt;
369c4762a1bSJed Brown     DMDACoor2d **LA_coor;
370c4762a1bSJed Brown     Vec          coor;
371c4762a1bSJed Brown     DM           dmcellcdm;
372c4762a1bSJed Brown 
3739566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell, &dmcellcdm));
3749566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell, &coor));
3759566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm, coor, &LA_coor));
3769566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
3779566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
3789566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
379c4762a1bSJed Brown     cnt = 0;
380c4762a1bSJed Brown     for (j = js; j < js + nj; j++) {
381c4762a1bSJed Brown       for (i = is; i < is + ni; i++) {
382c4762a1bSJed Brown         PetscReal xp, yp;
383c4762a1bSJed Brown 
384c4762a1bSJed Brown         xp                   = PetscRealPart(LA_coor[j][i].x);
385c4762a1bSJed Brown         yp                   = PetscRealPart(LA_coor[j][i].y);
386c4762a1bSJed Brown         array_x[4 * cnt + 0] = xp - dx * 0.1; /*if (array_x[4*cnt+0] < -1.0) { array_x[4*cnt+0] = -1.0+1.0e-12; }*/
387c4762a1bSJed Brown         array_x[4 * cnt + 1] = xp + dx * 0.1; /*if (array_x[4*cnt+1] > 1.0)  { array_x[4*cnt+1] =  1.0-1.0e-12; }*/
388c4762a1bSJed Brown         array_x[4 * cnt + 2] = xp - dx * 0.1; /*if (array_x[4*cnt+2] < -1.0) { array_x[4*cnt+2] = -1.0+1.0e-12; }*/
389c4762a1bSJed Brown         array_x[4 * cnt + 3] = xp + dx * 0.1; /*if (array_x[4*cnt+3] > 1.0)  { array_x[4*cnt+3] =  1.0-1.0e-12; }*/
390c4762a1bSJed Brown         array_y[4 * cnt + 0] = yp - dx * 0.1; /*if (array_y[4*cnt+0] < -1.0) { array_y[4*cnt+0] = -1.0+1.0e-12; }*/
391c4762a1bSJed Brown         array_y[4 * cnt + 1] = yp - dx * 0.1; /*if (array_y[4*cnt+1] < -1.0) { array_y[4*cnt+1] = -1.0+1.0e-12; }*/
392c4762a1bSJed Brown         array_y[4 * cnt + 2] = yp + dx * 0.1; /*if (array_y[4*cnt+2] > 1.0)  { array_y[4*cnt+2] =  1.0-1.0e-12; }*/
393c4762a1bSJed Brown         array_y[4 * cnt + 3] = yp + dx * 0.1; /*if (array_y[4*cnt+3] > 1.0)  { array_y[4*cnt+3] =  1.0-1.0e-12; }*/
394c4762a1bSJed Brown         cnt++;
395c4762a1bSJed Brown       }
396c4762a1bSJed Brown     }
3979566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
3989566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
3999566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm, coor, &LA_coor));
400c4762a1bSJed Brown   }
4019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &zone));
402c4762a1bSJed Brown   if (size == 4) {
403c4762a1bSJed Brown     if (rank == 0) {
404c4762a1bSJed Brown       zone->cx[0]  = 0.5;
405c4762a1bSJed Brown       zone->cx[1]  = 0.5;
406c4762a1bSJed Brown       zone->radius = 0.3;
407c4762a1bSJed Brown     }
408c4762a1bSJed Brown     if (rank == 1) {
409c4762a1bSJed Brown       zone->cx[0]  = -0.5;
410c4762a1bSJed Brown       zone->cx[1]  = 0.5;
411c4762a1bSJed Brown       zone->radius = 0.25;
412c4762a1bSJed Brown     }
413c4762a1bSJed Brown     if (rank == 2) {
414c4762a1bSJed Brown       zone->cx[0]  = 0.5;
415c4762a1bSJed Brown       zone->cx[1]  = -0.5;
416c4762a1bSJed Brown       zone->radius = 0.2;
417c4762a1bSJed Brown     }
418c4762a1bSJed Brown     if (rank == 3) {
419c4762a1bSJed Brown       zone->cx[0]  = -0.5;
420c4762a1bSJed Brown       zone->cx[1]  = -0.5;
421c4762a1bSJed Brown       zone->radius = 0.1;
422c4762a1bSJed Brown     }
423c4762a1bSJed Brown   } else {
424c4762a1bSJed Brown     if (rank == 0) {
425c4762a1bSJed Brown       zone->cx[0]  = 0.5;
426c4762a1bSJed Brown       zone->cx[1]  = 0.5;
427c4762a1bSJed Brown       zone->radius = 0.8;
428c4762a1bSJed Brown     } else {
429c4762a1bSJed Brown       zone->cx[0]  = 10.0;
430c4762a1bSJed Brown       zone->cx[1]  = 10.0;
431c4762a1bSJed Brown       zone->radius = 0.0;
432c4762a1bSJed Brown     }
433c4762a1bSJed Brown   }
434c4762a1bSJed Brown   {
435c4762a1bSJed Brown     PetscInt npoints[2], npoints_orig[2], ng;
436c4762a1bSJed Brown 
4379566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
4389566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
4399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
44063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
4419566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
4439566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_General(dms, collect_zone, sizeof(CollectZoneCtx), zone, &ng));
4449566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
4459566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
4469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
44763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
4489566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
450c4762a1bSJed Brown   }
451c4762a1bSJed Brown   {
452c4762a1bSJed Brown     PetscReal *array_x, *array_y;
453c4762a1bSJed Brown     PetscInt   npoints, p;
454c4762a1bSJed Brown     FILE      *fp = NULL;
455c4762a1bSJed Brown     char       name[PETSC_MAX_PATH_LEN];
456c4762a1bSJed Brown 
4579566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "c-rank%d.gp", rank));
458c4762a1bSJed Brown     fp = fopen(name, "w");
45928b400f6SJacob Faibussowitsch     PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
4609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
4619566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
4629566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
463*9371c9d4SSatish Balay     for (p = 0; p < npoints; p++) { fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array_x[p], array_y[p], (double)rank); }
4649566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
4659566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
466c4762a1bSJed Brown     fclose(fp);
467c4762a1bSJed Brown   }
4689566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
4699566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
4709566063dSJacob Faibussowitsch   PetscCall(PetscFree(zone));
471c4762a1bSJed Brown   PetscFunctionReturn(0);
472c4762a1bSJed Brown }
473c4762a1bSJed Brown 
474*9371c9d4SSatish Balay int main(int argc, char **argv) {
475c4762a1bSJed Brown   PetscInt test_mode = 4;
476c4762a1bSJed Brown 
477327415f7SBarry Smith   PetscFunctionBeginUser;
4789566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
4799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_mode", &test_mode, NULL));
480c4762a1bSJed Brown   if (test_mode == 1) {
4819566063dSJacob Faibussowitsch     PetscCall(ex1_1());
482c4762a1bSJed Brown   } else if (test_mode == 2) {
4839566063dSJacob Faibussowitsch     PetscCall(ex1_2());
484c4762a1bSJed Brown   } else if (test_mode == 3) {
4859566063dSJacob Faibussowitsch     PetscCall(ex1_3());
486c4762a1bSJed Brown   } else if (test_mode == 4) {
4879566063dSJacob Faibussowitsch     PetscCall(ex1_4());
488c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown test_mode value, should be 1,2,3,4");
4899566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
490b122ec5aSJacob Faibussowitsch   return 0;
491c4762a1bSJed Brown }
492c4762a1bSJed Brown 
493c4762a1bSJed Brown /*TEST
494c4762a1bSJed Brown 
495c4762a1bSJed Brown    build:
496c4762a1bSJed Brown       requires: !complex double
497c4762a1bSJed Brown 
498c4762a1bSJed Brown    test:
499c4762a1bSJed Brown       args: -test_mode 1
500c4762a1bSJed Brown       filter: grep -v atomic
501c4762a1bSJed Brown       filter_output: grep -v atomic
502c4762a1bSJed Brown 
503c4762a1bSJed Brown    test:
504c4762a1bSJed Brown       suffix: 2
505c4762a1bSJed Brown       args: -test_mode 2
506c4762a1bSJed Brown       filter: grep -v atomic
507c4762a1bSJed Brown       filter_output: grep -v atomic
508c4762a1bSJed Brown 
509c4762a1bSJed Brown    test:
510c4762a1bSJed Brown       suffix: 3
511c4762a1bSJed Brown       args: -test_mode 3
512c4762a1bSJed Brown       filter: grep -v atomic
513c4762a1bSJed Brown       filter_output: grep -v atomic
514c4762a1bSJed Brown       TODO: broken
515c4762a1bSJed Brown 
516c4762a1bSJed Brown    test:
517c4762a1bSJed Brown       suffix: 4
518c4762a1bSJed Brown       args: -test_mode 4
519c4762a1bSJed Brown       filter: grep -v atomic
520c4762a1bSJed Brown       filter_output: grep -v atomic
521c4762a1bSJed Brown 
522c4762a1bSJed Brown    test:
523c4762a1bSJed Brown       suffix: 5
524c4762a1bSJed Brown       nsize: 4
525c4762a1bSJed Brown       args: -test_mode 1
526c4762a1bSJed Brown       filter: grep -v atomic
527c4762a1bSJed Brown       filter_output: grep -v atomic
528c4762a1bSJed Brown 
529c4762a1bSJed Brown    test:
530c4762a1bSJed Brown       suffix: 6
531c4762a1bSJed Brown       nsize: 4
532c4762a1bSJed Brown       args: -test_mode 2
533c4762a1bSJed Brown       filter: grep -v atomic
534c4762a1bSJed Brown       filter_output: grep -v atomic
535c4762a1bSJed Brown 
536c4762a1bSJed Brown    test:
537c4762a1bSJed Brown       suffix: 7
538c4762a1bSJed Brown       nsize: 4
539c4762a1bSJed Brown       args: -test_mode 3
540c4762a1bSJed Brown       filter: grep -v atomic
541c4762a1bSJed Brown       filter_output: grep -v atomic
542c4762a1bSJed Brown       TODO: broken
543c4762a1bSJed Brown 
544c4762a1bSJed Brown    test:
545c4762a1bSJed Brown       suffix: 8
546c4762a1bSJed Brown       nsize: 4
547c4762a1bSJed Brown       args: -test_mode 4
548c4762a1bSJed Brown       filter: grep -v atomic
549c4762a1bSJed Brown       filter_output: grep -v atomic
550c4762a1bSJed Brown 
551c4762a1bSJed Brown TEST*/
552