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 11*d083f849SBarry Smith Collective on da 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; 279a42bb27SBarry Smith DM cda; 28bff4a2f0SMatthew G. Knepley DMBoundaryType bx,by,bz; 2947c6ae99SBarry Smith Vec xcoor; 3047c6ae99SBarry Smith PetscScalar *coors; 3147c6ae99SBarry Smith PetscReal hx,hy,hz_; 3247c6ae99SBarry Smith PetscInt i,j,k,M,N,P,istart,isize,jstart,jsize,kstart,ksize,dim,cnt; 3347c6ae99SBarry Smith PetscErrorCode ierr; 3447c6ae99SBarry Smith 3547c6ae99SBarry Smith PetscFunctionBegin; 36a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 37821c30cfSPeter Brune ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,0,0,&bx,&by,&bz,0);CHKERRQ(ierr); 38c593f006SBarry Smith if (xmax < xmin) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_INCOMP,"xmax must be larger than xmin %g %g",(double)xmin,(double)xmax); 39abf1d30eSPatrick 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); 40abf1d30eSPatrick 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); 4147c6ae99SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 42aa219208SBarry Smith ierr = DMDAGetCorners(da,&istart,&jstart,&kstart,&isize,&jsize,&ksize);CHKERRQ(ierr); 436636e97aSMatthew G Knepley ierr = DMGetCoordinateDM(da, &cda);CHKERRQ(ierr); 44564755cdSBarry Smith ierr = DMCreateGlobalVector(cda, &xcoor);CHKERRQ(ierr); 4547c6ae99SBarry Smith if (dim == 1) { 46bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/M; 47ce00eea3SSatish Balay else hx = (xmax-xmin)/(M-1); 4847c6ae99SBarry Smith ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr); 4947c6ae99SBarry Smith for (i=0; i<isize; i++) { 5047c6ae99SBarry Smith coors[i] = xmin + hx*(i+istart); 5147c6ae99SBarry Smith } 5247c6ae99SBarry Smith ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr); 5347c6ae99SBarry Smith } else if (dim == 2) { 54bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M); 5547c6ae99SBarry Smith else hx = (xmax-xmin)/(M-1); 56bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N); 5747c6ae99SBarry Smith else hy = (ymax-ymin)/(N-1); 5847c6ae99SBarry Smith ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr); 5947c6ae99SBarry Smith cnt = 0; 6047c6ae99SBarry Smith for (j=0; j<jsize; j++) { 6147c6ae99SBarry Smith for (i=0; i<isize; i++) { 6247c6ae99SBarry Smith coors[cnt++] = xmin + hx*(i+istart); 6347c6ae99SBarry Smith coors[cnt++] = ymin + hy*(j+jstart); 6447c6ae99SBarry Smith } 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr); 6747c6ae99SBarry Smith } else if (dim == 3) { 68bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax-xmin)/(M); 6947c6ae99SBarry Smith else hx = (xmax-xmin)/(M-1); 70bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) hy = (ymax-ymin)/(N); 7147c6ae99SBarry Smith else hy = (ymax-ymin)/(N-1); 72bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax-zmin)/(P); 7347c6ae99SBarry Smith else hz_ = (zmax-zmin)/(P-1); 7447c6ae99SBarry Smith ierr = VecGetArray(xcoor,&coors);CHKERRQ(ierr); 7547c6ae99SBarry Smith cnt = 0; 7647c6ae99SBarry Smith for (k=0; k<ksize; k++) { 7747c6ae99SBarry Smith for (j=0; j<jsize; j++) { 7847c6ae99SBarry Smith for (i=0; i<isize; i++) { 7947c6ae99SBarry Smith coors[cnt++] = xmin + hx*(i+istart); 8047c6ae99SBarry Smith coors[cnt++] = ymin + hy*(j+jstart); 8147c6ae99SBarry Smith coors[cnt++] = zmin + hz_*(k+kstart); 8247c6ae99SBarry Smith } 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith ierr = VecRestoreArray(xcoor,&coors);CHKERRQ(ierr); 86ce94432eSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot create uniform coordinates for this dimension %D\n",dim); 876636e97aSMatthew G Knepley ierr = DMSetCoordinates(da,xcoor);CHKERRQ(ierr); 883bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)xcoor);CHKERRQ(ierr); 89fcfd50ebSBarry Smith ierr = VecDestroy(&xcoor);CHKERRQ(ierr); 9047c6ae99SBarry Smith PetscFunctionReturn(0); 9147c6ae99SBarry Smith } 9247c6ae99SBarry Smith 939ae14b6eSBarry Smith /* 949ae14b6eSBarry Smith Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA 959ae14b6eSBarry Smith */ 964e6118eeSBarry Smith PetscErrorCode DMDASelectFields(DM da,PetscInt *outfields,PetscInt **fields) 974e6118eeSBarry Smith { 984e6118eeSBarry Smith PetscErrorCode ierr; 994e6118eeSBarry Smith PetscInt step,ndisplayfields,*displayfields,k,j; 1004e6118eeSBarry Smith PetscBool flg; 1014e6118eeSBarry Smith 1024e6118eeSBarry Smith PetscFunctionBegin; 1034e6118eeSBarry Smith ierr = DMDAGetInfo(da,0,0,0,0,0,0,0,&step,0,0,0,0,0);CHKERRQ(ierr); 104785e854fSJed Brown ierr = PetscMalloc1(step,&displayfields);CHKERRQ(ierr); 1054e6118eeSBarry Smith for (k=0; k<step; k++) displayfields[k] = k; 1064e6118eeSBarry Smith ndisplayfields = step; 107c5929fdfSBarry Smith ierr = PetscOptionsGetIntArray(NULL,NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr); 1084e6118eeSBarry Smith if (!ndisplayfields) ndisplayfields = step; 1094e6118eeSBarry Smith if (!flg) { 1104e6118eeSBarry Smith char **fields; 1114e6118eeSBarry Smith const char *fieldname; 1124e6118eeSBarry Smith PetscInt nfields = step; 113785e854fSJed Brown ierr = PetscMalloc1(step,&fields);CHKERRQ(ierr); 114c5929fdfSBarry Smith ierr = PetscOptionsGetStringArray(NULL,NULL,"-draw_fields_by_name",fields,&nfields,&flg);CHKERRQ(ierr); 1154e6118eeSBarry Smith if (flg) { 1164e6118eeSBarry Smith ndisplayfields = 0; 1174e6118eeSBarry Smith for (k=0; k<nfields;k++) { 1184e6118eeSBarry Smith for (j=0; j<step; j++) { 1194e6118eeSBarry Smith ierr = DMDAGetFieldName(da,j,&fieldname);CHKERRQ(ierr); 1204e6118eeSBarry Smith ierr = PetscStrcmp(fieldname,fields[k],&flg);CHKERRQ(ierr); 1214e6118eeSBarry Smith if (flg) { 1224e6118eeSBarry Smith goto found; 1234e6118eeSBarry Smith } 1244e6118eeSBarry Smith } 125ce94432eSBarry Smith SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Unknown fieldname %s",fields[k]); 1264e6118eeSBarry Smith found: displayfields[ndisplayfields++] = j; 1274e6118eeSBarry Smith } 1284e6118eeSBarry Smith } 129e3f3e4b6SBarry Smith for (k=0; k<nfields; k++) { 130e3f3e4b6SBarry Smith ierr = PetscFree(fields[k]);CHKERRQ(ierr); 131e3f3e4b6SBarry Smith } 1324e6118eeSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 1334e6118eeSBarry Smith } 1344e6118eeSBarry Smith *fields = displayfields; 1354e6118eeSBarry Smith *outfields = ndisplayfields; 1364e6118eeSBarry Smith PetscFunctionReturn(0); 1374e6118eeSBarry Smith } 1384e6118eeSBarry Smith 1399804daf3SBarry Smith #include <petscdraw.h> 14073f7a4c5SBarry Smith 14147c6ae99SBarry Smith PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v) 14247c6ae99SBarry Smith { 1439a42bb27SBarry Smith DM da; 14447c6ae99SBarry Smith PetscErrorCode ierr; 145e5ab1681SLisandro Dalcin PetscMPIInt rank,size,tag; 146e5ab1681SLisandro Dalcin PetscInt i,n,N,dof,istart,isize,j,nbounds; 14747c6ae99SBarry Smith MPI_Status status; 148e5ab1681SLisandro Dalcin PetscReal min,max,xmin = 0.0,xmax = 0.0,tmp = 0.0,xgtmp = 0.0; 14947c6ae99SBarry Smith const PetscScalar *array,*xg; 15047c6ae99SBarry Smith PetscDraw draw; 151e5ab1681SLisandro Dalcin PetscBool isnull,useports = PETSC_FALSE,showmarkers = PETSC_FALSE; 15247c6ae99SBarry Smith MPI_Comm comm; 15347c6ae99SBarry Smith PetscDrawAxis axis; 15447c6ae99SBarry Smith Vec xcoor; 155bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 156e5ab1681SLisandro Dalcin const char *tlabel = NULL,*xlabel = NULL; 157a56202b9SBarry Smith const PetscReal *bounds; 15803193ff8SBarry Smith PetscInt *displayfields; 15903193ff8SBarry Smith PetscInt k,ndisplayfields; 1604e6118eeSBarry Smith PetscBool hold; 161e5ab1681SLisandro Dalcin PetscDrawViewPorts *ports = NULL; 162e5ab1681SLisandro Dalcin PetscViewerFormat format; 16347c6ae99SBarry Smith 16447c6ae99SBarry Smith PetscFunctionBegin; 16547c6ae99SBarry Smith ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr); 16645f3bb6eSLisandro Dalcin ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 16745f3bb6eSLisandro Dalcin if (isnull) PetscFunctionReturn(0); 168a56202b9SBarry Smith ierr = PetscViewerDrawGetBounds(v,&nbounds,&bounds);CHKERRQ(ierr); 16947c6ae99SBarry Smith 170c688c046SMatthew G Knepley ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 171ce94432eSBarry Smith if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 172e5ab1681SLisandro Dalcin ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 173e5ab1681SLisandro Dalcin ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 174e5ab1681SLisandro Dalcin ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 17547c6ae99SBarry Smith 176c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-draw_vec_use_markers",&showmarkers,NULL);CHKERRQ(ierr); 17747c6ae99SBarry Smith 178e5ab1681SLisandro Dalcin ierr = DMDAGetInfo(da,NULL,&N,NULL,NULL,NULL,NULL,NULL,&dof,NULL,&bx,NULL,NULL,NULL);CHKERRQ(ierr); 179e5ab1681SLisandro Dalcin ierr = DMDAGetCorners(da,&istart,NULL,NULL,&isize,NULL,NULL);CHKERRQ(ierr); 18047c6ae99SBarry Smith ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 18147c6ae99SBarry Smith ierr = VecGetLocalSize(xin,&n);CHKERRQ(ierr); 182e5ab1681SLisandro Dalcin n = n/dof; 18347c6ae99SBarry Smith 184832b7cebSLisandro Dalcin /* Get coordinates of nodes */ 1856636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 18647c6ae99SBarry Smith if (!xcoor) { 187aa219208SBarry Smith ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,0.0,0.0,0.0);CHKERRQ(ierr); 1886636e97aSMatthew G Knepley ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 18947c6ae99SBarry Smith } 19047c6ae99SBarry Smith ierr = VecGetArrayRead(xcoor,&xg);CHKERRQ(ierr); 191832b7cebSLisandro Dalcin ierr = DMDAGetCoordinateName(da,0,&xlabel);CHKERRQ(ierr); 192832b7cebSLisandro Dalcin 193832b7cebSLisandro Dalcin /* Determine the min and max coordinate in plot */ 194e5ab1681SLisandro Dalcin if (!rank) xmin = PetscRealPart(xg[0]); 195e5ab1681SLisandro Dalcin if (rank == size-1) xmax = PetscRealPart(xg[n-1]); 19647c6ae99SBarry Smith ierr = MPI_Bcast(&xmin,1,MPIU_REAL,0,comm);CHKERRQ(ierr); 19747c6ae99SBarry Smith ierr = MPI_Bcast(&xmax,1,MPIU_REAL,size-1,comm);CHKERRQ(ierr); 19847c6ae99SBarry Smith 1994e6118eeSBarry Smith ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr); 200e5ab1681SLisandro Dalcin ierr = PetscViewerGetFormat(v,&format);CHKERRQ(ierr); 201e5ab1681SLisandro Dalcin ierr = PetscOptionsGetBool(NULL,NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr); 202832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 203832b7cebSLisandro Dalcin if (useports) { 204e5ab1681SLisandro Dalcin ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr); 205e5ab1681SLisandro Dalcin ierr = PetscViewerDrawGetDrawAxis(v,0,&axis);CHKERRQ(ierr); 206e5ab1681SLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 207e5ab1681SLisandro Dalcin ierr = PetscDrawClear(draw);CHKERRQ(ierr); 208e5ab1681SLisandro Dalcin ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr); 20973f7a4c5SBarry Smith } 210e5ab1681SLisandro Dalcin 211832b7cebSLisandro Dalcin /* Loop over each field; drawing each in a different window */ 21203193ff8SBarry Smith for (k=0; k<ndisplayfields; k++) { 21303193ff8SBarry Smith j = displayfields[k]; 21447c6ae99SBarry Smith 215e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 216e5ab1681SLisandro Dalcin ierr = VecStrideMin(xin,j,NULL,&min);CHKERRQ(ierr); 217e5ab1681SLisandro Dalcin ierr = VecStrideMax(xin,j,NULL,&max);CHKERRQ(ierr); 218a56202b9SBarry Smith if (j < nbounds) { 219a56202b9SBarry Smith min = PetscMin(min,bounds[2*j]); 220a56202b9SBarry Smith max = PetscMax(max,bounds[2*j+1]); 221a56202b9SBarry Smith } 222e5ab1681SLisandro Dalcin if (min == max) { 223e5ab1681SLisandro Dalcin min -= 1.e-5; 224e5ab1681SLisandro Dalcin max += 1.e-5; 225e5ab1681SLisandro Dalcin } 22647c6ae99SBarry Smith 227832b7cebSLisandro Dalcin if (useports) { 228832b7cebSLisandro Dalcin ierr = PetscDrawViewPortsSet(ports,k);CHKERRQ(ierr); 229832b7cebSLisandro Dalcin ierr = DMDAGetFieldName(da,j,&tlabel);CHKERRQ(ierr); 230832b7cebSLisandro Dalcin } else { 231832b7cebSLisandro Dalcin const char *title; 2326083293cSBarry Smith ierr = PetscViewerDrawGetHold(v,&hold);CHKERRQ(ierr); 233832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDraw(v,k,&draw);CHKERRQ(ierr); 234832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDrawAxis(v,k,&axis);CHKERRQ(ierr); 235832b7cebSLisandro Dalcin ierr = DMDAGetFieldName(da,j,&title);CHKERRQ(ierr); 236832b7cebSLisandro Dalcin if (title) {ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr);} 237832b7cebSLisandro Dalcin ierr = PetscDrawCheckResizedWindow(draw);CHKERRQ(ierr); 238e5ab1681SLisandro Dalcin if (!hold) {ierr = PetscDrawClear(draw);CHKERRQ(ierr);} 239832b7cebSLisandro Dalcin } 240832b7cebSLisandro Dalcin ierr = PetscDrawAxisSetLabels(axis,tlabel,xlabel,NULL);CHKERRQ(ierr); 241e5ab1681SLisandro Dalcin ierr = PetscDrawAxisSetLimits(axis,xmin,xmax,min,max);CHKERRQ(ierr); 24247c6ae99SBarry Smith ierr = PetscDrawAxisDraw(axis);CHKERRQ(ierr); 24347c6ae99SBarry Smith 24447c6ae99SBarry Smith /* draw local part of vector */ 245e5ab1681SLisandro Dalcin ierr = PetscObjectGetNewTag((PetscObject)xin,&tag);CHKERRQ(ierr); 24647c6ae99SBarry Smith if (rank < size-1) { /*send value to right */ 247e5ab1681SLisandro Dalcin ierr = MPI_Send((void*)&xg[n-1],1,MPIU_REAL,rank+1,tag,comm);CHKERRQ(ierr); 248e5ab1681SLisandro Dalcin ierr = MPI_Send((void*)&array[j+(n-1)*dof],1,MPIU_REAL,rank+1,tag,comm);CHKERRQ(ierr); 24947c6ae99SBarry Smith } 25047c6ae99SBarry Smith if (rank) { /* receive value from left */ 251e5ab1681SLisandro Dalcin ierr = MPI_Recv(&xgtmp,1,MPIU_REAL,rank-1,tag,comm,&status);CHKERRQ(ierr); 252e5ab1681SLisandro Dalcin ierr = MPI_Recv(&tmp,1,MPIU_REAL,rank-1,tag,comm,&status);CHKERRQ(ierr); 253e5ab1681SLisandro Dalcin } 254e5ab1681SLisandro Dalcin ierr = PetscDrawCollectiveBegin(draw);CHKERRQ(ierr); 255e5ab1681SLisandro Dalcin if (rank) { 25647c6ae99SBarry Smith ierr = PetscDrawLine(draw,xgtmp,tmp,PetscRealPart(xg[0]),PetscRealPart(array[j]),PETSC_DRAW_RED);CHKERRQ(ierr); 257e5ab1681SLisandro Dalcin if (showmarkers) {ierr = PetscDrawPoint(draw,xgtmp,tmp,PETSC_DRAW_BLACK);CHKERRQ(ierr);} 25847c6ae99SBarry Smith } 259e5ab1681SLisandro Dalcin for (i=1; i<n; i++) { 260e5ab1681SLisandro 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); 261e5ab1681SLisandro Dalcin if (showmarkers) {ierr = PetscDrawMarker(draw,PetscRealPart(xg[i-1]),PetscRealPart(array[j+dof*(i-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);} 26247c6ae99SBarry Smith } 263e5ab1681SLisandro Dalcin if (rank == size-1) { 264e5ab1681SLisandro Dalcin if (showmarkers) {ierr = PetscDrawMarker(draw,PetscRealPart(xg[n-1]),PetscRealPart(array[j+dof*(n-1)]),PETSC_DRAW_BLACK);CHKERRQ(ierr);} 26547c6ae99SBarry Smith } 266e5ab1681SLisandro Dalcin ierr = PetscDrawCollectiveEnd(draw);CHKERRQ(ierr); 2675b399a63SLisandro Dalcin ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 26847c6ae99SBarry Smith ierr = PetscDrawPause(draw);CHKERRQ(ierr); 269832b7cebSLisandro Dalcin if (!useports) {ierr = PetscDrawSave(draw);CHKERRQ(ierr);} 27047c6ae99SBarry Smith } 271832b7cebSLisandro Dalcin if (useports) { 272832b7cebSLisandro Dalcin ierr = PetscViewerDrawGetDraw(v,0,&draw);CHKERRQ(ierr); 273832b7cebSLisandro Dalcin ierr = PetscDrawSave(draw);CHKERRQ(ierr); 274832b7cebSLisandro Dalcin } 275832b7cebSLisandro Dalcin 276e5ab1681SLisandro Dalcin ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr); 277cd723cd1SBarry Smith ierr = PetscFree(displayfields);CHKERRQ(ierr); 27847c6ae99SBarry Smith ierr = VecRestoreArrayRead(xcoor,&xg);CHKERRQ(ierr); 27947c6ae99SBarry Smith ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 28047c6ae99SBarry Smith PetscFunctionReturn(0); 28147c6ae99SBarry Smith } 28247c6ae99SBarry Smith 283