xref: /petsc/src/dm/impls/da/gr1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
11d083f849SBarry 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 
21db781477SPatrick Sanan .seealso: `DMSetCoordinates()`, `DMGetCoordinates()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMStagSetUniformCoordinates()`
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith @*/
249371c9d4SSatish Balay PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax) {
2547c6ae99SBarry Smith   MPI_Comm       comm;
269a42bb27SBarry Smith   DM             cda;
27a5bc1bf3SBarry Smith   DM_DA         *dd = (DM_DA *)da->data;
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 
3447c6ae99SBarry Smith   PetscFunctionBegin;
35a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
367a8be351SBarry Smith   PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup");
379566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
3808401ef6SPierre Jolivet   PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax);
3908401ef6SPierre Jolivet   PetscCheck(!(dim > 1) || !(ymax < ymin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "ymax must be larger than ymin %g %g", (double)ymin, (double)ymax);
4008401ef6SPierre Jolivet   PetscCheck(!(dim > 2) || !(zmax < zmin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "zmax must be larger than zmin %g %g", (double)zmin, (double)zmax);
419566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
429566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize));
439566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
449566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(cda, &xcoor));
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);
489566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xcoor, &coors));
499371c9d4SSatish Balay     for (i = 0; i < isize; i++) { coors[i] = xmin + hx * (i + istart); }
509566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xcoor, &coors));
5147c6ae99SBarry Smith   } else if (dim == 2) {
52bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
5347c6ae99SBarry Smith     else hx = (xmax - xmin) / (M - 1);
54bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
5547c6ae99SBarry Smith     else hy = (ymax - ymin) / (N - 1);
569566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xcoor, &coors));
5747c6ae99SBarry Smith     cnt = 0;
5847c6ae99SBarry Smith     for (j = 0; j < jsize; j++) {
5947c6ae99SBarry Smith       for (i = 0; i < isize; i++) {
6047c6ae99SBarry Smith         coors[cnt++] = xmin + hx * (i + istart);
6147c6ae99SBarry Smith         coors[cnt++] = ymin + hy * (j + jstart);
6247c6ae99SBarry Smith       }
6347c6ae99SBarry Smith     }
649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xcoor, &coors));
6547c6ae99SBarry Smith   } else if (dim == 3) {
66bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
6747c6ae99SBarry Smith     else hx = (xmax - xmin) / (M - 1);
68bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
6947c6ae99SBarry Smith     else hy = (ymax - ymin) / (N - 1);
70bff4a2f0SMatthew G. Knepley     if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax - zmin) / (P);
7147c6ae99SBarry Smith     else hz_ = (zmax - zmin) / (P - 1);
729566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xcoor, &coors));
7347c6ae99SBarry Smith     cnt = 0;
7447c6ae99SBarry Smith     for (k = 0; k < ksize; k++) {
7547c6ae99SBarry Smith       for (j = 0; j < jsize; j++) {
7647c6ae99SBarry Smith         for (i = 0; i < isize; i++) {
7747c6ae99SBarry Smith           coors[cnt++] = xmin + hx * (i + istart);
7847c6ae99SBarry Smith           coors[cnt++] = ymin + hy * (j + jstart);
7947c6ae99SBarry Smith           coors[cnt++] = zmin + hz_ * (k + kstart);
8047c6ae99SBarry Smith         }
8147c6ae99SBarry Smith       }
8247c6ae99SBarry Smith     }
839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xcoor, &coors));
8463a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim);
859566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(da, xcoor));
869566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)da, (PetscObject)xcoor));
879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xcoor));
8847c6ae99SBarry Smith   PetscFunctionReturn(0);
8947c6ae99SBarry Smith }
9047c6ae99SBarry Smith 
919ae14b6eSBarry Smith /*
929ae14b6eSBarry Smith     Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA
939ae14b6eSBarry Smith */
949371c9d4SSatish Balay PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields) {
954e6118eeSBarry Smith   PetscInt  step, ndisplayfields, *displayfields, k, j;
964e6118eeSBarry Smith   PetscBool flg;
974e6118eeSBarry Smith 
984e6118eeSBarry Smith   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL));
1009566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(step, &displayfields));
1014e6118eeSBarry Smith   for (k = 0; k < step; k++) displayfields[k] = k;
1024e6118eeSBarry Smith   ndisplayfields = step;
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg));
1044e6118eeSBarry Smith   if (!ndisplayfields) ndisplayfields = step;
1054e6118eeSBarry Smith   if (!flg) {
1064e6118eeSBarry Smith     char      **fields;
1074e6118eeSBarry Smith     const char *fieldname;
1084e6118eeSBarry Smith     PetscInt    nfields = step;
1099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(step, &fields));
1109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg));
1114e6118eeSBarry Smith     if (flg) {
1124e6118eeSBarry Smith       ndisplayfields = 0;
1134e6118eeSBarry Smith       for (k = 0; k < nfields; k++) {
1144e6118eeSBarry Smith         for (j = 0; j < step; j++) {
1159566063dSJacob Faibussowitsch           PetscCall(DMDAGetFieldName(da, j, &fieldname));
1169566063dSJacob Faibussowitsch           PetscCall(PetscStrcmp(fieldname, fields[k], &flg));
1179371c9d4SSatish Balay           if (flg) { goto found; }
1184e6118eeSBarry Smith         }
11998921bdaSJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]);
1209371c9d4SSatish Balay       found:
1219371c9d4SSatish Balay         displayfields[ndisplayfields++] = j;
1224e6118eeSBarry Smith       }
1234e6118eeSBarry Smith     }
124*48a46eb9SPierre Jolivet     for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k]));
1259566063dSJacob Faibussowitsch     PetscCall(PetscFree(fields));
1264e6118eeSBarry Smith   }
1274e6118eeSBarry Smith   *fields    = displayfields;
1284e6118eeSBarry Smith   *outfields = ndisplayfields;
1294e6118eeSBarry Smith   PetscFunctionReturn(0);
1304e6118eeSBarry Smith }
1314e6118eeSBarry Smith 
1329804daf3SBarry Smith #include <petscdraw.h>
13373f7a4c5SBarry Smith 
1349371c9d4SSatish Balay PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v) {
1359a42bb27SBarry Smith   DM                  da;
136e5ab1681SLisandro Dalcin   PetscMPIInt         rank, size, tag;
137e5ab1681SLisandro Dalcin   PetscInt            i, n, N, dof, istart, isize, j, nbounds;
13847c6ae99SBarry Smith   MPI_Status          status;
139e5ab1681SLisandro Dalcin   PetscReal           min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0;
14047c6ae99SBarry Smith   const PetscScalar  *array, *xg;
14147c6ae99SBarry Smith   PetscDraw           draw;
142e5ab1681SLisandro Dalcin   PetscBool           isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE;
14347c6ae99SBarry Smith   MPI_Comm            comm;
14447c6ae99SBarry Smith   PetscDrawAxis       axis;
14547c6ae99SBarry Smith   Vec                 xcoor;
146bff4a2f0SMatthew G. Knepley   DMBoundaryType      bx;
147e5ab1681SLisandro Dalcin   const char         *tlabel = NULL, *xlabel = NULL;
148a56202b9SBarry Smith   const PetscReal    *bounds;
14903193ff8SBarry Smith   PetscInt           *displayfields;
15003193ff8SBarry Smith   PetscInt            k, ndisplayfields;
1514e6118eeSBarry Smith   PetscBool           hold;
152e5ab1681SLisandro Dalcin   PetscDrawViewPorts *ports = NULL;
153e5ab1681SLisandro Dalcin   PetscViewerFormat   format;
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
1579566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
15845f3bb6eSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1599566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds));
16047c6ae99SBarry Smith 
1619566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
1627a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
1639566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)xin, &comm));
1649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
16647c6ae99SBarry Smith 
1679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL));
16847c6ae99SBarry Smith 
1699566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL));
1709566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL));
1719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin, &array));
1729566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xin, &n));
173e5ab1681SLisandro Dalcin   n = n / dof;
17447c6ae99SBarry Smith 
175832b7cebSLisandro Dalcin   /* Get coordinates of nodes */
1769566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &xcoor));
17747c6ae99SBarry Smith   if (!xcoor) {
1789566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0));
1799566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da, &xcoor));
18047c6ae99SBarry Smith   }
1819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xcoor, &xg));
1829566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da, 0, &xlabel));
183832b7cebSLisandro Dalcin 
184832b7cebSLisandro Dalcin   /* Determine the min and max coordinate in plot */
185dd400576SPatrick Sanan   if (rank == 0) xmin = PetscRealPart(xg[0]);
186e5ab1681SLisandro Dalcin   if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]);
1879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm));
1889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm));
18947c6ae99SBarry Smith 
1909566063dSJacob Faibussowitsch   PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields));
1919566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(v, &format));
1929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL));
193832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
194832b7cebSLisandro Dalcin   if (useports) {
1959566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
1969566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis));
1979566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
1989566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
1999566063dSJacob Faibussowitsch     PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports));
20073f7a4c5SBarry Smith   }
201e5ab1681SLisandro Dalcin 
202832b7cebSLisandro Dalcin   /* Loop over each field; drawing each in a different window */
20303193ff8SBarry Smith   for (k = 0; k < ndisplayfields; k++) {
20403193ff8SBarry Smith     j = displayfields[k];
20547c6ae99SBarry Smith 
206e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2079566063dSJacob Faibussowitsch     PetscCall(VecStrideMin(xin, j, NULL, &min));
2089566063dSJacob Faibussowitsch     PetscCall(VecStrideMax(xin, j, NULL, &max));
209a56202b9SBarry Smith     if (j < nbounds) {
210a56202b9SBarry Smith       min = PetscMin(min, bounds[2 * j]);
211a56202b9SBarry Smith       max = PetscMax(max, bounds[2 * j + 1]);
212a56202b9SBarry Smith     }
213e5ab1681SLisandro Dalcin     if (min == max) {
214e5ab1681SLisandro Dalcin       min -= 1.e-5;
215e5ab1681SLisandro Dalcin       max += 1.e-5;
216e5ab1681SLisandro Dalcin     }
21747c6ae99SBarry Smith 
218832b7cebSLisandro Dalcin     if (useports) {
2199566063dSJacob Faibussowitsch       PetscCall(PetscDrawViewPortsSet(ports, k));
2209566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da, j, &tlabel));
221832b7cebSLisandro Dalcin     } else {
222832b7cebSLisandro Dalcin       const char *title;
2239566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetHold(v, &hold));
2249566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(v, k, &draw));
2259566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis));
2269566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da, j, &title));
2279566063dSJacob Faibussowitsch       if (title) PetscCall(PetscDrawSetTitle(draw, title));
2289566063dSJacob Faibussowitsch       PetscCall(PetscDrawCheckResizedWindow(draw));
2299566063dSJacob Faibussowitsch       if (!hold) PetscCall(PetscDrawClear(draw));
230832b7cebSLisandro Dalcin     }
2319566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL));
2329566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max));
2339566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
23447c6ae99SBarry Smith 
23547c6ae99SBarry Smith     /* draw local part of vector */
2369566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag));
23747c6ae99SBarry Smith     if (rank < size - 1) { /*send value to right */
2389566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm));
2399566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm));
24047c6ae99SBarry Smith     }
24147c6ae99SBarry Smith     if (rank) { /* receive value from left */
2429566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
2439566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
244e5ab1681SLisandro Dalcin     }
245d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
246e5ab1681SLisandro Dalcin     if (rank) {
2479566063dSJacob Faibussowitsch       PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED));
2489566063dSJacob Faibussowitsch       if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK));
24947c6ae99SBarry Smith     }
250e5ab1681SLisandro Dalcin     for (i = 1; i < n; i++) {
2519566063dSJacob Faibussowitsch       PetscCall(PetscDrawLine(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PetscRealPart(xg[i]), PetscRealPart(array[j + dof * i]), PETSC_DRAW_RED));
2529566063dSJacob Faibussowitsch       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK));
25347c6ae99SBarry Smith     }
254e5ab1681SLisandro Dalcin     if (rank == size - 1) {
2559566063dSJacob Faibussowitsch       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK));
25647c6ae99SBarry Smith     }
257d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
2589566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
2599566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
2609566063dSJacob Faibussowitsch     if (!useports) PetscCall(PetscDrawSave(draw));
26147c6ae99SBarry Smith   }
262832b7cebSLisandro Dalcin   if (useports) {
2639566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
2649566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
265832b7cebSLisandro Dalcin   }
266832b7cebSLisandro Dalcin 
2679566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsDestroy(ports));
2689566063dSJacob Faibussowitsch   PetscCall(PetscFree(displayfields));
2699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xcoor, &xg));
2709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin, &array));
27147c6ae99SBarry Smith   PetscFunctionReturn(0);
27247c6ae99SBarry Smith }
273