xref: /petsc/src/dm/tutorials/swarm_ex1.c (revision 6a5217c03994f2d95bb2e6dbd8bed42381aeb015)
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 
11c4762a1bSJed Brown PetscErrorCode ex1_1(void)
12c4762a1bSJed Brown {
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));
2108401ef6SPierre Jolivet   PetscCheck(!(size > 1) || !(size != 4),PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must be run wuth 4 MPI ranks");
22c4762a1bSJed Brown 
239566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&dms));
249566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms,DMSWARM));
25*6a5217c0SMatthew 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));
37c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
38c4762a1bSJed Brown       array[p] = 11.1 + p*0.1 + rank*100.0;
39c4762a1bSJed Brown     }
409566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array));
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   {
44c4762a1bSJed Brown     PetscReal *array;
459566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array));
46c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
47c4762a1bSJed Brown       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
48c4762a1bSJed Brown     }
499566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array));
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x));
539566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x));
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(dms,"strain"));
569566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dms,&x));
579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   {
60c4762a1bSJed Brown     PetscInt    *rankval;
61c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2];
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0]));
649566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints_orig[1]));
659566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval));
66c4762a1bSJed Brown     if ((rank == 0) && (size > 1)) {
67c4762a1bSJed Brown       rankval[0] = 1;
68c4762a1bSJed Brown       rankval[3] = 1;
69c4762a1bSJed Brown     }
70c4762a1bSJed Brown     if (rank == 3) {
71c4762a1bSJed Brown       rankval[2] = 1;
72c4762a1bSJed Brown     }
739566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval));
749566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(dms,PETSC_TRUE));
759566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints[0]));
769566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints[1]));
779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D) after(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1],npoints[0],npoints[1]));
799566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
809566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
81c4762a1bSJed Brown   }
82c4762a1bSJed Brown   {
839566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x));
849566063dSJacob Faibussowitsch     PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
859566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x));
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   {
889566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"strain",&x));
899566063dSJacob Faibussowitsch     PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
909566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown PetscErrorCode ex1_2(void)
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   DM             dms;
100c4762a1bSJed Brown   Vec            x;
101c4762a1bSJed Brown   PetscMPIInt    rank,size;
102c4762a1bSJed Brown   PetscInt       p;
103c4762a1bSJed Brown 
104362febeeSStefano Zampini   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1079566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&dms));
1089566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms,DMSWARM));
109*6a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) dms, "Particles"));
1109566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL));
1139566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL));
1149566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
1159566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms,5+rank,4));
1169566063dSJacob Faibussowitsch   PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD));
117c4762a1bSJed Brown   {
118c4762a1bSJed Brown     PetscReal *array;
1199566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array));
120c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
121c4762a1bSJed Brown       array[p] = 11.1 + p*0.1 + rank*100.0;
122c4762a1bSJed Brown     }
1239566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array));
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown   {
126c4762a1bSJed Brown     PetscReal *array;
1279566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array));
128c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
129c4762a1bSJed Brown       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
130c4762a1bSJed Brown     }
1319566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array));
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown   {
134c4762a1bSJed Brown     PetscInt    *rankval;
135c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2];
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0]));
1389566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints_orig[1]));
1399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
1409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]));
1419566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
143c4762a1bSJed Brown 
1449566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown     if (rank == 1) {
147c4762a1bSJed Brown       rankval[0] = -1;
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown     if (rank == 2) {
150c4762a1bSJed Brown       rankval[1] = -1;
151c4762a1bSJed Brown     }
152c4762a1bSJed Brown     if (rank == 3) {
153c4762a1bSJed Brown       rankval[3] = -1;
154c4762a1bSJed Brown       rankval[4] = -1;
155c4762a1bSJed Brown     }
1569566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval));
1579566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewCreate(dms));
1589566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints[0]));
1599566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints[1]));
1609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
1619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]));
1629566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x));
1669566063dSJacob Faibussowitsch     PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
1679566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x));
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollectViewDestroy(dms));
1709566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints[0]));
1719566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints[1]));
1729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
1739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after_v(%D,%D)\n",rank,npoints[0],npoints[1]));
1749566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x));
1789566063dSJacob Faibussowitsch     PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
1799566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x));
180c4762a1bSJed Brown   }
1819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /*
186c4762a1bSJed Brown  splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp"
187c4762a1bSJed Brown */
188c4762a1bSJed Brown PetscErrorCode ex1_3(void)
189c4762a1bSJed Brown {
190c4762a1bSJed Brown   DM             dms;
191c4762a1bSJed Brown   PetscMPIInt    rank,size;
192c4762a1bSJed Brown   PetscInt       is,js,ni,nj,overlap;
193c4762a1bSJed Brown   DM             dmcell;
194c4762a1bSJed Brown 
195362febeeSStefano Zampini   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
198c4762a1bSJed Brown   overlap = 2;
1999566063dSJacob 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));
2009566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
2019566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
2029566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0));
2039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL));
2049566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&dms));
2059566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms,DMSWARM));
206*6a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) dms, "Particles"));
2079566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(dms,dmcell));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* load in data types */
2109566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
2119566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL));
2129566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL));
2139566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL));
2149566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
2159566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms,ni*nj*4,4));
2169566063dSJacob Faibussowitsch   PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* set values within the swarm */
219c4762a1bSJed Brown   {
220c4762a1bSJed Brown     PetscReal  *array_x,*array_y;
221c4762a1bSJed Brown     PetscInt   npoints,i,j,cnt;
222c4762a1bSJed Brown     DMDACoor2d **LA_coor;
223c4762a1bSJed Brown     Vec        coor;
224c4762a1bSJed Brown     DM         dmcellcdm;
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell,&dmcellcdm));
2279566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell,&coor));
2289566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm,coor,&LA_coor));
2299566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints));
2309566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x));
2319566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y));
232c4762a1bSJed Brown     cnt = 0;
233c4762a1bSJed Brown     for (j=js; j<js+nj; j++) {
234c4762a1bSJed Brown       for (i=is; i<is+ni; i++) {
235c4762a1bSJed Brown         PetscReal xp,yp;
236c4762a1bSJed Brown         xp = PetscRealPart(LA_coor[j][i].x);
237c4762a1bSJed Brown         yp = PetscRealPart(LA_coor[j][i].y);
238c4762a1bSJed Brown         array_x[4*cnt+0] = xp - 0.05; if (array_x[4*cnt+0] < -1.0) { array_x[4*cnt+0] = -1.0+1.0e-12; }
239c4762a1bSJed Brown         array_x[4*cnt+1] = xp + 0.05; if (array_x[4*cnt+1] > 1.0)  { array_x[4*cnt+1] =  1.0-1.0e-12; }
240c4762a1bSJed Brown         array_x[4*cnt+2] = xp - 0.05; if (array_x[4*cnt+2] < -1.0) { array_x[4*cnt+2] = -1.0+1.0e-12; }
241c4762a1bSJed Brown         array_x[4*cnt+3] = xp + 0.05; if (array_x[4*cnt+3] > 1.0)  { array_x[4*cnt+3] =  1.0-1.0e-12; }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown         array_y[4*cnt+0] = yp - 0.05; if (array_y[4*cnt+0] < -1.0) { array_y[4*cnt+0] = -1.0+1.0e-12; }
244c4762a1bSJed Brown         array_y[4*cnt+1] = yp - 0.05; if (array_y[4*cnt+1] < -1.0) { array_y[4*cnt+1] = -1.0+1.0e-12; }
245c4762a1bSJed Brown         array_y[4*cnt+2] = yp + 0.05; if (array_y[4*cnt+2] > 1.0)  { array_y[4*cnt+2] =  1.0-1.0e-12; }
246c4762a1bSJed Brown         array_y[4*cnt+3] = yp + 0.05; if (array_y[4*cnt+3] > 1.0)  { array_y[4*cnt+3] =  1.0-1.0e-12; }
247c4762a1bSJed Brown         cnt++;
248c4762a1bSJed Brown       }
249c4762a1bSJed Brown     }
2509566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y));
2519566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x));
2529566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor));
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown   {
255c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2],ng;
256c4762a1bSJed Brown 
2579566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0]));
2589566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints_orig[1]));
2599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
2609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]));
2619566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
2639566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_DMDABoundingBox(dms,&ng));
264c4762a1bSJed Brown 
2659566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints[0]));
2669566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints[1]));
2679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
2689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]));
2699566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
2709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
271c4762a1bSJed Brown   }
272c4762a1bSJed Brown   {
273c4762a1bSJed Brown     PetscReal *array_x,*array_y;
274c4762a1bSJed Brown     PetscInt  npoints,p;
275c4762a1bSJed Brown     FILE      *fp = NULL;
276c4762a1bSJed Brown     char      name[PETSC_MAX_PATH_LEN];
277c4762a1bSJed Brown 
2789566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank));
279c4762a1bSJed Brown     fp = fopen(name,"w");
28028b400f6SJacob Faibussowitsch     PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
2819566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints));
2829566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x));
2839566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y));
284c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
285c4762a1bSJed Brown       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
286c4762a1bSJed Brown     }
2879566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y));
2889566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x));
289c4762a1bSJed Brown     fclose(fp);
290c4762a1bSJed Brown   }
2919566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
2929566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
293c4762a1bSJed Brown   PetscFunctionReturn(0);
294c4762a1bSJed Brown }
295c4762a1bSJed Brown 
296c4762a1bSJed Brown typedef struct {
297c4762a1bSJed Brown   PetscReal cx[2];
298c4762a1bSJed Brown   PetscReal radius;
299c4762a1bSJed Brown } CollectZoneCtx;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown PetscErrorCode collect_zone(DM dm,void *ctx,PetscInt *nfound,PetscInt **foundlist)
302c4762a1bSJed Brown {
303c4762a1bSJed Brown   CollectZoneCtx *zone = (CollectZoneCtx*)ctx;
304c4762a1bSJed Brown   PetscInt       p,npoints;
305c4762a1bSJed Brown   PetscReal      *array_x,*array_y,r2;
306c4762a1bSJed Brown   PetscInt       p2collect,*plist;
307c4762a1bSJed Brown   PetscMPIInt    rank;
308c4762a1bSJed Brown 
309362febeeSStefano Zampini   PetscFunctionBegin;
3109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
3119566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dm,&npoints));
3129566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x));
3139566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   r2 = zone->radius * zone->radius;
316c4762a1bSJed Brown   p2collect = 0;
317c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
318c4762a1bSJed Brown     PetscReal sep2;
319c4762a1bSJed Brown 
320c4762a1bSJed Brown     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
321c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
322c4762a1bSJed Brown     if (sep2 < r2) {
323c4762a1bSJed Brown       p2collect++;
324c4762a1bSJed Brown     }
325c4762a1bSJed Brown   }
326c4762a1bSJed Brown 
3279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(p2collect+1,&plist));
328c4762a1bSJed Brown   p2collect = 0;
329c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
330c4762a1bSJed Brown     PetscReal sep2;
331c4762a1bSJed Brown 
332c4762a1bSJed Brown     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
333c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
334c4762a1bSJed Brown     if (sep2 < r2) {
335c4762a1bSJed Brown       plist[p2collect] = p;
336c4762a1bSJed Brown       p2collect++;
337c4762a1bSJed Brown     }
338c4762a1bSJed Brown   }
3399566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y));
3409566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x));
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   *nfound = p2collect;
343c4762a1bSJed Brown   *foundlist = plist;
344c4762a1bSJed Brown   PetscFunctionReturn(0);
345c4762a1bSJed Brown }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown PetscErrorCode ex1_4(void)
348c4762a1bSJed Brown {
349c4762a1bSJed Brown   DM             dms;
350c4762a1bSJed Brown   PetscMPIInt    rank,size;
351c4762a1bSJed Brown   PetscInt       is,js,ni,nj,overlap,nn;
352c4762a1bSJed Brown   DM             dmcell;
353c4762a1bSJed Brown   CollectZoneCtx *zone;
354c4762a1bSJed Brown   PetscReal      dx;
355c4762a1bSJed Brown 
356362febeeSStefano Zampini   PetscFunctionBegin;
3579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
3589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
359c4762a1bSJed Brown   nn = 101;
360c4762a1bSJed Brown   dx = 2.0/ (PetscReal)(nn-1);
361c4762a1bSJed Brown   overlap = 0;
3629566063dSJacob 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));
3639566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dmcell));
3649566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmcell));
3659566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0));
3669566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL));
3679566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&dms));
3689566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms,DMSWARM));
369*6a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) dms, "Particles"));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /* load in data types */
3729566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
3739566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL));
3749566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL));
3759566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL));
3769566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
3779566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms,ni*nj*4,4));
3789566063dSJacob Faibussowitsch   PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD));
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   /* set values within the swarm */
381c4762a1bSJed Brown   {
382c4762a1bSJed Brown     PetscReal  *array_x,*array_y;
383c4762a1bSJed Brown     PetscInt   npoints,i,j,cnt;
384c4762a1bSJed Brown     DMDACoor2d **LA_coor;
385c4762a1bSJed Brown     Vec        coor;
386c4762a1bSJed Brown     DM         dmcellcdm;
387c4762a1bSJed Brown 
3889566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dmcell,&dmcellcdm));
3899566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(dmcell,&coor));
3909566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(dmcellcdm,coor,&LA_coor));
3919566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints));
3929566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x));
3939566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y));
394c4762a1bSJed Brown     cnt = 0;
395c4762a1bSJed Brown     for (j=js; j<js+nj; j++) {
396c4762a1bSJed Brown       for (i=is; i<is+ni; i++) {
397c4762a1bSJed Brown         PetscReal xp,yp;
398c4762a1bSJed Brown 
399c4762a1bSJed Brown         xp = PetscRealPart(LA_coor[j][i].x);
400c4762a1bSJed Brown         yp = PetscRealPart(LA_coor[j][i].y);
401c4762a1bSJed 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; }*/
402c4762a1bSJed 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; }*/
403c4762a1bSJed 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; }*/
404c4762a1bSJed 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; }*/
405c4762a1bSJed 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; }*/
406c4762a1bSJed 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; }*/
407c4762a1bSJed 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; }*/
408c4762a1bSJed 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; }*/
409c4762a1bSJed Brown         cnt++;
410c4762a1bSJed Brown       }
411c4762a1bSJed Brown     }
4129566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y));
4139566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x));
4149566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor));
415c4762a1bSJed Brown   }
4169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1,&zone));
417c4762a1bSJed Brown   if (size == 4) {
418c4762a1bSJed Brown     if (rank == 0) {
419c4762a1bSJed Brown       zone->cx[0] = 0.5;
420c4762a1bSJed Brown       zone->cx[1] = 0.5;
421c4762a1bSJed Brown       zone->radius = 0.3;
422c4762a1bSJed Brown     }
423c4762a1bSJed Brown     if (rank == 1) {
424c4762a1bSJed Brown       zone->cx[0] = -0.5;
425c4762a1bSJed Brown       zone->cx[1] = 0.5;
426c4762a1bSJed Brown       zone->radius = 0.25;
427c4762a1bSJed Brown     }
428c4762a1bSJed Brown     if (rank == 2) {
429c4762a1bSJed Brown       zone->cx[0] = 0.5;
430c4762a1bSJed Brown       zone->cx[1] = -0.5;
431c4762a1bSJed Brown       zone->radius = 0.2;
432c4762a1bSJed Brown     }
433c4762a1bSJed Brown     if (rank == 3) {
434c4762a1bSJed Brown       zone->cx[0] = -0.5;
435c4762a1bSJed Brown       zone->cx[1] = -0.5;
436c4762a1bSJed Brown       zone->radius = 0.1;
437c4762a1bSJed Brown     }
438c4762a1bSJed Brown   } else {
439c4762a1bSJed Brown     if (rank == 0) {
440c4762a1bSJed Brown       zone->cx[0] = 0.5;
441c4762a1bSJed Brown       zone->cx[1] = 0.5;
442c4762a1bSJed Brown       zone->radius = 0.8;
443c4762a1bSJed Brown     } else {
444c4762a1bSJed Brown       zone->cx[0] = 10.0;
445c4762a1bSJed Brown       zone->cx[1] = 10.0;
446c4762a1bSJed Brown       zone->radius = 0.0;
447c4762a1bSJed Brown     }
448c4762a1bSJed Brown   }
449c4762a1bSJed Brown   {
450c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2],ng;
451c4762a1bSJed Brown 
4529566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints_orig[0]));
4539566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints_orig[1]));
4549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
4559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]));
4569566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
4589566063dSJacob Faibussowitsch     PetscCall(DMSwarmCollect_General(dms,collect_zone,sizeof(CollectZoneCtx),zone,&ng));
4599566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints[0]));
4609566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetSize(dms,&npoints[1]));
4619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD));
4629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]));
4639566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
4649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD));
465c4762a1bSJed Brown   }
466c4762a1bSJed Brown   {
467c4762a1bSJed Brown     PetscReal *array_x,*array_y;
468c4762a1bSJed Brown     PetscInt  npoints,p;
469c4762a1bSJed Brown     FILE      *fp = NULL;
470c4762a1bSJed Brown     char      name[PETSC_MAX_PATH_LEN];
471c4762a1bSJed Brown 
4729566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank));
473c4762a1bSJed Brown     fp = fopen(name,"w");
47428b400f6SJacob Faibussowitsch     PetscCheck(fp,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
4759566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(dms,&npoints));
4769566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x));
4779566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y));
478c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
479c4762a1bSJed Brown       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
480c4762a1bSJed Brown     }
4819566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y));
4829566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x));
483c4762a1bSJed Brown     fclose(fp);
484c4762a1bSJed Brown   }
4859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmcell));
4869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
4879566063dSJacob Faibussowitsch   PetscCall(PetscFree(zone));
488c4762a1bSJed Brown   PetscFunctionReturn(0);
489c4762a1bSJed Brown }
490c4762a1bSJed Brown 
491c4762a1bSJed Brown int main(int argc,char **argv)
492c4762a1bSJed Brown {
493c4762a1bSJed Brown   PetscInt       test_mode = 4;
494c4762a1bSJed Brown 
4959566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
4969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL));
497c4762a1bSJed Brown   if (test_mode == 1) {
4989566063dSJacob Faibussowitsch     PetscCall(ex1_1());
499c4762a1bSJed Brown   } else if (test_mode == 2) {
5009566063dSJacob Faibussowitsch     PetscCall(ex1_2());
501c4762a1bSJed Brown   } else if (test_mode == 3) {
5029566063dSJacob Faibussowitsch     PetscCall(ex1_3());
503c4762a1bSJed Brown   } else if (test_mode == 4) {
5049566063dSJacob Faibussowitsch     PetscCall(ex1_4());
505c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown test_mode value, should be 1,2,3,4");
5069566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
507b122ec5aSJacob Faibussowitsch   return 0;
508c4762a1bSJed Brown }
509c4762a1bSJed Brown 
510c4762a1bSJed Brown /*TEST
511c4762a1bSJed Brown 
512c4762a1bSJed Brown    build:
513c4762a1bSJed Brown       requires: !complex double
514c4762a1bSJed Brown 
515c4762a1bSJed Brown    test:
516c4762a1bSJed Brown       args: -test_mode 1
517c4762a1bSJed Brown       filter: grep -v atomic
518c4762a1bSJed Brown       filter_output: grep -v atomic
519c4762a1bSJed Brown 
520c4762a1bSJed Brown    test:
521c4762a1bSJed Brown       suffix: 2
522c4762a1bSJed Brown       args: -test_mode 2
523c4762a1bSJed Brown       filter: grep -v atomic
524c4762a1bSJed Brown       filter_output: grep -v atomic
525c4762a1bSJed Brown 
526c4762a1bSJed Brown    test:
527c4762a1bSJed Brown       suffix: 3
528c4762a1bSJed Brown       args: -test_mode 3
529c4762a1bSJed Brown       filter: grep -v atomic
530c4762a1bSJed Brown       filter_output: grep -v atomic
531c4762a1bSJed Brown       TODO: broken
532c4762a1bSJed Brown 
533c4762a1bSJed Brown    test:
534c4762a1bSJed Brown       suffix: 4
535c4762a1bSJed Brown       args: -test_mode 4
536c4762a1bSJed Brown       filter: grep -v atomic
537c4762a1bSJed Brown       filter_output: grep -v atomic
538c4762a1bSJed Brown 
539c4762a1bSJed Brown    test:
540c4762a1bSJed Brown       suffix: 5
541c4762a1bSJed Brown       nsize: 4
542c4762a1bSJed Brown       args: -test_mode 1
543c4762a1bSJed Brown       filter: grep -v atomic
544c4762a1bSJed Brown       filter_output: grep -v atomic
545c4762a1bSJed Brown 
546c4762a1bSJed Brown    test:
547c4762a1bSJed Brown       suffix: 6
548c4762a1bSJed Brown       nsize: 4
549c4762a1bSJed Brown       args: -test_mode 2
550c4762a1bSJed Brown       filter: grep -v atomic
551c4762a1bSJed Brown       filter_output: grep -v atomic
552c4762a1bSJed Brown 
553c4762a1bSJed Brown    test:
554c4762a1bSJed Brown       suffix: 7
555c4762a1bSJed Brown       nsize: 4
556c4762a1bSJed Brown       args: -test_mode 3
557c4762a1bSJed Brown       filter: grep -v atomic
558c4762a1bSJed Brown       filter_output: grep -v atomic
559c4762a1bSJed Brown       TODO: broken
560c4762a1bSJed Brown 
561c4762a1bSJed Brown    test:
562c4762a1bSJed Brown       suffix: 8
563c4762a1bSJed Brown       nsize: 4
564c4762a1bSJed Brown       args: -test_mode 4
565c4762a1bSJed Brown       filter: grep -v atomic
566c4762a1bSJed Brown       filter_output: grep -v atomic
567c4762a1bSJed Brown 
568c4762a1bSJed Brown TEST*/
569