Actual source code: gr1.c

  1: /*$Id: gr1.c,v 1.30 2001/08/10 03:34:34 bsmith Exp $*/

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

 7:  #include petscda.h

  9: /*@
 10:     DASetUniformCoordinates - Sets a DA coordinates to be a uniform grid

 12:   Collective on DA

 14:   Input Parameters:
 15: +  da - the distributed array object
 16: .  xmin,xmax - extremes in the x direction
 17: .  xmin,xmax - extremes in the y direction
 18: -  xmin,xmax - extremes in the z direction

 20:   Level: beginner

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

 24: @*/
 25: int DASetUniformCoordinates(DA da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
 26: {
 27:   int            i,j,k,ierr,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;
 28:   PetscReal      hx,hy,hz_;
 29:   Vec            xcoor;
 30:   DAPeriodicType periodic;
 31:   PetscScalar    *coors;

 34:   if (xmax <= xmin) SETERRQ2(1,"Xmax must be larger than xmin %g %g",xmin,xmax);

 36:   DAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&periodic,0);
 37:   DAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);

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

 94:   return(0);
 95: }

 97: int VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
 98: {
 99:   DA             da;
100:   int            i,rank,size,ierr,n,tag1,tag2,N,step;
101:   int            istart,isize,j;
102:   MPI_Status     status;
103:   PetscReal      coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
104:   PetscScalar    *array,*xg;
105:   PetscDraw      draw;
106:   PetscTruth     isnull,showpoints;
107:   MPI_Comm       comm;
108:   PetscDrawAxis  axis;
109:   Vec            xcoor;
110:   DAPeriodicType periodic;

113:   PetscViewerDrawGetDraw(v,0,&draw);
114:   PetscDrawIsNull(draw,&isnull); if (isnull) return(0);

116:   PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
117:   if (!da) SETERRQ(1,"Vector not generated from a DA");

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

121:   DAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&periodic,0);
122:   DAGetCorners(da,&istart,0,0,&isize,0,0);
123:   VecGetArray(xin,&array);
124:   VecGetLocalSize(xin,&n);
125:   n    = n/step;

127:   /* get coordinates of nodes */
128:   DAGetCoordinates(da,&xcoor);
129:   if (!xcoor) {
130:     DASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);
131:     DAGetCoordinates(da,&xcoor);
132:   }
133:   VecGetArray(xcoor,&xg);

135:   PetscObjectGetComm((PetscObject)xin,&comm);
136:   MPI_Comm_size(comm,&size);
137:   MPI_Comm_rank(comm,&rank);

139:   /*
140:       Determine the min and max x coordinate in plot 
141:   */
142:   if (!rank) {
143:     xmin = PetscRealPart(xg[0]);
144:   }
145:   if (rank == size-1) {
146:     xmax = PetscRealPart(xg[n-1]);
147:   }
148:   MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);
149:   MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);

151:   for (j=0; j<step; j++) {
152:     PetscViewerDrawGetDraw(v,j,&draw);
153:     PetscDrawCheckResizedWindow(draw);

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

175:     PetscDrawSynchronizedClear(draw);
176:     PetscViewerDrawGetDrawAxis(v,j,&axis);
177:     PetscLogObjectParent(draw,axis);
178:     if (!rank) {
179:       char *title;

181:       PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);
182:       PetscDrawAxisDraw(axis);
183:       PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);
184:       DAGetFieldName(da,j,&title);
185:       if (title) {PetscDrawSetTitle(draw,title);}
186:     }
187:     MPI_Bcast(coors,4,MPIU_REAL,0,comm);
188:     if (rank) {
189:       PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
190:     }

192:     /* draw local part of vector */
193:     PetscObjectGetNewTag((PetscObject)xin,&tag1);
194:     PetscObjectGetNewTag((PetscObject)xin,&tag2);
195:     if (rank < size-1) { /*send value to right */
196:       MPI_Send(&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);
197:       MPI_Send(&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);
198:     }
199:     if (!rank && periodic) { /* first processor sends first value to last */
200:       MPI_Send(&array[j],1,MPIU_REAL,size-1,tag2,comm);
201:     }

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