xref: /petsc/src/dm/impls/da/gr1.c (revision abf1d30e7f103cb71b51198ec84d74c0c72b5d6b)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith    Plots vectors obtained with DMDACreate1d()
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>      /*I  "petscdmda.h"   I*/
747c6ae99SBarry Smith 
847c6ae99SBarry Smith /*@
9aa219208SBarry Smith     DMDASetUniformCoordinates - Sets a DMDA coordinates to be a uniform grid
1047c6ae99SBarry Smith 
11aa219208SBarry Smith   Collective on DMDA
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith   Input Parameters:
1447c6ae99SBarry Smith +  da - the distributed array object
1547c6ae99SBarry Smith .  xmin,xmax - extremes in the x direction
16362c83fbSPatrick Sanan .  ymin,ymax - extremes in the y direction (value ignored for 1 dimensional problems)
17362c83fbSPatrick Sanan -  zmin,zmax - extremes in the z direction (value ignored for 1 or 2 dimensional problems)
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith   Level: beginner
2047c6ae99SBarry Smith 
212150357eSBarry Smith .seealso: DMSetCoordinates(), DMGetCoordinates(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith @*/
247087cfbeSBarry Smith PetscErrorCode  DMDASetUniformCoordinates(DM da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
2547c6ae99SBarry Smith {
2647c6ae99SBarry Smith   MPI_Comm         comm;
2757459e9aSMatthew G Knepley   PetscSection     section;
289a42bb27SBarry Smith   DM               cda;
29bff4a2f0SMatthew G. Knepley   DMBoundaryType   bx,by,bz;
3047c6ae99SBarry Smith   Vec              xcoor;
3147c6ae99SBarry Smith   PetscScalar      *coors;
3247c6ae99SBarry Smith   PetscReal        hx,hy,hz_;
3347c6ae99SBarry Smith   PetscInt         i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;
3447c6ae99SBarry Smith   PetscErrorCode   ierr;
3547c6ae99SBarry Smith 
3647c6ae99SBarry Smith   PetscFunctionBegin;
375b990481SMatthew G Knepley   PetscValidHeaderSpecific(da,DM_CLASSID,1);
38821c30cfSPeter Brune   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
39c593f006SBarry Smith   if (xmax < xmin) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %g %g",(double)xmin,(double)xmax);
40*abf1d30eSPatrick Sanan   if ((dim > 1) && (ymax < ymin)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %g %g",(double)ymin,(double)ymax);
41*abf1d30eSPatrick Sanan   if ((dim > 2) && (zmax < zmin)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %g %g",(double)zmin,(double)zmax);
4247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
4357459e9aSMatthew G Knepley   ierr = DMGetDefaultSection(da,&section);CHKERRQ(ierr);
44aa219208SBarry Smith   ierr = DMDAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);CHKERRQ(ierr);
456636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(da, &cda);CHKERRQ(ierr);
4657459e9aSMatthew G Knepley   if (section) {
4757459e9aSMatthew G Knepley     /* This would be better as a vector, but this is compatible */
4857459e9aSMatthew G Knepley     PetscInt numComp[3]      = {1, 1, 1};
4957459e9aSMatthew G Knepley     PetscInt numVertexDof[3] = {1, 1, 1};
5057459e9aSMatthew G Knepley 
5157459e9aSMatthew G Knepley     ierr = DMDASetFieldName(cda, 0, "x");CHKERRQ(ierr);
5257459e9aSMatthew G Knepley     if (dim > 1) {ierr = DMDASetFieldName(cda, 1, "y");CHKERRQ(ierr);}
5357459e9aSMatthew G Knepley     if (dim > 2) {ierr = DMDASetFieldName(cda, 2, "z");CHKERRQ(ierr);}
540298fd71SBarry Smith     ierr = DMDACreateSection(cda, numComp, numVertexDof, NULL, NULL);CHKERRQ(ierr);
5557459e9aSMatthew G Knepley   }
56564755cdSBarry Smith   ierr = DMCreateGlobalVector(cda, &xcoor);CHKERRQ(ierr);
5757459e9aSMatthew G Knepley   if (section) {
5857459e9aSMatthew G Knepley     PetscSection csection;
5957459e9aSMatthew G Knepley     PetscInt     vStart, vEnd;
6057459e9aSMatthew G Knepley 
6157459e9aSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(cda,&csection);CHKERRQ(ierr);
6257459e9aSMatthew G Knepley     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
6357459e9aSMatthew G Knepley     ierr = DMDAGetHeightStratum(da, dim, &vStart, &vEnd);CHKERRQ(ierr);
64bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx  = (xmax-xmin)/(M+1);
6557459e9aSMatthew G Knepley     else                              hx  = (xmax-xmin)/(M ? M : 1);
66bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy  = (ymax-ymin)/(N+1);
6757459e9aSMatthew G Knepley     else                              hy  = (ymax-ymin)/(N ? N : 1);
68bff4a2f0SMatthew G. Knepley     if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P+1);
6957459e9aSMatthew G Knepley     else                              hz_ = (zmax-zmin)/(P ? P : 1);
7057459e9aSMatthew G Knepley     switch (dim) {
7157459e9aSMatthew G Knepley     case 1:
7257459e9aSMatthew G Knepley       for (i = 0; i < isize+1; ++i) {
7357459e9aSMatthew G Knepley         PetscInt v = i+vStart, dof, off;
7457459e9aSMatthew G Knepley 
7557459e9aSMatthew G Knepley         ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
7657459e9aSMatthew G Knepley         ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
7757459e9aSMatthew G Knepley         if (off >= 0) {
7857459e9aSMatthew G Knepley           coors[off] = xmin + hx*(i+istart);
7957459e9aSMatthew G Knepley         }
8057459e9aSMatthew G Knepley       }
8157459e9aSMatthew G Knepley       break;
8257459e9aSMatthew G Knepley     case 2:
8357459e9aSMatthew G Knepley       for (j = 0; j < jsize+1; ++j) {
8457459e9aSMatthew G Knepley         for (i = 0; i < isize+1; ++i) {
8557459e9aSMatthew G Knepley           PetscInt v = j*(isize+1)+i+vStart, dof, off;
8657459e9aSMatthew G Knepley 
8757459e9aSMatthew G Knepley           ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
8857459e9aSMatthew G Knepley           ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
8957459e9aSMatthew G Knepley           if (off >= 0) {
9057459e9aSMatthew G Knepley             coors[off+0] = xmin + hx*(i+istart);
9157459e9aSMatthew G Knepley             coors[off+1] = ymin + hy*(j+jstart);
9257459e9aSMatthew G Knepley           }
9357459e9aSMatthew G Knepley         }
9457459e9aSMatthew G Knepley       }
9557459e9aSMatthew G Knepley       break;
9657459e9aSMatthew G Knepley     case 3:
9757459e9aSMatthew G Knepley       for (k = 0; k < ksize+1; ++k) {
9857459e9aSMatthew G Knepley         for (j = 0; j < jsize+1; ++j) {
9957459e9aSMatthew G Knepley           for (i = 0; i < isize+1; ++i) {
10057459e9aSMatthew G Knepley             PetscInt v = (k*(jsize+1)+j)*(isize+1)+i+vStart, dof, off;
10157459e9aSMatthew G Knepley 
10257459e9aSMatthew G Knepley             ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
10357459e9aSMatthew G Knepley             ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
10457459e9aSMatthew G Knepley             if (off >= 0) {
10557459e9aSMatthew G Knepley               coors[off+0] = xmin + hx*(i+istart);
10657459e9aSMatthew G Knepley               coors[off+1] = ymin + hy*(j+jstart);
10757459e9aSMatthew G Knepley               coors[off+2] = zmin + hz_*(k+kstart);
10857459e9aSMatthew G Knepley             }
10957459e9aSMatthew G Knepley           }
11057459e9aSMatthew G Knepley         }
11157459e9aSMatthew G Knepley       }
11257459e9aSMatthew G Knepley       break;
11357459e9aSMatthew G Knepley     default:
114ce94432eSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
11557459e9aSMatthew G Knepley     }
11657459e9aSMatthew G Knepley     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
1176636e97aSMatthew G Knepley     ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr);
1183bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)xcoor);CHKERRQ(ierr);
11957459e9aSMatthew G Knepley     ierr = VecDestroy(&xcoor);CHKERRQ(ierr);
12057459e9aSMatthew G Knepley     PetscFunctionReturn(0);
12157459e9aSMatthew G Knepley   }
12247c6ae99SBarry Smith   if (dim == 1) {
123bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/M;
124ce00eea3SSatish Balay     else hx = (xmax-xmin)/(M-1);
12547c6ae99SBarry Smith     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
12647c6ae99SBarry Smith     for (i=0; i<isize; i++) {
12747c6ae99SBarry Smith       coors[i] = xmin + hx*(i+istart);
12847c6ae99SBarry Smith     }
12947c6ae99SBarry Smith     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
13047c6ae99SBarry Smith   } else if (dim == 2) {
131bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
13247c6ae99SBarry Smith     else hx = (xmax-xmin)/(M-1);
133bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
13447c6ae99SBarry Smith     else hy = (ymax-ymin)/(N-1);
13547c6ae99SBarry Smith     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
13647c6ae99SBarry Smith     cnt  = 0;
13747c6ae99SBarry Smith     for (j=0; j<jsize; j++) {
13847c6ae99SBarry Smith       for (i=0; i<isize; i++) {
13947c6ae99SBarry Smith         coors[cnt++] = xmin + hx*(i+istart);
14047c6ae99SBarry Smith         coors[cnt++] = ymin + hy*(j+jstart);
14147c6ae99SBarry Smith       }
14247c6ae99SBarry Smith     }
14347c6ae99SBarry Smith     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
14447c6ae99SBarry Smith   } else if (dim == 3) {
145bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
14647c6ae99SBarry Smith     else hx = (xmax-xmin)/(M-1);
147bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
14847c6ae99SBarry Smith     else hy = (ymax-ymin)/(N-1);
149bff4a2f0SMatthew G. Knepley     if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P);
15047c6ae99SBarry Smith     else hz_ = (zmax-zmin)/(P-1);
15147c6ae99SBarry Smith     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
15247c6ae99SBarry Smith     cnt  = 0;
15347c6ae99SBarry Smith     for (k=0; k<ksize; k++) {
15447c6ae99SBarry Smith       for (j=0; j<jsize; j++) {
15547c6ae99SBarry Smith         for (i=0; i<isize; i++) {
15647c6ae99SBarry Smith           coors[cnt++] = xmin + hx*(i+istart);
15747c6ae99SBarry Smith           coors[cnt++] = ymin + hy*(j+jstart);
15847c6ae99SBarry Smith           coors[cnt++] = zmin + hz_*(k+kstart);
15947c6ae99SBarry Smith         }
16047c6ae99SBarry Smith       }
16147c6ae99SBarry Smith     }
16247c6ae99SBarry Smith     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
163ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
1646636e97aSMatthew G Knepley   ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr);
1653bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)xcoor);CHKERRQ(ierr);
166fcfd50ebSBarry Smith   ierr = VecDestroy(&xcoor);CHKERRQ(ierr);
16747c6ae99SBarry Smith   PetscFunctionReturn(0);
16847c6ae99SBarry Smith }
16947c6ae99SBarry Smith 
1709ae14b6eSBarry Smith /*
1719ae14b6eSBarry Smith     Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA
1729ae14b6eSBarry Smith */
1734e6118eeSBarry Smith PetscErrorCode DMDASelectFields(DM da,PetscInt *outfields,PetscInt **fields)
1744e6118eeSBarry Smith {
1754e6118eeSBarry Smith   PetscErrorCode ierr;
1764e6118eeSBarry Smith   PetscInt       step,ndisplayfields,*displayfields,k,j;
1774e6118eeSBarry Smith   PetscBool      flg;
1784e6118eeSBarry Smith 
1794e6118eeSBarry Smith   PetscFunctionBegin;
1804e6118eeSBarry Smith   ierr = DMDAGetInfo(da,0,0,0,0,0,0,0,&step,0,0,0,0,0);CHKERRQ(ierr);
181785e854fSJed Brown   ierr = PetscMalloc1(step,&displayfields);CHKERRQ(ierr);
1824e6118eeSBarry Smith   for (k=0; k<step; k++) displayfields[k] = k;
1834e6118eeSBarry Smith   ndisplayfields = step;
184c5929fdfSBarry Smith   ierr           = PetscOptionsGetIntArray(NULL,NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr);
1854e6118eeSBarry Smith   if (!ndisplayfields) ndisplayfields = step;
1864e6118eeSBarry Smith   if (!flg) {
1874e6118eeSBarry Smith     char       **fields;
1884e6118eeSBarry Smith     const char *fieldname;
1894e6118eeSBarry Smith     PetscInt   nfields = step;
190785e854fSJed Brown     ierr = PetscMalloc1(step,&fields);CHKERRQ(ierr);
191c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-draw_fields_by_name",fields,&nfields,&flg);CHKERRQ(ierr);
1924e6118eeSBarry Smith     if (flg) {
1934e6118eeSBarry Smith       ndisplayfields = 0;
1944e6118eeSBarry Smith       for (k=0; k<nfields;k++) {
1954e6118eeSBarry Smith         for (j=0; j<step; j++) {
1964e6118eeSBarry Smith           ierr = DMDAGetFieldName(da,j,&fieldname);CHKERRQ(ierr);
1974e6118eeSBarry Smith           ierr = PetscStrcmp(fieldname,fields[k],&flg);CHKERRQ(ierr);
1984e6118eeSBarry Smith           if (flg) {
1994e6118eeSBarry Smith             goto found;
2004e6118eeSBarry Smith           }
2014e6118eeSBarry Smith         }
202ce94432eSBarry Smith         SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Unknown fieldname %s",fields[k]);
2034e6118eeSBarry Smith found:  displayfields[ndisplayfields++] = j;
2044e6118eeSBarry Smith       }
2054e6118eeSBarry Smith     }
206e3f3e4b6SBarry Smith     for (k=0; k<nfields; k++) {
207e3f3e4b6SBarry Smith       ierr = PetscFree(fields[k]);CHKERRQ(ierr);
208e3f3e4b6SBarry Smith     }
2094e6118eeSBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
2104e6118eeSBarry Smith   }
2114e6118eeSBarry Smith   *fields    = displayfields;
2124e6118eeSBarry Smith   *outfields = ndisplayfields;
2134e6118eeSBarry Smith   PetscFunctionReturn(0);
2144e6118eeSBarry Smith }
2154e6118eeSBarry Smith 
2169804daf3SBarry Smith #include <petscdraw.h>
21773f7a4c5SBarry Smith 
21847c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
21947c6ae99SBarry Smith {
2209a42bb27SBarry Smith   DM                da;
22147c6ae99SBarry Smith   PetscErrorCode    ierr;
222e5ab1681SLisandro Dalcin   PetscMPIInt       rank,size,tag;
223e5ab1681SLisandro Dalcin   PetscInt          i,n,N,dof,istart,isize,j,nbounds;
22447c6ae99SBarry Smith   MPI_Status        status;
225e5ab1681SLisandro Dalcin   PetscReal         min,max,xmin = 0.0,xmax = 0.0,tmp = 0.0,xgtmp = 0.0;
22647c6ae99SBarry Smith   const PetscScalar *array,*xg;
22747c6ae99SBarry Smith   PetscDraw         draw;
228e5ab1681SLisandro Dalcin   PetscBool         isnull,useports = PETSC_FALSE,showmarkers = PETSC_FALSE;
22947c6ae99SBarry Smith   MPI_Comm          comm;
23047c6ae99SBarry Smith   PetscDrawAxis     axis;
23147c6ae99SBarry Smith   Vec               xcoor;
232bff4a2f0SMatthew G. Knepley   DMBoundaryType    bx;
233e5ab1681SLisandro Dalcin   const char        *tlabel = NULL,*xlabel = NULL;
234a56202b9SBarry Smith   const PetscReal   *bounds;
23503193ff8SBarry Smith   PetscInt          *displayfields;
23603193ff8SBarry Smith   PetscInt          k,ndisplayfields;
2374e6118eeSBarry Smith   PetscBool         hold;
238e5ab1681SLisandro Dalcin   PetscDrawViewPorts *ports = NULL;
239e5ab1681SLisandro Dalcin   PetscViewerFormat  format;
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith   PetscFunctionBegin;
24247c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
24345f3bb6eSLisandro Dalcin   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr);
24445f3bb6eSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
245a56202b9SBarry Smith   ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr);
24647c6ae99SBarry Smith 
247c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
248ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
249e5ab1681SLisandro Dalcin   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
250e5ab1681SLisandro Dalcin   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
251e5ab1681SLisandro Dalcin   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
25247c6ae99SBarry Smith 
253c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_vec_use_markers",&showmarkers,NULL);CHKERRQ(ierr);
25447c6ae99SBarry Smith 
255e5ab1681SLisandro Dalcin   ierr = DMDAGetInfo(da,NULL,&N,NULL,NULL,NULL,NULL,NULL,&dof,NULL,&bx,NULL,NULL,NULL);CHKERRQ(ierr);
256e5ab1681SLisandro Dalcin   ierr = DMDAGetCorners(da,&istart,NULL,NULL,&isize,NULL,NULL);CHKERRQ(ierr);
25747c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
25847c6ae99SBarry Smith   ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr);
259e5ab1681SLisandro Dalcin   n    = n/dof;
26047c6ae99SBarry Smith 
261832b7cebSLisandro Dalcin   /* Get coordinates of nodes */
2626636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
26347c6ae99SBarry Smith   if (!xcoor) {
264aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr);
2656636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
26647c6ae99SBarry Smith   }
26747c6ae99SBarry Smith   ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr);
268832b7cebSLisandro Dalcin   ierr = DMDAGetCoordinateName(da,0,&xlabel);CHKERRQ(ierr);
269832b7cebSLisandro Dalcin 
270832b7cebSLisandro Dalcin   /* Determine the min and max coordinate in plot */
271e5ab1681SLisandro Dalcin   if (!rank) xmin = PetscRealPart(xg[0]);
272e5ab1681SLisandro Dalcin   if (rank == size-1) xmax = PetscRealPart(xg[n-1]);
27347c6ae99SBarry Smith   ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr);
27447c6ae99SBarry Smith   ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr);
27547c6ae99SBarry Smith 
2764e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
277e5ab1681SLisandro Dalcin   ierr = PetscViewerGetFormat(v,&format);CHKERRQ(ierr);
278e5ab1681SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
279832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
280832b7cebSLisandro Dalcin   if (useports) {
281e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
282e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDrawAxis(v,0,&axis);CHKERRQ(ierr);
283e5ab1681SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
284e5ab1681SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
285e5ab1681SLisandro Dalcin     ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
28673f7a4c5SBarry Smith   }
287e5ab1681SLisandro Dalcin 
288832b7cebSLisandro Dalcin   /* Loop over each field; drawing each in a different window */
28903193ff8SBarry Smith   for (k=0; k<ndisplayfields; k++) {
29003193ff8SBarry Smith     j = displayfields[k];
29147c6ae99SBarry Smith 
292e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
293e5ab1681SLisandro Dalcin     ierr = VecStrideMin(xin,j,NULL,&min);CHKERRQ(ierr);
294e5ab1681SLisandro Dalcin     ierr = VecStrideMax(xin,j,NULL,&max);CHKERRQ(ierr);
295a56202b9SBarry Smith     if (j < nbounds) {
296a56202b9SBarry Smith       min = PetscMin(min,bounds[2*j]);
297a56202b9SBarry Smith       max = PetscMax(max,bounds[2*j+1]);
298a56202b9SBarry Smith     }
299e5ab1681SLisandro Dalcin     if (min == max) {
300e5ab1681SLisandro Dalcin       min -= 1.e-5;
301e5ab1681SLisandro Dalcin       max += 1.e-5;
302e5ab1681SLisandro Dalcin     }
30347c6ae99SBarry Smith 
304832b7cebSLisandro Dalcin     if (useports) {
305832b7cebSLisandro Dalcin       ierr = PetscDrawViewPortsSet(ports,k);CHKERRQ(ierr);
306832b7cebSLisandro Dalcin       ierr = DMDAGetFieldName(da,j,&tlabel);CHKERRQ(ierr);
307832b7cebSLisandro Dalcin     } else {
308832b7cebSLisandro Dalcin       const char *title;
3096083293cSBarry Smith       ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr);
310832b7cebSLisandro Dalcin       ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr);
311832b7cebSLisandro Dalcin       ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr);
312832b7cebSLisandro Dalcin       ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr);
313832b7cebSLisandro Dalcin       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
314832b7cebSLisandro Dalcin       ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
315e5ab1681SLisandro Dalcin       if (!hold) {ierr = PetscDrawClear(draw);CHKERRQ(ierr);}
316832b7cebSLisandro Dalcin     }
317832b7cebSLisandro Dalcin     ierr = PetscDrawAxisSetLabels(axis,tlabel,xlabel,NULL);CHKERRQ(ierr);
318e5ab1681SLisandro Dalcin     ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,min,max);CHKERRQ(ierr);
31947c6ae99SBarry Smith     ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr);
32047c6ae99SBarry Smith 
32147c6ae99SBarry Smith     /* draw local part of vector */
322e5ab1681SLisandro Dalcin     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag);CHKERRQ(ierr);
32347c6ae99SBarry Smith     if (rank < size-1) { /*send value to right */
324e5ab1681SLisandro Dalcin       ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag,comm);CHKERRQ(ierr);
325e5ab1681SLisandro Dalcin       ierr = MPI_Send((void*)&array[j+(n-1)*dof],1,MPIU_REAL,rank+1,tag,comm);CHKERRQ(ierr);
32647c6ae99SBarry Smith     }
32747c6ae99SBarry Smith     if (rank) { /* receive value from left */
328e5ab1681SLisandro Dalcin       ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag,comm,&status);CHKERRQ(ierr);
329e5ab1681SLisandro Dalcin       ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,comm,&status);CHKERRQ(ierr);
330e5ab1681SLisandro Dalcin     }
331e5ab1681SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
332e5ab1681SLisandro Dalcin     if (rank) {
33347c6ae99SBarry Smith       ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr);
334e5ab1681SLisandro Dalcin       if (showmarkers) {ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr);}
33547c6ae99SBarry Smith     }
336e5ab1681SLisandro Dalcin     for (i=1; i<n; i++) {
337e5ab1681SLisandro Dalcin       ierr = PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+dof*(i-1)]),PetscRealPart(xg[i]),PetscRealPart(array[j+dof*i]),PETSC_DRAW_RED);CHKERRQ(ierr);
338e5ab1681SLisandro Dalcin       if (showmarkers) {ierr = PetscDrawMarker(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+dof*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);}
33947c6ae99SBarry Smith     }
340e5ab1681SLisandro Dalcin     if (rank == size-1) {
341e5ab1681SLisandro Dalcin       if (showmarkers) {ierr = PetscDrawMarker(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+dof*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);}
34247c6ae99SBarry Smith     }
343e5ab1681SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
3445b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
34547c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
346832b7cebSLisandro Dalcin     if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);}
34747c6ae99SBarry Smith   }
348832b7cebSLisandro Dalcin   if (useports) {
349832b7cebSLisandro Dalcin     ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
350832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
351832b7cebSLisandro Dalcin   }
352832b7cebSLisandro Dalcin 
353e5ab1681SLisandro Dalcin   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
354cd723cd1SBarry Smith   ierr = PetscFree(displayfields);CHKERRQ(ierr);
35547c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr);
35647c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
35747c6ae99SBarry Smith   PetscFunctionReturn(0);
35847c6ae99SBarry Smith }
35947c6ae99SBarry Smith 
360