Actual source code: swarm_ex1.c

petsc-3.8.4 2018-03-24
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 dm,PetscErrorCode (*collect)(DM,void*,PetscInt*,PetscInt**),size_t ctx_size,void *ctx,PetscInt *globalsize);
  9: PETSC_EXTERN PetscErrorCode DMSwarmCollect_DMDABoundingBox(DM dm,PetscInt *globalsize);

 11: PetscErrorCode ex1_1(void)
 12: {
 13:   DM dms;
 15:   Vec x;
 16:   PetscMPIInt rank,commsize;
 17:   PetscInt p;
 18: 
 19:   MPI_Comm_size(PETSC_COMM_WORLD,&commsize);
 20:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 22:   if ((commsize > 1) && (commsize != 4)) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must be run wuth 4 MPI ranks");
 23: 
 24:   DMCreate(PETSC_COMM_WORLD,&dms);
 25:   DMSetType(dms,DMSWARM);

 27:   DMSwarmInitializeFieldRegister(dms);
 28: 
 29:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
 30:   DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);
 31: 
 32:   DMSwarmFinalizeFieldRegister(dms);
 33: 
 34:   DMSwarmSetLocalSizes(dms,5+rank,4);
 35:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);
 36: 
 37:   {
 38:     PetscReal *array;
 39:     DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);
 40:     for (p=0; p<5+rank; p++) {
 41:       array[p] = 11.1 + p*0.1 + rank*100.0;
 42:     }
 43:     DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);
 44:   }

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

 55:   DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
 56:   /*VecView(x,PETSC_VIEWER_STDOUT_WORLD);*/
 57:   DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);
 58: 
 59:   DMSwarmVectorDefineField(dms,"strain");
 60:   DMCreateGlobalVector(dms,&x);
 61:   /*VecView(x,PETSC_VIEWER_STDOUT_WORLD);*/
 62:   VecDestroy(&x);
 63: 
 64:   {
 65:     PetscInt *rankval;
 66:     PetscInt npoints[2],npoints_orig[2];
 67: 
 68:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
 69:     DMSwarmGetSize(dms,&npoints_orig[1]);
 70:     DMSwarmGetField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);
 71:     if ((rank == 0) && (commsize > 1)) {
 72:       rankval[0] = 1;
 73:       rankval[3] = 1;
 74:     }
 75:     if (rank == 3) {
 76:       rankval[2] = 1;
 77:     }
 78:     DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);
 79: 
 80:     DMSwarmMigrate(dms,PETSC_TRUE);
 81:     DMSwarmGetLocalSize(dms,&npoints[0]);
 82:     DMSwarmGetSize(dms,&npoints[1]);

 84:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] before(%D,%D) after(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1],npoints[0],npoints[1]);
 85:   }
 86: 
 87:   {
 88:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
 89:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 90:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);
 91:   }
 92:   {
 93:     DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);
 94:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
 95:     DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);
 96:   }

 98:   DMDestroy(&dms);
 99: 
100:   return(0);
101: }

103: PetscErrorCode ex1_2(void)
104: {
105:   DM dms;
107:   Vec x;
108:   PetscMPIInt rank,commsize;
109:   PetscInt p;
110: 
111:   MPI_Comm_size(PETSC_COMM_WORLD,&commsize);
112:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
113: 
114:   DMCreate(PETSC_COMM_WORLD,&dms);
115:   DMSetType(dms,DMSWARM);
116: 
117:   DMSwarmInitializeFieldRegister(dms);
118: 
119:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
120:   DMSwarmRegisterPetscDatatypeField(dms,"strain",1,PETSC_REAL);
121: 
122:   DMSwarmFinalizeFieldRegister(dms);
123: 
124:   DMSwarmSetLocalSizes(dms,5+rank,4);
125:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);
126: 
127:   {
128:     PetscReal *array;
129:     DMSwarmGetField(dms,"viscosity",NULL,NULL,(void**)&array);
130:     for (p=0; p<5+rank; p++) {
131:       array[p] = 11.1 + p*0.1 + rank*100.0;
132:     }
133:     DMSwarmRestoreField(dms,"viscosity",NULL,NULL,(void**)&array);
134:   }
135: 
136:   {
137:     PetscReal *array;
138:     DMSwarmGetField(dms,"strain",NULL,NULL,(void**)&array);
139:     for (p=0; p<5+rank; p++) {
140:       array[p] = 2.0e-2 + p*0.001 + rank*1.0;
141:     }
142:     DMSwarmRestoreField(dms,"strain",NULL,NULL,(void**)&array);
143:   }
144: 
145:   {
146:     PetscInt *rankval;
147:     PetscInt npoints[2],npoints_orig[2];
148: 
149:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
150:     DMSwarmGetSize(dms,&npoints_orig[1]);
151:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);

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

