xref: /petsc/src/dm/tutorials/swarm_ex1.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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 
11d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_1(void)
12d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown   DM          dms;
14c4762a1bSJed Brown   Vec         x;
15c4762a1bSJed Brown   PetscMPIInt rank, size;
16c4762a1bSJed Brown   PetscInt    p;
17c4762a1bSJed Brown 
18362febeeSStefano Zampini   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
21*d5b43468SJose E. Roman   PetscCheck(!(size > 1) || !(size != 4), PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must be run with 4 MPI ranks");
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
249566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
256a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
289566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
299566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 1, PETSC_REAL));
309566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
319566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
329566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   {
35c4762a1bSJed Brown     PetscReal *array;
369566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", NULL, NULL, (void **)&array));
37ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 11.1 + p * 0.1 + rank * 100.0;
389566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", NULL, NULL, (void **)&array));
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   {
42c4762a1bSJed Brown     PetscReal *array;
439566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", NULL, NULL, (void **)&array));
44ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 2.0e-2 + p * 0.001 + rank * 1.0;
459566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", NULL, NULL, (void **)&array));
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
499566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(dms, "strain"));
529566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dms, &x));
539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   {
56c4762a1bSJed Brown     PetscInt *rankval;
57c4762a1bSJed Brown     PetscInt  npoints[2], npoints_orig[2];
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
619566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
62c4762a1bSJed Brown     if ((rank == 0) && (size > 1)) {
63c4762a1bSJed Brown       rankval[0] = 1;
64c4762a1bSJed Brown       rankval[3] = 1;
65c4762a1bSJed Brown     }
66ad540459SPierre Jolivet     if (rank == 3) rankval[2] = 1;
679566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
689566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
699566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
709566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
7263a3b9bcSJacob 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]));
739566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown   {
779566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
789566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
799566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown   {
829566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "strain", &x));
839566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
849566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "strain", &x));
85c4762a1bSJed Brown   }
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
88c4762a1bSJed Brown   PetscFunctionReturn(0);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_2(void)
92d71ae5a4SJacob Faibussowitsch {
93c4762a1bSJed Brown   DM          dms;
94c4762a1bSJed Brown   Vec         x;
95c4762a1bSJed Brown   PetscMPIInt rank, size;
96c4762a1bSJed Brown   PetscInt    p;
97c4762a1bSJed Brown 
98362febeeSStefano Zampini   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1019566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1029566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1036a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
1049566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
1079566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 1, PETSC_REAL));
1089566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
1099566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
1109566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
111c4762a1bSJed Brown   {
112c4762a1bSJed Brown     PetscReal *array;
1139566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", NULL, NULL, (void **)&array));
114ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 11.1 + p * 0.1 + rank * 100.0;
1159566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", NULL, NULL, (void **)&array));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown   {
118c4762a1bSJed Brown     PetscReal *array;
1199566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", NULL, NULL, (void **)&array));
120ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 2.0e-2 + p * 0.001 + rank * 1.0;
1219566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", NULL, NULL, (void **)&array));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown   {
124c4762a1bSJed Brown     PetscInt *rankval;
125c4762a1bSJed Brown     PetscInt  npoints[2], npoints_orig[2];
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
1289566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
1299566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
13063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
1319566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
135c4762a1bSJed Brown 
136ad540459SPierre Jolivet     if (rank == 1) rankval[0] = -1;
137ad540459SPierre Jolivet     if (rank == 2) rankval[1] = -1;
138c4762a1bSJed Brown     if (rank == 3) {
139c4762a1bSJed Brown       rankval[3] = -1;
140c4762a1bSJed Brown       rankval[4] = -1;
141c4762a1bSJed Brown     }
1429566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
1439566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewCreate(dms));
1449566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
1459566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
1469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
14763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
1489566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
1529566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1539566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
154c4762a1bSJed Brown 
1559566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewDestroy(dms));
1569566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
1579566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
1589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
15963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after_v(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
1609566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
162c4762a1bSJed Brown 
1639566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
1649566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1659566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
166c4762a1bSJed Brown   }
1679566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
168c4762a1bSJed Brown   PetscFunctionReturn(0);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown /*
172c4762a1bSJed Brown  splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp"
173c4762a1bSJed Brown */
174d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_3(void)
175d71ae5a4SJacob Faibussowitsch {
176c4762a1bSJed Brown   DM          dms;
177c4762a1bSJed Brown   PetscMPIInt rank, size;
178c4762a1bSJed Brown   PetscInt    is, js, ni, nj, overlap;
179c4762a1bSJed Brown   DM          dmcell;
180c4762a1bSJed Brown 
181362febeeSStefano Zampini   PetscFunctionBegin;
1829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
184c4762a1bSJed Brown   overlap = 2;
1859566063dSJacob 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));
1869566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
1879566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
1889566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
1899566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell, &is, &js, NULL, &ni, &nj, NULL));
1909566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1919566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1926a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
1939566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(dms, dmcell));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* load in data types */
1969566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
1979566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
1989566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coorx", 1, PETSC_REAL));
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coory", 1, PETSC_REAL));
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
2019566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, ni * nj * 4, 4));
2029566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* set values within the swarm */
205c4762a1bSJed Brown   {
206c4762a1bSJed Brown     PetscReal   *array_x, *array_y;
207c4762a1bSJed Brown     PetscInt     npoints, i, j, cnt;
208c4762a1bSJed Brown     DMDACoor2d **LA_coor;
209c4762a1bSJed Brown     Vec          coor;
210c4762a1bSJed Brown     DM           dmcellcdm;
211c4762a1bSJed Brown 
2129566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell, &dmcellcdm));
2139566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell, &coor));
2149566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm, coor, &LA_coor));
2159566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
2169566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
2179566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
218c4762a1bSJed Brown     cnt = 0;
219c4762a1bSJed Brown     for (j = js; j < js + nj; j++) {
220c4762a1bSJed Brown       for (i = is; i < is + ni; i++) {
221c4762a1bSJed Brown         PetscReal xp, yp;
222c4762a1bSJed Brown         xp                   = PetscRealPart(LA_coor[j][i].x);
223c4762a1bSJed Brown         yp                   = PetscRealPart(LA_coor[j][i].y);
2249371c9d4SSatish Balay         array_x[4 * cnt + 0] = xp - 0.05;
225ad540459SPierre Jolivet         if (array_x[4 * cnt + 0] < -1.0) array_x[4 * cnt + 0] = -1.0 + 1.0e-12;
2269371c9d4SSatish Balay         array_x[4 * cnt + 1] = xp + 0.05;
227ad540459SPierre Jolivet         if (array_x[4 * cnt + 1] > 1.0) array_x[4 * cnt + 1] = 1.0 - 1.0e-12;
2289371c9d4SSatish Balay         array_x[4 * cnt + 2] = xp - 0.05;
229ad540459SPierre Jolivet         if (array_x[4 * cnt + 2] < -1.0) array_x[4 * cnt + 2] = -1.0 + 1.0e-12;
2309371c9d4SSatish Balay         array_x[4 * cnt + 3] = xp + 0.05;
231ad540459SPierre Jolivet         if (array_x[4 * cnt + 3] > 1.0) array_x[4 * cnt + 3] = 1.0 - 1.0e-12;
232c4762a1bSJed Brown 
2339371c9d4SSatish Balay         array_y[4 * cnt + 0] = yp - 0.05;
234ad540459SPierre Jolivet         if (array_y[4 * cnt + 0] < -1.0) array_y[4 * cnt + 0] = -1.0 + 1.0e-12;
2359371c9d4SSatish Balay         array_y[4 * cnt + 1] = yp - 0.05;
236ad540459SPierre Jolivet         if (array_y[4 * cnt + 1] < -1.0) array_y[4 * cnt + 1] = -1.0 + 1.0e-12;
2379371c9d4SSatish Balay         array_y[4 * cnt + 2] = yp + 0.05;
238ad540459SPierre Jolivet         if (array_y[4 * cnt + 2] > 1.0) array_y[4 * cnt + 2] = 1.0 - 1.0e-12;
2399371c9d4SSatish Balay         array_y[4 * cnt + 3] = yp + 0.05;
240ad540459SPierre Jolivet         if (array_y[4 * cnt + 3] > 1.0) array_y[4 * cnt + 3] = 1.0 - 1.0e-12;
241c4762a1bSJed Brown         cnt++;
242c4762a1bSJed Brown       }
243c4762a1bSJed Brown     }
2449566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
2459566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
2469566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm, coor, &LA_coor));
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown   {
249c4762a1bSJed Brown     PetscInt npoints[2], npoints_orig[2], ng;
250c4762a1bSJed Brown 
2519566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
2529566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
25463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
2559566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
2579566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_DMDABoundingBox(dms, &ng));
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
2609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
2619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
26263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
2639566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
265c4762a1bSJed Brown   }
266c4762a1bSJed Brown   {
267c4762a1bSJed Brown     PetscReal *array_x, *array_y;
268c4762a1bSJed Brown     PetscInt   npoints, p;
269c4762a1bSJed Brown     FILE      *fp = NULL;
270c4762a1bSJed Brown     char       name[PETSC_MAX_PATH_LEN];
271c4762a1bSJed Brown 
2729566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "c-rank%d.gp", rank));
273c4762a1bSJed Brown     fp = fopen(name, "w");
27428b400f6SJacob Faibussowitsch     PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
2759566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
2769566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
2779566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
278ad540459SPierre Jolivet     for (p = 0; p < npoints; p++) fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array_x[p], array_y[p], (double)rank);
2799566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
2809566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
281c4762a1bSJed Brown     fclose(fp);
282c4762a1bSJed Brown   }
2839566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
2849566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown typedef struct {
289c4762a1bSJed Brown   PetscReal cx[2];
290c4762a1bSJed Brown   PetscReal radius;
291c4762a1bSJed Brown } CollectZoneCtx;
292c4762a1bSJed Brown 
293d71ae5a4SJacob Faibussowitsch PetscErrorCode collect_zone(DM dm, void *ctx, PetscInt *nfound, PetscInt **foundlist)
294d71ae5a4SJacob Faibussowitsch {
295c4762a1bSJed Brown   CollectZoneCtx *zone = (CollectZoneCtx *)ctx;
296c4762a1bSJed Brown   PetscInt        p, npoints;
297c4762a1bSJed Brown   PetscReal      *array_x, *array_y, r2;
298c4762a1bSJed Brown   PetscInt        p2collect, *plist;
299c4762a1bSJed Brown   PetscMPIInt     rank;
300c4762a1bSJed Brown 
301362febeeSStefano Zampini   PetscFunctionBegin;
3029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
3039566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &npoints));
3049566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
3059566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   r2        = zone->radius * zone->radius;
308c4762a1bSJed Brown   p2collect = 0;
309c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
310c4762a1bSJed Brown     PetscReal sep2;
311c4762a1bSJed Brown 
312c4762a1bSJed Brown     sep2 = (array_x[p] - zone->cx[0]) * (array_x[p] - zone->cx[0]);
313c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1]) * (array_y[p] - zone->cx[1]);
314ad540459SPierre Jolivet     if (sep2 < r2) p2collect++;
315c4762a1bSJed Brown   }
316c4762a1bSJed Brown 
3179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(p2collect + 1, &plist));
318c4762a1bSJed Brown   p2collect = 0;
319c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
320c4762a1bSJed Brown     PetscReal sep2;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown     sep2 = (array_x[p] - zone->cx[0]) * (array_x[p] - zone->cx[0]);
323c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1]) * (array_y[p] - zone->cx[1]);
324c4762a1bSJed Brown     if (sep2 < r2) {
325c4762a1bSJed Brown       plist[p2collect] = p;
326c4762a1bSJed Brown       p2collect++;
327c4762a1bSJed Brown     }
328c4762a1bSJed Brown   }
3299566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
3309566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   *nfound    = p2collect;
333c4762a1bSJed Brown   *foundlist = plist;
334c4762a1bSJed Brown   PetscFunctionReturn(0);
335c4762a1bSJed Brown }
336c4762a1bSJed Brown 
337d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_4(void)
338d71ae5a4SJacob Faibussowitsch {
339c4762a1bSJed Brown   DM              dms;
340c4762a1bSJed Brown   PetscMPIInt     rank, size;
341c4762a1bSJed Brown   PetscInt        is, js, ni, nj, overlap, nn;
342c4762a1bSJed Brown   DM              dmcell;
343c4762a1bSJed Brown   CollectZoneCtx *zone;
344c4762a1bSJed Brown   PetscReal       dx;
345c4762a1bSJed Brown 
346362febeeSStefano Zampini   PetscFunctionBegin;
3479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
349c4762a1bSJed Brown   nn      = 101;
350c4762a1bSJed Brown   dx      = 2.0 / (PetscReal)(nn - 1);
351c4762a1bSJed Brown   overlap = 0;
3529566063dSJacob 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));
3539566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
3549566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
3559566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
3569566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell, &is, &js, NULL, &ni, &nj, NULL));
3579566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
3589566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
3596a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* load in data types */
3629566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
3639566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
3649566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coorx", 1, PETSC_REAL));
3659566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coory", 1, PETSC_REAL));
3669566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
3679566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, ni * nj * 4, 4));
3689566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   /* set values within the swarm */
371c4762a1bSJed Brown   {
372c4762a1bSJed Brown     PetscReal   *array_x, *array_y;
373c4762a1bSJed Brown     PetscInt     npoints, i, j, cnt;
374c4762a1bSJed Brown     DMDACoor2d **LA_coor;
375c4762a1bSJed Brown     Vec          coor;
376c4762a1bSJed Brown     DM           dmcellcdm;
377c4762a1bSJed Brown 
3789566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell, &dmcellcdm));
3799566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell, &coor));
3809566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm, coor, &LA_coor));
3819566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
3829566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
3839566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
384c4762a1bSJed Brown     cnt = 0;
385c4762a1bSJed Brown     for (j = js; j < js + nj; j++) {
386c4762a1bSJed Brown       for (i = is; i < is + ni; i++) {
387c4762a1bSJed Brown         PetscReal xp, yp;
388c4762a1bSJed Brown 
389c4762a1bSJed Brown         xp                   = PetscRealPart(LA_coor[j][i].x);
390c4762a1bSJed Brown         yp                   = PetscRealPart(LA_coor[j][i].y);
391ad540459SPierre Jolivet         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;*/
392c4762a1bSJed 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; }*/
393ad540459SPierre Jolivet         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;*/
394c4762a1bSJed 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; }*/
395ad540459SPierre Jolivet         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;*/
396ad540459SPierre Jolivet         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;*/
397c4762a1bSJed 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; }*/
398c4762a1bSJed 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; }*/
399c4762a1bSJed Brown         cnt++;
400c4762a1bSJed Brown       }
401c4762a1bSJed Brown     }
4029566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
4039566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
4049566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm, coor, &LA_coor));
405c4762a1bSJed Brown   }
4069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &zone));
407c4762a1bSJed Brown   if (size == 4) {
408c4762a1bSJed Brown     if (rank == 0) {
409c4762a1bSJed Brown       zone->cx[0]  = 0.5;
410c4762a1bSJed Brown       zone->cx[1]  = 0.5;
411c4762a1bSJed Brown       zone->radius = 0.3;
412c4762a1bSJed Brown     }
413c4762a1bSJed Brown     if (rank == 1) {
414c4762a1bSJed Brown       zone->cx[0]  = -0.5;
415c4762a1bSJed Brown       zone->cx[1]  = 0.5;
416c4762a1bSJed Brown       zone->radius = 0.25;
417c4762a1bSJed Brown     }
418c4762a1bSJed Brown     if (rank == 2) {
419c4762a1bSJed Brown       zone->cx[0]  = 0.5;
420c4762a1bSJed Brown       zone->cx[1]  = -0.5;
421c4762a1bSJed Brown       zone->radius = 0.2;
422c4762a1bSJed Brown     }
423c4762a1bSJed Brown     if (rank == 3) {
424c4762a1bSJed Brown       zone->cx[0]  = -0.5;
425c4762a1bSJed Brown       zone->cx[1]  = -0.5;
426c4762a1bSJed Brown       zone->radius = 0.1;
427c4762a1bSJed Brown     }
428c4762a1bSJed Brown   } else {
429c4762a1bSJed Brown     if (rank == 0) {
430c4762a1bSJed Brown       zone->cx[0]  = 0.5;
431c4762a1bSJed Brown       zone->cx[1]  = 0.5;
432c4762a1bSJed Brown       zone->radius = 0.8;
433c4762a1bSJed Brown     } else {
434c4762a1bSJed Brown       zone->cx[0]  = 10.0;
435c4762a1bSJed Brown       zone->cx[1]  = 10.0;
436c4762a1bSJed Brown       zone->radius = 0.0;
437c4762a1bSJed Brown     }
438c4762a1bSJed Brown   }
439c4762a1bSJed Brown   {
440c4762a1bSJed Brown     PetscInt npoints[2], npoints_orig[2], ng;
441c4762a1bSJed Brown 
4429566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
4439566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
4449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
44563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
4469566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
4489566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_General(dms, collect_zone, sizeof(CollectZoneCtx), zone, &ng));
4499566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
4509566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
4519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
45263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
4539566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
455c4762a1bSJed Brown   }
456c4762a1bSJed Brown   {
457c4762a1bSJed Brown     PetscReal *array_x, *array_y;
458c4762a1bSJed Brown     PetscInt   npoints, p;
459c4762a1bSJed Brown     FILE      *fp = NULL;
460c4762a1bSJed Brown     char       name[PETSC_MAX_PATH_LEN];
461c4762a1bSJed Brown 
4629566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "c-rank%d.gp", rank));
463c4762a1bSJed Brown     fp = fopen(name, "w");
46428b400f6SJacob Faibussowitsch     PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
4659566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
4669566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
4679566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
468ad540459SPierre Jolivet     for (p = 0; p < npoints; p++) fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array_x[p], array_y[p], (double)rank);
4699566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
4709566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
471c4762a1bSJed Brown     fclose(fp);
472c4762a1bSJed Brown   }
4739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
4749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
4759566063dSJacob Faibussowitsch   PetscCall(PetscFree(zone));
476c4762a1bSJed Brown   PetscFunctionReturn(0);
477c4762a1bSJed Brown }
478c4762a1bSJed Brown 
479d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
480d71ae5a4SJacob Faibussowitsch {
481c4762a1bSJed Brown   PetscInt test_mode = 4;
482c4762a1bSJed Brown 
483327415f7SBarry Smith   PetscFunctionBeginUser;
4849566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
4859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_mode", &test_mode, NULL));
486c4762a1bSJed Brown   if (test_mode == 1) {
4879566063dSJacob Faibussowitsch     PetscCall(ex1_1());
488c4762a1bSJed Brown   } else if (test_mode == 2) {
4899566063dSJacob Faibussowitsch     PetscCall(ex1_2());
490c4762a1bSJed Brown   } else if (test_mode == 3) {
4919566063dSJacob Faibussowitsch     PetscCall(ex1_3());
492c4762a1bSJed Brown   } else if (test_mode == 4) {
4939566063dSJacob Faibussowitsch     PetscCall(ex1_4());
494c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown test_mode value, should be 1,2,3,4");
4959566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
496b122ec5aSJacob Faibussowitsch   return 0;
497c4762a1bSJed Brown }
498c4762a1bSJed Brown 
499c4762a1bSJed Brown /*TEST
500c4762a1bSJed Brown 
501c4762a1bSJed Brown    build:
502c4762a1bSJed Brown       requires: !complex double
503c4762a1bSJed Brown 
504c4762a1bSJed Brown    test:
505c4762a1bSJed Brown       args: -test_mode 1
506c4762a1bSJed Brown       filter: grep -v atomic
507c4762a1bSJed Brown       filter_output: grep -v atomic
508c4762a1bSJed Brown 
509c4762a1bSJed Brown    test:
510c4762a1bSJed Brown       suffix: 2
511c4762a1bSJed Brown       args: -test_mode 2
512c4762a1bSJed Brown       filter: grep -v atomic
513c4762a1bSJed Brown       filter_output: grep -v atomic
514c4762a1bSJed Brown 
515c4762a1bSJed Brown    test:
516c4762a1bSJed Brown       suffix: 3
517c4762a1bSJed Brown       args: -test_mode 3
518c4762a1bSJed Brown       filter: grep -v atomic
519c4762a1bSJed Brown       filter_output: grep -v atomic
520c4762a1bSJed Brown       TODO: broken
521c4762a1bSJed Brown 
522c4762a1bSJed Brown    test:
523c4762a1bSJed Brown       suffix: 4
524c4762a1bSJed Brown       args: -test_mode 4
525c4762a1bSJed Brown       filter: grep -v atomic
526c4762a1bSJed Brown       filter_output: grep -v atomic
527c4762a1bSJed Brown 
528c4762a1bSJed Brown    test:
529c4762a1bSJed Brown       suffix: 5
530c4762a1bSJed Brown       nsize: 4
531c4762a1bSJed Brown       args: -test_mode 1
532c4762a1bSJed Brown       filter: grep -v atomic
533c4762a1bSJed Brown       filter_output: grep -v atomic
534c4762a1bSJed Brown 
535c4762a1bSJed Brown    test:
536c4762a1bSJed Brown       suffix: 6
537c4762a1bSJed Brown       nsize: 4
538c4762a1bSJed Brown       args: -test_mode 2
539c4762a1bSJed Brown       filter: grep -v atomic
540c4762a1bSJed Brown       filter_output: grep -v atomic
541c4762a1bSJed Brown 
542c4762a1bSJed Brown    test:
543c4762a1bSJed Brown       suffix: 7
544c4762a1bSJed Brown       nsize: 4
545c4762a1bSJed Brown       args: -test_mode 3
546c4762a1bSJed Brown       filter: grep -v atomic
547c4762a1bSJed Brown       filter_output: grep -v atomic
548c4762a1bSJed Brown       TODO: broken
549c4762a1bSJed Brown 
550c4762a1bSJed Brown    test:
551c4762a1bSJed Brown       suffix: 8
552c4762a1bSJed Brown       nsize: 4
553c4762a1bSJed Brown       args: -test_mode 4
554c4762a1bSJed Brown       filter: grep -v atomic
555c4762a1bSJed Brown       filter_output: grep -v atomic
556c4762a1bSJed Brown 
557c4762a1bSJed Brown TEST*/
558