xref: /petsc/src/dm/impls/da/grglvis.c (revision 1ef7c9b6d69c78c0e1357aee83c8affecef1b942)
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);CHKERRQ(ierr);
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,*dims;
112   PetscInt           dim,M,N,P,m,n,p,dof,s,i,nf;
113   PetscBool          bsset;
114   PetscErrorCode     ierr;
115 
116   PetscFunctionBegin;
117   /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
118   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
119   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
120   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
121   ierr = PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");CHKERRQ(ierr);
122   switch (dim) {
123   case 1:
124     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);CHKERRQ(ierr);
125     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);CHKERRQ(ierr);
126     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_1D_P1");CHKERRQ(ierr);
127     break;
128   case 2:
129     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);
130     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);
131     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_2D_P1");CHKERRQ(ierr);
132     break;
133   case 3:
134     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);
135     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);
136     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_3D_P1");CHKERRQ(ierr);
137     break;
138   default:
139     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
140     break;
141   }
142   ierr = DMSetUp(daview);CHKERRQ(ierr);
143   ierr = DMSetUp(dacoord);CHKERRQ(ierr);
144   if (!xcoor) {
145     ierr = DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
146     ierr = DMGetCoordinates(daview,&xcoor);CHKERRQ(ierr);
147   }
148   ierr = DMCreateLocalVector(daview,&xlocal);CHKERRQ(ierr);
149   ierr = DMCreateLocalVector(dacoord,&xcoorl);CHKERRQ(ierr);
150   ierr = DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
151   ierr = DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
152 
153   /* xcoorl is composed with the original DMDA, ghosted coordinate DMDA is only available through this vector */
154   ierr = PetscObjectCompose(oda,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
155   ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
156 
157   /* customize the viewer */
158   ierr = DMDAGetFieldNames(da,(const char * const **)&dafieldname);CHKERRQ(ierr);
159   ierr = DMDAGetNumVerticesGhosted(daview,&M,&N,&P);CHKERRQ(ierr);
160   ierr = PetscMalloc5(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname);CHKERRQ(ierr);
161   for (i=0;i<dof;i++) bss[i] = 1;
162   nf = dof;
163 
164   ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");CHKERRQ(ierr);
165   ierr = PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);CHKERRQ(ierr);
166   ierr = PetscOptionsEnd();CHKERRQ(ierr);
167   if (bsset) {
168     PetscInt t;
169     for (i=0,t=0;i<nf;i++) t += bss[i];
170     if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
171   } else nf = dof;
172 
173   for (i=0,s=0;i<nf;i++) {
174     ierr = PetscStrallocpy(fec,&fec_type[i]);CHKERRQ(ierr);
175     if (bss[i] == 1) {
176       ierr = PetscStrallocpy(dafieldname[s],&fieldname[i]);CHKERRQ(ierr);
177     } else {
178       PetscInt b;
179       size_t tlen = 9; /* "Vector-" + end */
180       for (b=0;b<bss[i];b++) {
181         size_t len;
182         ierr = PetscStrlen(dafieldname[s+b],&len);CHKERRQ(ierr);
183         tlen += len + 1; /* field + "-" */
184       }
185       ierr = PetscMalloc1(tlen,&fieldname[i]);CHKERRQ(ierr);
186       ierr = PetscStrcpy(fieldname[i],"Vector-");CHKERRQ(ierr);
187       for (b=0;b<bss[i]-1;b++) {
188         ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
189         ierr = PetscStrcat(fieldname[i],"-");CHKERRQ(ierr);
190       }
191       ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
192     }
193     dims[i] = dim;
194     nlocal[i] = M*N*P*bss[i];
195     s += bss[i];
196   }
197 
198   /* the viewer context takes ownership of xlocal (and the the properly ghosted DMDA associated with it) */
199   ierr = PetscNew(&ctx);CHKERRQ(ierr);
200   ctx->xlocal = xlocal;
201 
202   ierr = PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,dims,DMDASampleGLVisFields_Private,ctx,DMDADestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
203   for (i=0;i<nf;i++) {
204     ierr = PetscFree(fec_type[i]);CHKERRQ(ierr);
205     ierr = PetscFree(fieldname[i]);CHKERRQ(ierr);
206   }
207   ierr = PetscFree5(fec_type,nlocal,bss,dims,fieldname);CHKERRQ(ierr);
208   ierr = DMDestroy(&dacoord);CHKERRQ(ierr);
209   ierr = DMDestroy(&daview);CHKERRQ(ierr);
210   PetscFunctionReturn(0);
211 }
212 
213 static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
214 {
215   DM                da;
216   Vec               xcoorl;
217   PetscMPIInt       commsize;
218   const PetscScalar *array;
219   PetscContainer    glvis_container;
220   PetscInt          dim,sdim,i,vid[8],mid,cid;
221   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken;
222   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
223   PetscBool         enabled = PETSC_TRUE, isascii;
224   PetscErrorCode    ierr;
225 
226   PetscFunctionBegin;
227   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
228   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
229   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
230   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
231   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
232   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
233   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
234 
235   /* get container: determines if a process visualizes is portion of the data or not */
236   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
237   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
238   {
239     PetscViewerGLVisInfo glvis_info;
240     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
241     enabled = glvis_info->enabled;
242   }
243   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
244   if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
245   ierr = VecGetDM(xcoorl,&da);CHKERRQ(ierr);
246   if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
247   ierr = DMGetCoordinateDim(da,&sdim);CHKERRQ(ierr);
248 
249   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
250   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
251   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
252 
253   if (!enabled) {
254     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
255     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
256     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
257     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
258     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
259     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
260     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
261     PetscFunctionReturn(0);
262   }
263 
264   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
265   i    = ien;
266   if (dim > 1) i *= jen;
267   if (dim > 2) i *= ken;
268   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
269   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",i);CHKERRQ(ierr);
270   switch (dim) {
271   case 1:
272     for (ie = 0; ie < ien; ie++) {
273       vid[0] = ie;
274       vid[1] = ie+1;
275       mid    = 1; /* material id */
276       cid    = 1; /* segment */
277       ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);CHKERRQ(ierr);
278     }
279     break;
280   case 2:
281     for (je = 0; je < jen; je++) {
282       for (ie = 0; ie < ien; ie++) {
283         vid[0] =     je*(ien+1) + ie;
284         vid[1] =     je*(ien+1) + ie+1;
285         vid[2] = (je+1)*(ien+1) + ie+1;
286         vid[3] = (je+1)*(ien+1) + ie;
287         mid    = 1; /* material id */
288         cid    = 3; /* quad */
289         ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);CHKERRQ(ierr);
290       }
291     }
292     break;
293   case 3:
294     for (ke = 0; ke < ken; ke++) {
295       for (je = 0; je < jen; je++) {
296         for (ie = 0; ie < ien; ie++) {
297           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
298           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
299           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
300           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
301           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
302           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
303           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
304           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
305           mid    = 1; /* material id */
306           cid    = 5; /* hex */
307           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);
308         }
309       }
310     }
311     break;
312   default:
313     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
314     break;
315   }
316   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
317   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
318 
319   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
320   ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
321   ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
322   kst  = gsz != sz ? 1 : 0;
323   jst  = gsy != sy ? 1 : 0;
324   ist  = gsx != sx ? 1 : 0;
325   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
326   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);CHKERRQ(ierr);
327   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
328   ierr = VecGetArrayRead(xcoorl,&array);CHKERRQ(ierr);
329   for (ke = kst; ke < kst + ken; ke++) {
330     for (je = jst; je < jst + jen; je++) {
331       for (ie = ist; ie < ist + ien; ie++) {
332         PetscInt d;
333 
334         i = ke * gm * gn + je * gm + ie;
335         for (d=0;d<sdim-1;d++) {
336           ierr = PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[sdim*i+d]));CHKERRQ(ierr);
337         }
338         ierr = PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[sdim*i+d]));CHKERRQ(ierr);
339       }
340     }
341   }
342   ierr = VecRestoreArrayRead(xcoorl,&array);CHKERRQ(ierr);
343   PetscFunctionReturn(0);
344 }
345 
346 /* 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 */
347 PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
348 {
349   PetscErrorCode ierr;
350   PetscBool      isglvis,isascii;
351 
352   PetscFunctionBegin;
353   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
354   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
355   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
356   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
357   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
358   if (isglvis) {
359     PetscViewer          view;
360     PetscViewerGLVisType type;
361 
362     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
363     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
364     if (view) { /* in the socket case, it may happen that the connection failed */
365       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
366         PetscMPIInt size,rank;
367         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
368         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
369         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
370       }
371       ierr = DMDAView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
372       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
373       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
374         PetscInt    dim;
375         const char* name;
376 
377         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
378         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
379         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
380         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
381       }
382     }
383   } else {
384     ierr = DMDAView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
385   }
386   PetscFunctionReturn(0);
387 }
388