xref: /petsc/src/dm/impls/da/grglvis.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(PetscFree(*vctx));
1446900a7fSStefano Zampini   PetscFunctionReturn(0);
1546900a7fSStefano Zampini }
1646900a7fSStefano Zampini 
1746900a7fSStefano Zampini typedef struct {
1846900a7fSStefano Zampini   Vec xlocal;
1946900a7fSStefano Zampini } DMDAFieldGLVisViewerCtx;
2046900a7fSStefano Zampini 
2146900a7fSStefano Zampini static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx)
2246900a7fSStefano Zampini {
2346900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
248135c375SStefano Zampini 
258135c375SStefano Zampini   PetscFunctionBegin;
269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->xlocal));
279566063dSJacob Faibussowitsch   PetscCall(PetscFree(vctx));
288135c375SStefano Zampini   PetscFunctionReturn(0);
298135c375SStefano Zampini }
308135c375SStefano Zampini 
3146900a7fSStefano Zampini /*
3246900a7fSStefano Zampini    dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right
3346900a7fSStefano Zampini    dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left
3446900a7fSStefano Zampini */
358135c375SStefano Zampini static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
368135c375SStefano Zampini {
3746900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
3846900a7fSStefano Zampini   PetscInt                  sx,sy,sz,ien,jen,ken;
398135c375SStefano Zampini 
408135c375SStefano Zampini   PetscFunctionBegin;
41aab5bcd8SJed Brown   /* Appease -Wmaybe-uninitialized */
42aab5bcd8SJed Brown   if (nex) *nex = -1;
43aab5bcd8SJed Brown   if (ney) *ney = -1;
44aab5bcd8SJed Brown   if (nez) *nez = -1;
459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&sx,&sy,&sz,&ien,&jen,&ken));
469566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da,&dactx));
4746900a7fSStefano Zampini   if (dactx->ll) {
4846900a7fSStefano Zampini     PetscInt dim;
4946900a7fSStefano Zampini 
509566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(da,&dim));
5146900a7fSStefano Zampini     if (!sx) ien--;
5246900a7fSStefano Zampini     if (!sy && dim > 1) jen--;
5346900a7fSStefano Zampini     if (!sz && dim > 2) ken--;
5446900a7fSStefano Zampini   } else {
5546900a7fSStefano Zampini     PetscInt M,N,P;
5646900a7fSStefano Zampini 
579566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
588135c375SStefano Zampini     if (sx+ien == M) ien--;
598135c375SStefano Zampini     if (sy+jen == N) jen--;
608135c375SStefano Zampini     if (sz+ken == P) ken--;
6146900a7fSStefano Zampini   }
628135c375SStefano Zampini   if (nex) *nex = ien;
638135c375SStefano Zampini   if (ney) *ney = jen;
648135c375SStefano Zampini   if (nez) *nez = ken;
658135c375SStefano Zampini   PetscFunctionReturn(0);
668135c375SStefano Zampini }
678135c375SStefano Zampini 
688135c375SStefano Zampini /* inherits number of vertices from DMDAGetNumElementsGhosted */
698135c375SStefano Zampini static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
708135c375SStefano Zampini {
7121414b21SStefano Zampini   PetscInt       ien = 0,jen = 0,ken = 0,dim;
7221414b21SStefano Zampini   PetscInt       tote;
738135c375SStefano Zampini 
748135c375SStefano Zampini   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da,&dim));
769566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumElementsGhosted(da,&ien,&jen,&ken));
7721414b21SStefano Zampini   tote = ien * (dim > 1 ? jen : 1) * (dim > 2 ? ken : 1);
7821414b21SStefano Zampini   if (tote) {
798135c375SStefano Zampini     ien = ien+1;
808135c375SStefano Zampini     jen = dim > 1 ? jen+1 : 1;
818135c375SStefano Zampini     ken = dim > 2 ? ken+1 : 1;
8221414b21SStefano Zampini   }
838135c375SStefano Zampini   if (nvx) *nvx = ien;
848135c375SStefano Zampini   if (nvy) *nvy = jen;
858135c375SStefano Zampini   if (nvz) *nvz = ken;
868135c375SStefano Zampini   PetscFunctionReturn(0);
878135c375SStefano Zampini }
888135c375SStefano Zampini 
898135c375SStefano Zampini static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
908135c375SStefano Zampini {
918135c375SStefano Zampini   DM                        da;
9246900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx   *ctx = (DMDAFieldGLVisViewerCtx*)vctx;
9346900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
948135c375SStefano Zampini   const PetscScalar         *array;
958135c375SStefano Zampini   PetscScalar               **arrayf;
968135c375SStefano Zampini   PetscInt                  i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss;
978135c375SStefano Zampini   PetscInt                  sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp;
988135c375SStefano Zampini 
998135c375SStefano Zampini   PetscFunctionBegin;
1009566063dSJacob Faibussowitsch   PetscCall(VecGetDM(ctx->xlocal,&da));
1017a8be351SBarry Smith   PetscCheck(da,PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
1029566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da,&dactx));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(ctx->xlocal,&bs));
1049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal));
1059566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal));
1069566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken));
1079566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp));
1089566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL));
10946900a7fSStefano Zampini   if (dactx->ll) {
11046900a7fSStefano Zampini     kst = jst = ist = 0;
11146900a7fSStefano Zampini   } else {
1128135c375SStefano Zampini     kst  = gsz != sz ? 1 : 0;
1138135c375SStefano Zampini     jst  = gsy != sy ? 1 : 0;
1148135c375SStefano Zampini     ist  = gsx != sx ? 1 : 0;
11546900a7fSStefano Zampini   }
1169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nf,&arrayf,nf,&bss));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->xlocal,&array));
1188135c375SStefano Zampini   for (f=0;f<nf;f++) {
1199566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize((Vec)oXf[f],&bss[f]));
1209566063dSJacob Faibussowitsch     PetscCall(VecGetArray((Vec)oXf[f],&arrayf[f]));
1218135c375SStefano Zampini   }
1228135c375SStefano Zampini   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
1238135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
1248135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
1258135c375SStefano Zampini         PetscInt cf,b;
1268135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
1278135c375SStefano Zampini         for (f=0,cf=0;f<nf;f++)
1288135c375SStefano Zampini           for (b=0;b<bss[f];b++)
1298135c375SStefano Zampini             arrayf[f][bss[f]*ii+b] = array[i*bs+cf++];
1308135c375SStefano Zampini         ii++;
1318135c375SStefano Zampini       }
1328135c375SStefano Zampini     }
1338135c375SStefano Zampini   }
1349566063dSJacob Faibussowitsch   for (f=0;f<nf;f++) PetscCall(VecRestoreArray((Vec)oXf[f],&arrayf[f]));
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->xlocal,&array));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFree2(arrayf,bss));
1378135c375SStefano Zampini   PetscFunctionReturn(0);
1388135c375SStefano Zampini }
1398135c375SStefano Zampini 
1408135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
1418135c375SStefano Zampini {
14246900a7fSStefano Zampini   DM             da = (DM)oda,daview;
1438135c375SStefano Zampini 
1448135c375SStefano Zampini   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery(oda,"GLVisGraphicsDMDAGhosted",(PetscObject*)&daview));
14646900a7fSStefano Zampini   if (!daview) {
14746900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
14846900a7fSStefano Zampini     DM                        dacoord = NULL;
14946900a7fSStefano Zampini     Vec                       xcoor,xcoorl;
15046900a7fSStefano Zampini     PetscBool                 hashocoord = PETSC_FALSE;
15146900a7fSStefano Zampini     const PetscInt            *lx,*ly,*lz;
15246900a7fSStefano Zampini     PetscInt                  dim,M,N,P,m,n,p,dof,s,i;
15346900a7fSStefano Zampini 
1549566063dSJacob Faibussowitsch     PetscCall(PetscNew(&dactx));
155a7bda6aaSStefano Zampini     dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
156d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");
1579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-viewer_glvis_dm_da_ll","Left-looking subdomain view",NULL,dactx->ll,&dactx->ll,NULL));
158d0609cedSBarry Smith     PetscOptionsEnd();
1598135c375SStefano Zampini     /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1609566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL));
1619566063dSJacob Faibussowitsch     PetscCall(DMDAGetOwnershipRanges(da,&lx,&ly,&lz));
1629566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor));
16377eacf09SStefano Zampini     if (!xcoor) {
1649566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(da,&xcoor));
16577eacf09SStefano Zampini     } else {
1668f07809aSStefano Zampini       hashocoord = PETSC_TRUE;
16777eacf09SStefano Zampini     }
1689566063dSJacob Faibussowitsch     PetscCall(PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n"));
1698135c375SStefano Zampini     switch (dim) {
1708135c375SStefano Zampini     case 1:
1719566063dSJacob Faibussowitsch       PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview));
1728f07809aSStefano Zampini       if (!hashocoord) {
1739566063dSJacob Faibussowitsch         PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord));
1748f07809aSStefano Zampini       }
1758135c375SStefano Zampini       break;
1768135c375SStefano Zampini     case 2:
1779566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,lx,ly,&daview));
1788f07809aSStefano Zampini       if (!hashocoord) {
1799566063dSJacob Faibussowitsch         PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,2,1,lx,ly,&dacoord));
1808f07809aSStefano Zampini       }
1818135c375SStefano Zampini       break;
1828135c375SStefano Zampini     case 3:
1839566063dSJacob Faibussowitsch       PetscCall(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));
1848f07809aSStefano Zampini       if (!hashocoord) {
1859566063dSJacob Faibussowitsch         PetscCall(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));
1868f07809aSStefano Zampini       }
1878135c375SStefano Zampini       break;
1888135c375SStefano Zampini     default:
18963a3b9bcSJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %" PetscInt_FMT,dim);
1908135c375SStefano Zampini     }
1919566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContext(daview,dactx));
1929566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private));
1939566063dSJacob Faibussowitsch     PetscCall(DMSetUp(daview));
1948135c375SStefano Zampini     if (!xcoor) {
1959566063dSJacob Faibussowitsch       PetscCall(DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0));
1969566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(daview,&xcoor));
1978135c375SStefano Zampini     }
1988f07809aSStefano Zampini     if (dacoord) {
1999566063dSJacob Faibussowitsch       PetscCall(DMSetUp(dacoord));
2009566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(dacoord,&xcoorl));
2019566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl));
2029566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl));
2039566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dacoord));
2048f07809aSStefano Zampini     } else {
2058f07809aSStefano Zampini       PetscInt   ien,jen,ken,nc,nl,cdof,deg;
2063e7cab22SLisandro Dalcin       char       fecmesh[64];
20721414b21SStefano Zampini       const char *name;
20821414b21SStefano Zampini       PetscBool  flg;
2098135c375SStefano Zampini 
2109566063dSJacob Faibussowitsch       PetscCall(DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken));
2118f07809aSStefano Zampini       nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
2128f07809aSStefano Zampini 
2139566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xcoor,&nl));
214*1dca8a05SBarry Smith       PetscCheck(!nc || (nl % nc) == 0,PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %" PetscInt_FMT " and number of cells %" PetscInt_FMT,nl,nc);
2159566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(xcoor,&xcoorl));
2169566063dSJacob Faibussowitsch       PetscCall(VecCopy(xcoor,xcoorl));
2179566063dSJacob Faibussowitsch       PetscCall(VecSetDM(xcoorl,NULL));
2189566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)xcoor,&name));
2199566063dSJacob Faibussowitsch       PetscCall(PetscStrbeginswith(name,"FiniteElementCollection:",&flg));
22021414b21SStefano Zampini       if (!flg) {
22121414b21SStefano Zampini         deg = 0;
22221414b21SStefano Zampini         if (nc && nl) {
2238f07809aSStefano Zampini           cdof = nl/(nc*dim);
2248f07809aSStefano Zampini           deg  = 1;
2258f07809aSStefano Zampini           while (1) {
2268f07809aSStefano Zampini             PetscInt degd = 1;
2278f07809aSStefano Zampini             for (i=0;i<dim;i++) degd *= (deg+1);
22863a3b9bcSJacob Faibussowitsch             PetscCheck(degd <= cdof,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %" PetscInt_FMT,cdof);
2298f07809aSStefano Zampini             if (degd == cdof) break;
2308f07809aSStefano Zampini             deg++;
2318f07809aSStefano Zampini           }
23221414b21SStefano Zampini         }
23363a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fecmesh,sizeof(fecmesh),"FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT,dim,deg));
2349566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)xcoorl,fecmesh));
23521414b21SStefano Zampini       } else {
2369566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)xcoorl,name));
23721414b21SStefano Zampini       }
2388f07809aSStefano Zampini     }
2398135c375SStefano Zampini 
24046900a7fSStefano Zampini     /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
2419566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl));
2429566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
2438135c375SStefano Zampini 
24477eacf09SStefano Zampini     /* daview is composed with the original DMDA */
2459566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview));
2469566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)daview));
24746900a7fSStefano Zampini   }
24877eacf09SStefano Zampini 
249c9493c35SStefano Zampini   /* customize the viewer if present */
250c9493c35SStefano Zampini   if (viewer) {
25146900a7fSStefano Zampini     DMDAFieldGLVisViewerCtx   *ctx;
25246900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
25346900a7fSStefano Zampini     char                      fec[64];
254ab166b77SStefano Zampini     Vec                       xlocal,*Ufield;
25546900a7fSStefano Zampini     const char                **dafieldname;
25646900a7fSStefano Zampini     char                      **fec_type,**fieldname;
25746900a7fSStefano Zampini     PetscInt                  *nlocal,*bss,*dims;
25846900a7fSStefano Zampini     PetscInt                  dim,M,N,P,dof,s,i,nf;
25946900a7fSStefano Zampini     PetscBool                 bsset;
2608f07809aSStefano Zampini 
2619566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL));
2629566063dSJacob Faibussowitsch     PetscCall(DMGetApplicationContext(daview,&dactx));
2639566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(daview,&xlocal));
2649566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldNames(da,(const char * const **)&dafieldname));
2659566063dSJacob Faibussowitsch     PetscCall(DMDAGetNumVerticesGhosted(daview,&M,&N,&P));
26663a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%" PetscInt_FMT "D_P1",dim));
2679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc6(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname,dof,&Ufield));
2688135c375SStefano Zampini     for (i=0;i<dof;i++) bss[i] = 1;
2698135c375SStefano Zampini     nf = dof;
2708135c375SStefano Zampini 
271d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");
2729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset));
273d0609cedSBarry Smith     PetscOptionsEnd();
2748135c375SStefano Zampini     if (bsset) {
2758135c375SStefano Zampini       PetscInt t;
2768135c375SStefano Zampini       for (i=0,t=0;i<nf;i++) t += bss[i];
27763a3b9bcSJacob Faibussowitsch       PetscCheck(t == dof,PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %" PetscInt_FMT " should equal %" PetscInt_FMT,t,dof);
2788135c375SStefano Zampini     } else nf = dof;
2798135c375SStefano Zampini 
2808135c375SStefano Zampini     for (i=0,s=0;i<nf;i++) {
2819566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec,&fec_type[i]));
2828135c375SStefano Zampini       if (bss[i] == 1) {
2839566063dSJacob Faibussowitsch         PetscCall(PetscStrallocpy(dafieldname[s],&fieldname[i]));
2848135c375SStefano Zampini       } else {
2858135c375SStefano Zampini         PetscInt b;
2868135c375SStefano Zampini         size_t tlen = 9; /* "Vector-" + end */
2878135c375SStefano Zampini         for (b=0;b<bss[i];b++) {
2888135c375SStefano Zampini           size_t len;
2899566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(dafieldname[s+b],&len));
2908135c375SStefano Zampini           tlen += len + 1; /* field + "-" */
2918135c375SStefano Zampini         }
2929566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(tlen,&fieldname[i]));
2939566063dSJacob Faibussowitsch         PetscCall(PetscStrcpy(fieldname[i],"Vector-"));
2948135c375SStefano Zampini         for (b=0;b<bss[i]-1;b++) {
2959566063dSJacob Faibussowitsch           PetscCall(PetscStrcat(fieldname[i],dafieldname[s+b]));
2969566063dSJacob Faibussowitsch           PetscCall(PetscStrcat(fieldname[i],"-"));
2978135c375SStefano Zampini         }
2989566063dSJacob Faibussowitsch         PetscCall(PetscStrcat(fieldname[i],dafieldname[s+b]));
2998135c375SStefano Zampini       }
300bb77a09fSStefano Zampini       dims[i] = dim;
3018135c375SStefano Zampini       nlocal[i] = M*N*P*bss[i];
3028135c375SStefano Zampini       s += bss[i];
3038135c375SStefano Zampini     }
3048135c375SStefano Zampini 
30546900a7fSStefano Zampini     /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
3069566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
3078135c375SStefano Zampini     ctx->xlocal = xlocal;
3088135c375SStefano Zampini 
3094cac2994SStefano Zampini     /* create work vectors */
3104cac2994SStefano Zampini     for (i=0;i<nf;i++) {
3119566063dSJacob Faibussowitsch       PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da),nlocal[i],PETSC_DECIDE,&Ufield[i]));
3129566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)Ufield[i],fieldname[i]));
3139566063dSJacob Faibussowitsch       PetscCall(VecSetBlockSize(Ufield[i],bss[i]));
3149566063dSJacob Faibussowitsch       PetscCall(VecSetDM(Ufield[i],da));
3154cac2994SStefano Zampini     }
3164cac2994SStefano Zampini 
3179566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetFields(viewer,nf,(const char**)fec_type,dims,DMDASampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DMDAFieldDestroyGLVisViewerCtx_Private));
3188135c375SStefano Zampini     for (i=0;i<nf;i++) {
3199566063dSJacob Faibussowitsch       PetscCall(PetscFree(fec_type[i]));
3209566063dSJacob Faibussowitsch       PetscCall(PetscFree(fieldname[i]));
3219566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&Ufield[i]));
3228135c375SStefano Zampini     }
3239566063dSJacob Faibussowitsch     PetscCall(PetscFree6(fec_type,nlocal,bss,dims,fieldname,Ufield));
324c9493c35SStefano Zampini   }
3258135c375SStefano Zampini   PetscFunctionReturn(0);
3268135c375SStefano Zampini }
3278135c375SStefano Zampini 
3288135c375SStefano Zampini static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
3298135c375SStefano Zampini {
33077eacf09SStefano Zampini   DM                da,cda;
3318135c375SStefano Zampini   Vec               xcoorl;
3323924b612SStefano Zampini   PetscMPIInt       size;
3338135c375SStefano Zampini   const PetscScalar *array;
3348135c375SStefano Zampini   PetscContainer    glvis_container;
33577eacf09SStefano Zampini   PetscInt          dim,sdim,i,vid[8],mid,cid,cdof;
33621414b21SStefano Zampini   PetscInt          sx,sy,sz,ie,je,ke,ien,jen,ken,nel;
3378135c375SStefano Zampini   PetscInt          gsx,gsy,gsz,gm,gn,gp,kst,jst,ist;
3388135c375SStefano Zampini   PetscBool         enabled = PETSC_TRUE, isascii;
3398f07809aSStefano Zampini   const char        *fmt;
3408135c375SStefano Zampini 
3418135c375SStefano Zampini   PetscFunctionBegin;
3428135c375SStefano Zampini   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3438135c375SStefano Zampini   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
3457a8be351SBarry Smith   PetscCheck(isascii,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII");
3469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size));
34708401ef6SPierre Jolivet   PetscCheck(size <= 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization");
3489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
3498135c375SStefano Zampini 
3508135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
3519566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container));
3527a8be351SBarry Smith   PetscCheck(glvis_container,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container");
3538135c375SStefano Zampini   {
3548135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
3559566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(glvis_container,(void**)&glvis_info));
3568135c375SStefano Zampini     enabled = glvis_info->enabled;
35777eacf09SStefano Zampini     fmt = glvis_info->fmt;
3588135c375SStefano Zampini   }
359c9493c35SStefano Zampini   /* this can happen if we are calling DMView outside of VecView_GLVis */
3609566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da));
3619566063dSJacob Faibussowitsch   if (!da) PetscCall(DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL));
3629566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da));
3637a8be351SBarry Smith   PetscCheck(da,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA");
3649566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(da,&sdim));
3658135c375SStefano Zampini 
3669566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.0\n"));
3679566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\ndimension\n"));
36863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",dim));
3698135c375SStefano Zampini 
3708135c375SStefano Zampini   if (!enabled) {
3719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nelements\n"));
37263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
3739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nboundary\n"));
37463a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
3759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
37663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
37763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",sdim));
3788135c375SStefano Zampini     PetscFunctionReturn(0);
3798135c375SStefano Zampini   }
3808135c375SStefano Zampini 
3819566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumElementsGhosted(da,&ien,&jen,&ken));
38221414b21SStefano Zampini   nel  = ien;
38321414b21SStefano Zampini   if (dim > 1) nel *= jen;
38421414b21SStefano Zampini   if (dim > 2) nel *= ken;
3859566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\nelements\n"));
38663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",nel));
3878135c375SStefano Zampini   switch (dim) {
3888135c375SStefano Zampini   case 1:
3898135c375SStefano Zampini     for (ie = 0; ie < ien; ie++) {
3908135c375SStefano Zampini       vid[0] = ie;
3918135c375SStefano Zampini       vid[1] = ie+1;
3928135c375SStefano Zampini       mid    = 1; /* material id */
3938135c375SStefano Zampini       cid    = 1; /* segment */
39463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",mid,cid,vid[0],vid[1]));
3958135c375SStefano Zampini     }
3968135c375SStefano Zampini     break;
3978135c375SStefano Zampini   case 2:
3988135c375SStefano Zampini     for (je = 0; je < jen; je++) {
3998135c375SStefano Zampini       for (ie = 0; ie < ien; ie++) {
4008135c375SStefano Zampini         vid[0] =     je*(ien+1) + ie;
4018135c375SStefano Zampini         vid[1] =     je*(ien+1) + ie+1;
4028135c375SStefano Zampini         vid[2] = (je+1)*(ien+1) + ie+1;
4038135c375SStefano Zampini         vid[3] = (je+1)*(ien+1) + ie;
4048135c375SStefano Zampini         mid    = 1; /* material id */
4058135c375SStefano Zampini         cid    = 3; /* quad */
40663a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",mid,cid,vid[0],vid[1],vid[2],vid[3]));
4078135c375SStefano Zampini       }
4088135c375SStefano Zampini     }
4098135c375SStefano Zampini     break;
4108135c375SStefano Zampini   case 3:
4118135c375SStefano Zampini     for (ke = 0; ke < ken; ke++) {
4128135c375SStefano Zampini       for (je = 0; je < jen; je++) {
4138135c375SStefano Zampini         for (ie = 0; ie < ien; ie++) {
4148135c375SStefano Zampini           vid[0] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie;
4158135c375SStefano Zampini           vid[1] =     ke*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
4168135c375SStefano Zampini           vid[2] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
4178135c375SStefano Zampini           vid[3] =     ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
4188135c375SStefano Zampini           vid[4] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie;
4198135c375SStefano Zampini           vid[5] = (ke+1)*(jen+1)*(ien+1) +     je*(ien+1) + ie+1;
4208135c375SStefano Zampini           vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1;
4218135c375SStefano Zampini           vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie;
4228135c375SStefano Zampini           mid    = 1; /* material id */
4238135c375SStefano Zampini           cid    = 5; /* hex */
42463a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",mid,cid,vid[0],vid[1],vid[2],vid[3],vid[4],vid[5],vid[6],vid[7]));
4258135c375SStefano Zampini         }
4268135c375SStefano Zampini       }
4278135c375SStefano Zampini     }
4288135c375SStefano Zampini     break;
4298135c375SStefano Zampini   default:
43063a3b9bcSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %" PetscInt_FMT,dim);
4318135c375SStefano Zampini   }
4329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\nboundary\n"));
43363a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"0\n"));
4348135c375SStefano Zampini 
43577eacf09SStefano Zampini   /* vertex coordinates */
4369566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl));
4377a8be351SBarry Smith   PetscCheck(xcoorl,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords");
4389566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken));
4399566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n"));
44063a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",ien*jen*ken));
44121414b21SStefano Zampini   if (nel) {
4429566063dSJacob Faibussowitsch     PetscCall(VecGetDM(xcoorl,&cda));
4439566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xcoorl,&array));
4448f07809aSStefano Zampini     if (!cda) { /* HO viz */
4458f07809aSStefano Zampini       const char *fecname;
4468f07809aSStefano Zampini       PetscInt   nc,nl;
4478f07809aSStefano Zampini 
4489566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)xcoorl,&fecname));
4499566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"nodes\n"));
4509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n"));
4519566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",fecname));
45263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"VDim: %" PetscInt_FMT "\n",sdim));
4539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n")); /*Ordering::byVDIM*/
4548f07809aSStefano Zampini       /* L2 coordinates */
4559566063dSJacob Faibussowitsch       PetscCall(DMDAGetNumElementsGhosted(da,&ien,&jen,&ken));
4569566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xcoorl,&nl));
4578f07809aSStefano Zampini       nc   = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1);
45821414b21SStefano Zampini       cdof = nc ? nl/nc : 0;
4598f07809aSStefano Zampini       if (!ien) ien++;
4608f07809aSStefano Zampini       if (!jen) jen++;
4618f07809aSStefano Zampini       if (!ken) ken++;
4628f07809aSStefano Zampini       ist = jst = kst = 0;
4638f07809aSStefano Zampini       gm = ien;
4648f07809aSStefano Zampini       gn = jen;
4658f07809aSStefano Zampini       gp = ken;
46677eacf09SStefano Zampini     } else {
46746900a7fSStefano Zampini       DMDAGhostedGLVisViewerCtx *dactx;
46846900a7fSStefano Zampini 
4699566063dSJacob Faibussowitsch       PetscCall(DMGetApplicationContext(da,&dactx));
47063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT "\n",sdim));
4718f07809aSStefano Zampini       cdof = sdim;
4729566063dSJacob Faibussowitsch       PetscCall(DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL));
4739566063dSJacob Faibussowitsch       PetscCall(DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp));
47446900a7fSStefano Zampini       if (dactx->ll) {
47546900a7fSStefano Zampini         kst = jst = ist = 0;
47646900a7fSStefano Zampini       } else {
4778135c375SStefano Zampini         kst  = gsz != sz ? 1 : 0;
4788135c375SStefano Zampini         jst  = gsy != sy ? 1 : 0;
4798135c375SStefano Zampini         ist  = gsx != sx ? 1 : 0;
4808f07809aSStefano Zampini       }
48146900a7fSStefano Zampini     }
4828135c375SStefano Zampini     for (ke = kst; ke < kst + ken; ke++) {
4838135c375SStefano Zampini       for (je = jst; je < jst + jen; je++) {
4848135c375SStefano Zampini         for (ie = ist; ie < ist + ien; ie++) {
48577eacf09SStefano Zampini           PetscInt c;
4868135c375SStefano Zampini 
4878135c375SStefano Zampini           i = ke * gm * gn + je * gm + ie;
48877eacf09SStefano Zampini           for (c=0;c<cdof/sdim;c++) {
48977eacf09SStefano Zampini             PetscInt d;
49077eacf09SStefano Zampini             for (d=0;d<sdim;d++) {
4919566063dSJacob Faibussowitsch               PetscCall(PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d])));
4928135c375SStefano Zampini             }
4939566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
49477eacf09SStefano Zampini           }
4958135c375SStefano Zampini         }
4968135c375SStefano Zampini       }
4978135c375SStefano Zampini     }
4989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xcoorl,&array));
49921414b21SStefano Zampini   }
5008135c375SStefano Zampini   PetscFunctionReturn(0);
5018135c375SStefano Zampini }
5028135c375SStefano Zampini 
5030286d493SLisandro Dalcin PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
5048135c375SStefano Zampini {
5058135c375SStefano Zampini   PetscFunctionBegin;
5069566063dSJacob Faibussowitsch   PetscCall(DMView_GLVis(dm,viewer,DMDAView_GLVis_ASCII));
5078135c375SStefano Zampini   PetscFunctionReturn(0);
5088135c375SStefano Zampini }
509