155:     if (rank == 1) {
156:       rankval[0] = -1;
157:     }
158:     if (rank == 2) {
159:       rankval[1] = -1;
160:     }
161:     if (rank == 3) {
162:       rankval[3] = -1;
163:       rankval[4] = -1;
164:     }
165:     DMSwarmRestoreField(dms,"DMSwarm_rank",NULL,NULL,(void**)&rankval);
166: 
167:     DMSwarmCollectViewCreate(dms);
168:     DMSwarmGetLocalSize(dms,&npoints[0]);
169:     DMSwarmGetSize(dms,&npoints[1]);
170:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);

172:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
173:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
174:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);

176:     DMSwarmCollectViewDestroy(dms);
177:     DMSwarmGetLocalSize(dms,&npoints[0]);
178:     DMSwarmGetSize(dms,&npoints[1]);
179:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] after_v(%D,%D)\n",rank,npoints[0],npoints[1]);

181:     DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);
182:     VecView(x,PETSC_VIEWER_STDOUT_WORLD);
183:     DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);
184:   }
185:   DMDestroy(&dms);
186: 
187:   return(0);
188: }

190: /*
191:  splot "c-rank0.gp","c-rank1.gp","c-rank2.gp","c-rank3.gp"
192: */
193: PetscErrorCode ex1_3(void)
194: {
195:   DM dms;
197:   PetscMPIInt rank,commsize;
198:   PetscInt is,js,ni,nj,overlap;
199:   DM dmcell;
200: 
201:   MPI_Comm_size(PETSC_COMM_WORLD,&commsize);
202:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
203: 
204:   overlap = 2;
205:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,13,13,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);
206:   DMSetFromOptions(dmcell);
207:   DMSetUp(dmcell);
208:   DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);
209: 
210:   DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);
211: 
212:   DMCreate(PETSC_COMM_WORLD,&dms);
213:   DMSetType(dms,DMSWARM);
214:   DMSwarmSetCellDM(dms,dmcell);
215: 
216:   /* load in data types */
217:   DMSwarmInitializeFieldRegister(dms);
218: 
219:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
220:   DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);
221:   DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);
222: 
223:   DMSwarmFinalizeFieldRegister(dms);
224: 
225:   DMSwarmSetLocalSizes(dms,ni*nj*4,4);
226:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);
227: 
228:   /* set values within the swarm */
229:   {
230:     PetscReal *array_x,*array_y;
231:     PetscInt npoints,i,j,cnt;
232:     DMDACoor2d **LA_coor;
233:     Vec coor;
234:     DM dmcellcdm;
235: 
236:     DMGetCoordinateDM(dmcell,&dmcellcdm);
237:     DMGetCoordinates(dmcell,&coor);
238:     DMDAVecGetArray(dmcellcdm,coor,&LA_coor);
239: 
240:     DMSwarmGetLocalSize(dms,&npoints);
241:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
242:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
243:     cnt = 0;
244:     for (j=js; j<js+nj; j++) {
245:       for (i=is; i<is+ni; i++) {
246:         PetscReal xp,yp;
247: 
248:         xp = PetscRealPart(LA_coor[j][i].x);
249:         yp = PetscRealPart(LA_coor[j][i].y);
250: 
251:         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; }
252:         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; }
253:         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; }
254:         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; }

256:         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; }
257:         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; }
258:         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; }
259:         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; }

261:         cnt++;
262:       }
263:     }
264: 
265:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
266:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
267:     DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);
268:   }
269: 
270:   {
271:     PetscInt npoints[2],npoints_orig[2],ng;
272: 
273:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
274:     DMSwarmGetSize(dms,&npoints_orig[1]);
275:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);
276: 
277:     DMSwarmCollect_DMDABoundingBox(dms,&ng);
278: 
279:     DMSwarmGetLocalSize(dms,&npoints[0]);
280:     DMSwarmGetSize(dms,&npoints[1]);
281:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);
282:   }
283: 
284:   {
285:     PetscReal *array_x,*array_y;
286:     PetscInt npoints,p;
287:     FILE *fp = NULL;
288:     char name[PETSC_MAX_PATH_LEN];
289: 
290:     PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);
291:     fp = fopen(name,"w");
292:     if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
293:     DMSwarmGetLocalSize(dms,&npoints);
294:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
295:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
296:     for (p=0; p<npoints; p++) {
297:       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
298:     }
299:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
300:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
301:     fclose(fp);
302:   }
303: 
304:   DMDestroy(&dmcell);
305:   DMDestroy(&dms);
306: 
307:   return(0);
308: }

