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 /*@ 9dce8aebaSBarry Smith DMDASetUniformCoordinates - Sets a `DMDA` coordinates to be a uniform grid 1047c6ae99SBarry Smith 1120f4b53cSBarry Smith Collective 1247c6ae99SBarry Smith 1347c6ae99SBarry Smith Input Parameters: 1447c6ae99SBarry Smith + da - the distributed array object 15*a4e35b19SJacob Faibussowitsch . xmin - min extreme in the x direction 16*a4e35b19SJacob Faibussowitsch . xmax - max extreme in the x direction 17*a4e35b19SJacob Faibussowitsch . ymin - min extreme in the y direction (value ignored for 1 dimensional problems) 18*a4e35b19SJacob Faibussowitsch . ymax - max extreme in the y direction (value ignored for 1 dimensional problems) 19*a4e35b19SJacob Faibussowitsch . zmin - min extreme in the z direction (value ignored for 1 or 2 dimensional problems) 20*a4e35b19SJacob Faibussowitsch - zmax - max extreme in the z direction (value ignored for 1 or 2 dimensional problems) 2147c6ae99SBarry Smith 2247c6ae99SBarry Smith Level: beginner 2347c6ae99SBarry Smith 24dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMStagSetUniformCoordinates()` 2547c6ae99SBarry Smith @*/ 26d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax) 27d71ae5a4SJacob Faibussowitsch { 2847c6ae99SBarry Smith MPI_Comm comm; 299a42bb27SBarry Smith DM cda; 30a5bc1bf3SBarry Smith DM_DA *dd = (DM_DA *)da->data; 31bff4a2f0SMatthew G. Knepley DMBoundaryType 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 3747c6ae99SBarry Smith PetscFunctionBegin; 38a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 397a8be351SBarry Smith PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup"); 409566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL)); 4108401ef6SPierre Jolivet PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax); 4208401ef6SPierre 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); 4308401ef6SPierre 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); 449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 459566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize)); 469566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 479566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(cda, &xcoor)); 4847c6ae99SBarry Smith if (dim == 1) { 49bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / M; 50ce00eea3SSatish Balay else hx = (xmax - xmin) / (M - 1); 519566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor, &coors)); 52ad540459SPierre Jolivet for (i = 0; i < isize; i++) coors[i] = xmin + hx * (i + istart); 539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor, &coors)); 5447c6ae99SBarry Smith } else if (dim == 2) { 55bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M); 5647c6ae99SBarry Smith else hx = (xmax - xmin) / (M - 1); 57bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N); 5847c6ae99SBarry Smith else hy = (ymax - ymin) / (N - 1); 599566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor, &coors)); 6047c6ae99SBarry Smith cnt = 0; 6147c6ae99SBarry Smith for (j = 0; j < jsize; j++) { 6247c6ae99SBarry Smith for (i = 0; i < isize; i++) { 6347c6ae99SBarry Smith coors[cnt++] = xmin + hx * (i + istart); 6447c6ae99SBarry Smith coors[cnt++] = ymin + hy * (j + jstart); 6547c6ae99SBarry Smith } 6647c6ae99SBarry Smith } 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor, &coors)); 6847c6ae99SBarry Smith } else if (dim == 3) { 69bff4a2f0SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M); 7047c6ae99SBarry Smith else hx = (xmax - xmin) / (M - 1); 71bff4a2f0SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N); 7247c6ae99SBarry Smith else hy = (ymax - ymin) / (N - 1); 73bff4a2f0SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax - zmin) / (P); 7447c6ae99SBarry Smith else hz_ = (zmax - zmin) / (P - 1); 759566063dSJacob Faibussowitsch PetscCall(VecGetArray(xcoor, &coors)); 7647c6ae99SBarry Smith cnt = 0; 7747c6ae99SBarry Smith for (k = 0; k < ksize; k++) { 7847c6ae99SBarry Smith for (j = 0; j < jsize; j++) { 7947c6ae99SBarry Smith for (i = 0; i < isize; i++) { 8047c6ae99SBarry Smith coors[cnt++] = xmin + hx * (i + istart); 8147c6ae99SBarry Smith coors[cnt++] = ymin + hy * (j + jstart); 8247c6ae99SBarry Smith coors[cnt++] = zmin + hz_ * (k + kstart); 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith } 8547c6ae99SBarry Smith } 869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xcoor, &coors)); 8763a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim); 889566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, xcoor)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xcoor)); 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 */ 96d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields) 97d71ae5a4SJacob Faibussowitsch { 984e6118eeSBarry Smith PetscInt step, ndisplayfields, *displayfields, k, j; 994e6118eeSBarry Smith PetscBool flg; 1004e6118eeSBarry Smith 1014e6118eeSBarry Smith PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL)); 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(step, &displayfields)); 1044e6118eeSBarry Smith for (k = 0; k < step; k++) displayfields[k] = k; 1054e6118eeSBarry Smith ndisplayfields = step; 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg)); 1074e6118eeSBarry Smith if (!ndisplayfields) ndisplayfields = step; 1084e6118eeSBarry Smith if (!flg) { 1094e6118eeSBarry Smith char **fields; 1104e6118eeSBarry Smith const char *fieldname; 1114e6118eeSBarry Smith PetscInt nfields = step; 1129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(step, &fields)); 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg)); 1144e6118eeSBarry Smith if (flg) { 1154e6118eeSBarry Smith ndisplayfields = 0; 1164e6118eeSBarry Smith for (k = 0; k < nfields; k++) { 1174e6118eeSBarry Smith for (j = 0; j < step; j++) { 1189566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, j, &fieldname)); 1199566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fieldname, fields[k], &flg)); 120ad540459SPierre Jolivet if (flg) goto found; 1214e6118eeSBarry Smith } 12298921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]); 1239371c9d4SSatish Balay found: 1249371c9d4SSatish Balay displayfields[ndisplayfields++] = j; 1254e6118eeSBarry Smith } 1264e6118eeSBarry Smith } 12748a46eb9SPierre Jolivet for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k])); 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(fields)); 1294e6118eeSBarry Smith } 1304e6118eeSBarry Smith *fields = displayfields; 1314e6118eeSBarry Smith *outfields = ndisplayfields; 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1334e6118eeSBarry Smith } 1344e6118eeSBarry Smith 1359804daf3SBarry Smith #include <petscdraw.h> 13673f7a4c5SBarry Smith 137d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v) 138d71ae5a4SJacob Faibussowitsch { 1399a42bb27SBarry Smith DM da; 140e5ab1681SLisandro Dalcin PetscMPIInt rank, size, tag; 141e5ab1681SLisandro Dalcin PetscInt i, n, N, dof, istart, isize, j, nbounds; 14247c6ae99SBarry Smith MPI_Status status; 143e5ab1681SLisandro Dalcin PetscReal min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0; 14447c6ae99SBarry Smith const PetscScalar *array, *xg; 14547c6ae99SBarry Smith PetscDraw draw; 146e5ab1681SLisandro Dalcin PetscBool isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE; 14747c6ae99SBarry Smith MPI_Comm comm; 14847c6ae99SBarry Smith PetscDrawAxis axis; 14947c6ae99SBarry Smith Vec xcoor; 150bff4a2f0SMatthew G. Knepley DMBoundaryType bx; 151e5ab1681SLisandro Dalcin const char *tlabel = NULL, *xlabel = NULL; 152a56202b9SBarry Smith const PetscReal *bounds; 15303193ff8SBarry Smith PetscInt *displayfields; 15403193ff8SBarry Smith PetscInt k, ndisplayfields; 1554e6118eeSBarry Smith PetscBool hold; 156e5ab1681SLisandro Dalcin PetscDrawViewPorts *ports = NULL; 157e5ab1681SLisandro Dalcin PetscViewerFormat format; 15847c6ae99SBarry Smith 15947c6ae99SBarry Smith PetscFunctionBegin; 1609566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, 0, &draw)); 1619566063dSJacob Faibussowitsch PetscCall(PetscDrawIsNull(draw, &isnull)); 1623ba16761SJacob Faibussowitsch if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 1639566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds)); 16447c6ae99SBarry Smith 1659566063dSJacob Faibussowitsch PetscCall(VecGetDM(xin, &da)); 1667a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 1679566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)xin, &comm)); 1689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 17047c6ae99SBarry Smith 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL)); 17247c6ae99SBarry Smith 1739566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL)); 1749566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL)); 1759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xin, &array)); 1769566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xin, &n)); 177e5ab1681SLisandro Dalcin n = n / dof; 17847c6ae99SBarry Smith 179832b7cebSLisandro Dalcin /* Get coordinates of nodes */ 1809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor)); 18147c6ae99SBarry Smith if (!xcoor) { 1829566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0)); 1839566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor)); 18447c6ae99SBarry Smith } 1859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoor, &xg)); 1869566063dSJacob Faibussowitsch PetscCall(DMDAGetCoordinateName(da, 0, &xlabel)); 187832b7cebSLisandro Dalcin 188832b7cebSLisandro Dalcin /* Determine the min and max coordinate in plot */ 189dd400576SPatrick Sanan if (rank == 0) xmin = PetscRealPart(xg[0]); 190e5ab1681SLisandro Dalcin if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]); 1919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm)); 1929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm)); 19347c6ae99SBarry Smith 1949566063dSJacob Faibussowitsch PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields)); 1959566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(v, &format)); 1969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL)); 197832b7cebSLisandro Dalcin if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 198832b7cebSLisandro Dalcin if (useports) { 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, 0, &draw)); 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 2029566063dSJacob Faibussowitsch PetscCall(PetscDrawClear(draw)); 2039566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports)); 20473f7a4c5SBarry Smith } 205e5ab1681SLisandro Dalcin 206832b7cebSLisandro Dalcin /* Loop over each field; drawing each in a different window */ 20703193ff8SBarry Smith for (k = 0; k < ndisplayfields; k++) { 20803193ff8SBarry Smith j = displayfields[k]; 20947c6ae99SBarry Smith 210e5ab1681SLisandro Dalcin /* determine the min and max value in plot */ 2119566063dSJacob Faibussowitsch PetscCall(VecStrideMin(xin, j, NULL, &min)); 2129566063dSJacob Faibussowitsch PetscCall(VecStrideMax(xin, j, NULL, &max)); 213a56202b9SBarry Smith if (j < nbounds) { 214a56202b9SBarry Smith min = PetscMin(min, bounds[2 * j]); 215a56202b9SBarry Smith max = PetscMax(max, bounds[2 * j + 1]); 216a56202b9SBarry Smith } 217e5ab1681SLisandro Dalcin if (min == max) { 218e5ab1681SLisandro Dalcin min -= 1.e-5; 219e5ab1681SLisandro Dalcin max += 1.e-5; 220e5ab1681SLisandro Dalcin } 22147c6ae99SBarry Smith 222832b7cebSLisandro Dalcin if (useports) { 2239566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, k)); 2249566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, j, &tlabel)); 225832b7cebSLisandro Dalcin } else { 226832b7cebSLisandro Dalcin const char *title; 2279566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetHold(v, &hold)); 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, k, &draw)); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis)); 2309566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da, j, &title)); 2319566063dSJacob Faibussowitsch if (title) PetscCall(PetscDrawSetTitle(draw, title)); 2329566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 2339566063dSJacob Faibussowitsch if (!hold) PetscCall(PetscDrawClear(draw)); 234832b7cebSLisandro Dalcin } 2359566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL)); 2369566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max)); 2379566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisDraw(axis)); 23847c6ae99SBarry Smith 23947c6ae99SBarry Smith /* draw local part of vector */ 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag)); 24147c6ae99SBarry Smith if (rank < size - 1) { /*send value to right */ 2429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm)); 2439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm)); 24447c6ae99SBarry Smith } 24547c6ae99SBarry Smith if (rank) { /* receive value from left */ 2469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status)); 2479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status)); 248e5ab1681SLisandro Dalcin } 249d0609cedSBarry Smith PetscDrawCollectiveBegin(draw); 250e5ab1681SLisandro Dalcin if (rank) { 2519566063dSJacob Faibussowitsch PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED)); 2529566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK)); 25347c6ae99SBarry Smith } 254e5ab1681SLisandro Dalcin for (i = 1; i < n; i++) { 2559566063dSJacob 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)); 2569566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK)); 25747c6ae99SBarry Smith } 258e5ab1681SLisandro Dalcin if (rank == size - 1) { 2599566063dSJacob Faibussowitsch if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK)); 26047c6ae99SBarry Smith } 261d0609cedSBarry Smith PetscDrawCollectiveEnd(draw); 2629566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 2639566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 2649566063dSJacob Faibussowitsch if (!useports) PetscCall(PetscDrawSave(draw)); 26547c6ae99SBarry Smith } 266832b7cebSLisandro Dalcin if (useports) { 2679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(v, 0, &draw)); 2689566063dSJacob Faibussowitsch PetscCall(PetscDrawSave(draw)); 269832b7cebSLisandro Dalcin } 270832b7cebSLisandro Dalcin 2719566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ports)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(displayfields)); 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoor, &xg)); 2749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xin, &array)); 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27647c6ae99SBarry Smith } 277