xref: /petsc/src/dm/impls/da/grglvis.c (revision 8f07809a058053454ca480c7b5a515cfac7e4450)
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;
28aab5bcd8SJed Brown   /* Appease -Wmaybe-uninitialized */
29aab5bcd8SJed Brown   if (nex) *nex = -1;
30aab5bcd8SJed Brown   if (ney) *ney = -1;
31aab5bcd8SJed 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 {
109*8f07809aSStefano Zampini   DM                 da = (DM)oda,daview,dacoord = NULL;
110*8f07809aSStefano Zampini   Vec                xcoor,xcoorl;
1118135c375SStefano Zampini   DMDAGLVisViewerCtx *ctx;
1128135c375SStefano Zampini   const char         **dafieldname;
113*8f07809aSStefano 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;
117*8f07809aSStefano Zampini   PetscBool          bsset,hashocoord = PETSC_FALSE;
1188135c375SStefano Zampini   PetscErrorCode     ierr;
1198135c375SStefano Zampini 
1208135c375SStefano Zampini   PetscFunctionBegin;
1218135c375SStefano Zampini   /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1228135c375SStefano Zampini   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1238135c375SStefano Zampini   ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
12477eacf09SStefano Zampini   ierr = PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor);CHKERRQ(ierr);
12577eacf09SStefano Zampini   if (!xcoor) {
1268135c375SStefano Zampini     ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
12777eacf09SStefano Zampini   } else {
128*8f07809aSStefano Zampini     hashocoord = PETSC_TRUE;
12977eacf09SStefano Zampini   }
1308135c375SStefano Zampini   ierr = PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");CHKERRQ(ierr);
1318135c375SStefano Zampini   switch (dim) {
1328135c375SStefano Zampini   case 1:
1338135c375SStefano Zampini     ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);CHKERRQ(ierr);
134*8f07809aSStefano Zampini     if (!hashocoord) {
135*8f07809aSStefano Zampini       ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);CHKERRQ(ierr);
136*8f07809aSStefano Zampini     }
1378135c375SStefano Zampini     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_1D_P1");CHKERRQ(ierr);
1388135c375SStefano Zampini     break;
1398135c375SStefano Zampini   case 2:
1408135c375SStefano 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);
141*8f07809aSStefano Zampini     if (!hashocoord) {
142*8f07809aSStefano 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);
143*8f07809aSStefano Zampini     }
1448135c375SStefano Zampini     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_2D_P1");CHKERRQ(ierr);
1458135c375SStefano Zampini     break;
1468135c375SStefano Zampini   case 3:
1478135c375SStefano 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);
148*8f07809aSStefano Zampini     if (!hashocoord) {
149*8f07809aSStefano 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);
150*8f07809aSStefano Zampini     }
1518135c375SStefano Zampini     ierr = PetscStrcpy(fec,"FiniteElementCollection: H1_3D_P1");CHKERRQ(ierr);
1528135c375SStefano Zampini     break;
1538135c375SStefano Zampini   default:
1548135c375SStefano Zampini     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
1558135c375SStefano Zampini     break;
1568135c375SStefano Zampini   }
1578135c375SStefano Zampini   ierr = DMSetUp(daview);CHKERRQ(ierr);
1588135c375SStefano Zampini   if (!xcoor) {
1598135c375SStefano Zampini     ierr = DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
1608135c375SStefano Zampini     ierr = DMGetCoordinates(daview,&xcoor);CHKERRQ(ierr);
1618135c375SStefano Zampini   }
162*8f07809aSStefano Zampini   if (dacoord) {
163*8f07809aSStefano Zampini     ierr = DMSetUp(dacoord);CHKERRQ(ierr);
1648135c375SStefano Zampini     ierr = DMCreateLocalVector(dacoord,&xcoorl);CHKERRQ(ierr);
1658135c375SStefano Zampini     ierr = DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
1668135c375SStefano Zampini     ierr = DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
167*8f07809aSStefano Zampini   } else {
168*8f07809aSStefano Zampini     PetscInt ien,jen,ken,nc,nl,cdof,deg;
169*8f07809aSStefano Zampini     char     fecmesh[32];
170*8f07809aSStefano Zampini 
171*8f07809aSStefano Zampini     ierr = DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken);CHKERRQ(ierr);
172*8f07809aSStefano Zampini     nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
173*8f07809aSStefano Zampini 
174*8f07809aSStefano Zampini     ierr = VecGetLocalSize(xcoor,&nl);CHKERRQ(ierr);
175*8f07809aSStefano Zampini     if (nc && nl % nc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc);
176*8f07809aSStefano Zampini     ierr = VecDuplicate(xcoor,&xcoorl);CHKERRQ(ierr);
177*8f07809aSStefano Zampini     ierr = VecCopy(xcoor,xcoorl);CHKERRQ(ierr);
178*8f07809aSStefano Zampini     ierr = VecSetDM(xcoorl,NULL);CHKERRQ(ierr);
179*8f07809aSStefano Zampini     cdof = nl/(nc*dim);
180*8f07809aSStefano Zampini     deg  = 1;
181*8f07809aSStefano Zampini     while (1) {
182*8f07809aSStefano Zampini       PetscInt degd = 1;
183*8f07809aSStefano Zampini       for (i=0;i<dim;i++) degd *= (deg+1);
184*8f07809aSStefano Zampini       if (degd > cdof) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof);
185*8f07809aSStefano Zampini       if (degd == cdof) break;
186*8f07809aSStefano Zampini       deg++;
187*8f07809aSStefano Zampini     }
188*8f07809aSStefano Zampini     ierr = PetscSNPrintf(fecmesh,sizeof(fecmesh),"L2_T1_%dD_P%d",dim,deg);CHKERRQ(ierr);
189*8f07809aSStefano Zampini     ierr = PetscObjectSetName((PetscObject)xcoorl,fecmesh);CHKERRQ(ierr);
190*8f07809aSStefano Zampini   }
1918135c375SStefano Zampini 
1928135c375SStefano Zampini   /* xcoorl is composed with the original DMDA, ghosted coordinate DMDA is only available through this vector */
1938135c375SStefano Zampini   ierr = PetscObjectCompose(oda,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
1948135c375SStefano Zampini   ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
1958135c375SStefano Zampini 
19677eacf09SStefano Zampini   /* daview is composed with the original DMDA */
19777eacf09SStefano Zampini   ierr = PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview);CHKERRQ(ierr);
19877eacf09SStefano Zampini   ierr = PetscObjectDereference((PetscObject)daview);CHKERRQ(ierr);
19977eacf09SStefano Zampini 
200c9493c35SStefano Zampini   /* customize the viewer if present */
201c9493c35SStefano Zampini   if (viewer) {
202*8f07809aSStefano Zampini     Vec xlocal;
203*8f07809aSStefano Zampini 
204*8f07809aSStefano Zampini     ierr = DMCreateLocalVector(daview,&xlocal);CHKERRQ(ierr);
2058135c375SStefano Zampini     ierr = DMDAGetFieldNames(da,(const char * const **)&dafieldname);CHKERRQ(ierr);
2068135c375SStefano Zampini     ierr = DMDAGetNumVerticesGhosted(daview,&M,&N,&P);CHKERRQ(ierr);
207bb77a09fSStefano Zampini     ierr = PetscMalloc5(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname);CHKERRQ(ierr);
2088135c375SStefano Zampini     for (i=0;i<dof;i++) bss[i] = 1;
2098135c375SStefano Zampini     nf = dof;
2108135c375SStefano Zampini 
2118135c375SStefano Zampini     ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");CHKERRQ(ierr);
2128135c375SStefano Zampini     ierr = PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);CHKERRQ(ierr);
2138135c375SStefano Zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2148135c375SStefano Zampini     if (bsset) {
2158135c375SStefano Zampini       PetscInt t;
2168135c375SStefano Zampini       for (i=0,t=0;i<nf;i++) t += bss[i];
2178135c375SStefano Zampini       if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
2188135c375SStefano Zampini     } else nf = dof;
2198135c375SStefano Zampini 
2208135c375SStefano Zampini     for (i=0,s=0;i<nf;i++) {
2218135c375SStefano Zampini       ierr = PetscStrallocpy(fec,&fec_type[i]);CHKERRQ(ierr);
2228135c375SStefano Zampini       if (bss[i] == 1) {
2238135c375SStefano Zampini         ierr = PetscStrallocpy(dafieldname[s],&fieldname[i]);CHKERRQ(ierr);
2248135c375SStefano Zampini       } else {
2258135c375SStefano Zampini         PetscInt b;
2268135c375SStefano Zampini         size_t tlen = 9; /* "Vector-" + end */
2278135c375SStefano Zampini         for (b=0;b<bss[i];b++) {
2288135c375SStefano Zampini           size_t len;
2298135c375SStefano Zampini           ierr = PetscStrlen(dafieldname[s+b],&len);CHKERRQ(ierr);
2308135c375SStefano Zampini           tlen += len + 1; /* field + "-" */
2318135c375SStefano Zampini         }
2328135c375SStefano Zampini         ierr = PetscMalloc1(tlen,&fieldname[i]);CHKERRQ(ierr);
2338135c375SStefano Zampini         ierr = PetscStrcpy(fieldname[i],"Vector-");CHKERRQ(ierr);
2348135c375SStefano Zampini         for (b=0;b<bss[i]-1;b++) {
2358135c375SStefano Zampini           ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
2368135c375SStefano Zampini           ierr = PetscStrcat(fieldname[i],"-");CHKERRQ(ierr);
2378135c375SStefano Zampini         }
2388135c375SStefano Zampini         ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
2398135c375SStefano Zampini       }
240bb77a09fSStefano Zampini       dims[i] = dim;
2418135c375SStefano Zampini       nlocal[i] = M*N*P*bss[i];
2428135c375SStefano Zampini       s += bss[i];
2438135c375SStefano Zampini     }
2448135c375SStefano Zampini 
24577eacf09SStefano Zampini     /* the viewer context takes ownership of xlocal */
2468135c375SStefano Zampini     ierr = PetscNew(&ctx);CHKERRQ(ierr);
2478135c375SStefano Zampini     ctx->xlocal = xlocal;
2488135c375SStefano Zampini 
249bb77a09fSStefano Zampini     ierr = PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,dims,DMDASampleGLVisFields_Private,ctx,DMDADestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
2508135c375SStefano Zampini     for (i=0;i<nf;i++) {
2518135c375SStefano Zampini       ierr = PetscFree(fec_type[i]);CHKERRQ(ierr);
2528135c375SStefano Zampini       ierr = PetscFree(fieldname[i]);CHKERRQ(ierr);
2538135c375SStefano Zampini     }
254bb77a09fSStefano Zampini     ierr = PetscFree5(fec_type,nlocal,bss,dims,fieldname);CHKERRQ(ierr);
255c9493c35SStefano Zampini   }
2568135c375SStefano Zampini   ierr = DMDestroy(&dacoord);CHKERRQ(ierr);
2578135c375SStefano Zampini   PetscFunctionReturn(0);
2588135c375SStefano Zampini }
2598135c375SStefano Zampini 
2608135c375SStefano Zampini static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
2618135c375SStefano Zampini {
26277eacf09SStefano Zampini   DM                da,cda;
2638135c375SStefano Zampini   Vec               xcoorl;
2648135c375SStefano Zampini   PetscMPIInt       commsize;
2658135c375SStefano Zampini   const PetscScalar *array;
2668135c375SStefano Zampini   PetscContainer    glvis_container;
26777eacf09SStefano Zampini   PetscInt          dim,sdim,i,vid[8],mid,cid,cdof;
2688135c375SStefano Zampini   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken;
2698135c375SStefano Zampini   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
2708135c375SStefano Zampini   PetscBool         enabled = PETSC_TRUE, isascii;
2718135c375SStefano Zampini   PetscErrorCode    ierr;
272*8f07809aSStefano Zampini   const char        *fmt;
2738135c375SStefano Zampini 
2748135c375SStefano Zampini   PetscFunctionBegin;
2758135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2768135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2778135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
2788135c375SStefano Zampini   if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
2798135c375SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr);
2808135c375SStefano Zampini   if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
2818135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
2828135c375SStefano Zampini 
2838135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
2848135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
2858135c375SStefano Zampini   if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
2868135c375SStefano Zampini   {
2878135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
2888135c375SStefano Zampini     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
2898135c375SStefano Zampini     enabled = glvis_info->enabled;
29077eacf09SStefano Zampini     fmt = glvis_info->fmt;
2918135c375SStefano Zampini   }
2928135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
293c9493c35SStefano Zampini   /* this can happen if we are calling DMView outside of VecView_GLVis */
294c9493c35SStefano Zampini   if (!xcoorl) {ierr = DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);CHKERRQ(ierr);}
29577eacf09SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);CHKERRQ(ierr);
2968135c375SStefano Zampini   if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
2978135c375SStefano Zampini   ierr = DMGetCoordinateDim(da,&sdim);CHKERRQ(ierr);
2988135c375SStefano Zampini 
2998135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr);
3008135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
3018135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
3028135c375SStefano Zampini 
3038135c375SStefano Zampini   if (!enabled) {
3048135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
3058135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3068135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
3078135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3088135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
3098135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3108135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
3118135c375SStefano Zampini     PetscFunctionReturn(0);
3128135c375SStefano Zampini   }
3138135c375SStefano Zampini 
3148135c375SStefano Zampini   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
3158135c375SStefano Zampini   i    = ien;
3168135c375SStefano Zampini   if (dim > 1) i *= jen;
3178135c375SStefano Zampini   if (dim > 2) i *= ken;
3188135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
3198135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",i);CHKERRQ(ierr);
3208135c375SStefano Zampini   switch (dim) {
3218135c375SStefano Zampini   case 1:
3228135c375SStefano Zampini     for (ie = 0; ie < ien; ie++) {
3238135c375SStefano Zampini       vid[0] = ie;
3248135c375SStefano Zampini       vid[1] = ie+1;
3258135c375SStefano Zampini       mid    = 1; /* material id */
3268135c375SStefano Zampini       cid    = 1; /* segment */
3278135c375SStefano Zampini       ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);CHKERRQ(ierr);
3288135c375SStefano Zampini     }
3298135c375SStefano Zampini     break;
3308135c375SStefano Zampini   case 2:
3318135c375SStefano Zampini     for (je = 0; je < jen; je++) {
3328135c375SStefano Zampini       for (ie = 0; ie < ien; ie++) {
3338135c375SStefano Zampini         vid[0] =     je*(ien+1) + ie;
3348135c375SStefano Zampini         vid[1] =     je*(ien+1) + ie+1;
3358135c375SStefano Zampini         vid[2] = (je+1)*(ien+1) + ie+1;
3368135c375SStefano Zampini         vid[3] = (je+1)*(ien+1) + ie;
3378135c375SStefano Zampini         mid    = 1; /* material id */
3388135c375SStefano Zampini         cid    = 3; /* quad */
3398135c375SStefano Zampini         ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);CHKERRQ(ierr);
3408135c375SStefano Zampini       }
3418135c375SStefano Zampini     }
3428135c375SStefano Zampini     break;
3438135c375SStefano Zampini   case 3:
3448135c375SStefano Zampini     for (ke = 0; ke < ken; ke++) {
3458135c375SStefano Zampini       for (je = 0; je < jen; je++) {
3468135c375SStefano Zampini         for (ie = 0; ie < ien; ie++) {
3478135c375SStefano Zampini           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
3488135c375SStefano Zampini           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
3498135c375SStefano Zampini           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
3508135c375SStefano Zampini           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
3518135c375SStefano Zampini           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
3528135c375SStefano Zampini           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
3538135c375SStefano Zampini           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
3548135c375SStefano Zampini           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
3558135c375SStefano Zampini           mid    = 1; /* material id */
3568135c375SStefano Zampini           cid    = 5; /* hex */
3578135c375SStefano 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);
3588135c375SStefano Zampini         }
3598135c375SStefano Zampini       }
3608135c375SStefano Zampini     }
3618135c375SStefano Zampini     break;
3628135c375SStefano Zampini   default:
3638135c375SStefano Zampini     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
3648135c375SStefano Zampini     break;
3658135c375SStefano Zampini   }
3668135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
3678135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3688135c375SStefano Zampini 
36977eacf09SStefano Zampini   /* vertex coordinates */
37077eacf09SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
37177eacf09SStefano Zampini   if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
372*8f07809aSStefano Zampini   ierr = VecGetDM(xcoorl,&cda);CHKERRQ(ierr);
3738135c375SStefano Zampini   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
37477eacf09SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
37577eacf09SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);CHKERRQ(ierr);
37677eacf09SStefano Zampini   ierr = VecGetArrayRead(xcoorl,&array);CHKERRQ(ierr);
377*8f07809aSStefano Zampini   if (!cda) { /* HO viz */
378*8f07809aSStefano Zampini     const char *fecname;
379*8f07809aSStefano Zampini     PetscInt   nc,nl;
380*8f07809aSStefano Zampini 
381*8f07809aSStefano Zampini     ierr = PetscObjectGetName((PetscObject)xcoorl,&fecname);CHKERRQ(ierr);
38277eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
38377eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
38477eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementCollection: %s\n",fecname);CHKERRQ(ierr);
38577eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
38677eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
387*8f07809aSStefano Zampini     /* L2 coordinates */
388*8f07809aSStefano Zampini     ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
389*8f07809aSStefano Zampini     ierr = VecGetLocalSize(xcoorl,&nl);CHKERRQ(ierr);
390*8f07809aSStefano Zampini     nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
391*8f07809aSStefano Zampini     cdof = nl/nc;
392*8f07809aSStefano Zampini     if (!ien) ien++;
393*8f07809aSStefano Zampini     if (!jen) jen++;
394*8f07809aSStefano Zampini     if (!ken) ken++;
395*8f07809aSStefano Zampini     ist = jst = kst = 0;
396*8f07809aSStefano Zampini     gm = ien;
397*8f07809aSStefano Zampini     gn = jen;
398*8f07809aSStefano Zampini     gp = ken;
39977eacf09SStefano Zampini   } else {
40077eacf09SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
401*8f07809aSStefano Zampini     cdof = sdim;
4028135c375SStefano Zampini     ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
4038135c375SStefano Zampini     ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
4048135c375SStefano Zampini     kst  = gsz != sz ? 1 : 0;
4058135c375SStefano Zampini     jst  = gsy != sy ? 1 : 0;
4068135c375SStefano Zampini     ist  = gsx != sx ? 1 : 0;
407*8f07809aSStefano Zampini   }
4088135c375SStefano Zampini   for (ke = kst; ke < kst + ken; ke++) {
4098135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
4108135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
41177eacf09SStefano Zampini         PetscInt c;
4128135c375SStefano Zampini 
4138135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
41477eacf09SStefano Zampini         for (c=0;c<cdof/sdim;c++) {
41577eacf09SStefano Zampini           PetscInt d;
41677eacf09SStefano Zampini           for (d=0;d<sdim;d++) {
417*8f07809aSStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]));CHKERRQ(ierr);
4188135c375SStefano Zampini           }
419*8f07809aSStefano Zampini           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
42077eacf09SStefano Zampini         }
4218135c375SStefano Zampini       }
4228135c375SStefano Zampini     }
4238135c375SStefano Zampini   }
4248135c375SStefano Zampini   ierr = VecRestoreArrayRead(xcoorl,&array);CHKERRQ(ierr);
4258135c375SStefano Zampini   PetscFunctionReturn(0);
4268135c375SStefano Zampini }
4278135c375SStefano Zampini 
4288135c375SStefano 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 */
4298135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
4308135c375SStefano Zampini {
4318135c375SStefano Zampini   PetscErrorCode ierr;
4328135c375SStefano Zampini   PetscBool      isglvis,isascii;
4338135c375SStefano Zampini 
4348135c375SStefano Zampini   PetscFunctionBegin;
4358135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4368135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4378135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
4388135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
4398135c375SStefano Zampini   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
4408135c375SStefano Zampini   if (isglvis) {
4418135c375SStefano Zampini     PetscViewer          view;
4428135c375SStefano Zampini     PetscViewerGLVisType type;
4438135c375SStefano Zampini 
4448135c375SStefano Zampini     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
4458135c375SStefano Zampini     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
4468135c375SStefano Zampini     if (view) { /* in the socket case, it may happen that the connection failed */
4478135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
4488135c375SStefano Zampini         PetscMPIInt size,rank;
4498135c375SStefano Zampini         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
4508135c375SStefano Zampini         ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
4518135c375SStefano Zampini         ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr);
4528135c375SStefano Zampini       }
4538135c375SStefano Zampini       ierr = DMDAView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
4548135c375SStefano Zampini       ierr = PetscViewerFlush(view);CHKERRQ(ierr);
4558135c375SStefano Zampini       if (type == PETSC_VIEWER_GLVIS_SOCKET) {
4568135c375SStefano Zampini         PetscInt    dim;
4578135c375SStefano Zampini         const char* name;
4588135c375SStefano Zampini 
4598135c375SStefano Zampini         ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
4608135c375SStefano Zampini         ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
4618135c375SStefano Zampini         ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr);
4628135c375SStefano Zampini         ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr);
4638135c375SStefano Zampini       }
4648135c375SStefano Zampini     }
4658135c375SStefano Zampini   } else {
4668135c375SStefano Zampini     ierr = DMDAView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
4678135c375SStefano Zampini   }
4688135c375SStefano Zampini   PetscFunctionReturn(0);
4698135c375SStefano Zampini }
470