xref: /petsc/src/dm/impls/da/grglvis.c (revision 8135c375e96a7bf70d67c02f02cd826154817056)
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   Vec xlocal;
8 } DMDAGLVisViewerCtx;
9 
10 static PetscErrorCode DMDADestroyGLVisViewerCtx_Private(void *vctx)
11 {
12   DMDAGLVisViewerCtx *ctx = (DMDAGLVisViewerCtx*)vctx;
13   PetscErrorCode     ierr;
14 
15   PetscFunctionBegin;
16   ierr = VecDestroy(&ctx->xlocal);CHKERRQ(ierr);
17   ierr = PetscFree(vctx);
18   PetscFunctionReturn(0);
19 }
20 
21 /* all but the last proc per dimension claim the ghosted node */
22 static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
23 {
24   PetscInt       M,N,P,sx,sy,sz,ien,jen,ken;
25   PetscErrorCode ierr;
26 
27   PetscFunctionBegin;
28   ierr = DMDAGetInfo(da,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
29   ierr = DMDAGetCorners(da,&sx,&sy,&sz,&ien,&jen,&ken);CHKERRQ(ierr);
30   if (sx+ien == M) ien--;
31   if (sy+jen == N) jen--;
32   if (sz+ken == P) ken--;
33   if (nex) *nex = ien;
34   if (ney) *ney = jen;
35   if (nez) *nez = ken;
36   PetscFunctionReturn(0);
37 }
38 
39 /* inherits number of vertices from DMDAGetNumElementsGhosted */
40 static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
41 {
42   PetscInt       ien,jen,ken,dim;
43   PetscErrorCode ierr;
44 
45   PetscFunctionBegin;
46   ierr = DMGetDimension(da,&dim);CHKERRQ(ierr);
47   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
48   ien  = ien+1;
49   jen  = dim > 1 ? jen+1 : 1;
50   ken  = dim > 2 ? ken+1 : 1;
51   if (nvx) *nvx = ien;
52   if (nvy) *nvy = jen;
53   if (nvz) *nvz = ken;
54   PetscFunctionReturn(0);
55 }
56 
57 static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
58 {
59   DM                 da;
60   DMDAGLVisViewerCtx *ctx = (DMDAGLVisViewerCtx*)vctx;
61   const PetscScalar  *array;
62   PetscScalar        **arrayf;
63   PetscInt           i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
64   PetscInt           sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
65   PetscErrorCode     ierr;
66 
67   PetscFunctionBegin;
68   ierr = VecGetDM(ctx->xlocal,&da);CHKERRQ(ierr);
69   if (!da) SETERRQ(PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
70   ierr = VecGetBlockSize(ctx->xlocal,&bs);CHKERRQ(ierr);
71   ierr = DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr);
72   ierr = DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr);
73   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
74   ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
75   ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
76   kst  = gsz != sz ? 1 : 0;
77   jst  = gsy != sy ? 1 : 0;
78   ist  = gsx != sx ? 1 : 0;
79   ierr = PetscMalloc2(nf,&arrayf,nf,&bss);CHKERRQ(ierr);
80   ierr = VecGetArrayRead(ctx->xlocal,&array);CHKERRQ(ierr);
81   for (f=0;f<nf;f++) {
82     ierr = VecGetBlockSize((Vec)oXf[f],&bss[f]);CHKERRQ(ierr);
83     ierr = VecGetArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr);
84   }
85   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
86     for (je = jst; je < jst + jen; je++) {
87       for (ie = ist; ie < ist + ien; ie++) {
88         PetscInt cf,b;
89         i = ke * gm * gn + je * gm + ie;
90         for (f=0,cf=0;f<nf;f++)
91           for (b=0;b<bss[f];b++)
92             arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
93         ii++;
94       }
95     }
96   }
97   for (f=0;f<nf;f++) { ierr = VecRestoreArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr); }
98   ierr = VecRestoreArrayRead(ctx->xlocal,&array);CHKERRQ(ierr);
99   ierr = PetscFree2(arrayf,bss);CHKERRQ(ierr);
100   PetscFunctionReturn(0);
101 }
102 
103 PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
104 {
105   DM                 da = (DM)oda,daview,dacoord;
106   Vec                xcoor,xcoorl,xlocal;
107   DMDAGLVisViewerCtx *ctx;
108   const char         **dafieldname;
109   char               **fec_type,**fieldname,fec[64];
110   const PetscInt     *lx,*ly,*lz;
111   PetscInt           *nlocal,*bss;
112   PetscInt           dim,M,N,P,m,n,p,dof,s,i,nf;
113   PetscBool          bsset;
114   PetscErrorCode     ierr;
115 
116   PetscFunctionBegin;
117   ierr = PetscObjectQuery(oda,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
118   if (xcoorl) PetscFunctionReturn(0);
119   /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
120   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
121   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
122   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
123   ierr = PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");CHKERRQ(ierr);
124   switch (dim) {
125   case 1:
126     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);CHKERRQ(ierr);
127     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);CHKERRQ(ierr);
128     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_1D_P1");CHKERRQ(ierr);
129     break;
130   case 2:
131     ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,lx,ly,&daview);CHKERRQ(ierr);
132     ierr = DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,2,1,lx,ly,&dacoord);CHKERRQ(ierr);
133     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_2D_P1");CHKERRQ(ierr);
134     break;
135   case 3:
136     ierr = 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);CHKERRQ(ierr);
137     ierr = 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);CHKERRQ(ierr);
138     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_3D_P1");CHKERRQ(ierr);
139     break;
140   default:
141     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
142     break;
143   }
144   ierr = DMSetUp(daview);CHKERRQ(ierr);
145   ierr = DMSetUp(dacoord);CHKERRQ(ierr);
146   if (!xcoor) {
147     ierr = DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
148     ierr = DMGetCoordinates(daview,&xcoor);CHKERRQ(ierr);
149   }
150   ierr = DMCreateLocalVector(daview,&xlocal);CHKERRQ(ierr);
151   ierr = DMCreateLocalVector(dacoord,&xcoorl);CHKERRQ(ierr);
152   ierr = DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
153   ierr = DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
154 
155   /* xcoorl is composed with the original DMDA, ghosted coordinate DMDA is only available through this vector */
156   ierr = PetscObjectCompose(oda,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
157   ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
158 
159   /* customize the viewer */
160   ierr = DMDAGetFieldNames(da,(const char * const **)&dafieldname);CHKERRQ(ierr);
161   ierr = DMDAGetNumVerticesGhosted(daview,&M,&N,&P);CHKERRQ(ierr);
162   ierr = PetscMalloc4(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&fieldname);CHKERRQ(ierr);
163   for (i=0;i<dof;i++) bss[i] = 1;
164   nf = dof;
165 
166   ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");CHKERRQ(ierr);
167   ierr = PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);CHKERRQ(ierr);
168   ierr = PetscOptionsEnd();CHKERRQ(ierr);
169   if (bsset) {
170     PetscInt t;
171     for (i=0,t=0;i<nf;i++) t += bss[i];
172     if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
173   } else nf = dof;
174 
175   for (i=0,s=0;i<nf;i++) {
176     ierr = PetscStrallocpy(fec,&fec_type[i]);CHKERRQ(ierr);
177     if (bss[i] == 1) {
178       ierr = PetscStrallocpy(dafieldname[s],&fieldname[i]);CHKERRQ(ierr);
179     } else {
180       PetscInt b;
181       size_t tlen = 9; /* "Vector-" + end */
182       for (b=0;b<bss[i];b++) {
183         size_t len;
184         ierr = PetscStrlen(dafieldname[s+b],&len);CHKERRQ(ierr);
185         tlen += len + 1; /* field + "-" */
186       }
187       ierr = PetscMalloc1(tlen,&fieldname[i]);CHKERRQ(ierr);
188       ierr = PetscStrcpy(fieldname[i],"Vector-");CHKERRQ(ierr);
189       for (b=0;b<bss[i]-1;b++) {
190         ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
191         ierr = PetscStrcat(fieldname[i],"-");CHKERRQ(ierr);
192       }
193       ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
194     }
195     nlocal[i] = M*N*P*bss[i];
196     s += bss[i];
197   }
198 
199   /* the viewer context takes ownership of xlocal (and the the properly ghosted DMDA associated with it) */
200   ierr = PetscNew(&ctx);CHKERRQ(ierr);
201   ctx->xlocal = xlocal;
202 
203   ierr = PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,DMDASampleGLVisFields_Private,ctx,DMDADestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
204   for (i=0;i<nf;i++) {
205     ierr = PetscFree(fec_type[i]);CHKERRQ(ierr);
206     ierr = PetscFree(fieldname[i]);CHKERRQ(ierr);
207   }
208   ierr = PetscFree4(fec_type,nlocal,bss,fieldname);CHKERRQ(ierr);
209   ierr = DMDestroy(&dacoord);CHKERRQ(ierr);
210   ierr = DMDestroy(&daview);CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
215 {
216   DM                da;
217   Vec               xcoorl;
218   PetscMPIInt       commsize;
219   const PetscScalar *array;
220   PetscContainer    glvis_container;
221   PetscInt          dim,sdim,i,vid[8],mid,cid;
222   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken;
223   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
224   PetscBool         enabled = PETSC_TRUE, isascii;
225   PetscErrorCode    ierr;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
229   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
230   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
231   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
232   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
233   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
234   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
235 
236   /* get container: determines if a process visualizes is portion of the data or not */
237   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
238   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
239   {
240     PetscViewerGLVisInfo glvis_info;
241     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
242     enabled = glvis_info->enabled;
243   }
244   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
245   if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
246   ierr = VecGetDM(xcoorl,&da);CHKERRQ(ierr);
247   if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
248   ierr = DMGetCoordinateDim(da,&sdim);CHKERRQ(ierr);
249 
250   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
251   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
252   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
253 
254   if (!enabled) {
255     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
256     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
257     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
258     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
259     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
260     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
261     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
262     PetscFunctionReturn(0);
263   }
264 
265   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
266   i    = ien;
267   if (dim > 1) i *= jen;
268   if (dim > 2) i *= ken;
269   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
270   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",i);CHKERRQ(ierr);
271   switch (dim) {
272   case 1:
273     for (ie = 0; ie < ien; ie++) {
274       vid[0] = ie;
275       vid[1] = ie+1;
276       mid    = 1; /* material id */
277       cid    = 1; /* segment */
278       ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);CHKERRQ(ierr);
279     }
280     break;
281   case 2:
282     for (je = 0; je < jen; je++) {
283       for (ie = 0; ie < ien; ie++) {
284         vid[0] =     je*(ien+1) + ie;
285         vid[1] =     je*(ien+1) + ie+1;
286         vid[2] = (je+1)*(ien+1) + ie+1;
287         vid[3] = (je+1)*(ien+1) + ie;
288         mid    = 1; /* material id */
289         cid    = 3; /* quad */
290         ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);CHKERRQ(ierr);
291       }
292     }
293     break;
294   case 3:
295     for (ke = 0; ke < ken; ke++) {
296       for (je = 0; je < jen; je++) {
297         for (ie = 0; ie < ien; ie++) {
298           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
299           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
300           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
301           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
302           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
303           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
304           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
305           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
306           mid    = 1; /* material id */
307           cid    = 5; /* hex */
308           ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3],vid[4],vid[5],vid[6],vid[7]);CHKERRQ(ierr);
309         }
310       }
311     }
312     break;
313   default:
314     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
315     break;
316   }
317   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
318   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
319 
320   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
321   ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
322   ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
323   kst  = gsz != sz ? 1 : 0;
324   jst  = gsy != sy ? 1 : 0;
325   ist  = gsx != sx ? 1 : 0;
326   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
327   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);CHKERRQ(ierr);
328   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
329   ierr = VecGetArrayRead(xcoorl,&array);CHKERRQ(ierr);
330   for (ke = kst; ke < kst + ken; ke++) {
331     for (je = jst; je < jst + jen; je++) {
332       for (ie = ist; ie < ist + ien; ie++) {
333         PetscInt d;
334 
335         i = ke * gm * gn + je * gm + ie;
336         for (d=0;d<sdim-1;d++) {
337           ierr = PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[sdim*i+d]));CHKERRQ(ierr);
338         }
339         ierr = PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[sdim*i+d]));CHKERRQ(ierr);
340       }
341     }
342   }
343   ierr = VecRestoreArrayRead(xcoorl,&array);CHKERRQ(ierr);
344   PetscFunctionReturn(0);
345 }
346 
347 /* dispatching, prints through the socket by prepending the mesh keyword to the usual ASCII dump: duplicated code as in plexglvis.c, should be merged together */
348 PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
349 {
350   PetscErrorCode ierr;
351   PetscBool      isglvis,isascii;
352 
353   PetscFunctionBegin;
354   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
355   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
356   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
357   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
358   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
359   if (isglvis) {
360     PetscViewer          view;
361     PetscViewerGLVisType type;
362 
363     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
364     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
365     if (view) { /* in the socket case, it may happen that the connection failed */
366       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
367         PetscMPIInt size,rank;
368         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
369         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
370         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
371       }
372       ierr = DMDAView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
373       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
374       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
375         PetscInt    dim;
376         const char* name;
377 
378         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
379         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
380         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
381         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
382       }
383     }
384   } else {
385     ierr = DMDAView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
386   }
387   PetscFunctionReturn(0);
388 }
389