310: typedef struct {
311:   PetscReal cx[2];
312:   PetscReal radius;
313: } CollectZoneCtx;

315: PetscErrorCode collect_zone(DM dm,void *ctx,PetscInt *nfound,PetscInt **foundlist)
316: {
317:   CollectZoneCtx *zone = (CollectZoneCtx*)ctx;
318:   PetscInt p,npoints;
319:   PetscReal *array_x,*array_y,r2;
320:   PetscInt p2collect,*plist;
321:   PetscMPIInt rank;
323: 
324:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

326:   /*PetscPrintf(PETSC_COMM_WORLD,"z %1.4e,%1.4e --> %1.4e\n",zone->cx[0],zone->cx[1],zone->radius);*/
327: 
328:   DMSwarmGetLocalSize(dm,&npoints);
329:   DMSwarmGetField(dm,"coorx",NULL,NULL,(void**)&array_x);
330:   DMSwarmGetField(dm,"coory",NULL,NULL,(void**)&array_y);

332:   r2 = zone->radius * zone->radius;
333:   p2collect = 0;
334:   for (p=0; p<npoints; p++) {
335:     PetscReal sep2;
336: 
337:     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
338:     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
339:     if (sep2 < r2) {
340:       p2collect++;
341:     }
342:   }
343: 
344:   PetscMalloc1(p2collect+1,&plist);
345:   p2collect = 0;
346:   for (p=0; p<npoints; p++) {
347:     PetscReal sep2;
348: 
349:     sep2  = (array_x[p] - zone->cx[0])*(array_x[p] - zone->cx[0]);
350:     sep2 += (array_y[p] - zone->cx[1])*(array_y[p] - zone->cx[1]);
351:     if (sep2 < r2) {
352:       plist[p2collect] = p;
353:       p2collect++;
354:     }
355:   }
356:   DMSwarmRestoreField(dm,"coory",NULL,NULL,(void**)&array_y);
357:   DMSwarmRestoreField(dm,"coorx",NULL,NULL,(void**)&array_x);
358: 
359:   /*PetscPrintf(PETSC_COMM_WORLD,"rank[%d]: p2collect = %d\n",rank,p2collect);*/

361:   *nfound = p2collect;
362:   *foundlist = plist;
363:   return(0);
364: }

