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