xref: /petsc/src/dm/tutorials/swarm_ex1.c (revision 362febeeeb69b91ebadcb4b2dc0a22cb6dfc4097)
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   PetscErrorCode ierr;
15c4762a1bSJed Brown   Vec            x;
16c4762a1bSJed Brown   PetscMPIInt    rank,size;
17c4762a1bSJed Brown   PetscInt       p;
18c4762a1bSJed Brown 
19*362febeeSStefano Zampini   PetscFunctionBegin;
20ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
21ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
22c4762a1bSJed Brown   if ((size > 1) && (size != 4)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must be run wuth 4 MPI ranks");
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr);
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(dms,5+rank,4);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   {
35c4762a1bSJed Brown     PetscReal *array;
36c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr);
37c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
38c4762a1bSJed Brown       array[p] = 11.1 + p*0.1 + rank*100.0;
39c4762a1bSJed Brown     }
40c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr);
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   {
44c4762a1bSJed Brown     PetscReal *array;
45c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr);
46c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
47c4762a1bSJed Brown       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr);
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   ierr = DMSwarmVectorDefineField(dms,"strain");CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dms,&x);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   {
60c4762a1bSJed Brown     PetscInt    *rankval;
61c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2];
62c4762a1bSJed Brown 
63c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr);
64c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr);
65c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
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     }
73c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
74c4762a1bSJed Brown     ierr = DMSwarmMigrate(dms,PETSC_TRUE);CHKERRQ(ierr);
75c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr);
76c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr);
77c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
78c4762a1bSJed Brown     ierr = 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]);CHKERRQ(ierr);
79c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
80c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
81c4762a1bSJed Brown   }
82c4762a1bSJed Brown   {
83c4762a1bSJed Brown     ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
84c4762a1bSJed Brown     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
85c4762a1bSJed Brown     ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   {
88c4762a1bSJed Brown     ierr = DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr);
89c4762a1bSJed Brown     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
90c4762a1bSJed Brown     ierr = DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr);
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   ierr = DMDestroy(&dms);CHKERRQ(ierr);
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown PetscErrorCode ex1_2(void)
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   DM             dms;
100c4762a1bSJed Brown   PetscErrorCode ierr;
101c4762a1bSJed Brown   Vec            x;
102c4762a1bSJed Brown   PetscMPIInt    rank,size;
103c4762a1bSJed Brown   PetscInt       p;
104c4762a1bSJed Brown 
105*362febeeSStefano Zampini   PetscFunctionBegin;
106ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
107ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
108c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr);
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(dms,5+rank,4);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
117c4762a1bSJed Brown   {
118c4762a1bSJed Brown     PetscReal *array;
119c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr);
120c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
121c4762a1bSJed Brown       array[p] = 11.1 + p*0.1 + rank*100.0;
122c4762a1bSJed Brown     }
123c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);CHKERRQ(ierr);
124c4762a1bSJed Brown   }
125c4762a1bSJed Brown   {
126c4762a1bSJed Brown     PetscReal *array;
127c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr);
128c4762a1bSJed Brown     for (p=0; p<5+rank; p++) {
129c4762a1bSJed Brown       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
130c4762a1bSJed Brown     }
131c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);CHKERRQ(ierr);
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown   {
134c4762a1bSJed Brown     PetscInt    *rankval;
135c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2];
136c4762a1bSJed Brown 
137c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr);
138c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr);
139c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
140c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);CHKERRQ(ierr);
141c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
142c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
143c4762a1bSJed Brown 
144c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
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     }
156c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
157c4762a1bSJed Brown     ierr = DMSwarmCollectViewCreate(dms);CHKERRQ(ierr);
158c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr);
159c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr);
160c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
161c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr);
162c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
163c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
164c4762a1bSJed Brown 
165c4762a1bSJed Brown     ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
166c4762a1bSJed Brown     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
167c4762a1bSJed Brown     ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
168c4762a1bSJed Brown 
169c4762a1bSJed Brown     ierr = DMSwarmCollectViewDestroy(dms);CHKERRQ(ierr);
170c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr);
171c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr);
172c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
173c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after_v(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr);
174c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
175c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
176c4762a1bSJed Brown 
177c4762a1bSJed Brown     ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
178c4762a1bSJed Brown     ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
179c4762a1bSJed Brown     ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown   ierr = DMDestroy(&dms);CHKERRQ(ierr);
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   PetscErrorCode ierr;
192c4762a1bSJed Brown   PetscMPIInt    rank,size;
193c4762a1bSJed Brown   PetscInt       is,js,ni,nj,overlap;
194c4762a1bSJed Brown   DM             dmcell;
195c4762a1bSJed Brown 
196*362febeeSStefano Zampini   PetscFunctionBegin;
197ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
198ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
199c4762a1bSJed Brown   overlap = 2;
200c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,13,13,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = DMSetFromOptions(dmcell);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = DMSetUp(dmcell);CHKERRQ(ierr);
203c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);CHKERRQ(ierr);
204c4762a1bSJed Brown   ierr = DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);CHKERRQ(ierr);
205c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr);
206c4762a1bSJed Brown   ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr);
207c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(dms,dmcell);CHKERRQ(ierr);
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* load in data types */
210c4762a1bSJed Brown   ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr);
211c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr);
212c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);CHKERRQ(ierr);
214c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr);
215c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(dms,ni*nj*4,4);CHKERRQ(ierr);
216c4762a1bSJed Brown   ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
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 
226c4762a1bSJed Brown     ierr = DMGetCoordinateDM(dmcell,&dmcellcdm);CHKERRQ(ierr);
227c4762a1bSJed Brown     ierr = DMGetCoordinates(dmcell,&coor);CHKERRQ(ierr);
228c4762a1bSJed Brown     ierr = DMDAVecGetArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr);
229c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr);
230c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
231c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
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     }
250c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
251c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
252c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr);
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown   {
255c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2],ng;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr);
258c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr);
259c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
260c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);CHKERRQ(ierr);
261c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
262c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
263c4762a1bSJed Brown     ierr = DMSwarmCollect_DMDABoundingBox(dms,&ng);CHKERRQ(ierr);
264c4762a1bSJed Brown 
265c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr);
266c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr);
267c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
268c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr);
269c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
270c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
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 
278c4762a1bSJed Brown     ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);CHKERRQ(ierr);
279c4762a1bSJed Brown     fp = fopen(name,"w");
280c4762a1bSJed Brown     if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
281c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr);
282c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
283c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
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     }
287c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
288c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
289c4762a1bSJed Brown     fclose(fp);
290c4762a1bSJed Brown   }
291c4762a1bSJed Brown   ierr = DMDestroy(&dmcell);CHKERRQ(ierr);
292c4762a1bSJed Brown   ierr = DMDestroy(&dms);CHKERRQ(ierr);
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   PetscErrorCode ierr;
309c4762a1bSJed Brown 
310*362febeeSStefano Zampini   PetscFunctionBegin;
311ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
312c4762a1bSJed Brown   ierr = DMSwarmGetLocalSize(dm,&npoints);CHKERRQ(ierr);
313c4762a1bSJed Brown   ierr = DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
314c4762a1bSJed Brown   ierr = DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   r2 = zone->radius * zone->radius;
317c4762a1bSJed Brown   p2collect = 0;
318c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
319c4762a1bSJed Brown     PetscReal sep2;
320c4762a1bSJed Brown 
321c4762a1bSJed Brown     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
322c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
323c4762a1bSJed Brown     if (sep2 < r2) {
324c4762a1bSJed Brown       p2collect++;
325c4762a1bSJed Brown     }
326c4762a1bSJed Brown   }
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   ierr = PetscMalloc1(p2collect+1,&plist);CHKERRQ(ierr);
329c4762a1bSJed Brown   p2collect = 0;
330c4762a1bSJed Brown   for (p=0; p<npoints; p++) {
331c4762a1bSJed Brown     PetscReal sep2;
332c4762a1bSJed Brown 
333c4762a1bSJed Brown     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
334c4762a1bSJed Brown     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
335c4762a1bSJed Brown     if (sep2 < r2) {
336c4762a1bSJed Brown       plist[p2collect] = p;
337c4762a1bSJed Brown       p2collect++;
338c4762a1bSJed Brown     }
339c4762a1bSJed Brown   }
340c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
341c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   *nfound = p2collect;
344c4762a1bSJed Brown   *foundlist = plist;
345c4762a1bSJed Brown   PetscFunctionReturn(0);
346c4762a1bSJed Brown }
347c4762a1bSJed Brown 
348c4762a1bSJed Brown PetscErrorCode ex1_4(void)
349c4762a1bSJed Brown {
350c4762a1bSJed Brown   DM             dms;
351c4762a1bSJed Brown   PetscErrorCode ierr;
352c4762a1bSJed Brown   PetscMPIInt    rank,size;
353c4762a1bSJed Brown   PetscInt       is,js,ni,nj,overlap,nn;
354c4762a1bSJed Brown   DM             dmcell;
355c4762a1bSJed Brown   CollectZoneCtx *zone;
356c4762a1bSJed Brown   PetscReal      dx;
357c4762a1bSJed Brown 
358*362febeeSStefano Zampini   PetscFunctionBegin;
359ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
360ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
361c4762a1bSJed Brown   nn = 101;
362c4762a1bSJed Brown   dx = 2.0/ (PetscReal)(nn-1);
363c4762a1bSJed Brown   overlap = 0;
364c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nn,nn,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);CHKERRQ(ierr);
365c4762a1bSJed Brown   ierr = DMSetFromOptions(dmcell);CHKERRQ(ierr);
366c4762a1bSJed Brown   ierr = DMSetUp(dmcell);CHKERRQ(ierr);
367c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);CHKERRQ(ierr);
368c4762a1bSJed Brown   ierr = DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr);
370c4762a1bSJed Brown   ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr);
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /* load in data types */
373c4762a1bSJed Brown   ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr);
374c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr);
375c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);CHKERRQ(ierr);
377c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr);
378c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(dms,ni*nj*4,4);CHKERRQ(ierr);
379c4762a1bSJed Brown   ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
380c4762a1bSJed Brown 
381c4762a1bSJed Brown   /* set values within the swarm */
382c4762a1bSJed Brown   {
383c4762a1bSJed Brown     PetscReal  *array_x,*array_y;
384c4762a1bSJed Brown     PetscInt   npoints,i,j,cnt;
385c4762a1bSJed Brown     DMDACoor2d **LA_coor;
386c4762a1bSJed Brown     Vec        coor;
387c4762a1bSJed Brown     DM         dmcellcdm;
388c4762a1bSJed Brown 
389c4762a1bSJed Brown     ierr = DMGetCoordinateDM(dmcell,&dmcellcdm);CHKERRQ(ierr);
390c4762a1bSJed Brown     ierr = DMGetCoordinates(dmcell,&coor);CHKERRQ(ierr);
391c4762a1bSJed Brown     ierr = DMDAVecGetArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr);
392c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr);
393c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
394c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
395c4762a1bSJed Brown     cnt = 0;
396c4762a1bSJed Brown     for (j=js; j<js+nj; j++) {
397c4762a1bSJed Brown       for (i=is; i<is+ni; i++) {
398c4762a1bSJed Brown         PetscReal xp,yp;
399c4762a1bSJed Brown 
400c4762a1bSJed Brown         xp = PetscRealPart(LA_coor[j][i].x);
401c4762a1bSJed Brown         yp = PetscRealPart(LA_coor[j][i].y);
402c4762a1bSJed 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; }*/
403c4762a1bSJed 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; }*/
404c4762a1bSJed 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; }*/
405c4762a1bSJed 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; }*/
406c4762a1bSJed 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; }*/
407c4762a1bSJed 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; }*/
408c4762a1bSJed 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; }*/
409c4762a1bSJed 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; }*/
410c4762a1bSJed Brown         cnt++;
411c4762a1bSJed Brown       }
412c4762a1bSJed Brown     }
413c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
414c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
415c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);CHKERRQ(ierr);
416c4762a1bSJed Brown   }
417c4762a1bSJed Brown   ierr = PetscMalloc1(1,&zone);CHKERRQ(ierr);
418c4762a1bSJed Brown   if (size == 4) {
419c4762a1bSJed Brown     if (rank == 0) {
420c4762a1bSJed Brown       zone->cx[0] = 0.5;
421c4762a1bSJed Brown       zone->cx[1] = 0.5;
422c4762a1bSJed Brown       zone->radius = 0.3;
423c4762a1bSJed Brown     }
424c4762a1bSJed Brown     if (rank == 1) {
425c4762a1bSJed Brown       zone->cx[0] = -0.5;
426c4762a1bSJed Brown       zone->cx[1] = 0.5;
427c4762a1bSJed Brown       zone->radius = 0.25;
428c4762a1bSJed Brown     }
429c4762a1bSJed Brown     if (rank == 2) {
430c4762a1bSJed Brown       zone->cx[0] = 0.5;
431c4762a1bSJed Brown       zone->cx[1] = -0.5;
432c4762a1bSJed Brown       zone->radius = 0.2;
433c4762a1bSJed Brown     }
434c4762a1bSJed Brown     if (rank == 3) {
435c4762a1bSJed Brown       zone->cx[0] = -0.5;
436c4762a1bSJed Brown       zone->cx[1] = -0.5;
437c4762a1bSJed Brown       zone->radius = 0.1;
438c4762a1bSJed Brown     }
439c4762a1bSJed Brown   } else {
440c4762a1bSJed Brown     if (rank == 0) {
441c4762a1bSJed Brown       zone->cx[0] = 0.5;
442c4762a1bSJed Brown       zone->cx[1] = 0.5;
443c4762a1bSJed Brown       zone->radius = 0.8;
444c4762a1bSJed Brown     } else {
445c4762a1bSJed Brown       zone->cx[0] = 10.0;
446c4762a1bSJed Brown       zone->cx[1] = 10.0;
447c4762a1bSJed Brown       zone->radius = 0.0;
448c4762a1bSJed Brown     }
449c4762a1bSJed Brown   }
450c4762a1bSJed Brown   {
451c4762a1bSJed Brown     PetscInt    npoints[2],npoints_orig[2],ng;
452c4762a1bSJed Brown 
453c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints_orig[0]);CHKERRQ(ierr);
454c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints_orig[1]);CHKERRQ(ierr);
455c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
456c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);CHKERRQ(ierr);
457c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
458c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
459c4762a1bSJed Brown     ierr = DMSwarmCollect_General(dms,collect_zone,sizeof(CollectZoneCtx),zone,&ng);CHKERRQ(ierr);
460c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints[0]);CHKERRQ(ierr);
461c4762a1bSJed Brown     ierr = DMSwarmGetSize(dms,&npoints[1]);CHKERRQ(ierr);
462c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
463c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);CHKERRQ(ierr);
464c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
465c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
466c4762a1bSJed Brown   }
467c4762a1bSJed Brown   {
468c4762a1bSJed Brown     PetscReal *array_x,*array_y;
469c4762a1bSJed Brown     PetscInt  npoints,p;
470c4762a1bSJed Brown     FILE      *fp = NULL;
471c4762a1bSJed Brown     char      name[PETSC_MAX_PATH_LEN];
472c4762a1bSJed Brown 
473c4762a1bSJed Brown     ierr = PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);CHKERRQ(ierr);
474c4762a1bSJed Brown     fp = fopen(name,"w");
475c4762a1bSJed Brown     if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
476c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(dms,&npoints);CHKERRQ(ierr);
477c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
478c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
479c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
480c4762a1bSJed Brown       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
481c4762a1bSJed Brown     }
482c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);CHKERRQ(ierr);
483c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);CHKERRQ(ierr);
484c4762a1bSJed Brown     fclose(fp);
485c4762a1bSJed Brown   }
486c4762a1bSJed Brown   ierr = DMDestroy(&dmcell);CHKERRQ(ierr);
487c4762a1bSJed Brown   ierr = DMDestroy(&dms);CHKERRQ(ierr);
488c4762a1bSJed Brown   ierr = PetscFree(zone);CHKERRQ(ierr);
489c4762a1bSJed Brown   PetscFunctionReturn(0);
490c4762a1bSJed Brown }
491c4762a1bSJed Brown 
492c4762a1bSJed Brown int main(int argc,char **argv)
493c4762a1bSJed Brown {
494c4762a1bSJed Brown   PetscErrorCode ierr;
495c4762a1bSJed Brown   PetscInt       test_mode = 4;
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
498c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL);CHKERRQ(ierr);
499c4762a1bSJed Brown   if (test_mode == 1) {
500c4762a1bSJed Brown     ierr = ex1_1();CHKERRQ(ierr);
501c4762a1bSJed Brown   } else if (test_mode == 2) {
502c4762a1bSJed Brown     ierr = ex1_2();CHKERRQ(ierr);
503c4762a1bSJed Brown   } else if (test_mode == 3) {
504c4762a1bSJed Brown     ierr = ex1_3();CHKERRQ(ierr);
505c4762a1bSJed Brown   } else if (test_mode == 4) {
506c4762a1bSJed Brown     ierr = ex1_4();CHKERRQ(ierr);
507c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown test_mode value, should be 1,2,3,4");
508c4762a1bSJed Brown   ierr = PetscFinalize();
509c4762a1bSJed Brown   return ierr;
510c4762a1bSJed Brown }
511c4762a1bSJed Brown 
512c4762a1bSJed Brown /*TEST
513c4762a1bSJed Brown 
514c4762a1bSJed Brown    build:
515c4762a1bSJed Brown       requires: !complex double
516c4762a1bSJed Brown 
517c4762a1bSJed Brown    test:
518c4762a1bSJed Brown       args: -test_mode 1
519c4762a1bSJed Brown       filter: grep -v atomic
520c4762a1bSJed Brown       filter_output: grep -v atomic
521c4762a1bSJed Brown 
522c4762a1bSJed Brown    test:
523c4762a1bSJed Brown       suffix: 2
524c4762a1bSJed Brown       args: -test_mode 2
525c4762a1bSJed Brown       filter: grep -v atomic
526c4762a1bSJed Brown       filter_output: grep -v atomic
527c4762a1bSJed Brown 
528c4762a1bSJed Brown    test:
529c4762a1bSJed Brown       suffix: 3
530c4762a1bSJed Brown       args: -test_mode 3
531c4762a1bSJed Brown       filter: grep -v atomic
532c4762a1bSJed Brown       filter_output: grep -v atomic
533c4762a1bSJed Brown       TODO: broken
534c4762a1bSJed Brown 
535c4762a1bSJed Brown    test:
536c4762a1bSJed Brown       suffix: 4
537c4762a1bSJed Brown       args: -test_mode 4
538c4762a1bSJed Brown       filter: grep -v atomic
539c4762a1bSJed Brown       filter_output: grep -v atomic
540c4762a1bSJed Brown 
541c4762a1bSJed Brown    test:
542c4762a1bSJed Brown       suffix: 5
543c4762a1bSJed Brown       nsize: 4
544c4762a1bSJed Brown       args: -test_mode 1
545c4762a1bSJed Brown       filter: grep -v atomic
546c4762a1bSJed Brown       filter_output: grep -v atomic
547c4762a1bSJed Brown 
548c4762a1bSJed Brown    test:
549c4762a1bSJed Brown       suffix: 6
550c4762a1bSJed Brown       nsize: 4
551c4762a1bSJed Brown       args: -test_mode 2
552c4762a1bSJed Brown       filter: grep -v atomic
553c4762a1bSJed Brown       filter_output: grep -v atomic
554c4762a1bSJed Brown 
555c4762a1bSJed Brown    test:
556c4762a1bSJed Brown       suffix: 7
557c4762a1bSJed Brown       nsize: 4
558c4762a1bSJed Brown       args: -test_mode 3
559c4762a1bSJed Brown       filter: grep -v atomic
560c4762a1bSJed Brown       filter_output: grep -v atomic
561c4762a1bSJed Brown       TODO: broken
562c4762a1bSJed Brown 
563c4762a1bSJed Brown    test:
564c4762a1bSJed Brown       suffix: 8
565c4762a1bSJed Brown       nsize: 4
566c4762a1bSJed Brown       args: -test_mode 4
567c4762a1bSJed Brown       filter: grep -v atomic
568c4762a1bSJed Brown       filter_output: grep -v atomic
569c4762a1bSJed Brown 
570c4762a1bSJed Brown TEST*/
571