147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 3aa219208SBarry Smith Plots vectors obtained with DMDACreate1d() 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6c6db04a5SJed Brown #include <petscdmda.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 1847c6ae99SBarry Smith . ymin,ymax - extremes in the y direction (use PETSC_NULL for 1 dimensional problems) 1947c6ae99SBarry Smith - zmin,zmax - extremes in the z direction (use PETSC_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); 4147c6ae99SBarry Smith if (xmax <= xmin) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %G %G",xmin,xmax); 4257459e9aSMatthew G Knepley if ((ymax <= ymin) && (dim > 1)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_ARG_INCOMP,"ymax must be larger than ymin %G %G",ymin,ymax); 4357459e9aSMatthew G Knepley if ((zmax <= zmin) && (dim > 2)) SETERRQ2(((PetscObject)da)->comm,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,§ion);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);} 5657459e9aSMatthew G Knepley ierr = DMDACreateSection(cda, numComp, numVertexDof, PETSC_NULL, PETSC_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: 11657459e9aSMatthew G Knepley SETERRQ1(((PetscObject)da)->comm,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); 12057459e9aSMatthew G Knepley ierr = PetscLogObjectParent(da,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); 165bfec8eecSBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim); 1666636e97aSMatthew G Knepley ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr); 16747c6ae99SBarry Smith ierr = PetscLogObjectParent(da,xcoor);CHKERRQ(ierr); 168fcfd50ebSBarry Smith ierr = VecDestroy(&xcoor);CHKERRQ(ierr); 16947c6ae99SBarry Smith PetscFunctionReturn(0); 17047c6ae99SBarry Smith } 17147c6ae99SBarry Smith 17247c6ae99SBarry Smith #undef __FUNCT__ 17347c6ae99SBarry Smith #define __FUNCT__ "VecView_MPI_Draw_DA1d" 17447c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v) 17547c6ae99SBarry Smith { 1769a42bb27SBarry Smith DM da; 17747c6ae99SBarry Smith PetscErrorCode ierr; 17847c6ae99SBarry Smith PetscMPIInt rank,size,tag1,tag2; 179a56202b9SBarry Smith PetscInt i,n,N,step,istart,isize,j,nbounds; 18047c6ae99SBarry Smith MPI_Status status; 18147c6ae99SBarry Smith PetscReal coors[4],ymin,ymax,min,max,xmin,xmax,tmp,xgtmp; 18247c6ae99SBarry Smith const PetscScalar *array,*xg; 18347c6ae99SBarry Smith PetscDraw draw; 18447c6ae99SBarry Smith PetscBool isnull,showpoints = PETSC_FALSE; 18547c6ae99SBarry Smith MPI_Comm comm; 18647c6ae99SBarry Smith PetscDrawAxis axis; 18747c6ae99SBarry Smith Vec xcoor; 1881321219cSEthan Coon DMDABoundaryType bx; 189a56202b9SBarry Smith const PetscReal *bounds; 19003193ff8SBarry Smith PetscInt *displayfields; 19103193ff8SBarry Smith PetscInt k,ndisplayfields; 1926083293cSBarry Smith PetscBool flg,hold; 19347c6ae99SBarry Smith 19447c6ae99SBarry Smith PetscFunctionBegin; 19547c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr); 19647c6ae99SBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 197a56202b9SBarry Smith ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr); 19847c6ae99SBarry Smith 199*c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 200aa219208SBarry Smith if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 20147c6ae99SBarry Smith 202671f6225SBarry Smith ierr = PetscOptionsGetBool(PETSC_NULL,"-draw_vec_mark_points",&showpoints,PETSC_NULL);CHKERRQ(ierr); 20347c6ae99SBarry Smith 2041321219cSEthan Coon ierr = DMDAGetInfo(da,0,&N,0,0,0,0,0,&step,0,&bx,0,0,0);CHKERRQ(ierr); 205aa219208SBarry Smith ierr = DMDAGetCorners(da,&istart,0,0,&isize,0,0);CHKERRQ(ierr); 20647c6ae99SBarry Smith ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 20747c6ae99SBarry Smith ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr); 20847c6ae99SBarry Smith n = n/step; 20947c6ae99SBarry Smith 21047c6ae99SBarry Smith /* get coordinates of nodes */ 2116636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 21247c6ae99SBarry Smith if (!xcoor) { 213aa219208SBarry Smith ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr); 2146636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 21547c6ae99SBarry Smith } 21647c6ae99SBarry Smith ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr); 21747c6ae99SBarry Smith 21847c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 21947c6ae99SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 22047c6ae99SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 22147c6ae99SBarry Smith 22247c6ae99SBarry Smith /* 22347c6ae99SBarry Smith Determine the min and max x coordinate in plot 22447c6ae99SBarry Smith */ 22547c6ae99SBarry Smith if (!rank) { 22647c6ae99SBarry Smith xmin = PetscRealPart(xg[0]); 22747c6ae99SBarry Smith } 22847c6ae99SBarry Smith if (rank == size-1) { 22947c6ae99SBarry Smith xmax = PetscRealPart(xg[n-1]); 23047c6ae99SBarry Smith } 23147c6ae99SBarry Smith ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 23247c6ae99SBarry Smith ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr); 23347c6ae99SBarry Smith 23403193ff8SBarry Smith ierr = PetscMalloc(step*sizeof(PetscInt),&displayfields);CHKERRQ(ierr); 23503193ff8SBarry Smith for (i=0; i<step; i++) displayfields[i] = i; 23603193ff8SBarry Smith ndisplayfields = step; 23703193ff8SBarry Smith ierr = PetscOptionsGetIntArray(PETSC_NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr); 23803193ff8SBarry Smith if (!flg) ndisplayfields = step; 23903193ff8SBarry Smith for (k=0; k<ndisplayfields; k++) { 24003193ff8SBarry Smith j = displayfields[k]; 24103193ff8SBarry Smith ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr); 24247c6ae99SBarry Smith ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 24347c6ae99SBarry Smith 24447c6ae99SBarry Smith /* 24547c6ae99SBarry Smith Determine the min and max y coordinate in plot 24647c6ae99SBarry Smith */ 24747c6ae99SBarry Smith min = 1.e20; max = -1.e20; 24847c6ae99SBarry Smith for (i=0; i<n; i++) { 24947c6ae99SBarry Smith if (PetscRealPart(array[j+i*step]) < min) min = PetscRealPart(array[j+i*step]); 25047c6ae99SBarry Smith if (PetscRealPart(array[j+i*step]) > max) max = PetscRealPart(array[j+i*step]); 25147c6ae99SBarry Smith } 25247c6ae99SBarry Smith if (min + 1.e-10 > max) { 25347c6ae99SBarry Smith min -= 1.e-5; 25447c6ae99SBarry Smith max += 1.e-5; 25547c6ae99SBarry Smith } 256a56202b9SBarry Smith if (j < nbounds) { 257a56202b9SBarry Smith min = PetscMin(min,bounds[2*j]); 258a56202b9SBarry Smith max = PetscMax(max,bounds[2*j+1]); 259a56202b9SBarry Smith } 260a56202b9SBarry Smith 261d9822059SBarry Smith ierr = MPI_Reduce(&min,&ymin,1,MPIU_REAL,MPIU_MIN,0,comm);CHKERRQ(ierr); 262d9822059SBarry Smith ierr = MPI_Reduce(&max,&ymax,1,MPIU_REAL,MPIU_MAX,0,comm);CHKERRQ(ierr); 26347c6ae99SBarry Smith 2646083293cSBarry Smith ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr); 2656083293cSBarry Smith if (!hold) { 26647c6ae99SBarry Smith ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 2676083293cSBarry Smith } 26803193ff8SBarry Smith ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr); 26947c6ae99SBarry Smith ierr = PetscLogObjectParent(draw,axis);CHKERRQ(ierr); 27047c6ae99SBarry Smith if (!rank) { 27147c6ae99SBarry Smith const char *title; 27247c6ae99SBarry Smith 27347c6ae99SBarry Smith ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,ymin,ymax);CHKERRQ(ierr); 27447c6ae99SBarry Smith ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr); 27547c6ae99SBarry Smith ierr = PetscDrawGetCoordinates(draw,coors,coors+1,coors+2,coors+3);CHKERRQ(ierr); 276aa219208SBarry Smith ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr); 27747c6ae99SBarry Smith if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);} 27847c6ae99SBarry Smith } 27947c6ae99SBarry Smith ierr = MPI_Bcast(coors,4,MPIU_REAL,0,comm);CHKERRQ(ierr); 28047c6ae99SBarry Smith if (rank) { 28147c6ae99SBarry Smith ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 28247c6ae99SBarry Smith } 28347c6ae99SBarry Smith 28447c6ae99SBarry Smith /* draw local part of vector */ 2857d1ecd34SBarry Smith ierr = PetscObjectGetNewTag((PetscObject)xin,&tag1);CHKERRQ(ierr); 2867d1ecd34SBarry Smith ierr = PetscObjectGetNewTag((PetscObject)xin,&tag2);CHKERRQ(ierr); 28747c6ae99SBarry Smith if (rank < size-1) { /*send value to right */ 28847c6ae99SBarry Smith ierr = MPI_Send((void*)&array[j+(n-1)*step],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr); 28947c6ae99SBarry Smith ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag1,comm);CHKERRQ(ierr); 29047c6ae99SBarry Smith } 2911321219cSEthan Coon if (!rank && bx == DMDA_BOUNDARY_PERIODIC && size > 1) { /* first processor sends first value to last */ 29247c6ae99SBarry Smith ierr = MPI_Send((void*)&array[j],1,MPIU_REAL,size-1,tag2,comm);CHKERRQ(ierr); 29347c6ae99SBarry Smith } 29447c6ae99SBarry Smith 29547c6ae99SBarry Smith for (i=1; i<n; i++) { 29647c6ae99SBarry 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); 29747c6ae99SBarry Smith if (showpoints) { 29847c6ae99SBarry Smith ierr = PetscDrawPoint(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+step*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr); 29947c6ae99SBarry Smith } 30047c6ae99SBarry Smith } 30147c6ae99SBarry Smith if (rank) { /* receive value from left */ 30247c6ae99SBarry Smith ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr); 30347c6ae99SBarry Smith ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag1,comm,&status);CHKERRQ(ierr); 30447c6ae99SBarry Smith ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr); 30547c6ae99SBarry Smith if (showpoints) { 30647c6ae99SBarry Smith ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr); 30747c6ae99SBarry Smith } 30847c6ae99SBarry Smith } 3091321219cSEthan Coon if (rank == size-1 && bx == DMDA_BOUNDARY_PERIODIC && size > 1) { 31047c6ae99SBarry Smith ierr = MPI_Recv(&tmp,1,MPIU_REAL,0,tag2,comm,&status);CHKERRQ(ierr); 31117eed501SBarry 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 */ 31217eed501SBarry 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); 31347c6ae99SBarry Smith if (showpoints) { 31447c6ae99SBarry Smith ierr = PetscDrawPoint(draw,PetscRealPart(xg[n-2]),PetscRealPart(array[j+step*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr); 31547c6ae99SBarry Smith } 31647c6ae99SBarry Smith } 31747c6ae99SBarry Smith ierr = PetscDrawSynchronizedFlush(draw);CHKERRQ(ierr); 31847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 31947c6ae99SBarry Smith } 320cd723cd1SBarry Smith ierr = PetscFree(displayfields);CHKERRQ(ierr); 32147c6ae99SBarry Smith ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr); 32247c6ae99SBarry Smith ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 32347c6ae99SBarry Smith PetscFunctionReturn(0); 32447c6ae99SBarry Smith } 32547c6ae99SBarry Smith 326