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