Actual source code: swarm_ex1.c

petsc-3.12.5 2020-03-29
Report Typos and Errors

  2: static char help[] = "Tests DMSwarm\n\n";

  4:  #include <petscdm.h>
  5:  #include <petscdmda.h>
  6:  #include <petscdmswarm.h>

  8: PETSC_EXTERN PetscErrorCode DMSwarmCollect_General(DM,PetscErrorCode (*)(DM,void*,PetscInt*,PetscInt**),size_t,void*,PetscInt*);
  9: PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM,PetscInt*);

 11: PetscErrorCode ex1_1(void)
 12: {
 13:   DM             dms;
 15:   Vec            x;
 16:   PetscMPIInt    rank,size;
 17:   PetscInt       p;

 19:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 20:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 21:   if ((size > 1) && (size != 4)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must be run wuth 4 MPI ranks");

 23:   DMCreate(PETSC_COMM_WORLD,&dms);
 24:   DMSetType(dms,DMSWARM);

 26:   DMSwarmInitializeFieldRegister(dms);
 27:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
 28:   DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);
 29:   DMSwarmFinalizeFieldRegister(dms);
 30:   DMSwarmSetLocalSizes(dms,5+rank,4);
 31:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);

 33:   {
 34:     PetscReal *array;
 35:     DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);
 36:     for (p=0; p<5+rank; p++) {
 37:       array[p] = 11.1 + p*0.1 + rank*100.0;
 38:     }
 39:     DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);
 40:   }

 42:   {
 43:     PetscReal *array;
 44:     DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);
 45:     for (p=0; p<5+rank; p++) {
 46:       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
 47:     }
 48:     DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);
 49:   }

 51:   DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
 52:   DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);

 54:   DMSwarmVectorDefineField(dms,"strain");
 55:   DMCreateGlobalVector(dms,&x);
 56:   VecDestroy(&x);

 58:   {
 59:     PetscInt    *rankval;
 60:     PetscInt    npoints[2],npoints_orig[2];

 62:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
 63:     DMSwarmGetSize(dms,&npoints_orig[1]);
 64:     DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);
 65:     if ((rank == 0) && (size > 1)) {
 66:       rankval[0] = 1;
 67:       rankval[3] = 1;
 68:     }
 69:     if (rank == 3) {
 70:       rankval[2] = 1;
 71:     }
 72:     DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);
 73:     DMSwarmMigrate(dms,PETSC_TRUE);
 74:     DMSwarmGetLocalSize(dms,&npoints[0]);
 75:     DMSwarmGetSize(dms,&npoints[1]);
 76:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
 77:     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]);
 78:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
 79:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
 80:   }
 81:   {
 82:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
 83:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 84:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);
 85:   }
 86:   {
 87:     DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);
 88:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 89:     DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);
 90:   }

 92:   DMDestroy(&dms);
 93:   return(0);
 94: }

 96: PetscErrorCode ex1_2(void)
 97: {
 98:   DM             dms;
100:   Vec            x;
101:   PetscMPIInt    rank,size;
102:   PetscInt       p;

104:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
105:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
106:   DMCreate(PETSC_COMM_WORLD,&dms);
107:   DMSetType(dms,DMSWARM);
108:   DMSwarmInitializeFieldRegister(dms);

110:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
111:   DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);
112:   DMSwarmFinalizeFieldRegister(dms);
113:   DMSwarmSetLocalSizes(dms,5+rank,4);
114:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);
115:   {
116:     PetscReal *array;
117:     DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);
118:     for (p=0; p<5+rank; p++) {
119:       array[p] = 11.1 + p*0.1 + rank*100.0;
120:     }
121:     DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);
122:   }
123:   {
124:     PetscReal *array;
125:     DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);
126:     for (p=0; p<5+rank; p++) {
127:       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
128:     }
129:     DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);
130:   }
131:   {
132:     PetscInt    *rankval;
133:     PetscInt    npoints[2],npoints_orig[2];

135:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
136:     DMSwarmGetSize(dms,&npoints_orig[1]);
137:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
138:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);
139:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
140:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);

