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