1 2 /* 3 Plots vectors obtained with DMDACreate1d() 4 */ 5 6 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7 8 /*@ 9 DMDASetUniformCoordinates - Sets a `DMDA` coordinates to be a uniform grid 10 11 Collective 12 13 Input Parameters: 14 + da - the distributed array object 15 . xmin - min extreme in the x direction 16 . xmax - max extreme in the x direction 17 . ymin - min extreme in the y direction (value ignored for 1 dimensional problems) 18 . ymax - max extreme in the y direction (value ignored for 1 dimensional problems) 19 . zmin - min extreme in the z direction (value ignored for 1 or 2 dimensional problems) 20 - zmax - max extreme in the z direction (value ignored for 1 or 2 dimensional problems) 21 22 Level: beginner 23 24 .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMStagSetUniformCoordinates()` 25 @*/ 26 PetscErrorCode DMDASetUniformCoordinates(DM da, PetscReal xmin, PetscReal xmax, PetscReal ymin, PetscReal ymax, PetscReal zmin, PetscReal zmax) 27 { 28 MPI_Comm comm; 29 DM cda; 30 DM_DA *dd = (DM_DA *)da->data; 31 DMBoundaryType bx, by, bz; 32 Vec xcoor; 33 PetscScalar *coors; 34 PetscReal hx, hy, hz_; 35 PetscInt i, j, k, M, N, P, istart, isize, jstart, jsize, kstart, ksize, dim, cnt; 36 37 PetscFunctionBegin; 38 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 39 PetscCheck(dd->gtol, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONGSTATE, "Cannot set coordinates until after DMDA has been setup"); 40 PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL)); 41 PetscCheck(xmax >= xmin, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "xmax must be larger than xmin %g %g", (double)xmin, (double)xmax); 42 PetscCheck(!(dim > 1) || !(ymax < ymin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "ymax must be larger than ymin %g %g", (double)ymin, (double)ymax); 43 PetscCheck(!(dim > 2) || !(zmax < zmin), PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_INCOMP, "zmax must be larger than zmin %g %g", (double)zmin, (double)zmax); 44 PetscCall(PetscObjectGetComm((PetscObject)da, &comm)); 45 PetscCall(DMDAGetCorners(da, &istart, &jstart, &kstart, &isize, &jsize, &ksize)); 46 PetscCall(DMGetCoordinateDM(da, &cda)); 47 PetscCall(DMCreateGlobalVector(cda, &xcoor)); 48 if (dim == 1) { 49 if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / M; 50 else hx = (xmax - xmin) / (M - 1); 51 PetscCall(VecGetArray(xcoor, &coors)); 52 for (i = 0; i < isize; i++) coors[i] = xmin + hx * (i + istart); 53 PetscCall(VecRestoreArray(xcoor, &coors)); 54 } else if (dim == 2) { 55 if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M); 56 else hx = (xmax - xmin) / (M - 1); 57 if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N); 58 else hy = (ymax - ymin) / (N - 1); 59 PetscCall(VecGetArray(xcoor, &coors)); 60 cnt = 0; 61 for (j = 0; j < jsize; j++) { 62 for (i = 0; i < isize; i++) { 63 coors[cnt++] = xmin + hx * (i + istart); 64 coors[cnt++] = ymin + hy * (j + jstart); 65 } 66 } 67 PetscCall(VecRestoreArray(xcoor, &coors)); 68 } else if (dim == 3) { 69 if (bx == DM_BOUNDARY_PERIODIC) hx = (xmax - xmin) / (M); 70 else hx = (xmax - xmin) / (M - 1); 71 if (by == DM_BOUNDARY_PERIODIC) hy = (ymax - ymin) / (N); 72 else hy = (ymax - ymin) / (N - 1); 73 if (bz == DM_BOUNDARY_PERIODIC) hz_ = (zmax - zmin) / (P); 74 else hz_ = (zmax - zmin) / (P - 1); 75 PetscCall(VecGetArray(xcoor, &coors)); 76 cnt = 0; 77 for (k = 0; k < ksize; k++) { 78 for (j = 0; j < jsize; j++) { 79 for (i = 0; i < isize; i++) { 80 coors[cnt++] = xmin + hx * (i + istart); 81 coors[cnt++] = ymin + hy * (j + jstart); 82 coors[cnt++] = zmin + hz_ * (k + kstart); 83 } 84 } 85 } 86 PetscCall(VecRestoreArray(xcoor, &coors)); 87 } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Cannot create uniform coordinates for this dimension %" PetscInt_FMT, dim); 88 PetscCall(DMSetCoordinates(da, xcoor)); 89 PetscCall(VecDestroy(&xcoor)); 90 PetscFunctionReturn(PETSC_SUCCESS); 91 } 92 93 /* 94 Allows a user to select a subset of the fields to be drawn by VecView() when the vector comes from a DMDA 95 */ 96 PetscErrorCode DMDASelectFields(DM da, PetscInt *outfields, PetscInt **fields) 97 { 98 PetscInt step, ndisplayfields, *displayfields, k, j; 99 PetscBool flg; 100 101 PetscFunctionBegin; 102 PetscCall(DMDAGetInfo(da, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &step, NULL, NULL, NULL, NULL, NULL)); 103 PetscCall(PetscMalloc1(step, &displayfields)); 104 for (k = 0; k < step; k++) displayfields[k] = k; 105 ndisplayfields = step; 106 PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-draw_fields", displayfields, &ndisplayfields, &flg)); 107 if (!ndisplayfields) ndisplayfields = step; 108 if (!flg) { 109 char **fields; 110 const char *fieldname; 111 PetscInt nfields = step; 112 PetscCall(PetscMalloc1(step, &fields)); 113 PetscCall(PetscOptionsGetStringArray(NULL, NULL, "-draw_fields_by_name", fields, &nfields, &flg)); 114 if (flg) { 115 ndisplayfields = 0; 116 for (k = 0; k < nfields; k++) { 117 for (j = 0; j < step; j++) { 118 PetscCall(DMDAGetFieldName(da, j, &fieldname)); 119 PetscCall(PetscStrcmp(fieldname, fields[k], &flg)); 120 if (flg) goto found; 121 } 122 SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_USER, "Unknown fieldname %s", fields[k]); 123 found: 124 displayfields[ndisplayfields++] = j; 125 } 126 } 127 for (k = 0; k < nfields; k++) PetscCall(PetscFree(fields[k])); 128 PetscCall(PetscFree(fields)); 129 } 130 *fields = displayfields; 131 *outfields = ndisplayfields; 132 PetscFunctionReturn(PETSC_SUCCESS); 133 } 134 135 #include <petscdraw.h> 136 137 PetscErrorCode VecView_MPI_Draw_DA1d(Vec xin, PetscViewer v) 138 { 139 DM da; 140 PetscMPIInt rank, size, tag; 141 PetscInt i, n, N, dof, istart, isize, j, nbounds; 142 MPI_Status status; 143 PetscReal min, max, xmin = 0.0, xmax = 0.0, tmp = 0.0, xgtmp = 0.0; 144 const PetscScalar *array, *xg; 145 PetscDraw draw; 146 PetscBool isnull, useports = PETSC_FALSE, showmarkers = PETSC_FALSE; 147 MPI_Comm comm; 148 PetscDrawAxis axis; 149 Vec xcoor; 150 DMBoundaryType bx; 151 const char *tlabel = NULL, *xlabel = NULL; 152 const PetscReal *bounds; 153 PetscInt *displayfields; 154 PetscInt k, ndisplayfields; 155 PetscBool hold; 156 PetscDrawViewPorts *ports = NULL; 157 PetscViewerFormat format; 158 159 PetscFunctionBegin; 160 PetscCall(PetscViewerDrawGetDraw(v, 0, &draw)); 161 PetscCall(PetscDrawIsNull(draw, &isnull)); 162 if (isnull) PetscFunctionReturn(PETSC_SUCCESS); 163 PetscCall(PetscViewerDrawGetBounds(v, &nbounds, &bounds)); 164 165 PetscCall(VecGetDM(xin, &da)); 166 PetscCheck(da, PetscObjectComm((PetscObject)xin), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA"); 167 PetscCall(PetscObjectGetComm((PetscObject)xin, &comm)); 168 PetscCallMPI(MPI_Comm_size(comm, &size)); 169 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 170 171 PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_vec_use_markers", &showmarkers, NULL)); 172 173 PetscCall(DMDAGetInfo(da, NULL, &N, NULL, NULL, NULL, NULL, NULL, &dof, NULL, &bx, NULL, NULL, NULL)); 174 PetscCall(DMDAGetCorners(da, &istart, NULL, NULL, &isize, NULL, NULL)); 175 PetscCall(VecGetArrayRead(xin, &array)); 176 PetscCall(VecGetLocalSize(xin, &n)); 177 n = n / dof; 178 179 /* Get coordinates of nodes */ 180 PetscCall(DMGetCoordinates(da, &xcoor)); 181 if (!xcoor) { 182 PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0)); 183 PetscCall(DMGetCoordinates(da, &xcoor)); 184 } 185 PetscCall(VecGetArrayRead(xcoor, &xg)); 186 PetscCall(DMDAGetCoordinateName(da, 0, &xlabel)); 187 188 /* Determine the min and max coordinate in plot */ 189 if (rank == 0) xmin = PetscRealPart(xg[0]); 190 if (rank == size - 1) xmax = PetscRealPart(xg[n - 1]); 191 PetscCallMPI(MPI_Bcast(&xmin, 1, MPIU_REAL, 0, comm)); 192 PetscCallMPI(MPI_Bcast(&xmax, 1, MPIU_REAL, size - 1, comm)); 193 194 PetscCall(DMDASelectFields(da, &ndisplayfields, &displayfields)); 195 PetscCall(PetscViewerGetFormat(v, &format)); 196 PetscCall(PetscOptionsGetBool(NULL, NULL, "-draw_ports", &useports, NULL)); 197 if (format == PETSC_VIEWER_DRAW_PORTS) useports = PETSC_TRUE; 198 if (useports) { 199 PetscCall(PetscViewerDrawGetDraw(v, 0, &draw)); 200 PetscCall(PetscViewerDrawGetDrawAxis(v, 0, &axis)); 201 PetscCall(PetscDrawCheckResizedWindow(draw)); 202 PetscCall(PetscDrawClear(draw)); 203 PetscCall(PetscDrawViewPortsCreate(draw, ndisplayfields, &ports)); 204 } 205 206 /* Loop over each field; drawing each in a different window */ 207 for (k = 0; k < ndisplayfields; k++) { 208 j = displayfields[k]; 209 210 /* determine the min and max value in plot */ 211 PetscCall(VecStrideMin(xin, j, NULL, &min)); 212 PetscCall(VecStrideMax(xin, j, NULL, &max)); 213 if (j < nbounds) { 214 min = PetscMin(min, bounds[2 * j]); 215 max = PetscMax(max, bounds[2 * j + 1]); 216 } 217 if (min == max) { 218 min -= 1.e-5; 219 max += 1.e-5; 220 } 221 222 if (useports) { 223 PetscCall(PetscDrawViewPortsSet(ports, k)); 224 PetscCall(DMDAGetFieldName(da, j, &tlabel)); 225 } else { 226 const char *title; 227 PetscCall(PetscViewerDrawGetHold(v, &hold)); 228 PetscCall(PetscViewerDrawGetDraw(v, k, &draw)); 229 PetscCall(PetscViewerDrawGetDrawAxis(v, k, &axis)); 230 PetscCall(DMDAGetFieldName(da, j, &title)); 231 if (title) PetscCall(PetscDrawSetTitle(draw, title)); 232 PetscCall(PetscDrawCheckResizedWindow(draw)); 233 if (!hold) PetscCall(PetscDrawClear(draw)); 234 } 235 PetscCall(PetscDrawAxisSetLabels(axis, tlabel, xlabel, NULL)); 236 PetscCall(PetscDrawAxisSetLimits(axis, xmin, xmax, min, max)); 237 PetscCall(PetscDrawAxisDraw(axis)); 238 239 /* draw local part of vector */ 240 PetscCall(PetscObjectGetNewTag((PetscObject)xin, &tag)); 241 if (rank < size - 1) { /*send value to right */ 242 PetscCallMPI(MPI_Send((void *)&xg[n - 1], 1, MPIU_REAL, rank + 1, tag, comm)); 243 PetscCallMPI(MPI_Send((void *)&array[j + (n - 1) * dof], 1, MPIU_REAL, rank + 1, tag, comm)); 244 } 245 if (rank) { /* receive value from left */ 246 PetscCallMPI(MPI_Recv(&xgtmp, 1, MPIU_REAL, rank - 1, tag, comm, &status)); 247 PetscCallMPI(MPI_Recv(&tmp, 1, MPIU_REAL, rank - 1, tag, comm, &status)); 248 } 249 PetscDrawCollectiveBegin(draw); 250 if (rank) { 251 PetscCall(PetscDrawLine(draw, xgtmp, tmp, PetscRealPart(xg[0]), PetscRealPart(array[j]), PETSC_DRAW_RED)); 252 if (showmarkers) PetscCall(PetscDrawPoint(draw, xgtmp, tmp, PETSC_DRAW_BLACK)); 253 } 254 for (i = 1; i < n; i++) { 255 PetscCall(PetscDrawLine(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PetscRealPart(xg[i]), PetscRealPart(array[j + dof * i]), PETSC_DRAW_RED)); 256 if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[i - 1]), PetscRealPart(array[j + dof * (i - 1)]), PETSC_DRAW_BLACK)); 257 } 258 if (rank == size - 1) { 259 if (showmarkers) PetscCall(PetscDrawMarker(draw, PetscRealPart(xg[n - 1]), PetscRealPart(array[j + dof * (n - 1)]), PETSC_DRAW_BLACK)); 260 } 261 PetscDrawCollectiveEnd(draw); 262 PetscCall(PetscDrawFlush(draw)); 263 PetscCall(PetscDrawPause(draw)); 264 if (!useports) PetscCall(PetscDrawSave(draw)); 265 } 266 if (useports) { 267 PetscCall(PetscViewerDrawGetDraw(v, 0, &draw)); 268 PetscCall(PetscDrawSave(draw)); 269 } 270 271 PetscCall(PetscDrawViewPortsDestroy(ports)); 272 PetscCall(PetscFree(displayfields)); 273 PetscCall(VecRestoreArrayRead(xcoor, &xg)); 274 PetscCall(VecRestoreArrayRead(xin, &array)); 275 PetscFunctionReturn(PETSC_SUCCESS); 276 } 277