142:     DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);

144:     if (rank == 1) {
145:       rankval[0] = -1;
146:     }
147:     if (rank == 2) {
148:       rankval[1] = -1;
149:     }
150:     if (rank == 3) {
151:       rankval[3] = -1;
152:       rankval[4] = -1;
153:     }
154:     DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);
155:     DMSwarmCollectViewCreate(dms);
156:     DMSwarmGetLocalSize(dms,&npoints[0]);
157:     DMSwarmGetSize(dms,&npoints[1]);
158:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
159:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);
160:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
161:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);

163:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
164:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
165:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);

167:     DMSwarmCollectViewDestroy(dms);
168:     DMSwarmGetLocalSize(dms,&npoints[0]);
169:     DMSwarmGetSize(dms,&npoints[1]);
170:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
171:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after_v(%D,%D)\n",rank,npoints[0],npoints[1]);
172:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
173:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);

175:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
176:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
177:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);
178:   }
179:   DMDestroy(&dms);
180:   return(0);
181: }

183: /*
184:  splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp"
185: */
186: PetscErrorCode ex1_3(void)
187: {
188:   DM             dms;
190:   PetscMPIInt    rank,size;
191:   PetscInt       is,js,ni,nj,overlap;
192:   DM             dmcell;

194:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
195:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
196:   overlap = 2;
197:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,13,13,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);
198:   DMSetFromOptions(dmcell);
199:   DMSetUp(dmcell);
200:   DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);
201:   DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);
202:   DMCreate(PETSC_COMM_WORLD,&dms);
203:   DMSetType(dms,DMSWARM);
204:   DMSwarmSetCellDM(dms,dmcell);

206:   /* load in data types */
207:   DMSwarmInitializeFieldRegister(dms);
208:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
209:   DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);
210:   DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);
211:   DMSwarmFinalizeFieldRegister(dms);
212:   DMSwarmSetLocalSizes(dms,ni*nj*4,4);
213:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);

215:   /* set values within the swarm */
216:   {
217:     PetscReal  *array_x,*array_y;
218:     PetscInt   npoints,i,j,cnt;
219:     DMDACoor2d **LA_coor;
220:     Vec        coor;
221:     DM         dmcellcdm;

223:     DMGetCoordinateDM(dmcell,&dmcellcdm);
224:     DMGetCoordinates(dmcell,&coor);
225:     DMDAVecGetArray(dmcellcdm,coor,&LA_coor);
226:     DMSwarmGetLocalSize(dms,&npoints);
227:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
228:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
229:     cnt = 0;
230:     for (j=js; j<js+nj; j++) {
231:       for (i=is; i<is+ni; i++) {
232:         PetscReal xp,yp;
233:         xp = PetscRealPart(LA_coor[j][i].x);
234:         yp = PetscRealPart(LA_coor[j][i].y);
235:         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; }
236:         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; }
237:         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; }
238:         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; }

240:         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; }
241:         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; }
242:         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; }
243:         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; }
244:         cnt++;
245:       }
246:     }
247:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
248:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
249:     DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);
250:   }
251:   {
252:     PetscInt    npoints[2],npoints_orig[2],ng;

254:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
255:     DMSwarmGetSize(dms,&npoints_orig[1]);
256:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
257:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);
258:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
259:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
260:     DMSwarmCollect_DMDABoundingBox(dms,&ng);

262:     DMSwarmGetLocalSize(dms,&npoints[0]);
263:     DMSwarmGetSize(dms,&npoints[1]);
264:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
265:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);
266:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
267:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
268:   }
269:   {
270:     PetscReal *array_x,*array_y;
271:     PetscInt  npoints,p;
272:     FILE      *fp = NULL;
273:     char      name[PETSC_MAX_PATH_LEN];

275:     PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);
276:     fp = fopen(name,"w");
277:     if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
278:     DMSwarmGetLocalSize(dms,&npoints);
279:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
280:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
281:     for (p=0; p<npoints; p++) {
282:       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
283:     }
284:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
285:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
286:     fclose(fp);
287:   }
288:   DMDestroy(&dmcell);
289:   DMDestroy(&dms);
290:   return(0);
291: }

