Actual source code: swarm_ex1.c


  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;
 14:   Vec            x;
 15:   PetscMPIInt    rank,size;
 16:   PetscInt       p;

 18:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

161:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
162:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
163:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);

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

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

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

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

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

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

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

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

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

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

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

290: typedef struct {
291:   PetscReal cx[2];
292:   PetscReal radius;
293: } CollectZoneCtx;

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

303:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
304:   DMSwarmGetLocalSize(dm,&npoints);
305:   DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);
306:   DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);

308:   r2 = zone->radius * zone->radius;
309:   p2collect = 0;
310:   for (p=0; p<npoints; p++) {
311:     PetscReal sep2;

313:     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
314:     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
315:     if (sep2 < r2) {
316:       p2collect++;
317:     }
318:   }

320:   PetscMalloc1(p2collect+1,&plist);
321:   p2collect = 0;
322:   for (p=0; p<npoints; p++) {
323:     PetscReal sep2;

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

335:   *nfound = p2collect;
336:   *foundlist = plist;
337:   return 0;
338: }

340: PetscErrorCode ex1_4(void)
341: {
342:   DM             dms;
343:   PetscMPIInt    rank,size;
344:   PetscInt       is,js,ni,nj,overlap,nn;
345:   DM             dmcell;
346:   CollectZoneCtx *zone;
347:   PetscReal      dx;

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

362:   /* load in data types */
363:   DMSwarmInitializeFieldRegister(dms);
364:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
365:   DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);
366:   DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);
367:   DMSwarmFinalizeFieldRegister(dms);
368:   DMSwarmSetLocalSizes(dms,ni*nj*4,4);
369:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);

371:   /* set values within the swarm */
372:   {
373:     PetscReal  *array_x,*array_y;
374:     PetscInt   npoints,i,j,cnt;
375:     DMDACoor2d **LA_coor;
376:     Vec        coor;
377:     DM         dmcellcdm;

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

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

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

463:     PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);
464:     fp = fopen(name,"w");
466:     DMSwarmGetLocalSize(dms,&npoints);
467:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
468:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
469:     for (p=0; p<npoints; p++) {
470:       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
471:     }
472:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
473:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
474:     fclose(fp);
475:   }
476:   DMDestroy(&dmcell);
477:   DMDestroy(&dms);
478:   PetscFree(zone);
479:   return 0;
480: }

482: int main(int argc,char **argv)
483: {
484:   PetscInt       test_mode = 4;

486:   PetscInitialize(&argc,&argv,(char*)0,help);
487:   PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL);
488:   if (test_mode == 1) {
489:     ex1_1();
490:   } else if (test_mode == 2) {
491:     ex1_2();
492:   } else if (test_mode == 3) {
493:     ex1_3();
494:   } else if (test_mode == 4) {
495:     ex1_4();
496:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown test_mode value, should be 1,2,3,4");
497:   PetscFinalize();
498:   return 0;
499: }

501: /*TEST

503:    build:
504:       requires: !complex double

506:    test:
507:       args: -test_mode 1
508:       filter: grep -v atomic
509:       filter_output: grep -v atomic

511:    test:
512:       suffix: 2
513:       args: -test_mode 2
514:       filter: grep -v atomic
515:       filter_output: grep -v atomic

517:    test:
518:       suffix: 3
519:       args: -test_mode 3
520:       filter: grep -v atomic
521:       filter_output: grep -v atomic
522:       TODO: broken

524:    test:
525:       suffix: 4
526:       args: -test_mode 4
527:       filter: grep -v atomic
528:       filter_output: grep -v atomic

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

537:    test:
538:       suffix: 6
539:       nsize: 4
540:       args: -test_mode 2
541:       filter: grep -v atomic
542:       filter_output: grep -v atomic

544:    test:
545:       suffix: 7
546:       nsize: 4
547:       args: -test_mode 3
548:       filter: grep -v atomic
549:       filter_output: grep -v atomic
550:       TODO: broken

552:    test:
553:       suffix: 8
554:       nsize: 4
555:       args: -test_mode 4
556:       filter: grep -v atomic
557:       filter_output: grep -v atomic

559: TEST*/