xref: /petsc/src/dm/impls/da/gr1.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 @*/
24*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax)
25*d71ae5a4SJacob Faibussowitsch {
2647c6ae99SBarry Smith   MPI_Comm       comm;
279a42bb27SBarry Smith   DM             cda;
28a5bc1bf3SBarry Smith   DM_DA         *dd = (DM_DA *)da->data;
29bff4a2f0SMatthew G. Knepley   DMBoundaryType bx, by, bz;
3047c6ae99SBarry Smith   Vec            xcoor;
3147c6ae99SBarry Smith   PetscScalar   *coors;
3247c6ae99SBarry Smith   PetscReal      hx, hy, hz_;
3347c6ae99SBarry Smith   PetscInt       i, j, k, M, N, P, istart, isize, jstart, jsize, kstart, ksize, dim, cnt;
3447c6ae99SBarry Smith 
3547c6ae99SBarry Smith   PetscFunctionBegin;
36a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
377a8be351SBarry Smith   PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup");
389566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
3908401ef6SPierre Jolivet   PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax);
4008401ef6SPierre 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);
4108401ef6SPierre 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);
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
439566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize));
449566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
459566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(cda, &xcoor));
4647c6ae99SBarry Smith   if (dim == 1) {
47bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / M;
48ce00eea3SSatish Balay     else hx = (xmax - xmin) / (M - 1);
499566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xcoor, &coors));
50ad540459SPierre Jolivet     for (i = 0; i < isize; i++) coors[i] = xmin + hx * (i + istart);
519566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xcoor, &coors));
5247c6ae99SBarry Smith   } else if (dim == 2) {
53bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
5447c6ae99SBarry Smith     else hx = (xmax - xmin) / (M - 1);
55bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
5647c6ae99SBarry Smith     else hy = (ymax - ymin) / (N - 1);
579566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xcoor, &coors));
5847c6ae99SBarry Smith     cnt = 0;
5947c6ae99SBarry Smith     for (j = 0; j < jsize; j++) {
6047c6ae99SBarry Smith       for (i = 0; i < isize; i++) {
6147c6ae99SBarry Smith         coors[cnt++] = xmin + hx * (i + istart);
6247c6ae99SBarry Smith         coors[cnt++] = ymin + hy * (j + jstart);
6347c6ae99SBarry Smith       }
6447c6ae99SBarry Smith     }
659566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xcoor, &coors));
6647c6ae99SBarry Smith   } else if (dim == 3) {
67bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M);
6847c6ae99SBarry Smith     else hx = (xmax - xmin) / (M - 1);
69bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N);
7047c6ae99SBarry Smith     else hy = (ymax - ymin) / (N - 1);
71bff4a2f0SMatthew G. Knepley     if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax - zmin) / (P);
7247c6ae99SBarry Smith     else hz_ = (zmax - zmin) / (P - 1);
739566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xcoor, &coors));
7447c6ae99SBarry Smith     cnt = 0;
7547c6ae99SBarry Smith     for (k = 0; k < ksize; k++) {
7647c6ae99SBarry Smith       for (j = 0; j < jsize; j++) {
7747c6ae99SBarry Smith         for (i = 0; i < isize; i++) {
7847c6ae99SBarry Smith           coors[cnt++] = xmin + hx * (i + istart);
7947c6ae99SBarry Smith           coors[cnt++] = ymin + hy * (j + jstart);
8047c6ae99SBarry Smith           coors[cnt++] = zmin + hz_ * (k + kstart);
8147c6ae99SBarry Smith         }
8247c6ae99SBarry Smith       }
8347c6ae99SBarry Smith     }
849566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xcoor, &coors));
8563a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim);
869566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(da, 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 */
94*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields)
95*d71ae5a4SJacob Faibussowitsch {
964e6118eeSBarry Smith   PetscInt  step, ndisplayfields, *displayfields, k, j;
974e6118eeSBarry Smith   PetscBool flg;
984e6118eeSBarry Smith 
994e6118eeSBarry Smith   PetscFunctionBegin;
1009566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL));
1019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(step, &displayfields));
1024e6118eeSBarry Smith   for (k = 0; k < step; k++) displayfields[k] = k;
1034e6118eeSBarry Smith   ndisplayfields = step;
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg));
1054e6118eeSBarry Smith   if (!ndisplayfields) ndisplayfields = step;
1064e6118eeSBarry Smith   if (!flg) {
1074e6118eeSBarry Smith     char      **fields;
1084e6118eeSBarry Smith     const char *fieldname;
1094e6118eeSBarry Smith     PetscInt    nfields = step;
1109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(step, &fields));
1119566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg));
1124e6118eeSBarry Smith     if (flg) {
1134e6118eeSBarry Smith       ndisplayfields = 0;
1144e6118eeSBarry Smith       for (k = 0; k < nfields; k++) {
1154e6118eeSBarry Smith         for (j = 0; j < step; j++) {
1169566063dSJacob Faibussowitsch           PetscCall(DMDAGetFieldName(da, j, &fieldname));
1179566063dSJacob Faibussowitsch           PetscCall(PetscStrcmp(fieldname, fields[k], &flg));
118ad540459SPierre Jolivet           if (flg) goto found;
1194e6118eeSBarry Smith         }
12098921bdaSJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]);
1219371c9d4SSatish Balay       found:
1229371c9d4SSatish Balay         displayfields[ndisplayfields++] = j;
1234e6118eeSBarry Smith       }
1244e6118eeSBarry Smith     }
12548a46eb9SPierre Jolivet     for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k]));
1269566063dSJacob Faibussowitsch     PetscCall(PetscFree(fields));
1274e6118eeSBarry Smith   }
1284e6118eeSBarry Smith   *fields    = displayfields;
1294e6118eeSBarry Smith   *outfields = ndisplayfields;
1304e6118eeSBarry Smith   PetscFunctionReturn(0);
1314e6118eeSBarry Smith }
1324e6118eeSBarry Smith 
1339804daf3SBarry Smith #include <petscdraw.h>
13473f7a4c5SBarry Smith 
135*d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v)
136*d71ae5a4SJacob Faibussowitsch {
1379a42bb27SBarry Smith   DM                  da;
138e5ab1681SLisandro Dalcin   PetscMPIInt         rank, size, tag;
139e5ab1681SLisandro Dalcin   PetscInt            i, n, N, dof, istart, isize, j, nbounds;
14047c6ae99SBarry Smith   MPI_Status          status;
141e5ab1681SLisandro Dalcin   PetscReal           min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0;
14247c6ae99SBarry Smith   const PetscScalar  *array, *xg;
14347c6ae99SBarry Smith   PetscDraw           draw;
144e5ab1681SLisandro Dalcin   PetscBool           isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE;
14547c6ae99SBarry Smith   MPI_Comm            comm;
14647c6ae99SBarry Smith   PetscDrawAxis       axis;
14747c6ae99SBarry Smith   Vec                 xcoor;
148bff4a2f0SMatthew G. Knepley   DMBoundaryType      bx;
149e5ab1681SLisandro Dalcin   const char         *tlabel = NULL, *xlabel = NULL;
150a56202b9SBarry Smith   const PetscReal    *bounds;
15103193ff8SBarry Smith   PetscInt           *displayfields;
15203193ff8SBarry Smith   PetscInt            k, ndisplayfields;
1534e6118eeSBarry Smith   PetscBool           hold;
154e5ab1681SLisandro Dalcin   PetscDrawViewPorts *ports = NULL;
155e5ab1681SLisandro Dalcin   PetscViewerFormat   format;
15647c6ae99SBarry Smith 
15747c6ae99SBarry Smith   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
1599566063dSJacob Faibussowitsch   PetscCall(PetscDrawIsNull(draw, &isnull));
16045f3bb6eSLisandro Dalcin   if (isnull) PetscFunctionReturn(0);
1619566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds));
16247c6ae99SBarry Smith 
1639566063dSJacob Faibussowitsch   PetscCall(VecGetDM(xin, &da));
1647a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
1659566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)xin, &comm));
1669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
16847c6ae99SBarry Smith 
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL));
17047c6ae99SBarry Smith 
1719566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL));
1729566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL));
1739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xin, &array));
1749566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(xin, &n));
175e5ab1681SLisandro Dalcin   n = n / dof;
17647c6ae99SBarry Smith 
177832b7cebSLisandro Dalcin   /* Get coordinates of nodes */
1789566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &xcoor));
17947c6ae99SBarry Smith   if (!xcoor) {
1809566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0));
1819566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da, &xcoor));
18247c6ae99SBarry Smith   }
1839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xcoor, &xg));
1849566063dSJacob Faibussowitsch   PetscCall(DMDAGetCoordinateName(da, 0, &xlabel));
185832b7cebSLisandro Dalcin 
186832b7cebSLisandro Dalcin   /* Determine the min and max coordinate in plot */
187dd400576SPatrick Sanan   if (rank == 0) xmin = PetscRealPart(xg[0]);
188e5ab1681SLisandro Dalcin   if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]);
1899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm));
1909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm));
19147c6ae99SBarry Smith 
1929566063dSJacob Faibussowitsch   PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields));
1939566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(v, &format));
1949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL));
195832b7cebSLisandro Dalcin   if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
196832b7cebSLisandro Dalcin   if (useports) {
1979566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
1989566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis));
1999566063dSJacob Faibussowitsch     PetscCall(PetscDrawCheckResizedWindow(draw));
2009566063dSJacob Faibussowitsch     PetscCall(PetscDrawClear(draw));
2019566063dSJacob Faibussowitsch     PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports));
20273f7a4c5SBarry Smith   }
203e5ab1681SLisandro Dalcin 
204832b7cebSLisandro Dalcin   /* Loop over each field; drawing each in a different window */
20503193ff8SBarry Smith   for (k = 0; k < ndisplayfields; k++) {
20603193ff8SBarry Smith     j = displayfields[k];
20747c6ae99SBarry Smith 
208e5ab1681SLisandro Dalcin     /* determine the min and max value in plot */
2099566063dSJacob Faibussowitsch     PetscCall(VecStrideMin(xin, j, NULL, &min));
2109566063dSJacob Faibussowitsch     PetscCall(VecStrideMax(xin, j, NULL, &max));
211a56202b9SBarry Smith     if (j < nbounds) {
212a56202b9SBarry Smith       min = PetscMin(min, bounds[2 * j]);
213a56202b9SBarry Smith       max = PetscMax(max, bounds[2 * j + 1]);
214a56202b9SBarry Smith     }
215e5ab1681SLisandro Dalcin     if (min == max) {
216e5ab1681SLisandro Dalcin       min -= 1.e-5;
217e5ab1681SLisandro Dalcin       max += 1.e-5;
218e5ab1681SLisandro Dalcin     }
21947c6ae99SBarry Smith 
220832b7cebSLisandro Dalcin     if (useports) {
2219566063dSJacob Faibussowitsch       PetscCall(PetscDrawViewPortsSet(ports, k));
2229566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da, j, &tlabel));
223832b7cebSLisandro Dalcin     } else {
224832b7cebSLisandro Dalcin       const char *title;
2259566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetHold(v, &hold));
2269566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(v, k, &draw));
2279566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis));
2289566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldName(da, j, &title));
2299566063dSJacob Faibussowitsch       if (title) PetscCall(PetscDrawSetTitle(draw, title));
2309566063dSJacob Faibussowitsch       PetscCall(PetscDrawCheckResizedWindow(draw));
2319566063dSJacob Faibussowitsch       if (!hold) PetscCall(PetscDrawClear(draw));
232832b7cebSLisandro Dalcin     }
2339566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL));
2349566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max));
2359566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisDraw(axis));
23647c6ae99SBarry Smith 
23747c6ae99SBarry Smith     /* draw local part of vector */
2389566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag));
23947c6ae99SBarry Smith     if (rank < size - 1) { /*send value to right */
2409566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm));
2419566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm));
24247c6ae99SBarry Smith     }
24347c6ae99SBarry Smith     if (rank) { /* receive value from left */
2449566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
2459566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
246e5ab1681SLisandro Dalcin     }
247d0609cedSBarry Smith     PetscDrawCollectiveBegin(draw);
248e5ab1681SLisandro Dalcin     if (rank) {
2499566063dSJacob Faibussowitsch       PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED));
2509566063dSJacob Faibussowitsch       if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK));
25147c6ae99SBarry Smith     }
252e5ab1681SLisandro Dalcin     for (i = 1; i < n; i++) {
2539566063dSJacob 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));
2549566063dSJacob Faibussowitsch       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK));
25547c6ae99SBarry Smith     }
256e5ab1681SLisandro Dalcin     if (rank == size - 1) {
2579566063dSJacob Faibussowitsch       if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK));
25847c6ae99SBarry Smith     }
259d0609cedSBarry Smith     PetscDrawCollectiveEnd(draw);
2609566063dSJacob Faibussowitsch     PetscCall(PetscDrawFlush(draw));
2619566063dSJacob Faibussowitsch     PetscCall(PetscDrawPause(draw));
2629566063dSJacob Faibussowitsch     if (!useports) PetscCall(PetscDrawSave(draw));
26347c6ae99SBarry Smith   }
264832b7cebSLisandro Dalcin   if (useports) {
2659566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
2669566063dSJacob Faibussowitsch     PetscCall(PetscDrawSave(draw));
267832b7cebSLisandro Dalcin   }
268832b7cebSLisandro Dalcin 
2699566063dSJacob Faibussowitsch   PetscCall(PetscDrawViewPortsDestroy(ports));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFree(displayfields));
2719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xcoor, &xg));
2729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xin, &array));
27347c6ae99SBarry Smith   PetscFunctionReturn(0);
27447c6ae99SBarry Smith }
275