293: typedef struct {
294:   PetscReal cx[2];
295:   PetscReal radius;
296: } CollectZoneCtx;

298: PetscErrorCode collect_zone(DM dm,void *ctx,PetscInt *nfound,PetscInt **foundlist)
299: {
300:   CollectZoneCtx *zone = (CollectZoneCtx*)ctx;
301:   PetscInt       p,npoints;
302:   PetscReal      *array_x,*array_y,r2;
303:   PetscInt       p2collect,*plist;
304:   PetscMPIInt    rank;

307:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
308:   DMSwarmGetLocalSize(dm,&npoints);
309:   DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);
310:   DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);

312:   r2 = zone->radius * zone->radius;
313:   p2collect = 0;
314:   for (p=0; p<npoints; p++) {
315:     PetscReal sep2;

317:     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
318:     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
319:     if (sep2 < r2) {
320:       p2collect++;
321:     }
322:   }

324:   PetscMalloc1(p2collect+1,&plist);
325:   p2collect = 0;
326:   for (p=0; p<npoints; p++) {
327:     PetscReal sep2;

329:     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
330:     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
331:     if (sep2 < r2) {
332:       plist[p2collect] = p;
333:       p2collect++;
334:     }
335:   }
336:   DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);
337:   DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);

339:   *nfound = p2collect;
340:   *foundlist = plist;
341:   return(0);
342: }

344: PetscErrorCode ex1_4(void)
345: {
346:   DM             dms;
348:   PetscMPIInt    rank,size;
349:   PetscInt       is,js,ni,nj,overlap,nn;
350:   DM             dmcell;
351:   CollectZoneCtx *zone;
352:   PetscReal      dx;

354:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
355:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
356:   nn = 101;
357:   dx = 2.0/ (PetscReal)(nn-1);
358:   overlap = 0;
359:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nn,nn,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);
360:   DMSetFromOptions(dmcell);
361:   DMSetUp(dmcell);
362:   DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);
363:   DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);
364:   DMCreate(PETSC_COMM_WORLD,&dms);
365:   DMSetType(dms,DMSWARM);

367:   /* load in data types */
368:   DMSwarmInitializeFieldRegister(dms);
369:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
370:   DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);
371:   DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);
372:   DMSwarmFinalizeFieldRegister(dms);
373:   DMSwarmSetLocalSizes(dms,ni*nj*4,4);
374:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);

