xref: /petsc/src/dm/impls/da/gr1.c (revision 3bb1ff401821b9e2ae019d3e61bc8ab4bd4e59d5)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
3aa219208SBarry Smith    Plots vectors obtained with DMDACreate1d()
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>      /*I  "petscdmda.h"   I*/
747c6ae99SBarry Smith 
847c6ae99SBarry Smith #undef __FUNCT__
9aa219208SBarry Smith #define __FUNCT__ "DMDASetUniformCoordinates"
1047c6ae99SBarry Smith /*@
11aa219208SBarry Smith     DMDASetUniformCoordinates - Sets a DMDA coordinates to be a uniform grid
1247c6ae99SBarry Smith 
13aa219208SBarry Smith   Collective on DMDA
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith   Input Parameters:
1647c6ae99SBarry Smith +  da - the distributed array object
1747c6ae99SBarry Smith .  xmin,xmax - extremes in the x direction
180298fd71SBarry Smith .  ymin,ymax - extremes in the y direction (use NULL for 1 dimensional problems)
190298fd71SBarry Smith -  zmin,zmax - extremes in the z direction (use NULL for 1 or 2 dimensional problems)
2047c6ae99SBarry Smith 
2147c6ae99SBarry Smith   Level: beginner
2247c6ae99SBarry Smith 
232150357eSBarry Smith .seealso: DMSetCoordinates(), DMGetCoordinates(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith @*/
267087cfbeSBarry Smith PetscErrorCode  DMDASetUniformCoordinates(DM da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)
2747c6ae99SBarry Smith {
2847c6ae99SBarry Smith   MPI_Comm         comm;
2957459e9aSMatthew G Knepley   PetscSection     section;
309a42bb27SBarry Smith   DM               cda;
311321219cSEthan Coon   DMDABoundaryType bx,by,bz;
3247c6ae99SBarry Smith   Vec              xcoor;
3347c6ae99SBarry Smith   PetscScalar      *coors;
3447c6ae99SBarry Smith   PetscReal        hx,hy,hz_;
3547c6ae99SBarry Smith   PetscInt         i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt;
3647c6ae99SBarry Smith   PetscErrorCode   ierr;
3747c6ae99SBarry Smith 
3847c6ae99SBarry Smith   PetscFunctionBegin;
395b990481SMatthew G Knepley   PetscValidHeaderSpecific(da,DM_CLASSID,1);
40821c30cfSPeter Brune   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr);
41ce94432eSBarry Smith   if (xmax <= xmin) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %G %G",xmin,xmax);
42ce94432eSBarry Smith   if ((ymax <= ymin) && (dim > 1)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax);
43ce94432eSBarry Smith   if ((zmax <= zmin) && (dim > 2)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %G %G",zmin,zmax);
4447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
4557459e9aSMatthew G Knepley   ierr = DMGetDefaultSection(da,&section);CHKERRQ(ierr);
46aa219208SBarry Smith   ierr = DMDAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);CHKERRQ(ierr);
476636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(da, &cda);CHKERRQ(ierr);
4857459e9aSMatthew G Knepley   if (section) {
4957459e9aSMatthew G Knepley     /* This would be better as a vector, but this is compatible */
5057459e9aSMatthew G Knepley     PetscInt numComp[3]      = {1, 1, 1};
5157459e9aSMatthew G Knepley     PetscInt numVertexDof[3] = {1, 1, 1};
5257459e9aSMatthew G Knepley 
5357459e9aSMatthew G Knepley     ierr = DMDASetFieldName(cda, 0, "x");CHKERRQ(ierr);
5457459e9aSMatthew G Knepley     if (dim > 1) {ierr = DMDASetFieldName(cda, 1, "y");CHKERRQ(ierr);}
5557459e9aSMatthew G Knepley     if (dim > 2) {ierr = DMDASetFieldName(cda, 2, "z");CHKERRQ(ierr);}
560298fd71SBarry Smith     ierr = DMDACreateSection(cda, numComp, numVertexDof, NULL, NULL);CHKERRQ(ierr);
5757459e9aSMatthew G Knepley   }
58564755cdSBarry Smith   ierr = DMCreateGlobalVector(cda, &xcoor);CHKERRQ(ierr);
5957459e9aSMatthew G Knepley   if (section) {
6057459e9aSMatthew G Knepley     PetscSection csection;
6157459e9aSMatthew G Knepley     PetscInt     vStart, vEnd;
6257459e9aSMatthew G Knepley 
6357459e9aSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(cda,&csection);CHKERRQ(ierr);
6457459e9aSMatthew G Knepley     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
6557459e9aSMatthew G Knepley     ierr = DMDAGetHeightStratum(da, dim, &vStart, &vEnd);CHKERRQ(ierr);
6657459e9aSMatthew G Knepley     if (bx == DMDA_BOUNDARY_PERIODIC) hx  = (xmax-xmin)/(M+1);
6757459e9aSMatthew G Knepley     else                              hx  = (xmax-xmin)/(M ? M : 1);
6857459e9aSMatthew G Knepley     if (by == DMDA_BOUNDARY_PERIODIC) hy  = (ymax-ymin)/(N+1);
6957459e9aSMatthew G Knepley     else                              hy  = (ymax-ymin)/(N ? N : 1);
7057459e9aSMatthew G Knepley     if (bz == DMDA_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P+1);
7157459e9aSMatthew G Knepley     else                              hz_ = (zmax-zmin)/(P ? P : 1);
7257459e9aSMatthew G Knepley     switch (dim) {
7357459e9aSMatthew G Knepley     case 1:
7457459e9aSMatthew G Knepley       for (i = 0; i < isize+1; ++i) {
7557459e9aSMatthew G Knepley         PetscInt v = i+vStart, dof, off;
7657459e9aSMatthew G Knepley 
7757459e9aSMatthew G Knepley         ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
7857459e9aSMatthew G Knepley         ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
7957459e9aSMatthew G Knepley         if (off >= 0) {
8057459e9aSMatthew G Knepley           coors[off] = xmin + hx*(i+istart);
8157459e9aSMatthew G Knepley         }
8257459e9aSMatthew G Knepley       }
8357459e9aSMatthew G Knepley       break;
8457459e9aSMatthew G Knepley     case 2:
8557459e9aSMatthew G Knepley       for (j = 0; j < jsize+1; ++j) {
8657459e9aSMatthew G Knepley         for (i = 0; i < isize+1; ++i) {
8757459e9aSMatthew G Knepley           PetscInt v = j*(isize+1)+i+vStart, dof, off;
8857459e9aSMatthew G Knepley 
8957459e9aSMatthew G Knepley           ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
9057459e9aSMatthew G Knepley           ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
9157459e9aSMatthew G Knepley           if (off >= 0) {
9257459e9aSMatthew G Knepley             coors[off+0] = xmin + hx*(i+istart);
9357459e9aSMatthew G Knepley             coors[off+1] = ymin + hy*(j+jstart);
9457459e9aSMatthew G Knepley           }
9557459e9aSMatthew G Knepley         }
9657459e9aSMatthew G Knepley       }
9757459e9aSMatthew G Knepley       break;
9857459e9aSMatthew G Knepley     case 3:
9957459e9aSMatthew G Knepley       for (k = 0; k < ksize+1; ++k) {
10057459e9aSMatthew G Knepley         for (j = 0; j < jsize+1; ++j) {
10157459e9aSMatthew G Knepley           for (i = 0; i < isize+1; ++i) {
10257459e9aSMatthew G Knepley             PetscInt v = (k*(jsize+1)+j)*(isize+1)+i+vStart, dof, off;
10357459e9aSMatthew G Knepley 
10457459e9aSMatthew G Knepley             ierr = PetscSectionGetDof(csection, v, &dof);CHKERRQ(ierr);
10557459e9aSMatthew G Knepley             ierr = PetscSectionGetOffset(csection, v, &off);CHKERRQ(ierr);
10657459e9aSMatthew G Knepley             if (off >= 0) {
10757459e9aSMatthew G Knepley               coors[off+0] = xmin + hx*(i+istart);
10857459e9aSMatthew G Knepley               coors[off+1] = ymin + hy*(j+jstart);
10957459e9aSMatthew G Knepley               coors[off+2] = zmin + hz_*(k+kstart);
11057459e9aSMatthew G Knepley             }
11157459e9aSMatthew G Knepley           }
11257459e9aSMatthew G Knepley         }
11357459e9aSMatthew G Knepley       }
11457459e9aSMatthew G Knepley       break;
11557459e9aSMatthew G Knepley     default:
116ce94432eSBarry Smith       SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
11757459e9aSMatthew G Knepley     }
11857459e9aSMatthew G Knepley     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
1196636e97aSMatthew G Knepley     ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr);
120*3bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)xcoor);CHKERRQ(ierr);
12157459e9aSMatthew G Knepley     ierr = VecDestroy(&xcoor);CHKERRQ(ierr);
12257459e9aSMatthew G Knepley     PetscFunctionReturn(0);
12357459e9aSMatthew G Knepley   }
12447c6ae99SBarry Smith   if (dim == 1) {
1251321219cSEthan Coon     if (bx == DMDA_BOUNDARY_PERIODIC) hx = (xmax-xmin)/M;
126ce00eea3SSatish Balay     else hx = (xmax-xmin)/(M-1);
12747c6ae99SBarry Smith     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
12847c6ae99SBarry Smith     for (i=0; i<isize; i++) {
12947c6ae99SBarry Smith       coors[i] = xmin + hx*(i+istart);
13047c6ae99SBarry Smith     }
13147c6ae99SBarry Smith     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
13247c6ae99SBarry Smith   } else if (dim == 2) {
1331321219cSEthan Coon     if (bx == DMDA_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
13447c6ae99SBarry Smith     else hx = (xmax-xmin)/(M-1);
1351321219cSEthan Coon     if (by == DMDA_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
13647c6ae99SBarry Smith     else hy = (ymax-ymin)/(N-1);
13747c6ae99SBarry Smith     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
13847c6ae99SBarry Smith     cnt  = 0;
13947c6ae99SBarry Smith     for (j=0; j<jsize; j++) {
14047c6ae99SBarry Smith       for (i=0; i<isize; i++) {
14147c6ae99SBarry Smith         coors[cnt++] = xmin + hx*(i+istart);
14247c6ae99SBarry Smith         coors[cnt++] = ymin + hy*(j+jstart);
14347c6ae99SBarry Smith       }
14447c6ae99SBarry Smith     }
14547c6ae99SBarry Smith     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
14647c6ae99SBarry Smith   } else if (dim == 3) {
1471321219cSEthan Coon     if (bx == DMDA_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
14847c6ae99SBarry Smith     else hx = (xmax-xmin)/(M-1);
1491321219cSEthan Coon     if (by == DMDA_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
15047c6ae99SBarry Smith     else hy = (ymax-ymin)/(N-1);
1511321219cSEthan Coon     if (bz == DMDA_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P);
15247c6ae99SBarry Smith     else hz_ = (zmax-zmin)/(P-1);
15347c6ae99SBarry Smith     ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr);
15447c6ae99SBarry Smith     cnt  = 0;
15547c6ae99SBarry Smith     for (k=0; k<ksize; k++) {
15647c6ae99SBarry Smith       for (j=0; j<jsize; j++) {
15747c6ae99SBarry Smith         for (i=0; i<isize; i++) {
15847c6ae99SBarry Smith           coors[cnt++] = xmin + hx*(i+istart);
15947c6ae99SBarry Smith           coors[cnt++] = ymin + hy*(j+jstart);
16047c6ae99SBarry Smith           coors[cnt++] = zmin + hz_*(k+kstart);
16147c6ae99SBarry Smith         }
16247c6ae99SBarry Smith       }
16347c6ae99SBarry Smith     }
16447c6ae99SBarry Smith     ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr);
165ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim);
1666636e97aSMatthew G Knepley   ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr);
167*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)xcoor);CHKERRQ(ierr);
168fcfd50ebSBarry Smith   ierr = VecDestroy(&xcoor);CHKERRQ(ierr);
16947c6ae99SBarry Smith   PetscFunctionReturn(0);
17047c6ae99SBarry Smith }
17147c6ae99SBarry Smith 
17247c6ae99SBarry Smith #undef __FUNCT__
1734e6118eeSBarry Smith #define __FUNCT__ "DMDASelectFields"
1744e6118eeSBarry Smith PetscErrorCode DMDASelectFields(DM da,PetscInt *outfields,PetscInt **fields)
1754e6118eeSBarry Smith {
1764e6118eeSBarry Smith   PetscErrorCode ierr;
1774e6118eeSBarry Smith   PetscInt       step,ndisplayfields,*displayfields,k,j;
1784e6118eeSBarry Smith   PetscBool      flg;
1794e6118eeSBarry Smith 
1804e6118eeSBarry Smith   PetscFunctionBegin;
1814e6118eeSBarry Smith   ierr = DMDAGetInfo(da,0,0,0,0,0,0,0,&step,0,0,0,0,0);CHKERRQ(ierr);
1824e6118eeSBarry Smith   ierr = PetscMalloc(step*sizeof(PetscInt),&displayfields);CHKERRQ(ierr);
1834e6118eeSBarry Smith   for (k=0; k<step; k++) displayfields[k] = k;
1844e6118eeSBarry Smith   ndisplayfields = step;
1850298fd71SBarry Smith   ierr           = PetscOptionsGetIntArray(NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr);
1864e6118eeSBarry Smith   if (!ndisplayfields) ndisplayfields = step;
1874e6118eeSBarry Smith   if (!flg) {
1884e6118eeSBarry Smith     char       **fields;
1894e6118eeSBarry Smith     const char *fieldname;
1904e6118eeSBarry Smith     PetscInt   nfields = step;
1914e6118eeSBarry Smith     ierr = PetscMalloc(step*sizeof(char*),&fields);CHKERRQ(ierr);
1920298fd71SBarry Smith     ierr = PetscOptionsGetStringArray(NULL,"-draw_fields_by_name",fields,&nfields,&flg);CHKERRQ(ierr);
1934e6118eeSBarry Smith     if (flg) {
1944e6118eeSBarry Smith       ndisplayfields = 0;
1954e6118eeSBarry Smith       for (k=0; k<nfields;k++) {
1964e6118eeSBarry Smith         for (j=0; j<step; j++) {
1974e6118eeSBarry Smith           ierr = DMDAGetFieldName(da,j,&fieldname);CHKERRQ(ierr);
1984e6118eeSBarry Smith           ierr = PetscStrcmp(fieldname,fields[k],&flg);CHKERRQ(ierr);
1994e6118eeSBarry Smith           if (flg) {
2004e6118eeSBarry Smith             goto found;
2014e6118eeSBarry Smith           }
2024e6118eeSBarry Smith         }
203ce94432eSBarry Smith         SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Unknown fieldname %s",fields[k]);
2044e6118eeSBarry Smith found:  displayfields[ndisplayfields++] = j;
2054e6118eeSBarry Smith       }
2064e6118eeSBarry Smith     }
207e3f3e4b6SBarry Smith     for (k=0; k<nfields; k++) {
208e3f3e4b6SBarry Smith       ierr = PetscFree(fields[k]);CHKERRQ(ierr);
209e3f3e4b6SBarry Smith     }
2104e6118eeSBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
2114e6118eeSBarry Smith   }
2124e6118eeSBarry Smith   *fields    = displayfields;
2134e6118eeSBarry Smith   *outfields = ndisplayfields;
2144e6118eeSBarry Smith   PetscFunctionReturn(0);
2154e6118eeSBarry Smith }
2164e6118eeSBarry Smith 
2179804daf3SBarry Smith #include <petscdraw.h>
2184e6118eeSBarry Smith #undef __FUNCT__
21947c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA1d"
22047c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
22147c6ae99SBarry Smith {
2229a42bb27SBarry Smith   DM                da;
22347c6ae99SBarry Smith   PetscErrorCode    ierr;
22447c6ae99SBarry Smith   PetscMPIInt       rank,size,tag1,tag2;
225a56202b9SBarry Smith   PetscInt          i,n,N,step,istart,isize,j,nbounds;
22647c6ae99SBarry Smith   MPI_Status        status;
22747c6ae99SBarry Smith   PetscReal         coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp;
22847c6ae99SBarry Smith   const PetscScalar *array,*xg;
22947c6ae99SBarry Smith   PetscDraw         draw;
23047c6ae99SBarry Smith   PetscBool         isnull,showpoints = PETSC_FALSE;
23147c6ae99SBarry Smith   MPI_Comm          comm;
23247c6ae99SBarry Smith   PetscDrawAxis     axis;
23347c6ae99SBarry Smith   Vec               xcoor;
2341321219cSEthan Coon   DMDABoundaryType  bx;
235a56202b9SBarry Smith   const PetscReal   *bounds;
23603193ff8SBarry Smith   PetscInt          *displayfields;
23703193ff8SBarry Smith   PetscInt          k,ndisplayfields;
2384e6118eeSBarry Smith   PetscBool         hold;
23947c6ae99SBarry Smith 
24047c6ae99SBarry Smith   PetscFunctionBegin;
24147c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
24247c6ae99SBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
243a56202b9SBarry Smith   ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr);
24447c6ae99SBarry Smith 
245c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
246ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
24747c6ae99SBarry Smith 
2480298fd71SBarry Smith   ierr = PetscOptionsGetBool(NULL,"-draw_vec_mark_points",&showpoints,NULL);CHKERRQ(ierr);
24947c6ae99SBarry Smith 
2501321219cSEthan Coon   ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&bx,0,0,0);CHKERRQ(ierr);
251aa219208SBarry Smith   ierr = DMDAGetCorners(da,&istart,0,0,&isize,0,0);CHKERRQ(ierr);
25247c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
25347c6ae99SBarry Smith   ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr);
25447c6ae99SBarry Smith   n    = n/step;
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith   /* get coordinates of nodes */
2576636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
25847c6ae99SBarry Smith   if (!xcoor) {
259aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr);
2606636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
26147c6ae99SBarry Smith   }
26247c6ae99SBarry Smith   ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr);
26347c6ae99SBarry Smith 
26447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
26547c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
26647c6ae99SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith   /*
26947c6ae99SBarry Smith       Determine the min and max x coordinate in plot
27047c6ae99SBarry Smith   */
27147c6ae99SBarry Smith   if (!rank) {
27247c6ae99SBarry Smith     xmin = PetscRealPart(xg[0]);
27347c6ae99SBarry Smith   }
27447c6ae99SBarry Smith   if (rank == size-1) {
27547c6ae99SBarry Smith     xmax = PetscRealPart(xg[n-1]);
27647c6ae99SBarry Smith   }
27747c6ae99SBarry Smith   ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr);
27847c6ae99SBarry Smith   ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr);
27947c6ae99SBarry Smith 
2804e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
28103193ff8SBarry Smith   for (k=0; k<ndisplayfields; k++) {
28203193ff8SBarry Smith     j    = displayfields[k];
28303193ff8SBarry Smith     ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr);
28447c6ae99SBarry Smith     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
28547c6ae99SBarry Smith 
28647c6ae99SBarry Smith     /*
28747c6ae99SBarry Smith         Determine the min and max y coordinate in plot
28847c6ae99SBarry Smith     */
28947c6ae99SBarry Smith     min = 1.e20; max = -1.e20;
29047c6ae99SBarry Smith     for (i=0; i<n; i++) {
29147c6ae99SBarry Smith       if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]);
29247c6ae99SBarry Smith       if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]);
29347c6ae99SBarry Smith     }
29447c6ae99SBarry Smith     if (min + 1.e-10 > max) {
29547c6ae99SBarry Smith       min -= 1.e-5;
29647c6ae99SBarry Smith       max += 1.e-5;
29747c6ae99SBarry Smith     }
298a56202b9SBarry Smith     if (j < nbounds) {
299a56202b9SBarry Smith       min = PetscMin(min,bounds[2*j]);
300a56202b9SBarry Smith       max = PetscMax(max,bounds[2*j+1]);
301a56202b9SBarry Smith     }
302a56202b9SBarry Smith 
303d9822059SBarry Smith     ierr = MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPIU_MIN,0,comm);CHKERRQ(ierr);
304d9822059SBarry Smith     ierr = MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPIU_MAX,0,comm);CHKERRQ(ierr);
30547c6ae99SBarry Smith 
3066083293cSBarry Smith     ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr);
3076083293cSBarry Smith     if (!hold) {
30847c6ae99SBarry Smith       ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr);
3096083293cSBarry Smith     }
31003193ff8SBarry Smith     ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr);
311*3bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)draw,(PetscObject)axis);CHKERRQ(ierr);
31247c6ae99SBarry Smith     if (!rank) {
31347c6ae99SBarry Smith       const char *title;
31447c6ae99SBarry Smith 
31547c6ae99SBarry Smith       ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr);
31647c6ae99SBarry Smith       ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr);
31747c6ae99SBarry Smith       ierr = PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);CHKERRQ(ierr);
318aa219208SBarry Smith       ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr);
31947c6ae99SBarry Smith       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
32047c6ae99SBarry Smith     }
32147c6ae99SBarry Smith     ierr = MPI_Bcast(coors,4,MPIU_REAL,0,comm);CHKERRQ(ierr);
32247c6ae99SBarry Smith     if (rank) {
32347c6ae99SBarry Smith       ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr);
32447c6ae99SBarry Smith     }
32547c6ae99SBarry Smith 
32647c6ae99SBarry Smith     /* draw local part of vector */
3277d1ecd34SBarry Smith     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag1);CHKERRQ(ierr);
3287d1ecd34SBarry Smith     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag2);CHKERRQ(ierr);
32947c6ae99SBarry Smith     if (rank < size-1) { /*send value to right */
33047c6ae99SBarry Smith       ierr = MPI_Send((void*)&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr);
33147c6ae99SBarry Smith       ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr);
33247c6ae99SBarry Smith     }
3331321219cSEthan Coon     if (!rank && bx == DMDA_BOUNDARY_PERIODIC && size > 1) { /* first processor sends first value to last */
33447c6ae99SBarry Smith       ierr = MPI_Send((void*)&array[j],1,MPIU_REAL,size-1,tag2,comm);CHKERRQ(ierr);
33547c6ae99SBarry Smith     }
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith     for (i=1; i<n; i++) {
33847c6ae99SBarry Smith       ierr = PetscDrawLine(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PetscRealPart(xg[i]),PetscRealPart(array[j+step*i]),PETSC_DRAW_RED);CHKERRQ(ierr);
33947c6ae99SBarry Smith       if (showpoints) {
34047c6ae99SBarry Smith         ierr = PetscDrawPoint(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);
34147c6ae99SBarry Smith       }
34247c6ae99SBarry Smith     }
34347c6ae99SBarry Smith     if (rank) { /* receive value from left */
34447c6ae99SBarry Smith       ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr);
34547c6ae99SBarry Smith       ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr);
34647c6ae99SBarry Smith       ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr);
34747c6ae99SBarry Smith       if (showpoints) {
34847c6ae99SBarry Smith         ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr);
34947c6ae99SBarry Smith       }
35047c6ae99SBarry Smith     }
3511321219cSEthan Coon     if (rank == size-1 && bx == DMDA_BOUNDARY_PERIODIC && size > 1) {
35247c6ae99SBarry Smith       ierr = MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);CHKERRQ(ierr);
35317eed501SBarry Smith       /* If the mesh is not uniform we do not know the mesh spacing between the last point on the right and the first ghost point */
35417eed501SBarry Smith       ierr = PetscDrawLine(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+step*(n-1)]),PetscRealPart(xg[n-1]+(xg[n-1]-xg[n-2])),tmp,PETSC_DRAW_RED);CHKERRQ(ierr);
35547c6ae99SBarry Smith       if (showpoints) {
35647c6ae99SBarry Smith         ierr = PetscDrawPoint(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);
35747c6ae99SBarry Smith       }
35847c6ae99SBarry Smith     }
35947c6ae99SBarry Smith     ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr);
36047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
36147c6ae99SBarry Smith   }
362cd723cd1SBarry Smith   ierr = PetscFree(displayfields);CHKERRQ(ierr);
36347c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr);
36447c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
36547c6ae99SBarry Smith   PetscFunctionReturn(0);
36647c6ae99SBarry Smith }
36747c6ae99SBarry Smith 
368