147c6ae99SBarry Smith /*
2aa219208SBarry Smith Plots vectors obtained with DMDACreate1d()
347c6ae99SBarry Smith */
447c6ae99SBarry Smith
5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
647c6ae99SBarry Smith
747c6ae99SBarry Smith /*@
8dce8aebaSBarry Smith DMDASetUniformCoordinates - Sets a `DMDA` coordinates to be a uniform grid
947c6ae99SBarry Smith
1020f4b53cSBarry Smith Collective
1147c6ae99SBarry Smith
1247c6ae99SBarry Smith Input Parameters:
13*72fd0fbdSBarry Smith + da - the `DMDA` object
14a4e35b19SJacob Faibussowitsch . xmin - min extreme in the x direction
15a4e35b19SJacob Faibussowitsch . xmax - max extreme in the x direction
16a4e35b19SJacob Faibussowitsch . ymin - min extreme in the y direction (value ignored for 1 dimensional problems)
17a4e35b19SJacob Faibussowitsch . ymax - max extreme in the y direction (value ignored for 1 dimensional problems)
18a4e35b19SJacob Faibussowitsch . zmin - min extreme in the z direction (value ignored for 1 or 2 dimensional problems)
19a4e35b19SJacob Faibussowitsch - zmax - max extreme in the z direction (value ignored for 1 or 2 dimensional problems)
2047c6ae99SBarry Smith
2147c6ae99SBarry Smith Level: beginner
2247c6ae99SBarry Smith
2312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMStagSetUniformCoordinates()`
2447c6ae99SBarry Smith @*/
DMDASetUniformCoordinates(DM da,PetscReal xmin,PetscReal xmax,PetscReal ymin,PetscReal ymax,PetscReal zmin,PetscReal zmax)25d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax)
26d71ae5a4SJacob Faibussowitsch {
2747c6ae99SBarry Smith MPI_Comm comm;
289a42bb27SBarry Smith DM cda;
29a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA *)da->data;
30bff4a2f0SMatthew G. Knepley DMBoundaryType bx, by, bz;
3147c6ae99SBarry Smith Vec xcoor;
3247c6ae99SBarry Smith PetscScalar *coors;
336c6a6b79SMatthew G. Knepley PetscReal hx = 0., hy = 0., hz_ = 0.;
3447c6ae99SBarry Smith PetscInt i, j, k, M, N, P, istart, isize, jstart, jsize, kstart, ksize, dim, cnt;
3547c6ae99SBarry Smith
3647c6ae99SBarry Smith PetscFunctionBegin;
37a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
387a8be351SBarry Smith PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup");
399566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
4008401ef6SPierre Jolivet PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax);
4108401ef6SPierre 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);
4208401ef6SPierre 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);
439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize));
459566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda));
469566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(cda, &xcoor));
4747c6ae99SBarry Smith if (dim == 1) {
48bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / M;
49ce00eea3SSatish Balay else hx = (xmax - xmin) / (M - 1);
509566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor, &coors));
51ad540459SPierre Jolivet for (i = 0; i < isize; i++) coors[i] = xmin + hx * (i + istart);
529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor, &coors));
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);
589566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor, &coors));
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 }
669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor, &coors));
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);
749566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor, &coors));
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 }
859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor, &coors));
8663a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim);
879566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, xcoor));
889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xcoor));
899dc49dd9SMatthew G. Knepley // Handle periodicity
909dc49dd9SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC || by == DM_BOUNDARY_PERIODIC || bz == DM_BOUNDARY_PERIODIC) {
919dc49dd9SMatthew G. Knepley PetscReal maxCell[3] = {hx, hy, hz_};
929dc49dd9SMatthew G. Knepley PetscReal Lstart[3] = {xmin, ymin, zmin};
939dc49dd9SMatthew G. Knepley PetscReal L[3] = {xmax - xmin, ymax - ymin, zmax - zmin};
949dc49dd9SMatthew G. Knepley
959dc49dd9SMatthew G. Knepley PetscCall(DMSetPeriodicity(da, maxCell, Lstart, L));
969dc49dd9SMatthew G. Knepley }
973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9847c6ae99SBarry Smith }
9947c6ae99SBarry Smith
1009ae14b6eSBarry Smith /*
1019ae14b6eSBarry Smith Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA
1029ae14b6eSBarry Smith */
DMDASelectFields(DM da,PetscInt * outfields,PetscInt ** fields)103d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields)
104d71ae5a4SJacob Faibussowitsch {
1054e6118eeSBarry Smith PetscInt step, ndisplayfields, *displayfields, k, j;
1064e6118eeSBarry Smith PetscBool flg;
1074e6118eeSBarry Smith
1084e6118eeSBarry Smith PetscFunctionBegin;
1099566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL));
1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(step, &displayfields));
1114e6118eeSBarry Smith for (k = 0; k < step; k++) displayfields[k] = k;
1124e6118eeSBarry Smith ndisplayfields = step;
1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg));
1144e6118eeSBarry Smith if (!ndisplayfields) ndisplayfields = step;
1154e6118eeSBarry Smith if (!flg) {
1164e6118eeSBarry Smith char **fields;
1174e6118eeSBarry Smith const char *fieldname;
1184e6118eeSBarry Smith PetscInt nfields = step;
1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(step, &fields));
1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg));
1214e6118eeSBarry Smith if (flg) {
1224e6118eeSBarry Smith ndisplayfields = 0;
1234e6118eeSBarry Smith for (k = 0; k < nfields; k++) {
1244e6118eeSBarry Smith for (j = 0; j < step; j++) {
1259566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, j, &fieldname));
1269566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fieldname, fields[k], &flg));
127ad540459SPierre Jolivet if (flg) goto found;
1284e6118eeSBarry Smith }
12998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]);
1309371c9d4SSatish Balay found:
1319371c9d4SSatish Balay displayfields[ndisplayfields++] = j;
1324e6118eeSBarry Smith }
1334e6118eeSBarry Smith }
13448a46eb9SPierre Jolivet for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k]));
1359566063dSJacob Faibussowitsch PetscCall(PetscFree(fields));
1364e6118eeSBarry Smith }
1374e6118eeSBarry Smith *fields = displayfields;
1384e6118eeSBarry Smith *outfields = ndisplayfields;
1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1404e6118eeSBarry Smith }
1414e6118eeSBarry Smith
1429804daf3SBarry Smith #include <petscdraw.h>
14373f7a4c5SBarry Smith
VecView_MPI_Draw_DA1d(Vec xin,PetscViewer v)144d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v)
145d71ae5a4SJacob Faibussowitsch {
1469a42bb27SBarry Smith DM da;
147e5ab1681SLisandro Dalcin PetscMPIInt rank, size, tag;
148e5ab1681SLisandro Dalcin PetscInt i, n, N, dof, istart, isize, j, nbounds;
14947c6ae99SBarry Smith MPI_Status status;
150e5ab1681SLisandro Dalcin PetscReal min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0;
15147c6ae99SBarry Smith const PetscScalar *array, *xg;
15247c6ae99SBarry Smith PetscDraw draw;
153e5ab1681SLisandro Dalcin PetscBool isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE;
15447c6ae99SBarry Smith MPI_Comm comm;
15547c6ae99SBarry Smith PetscDrawAxis axis;
15647c6ae99SBarry Smith Vec xcoor;
157bff4a2f0SMatthew G. Knepley DMBoundaryType bx;
158e5ab1681SLisandro Dalcin const char *tlabel = NULL, *xlabel = NULL;
159a56202b9SBarry Smith const PetscReal *bounds;
16003193ff8SBarry Smith PetscInt *displayfields;
16103193ff8SBarry Smith PetscInt k, ndisplayfields;
1624e6118eeSBarry Smith PetscBool hold;
163e5ab1681SLisandro Dalcin PetscDrawViewPorts *ports = NULL;
164e5ab1681SLisandro Dalcin PetscViewerFormat format;
16547c6ae99SBarry Smith
16647c6ae99SBarry Smith PetscFunctionBegin;
1679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
1689566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull));
1693ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS);
1709566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds));
17147c6ae99SBarry Smith
1729566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da));
1737a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
1749566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)xin, &comm));
1759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
1769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
17747c6ae99SBarry Smith
1789566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL));
17947c6ae99SBarry Smith
1809566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL));
1819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL));
1829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin, &array));
1839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xin, &n));
184e5ab1681SLisandro Dalcin n = n / dof;
18547c6ae99SBarry Smith
186832b7cebSLisandro Dalcin /* Get coordinates of nodes */
1879566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor));
18847c6ae99SBarry Smith if (!xcoor) {
1899566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0));
1909566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor));
19147c6ae99SBarry Smith }
1929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoor, &xg));
1939566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da, 0, &xlabel));
194832b7cebSLisandro Dalcin
195832b7cebSLisandro Dalcin /* Determine the min and max coordinate in plot */
196dd400576SPatrick Sanan if (rank == 0) xmin = PetscRealPart(xg[0]);
197e5ab1681SLisandro Dalcin if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]);
1989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm));
1999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm));
20047c6ae99SBarry Smith
2019566063dSJacob Faibussowitsch PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields));
2029566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(v, &format));
2039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL));
204832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE;
205832b7cebSLisandro Dalcin if (useports) {
2069566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
2079566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis));
2089566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw));
2099566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw));
2109566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports));
21173f7a4c5SBarry Smith }
212e5ab1681SLisandro Dalcin
213832b7cebSLisandro Dalcin /* Loop over each field; drawing each in a different window */
21403193ff8SBarry Smith for (k = 0; k < ndisplayfields; k++) {
21503193ff8SBarry Smith j = displayfields[k];
21647c6ae99SBarry Smith
217e5ab1681SLisandro Dalcin /* determine the min and max value in plot */
2189566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xin, j, NULL, &min));
2199566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xin, j, NULL, &max));
220a56202b9SBarry Smith if (j < nbounds) {
221a56202b9SBarry Smith min = PetscMin(min, bounds[2 * j]);
222a56202b9SBarry Smith max = PetscMax(max, bounds[2 * j + 1]);
223a56202b9SBarry Smith }
224e5ab1681SLisandro Dalcin if (min == max) {
225e5ab1681SLisandro Dalcin min -= 1.e-5;
226e5ab1681SLisandro Dalcin max += 1.e-5;
227e5ab1681SLisandro Dalcin }
22847c6ae99SBarry Smith
229832b7cebSLisandro Dalcin if (useports) {
2309566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, k));
2319566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, j, &tlabel));
232832b7cebSLisandro Dalcin } else {
233832b7cebSLisandro Dalcin const char *title;
2349566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetHold(v, &hold));
2359566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, k, &draw));
2369566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis));
2379566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, j, &title));
2389566063dSJacob Faibussowitsch if (title) PetscCall(PetscDrawSetTitle(draw, title));
2399566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw));
2409566063dSJacob Faibussowitsch if (!hold) PetscCall(PetscDrawClear(draw));
241832b7cebSLisandro Dalcin }
2429566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL));
2439566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max));
2449566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis));
24547c6ae99SBarry Smith
24647c6ae99SBarry Smith /* draw local part of vector */
2479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag));
24847c6ae99SBarry Smith if (rank < size - 1) { /*send value to right */
2499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm));
2509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm));
25147c6ae99SBarry Smith }
25247c6ae99SBarry Smith if (rank) { /* receive value from left */
2539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
2549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status));
255e5ab1681SLisandro Dalcin }
256d0609cedSBarry Smith PetscDrawCollectiveBegin(draw);
257e5ab1681SLisandro Dalcin if (rank) {
2589566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED));
2599566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK));
26047c6ae99SBarry Smith }
261e5ab1681SLisandro Dalcin for (i = 1; i < n; i++) {
2629566063dSJacob 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));
2639566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK));
26447c6ae99SBarry Smith }
265e5ab1681SLisandro Dalcin if (rank == size - 1) {
2669566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK));
26747c6ae99SBarry Smith }
268d0609cedSBarry Smith PetscDrawCollectiveEnd(draw);
2699566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw));
2709566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw));
2719566063dSJacob Faibussowitsch if (!useports) PetscCall(PetscDrawSave(draw));
27247c6ae99SBarry Smith }
273832b7cebSLisandro Dalcin if (useports) {
2749566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
2759566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw));
276832b7cebSLisandro Dalcin }
277832b7cebSLisandro Dalcin
2789566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ports));
2799566063dSJacob Faibussowitsch PetscCall(PetscFree(displayfields));
2809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoor, &xg));
2819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin, &array));
2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28347c6ae99SBarry Smith }
284