376:   /* set values within the swarm */
377:   {
378:     PetscReal  *array_x,*array_y;
379:     PetscInt   npoints,i,j,cnt;
380:     DMDACoor2d **LA_coor;
381:     Vec        coor;
382:     DM         dmcellcdm;

384:     DMGetCoordinateDM(dmcell,&dmcellcdm);
385:     DMGetCoordinates(dmcell,&coor);
386:     DMDAVecGetArray(dmcellcdm,coor,&LA_coor);
387:     DMSwarmGetLocalSize(dms,&npoints);
388:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
389:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
390:     cnt = 0;
391:     for (j=js; j<js+nj; j++) {
392:       for (i=is; i<is+ni; i++) {
393:         PetscReal xp,yp;

395:         xp = PetscRealPart(LA_coor[j][i].x);
396:         yp = PetscRealPart(LA_coor[j][i].y);
397:         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; }*/
398:         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; }*/
399:         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; }*/
400:         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; }*/
401:         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; }*/
402:         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; }*/
403:         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; }*/
404:         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; }*/
405:         cnt++;
406:       }
407:     }
408:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
409:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
410:     DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);
411:   }
412:   PetscMalloc1(1,&zone);
413:   if (size == 4) {
414:     if (rank == 0) {
415:       zone->cx[0] = 0.5;
416:       zone->cx[1] = 0.5;
417:       zone->radius = 0.3;
418:     }
419:     if (rank == 1) {
420:       zone->cx[0] = -0.5;
421:       zone->cx[1] = 0.5;
422:       zone->radius = 0.25;
423:     }
424:     if (rank == 2) {
425:       zone->cx[0] = 0.5;
426:       zone->cx[1] = -0.5;
427:       zone->radius = 0.2;
428:     }
429:     if (rank == 3) {
430:       zone->cx[0] = -0.5;
431:       zone->cx[1] = -0.5;
432:       zone->radius = 0.1;
433:     }
434:   } else {
435:     if (rank == 0) {
436:       zone->cx[0] = 0.5;
437:       zone->cx[1] = 0.5;
438:       zone->radius = 0.8;
439:     } else {
440:       zone->cx[0] = 10.0;
441:       zone->cx[1] = 10.0;
442:       zone->radius = 0.0;
443:     }
444:   }
445:   {
446:     PetscInt    npoints[2],npoints_orig[2],ng;

448:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
449:     DMSwarmGetSize(dms,&npoints_orig[1]);
450:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
451:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);
452:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
453:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
454:     DMSwarmCollect_General(dms,collect_zone,sizeof(CollectZoneCtx),zone,&ng);
455:     DMSwarmGetLocalSize(dms,&npoints[0]);
456:     DMSwarmGetSize(dms,&npoints[1]);
457:     PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
458:     PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);
459:     PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
460:     PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
461:   }
462:   {
463:     PetscReal *array_x,*array_y;
464:     PetscInt  npoints,p;
465:     FILE      *fp = NULL;
466:     char      name[PETSC_MAX_PATH_LEN];

468:     PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);
469:     fp = fopen(name,"w");
470:     if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
471:     DMSwarmGetLocalSize(dms,&npoints);
472:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
473:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
474:     for (p=0; p<npoints; p++) {
475:       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
476:     }
477:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
478:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
479:     fclose(fp);
480:   }
481:   DMDestroy(&dmcell);
482:   DMDestroy(&dms);
483:   PetscFree(zone);
484:   return(0);
485: }

487: int main(int argc,char **argv)
488: {
490:   PetscInt       test_mode = 4;

492:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
493:   PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL);
494:   if (test_mode == 1) {
495:     ex1_1();
496:   } else if (test_mode == 2) {
497:     ex1_2();
498:   } else if (test_mode == 3) {
499:     ex1_3();
500:   } else if (test_mode == 4) {
501:     ex1_4();
502:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown test_mode value, should be 1,2,3,4");
503:   PetscFinalize();
504:   return ierr;
505: }

507: /*TEST

509:    build:
510:       requires: !complex double

512:    test:
513:       args: -test_mode 1
514:       filter: grep -v atomic
515:       filter_output: grep -v atomic

517:    test:
518:       suffix: 2
519:       args: -test_mode 2
520:       filter: grep -v atomic
521:       filter_output: grep -v atomic

523:    test:
524:       suffix: 3
525:       args: -test_mode 3
526:       filter: grep -v atomic
527:       filter_output: grep -v atomic
528:       TODO: broken

530:    test:
531:       suffix: 4
532:       args: -test_mode 4
533:       filter: grep -v atomic
534:       filter_output: grep -v atomic

536:    test:
537:       suffix: 5
538:       nsize: 4
539:       args: -test_mode 1
540:       filter: grep -v atomic
541:       filter_output: grep -v atomic

543:    test:
544:       suffix: 6
545:       nsize: 4
546:       args: -test_mode 2
547:       filter: grep -v atomic
548:       filter_output: grep -v atomic

550:    test:
551:       suffix: 7
552:       nsize: 4
553:       args: -test_mode 3
554:       filter: grep -v atomic
555:       filter_output: grep -v atomic
556:       TODO: broken

558:    test:
559:       suffix: 8
560:       nsize: 4
561:       args: -test_mode 4
562:       filter: grep -v atomic
563:       filter_output: grep -v atomic

565: TEST*/