Actual source code: gr1.c

  1: /*$Id: gr1.c,v 1.27 2001/04/10 19:37:23 bsmith Exp $*/

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

  7: #include "petscda.h"      /*I  "petscda.h"   I*/

  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,double xmin,double xmax,double ymin,double ymax,double zmin,double zmax)
 26: {
 27:   int            i,j,k,ierr,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;
 28:   double         hx,hy,hz_;
 29:   Vec            xcoor;
 30:   DAPeriodicType periodic;
 31:   Scalar         *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:   double         coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
104:   Scalar         *array,*xg;
105:   PetscDraw      draw;
106:   PetscTruth     isnull;
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:   DAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&periodic,0);
120:   DAGetCorners(da,&istart,0,0,&isize,0,0);
121:   VecGetArray(xin,&array);
122:   VecGetLocalSize(xin,&n);
123:   n    = n/step;

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

133:   PetscObjectGetComm((PetscObject)xin,&comm);
134:   MPI_Comm_size(comm,&size);
135:   MPI_Comm_rank(comm,&rank);

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

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

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

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

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

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

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