Actual source code: gr1.c

  1: #define PETSCDM_DLL

  3: /* 
  4:    Plots vectors obtained with DACreate1d()
  5: */

 7:  #include petscda.h

 11: /*@
 12:     DASetUniformCoordinates - Sets a DA coordinates to be a uniform grid

 14:   Collective on DA

 16:   Input Parameters:
 17: +  da - the distributed array object
 18: .  xmin,xmax - extremes in the x direction
 19: .  ymin,ymax - extremes in the y direction (use PETSC_NULL for 1 dimensional problems)
 20: -  zmin,zmax - extremes in the z direction (use PETSC_NULL for 1 or 2 dimensional problems)

 22:   Level: beginner

 24: .seealso: DASetCoordinates(), DAGetCoordinates(), DACreate1d(), DACreate2d(), DACreate3d()

 26: @*/
 27: PetscErrorCode  DASetUniformCoordinates(DA da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
 28: {
 30:   PetscInt       i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;
 31:   PetscReal      hx,hy,hz_;
 32:   Vec            xcoor;
 33:   DAPeriodicType periodic;
 34:   PetscScalar    *coors;
 35:   MPI_Comm       comm;

 38:   if (xmax <= xmin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %G %G",xmin,xmax);

 40:   PetscObjectGetComm((PetscObject)da,&comm);
 41:   DAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&periodic,0);
 42:   DAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);

 44:   if (dim == 1) {
 45:     VecCreateMPI(comm,isize,PETSC_DETERMINE,&xcoor);
 46:     if (periodic == DA_NONPERIODIC) hx = (xmax-xmin)/(M-1);
 47:     else                            hx = (xmax-xmin)/M;
 48:     VecGetArray(xcoor,&coors);
 49:     for (i=0; i<isize; i++) {
 50:       coors[i] = xmin + hx*(i+istart);
 51:     }
 52:     VecRestoreArray(xcoor,&coors);
 53:   } else if (dim == 2) {
 54:     if (ymax <= ymin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax);
 55:     VecCreateMPI(comm,2*isize*jsize,PETSC_DETERMINE,&xcoor);
 56:     VecSetBlockSize(xcoor,2);
 57:     if (DAXPeriodic(periodic)) hx = (xmax-xmin)/(M);
 58:     else                       hx = (xmax-xmin)/(M-1);
 59:     if (DAYPeriodic(periodic)) hy = (ymax-ymin)/(N);
 60:     else                       hy = (ymax-ymin)/(N-1);
 61:     VecGetArray(xcoor,&coors);
 62:     cnt  = 0;
 63:     for (j=0; j<jsize; j++) {
 64:       for (i=0; i<isize; i++) {
 65:         coors[cnt++] = xmin + hx*(i+istart);
 66:         coors[cnt++] = ymin + hy*(j+jstart);
 67:       }
 68:     }
 69:     VecRestoreArray(xcoor,&coors);
 70:   } else if (dim == 3) {
 71:     if (ymax <= ymin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax);
 72:     if (zmax <= zmin) SETERRQ2(PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %G %G",zmin,zmax);
 73:     VecCreateMPI(comm,3*isize*jsize*ksize,PETSC_DETERMINE,&xcoor);
 74:     VecSetBlockSize(xcoor,3);
 75:     if (DAXPeriodic(periodic)) hx = (xmax-xmin)/(M);
 76:     else                       hx = (xmax-xmin)/(M-1);
 77:     if (DAYPeriodic(periodic)) hy = (ymax-ymin)/(N);
 78:     else                       hy = (ymax-ymin)/(N-1);
 79:     if (DAZPeriodic(periodic)) hz_ = (zmax-zmin)/(P);
 80:     else                       hz_ = (zmax-zmin)/(P-1);
 81:     VecGetArray(xcoor,&coors);
 82:     cnt  = 0;
 83:     for (k=0; k<ksize; k++) {
 84:       for (j=0; j<jsize; j++) {
 85:         for (i=0; i<isize; i++) {
 86:           coors[cnt++] = xmin + hx*(i+istart);
 87:           coors[cnt++] = ymin + hy*(j+jstart);
 88:           coors[cnt++] = zmin + hz_*(k+kstart);
 89:         }
 90:       }
 91:     }
 92:     VecRestoreArray(xcoor,&coors);
 93:   } else {
 94:     SETERRQ1(PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
 95:   }
 96:   DASetCoordinates(da,xcoor);
 97:   PetscLogObjectParent(da,xcoor);

 99:   return(0);
100: }

104: PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
105: {
106:   DA             da;
108:   PetscMPIInt    rank,size,tag1,tag2;
109:   PetscInt       i,n,N,step,istart,isize,j;
110:   MPI_Status     status;
111:   PetscReal      coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
112:   PetscScalar    *array,*xg;
113:   PetscDraw      draw;
114:   PetscTruth     isnull,showpoints;
115:   MPI_Comm       comm;
116:   PetscDrawAxis  axis;
117:   Vec            xcoor;
118:   DAPeriodicType periodic;

121:   PetscViewerDrawGetDraw(v,0,&draw);
122:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

124:   PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
125:   if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");

127:   PetscOptionsHasName(PETSC_NULL,"-draw_vec_mark_points",&showpoints);

129:   DAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&periodic,0);
130:   DAGetCorners(da,&istart,0,0,&isize,0,0);
131:   VecGetArray(xin,&array);
132:   VecGetLocalSize(xin,&n);
133:   n    = n/step;

135:   /* get coordinates of nodes */
136:   DAGetCoordinates(da,&xcoor);
137:   if (!xcoor) {
138:     DASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);
139:     DAGetCoordinates(da,&xcoor);
140:   }
141:   VecGetArray(xcoor,&xg);

143:   PetscObjectGetComm((PetscObject)xin,&comm);
144:   MPI_Comm_size(comm,&size);
145:   MPI_Comm_rank(comm,&rank);

147:   /*
148:       Determine the min and max x coordinate in plot 
149:   */
150:   if (!rank) {
151:     xmin = PetscRealPart(xg[0]);
152:   }
153:   if (rank == size-1) {
154:     xmax = PetscRealPart(xg[n-1]);
155:   }
156:   MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);
157:   MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);

