xref: /petsc/src/dm/impls/da/gr1.c (revision 832b7ceb95bd65cfd99922d89fc8cdbc69c04de1)
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 #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
18362c83fbSPatrick Sanan .  ymin,ymax - extremes in the y direction (value ignored for 1 dimensional problems)
19362c83fbSPatrick Sanan -  zmin,zmax - extremes in the z direction (value ignored 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;
31bff4a2f0SMatthew G. Knepley   DMBoundaryType   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);
41c593f006SBarry Smith   if (xmax < xmin) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %g %g",(double)xmin,(double)xmax);
42c593f006SBarry Smith   if ((ymax < ymin) && (dim > 1)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %g %g",(double)ymin,(double)ymax);
43c593f006SBarry Smith   if ((zmax < zmin) && (dim > 2)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"zmax must be larger than zmin %g %g",(double)zmin,(double)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);
66bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx  = (xmax-xmin)/(M+1);
6757459e9aSMatthew G Knepley     else                              hx  = (xmax-xmin)/(M ? M : 1);
68bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy  = (ymax-ymin)/(N+1);
6957459e9aSMatthew G Knepley     else                              hy  = (ymax-ymin)/(N ? N : 1);
70bff4a2f0SMatthew G. Knepley     if (bz == DM_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);
1203bb1ff40SBarry 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) {
125bff4a2f0SMatthew G. Knepley     if (bx == DM_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) {
133bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
13447c6ae99SBarry Smith     else hx = (xmax-xmin)/(M-1);
135bff4a2f0SMatthew G. Knepley     if (by == DM_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) {
147bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M);
14847c6ae99SBarry Smith     else hx = (xmax-xmin)/(M-1);
149bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N);
15047c6ae99SBarry Smith     else hy = (ymax-ymin)/(N-1);
151bff4a2f0SMatthew G. Knepley     if (bz == DM_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);
1673bb1ff40SBarry 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"
1749ae14b6eSBarry Smith /*
1759ae14b6eSBarry Smith     Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA
1769ae14b6eSBarry Smith */
1774e6118eeSBarry Smith PetscErrorCode DMDASelectFields(DM da,PetscInt *outfields,PetscInt **fields)
1784e6118eeSBarry Smith {
1794e6118eeSBarry Smith   PetscErrorCode ierr;
1804e6118eeSBarry Smith   PetscInt       step,ndisplayfields,*displayfields,k,j;
1814e6118eeSBarry Smith   PetscBool      flg;
1824e6118eeSBarry Smith 
1834e6118eeSBarry Smith   PetscFunctionBegin;
1844e6118eeSBarry Smith   ierr = DMDAGetInfo(da,0,0,0,0,0,0,0,&step,0,0,0,0,0);CHKERRQ(ierr);
185785e854fSJed Brown   ierr = PetscMalloc1(step,&displayfields);CHKERRQ(ierr);
1864e6118eeSBarry Smith   for (k=0; k<step; k++) displayfields[k] = k;
1874e6118eeSBarry Smith   ndisplayfields = step;
188c5929fdfSBarry Smith   ierr           = PetscOptionsGetIntArray(NULL,NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr);
1894e6118eeSBarry Smith   if (!ndisplayfields) ndisplayfields = step;
1904e6118eeSBarry Smith   if (!flg) {
1914e6118eeSBarry Smith     char       **fields;
1924e6118eeSBarry Smith     const char *fieldname;
1934e6118eeSBarry Smith     PetscInt   nfields = step;
194785e854fSJed Brown     ierr = PetscMalloc1(step,&fields);CHKERRQ(ierr);
195c5929fdfSBarry Smith     ierr = PetscOptionsGetStringArray(NULL,NULL,"-draw_fields_by_name",fields,&nfields,&flg);CHKERRQ(ierr);
1964e6118eeSBarry Smith     if (flg) {
1974e6118eeSBarry Smith       ndisplayfields = 0;
1984e6118eeSBarry Smith       for (k=0; k<nfields;k++) {
1994e6118eeSBarry Smith         for (j=0; j<step; j++) {
2004e6118eeSBarry Smith           ierr = DMDAGetFieldName(da,j,&fieldname);CHKERRQ(ierr);
2014e6118eeSBarry Smith           ierr = PetscStrcmp(fieldname,fields[k],&flg);CHKERRQ(ierr);
2024e6118eeSBarry Smith           if (flg) {
2034e6118eeSBarry Smith             goto found;
2044e6118eeSBarry Smith           }
2054e6118eeSBarry Smith         }
206ce94432eSBarry Smith         SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Unknown fieldname %s",fields[k]);
2074e6118eeSBarry Smith found:  displayfields[ndisplayfields++] = j;
2084e6118eeSBarry Smith       }
2094e6118eeSBarry Smith     }
210e3f3e4b6SBarry Smith     for (k=0; k<nfields; k++) {
211e3f3e4b6SBarry Smith       ierr = PetscFree(fields[k]);CHKERRQ(ierr);
212e3f3e4b6SBarry Smith     }
2134e6118eeSBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
2144e6118eeSBarry Smith   }
2154e6118eeSBarry Smith   *fields    = displayfields;
2164e6118eeSBarry Smith   *outfields = ndisplayfields;
2174e6118eeSBarry Smith   PetscFunctionReturn(0);
2184e6118eeSBarry Smith }
2194e6118eeSBarry Smith 
2209804daf3SBarry Smith #include <petscdraw.h>
22173f7a4c5SBarry Smith 
2224e6118eeSBarry Smith #undef __FUNCT__
22347c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA1d"
22447c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)
22547c6ae99SBarry Smith {
2269a42bb27SBarry Smith   DM                da;
22747c6ae99SBarry Smith   PetscErrorCode    ierr;
228e5ab1681SLisandro Dalcin   PetscMPIInt       rank,size,tag;
229e5ab1681SLisandro Dalcin   PetscInt          i,n,N,dof,istart,isize,j,nbounds;
23047c6ae99SBarry Smith   MPI_Status        status;
231e5ab1681SLisandro Dalcin   PetscReal         min,max,xmin = 0.0,xmax = 0.0,tmp = 0.0,xgtmp = 0.0;
23247c6ae99SBarry Smith   const PetscScalar *array,*xg;
23347c6ae99SBarry Smith   PetscDraw         draw;
234e5ab1681SLisandro Dalcin   PetscBool         isnull,useports = PETSC_FALSE,showmarkers = PETSC_FALSE;
23547c6ae99SBarry Smith   MPI_Comm          comm;
23647c6ae99SBarry Smith   PetscDrawAxis     axis;
23747c6ae99SBarry Smith   Vec               xcoor;
238bff4a2f0SMatthew G. Knepley   DMBoundaryType    bx;
239e5ab1681SLisandro Dalcin   const char        *tlabel = NULL,*xlabel = NULL;
240a56202b9SBarry Smith   const PetscReal   *bounds;
24103193ff8SBarry Smith   PetscInt          *displayfields;
24203193ff8SBarry Smith   PetscInt          k,ndisplayfields;
2434e6118eeSBarry Smith   PetscBool         hold;
244e5ab1681SLisandro Dalcin   PetscDrawViewPorts *ports = NULL;
245e5ab1681SLisandro Dalcin   PetscViewerFormat  format;
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith   PetscFunctionBegin;
24847c6ae99SBarry Smith   ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
24947c6ae99SBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
250a56202b9SBarry Smith   ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr);
25147c6ae99SBarry Smith 
252c688c046SMatthew G Knepley   ierr = VecGetDM(xin,&da);CHKERRQ(ierr);
253ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
254e5ab1681SLisandro Dalcin   ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr);
255e5ab1681SLisandro Dalcin   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
256e5ab1681SLisandro Dalcin   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
25747c6ae99SBarry Smith 
258c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_vec_use_markers",&showmarkers,NULL);CHKERRQ(ierr);
25947c6ae99SBarry Smith 
260e5ab1681SLisandro Dalcin   ierr = DMDAGetInfo(da,NULL,&N,NULL,NULL,NULL,NULL,NULL,&dof,NULL,&bx,NULL,NULL,NULL);CHKERRQ(ierr);
261e5ab1681SLisandro Dalcin   ierr = DMDAGetCorners(da,&istart,NULL,NULL,&isize,NULL,NULL);CHKERRQ(ierr);
26247c6ae99SBarry Smith   ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr);
26347c6ae99SBarry Smith   ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr);
264e5ab1681SLisandro Dalcin   n    = n/dof;
26547c6ae99SBarry Smith 
266*832b7cebSLisandro Dalcin   /* Get coordinates of nodes */
2676636e97aSMatthew G Knepley   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
26847c6ae99SBarry Smith   if (!xcoor) {
269aa219208SBarry Smith     ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr);
2706636e97aSMatthew G Knepley     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
27147c6ae99SBarry Smith   }
27247c6ae99SBarry Smith   ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr);
273*832b7cebSLisandro Dalcin   ierr = DMDAGetCoordinateName(da,0,&xlabel);CHKERRQ(ierr);
274*832b7cebSLisandro Dalcin 
275*832b7cebSLisandro Dalcin   /* Determine the min and max coordinate in plot */
276e5ab1681SLisandro Dalcin   if (!rank) xmin = PetscRealPart(xg[0]);
277e5ab1681SLisandro Dalcin   if (rank == size-1) xmax = PetscRealPart(xg[n-1]);
27847c6ae99SBarry Smith   ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr);
27947c6ae99SBarry Smith   ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr);
28047c6ae99SBarry Smith 
2814e6118eeSBarry Smith   ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr);
282e5ab1681SLisandro Dalcin   ierr = PetscViewerGetFormat(v,&format);CHKERRQ(ierr);
283e5ab1681SLisandro Dalcin   ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr);
284*832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
285*832b7cebSLisandro Dalcin   if (useports) {
286e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
287e5ab1681SLisandro Dalcin     ierr = PetscViewerDrawGetDrawAxis(v,0,&axis);CHKERRQ(ierr);
288e5ab1681SLisandro Dalcin     ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
289e5ab1681SLisandro Dalcin     ierr = PetscDrawClear(draw);CHKERRQ(ierr);
290e5ab1681SLisandro Dalcin     ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr);
29173f7a4c5SBarry Smith   }
292e5ab1681SLisandro Dalcin 
293*832b7cebSLisandro Dalcin   /* Loop over each field; drawing each in a different window */
29403193ff8SBarry Smith   for (k=0; k<ndisplayfields; k++) {
29503193ff8SBarry Smith     j = displayfields[k];
29647c6ae99SBarry Smith 
297e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
298e5ab1681SLisandro Dalcin     ierr = VecStrideMin(xin,j,NULL,&min);CHKERRQ(ierr);
299e5ab1681SLisandro Dalcin     ierr = VecStrideMax(xin,j,NULL,&max);CHKERRQ(ierr);
300a56202b9SBarry Smith     if (j < nbounds) {
301a56202b9SBarry Smith       min = PetscMin(min,bounds[2*j]);
302a56202b9SBarry Smith       max = PetscMax(max,bounds[2*j+1]);
303a56202b9SBarry Smith     }
304e5ab1681SLisandro Dalcin     if (min == max) {
305e5ab1681SLisandro Dalcin       min -= 1.e-5;
306e5ab1681SLisandro Dalcin       max += 1.e-5;
307e5ab1681SLisandro Dalcin     }
30847c6ae99SBarry Smith 
309*832b7cebSLisandro Dalcin     if (useports) {
310*832b7cebSLisandro Dalcin       ierr = PetscDrawViewPortsSet(ports,k);CHKERRQ(ierr);
311*832b7cebSLisandro Dalcin       ierr = DMDAGetFieldName(da,j,&tlabel);CHKERRQ(ierr);
312*832b7cebSLisandro Dalcin     } else {
313*832b7cebSLisandro Dalcin       const char *title;
3146083293cSBarry Smith       ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr);
315*832b7cebSLisandro Dalcin       ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr);
316*832b7cebSLisandro Dalcin       ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr);
317*832b7cebSLisandro Dalcin       ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr);
318*832b7cebSLisandro Dalcin       if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);}
319*832b7cebSLisandro Dalcin       ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr);
320e5ab1681SLisandro Dalcin       if (!hold) {ierr = PetscDrawClear(draw);CHKERRQ(ierr);}
321*832b7cebSLisandro Dalcin     }
322*832b7cebSLisandro Dalcin     ierr = PetscDrawAxisSetLabels(axis,tlabel,xlabel,NULL);CHKERRQ(ierr);
323e5ab1681SLisandro Dalcin     ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,min,max);CHKERRQ(ierr);
32447c6ae99SBarry Smith     ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr);
32547c6ae99SBarry Smith 
32647c6ae99SBarry Smith     /* draw local part of vector */
327e5ab1681SLisandro Dalcin     ierr = PetscObjectGetNewTag((PetscObject)xin,&tag);CHKERRQ(ierr);
32847c6ae99SBarry Smith     if (rank < size-1) { /*send value to right */
329e5ab1681SLisandro Dalcin       ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag,comm);CHKERRQ(ierr);
330e5ab1681SLisandro Dalcin       ierr = MPI_Send((void*)&array[j+(n-1)*dof],1,MPIU_REAL,rank+1,tag,comm);CHKERRQ(ierr);
33147c6ae99SBarry Smith     }
33247c6ae99SBarry Smith     if (rank) { /* receive value from left */
333e5ab1681SLisandro Dalcin       ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag,comm,&status);CHKERRQ(ierr);
334e5ab1681SLisandro Dalcin       ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,comm,&status);CHKERRQ(ierr);
335e5ab1681SLisandro Dalcin     }
336e5ab1681SLisandro Dalcin     ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr);
337e5ab1681SLisandro Dalcin     if (rank) {
33847c6ae99SBarry Smith       ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr);
339e5ab1681SLisandro Dalcin       if (showmarkers) {ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr);}
34047c6ae99SBarry Smith     }
341e5ab1681SLisandro Dalcin     for (i=1; i<n; i++) {
342e5ab1681SLisandro 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);
343e5ab1681SLisandro Dalcin       if (showmarkers) {ierr = PetscDrawMarker(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+dof*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);}
34447c6ae99SBarry Smith     }
345e5ab1681SLisandro Dalcin     if (rank == size-1) {
346e5ab1681SLisandro Dalcin       if (showmarkers) {ierr = PetscDrawMarker(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+dof*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);}
34747c6ae99SBarry Smith     }
348e5ab1681SLisandro Dalcin     ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr);
3495b399a63SLisandro Dalcin     ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
35047c6ae99SBarry Smith     ierr = PetscDrawPause(draw);CHKERRQ(ierr);
351*832b7cebSLisandro Dalcin     if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);}
35247c6ae99SBarry Smith   }
353*832b7cebSLisandro Dalcin   if (useports) {
354*832b7cebSLisandro Dalcin     ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr);
355*832b7cebSLisandro Dalcin     ierr = PetscDrawSave(draw);CHKERRQ(ierr);
356*832b7cebSLisandro Dalcin   }
357*832b7cebSLisandro Dalcin 
358e5ab1681SLisandro Dalcin   ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr);
359cd723cd1SBarry Smith   ierr = PetscFree(displayfields);CHKERRQ(ierr);
36047c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr);
36147c6ae99SBarry Smith   ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr);
36247c6ae99SBarry Smith   PetscFunctionReturn(0);
36347c6ae99SBarry Smith }
36447c6ae99SBarry Smith 
365