xref: /petsc/src/dm/impls/da/grglvis.c (revision aab5bcd8e18f0808071e10849b7b64c13bdb175c)
18135c375SStefano Zampini /* Routines to visualize DMDAs and fields through GLVis */
28135c375SStefano Zampini 
38135c375SStefano Zampini #include <petsc/private/dmdaimpl.h>
48135c375SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
58135c375SStefano Zampini 
68135c375SStefano Zampini typedef struct {
78135c375SStefano Zampini   Vec xlocal;
88135c375SStefano Zampini } DMDAGLVisViewerCtx;
98135c375SStefano Zampini 
108135c375SStefano Zampini static PetscErrorCode DMDADestroyGLVisViewerCtx_Private(void *vctx)
118135c375SStefano Zampini {
128135c375SStefano Zampini   DMDAGLVisViewerCtx *ctx = (DMDAGLVisViewerCtx*)vctx;
138135c375SStefano Zampini   PetscErrorCode     ierr;
148135c375SStefano Zampini 
158135c375SStefano Zampini   PetscFunctionBegin;
168135c375SStefano Zampini   ierr = VecDestroy(&ctx->xlocal);CHKERRQ(ierr);
1751f51421SSatish Balay   ierr = PetscFree(vctx);CHKERRQ(ierr);
188135c375SStefano Zampini   PetscFunctionReturn(0);
198135c375SStefano Zampini }
208135c375SStefano Zampini 
218135c375SStefano Zampini /* all but the last proc per dimension claim the ghosted node */
228135c375SStefano Zampini static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
238135c375SStefano Zampini {
248135c375SStefano Zampini   PetscInt       M,N,P,sx,sy,sz,ien,jen,ken;
258135c375SStefano Zampini   PetscErrorCode ierr;
268135c375SStefano Zampini 
278135c375SStefano Zampini   PetscFunctionBegin;
28*aab5bcd8SJed Brown   /* Appease -Wmaybe-uninitialized */
29*aab5bcd8SJed Brown   if (nex) *nex = -1;
30*aab5bcd8SJed Brown   if (ney) *ney = -1;
31*aab5bcd8SJed Brown   if (nez) *nez = -1;
328135c375SStefano Zampini   ierr = DMDAGetInfo(da,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
338135c375SStefano Zampini   ierr = DMDAGetCorners(da,&sx,&sy,&sz,&ien,&jen,&ken);CHKERRQ(ierr);
348135c375SStefano Zampini   if (sx+ien == M) ien--;
358135c375SStefano Zampini   if (sy+jen == N) jen--;
368135c375SStefano Zampini   if (sz+ken == P) ken--;
378135c375SStefano Zampini   if (nex) *nex = ien;
388135c375SStefano Zampini   if (ney) *ney = jen;
398135c375SStefano Zampini   if (nez) *nez = ken;
408135c375SStefano Zampini   PetscFunctionReturn(0);
418135c375SStefano Zampini }
428135c375SStefano Zampini 
438135c375SStefano Zampini /* inherits number of vertices from DMDAGetNumElementsGhosted */
448135c375SStefano Zampini static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
458135c375SStefano Zampini {
468135c375SStefano Zampini   PetscInt       ien,jen,ken,dim;
478135c375SStefano Zampini   PetscErrorCode ierr;
488135c375SStefano Zampini 
498135c375SStefano Zampini   PetscFunctionBegin;
508135c375SStefano Zampini   ierr = DMGetDimension(da,&dim);CHKERRQ(ierr);
518135c375SStefano Zampini   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
528135c375SStefano Zampini   ien  = ien+1;
538135c375SStefano Zampini   jen  = dim > 1 ? jen+1 : 1;
548135c375SStefano Zampini   ken  = dim > 2 ? ken+1 : 1;
558135c375SStefano Zampini   if (nvx) *nvx = ien;
568135c375SStefano Zampini   if (nvy) *nvy = jen;
578135c375SStefano Zampini   if (nvz) *nvz = ken;
588135c375SStefano Zampini   PetscFunctionReturn(0);
598135c375SStefano Zampini }
608135c375SStefano Zampini 
618135c375SStefano Zampini static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
628135c375SStefano Zampini {
638135c375SStefano Zampini   DM                 da;
648135c375SStefano Zampini   DMDAGLVisViewerCtx *ctx = (DMDAGLVisViewerCtx*)vctx;
658135c375SStefano Zampini   const PetscScalar  *array;
668135c375SStefano Zampini   PetscScalar        **arrayf;
678135c375SStefano Zampini   PetscInt           i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
688135c375SStefano Zampini   PetscInt           sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
698135c375SStefano Zampini   PetscErrorCode     ierr;
708135c375SStefano Zampini 
718135c375SStefano Zampini   PetscFunctionBegin;
728135c375SStefano Zampini   ierr = VecGetDM(ctx->xlocal,&da);CHKERRQ(ierr);
738135c375SStefano Zampini   if (!da) SETERRQ(PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
748135c375SStefano Zampini   ierr = VecGetBlockSize(ctx->xlocal,&bs);CHKERRQ(ierr);
758135c375SStefano Zampini   ierr = DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr);
768135c375SStefano Zampini   ierr = DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr);
778135c375SStefano Zampini   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
788135c375SStefano Zampini   ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
798135c375SStefano Zampini   ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
808135c375SStefano Zampini   kst  = gsz != sz ? 1 : 0;
818135c375SStefano Zampini   jst  = gsy != sy ? 1 : 0;
828135c375SStefano Zampini   ist  = gsx != sx ? 1 : 0;
838135c375SStefano Zampini   ierr = PetscMalloc2(nf,&arrayf,nf,&bss);CHKERRQ(ierr);
848135c375SStefano Zampini   ierr = VecGetArrayRead(ctx->xlocal,&array);CHKERRQ(ierr);
858135c375SStefano Zampini   for (f=0;f<nf;f++) {
868135c375SStefano Zampini     ierr = VecGetBlockSize((Vec)oXf[f],&bss[f]);CHKERRQ(ierr);
878135c375SStefano Zampini     ierr = VecGetArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr);
888135c375SStefano Zampini   }
898135c375SStefano Zampini   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
908135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
918135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
928135c375SStefano Zampini         PetscInt cf,b;
938135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
948135c375SStefano Zampini         for (f=0,cf=0;f<nf;f++)
958135c375SStefano Zampini           for (b=0;b<bss[f];b++)
968135c375SStefano Zampini             arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
978135c375SStefano Zampini         ii++;
988135c375SStefano Zampini       }
998135c375SStefano Zampini     }
1008135c375SStefano Zampini   }
1018135c375SStefano Zampini   for (f=0;f<nf;f++) { ierr = VecRestoreArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr); }
1028135c375SStefano Zampini   ierr = VecRestoreArrayRead(ctx->xlocal,&array);CHKERRQ(ierr);
1038135c375SStefano Zampini   ierr = PetscFree2(arrayf,bss);CHKERRQ(ierr);
1048135c375SStefano Zampini   PetscFunctionReturn(0);
1058135c375SStefano Zampini }
1068135c375SStefano Zampini 
1078135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
1088135c375SStefano Zampini {
1098135c375SStefano Zampini   DM                 da = (DM)oda,daview,dacoord;
1108135c375SStefano Zampini   Vec                xcoor,xcoorl,xlocal;
1118135c375SStefano Zampini   DMDAGLVisViewerCtx *ctx;
1128135c375SStefano Zampini   const char         **dafieldname;
1138135c375SStefano Zampini   char               **fec_type,**fieldname,fec[64];
1148135c375SStefano Zampini   const PetscInt     *lx,*ly,*lz;
115bb77a09fSStefano Zampini   PetscInt           *nlocal,*bss,*dims;
1168135c375SStefano Zampini   PetscInt           dim,M,N,P,m,n,p,dof,s,i,nf;
1178135c375SStefano Zampini   PetscBool          bsset;
1188135c375SStefano Zampini   PetscErrorCode     ierr;
1198135c375SStefano Zampini 
1208135c375SStefano Zampini   PetscFunctionBegin;
1218135c375SStefano Zampini   ierr = PetscObjectQuery(oda,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
1228135c375SStefano Zampini   if (xcoorl) PetscFunctionReturn(0);
1238135c375SStefano Zampini   /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1248135c375SStefano Zampini   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1258135c375SStefano Zampini   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
1268135c375SStefano Zampini   ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
1278135c375SStefano Zampini   ierr = PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");CHKERRQ(ierr);
1288135c375SStefano Zampini   switch (dim) {
1298135c375SStefano Zampini   case 1:
1308135c375SStefano Zampini     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);CHKERRQ(ierr);
1318135c375SStefano Zampini     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);CHKERRQ(ierr);
1328135c375SStefano Zampini     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_1D_P1");CHKERRQ(ierr);
1338135c375SStefano Zampini     break;
1348135c375SStefano Zampini   case 2:
1358135c375SStefano Zampini     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);
1368135c375SStefano Zampini     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);
1378135c375SStefano Zampini     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_2D_P1");CHKERRQ(ierr);
1388135c375SStefano Zampini     break;
1398135c375SStefano Zampini   case 3:
1408135c375SStefano Zampini     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);
1418135c375SStefano Zampini     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);
1428135c375SStefano Zampini     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_3D_P1");CHKERRQ(ierr);
1438135c375SStefano Zampini     break;
1448135c375SStefano Zampini   default:
1458135c375SStefano Zampini     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
1468135c375SStefano Zampini     break;
1478135c375SStefano Zampini   }
1488135c375SStefano Zampini   ierr = DMSetUp(daview);CHKERRQ(ierr);
1498135c375SStefano Zampini   ierr = DMSetUp(dacoord);CHKERRQ(ierr);
1508135c375SStefano Zampini   if (!xcoor) {
1518135c375SStefano Zampini     ierr = DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
1528135c375SStefano Zampini     ierr = DMGetCoordinates(daview,&xcoor);CHKERRQ(ierr);
1538135c375SStefano Zampini   }
1548135c375SStefano Zampini   ierr = DMCreateLocalVector(daview,&xlocal);CHKERRQ(ierr);
1558135c375SStefano Zampini   ierr = DMCreateLocalVector(dacoord,&xcoorl);CHKERRQ(ierr);
1568135c375SStefano Zampini   ierr = DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
1578135c375SStefano Zampini   ierr = DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
1588135c375SStefano Zampini 
1598135c375SStefano Zampini   /* xcoorl is composed with the original DMDA, ghosted coordinate DMDA is only available through this vector */
1608135c375SStefano Zampini   ierr = PetscObjectCompose(oda,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
1618135c375SStefano Zampini   ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
1628135c375SStefano Zampini 
1638135c375SStefano Zampini   /* customize the viewer */
1648135c375SStefano Zampini   ierr = DMDAGetFieldNames(da,(const char * const **)&dafieldname);CHKERRQ(ierr);
1658135c375SStefano Zampini   ierr = DMDAGetNumVerticesGhosted(daview,&M,&N,&P);CHKERRQ(ierr);
166bb77a09fSStefano Zampini   ierr = PetscMalloc5(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname);CHKERRQ(ierr);
1678135c375SStefano Zampini   for (i=0;i<dof;i++) bss[i] = 1;
1688135c375SStefano Zampini   nf = dof;
1698135c375SStefano Zampini 
1708135c375SStefano Zampini   ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");CHKERRQ(ierr);
1718135c375SStefano Zampini   ierr = PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);CHKERRQ(ierr);
1728135c375SStefano Zampini   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1738135c375SStefano Zampini   if (bsset) {
1748135c375SStefano Zampini     PetscInt t;
1758135c375SStefano Zampini     for (i=0,t=0;i<nf;i++) t += bss[i];
1768135c375SStefano Zampini     if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
1778135c375SStefano Zampini   } else nf = dof;
1788135c375SStefano Zampini 
1798135c375SStefano Zampini   for (i=0,s=0;i<nf;i++) {
1808135c375SStefano Zampini     ierr = PetscStrallocpy(fec,&fec_type[i]);CHKERRQ(ierr);
1818135c375SStefano Zampini     if (bss[i] == 1) {
1828135c375SStefano Zampini       ierr = PetscStrallocpy(dafieldname[s],&fieldname[i]);CHKERRQ(ierr);
1838135c375SStefano Zampini     } else {
1848135c375SStefano Zampini       PetscInt b;
1858135c375SStefano Zampini       size_t tlen = 9; /* "Vector-" + end */
1868135c375SStefano Zampini       for (b=0;b<bss[i];b++) {
1878135c375SStefano Zampini         size_t len;
1888135c375SStefano Zampini         ierr = PetscStrlen(dafieldname[s+b],&len);CHKERRQ(ierr);
1898135c375SStefano Zampini         tlen += len + 1; /* field + "-" */
1908135c375SStefano Zampini       }
1918135c375SStefano Zampini       ierr = PetscMalloc1(tlen,&fieldname[i]);CHKERRQ(ierr);
1928135c375SStefano Zampini       ierr = PetscStrcpy(fieldname[i],"Vector-");CHKERRQ(ierr);
1938135c375SStefano Zampini       for (b=0;b<bss[i]-1;b++) {
1948135c375SStefano Zampini         ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
1958135c375SStefano Zampini         ierr = PetscStrcat(fieldname[i],"-");CHKERRQ(ierr);
1968135c375SStefano Zampini       }
1978135c375SStefano Zampini       ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
1988135c375SStefano Zampini     }
199bb77a09fSStefano Zampini     dims[i] = dim;
2008135c375SStefano Zampini     nlocal[i] = M*N*P*bss[i];
2018135c375SStefano Zampini     s += bss[i];
2028135c375SStefano Zampini   }
2038135c375SStefano Zampini 
2048135c375SStefano Zampini   /* the viewer context takes ownership of xlocal (and the the properly ghosted DMDA associated with it) */
2058135c375SStefano Zampini   ierr = PetscNew(&ctx);CHKERRQ(ierr);
2068135c375SStefano Zampini   ctx->xlocal = xlocal;
2078135c375SStefano Zampini 
208bb77a09fSStefano Zampini   ierr = PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,dims,DMDASampleGLVisFields_Private,ctx,DMDADestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
2098135c375SStefano Zampini   for (i=0;i<nf;i++) {
2108135c375SStefano Zampini     ierr = PetscFree(fec_type[i]);CHKERRQ(ierr);
2118135c375SStefano Zampini     ierr = PetscFree(fieldname[i]);CHKERRQ(ierr);
2128135c375SStefano Zampini   }
213bb77a09fSStefano Zampini   ierr = PetscFree5(fec_type,nlocal,bss,dims,fieldname);CHKERRQ(ierr);
2148135c375SStefano Zampini   ierr = DMDestroy(&dacoord);CHKERRQ(ierr);
2158135c375SStefano Zampini   ierr = DMDestroy(&daview);CHKERRQ(ierr);
2168135c375SStefano Zampini   PetscFunctionReturn(0);
2178135c375SStefano Zampini }
2188135c375SStefano Zampini 
2198135c375SStefano Zampini static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
2208135c375SStefano Zampini {
2218135c375SStefano Zampini   DM                da;
2228135c375SStefano Zampini   Vec               xcoorl;
2238135c375SStefano Zampini   PetscMPIInt       commsize;
2248135c375SStefano Zampini   const PetscScalar *array;
2258135c375SStefano Zampini   PetscContainer    glvis_container;
2268135c375SStefano Zampini   PetscInt          dim,sdim,i,vid[8],mid,cid;
2278135c375SStefano Zampini   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken;
2288135c375SStefano Zampini   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
2298135c375SStefano Zampini   PetscBool         enabled = PETSC_TRUE, isascii;
2308135c375SStefano Zampini   PetscErrorCode    ierr;
2318135c375SStefano Zampini 
2328135c375SStefano Zampini   PetscFunctionBegin;
2338135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2348135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2358135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
2368135c375SStefano Zampini   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
2378135c375SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
2388135c375SStefano Zampini   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
2398135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
2408135c375SStefano Zampini 
2418135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
2428135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
2438135c375SStefano Zampini   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
2448135c375SStefano Zampini   {
2458135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
2468135c375SStefano Zampini     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
2478135c375SStefano Zampini     enabled = glvis_info->enabled;
2488135c375SStefano Zampini   }
2498135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
2508135c375SStefano Zampini   if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
2518135c375SStefano Zampini   ierr = VecGetDM(xcoorl,&da);CHKERRQ(ierr);
2528135c375SStefano Zampini   if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
2538135c375SStefano Zampini   ierr = DMGetCoordinateDim(da,&sdim);CHKERRQ(ierr);
2548135c375SStefano Zampini 
2558135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
2568135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
2578135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
2588135c375SStefano Zampini 
2598135c375SStefano Zampini   if (!enabled) {
2608135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
2618135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
2628135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
2638135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
2648135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
2658135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
2668135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
2678135c375SStefano Zampini     PetscFunctionReturn(0);
2688135c375SStefano Zampini   }
2698135c375SStefano Zampini 
2708135c375SStefano Zampini   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
2718135c375SStefano Zampini   i    = ien;
2728135c375SStefano Zampini   if (dim > 1) i *= jen;
2738135c375SStefano Zampini   if (dim > 2) i *= ken;
2748135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
2758135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",i);CHKERRQ(ierr);
2768135c375SStefano Zampini   switch (dim) {
2778135c375SStefano Zampini   case 1:
2788135c375SStefano Zampini     for (ie = 0; ie < ien; ie++) {
2798135c375SStefano Zampini       vid[0] = ie;
2808135c375SStefano Zampini       vid[1] = ie+1;
2818135c375SStefano Zampini       mid    = 1; /* material id */
2828135c375SStefano Zampini       cid    = 1; /* segment */
2838135c375SStefano Zampini       ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);CHKERRQ(ierr);
2848135c375SStefano Zampini     }
2858135c375SStefano Zampini     break;
2868135c375SStefano Zampini   case 2:
2878135c375SStefano Zampini     for (je = 0; je < jen; je++) {
2888135c375SStefano Zampini       for (ie = 0; ie < ien; ie++) {
2898135c375SStefano Zampini         vid[0] =     je*(ien+1) + ie;
2908135c375SStefano Zampini         vid[1] =     je*(ien+1) + ie+1;
2918135c375SStefano Zampini         vid[2] = (je+1)*(ien+1) + ie+1;
2928135c375SStefano Zampini         vid[3] = (je+1)*(ien+1) + ie;
2938135c375SStefano Zampini         mid    = 1; /* material id */
2948135c375SStefano Zampini         cid    = 3; /* quad */
2958135c375SStefano Zampini         ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);CHKERRQ(ierr);
2968135c375SStefano Zampini       }
2978135c375SStefano Zampini     }
2988135c375SStefano Zampini     break;
2998135c375SStefano Zampini   case 3:
3008135c375SStefano Zampini     for (ke = 0; ke < ken; ke++) {
3018135c375SStefano Zampini       for (je = 0; je < jen; je++) {
3028135c375SStefano Zampini         for (ie = 0; ie < ien; ie++) {
3038135c375SStefano Zampini           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
3048135c375SStefano Zampini           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
3058135c375SStefano Zampini           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
3068135c375SStefano Zampini           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
3078135c375SStefano Zampini           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
3088135c375SStefano Zampini           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
3098135c375SStefano Zampini           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
3108135c375SStefano Zampini           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
3118135c375SStefano Zampini           mid    = 1; /* material id */
3128135c375SStefano Zampini           cid    = 5; /* hex */
3138135c375SStefano Zampini           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);
3148135c375SStefano Zampini         }
3158135c375SStefano Zampini       }
3168135c375SStefano Zampini     }
3178135c375SStefano Zampini     break;
3188135c375SStefano Zampini   default:
3198135c375SStefano Zampini     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
3208135c375SStefano Zampini     break;
3218135c375SStefano Zampini   }
3228135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
3238135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3248135c375SStefano Zampini 
3258135c375SStefano Zampini   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
3268135c375SStefano Zampini   ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
3278135c375SStefano Zampini   ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
3288135c375SStefano Zampini   kst  = gsz != sz ? 1 : 0;
3298135c375SStefano Zampini   jst  = gsy != sy ? 1 : 0;
3308135c375SStefano Zampini   ist  = gsx != sx ? 1 : 0;
3318135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
3328135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);CHKERRQ(ierr);
3338135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
3348135c375SStefano Zampini   ierr = VecGetArrayRead(xcoorl,&array);CHKERRQ(ierr);
3358135c375SStefano Zampini   for (ke = kst; ke < kst + ken; ke++) {
3368135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
3378135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
3388135c375SStefano Zampini         PetscInt d;
3398135c375SStefano Zampini 
3408135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
3418135c375SStefano Zampini         for (d=0;d<sdim-1;d++) {
3428135c375SStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"%g ",PetscRealPart(array[sdim*i+d]));CHKERRQ(ierr);
3438135c375SStefano Zampini         }
3448135c375SStefano Zampini         ierr = PetscViewerASCIIPrintf(viewer,"%g\n",PetscRealPart(array[sdim*i+d]));CHKERRQ(ierr);
3458135c375SStefano Zampini       }
3468135c375SStefano Zampini     }
3478135c375SStefano Zampini   }
3488135c375SStefano Zampini   ierr = VecRestoreArrayRead(xcoorl,&array);CHKERRQ(ierr);
3498135c375SStefano Zampini   PetscFunctionReturn(0);
3508135c375SStefano Zampini }
3518135c375SStefano Zampini 
3528135c375SStefano Zampini /* 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 */
3538135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
3548135c375SStefano Zampini {
3558135c375SStefano Zampini   PetscErrorCode ierr;
3568135c375SStefano Zampini   PetscBool      isglvis,isascii;
3578135c375SStefano Zampini 
3588135c375SStefano Zampini   PetscFunctionBegin;
3598135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3608135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3618135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
3628135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
3638135c375SStefano Zampini   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
3648135c375SStefano Zampini   if (isglvis) {
3658135c375SStefano Zampini     PetscViewer          view;
3668135c375SStefano Zampini     PetscViewerGLVisType type;
3678135c375SStefano Zampini 
3688135c375SStefano Zampini     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
3698135c375SStefano Zampini     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
3708135c375SStefano Zampini     if (view) { /* in the socket case, it may happen that the connection failed */
3718135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
3728135c375SStefano Zampini         PetscMPIInt size,rank;
3738135c375SStefano Zampini         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
3748135c375SStefano Zampini         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
3758135c375SStefano Zampini         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
3768135c375SStefano Zampini       }
3778135c375SStefano Zampini       ierr = DMDAView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
3788135c375SStefano Zampini       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
3798135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
3808135c375SStefano Zampini         PetscInt    dim;
3818135c375SStefano Zampini         const char* name;
3828135c375SStefano Zampini 
3838135c375SStefano Zampini         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
3848135c375SStefano Zampini         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
3858135c375SStefano Zampini         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
3868135c375SStefano Zampini         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
3878135c375SStefano Zampini       }
3888135c375SStefano Zampini     }
3898135c375SStefano Zampini   } else {
3908135c375SStefano Zampini     ierr = DMDAView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
3918135c375SStefano Zampini   }
3928135c375SStefano Zampini   PetscFunctionReturn(0);
3938135c375SStefano Zampini }
394