159:   for (j=0; j<step; j++) {
160:     PetscViewerDrawGetDraw(v,j,&draw);
161:     PetscDrawCheckResizedWindow(draw);

163:     /*
164:         Determine the min and max y coordinate in plot 
165:     */
166:     min = 1.e20; max = -1.e20;
167:     for (i=0; i<n; i++) {
168: #if defined(PETSC_USE_COMPLEX)
169:       if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]);
170:       if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]);
171: #else
172:       if (array[j+i*step] < min) min = array[j+i*step];
173:       if (array[j+i*step] > max) max = array[j+i*step];
174: #endif
175:     }
176:     if (min + 1.e-10 > max) {
177:       min -= 1.e-5;
178:       max += 1.e-5;
179:     }
180:     MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPI_MIN,0,comm);
181:     MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPI_MAX,0,comm);

183:     PetscDrawSynchronizedClear(draw);
184:     PetscViewerDrawGetDrawAxis(v,j,&axis);
185:     PetscLogObjectParent(draw,axis);
186:     if (!rank) {
187:       char *title;

189:       PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);
190:       PetscDrawAxisDraw(axis);
191:       PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);
192:       DAGetFieldName(da,j,&title);
193:       if (title) {PetscDrawSetTitle(draw,title);}
194:     }
195:     MPI_Bcast(coors,4,MPIU_REAL,0,comm);
196:     if (rank) {
197:       PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
198:     }

200:     /* draw local part of vector */
201:     PetscObjectGetNewTag((PetscObject)xin,&tag1);
202:     PetscObjectGetNewTag((PetscObject)xin,&tag2);
203:     if (rank < size-1) { /*send value to right */
204:       MPI_Send(&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);
205:       MPI_Send(&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);
206:     }
207:     if (!rank && periodic && size > 1) { /* first processor sends first value to last */
208:       MPI_Send(&array[j],1,MPIU_REAL,size-1,tag2,comm);
209:     }

211:     for (i=1; i<n; i++) {
212: #if !defined(PETSC_USE_COMPLEX)
213:       PetscDrawLine(draw,xg[i-1],array[j+step*(i-1)],xg[i],array[j+step*i],
214:                       PETSC_DRAW_RED);
215: #else
216:       PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),
217:                       PetscRealPart(xg[i]),PetscRealPart(array[j+step*i]),PETSC_DRAW_RED);
218: #endif
219:       if (showpoints) {
220:         PetscDrawPoint(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);
221:       }
222:     }
223:     if (rank) { /* receive value from left */
224:       MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);
225:       MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);
226: #if !defined(PETSC_USE_COMPLEX)
227:       PetscDrawLine(draw,xgtmp,tmp,xg[0],array[j],PETSC_DRAW_RED);
228: #else
229:       PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),
230:                       PETSC_DRAW_RED);
231: #endif
232:       if (showpoints) {
233:         PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);
234:       }
235:     }
236:     if (rank == size-1 && periodic && size > 1) {
237:       MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);
238: #if !defined(PETSC_USE_COMPLEX)
239:       PetscDrawLine(draw,xg[n-2],array[j+step*(n-1)],xg[n-1],tmp,PETSC_DRAW_RED);
240: #else
241:       PetscDrawLine(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),
242:                       PetscRealPart(xg[n-1]),tmp,PETSC_DRAW_RED);
243: #endif
244:       if (showpoints) {
245:         PetscDrawPoint(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);
246:       }
247:     }
248:     PetscDrawSynchronizedFlush(draw);
249:     PetscDrawPause(draw);
250:   }
251:   VecRestoreArray(xcoor,&xg);
252:   VecRestoreArray(xin,&array);
253:   return(0);
254: }