xref: /petsc/src/dm/impls/da/grglvis.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1 /* Routines to visualize DMDAs and fields through GLVis */
2 
3 #include <petsc/private/dmdaimpl.h>
4 #include <petsc/private/glvisviewerimpl.h>
5 
6 typedef struct {
7   PetscBool ll;
8 } DMDAGhostedGLVisViewerCtx;
9 
10 static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(void **vctx)
11 {
12   PetscFunctionBegin;
13   PetscCall(PetscFree(*vctx));
14   PetscFunctionReturn(0);
15 }
16 
17 typedef struct {
18   Vec xlocal;
19 } DMDAFieldGLVisViewerCtx;
20 
21 static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx)
22 {
23   DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx *)vctx;
24 
25   PetscFunctionBegin;
26   PetscCall(VecDestroy(&ctx->xlocal));
27   PetscCall(PetscFree(vctx));
28   PetscFunctionReturn(0);
29 }
30 
31 /*
32    dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right
33    dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left
34 */
35 static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
36 {
37   DMDAGhostedGLVisViewerCtx *dactx;
38   PetscInt                   sx, sy, sz, ien, jen, ken;
39 
40   PetscFunctionBegin;
41   /* Appease -Wmaybe-uninitialized */
42   if (nex) *nex = -1;
43   if (ney) *ney = -1;
44   if (nez) *nez = -1;
45   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &ien, &jen, &ken));
46   PetscCall(DMGetApplicationContext(da, &dactx));
47   if (dactx->ll) {
48     PetscInt dim;
49 
50     PetscCall(DMGetDimension(da, &dim));
51     if (!sx) ien--;
52     if (!sy && dim > 1) jen--;
53     if (!sz && dim > 2) ken--;
54   } else {
55     PetscInt M, N, P;
56 
57     PetscCall(DMDAGetInfo(da, NULL, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
58     if (sx + ien == M) ien--;
59     if (sy + jen == N) jen--;
60     if (sz + ken == P) ken--;
61   }
62   if (nex) *nex = ien;
63   if (ney) *ney = jen;
64   if (nez) *nez = ken;
65   PetscFunctionReturn(0);
66 }
67 
68 /* inherits number of vertices from DMDAGetNumElementsGhosted */
69 static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
70 {
71   PetscInt ien = 0, jen = 0, ken = 0, dim;
72   PetscInt tote;
73 
74   PetscFunctionBegin;
75   PetscCall(DMGetDimension(da, &dim));
76   PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
77   tote = ien * (dim > 1 ? jen : 1) * (dim > 2 ? ken : 1);
78   if (tote) {
79     ien = ien + 1;
80     jen = dim > 1 ? jen + 1 : 1;
81     ken = dim > 2 ? ken + 1 : 1;
82   }
83   if (nvx) *nvx = ien;
84   if (nvy) *nvy = jen;
85   if (nvz) *nvz = ken;
86   PetscFunctionReturn(0);
87 }
88 
89 static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
90 {
91   DM                         da;
92   DMDAFieldGLVisViewerCtx   *ctx = (DMDAFieldGLVisViewerCtx *)vctx;
93   DMDAGhostedGLVisViewerCtx *dactx;
94   const PetscScalar         *array;
95   PetscScalar              **arrayf;
96   PetscInt                   i, f, ii, ien, jen, ken, ie, je, ke, bs, *bss;
97   PetscInt                   sx, sy, sz, gsx, gsy, gsz, ist, jst, kst, gm, gn, gp;
98 
99   PetscFunctionBegin;
100   PetscCall(VecGetDM(ctx->xlocal, &da));
101   PetscCheck(da, PetscObjectComm(oX), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
102   PetscCall(DMGetApplicationContext(da, &dactx));
103   PetscCall(VecGetBlockSize(ctx->xlocal, &bs));
104   PetscCall(DMGlobalToLocalBegin(da, (Vec)oX, INSERT_VALUES, ctx->xlocal));
105   PetscCall(DMGlobalToLocalEnd(da, (Vec)oX, INSERT_VALUES, ctx->xlocal));
106   PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
107   PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
108   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
109   if (dactx->ll) {
110     kst = jst = ist = 0;
111   } else {
112     kst = gsz != sz ? 1 : 0;
113     jst = gsy != sy ? 1 : 0;
114     ist = gsx != sx ? 1 : 0;
115   }
116   PetscCall(PetscMalloc2(nf, &arrayf, nf, &bss));
117   PetscCall(VecGetArrayRead(ctx->xlocal, &array));
118   for (f = 0; f < nf; f++) {
119     PetscCall(VecGetBlockSize((Vec)oXf[f], &bss[f]));
120     PetscCall(VecGetArray((Vec)oXf[f], &arrayf[f]));
121   }
122   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
123     for (je = jst; je < jst + jen; je++) {
124       for (ie = ist; ie < ist + ien; ie++) {
125         PetscInt cf, b;
126         i = ke * gm * gn + je * gm + ie;
127         for (f = 0, cf = 0; f < nf; f++)
128           for (b = 0; b < bss[f]; b++) arrayf[f][bss[f] * ii + b] = array[i * bs + cf++];
129         ii++;
130       }
131     }
132   }
133   for (f = 0; f < nf; f++) PetscCall(VecRestoreArray((Vec)oXf[f], &arrayf[f]));
134   PetscCall(VecRestoreArrayRead(ctx->xlocal, &array));
135   PetscCall(PetscFree2(arrayf, bss));
136   PetscFunctionReturn(0);
137 }
138 
139 PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
140 {
141   DM da = (DM)oda, daview;
142 
143   PetscFunctionBegin;
144   PetscCall(PetscObjectQuery(oda, "GLVisGraphicsDMDAGhosted", (PetscObject *)&daview));
145   if (!daview) {
146     DMDAGhostedGLVisViewerCtx *dactx;
147     DM                         dacoord = NULL;
148     Vec                        xcoor, xcoorl;
149     PetscBool                  hashocoord = PETSC_FALSE;
150     const PetscInt            *lx, *ly, *lz;
151     PetscInt                   dim, M, N, P, m, n, p, dof, s, i;
152 
153     PetscCall(PetscNew(&dactx));
154     dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
155     PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Options", "PetscViewer");
156     PetscCall(PetscOptionsBool("-viewer_glvis_dm_da_ll", "Left-looking subdomain view", NULL, dactx->ll, &dactx->ll, NULL));
157     PetscOptionsEnd();
158     /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
159     PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, &m, &n, &p, &dof, &s, NULL, NULL, NULL, NULL));
160     PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz));
161     PetscCall(PetscObjectQuery((PetscObject)da, "_glvis_mesh_coords", (PetscObject *)&xcoor));
162     if (!xcoor) {
163       PetscCall(DMGetCoordinates(da, &xcoor));
164     } else {
165       hashocoord = PETSC_TRUE;
166     }
167     PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing GLVis graphics\n"));
168     switch (dim) {
169     case 1:
170       PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, dof, 1, lx, &daview));
171       if (!hashocoord) PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, 1, 1, lx, &dacoord));
172       break;
173     case 2:
174       PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, dof, 1, lx, ly, &daview));
175       if (!hashocoord) PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, 2, 1, lx, ly, &dacoord));
176       break;
177     case 3:
178       PetscCall(DMDACreate3d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, P, m, n, p, dof, 1, lx, ly, lz, &daview));
179       if (!hashocoord) PetscCall(DMDACreate3d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, P, m, n, p, 3, 1, lx, ly, lz, &dacoord));
180       break;
181     default:
182       SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
183     }
184     PetscCall(DMSetApplicationContext(daview, dactx));
185     PetscCall(DMSetApplicationContextDestroy(daview, DMDAGhostedDestroyGLVisViewerCtx_Private));
186     PetscCall(DMSetUp(daview));
187     if (!xcoor) {
188       PetscCall(DMDASetUniformCoordinates(daview, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
189       PetscCall(DMGetCoordinates(daview, &xcoor));
190     }
191     if (dacoord) {
192       PetscCall(DMSetUp(dacoord));
193       PetscCall(DMCreateLocalVector(dacoord, &xcoorl));
194       PetscCall(DMGlobalToLocalBegin(dacoord, xcoor, INSERT_VALUES, xcoorl));
195       PetscCall(DMGlobalToLocalEnd(dacoord, xcoor, INSERT_VALUES, xcoorl));
196       PetscCall(DMDestroy(&dacoord));
197     } else {
198       PetscInt    ien, jen, ken, nc, nl, cdof, deg;
199       char        fecmesh[64];
200       const char *name;
201       PetscBool   flg;
202 
203       PetscCall(DMDAGetNumElementsGhosted(daview, &ien, &jen, &ken));
204       nc = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
205 
206       PetscCall(VecGetLocalSize(xcoor, &nl));
207       PetscCheck(!nc || (nl % nc) == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Incompatible local coordinate size %" PetscInt_FMT " and number of cells %" PetscInt_FMT, nl, nc);
208       PetscCall(VecDuplicate(xcoor, &xcoorl));
209       PetscCall(VecCopy(xcoor, xcoorl));
210       PetscCall(VecSetDM(xcoorl, NULL));
211       PetscCall(PetscObjectGetName((PetscObject)xcoor, &name));
212       PetscCall(PetscStrbeginswith(name, "FiniteElementCollection:", &flg));
213       if (!flg) {
214         deg = 0;
215         if (nc && nl) {
216           cdof = nl / (nc * dim);
217           deg  = 1;
218           while (1) {
219             PetscInt degd = 1;
220             for (i = 0; i < dim; i++) degd *= (deg + 1);
221             PetscCheck(degd <= cdof, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell dofs %" PetscInt_FMT, cdof);
222             if (degd == cdof) break;
223             deg++;
224           }
225         }
226         PetscCall(PetscSNPrintf(fecmesh, sizeof(fecmesh), "FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, deg));
227         PetscCall(PetscObjectSetName((PetscObject)xcoorl, fecmesh));
228       } else PetscCall(PetscObjectSetName((PetscObject)xcoorl, name));
229     }
230 
231     /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
232     PetscCall(PetscObjectCompose((PetscObject)daview, "GLVisGraphicsCoordsGhosted", (PetscObject)xcoorl));
233     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
234 
235     /* daview is composed with the original DMDA */
236     PetscCall(PetscObjectCompose(oda, "GLVisGraphicsDMDAGhosted", (PetscObject)daview));
237     PetscCall(PetscObjectDereference((PetscObject)daview));
238   }
239 
240   /* customize the viewer if present */
241   if (viewer) {
242     DMDAFieldGLVisViewerCtx   *ctx;
243     DMDAGhostedGLVisViewerCtx *dactx;
244     char                       fec[64];
245     Vec                        xlocal, *Ufield;
246     const char               **dafieldname;
247     char                     **fec_type, **fieldname;
248     PetscInt                  *nlocal, *bss, *dims;
249     PetscInt                   dim, M, N, P, dof, s, i, nf;
250     PetscBool                  bsset;
251 
252     PetscCall(DMDAGetInfo(daview, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
253     PetscCall(DMGetApplicationContext(daview, &dactx));
254     PetscCall(DMCreateLocalVector(daview, &xlocal));
255     PetscCall(DMDAGetFieldNames(da, (const char *const **)&dafieldname));
256     PetscCall(DMDAGetNumVerticesGhosted(daview, &M, &N, &P));
257     PetscCall(PetscSNPrintf(fec, sizeof(fec), "FiniteElementCollection: H1_%" PetscInt_FMT "D_P1", dim));
258     PetscCall(PetscMalloc6(dof, &fec_type, dof, &nlocal, dof, &bss, dof, &dims, dof, &fieldname, dof, &Ufield));
259     for (i = 0; i < dof; i++) bss[i] = 1;
260     nf = dof;
261 
262     PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Field options", "PetscViewer");
263     PetscCall(PetscOptionsIntArray("-viewer_glvis_dm_da_bs", "Block sizes for subfields; enable vector representation", NULL, bss, &nf, &bsset));
264     PetscOptionsEnd();
265     if (bsset) {
266       PetscInt t;
267       for (i = 0, t = 0; i < nf; i++) t += bss[i];
268       PetscCheck(t == dof, PetscObjectComm(oda), PETSC_ERR_USER, "Sum of block sizes %" PetscInt_FMT " should equal %" PetscInt_FMT, t, dof);
269     } else nf = dof;
270 
271     for (i = 0, s = 0; i < nf; i++) {
272       PetscCall(PetscStrallocpy(fec, &fec_type[i]));
273       if (bss[i] == 1) {
274         PetscCall(PetscStrallocpy(dafieldname[s], &fieldname[i]));
275       } else {
276         PetscInt b;
277         size_t   tlen = 9; /* "Vector-" + end */
278         for (b = 0; b < bss[i]; b++) {
279           size_t len;
280           PetscCall(PetscStrlen(dafieldname[s + b], &len));
281           tlen += len + 1; /* field + "-" */
282         }
283         PetscCall(PetscMalloc1(tlen, &fieldname[i]));
284         PetscCall(PetscStrcpy(fieldname[i], "Vector-"));
285         for (b = 0; b < bss[i] - 1; b++) {
286           PetscCall(PetscStrcat(fieldname[i], dafieldname[s + b]));
287           PetscCall(PetscStrcat(fieldname[i], "-"));
288         }
289         PetscCall(PetscStrcat(fieldname[i], dafieldname[s + b]));
290       }
291       dims[i]   = dim;
292       nlocal[i] = M * N * P * bss[i];
293       s += bss[i];
294     }
295 
296     /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
297     PetscCall(PetscNew(&ctx));
298     ctx->xlocal = xlocal;
299 
300     /* create work vectors */
301     for (i = 0; i < nf; i++) {
302       PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), nlocal[i], PETSC_DECIDE, &Ufield[i]));
303       PetscCall(PetscObjectSetName((PetscObject)Ufield[i], fieldname[i]));
304       PetscCall(VecSetBlockSize(Ufield[i], bss[i]));
305       PetscCall(VecSetDM(Ufield[i], da));
306     }
307 
308     PetscCall(PetscViewerGLVisSetFields(viewer, nf, (const char **)fec_type, dims, DMDASampleGLVisFields_Private, (PetscObject *)Ufield, ctx, DMDAFieldDestroyGLVisViewerCtx_Private));
309     for (i = 0; i < nf; i++) {
310       PetscCall(PetscFree(fec_type[i]));
311       PetscCall(PetscFree(fieldname[i]));
312       PetscCall(VecDestroy(&Ufield[i]));
313     }
314     PetscCall(PetscFree6(fec_type, nlocal, bss, dims, fieldname, Ufield));
315   }
316   PetscFunctionReturn(0);
317 }
318 
319 static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
320 {
321   DM                 da, cda;
322   Vec                xcoorl;
323   PetscMPIInt        size;
324   const PetscScalar *array;
325   PetscContainer     glvis_container;
326   PetscInt           dim, sdim, i, vid[8], mid, cid, cdof;
327   PetscInt           sx, sy, sz, ie, je, ke, ien, jen, ken, nel;
328   PetscInt           gsx, gsy, gsz, gm, gn, gp, kst, jst, ist;
329   PetscBool          enabled = PETSC_TRUE, isascii;
330   const char        *fmt;
331 
332   PetscFunctionBegin;
333   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
334   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
335   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
336   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer must be of type VIEWERASCII");
337   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
338   PetscCheck(size <= 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Use single sequential viewers for parallel visualization");
339   PetscCall(DMGetDimension(dm, &dim));
340 
341   /* get container: determines if a process visualizes is portion of the data or not */
342   PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
343   PetscCheck(glvis_container, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis container");
344   {
345     PetscViewerGLVisInfo glvis_info;
346     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info));
347     enabled = glvis_info->enabled;
348     fmt     = glvis_info->fmt;
349   }
350   /* this can happen if we are calling DMView outside of VecView_GLVis */
351   PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
352   if (!da) PetscCall(DMSetUpGLVisViewer_DMDA((PetscObject)dm, NULL));
353   PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
354   PetscCheck(da, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted DMDA");
355   PetscCall(DMGetCoordinateDim(da, &sdim));
356 
357   PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh v1.0\n"));
358   PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n"));
359   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim));
360 
361   if (!enabled) {
362     PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
363     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
364     PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
365     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
366     PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
367     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
368     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
369     PetscFunctionReturn(0);
370   }
371 
372   PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
373   nel = ien;
374   if (dim > 1) nel *= jen;
375   if (dim > 2) nel *= ken;
376   PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
377   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", nel));
378   switch (dim) {
379   case 1:
380     for (ie = 0; ie < ien; ie++) {
381       vid[0] = ie;
382       vid[1] = ie + 1;
383       mid    = 1; /* material id */
384       cid    = 1; /* segment */
385       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1]));
386     }
387     break;
388   case 2:
389     for (je = 0; je < jen; je++) {
390       for (ie = 0; ie < ien; ie++) {
391         vid[0] = je * (ien + 1) + ie;
392         vid[1] = je * (ien + 1) + ie + 1;
393         vid[2] = (je + 1) * (ien + 1) + ie + 1;
394         vid[3] = (je + 1) * (ien + 1) + ie;
395         mid    = 1; /* material id */
396         cid    = 3; /* quad */
397         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1], vid[2], vid[3]));
398       }
399     }
400     break;
401   case 3:
402     for (ke = 0; ke < ken; ke++) {
403       for (je = 0; je < jen; je++) {
404         for (ie = 0; ie < ien; ie++) {
405           vid[0] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
406           vid[1] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
407           vid[2] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
408           vid[3] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
409           vid[4] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
410           vid[5] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
411           vid[6] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
412           vid[7] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
413           mid    = 1; /* material id */
414           cid    = 5; /* hex */
415           PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1], vid[2], vid[3], vid[4], vid[5], vid[6], vid[7]));
416         }
417       }
418     }
419     break;
420   default:
421     SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
422   }
423   PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
424   PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
425 
426   /* vertex coordinates */
427   PetscCall(PetscObjectQuery((PetscObject)da, "GLVisGraphicsCoordsGhosted", (PetscObject *)&xcoorl));
428   PetscCheck(xcoorl, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted coords");
429   PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
430   PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
431   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", ien * jen * ken));
432   if (nel) {
433     PetscCall(VecGetDM(xcoorl, &cda));
434     PetscCall(VecGetArrayRead(xcoorl, &array));
435     if (!cda) { /* HO viz */
436       const char *fecname;
437       PetscInt    nc, nl;
438 
439       PetscCall(PetscObjectGetName((PetscObject)xcoorl, &fecname));
440       PetscCall(PetscViewerASCIIPrintf(viewer, "nodes\n"));
441       PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
442       PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", fecname));
443       PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", sdim));
444       PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: 1\n\n")); /*Ordering::byVDIM*/
445       /* L2 coordinates */
446       PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
447       PetscCall(VecGetLocalSize(xcoorl, &nl));
448       nc   = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
449       cdof = nc ? nl / nc : 0;
450       if (!ien) ien++;
451       if (!jen) jen++;
452       if (!ken) ken++;
453       ist = jst = kst = 0;
454       gm              = ien;
455       gn              = jen;
456       gp              = ken;
457     } else {
458       DMDAGhostedGLVisViewerCtx *dactx;
459 
460       PetscCall(DMGetApplicationContext(da, &dactx));
461       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
462       cdof = sdim;
463       PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
464       PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
465       if (dactx->ll) {
466         kst = jst = ist = 0;
467       } else {
468         kst = gsz != sz ? 1 : 0;
469         jst = gsy != sy ? 1 : 0;
470         ist = gsx != sx ? 1 : 0;
471       }
472     }
473     for (ke = kst; ke < kst + ken; ke++) {
474       for (je = jst; je < jst + jen; je++) {
475         for (ie = ist; ie < ist + ien; ie++) {
476           PetscInt c;
477 
478           i = ke * gm * gn + je * gm + ie;
479           for (c = 0; c < cdof / sdim; c++) {
480             PetscInt d;
481             for (d = 0; d < sdim; d++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, PetscRealPart(array[cdof * i + c * sdim + d])));
482             PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
483           }
484         }
485       }
486     }
487     PetscCall(VecRestoreArrayRead(xcoorl, &array));
488   }
489   PetscFunctionReturn(0);
490 }
491 
492 PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
493 {
494   PetscFunctionBegin;
495   PetscCall(DMView_GLVis(dm, viewer, DMDAView_GLVis_ASCII));
496   PetscFunctionReturn(0);
497 }
498