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