xref: /petsc/src/dm/impls/da/grglvis.c (revision 7a8be3513cf479fb6a672bd91de7eb6883fcfd02)
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 {
746900a7fSStefano Zampini   PetscBool ll;
846900a7fSStefano Zampini } DMDAGhostedGLVisViewerCtx;
98135c375SStefano Zampini 
1046900a7fSStefano Zampini static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(void **vctx)
118135c375SStefano Zampini {
1246900a7fSStefano Zampini   PetscErrorCode ierr;
1346900a7fSStefano Zampini 
1446900a7fSStefano Zampini   PetscFunctionBegin;
1546900a7fSStefano Zampini   ierr = PetscFree(*vctx);CHKERRQ(ierr);
1646900a7fSStefano Zampini   PetscFunctionReturn(0);
1746900a7fSStefano Zampini }
1846900a7fSStefano Zampini 
1946900a7fSStefano Zampini typedef struct {
2046900a7fSStefano Zampini   Vec xlocal;
2146900a7fSStefano Zampini } DMDAFieldGLVisViewerCtx;
2246900a7fSStefano Zampini 
2346900a7fSStefano Zampini static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx)
2446900a7fSStefano Zampini {
2546900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
268135c375SStefano Zampini   PetscErrorCode          ierr;
278135c375SStefano Zampini 
288135c375SStefano Zampini   PetscFunctionBegin;
298135c375SStefano Zampini   ierr = VecDestroy(&ctx->xlocal);CHKERRQ(ierr);
3051f51421SSatish Balay   ierr = PetscFree(vctx);CHKERRQ(ierr);
318135c375SStefano Zampini   PetscFunctionReturn(0);
328135c375SStefano Zampini }
338135c375SStefano Zampini 
3446900a7fSStefano Zampini /*
3546900a7fSStefano Zampini    dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right
3646900a7fSStefano Zampini    dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left
3746900a7fSStefano Zampini */
388135c375SStefano Zampini static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
398135c375SStefano Zampini {
4046900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
4146900a7fSStefano Zampini   PetscInt                  sx,sy,sz,ien,jen,ken;
428135c375SStefano Zampini   PetscErrorCode            ierr;
438135c375SStefano Zampini 
448135c375SStefano Zampini   PetscFunctionBegin;
45aab5bcd8SJed Brown   /* Appease -Wmaybe-uninitialized */
46aab5bcd8SJed Brown   if (nex) *nex = -1;
47aab5bcd8SJed Brown   if (ney) *ney = -1;
48aab5bcd8SJed Brown   if (nez) *nez = -1;
498135c375SStefano Zampini   ierr = DMDAGetCorners(da,&sx,&sy,&sz,&ien,&jen,&ken);CHKERRQ(ierr);
503ec1f749SStefano Zampini   ierr = DMGetApplicationContext(da,&dactx);CHKERRQ(ierr);
5146900a7fSStefano Zampini   if (dactx->ll) {
5246900a7fSStefano Zampini     PetscInt dim;
5346900a7fSStefano Zampini 
5446900a7fSStefano Zampini     ierr = DMGetDimension(da,&dim);CHKERRQ(ierr);
5546900a7fSStefano Zampini     if (!sx) ien--;
5646900a7fSStefano Zampini     if (!sy && dim > 1) jen--;
5746900a7fSStefano Zampini     if (!sz && dim > 2) ken--;
5846900a7fSStefano Zampini   } else {
5946900a7fSStefano Zampini     PetscInt M,N,P;
6046900a7fSStefano Zampini 
6146900a7fSStefano Zampini     ierr = DMDAGetInfo(da,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
628135c375SStefano Zampini     if (sx+ien == M) ien--;
638135c375SStefano Zampini     if (sy+jen == N) jen--;
648135c375SStefano Zampini     if (sz+ken == P) ken--;
6546900a7fSStefano Zampini   }
668135c375SStefano Zampini   if (nex) *nex = ien;
678135c375SStefano Zampini   if (ney) *ney = jen;
688135c375SStefano Zampini   if (nez) *nez = ken;
698135c375SStefano Zampini   PetscFunctionReturn(0);
708135c375SStefano Zampini }
718135c375SStefano Zampini 
728135c375SStefano Zampini /* inherits number of vertices from DMDAGetNumElementsGhosted */
738135c375SStefano Zampini static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
748135c375SStefano Zampini {
7521414b21SStefano Zampini   PetscInt       ien = 0,jen = 0,ken = 0,dim;
7621414b21SStefano Zampini   PetscInt       tote;
778135c375SStefano Zampini   PetscErrorCode ierr;
788135c375SStefano Zampini 
798135c375SStefano Zampini   PetscFunctionBegin;
808135c375SStefano Zampini   ierr = DMGetDimension(da,&dim);CHKERRQ(ierr);
818135c375SStefano Zampini   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
8221414b21SStefano Zampini   tote = ien * (dim > 1 ? jen : 1) * (dim > 2 ? ken : 1);
8321414b21SStefano Zampini   if (tote) {
848135c375SStefano Zampini     ien = ien+1;
858135c375SStefano Zampini     jen = dim > 1 ? jen+1 : 1;
868135c375SStefano Zampini     ken = dim > 2 ? ken+1 : 1;
8721414b21SStefano Zampini   }
888135c375SStefano Zampini   if (nvx) *nvx = ien;
898135c375SStefano Zampini   if (nvy) *nvy = jen;
908135c375SStefano Zampini   if (nvz) *nvz = ken;
918135c375SStefano Zampini   PetscFunctionReturn(0);
928135c375SStefano Zampini }
938135c375SStefano Zampini 
948135c375SStefano Zampini static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
958135c375SStefano Zampini {
968135c375SStefano Zampini   DM                        da;
9746900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx   *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
9846900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
998135c375SStefano Zampini   const PetscScalar         *array;
1008135c375SStefano Zampini   PetscScalar               **arrayf;
1018135c375SStefano Zampini   PetscInt                  i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
1028135c375SStefano Zampini   PetscInt                  sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
1038135c375SStefano Zampini   PetscErrorCode            ierr;
1048135c375SStefano Zampini 
1058135c375SStefano Zampini   PetscFunctionBegin;
1068135c375SStefano Zampini   ierr = VecGetDM(ctx->xlocal,&da);CHKERRQ(ierr);
107*7a8be351SBarry Smith   PetscCheck(da,PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
1083ec1f749SStefano Zampini   ierr = DMGetApplicationContext(da,&dactx);CHKERRQ(ierr);
1098135c375SStefano Zampini   ierr = VecGetBlockSize(ctx->xlocal,&bs);CHKERRQ(ierr);
1108135c375SStefano Zampini   ierr = DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr);
1118135c375SStefano Zampini   ierr = DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr);
1128135c375SStefano Zampini   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
1138135c375SStefano Zampini   ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
1148135c375SStefano Zampini   ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
11546900a7fSStefano Zampini   if (dactx->ll) {
11646900a7fSStefano Zampini     kst = jst = ist = 0;
11746900a7fSStefano Zampini   } else {
1188135c375SStefano Zampini     kst  = gsz != sz ? 1 : 0;
1198135c375SStefano Zampini     jst  = gsy != sy ? 1 : 0;
1208135c375SStefano Zampini     ist  = gsx != sx ? 1 : 0;
12146900a7fSStefano Zampini   }
1228135c375SStefano Zampini   ierr = PetscMalloc2(nf,&arrayf,nf,&bss);CHKERRQ(ierr);
1238135c375SStefano Zampini   ierr = VecGetArrayRead(ctx->xlocal,&array);CHKERRQ(ierr);
1248135c375SStefano Zampini   for (f=0;f<nf;f++) {
1258135c375SStefano Zampini     ierr = VecGetBlockSize((Vec)oXf[f],&bss[f]);CHKERRQ(ierr);
1268135c375SStefano Zampini     ierr = VecGetArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr);
1278135c375SStefano Zampini   }
1288135c375SStefano Zampini   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
1298135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
1308135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
1318135c375SStefano Zampini         PetscInt cf,b;
1328135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
1338135c375SStefano Zampini         for (f=0,cf=0;f<nf;f++)
1348135c375SStefano Zampini           for (b=0;b<bss[f];b++)
1358135c375SStefano Zampini             arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
1368135c375SStefano Zampini         ii++;
1378135c375SStefano Zampini       }
1388135c375SStefano Zampini     }
1398135c375SStefano Zampini   }
1408135c375SStefano Zampini   for (f=0;f<nf;f++) { ierr = VecRestoreArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr); }
1418135c375SStefano Zampini   ierr = VecRestoreArrayRead(ctx->xlocal,&array);CHKERRQ(ierr);
1428135c375SStefano Zampini   ierr = PetscFree2(arrayf,bss);CHKERRQ(ierr);
1438135c375SStefano Zampini   PetscFunctionReturn(0);
1448135c375SStefano Zampini }
1458135c375SStefano Zampini 
1468135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
1478135c375SStefano Zampini {
14846900a7fSStefano Zampini   DM             da = (DM)oda,daview;
1498135c375SStefano Zampini   PetscErrorCode ierr;
1508135c375SStefano Zampini 
1518135c375SStefano Zampini   PetscFunctionBegin;
15246900a7fSStefano Zampini   ierr = PetscObjectQuery(oda,"GLVisGraphicsDMDAGhosted",(PetscObject*)&daview);CHKERRQ(ierr);
15346900a7fSStefano Zampini   if (!daview) {
15446900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
15546900a7fSStefano Zampini     DM                        dacoord = NULL;
15646900a7fSStefano Zampini     Vec                       xcoor,xcoorl;
15746900a7fSStefano Zampini     PetscBool                 hashocoord = PETSC_FALSE;
15846900a7fSStefano Zampini     const PetscInt            *lx,*ly,*lz;
15946900a7fSStefano Zampini     PetscInt                  dim,M,N,P,m,n,p,dof,s,i;
16046900a7fSStefano Zampini 
16146900a7fSStefano Zampini     ierr = PetscNew(&dactx);CHKERRQ(ierr);
162a7bda6aaSStefano Zampini     dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
16346900a7fSStefano Zampini     ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");CHKERRQ(ierr);
16446900a7fSStefano Zampini     ierr = PetscOptionsBool("-viewer_glvis_dm_da_ll","Left-looking subdomain view",NULL,dactx->ll,&dactx->ll,NULL);CHKERRQ(ierr);
16546900a7fSStefano Zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1668135c375SStefano Zampini     /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1678135c375SStefano Zampini     ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
1688135c375SStefano Zampini     ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr);
16977eacf09SStefano Zampini     ierr = PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor);CHKERRQ(ierr);
17077eacf09SStefano Zampini     if (!xcoor) {
1718135c375SStefano Zampini       ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr);
17277eacf09SStefano Zampini     } else {
1738f07809aSStefano Zampini       hashocoord = PETSC_TRUE;
17477eacf09SStefano Zampini     }
1758135c375SStefano Zampini     ierr = PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");CHKERRQ(ierr);
1768135c375SStefano Zampini     switch (dim) {
1778135c375SStefano Zampini     case 1:
1788135c375SStefano Zampini       ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);CHKERRQ(ierr);
1798f07809aSStefano Zampini       if (!hashocoord) {
1808135c375SStefano Zampini         ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);CHKERRQ(ierr);
1818f07809aSStefano Zampini       }
1828135c375SStefano Zampini       break;
1838135c375SStefano Zampini     case 2:
1848135c375SStefano 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);
1858f07809aSStefano Zampini       if (!hashocoord) {
1868135c375SStefano 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);
1878f07809aSStefano Zampini       }
1888135c375SStefano Zampini       break;
1898135c375SStefano Zampini     case 3:
1908135c375SStefano 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);
1918f07809aSStefano Zampini       if (!hashocoord) {
1928135c375SStefano 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);
1938f07809aSStefano Zampini       }
1948135c375SStefano Zampini       break;
1958135c375SStefano Zampini     default:
19698921bdaSJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
1978135c375SStefano Zampini     }
19846900a7fSStefano Zampini     ierr = DMSetApplicationContext(daview,dactx);CHKERRQ(ierr);
19946900a7fSStefano Zampini     ierr = DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
2008135c375SStefano Zampini     ierr = DMSetUp(daview);CHKERRQ(ierr);
2018135c375SStefano Zampini     if (!xcoor) {
2028135c375SStefano Zampini       ierr = DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
2038135c375SStefano Zampini       ierr = DMGetCoordinates(daview,&xcoor);CHKERRQ(ierr);
2048135c375SStefano Zampini     }
2058f07809aSStefano Zampini     if (dacoord) {
2068f07809aSStefano Zampini       ierr = DMSetUp(dacoord);CHKERRQ(ierr);
2078135c375SStefano Zampini       ierr = DMCreateLocalVector(dacoord,&xcoorl);CHKERRQ(ierr);
2088135c375SStefano Zampini       ierr = DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
2098135c375SStefano Zampini       ierr = DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr);
21046900a7fSStefano Zampini       ierr = DMDestroy(&dacoord);CHKERRQ(ierr);
2118f07809aSStefano Zampini     } else {
2128f07809aSStefano Zampini       PetscInt   ien,jen,ken,nc,nl,cdof,deg;
2133e7cab22SLisandro Dalcin       char       fecmesh[64];
21421414b21SStefano Zampini       const char *name;
21521414b21SStefano Zampini       PetscBool  flg;
2168135c375SStefano Zampini 
2178f07809aSStefano Zampini       ierr = DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken);CHKERRQ(ierr);
2188f07809aSStefano Zampini       nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
2198f07809aSStefano Zampini 
2208f07809aSStefano Zampini       ierr = VecGetLocalSize(xcoor,&nl);CHKERRQ(ierr);
2212c71b3e2SJacob Faibussowitsch       PetscCheckFalse(nc && nl % nc,PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc);
2228f07809aSStefano Zampini       ierr = VecDuplicate(xcoor,&xcoorl);CHKERRQ(ierr);
2238f07809aSStefano Zampini       ierr = VecCopy(xcoor,xcoorl);CHKERRQ(ierr);
2248f07809aSStefano Zampini       ierr = VecSetDM(xcoorl,NULL);CHKERRQ(ierr);
22521414b21SStefano Zampini       ierr = PetscObjectGetName((PetscObject)xcoor,&name);CHKERRQ(ierr);
22621414b21SStefano Zampini       ierr = PetscStrbeginswith(name,"FiniteElementCollection:",&flg);CHKERRQ(ierr);
22721414b21SStefano Zampini       if (!flg) {
22821414b21SStefano Zampini         deg = 0;
22921414b21SStefano Zampini         if (nc && nl) {
2308f07809aSStefano Zampini           cdof = nl/(nc*dim);
2318f07809aSStefano Zampini           deg  = 1;
2328f07809aSStefano Zampini           while (1) {
2338f07809aSStefano Zampini             PetscInt degd = 1;
2348f07809aSStefano Zampini             for (i=0;i<dim;i++) degd *= (deg+1);
2352c71b3e2SJacob Faibussowitsch             PetscCheckFalse(degd > cdof,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof);
2368f07809aSStefano Zampini             if (degd == cdof) break;
2378f07809aSStefano Zampini             deg++;
2388f07809aSStefano Zampini           }
23921414b21SStefano Zampini         }
240bfadf1d8SStefano Zampini         ierr = PetscSNPrintf(fecmesh,sizeof(fecmesh),"FiniteElementCollection: L2_T1_%DD_P%D",dim,deg);CHKERRQ(ierr);
2418f07809aSStefano Zampini         ierr = PetscObjectSetName((PetscObject)xcoorl,fecmesh);CHKERRQ(ierr);
24221414b21SStefano Zampini       } else {
24321414b21SStefano Zampini         ierr = PetscObjectSetName((PetscObject)xcoorl,name);CHKERRQ(ierr);
24421414b21SStefano Zampini       }
2458f07809aSStefano Zampini     }
2468135c375SStefano Zampini 
24746900a7fSStefano Zampini     /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
24846900a7fSStefano Zampini     ierr = PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);CHKERRQ(ierr);
2498135c375SStefano Zampini     ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr);
2508135c375SStefano Zampini 
25177eacf09SStefano Zampini     /* daview is composed with the original DMDA */
25277eacf09SStefano Zampini     ierr = PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview);CHKERRQ(ierr);
25377eacf09SStefano Zampini     ierr = PetscObjectDereference((PetscObject)daview);CHKERRQ(ierr);
25446900a7fSStefano Zampini   }
25577eacf09SStefano Zampini 
256c9493c35SStefano Zampini   /* customize the viewer if present */
257c9493c35SStefano Zampini   if (viewer) {
25846900a7fSStefano Zampini     DMDAFieldGLVisViewerCtx   *ctx;
25946900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
26046900a7fSStefano Zampini     char                      fec[64];
261ab166b77SStefano Zampini     Vec                       xlocal,*Ufield;
26246900a7fSStefano Zampini     const char                **dafieldname;
26346900a7fSStefano Zampini     char                      **fec_type,**fieldname;
26446900a7fSStefano Zampini     PetscInt                  *nlocal,*bss,*dims;
26546900a7fSStefano Zampini     PetscInt                  dim,M,N,P,dof,s,i,nf;
26646900a7fSStefano Zampini     PetscBool                 bsset;
2678f07809aSStefano Zampini 
26821414b21SStefano Zampini     ierr = DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2693ec1f749SStefano Zampini     ierr = DMGetApplicationContext(daview,&dactx);CHKERRQ(ierr);
2708f07809aSStefano Zampini     ierr = DMCreateLocalVector(daview,&xlocal);CHKERRQ(ierr);
2718135c375SStefano Zampini     ierr = DMDAGetFieldNames(da,(const char * const **)&dafieldname);CHKERRQ(ierr);
2728135c375SStefano Zampini     ierr = DMDAGetNumVerticesGhosted(daview,&M,&N,&P);CHKERRQ(ierr);
273bfadf1d8SStefano Zampini     ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%DD_P1",dim);CHKERRQ(ierr);
2744cac2994SStefano Zampini     ierr = PetscMalloc6(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname,dof,&Ufield);CHKERRQ(ierr);
2758135c375SStefano Zampini     for (i=0;i<dof;i++) bss[i] = 1;
2768135c375SStefano Zampini     nf = dof;
2778135c375SStefano Zampini 
27846900a7fSStefano Zampini     ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");CHKERRQ(ierr);
2798135c375SStefano Zampini     ierr = PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);CHKERRQ(ierr);
2808135c375SStefano Zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2818135c375SStefano Zampini     if (bsset) {
2828135c375SStefano Zampini       PetscInt t;
2838135c375SStefano Zampini       for (i=0,t=0;i<nf;i++) t += bss[i];
2842c71b3e2SJacob Faibussowitsch       PetscCheckFalse(t != dof,PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof);
2858135c375SStefano Zampini     } else nf = dof;
2868135c375SStefano Zampini 
2878135c375SStefano Zampini     for (i=0,s=0;i<nf;i++) {
2888135c375SStefano Zampini       ierr = PetscStrallocpy(fec,&fec_type[i]);CHKERRQ(ierr);
2898135c375SStefano Zampini       if (bss[i] == 1) {
2908135c375SStefano Zampini         ierr = PetscStrallocpy(dafieldname[s],&fieldname[i]);CHKERRQ(ierr);
2918135c375SStefano Zampini       } else {
2928135c375SStefano Zampini         PetscInt b;
2938135c375SStefano Zampini         size_t tlen = 9; /* "Vector-" + end */
2948135c375SStefano Zampini         for (b=0;b<bss[i];b++) {
2958135c375SStefano Zampini           size_t len;
2968135c375SStefano Zampini           ierr = PetscStrlen(dafieldname[s+b],&len);CHKERRQ(ierr);
2978135c375SStefano Zampini           tlen += len + 1; /* field + "-" */
2988135c375SStefano Zampini         }
2998135c375SStefano Zampini         ierr = PetscMalloc1(tlen,&fieldname[i]);CHKERRQ(ierr);
3008135c375SStefano Zampini         ierr = PetscStrcpy(fieldname[i],"Vector-");CHKERRQ(ierr);
3018135c375SStefano Zampini         for (b=0;b<bss[i]-1;b++) {
3028135c375SStefano Zampini           ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
3038135c375SStefano Zampini           ierr = PetscStrcat(fieldname[i],"-");CHKERRQ(ierr);
3048135c375SStefano Zampini         }
3058135c375SStefano Zampini         ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr);
3068135c375SStefano Zampini       }
307bb77a09fSStefano Zampini       dims[i] = dim;
3088135c375SStefano Zampini       nlocal[i] = M*N*P*bss[i];
3098135c375SStefano Zampini       s += bss[i];
3108135c375SStefano Zampini     }
3118135c375SStefano Zampini 
31246900a7fSStefano Zampini     /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
3138135c375SStefano Zampini     ierr = PetscNew(&ctx);CHKERRQ(ierr);
3148135c375SStefano Zampini     ctx->xlocal = xlocal;
3158135c375SStefano Zampini 
3164cac2994SStefano Zampini     /* create work vectors */
3174cac2994SStefano Zampini     for (i=0;i<nf;i++) {
3184cac2994SStefano Zampini       ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),nlocal[i],PETSC_DECIDE,&Ufield[i]);CHKERRQ(ierr);
3194cac2994SStefano Zampini       ierr = PetscObjectSetName((PetscObject)Ufield[i],fieldname[i]);CHKERRQ(ierr);
3204cac2994SStefano Zampini       ierr = VecSetBlockSize(Ufield[i],bss[i]);CHKERRQ(ierr);
3214cac2994SStefano Zampini       ierr = VecSetDM(Ufield[i],da);CHKERRQ(ierr);
3224cac2994SStefano Zampini     }
3234cac2994SStefano Zampini 
324ab166b77SStefano Zampini     ierr = PetscViewerGLVisSetFields(viewer,nf,(const char**)fec_type,dims,DMDASampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DMDAFieldDestroyGLVisViewerCtx_Private);CHKERRQ(ierr);
3258135c375SStefano Zampini     for (i=0;i<nf;i++) {
3268135c375SStefano Zampini       ierr = PetscFree(fec_type[i]);CHKERRQ(ierr);
3278135c375SStefano Zampini       ierr = PetscFree(fieldname[i]);CHKERRQ(ierr);
3284cac2994SStefano Zampini       ierr = VecDestroy(&Ufield[i]);CHKERRQ(ierr);
3298135c375SStefano Zampini     }
3304cac2994SStefano Zampini     ierr = PetscFree6(fec_type,nlocal,bss,dims,fieldname,Ufield);CHKERRQ(ierr);
331c9493c35SStefano Zampini   }
3328135c375SStefano Zampini   PetscFunctionReturn(0);
3338135c375SStefano Zampini }
3348135c375SStefano Zampini 
3358135c375SStefano Zampini static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
3368135c375SStefano Zampini {
33777eacf09SStefano Zampini   DM                da,cda;
3388135c375SStefano Zampini   Vec               xcoorl;
3393924b612SStefano Zampini   PetscMPIInt       size;
3408135c375SStefano Zampini   const PetscScalar *array;
3418135c375SStefano Zampini   PetscContainer    glvis_container;
34277eacf09SStefano Zampini   PetscInt          dim,sdim,i,vid[8],mid,cid,cdof;
34321414b21SStefano Zampini   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken,nel;
3448135c375SStefano Zampini   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
3458135c375SStefano Zampini   PetscBool         enabled = PETSC_TRUE, isascii;
3468135c375SStefano Zampini   PetscErrorCode    ierr;
3478f07809aSStefano Zampini   const char        *fmt;
3488135c375SStefano Zampini 
3498135c375SStefano Zampini   PetscFunctionBegin;
3508135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3518135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3528135c375SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
353*7a8be351SBarry Smith   PetscCheck(isascii,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
354ffc4695bSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRMPI(ierr);
3552c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
3568135c375SStefano Zampini   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
3578135c375SStefano Zampini 
3588135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
3598135c375SStefano Zampini   ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr);
360*7a8be351SBarry Smith   PetscCheck(glvis_container,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
3618135c375SStefano Zampini   {
3628135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
3638135c375SStefano Zampini     ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr);
3648135c375SStefano Zampini     enabled = glvis_info->enabled;
36577eacf09SStefano Zampini     fmt = glvis_info->fmt;
3668135c375SStefano Zampini   }
367c9493c35SStefano Zampini   /* this can happen if we are calling DMView outside of VecView_GLVis */
36846900a7fSStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);CHKERRQ(ierr);
36946900a7fSStefano Zampini   if (!da) {ierr = DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);CHKERRQ(ierr);}
37077eacf09SStefano Zampini   ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);CHKERRQ(ierr);
371*7a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
3728135c375SStefano Zampini   ierr = DMGetCoordinateDim(da,&sdim);CHKERRQ(ierr);
3738135c375SStefano Zampini 
3741d4815f0SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.0\n");CHKERRQ(ierr);
3758135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr);
3768135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr);
3778135c375SStefano Zampini 
3788135c375SStefano Zampini   if (!enabled) {
3798135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
3808135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3818135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
3828135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3838135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
3848135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
3858135c375SStefano Zampini     ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
3868135c375SStefano Zampini     PetscFunctionReturn(0);
3878135c375SStefano Zampini   }
3888135c375SStefano Zampini 
3898135c375SStefano Zampini   ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
39021414b21SStefano Zampini   nel  = ien;
39121414b21SStefano Zampini   if (dim > 1) nel *= jen;
39221414b21SStefano Zampini   if (dim > 2) nel *= ken;
3938135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr);
39421414b21SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",nel);CHKERRQ(ierr);
3958135c375SStefano Zampini   switch (dim) {
3968135c375SStefano Zampini   case 1:
3978135c375SStefano Zampini     for (ie = 0; ie < ien; ie++) {
3988135c375SStefano Zampini       vid[0] = ie;
3998135c375SStefano Zampini       vid[1] = ie+1;
4008135c375SStefano Zampini       mid    = 1; /* material id */
4018135c375SStefano Zampini       cid    = 1; /* segment */
4028135c375SStefano Zampini       ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);CHKERRQ(ierr);
4038135c375SStefano Zampini     }
4048135c375SStefano Zampini     break;
4058135c375SStefano Zampini   case 2:
4068135c375SStefano Zampini     for (je = 0; je < jen; je++) {
4078135c375SStefano Zampini       for (ie = 0; ie < ien; ie++) {
4088135c375SStefano Zampini         vid[0] =     je*(ien+1) + ie;
4098135c375SStefano Zampini         vid[1] =     je*(ien+1) + ie+1;
4108135c375SStefano Zampini         vid[2] = (je+1)*(ien+1) + ie+1;
4118135c375SStefano Zampini         vid[3] = (je+1)*(ien+1) + ie;
4128135c375SStefano Zampini         mid    = 1; /* material id */
4138135c375SStefano Zampini         cid    = 3; /* quad */
4148135c375SStefano Zampini         ierr   = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);CHKERRQ(ierr);
4158135c375SStefano Zampini       }
4168135c375SStefano Zampini     }
4178135c375SStefano Zampini     break;
4188135c375SStefano Zampini   case 3:
4198135c375SStefano Zampini     for (ke = 0; ke < ken; ke++) {
4208135c375SStefano Zampini       for (je = 0; je < jen; je++) {
4218135c375SStefano Zampini         for (ie = 0; ie < ien; ie++) {
4228135c375SStefano Zampini           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
4238135c375SStefano Zampini           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
4248135c375SStefano Zampini           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
4258135c375SStefano Zampini           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
4268135c375SStefano Zampini           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
4278135c375SStefano Zampini           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
4288135c375SStefano Zampini           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
4298135c375SStefano Zampini           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
4308135c375SStefano Zampini           mid    = 1; /* material id */
4318135c375SStefano Zampini           cid    = 5; /* hex */
4328135c375SStefano 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);
4338135c375SStefano Zampini         }
4348135c375SStefano Zampini       }
4358135c375SStefano Zampini     }
4368135c375SStefano Zampini     break;
4378135c375SStefano Zampini   default:
43898921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim);
4398135c375SStefano Zampini   }
4408135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr);
4418135c375SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr);
4428135c375SStefano Zampini 
44377eacf09SStefano Zampini   /* vertex coordinates */
44446900a7fSStefano Zampini   ierr = PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr);
445*7a8be351SBarry Smith   PetscCheck(xcoorl,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
4468135c375SStefano Zampini   ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
44777eacf09SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr);
44877eacf09SStefano Zampini   ierr = PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);CHKERRQ(ierr);
44921414b21SStefano Zampini   if (nel) {
45046900a7fSStefano Zampini     ierr = VecGetDM(xcoorl,&cda);CHKERRQ(ierr);
45177eacf09SStefano Zampini     ierr = VecGetArrayRead(xcoorl,&array);CHKERRQ(ierr);
4528f07809aSStefano Zampini     if (!cda) { /* HO viz */
4538f07809aSStefano Zampini       const char *fecname;
4548f07809aSStefano Zampini       PetscInt   nc,nl;
4558f07809aSStefano Zampini 
4568f07809aSStefano Zampini       ierr = PetscObjectGetName((PetscObject)xcoorl,&fecname);CHKERRQ(ierr);
45777eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr);
45877eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr);
45921414b21SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%s\n",fecname);CHKERRQ(ierr);
46077eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr);
46177eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/
4628f07809aSStefano Zampini       /* L2 coordinates */
4638f07809aSStefano Zampini       ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr);
4648f07809aSStefano Zampini       ierr = VecGetLocalSize(xcoorl,&nl);CHKERRQ(ierr);
4658f07809aSStefano Zampini       nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
46621414b21SStefano Zampini       cdof = nc ? nl/nc : 0;
4678f07809aSStefano Zampini       if (!ien) ien++;
4688f07809aSStefano Zampini       if (!jen) jen++;
4698f07809aSStefano Zampini       if (!ken) ken++;
4708f07809aSStefano Zampini       ist = jst = kst = 0;
4718f07809aSStefano Zampini       gm = ien;
4728f07809aSStefano Zampini       gn = jen;
4738f07809aSStefano Zampini       gp = ken;
47477eacf09SStefano Zampini     } else {
47546900a7fSStefano Zampini       DMDAGhostedGLVisViewerCtx *dactx;
47646900a7fSStefano Zampini 
4773ec1f749SStefano Zampini       ierr = DMGetApplicationContext(da,&dactx);CHKERRQ(ierr);
47877eacf09SStefano Zampini       ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr);
4798f07809aSStefano Zampini       cdof = sdim;
4808135c375SStefano Zampini       ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr);
4818135c375SStefano Zampini       ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr);
48246900a7fSStefano Zampini       if (dactx->ll) {
48346900a7fSStefano Zampini         kst = jst = ist = 0;
48446900a7fSStefano Zampini       } else {
4858135c375SStefano Zampini         kst  = gsz != sz ? 1 : 0;
4868135c375SStefano Zampini         jst  = gsy != sy ? 1 : 0;
4878135c375SStefano Zampini         ist  = gsx != sx ? 1 : 0;
4888f07809aSStefano Zampini       }
48946900a7fSStefano Zampini     }
4908135c375SStefano Zampini     for (ke = kst; ke < kst + ken; ke++) {
4918135c375SStefano Zampini       for (je = jst; je < jst + jen; je++) {
4928135c375SStefano Zampini         for (ie = ist; ie < ist + ien; ie++) {
49377eacf09SStefano Zampini           PetscInt c;
4948135c375SStefano Zampini 
4958135c375SStefano Zampini           i = ke * gm * gn + je * gm + ie;
49677eacf09SStefano Zampini           for (c=0;c<cdof/sdim;c++) {
49777eacf09SStefano Zampini             PetscInt d;
49877eacf09SStefano Zampini             for (d=0;d<sdim;d++) {
4998f07809aSStefano Zampini               ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]));CHKERRQ(ierr);
5008135c375SStefano Zampini             }
5018f07809aSStefano Zampini             ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
50277eacf09SStefano Zampini           }
5038135c375SStefano Zampini         }
5048135c375SStefano Zampini       }
5058135c375SStefano Zampini     }
5068135c375SStefano Zampini     ierr = VecRestoreArrayRead(xcoorl,&array);CHKERRQ(ierr);
50721414b21SStefano Zampini   }
5088135c375SStefano Zampini   PetscFunctionReturn(0);
5098135c375SStefano Zampini }
5108135c375SStefano Zampini 
5110286d493SLisandro Dalcin PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
5128135c375SStefano Zampini {
5138135c375SStefano Zampini   PetscErrorCode ierr;
5148135c375SStefano Zampini   PetscFunctionBegin;
5150286d493SLisandro Dalcin   ierr = DMView_GLVis(dm,viewer,DMDAView_GLVis_ASCII);CHKERRQ(ierr);
5168135c375SStefano Zampini   PetscFunctionReturn(0);
5178135c375SStefano Zampini }
518