366: PetscErrorCode ex1_4(void)
367: {
368:   DM dms;
370:   PetscMPIInt rank,commsize;
371:   PetscInt is,js,ni,nj,overlap,nn;
372:   DM dmcell;
373:   CollectZoneCtx *zone;
374:   PetscReal dx;
375: 
376:   MPI_Comm_size(PETSC_COMM_WORLD,&commsize);
377:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
378: 
379:   nn = 101;
380:   dx = 2.0/ (PetscReal)(nn-1);
381:   overlap = 0;
382:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nn,nn,PETSC_DECIDE,PETSC_DECIDE,1,overlap,NULL,NULL,&dmcell);
383:   DMSetFromOptions(dmcell);
384:   DMSetUp(dmcell);
385:   DMDASetUniformCoordinates(dmcell,-1.0,1.0,-1.0,1.0,-1.0,1.0);
386: 
387:   DMDAGetCorners(dmcell,&is,&js,NULL,&ni,&nj,NULL);
388: 
389:   DMCreate(PETSC_COMM_WORLD,&dms);
390:   DMSetType(dms,DMSWARM);
391: 
392:   /* load in data types */
393:   DMSwarmInitializeFieldRegister(dms);
394: 
395:   DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);
396:   DMSwarmRegisterPetscDatatypeField(dms,"coorx",1,PETSC_REAL);
397:   DMSwarmRegisterPetscDatatypeField(dms,"coory",1,PETSC_REAL);
398: 
399:   DMSwarmFinalizeFieldRegister(dms);
400: 
401:   DMSwarmSetLocalSizes(dms,ni*nj*4,4);
402:   DMView(dms,PETSC_VIEWER_STDOUT_WORLD);
403: 
404:   /* set values within the swarm */
405:   {
406:     PetscReal *array_x,*array_y;
407:     PetscInt npoints,i,j,cnt;
408:     DMDACoor2d **LA_coor;
409:     Vec coor;
410:     DM dmcellcdm;
411: 
412:     DMGetCoordinateDM(dmcell,&dmcellcdm);
413:     DMGetCoordinates(dmcell,&coor);
414:     DMDAVecGetArray(dmcellcdm,coor,&LA_coor);
415: 
416:     DMSwarmGetLocalSize(dms,&npoints);
417:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
418:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
419:     cnt = 0;
420:     for (j=js; j<js+nj; j++) {
421:       for (i=is; i<is+ni; i++) {
422:         PetscReal xp,yp;
423: 
424:         xp = PetscRealPart(LA_coor[j][i].x);
425:         yp = PetscRealPart(LA_coor[j][i].y);
426: 
427:         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; }*/
428:         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; }*/
429:         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; }*/
430:         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; }*/
431: 
432:         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; }*/
433:         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; }*/
434:         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; }*/
435:         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; }*/
436: 
437:         cnt++;
438:       }
439:     }
440: 
441:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
442:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
443:     DMDAVecRestoreArray(dmcellcdm,coor,&LA_coor);
444:   }
445: 
446:   PetscMalloc1(1,&zone);
447:   if (commsize == 4) {
448:     if (rank == 0) {
449:       zone->cx[0] = 0.5;
450:       zone->cx[1] = 0.5;
451:       zone->radius = 0.3;
452:     }
453:     if (rank == 1) {
454:       zone->cx[0] = -0.5;
455:       zone->cx[1] = 0.5;
456:       zone->radius = 0.25;
457:     }
458:     if (rank == 2) {
459:       zone->cx[0] = 0.5;
460:       zone->cx[1] = -0.5;
461:       zone->radius = 0.2;
462:     }
463:     if (rank == 3) {
464:       zone->cx[0] = -0.5;
465:       zone->cx[1] = -0.5;
466:       zone->radius = 0.1;
467:     }
468:   } else {
469:     if (rank == 0) {
470:       zone->cx[0] = 0.5;
471:       zone->cx[1] = 0.5;
472:       zone->radius = 0.8;
473:     } else {
474:       zone->cx[0] = 10.0;
475:       zone->cx[1] = 10.0;
476:       zone->radius = 0.0;
477:     }
478:   }
479: 
480:   {
481:     PetscInt npoints[2],npoints_orig[2],ng;
482: 
483:     DMSwarmGetLocalSize(dms,&npoints_orig[0]);
484:     DMSwarmGetSize(dms,&npoints_orig[1]);
485:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] before(%D,%D)\n",rank,npoints_orig[0],npoints_orig[1]);
486: 
487:     DMSwarmCollect_General(dms,collect_zone,sizeof(CollectZoneCtx),zone,&ng);
488: 
489:     DMSwarmGetLocalSize(dms,&npoints[0]);
490:     DMSwarmGetSize(dms,&npoints[1]);
491:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] after(%D,%D)\n",rank,npoints[0],npoints[1]);
492:   }
493: 
494:   {
495:     PetscReal *array_x,*array_y;
496:     PetscInt npoints,p;
497:     FILE *fp = NULL;
498:     char name[PETSC_MAX_PATH_LEN];
499: 
500:     PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"c-rank%d.gp",rank);
501:     fp = fopen(name,"w");
502:     if (!fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"Cannot open file %s",name);
503:     DMSwarmGetLocalSize(dms,&npoints);
504:     DMSwarmGetField(dms,"coorx",NULL,NULL,(void**)&array_x);
505:     DMSwarmGetField(dms,"coory",NULL,NULL,(void**)&array_y);
506:     for (p=0; p<npoints; p++) {
507:       fprintf(fp,"%+1.4e %+1.4e %1.4e\n",array_x[p],array_y[p],(double)rank);
508:     }
509:     DMSwarmRestoreField(dms,"coory",NULL,NULL,(void**)&array_y);
510:     DMSwarmRestoreField(dms,"coorx",NULL,NULL,(void**)&array_x);
511:     fclose(fp);
512:   }
513: 
514:   DMDestroy(&dmcell);
515:   DMDestroy(&dms);
516:   PetscFree(zone);
517: 
518:   return(0);
519: }

521: int main(int argc,char **argv)
522: {
524:   PetscInt test_mode = 4;
525: 
526:   PetscInitialize(&argc,&argv,(char*)0,help);
527:   PetscOptionsGetInt(NULL,NULL,"-test_mode",&test_mode,NULL);
528:   if (test_mode == 1) {
529:     ex1_1();
530:   } else if (test_mode == 2) {
531:     ex1_2();
532:   } else if (test_mode == 3) {
533:     ex1_3();
534:   } else if (test_mode == 4) {
535:     ex1_4();
536:   } else {
537:     PetscPrintf(PETSC_COMM_WORLD,"Unknown test_mode value, should be 1,2,3,4\n");
538:   }
539:   PetscFinalize();
540:   return 0;
541: }