xref: /petsc/src/dm/tutorials/swarm_ex1.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1fca3fa51SMatthew G. Knepley #include "petscis.h"
2fca3fa51SMatthew G. Knepley #include "petscsys.h"
3fca3fa51SMatthew G. Knepley #include "petscsystypes.h"
4fca3fa51SMatthew G. Knepley #include "petscvec.h"
5c4762a1bSJed Brown static char help[] = "Tests DMSwarm\n\n";
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscdm.h>
8c4762a1bSJed Brown #include <petscdmda.h>
9c4762a1bSJed Brown #include <petscdmswarm.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM, PetscErrorCode (*)(DM, void *, PetscInt *, PetscInt **), size_t, void *, PetscInt *);
12c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM, PetscInt *);
13c4762a1bSJed Brown 
ex1_1(void)14d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_1(void)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown   DM          dms;
17c4762a1bSJed Brown   Vec         x;
18c4762a1bSJed Brown   PetscMPIInt rank, size;
19c4762a1bSJed Brown   PetscInt    p;
20c4762a1bSJed Brown 
21362febeeSStefano Zampini   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
24d5b43468SJose E. Roman   PetscCheck(!(size > 1) || !(size != 4), PETSC_COMM_WORLD, PETSC_ERR_SUP, "Must be run with 4 MPI ranks");
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
279566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
286a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
319566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
329566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 1, PETSC_REAL));
339566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
349566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
359566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   {
38c4762a1bSJed Brown     PetscReal *array;
399566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", NULL, NULL, (void **)&array));
40ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 11.1 + p * 0.1 + rank * 100.0;
419566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", NULL, NULL, (void **)&array));
42c4762a1bSJed Brown   }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   {
45c4762a1bSJed Brown     PetscReal *array;
469566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", NULL, NULL, (void **)&array));
47ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 2.0e-2 + p * 0.001 + rank * 1.0;
489566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", NULL, NULL, (void **)&array));
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
529566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(dms, "strain"));
559566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dms, &x));
569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   {
59c4762a1bSJed Brown     PetscInt *rankval;
60c4762a1bSJed Brown     PetscInt  npoints[2], npoints_orig[2];
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
639566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
64fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmGetField(dms, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
65c4762a1bSJed Brown     if ((rank == 0) && (size > 1)) {
66c4762a1bSJed Brown       rankval[0] = 1;
67c4762a1bSJed Brown       rankval[3] = 1;
68c4762a1bSJed Brown     }
69ad540459SPierre Jolivet     if (rank == 3) rankval[2] = 1;
70fca3fa51SMatthew G. Knepley     PetscCall(DMSwarmRestoreField(dms, DMSwarmField_rank, NULL, NULL, (void **)&rankval));
719566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
729566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
739566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
7563a3b9bcSJacob 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]));
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown   {
809566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
819566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
829566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown   {
859566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "strain", &x));
869566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
879566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "strain", &x));
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
ex1_2(void)94d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_2(void)
95d71ae5a4SJacob Faibussowitsch {
96c4762a1bSJed Brown   DM          dms;
97c4762a1bSJed Brown   Vec         x;
98c4762a1bSJed Brown   PetscMPIInt rank, size;
99c4762a1bSJed Brown   PetscInt    p;
100c4762a1bSJed Brown 
101362febeeSStefano Zampini   PetscFunctionBegin;
1029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1049566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1059566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1066a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
1079566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
108c4762a1bSJed Brown 
1099566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
1109566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 1, PETSC_REAL));
1119566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
1129566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
1139566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
114c4762a1bSJed Brown   {
115c4762a1bSJed Brown     PetscReal *array;
1169566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", NULL, NULL, (void **)&array));
117ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 11.1 + p * 0.1 + rank * 100.0;
1189566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", NULL, NULL, (void **)&array));
119c4762a1bSJed Brown   }
120c4762a1bSJed Brown   {
121c4762a1bSJed Brown     PetscReal *array;
1229566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", NULL, NULL, (void **)&array));
123ad540459SPierre Jolivet     for (p = 0; p < 5 + rank; p++) array[p] = 2.0e-2 + p * 0.001 + rank * 1.0;
1249566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", NULL, NULL, (void **)&array));
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown   {
127c4762a1bSJed Brown     PetscInt *rankval;
128c4762a1bSJed Brown     PetscInt  npoints[2], npoints_orig[2];
129c4762a1bSJed Brown 
1309566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
1319566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
1329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
13363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
1349566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
138c4762a1bSJed Brown 
139ad540459SPierre Jolivet     if (rank == 1) rankval[0] = -1;
140ad540459SPierre Jolivet     if (rank == 2) rankval[1] = -1;
141c4762a1bSJed Brown     if (rank == 3) {
142c4762a1bSJed Brown       rankval[3] = -1;
143c4762a1bSJed Brown       rankval[4] = -1;
144c4762a1bSJed Brown     }
1459566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "DMSwarm_rank", NULL, NULL, (void **)&rankval));
1469566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewCreate(dms));
1479566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
1489566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
15063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
1519566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
1559566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1569566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewDestroy(dms));
1599566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
1609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
1619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
16263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after_v(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
1639566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
165c4762a1bSJed Brown 
1669566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
1679566063dSJacob Faibussowitsch     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
1689566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
169c4762a1bSJed Brown   }
1709566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172c4762a1bSJed Brown }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown /*
175c4762a1bSJed Brown  splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp"
176c4762a1bSJed Brown */
ex1_3(void)177d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_3(void)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown   DM          dms;
180c4762a1bSJed Brown   PetscMPIInt rank, size;
181c4762a1bSJed Brown   PetscInt    is, js, ni, nj, overlap;
182c4762a1bSJed Brown   DM          dmcell;
183c4762a1bSJed Brown 
184362febeeSStefano Zampini   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
187c4762a1bSJed Brown   overlap = 2;
1889566063dSJacob 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));
1899566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
1909566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
1919566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
1929566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell, &is, &js, NULL, &ni, &nj, NULL));
1939566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
1949566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
1956a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
1969566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(dms, dmcell));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* load in data types */
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
2019566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coorx", 1, PETSC_REAL));
2029566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coory", 1, PETSC_REAL));
2039566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
2049566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, ni * nj * 4, 4));
2059566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /* set values within the swarm */
208c4762a1bSJed Brown   {
209c4762a1bSJed Brown     PetscReal   *array_x, *array_y;
210c4762a1bSJed Brown     PetscInt     npoints, i, j, cnt;
211c4762a1bSJed Brown     DMDACoor2d **LA_coor;
212c4762a1bSJed Brown     Vec          coor;
213c4762a1bSJed Brown     DM           dmcellcdm;
214c4762a1bSJed Brown 
2159566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell, &dmcellcdm));
2169566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell, &coor));
2179566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm, coor, &LA_coor));
2189566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
2199566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
2209566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
221c4762a1bSJed Brown     cnt = 0;
222c4762a1bSJed Brown     for (j = js; j < js + nj; j++) {
223c4762a1bSJed Brown       for (i = is; i < is + ni; i++) {
224c4762a1bSJed Brown         PetscReal xp, yp;
225c4762a1bSJed Brown         xp                   = PetscRealPart(LA_coor[j][i].x);
226c4762a1bSJed Brown         yp                   = PetscRealPart(LA_coor[j][i].y);
2279371c9d4SSatish Balay         array_x[4 * cnt + 0] = xp - 0.05;
228ad540459SPierre Jolivet         if (array_x[4 * cnt + 0] < -1.0) array_x[4 * cnt + 0] = -1.0 + 1.0e-12;
2299371c9d4SSatish Balay         array_x[4 * cnt + 1] = xp + 0.05;
230ad540459SPierre Jolivet         if (array_x[4 * cnt + 1] > 1.0) array_x[4 * cnt + 1] = 1.0 - 1.0e-12;
2319371c9d4SSatish Balay         array_x[4 * cnt + 2] = xp - 0.05;
232ad540459SPierre Jolivet         if (array_x[4 * cnt + 2] < -1.0) array_x[4 * cnt + 2] = -1.0 + 1.0e-12;
2339371c9d4SSatish Balay         array_x[4 * cnt + 3] = xp + 0.05;
234ad540459SPierre Jolivet         if (array_x[4 * cnt + 3] > 1.0) array_x[4 * cnt + 3] = 1.0 - 1.0e-12;
235c4762a1bSJed Brown 
2369371c9d4SSatish Balay         array_y[4 * cnt + 0] = yp - 0.05;
237ad540459SPierre Jolivet         if (array_y[4 * cnt + 0] < -1.0) array_y[4 * cnt + 0] = -1.0 + 1.0e-12;
2389371c9d4SSatish Balay         array_y[4 * cnt + 1] = yp - 0.05;
239ad540459SPierre Jolivet         if (array_y[4 * cnt + 1] < -1.0) array_y[4 * cnt + 1] = -1.0 + 1.0e-12;
2409371c9d4SSatish Balay         array_y[4 * cnt + 2] = yp + 0.05;
241ad540459SPierre Jolivet         if (array_y[4 * cnt + 2] > 1.0) array_y[4 * cnt + 2] = 1.0 - 1.0e-12;
2429371c9d4SSatish Balay         array_y[4 * cnt + 3] = yp + 0.05;
243ad540459SPierre Jolivet         if (array_y[4 * cnt + 3] > 1.0) array_y[4 * cnt + 3] = 1.0 - 1.0e-12;
244c4762a1bSJed Brown         cnt++;
245c4762a1bSJed Brown       }
246c4762a1bSJed Brown     }
2479566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
2489566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
2499566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm, coor, &LA_coor));
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown   {
252c4762a1bSJed Brown     PetscInt npoints[2], npoints_orig[2], ng;
253c4762a1bSJed Brown 
2549566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
2559566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
2569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
25763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
2589566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
2609566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_DMDABoundingBox(dms, &ng));
261c4762a1bSJed Brown 
2629566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
2639566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
2649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
26563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
2669566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
268c4762a1bSJed Brown   }
269c4762a1bSJed Brown   {
270c4762a1bSJed Brown     PetscReal *array_x, *array_y;
271c4762a1bSJed Brown     PetscInt   npoints, p;
272c4762a1bSJed Brown     FILE      *fp = NULL;
273c4762a1bSJed Brown     char       name[PETSC_MAX_PATH_LEN];
274c4762a1bSJed Brown 
2759566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "c-rank%d.gp", rank));
276c4762a1bSJed Brown     fp = fopen(name, "w");
27728b400f6SJacob Faibussowitsch     PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
2789566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
2799566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
2809566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
281ad540459SPierre Jolivet     for (p = 0; p < npoints; p++) fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array_x[p], array_y[p], (double)rank);
2829566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
2839566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
284c4762a1bSJed Brown     fclose(fp);
285c4762a1bSJed Brown   }
2869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
2879566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
291c4762a1bSJed Brown typedef struct {
292c4762a1bSJed Brown   PetscReal cx[2];
293c4762a1bSJed Brown   PetscReal radius;
294c4762a1bSJed Brown } CollectZoneCtx;
295c4762a1bSJed Brown 
collect_zone(DM dm,PetscCtx ctx,PetscInt * nfound,PetscInt ** foundlist)296*2a8381b2SBarry Smith PetscErrorCode collect_zone(DM dm, PetscCtx ctx, PetscInt *nfound, PetscInt **foundlist)
297d71ae5a4SJacob Faibussowitsch {
298c4762a1bSJed Brown   CollectZoneCtx *zone = (CollectZoneCtx *)ctx;
299c4762a1bSJed Brown   PetscInt        p, npoints;
300c4762a1bSJed Brown   PetscReal      *array_x, *array_y, r2;
301c4762a1bSJed Brown   PetscInt        p2collect, *plist;
302c4762a1bSJed Brown   PetscMPIInt     rank;
303c4762a1bSJed Brown 
304362febeeSStefano Zampini   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
3069566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm, &npoints));
3079566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, "coorx", NULL, NULL, (void **)&array_x));
3089566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, "coory", NULL, NULL, (void **)&array_y));
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   r2        = zone->radius * zone->radius;
311c4762a1bSJed Brown   p2collect = 0;
312c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
313c4762a1bSJed Brown     PetscReal sep2;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown     sep2 = (array_x[p] - zone->cx[0]) * (array_x[p] - zone->cx[0]);
316c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1]) * (array_y[p] - zone->cx[1]);
317ad540459SPierre Jolivet     if (sep2 < r2) p2collect++;
318c4762a1bSJed Brown   }
319c4762a1bSJed Brown 
3209566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(p2collect + 1, &plist));
321c4762a1bSJed Brown   p2collect = 0;
322c4762a1bSJed Brown   for (p = 0; p < npoints; p++) {
323c4762a1bSJed Brown     PetscReal sep2;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     sep2 = (array_x[p] - zone->cx[0]) * (array_x[p] - zone->cx[0]);
326c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1]) * (array_y[p] - zone->cx[1]);
327c4762a1bSJed Brown     if (sep2 < r2) {
328c4762a1bSJed Brown       plist[p2collect] = p;
329c4762a1bSJed Brown       p2collect++;
330c4762a1bSJed Brown     }
331c4762a1bSJed Brown   }
3329566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, "coory", NULL, NULL, (void **)&array_y));
3339566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, "coorx", NULL, NULL, (void **)&array_x));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   *nfound    = p2collect;
336c4762a1bSJed Brown   *foundlist = plist;
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
ex1_4(void)340d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_4(void)
341d71ae5a4SJacob Faibussowitsch {
342c4762a1bSJed Brown   DM              dms;
343c4762a1bSJed Brown   PetscMPIInt     rank, size;
344c4762a1bSJed Brown   PetscInt        is, js, ni, nj, overlap, nn;
345c4762a1bSJed Brown   DM              dmcell;
346c4762a1bSJed Brown   CollectZoneCtx *zone;
347c4762a1bSJed Brown   PetscReal       dx;
348c4762a1bSJed Brown 
349362febeeSStefano Zampini   PetscFunctionBegin;
3509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
352c4762a1bSJed Brown   nn      = 101;
353c4762a1bSJed Brown   dx      = 2.0 / (PetscReal)(nn - 1);
354c4762a1bSJed Brown   overlap = 0;
3559566063dSJacob 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));
3569566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
3579566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
3589566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0));
3599566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell, &is, &js, NULL, &ni, &nj, NULL));
3609566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
3619566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
3626a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /* load in data types */
3659566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
3669566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
3679566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coorx", 1, PETSC_REAL));
3689566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "coory", 1, PETSC_REAL));
3699566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
3709566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, ni * nj * 4, 4));
3719566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   /* set values within the swarm */
374c4762a1bSJed Brown   {
375c4762a1bSJed Brown     PetscReal   *array_x, *array_y;
376c4762a1bSJed Brown     PetscInt     npoints, i, j, cnt;
377c4762a1bSJed Brown     DMDACoor2d **LA_coor;
378c4762a1bSJed Brown     Vec          coor;
379c4762a1bSJed Brown     DM           dmcellcdm;
380c4762a1bSJed Brown 
3819566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell, &dmcellcdm));
3829566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell, &coor));
3839566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm, coor, &LA_coor));
3849566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
3859566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
3869566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
387c4762a1bSJed Brown     cnt = 0;
388c4762a1bSJed Brown     for (j = js; j < js + nj; j++) {
389c4762a1bSJed Brown       for (i = is; i < is + ni; i++) {
390c4762a1bSJed Brown         PetscReal xp, yp;
391c4762a1bSJed Brown 
392c4762a1bSJed Brown         xp                   = PetscRealPart(LA_coor[j][i].x);
393c4762a1bSJed Brown         yp                   = PetscRealPart(LA_coor[j][i].y);
394ad540459SPierre 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;*/
395c4762a1bSJed 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; }*/
396ad540459SPierre 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;*/
397c4762a1bSJed 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; }*/
398ad540459SPierre 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;*/
399ad540459SPierre 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;*/
400c4762a1bSJed 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; }*/
401c4762a1bSJed 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; }*/
402c4762a1bSJed Brown         cnt++;
403c4762a1bSJed Brown       }
404c4762a1bSJed Brown     }
4059566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
4069566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
4079566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm, coor, &LA_coor));
408c4762a1bSJed Brown   }
4099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, &zone));
410c4762a1bSJed Brown   if (size == 4) {
411c4762a1bSJed Brown     if (rank == 0) {
412c4762a1bSJed Brown       zone->cx[0]  = 0.5;
413c4762a1bSJed Brown       zone->cx[1]  = 0.5;
414c4762a1bSJed Brown       zone->radius = 0.3;
415c4762a1bSJed Brown     }
416c4762a1bSJed Brown     if (rank == 1) {
417c4762a1bSJed Brown       zone->cx[0]  = -0.5;
418c4762a1bSJed Brown       zone->cx[1]  = 0.5;
419c4762a1bSJed Brown       zone->radius = 0.25;
420c4762a1bSJed Brown     }
421c4762a1bSJed Brown     if (rank == 2) {
422c4762a1bSJed Brown       zone->cx[0]  = 0.5;
423c4762a1bSJed Brown       zone->cx[1]  = -0.5;
424c4762a1bSJed Brown       zone->radius = 0.2;
425c4762a1bSJed Brown     }
426c4762a1bSJed Brown     if (rank == 3) {
427c4762a1bSJed Brown       zone->cx[0]  = -0.5;
428c4762a1bSJed Brown       zone->cx[1]  = -0.5;
429c4762a1bSJed Brown       zone->radius = 0.1;
430c4762a1bSJed Brown     }
431c4762a1bSJed Brown   } else {
432c4762a1bSJed Brown     if (rank == 0) {
433c4762a1bSJed Brown       zone->cx[0]  = 0.5;
434c4762a1bSJed Brown       zone->cx[1]  = 0.5;
435c4762a1bSJed Brown       zone->radius = 0.8;
436c4762a1bSJed Brown     } else {
437c4762a1bSJed Brown       zone->cx[0]  = 10.0;
438c4762a1bSJed Brown       zone->cx[1]  = 10.0;
439c4762a1bSJed Brown       zone->radius = 0.0;
440c4762a1bSJed Brown     }
441c4762a1bSJed Brown   }
442c4762a1bSJed Brown   {
443c4762a1bSJed Brown     PetscInt npoints[2], npoints_orig[2], ng;
444c4762a1bSJed Brown 
4459566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints_orig[0]));
4469566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints_orig[1]));
4479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
44863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] before(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints_orig[0], npoints_orig[1]));
4499566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
4519566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_General(dms, collect_zone, sizeof(CollectZoneCtx), zone, &ng));
4529566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints[0]));
4539566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms, &npoints[1]));
4549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
45563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "rank[%d] after(%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, npoints[0], npoints[1]));
4569566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
458c4762a1bSJed Brown   }
459c4762a1bSJed Brown   {
460c4762a1bSJed Brown     PetscReal *array_x, *array_y;
461c4762a1bSJed Brown     PetscInt   npoints, p;
462c4762a1bSJed Brown     FILE      *fp = NULL;
463c4762a1bSJed Brown     char       name[PETSC_MAX_PATH_LEN];
464c4762a1bSJed Brown 
4659566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "c-rank%d.gp", rank));
466c4762a1bSJed Brown     fp = fopen(name, "w");
46728b400f6SJacob Faibussowitsch     PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
4689566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms, &npoints));
4699566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coorx", NULL, NULL, (void **)&array_x));
4709566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "coory", NULL, NULL, (void **)&array_y));
471ad540459SPierre Jolivet     for (p = 0; p < npoints; p++) fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array_x[p], array_y[p], (double)rank);
4729566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coory", NULL, NULL, (void **)&array_y));
4739566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "coorx", NULL, NULL, (void **)&array_x));
474c4762a1bSJed Brown     fclose(fp);
475c4762a1bSJed Brown   }
4769566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
4779566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
4789566063dSJacob Faibussowitsch   PetscCall(PetscFree(zone));
4793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
480c4762a1bSJed Brown }
481c4762a1bSJed Brown 
main(int argc,char ** argv)482d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
483d71ae5a4SJacob Faibussowitsch {
484c4762a1bSJed Brown   PetscInt test_mode = 4;
485c4762a1bSJed Brown 
486327415f7SBarry Smith   PetscFunctionBeginUser;
487c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_mode", &test_mode, NULL));
489c4762a1bSJed Brown   if (test_mode == 1) {
4909566063dSJacob Faibussowitsch     PetscCall(ex1_1());
491c4762a1bSJed Brown   } else if (test_mode == 2) {
4929566063dSJacob Faibussowitsch     PetscCall(ex1_2());
493c4762a1bSJed Brown   } else if (test_mode == 3) {
4949566063dSJacob Faibussowitsch     PetscCall(ex1_3());
495c4762a1bSJed Brown   } else if (test_mode == 4) {
4969566063dSJacob Faibussowitsch     PetscCall(ex1_4());
497c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_USER, "Unknown test_mode value, should be 1,2,3,4");
4989566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
499b122ec5aSJacob Faibussowitsch   return 0;
500c4762a1bSJed Brown }
501c4762a1bSJed Brown 
502c4762a1bSJed Brown /*TEST
503c4762a1bSJed Brown 
504c4762a1bSJed Brown    build:
505c4762a1bSJed Brown       requires: !complex double
506c4762a1bSJed Brown 
507c4762a1bSJed Brown    test:
508c4762a1bSJed Brown       args: -test_mode 1
509c4762a1bSJed Brown       filter: grep -v atomic
510c4762a1bSJed Brown       filter_output: grep -v atomic
511c4762a1bSJed Brown 
512c4762a1bSJed Brown    test:
513c4762a1bSJed Brown       suffix: 2
514c4762a1bSJed Brown       args: -test_mode 2
515c4762a1bSJed Brown       filter: grep -v atomic
516c4762a1bSJed Brown       filter_output: grep -v atomic
517c4762a1bSJed Brown 
518c4762a1bSJed Brown    test:
519c4762a1bSJed Brown       suffix: 3
520c4762a1bSJed Brown       args: -test_mode 3
521c4762a1bSJed Brown       filter: grep -v atomic
522c4762a1bSJed Brown       filter_output: grep -v atomic
523c4762a1bSJed Brown       TODO: broken
524c4762a1bSJed Brown 
525c4762a1bSJed Brown    test:
526c4762a1bSJed Brown       suffix: 4
527c4762a1bSJed Brown       args: -test_mode 4
528c4762a1bSJed Brown       filter: grep -v atomic
529c4762a1bSJed Brown       filter_output: grep -v atomic
530c4762a1bSJed Brown 
531c4762a1bSJed Brown    test:
532c4762a1bSJed Brown       suffix: 5
533c4762a1bSJed Brown       nsize: 4
534c4762a1bSJed Brown       args: -test_mode 1
535c4762a1bSJed Brown       filter: grep -v atomic
536c4762a1bSJed Brown       filter_output: grep -v atomic
537c4762a1bSJed Brown 
538c4762a1bSJed Brown    test:
539c4762a1bSJed Brown       suffix: 6
540c4762a1bSJed Brown       nsize: 4
541c4762a1bSJed Brown       args: -test_mode 2
542c4762a1bSJed Brown       filter: grep -v atomic
543c4762a1bSJed Brown       filter_output: grep -v atomic
544c4762a1bSJed Brown 
545c4762a1bSJed Brown    test:
546c4762a1bSJed Brown       suffix: 7
547c4762a1bSJed Brown       nsize: 4
548c4762a1bSJed Brown       args: -test_mode 3
549c4762a1bSJed Brown       filter: grep -v atomic
550c4762a1bSJed Brown       filter_output: grep -v atomic
551c4762a1bSJed Brown       TODO: broken
552c4762a1bSJed Brown 
553c4762a1bSJed Brown    test:
554c4762a1bSJed Brown       suffix: 8
555c4762a1bSJed Brown       nsize: 4
556c4762a1bSJed Brown       args: -test_mode 4
557c4762a1bSJed Brown       filter: grep -v atomic
558c4762a1bSJed Brown       filter_output: grep -v atomic
559c4762a1bSJed Brown 
560c4762a1bSJed Brown TEST*/
561