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 { 7*46900a7fSStefano Zampini PetscBool ll; 8*46900a7fSStefano Zampini } DMDAGhostedGLVisViewerCtx; 98135c375SStefano Zampini 10*46900a7fSStefano Zampini static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(void **vctx) 118135c375SStefano Zampini { 12*46900a7fSStefano Zampini PetscErrorCode ierr; 13*46900a7fSStefano Zampini 14*46900a7fSStefano Zampini PetscFunctionBegin; 15*46900a7fSStefano Zampini ierr = PetscFree(*vctx);CHKERRQ(ierr); 16*46900a7fSStefano Zampini PetscFunctionReturn(0); 17*46900a7fSStefano Zampini } 18*46900a7fSStefano Zampini 19*46900a7fSStefano Zampini typedef struct { 20*46900a7fSStefano Zampini Vec xlocal; 21*46900a7fSStefano Zampini } DMDAFieldGLVisViewerCtx; 22*46900a7fSStefano Zampini 23*46900a7fSStefano Zampini static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx) 24*46900a7fSStefano Zampini { 25*46900a7fSStefano 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 34*46900a7fSStefano Zampini /* 35*46900a7fSStefano Zampini dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right 36*46900a7fSStefano Zampini dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left 37*46900a7fSStefano Zampini */ 388135c375SStefano Zampini static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez) 398135c375SStefano Zampini { 40*46900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 41*46900a7fSStefano 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); 50*46900a7fSStefano Zampini ierr = DMGetApplicationContext(da,(void**)&dactx);CHKERRQ(ierr); 51*46900a7fSStefano Zampini if (dactx->ll) { 52*46900a7fSStefano Zampini PetscInt dim; 53*46900a7fSStefano Zampini 54*46900a7fSStefano Zampini ierr = DMGetDimension(da,&dim);CHKERRQ(ierr); 55*46900a7fSStefano Zampini if (!sx) ien--; 56*46900a7fSStefano Zampini if (!sy && dim > 1) jen--; 57*46900a7fSStefano Zampini if (!sz && dim > 2) ken--; 58*46900a7fSStefano Zampini } else { 59*46900a7fSStefano Zampini PetscInt M,N,P; 60*46900a7fSStefano Zampini 61*46900a7fSStefano 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--; 65*46900a7fSStefano 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 { 758135c375SStefano Zampini PetscInt ien,jen,ken,dim; 768135c375SStefano Zampini PetscErrorCode ierr; 778135c375SStefano Zampini 788135c375SStefano Zampini PetscFunctionBegin; 798135c375SStefano Zampini ierr = DMGetDimension(da,&dim);CHKERRQ(ierr); 808135c375SStefano Zampini ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr); 818135c375SStefano Zampini ien = ien+1; 828135c375SStefano Zampini jen = dim > 1 ? jen+1 : 1; 838135c375SStefano Zampini ken = dim > 2 ? ken+1 : 1; 848135c375SStefano Zampini if (nvx) *nvx = ien; 858135c375SStefano Zampini if (nvy) *nvy = jen; 868135c375SStefano Zampini if (nvz) *nvz = ken; 878135c375SStefano Zampini PetscFunctionReturn(0); 888135c375SStefano Zampini } 898135c375SStefano Zampini 908135c375SStefano Zampini static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx) 918135c375SStefano Zampini { 928135c375SStefano Zampini DM da; 93*46900a7fSStefano Zampini DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx*)vctx; 94*46900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 958135c375SStefano Zampini const PetscScalar *array; 968135c375SStefano Zampini PetscScalar **arrayf; 978135c375SStefano Zampini PetscInt i,f,ii,ien,jen,ken,ie,je,ke,bs,*bss; 988135c375SStefano Zampini PetscInt sx,sy,sz,gsx,gsy,gsz,ist,jst,kst,gm,gn,gp; 998135c375SStefano Zampini PetscErrorCode ierr; 1008135c375SStefano Zampini 1018135c375SStefano Zampini PetscFunctionBegin; 1028135c375SStefano Zampini ierr = VecGetDM(ctx->xlocal,&da);CHKERRQ(ierr); 1038135c375SStefano Zampini if (!da) SETERRQ(PetscObjectComm(oX),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 104*46900a7fSStefano Zampini ierr = DMGetApplicationContext(da,(void**)&dactx);CHKERRQ(ierr); 1058135c375SStefano Zampini ierr = VecGetBlockSize(ctx->xlocal,&bs);CHKERRQ(ierr); 1068135c375SStefano Zampini ierr = DMGlobalToLocalBegin(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr); 1078135c375SStefano Zampini ierr = DMGlobalToLocalEnd(da,(Vec)oX,INSERT_VALUES,ctx->xlocal);CHKERRQ(ierr); 1088135c375SStefano Zampini ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr); 1098135c375SStefano Zampini ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr); 1108135c375SStefano Zampini ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr); 111*46900a7fSStefano Zampini if (dactx->ll) { 112*46900a7fSStefano Zampini kst = jst = ist = 0; 113*46900a7fSStefano Zampini } else { 1148135c375SStefano Zampini kst = gsz != sz ? 1 : 0; 1158135c375SStefano Zampini jst = gsy != sy ? 1 : 0; 1168135c375SStefano Zampini ist = gsx != sx ? 1 : 0; 117*46900a7fSStefano Zampini } 1188135c375SStefano Zampini ierr = PetscMalloc2(nf,&arrayf,nf,&bss);CHKERRQ(ierr); 1198135c375SStefano Zampini ierr = VecGetArrayRead(ctx->xlocal,&array);CHKERRQ(ierr); 1208135c375SStefano Zampini for (f=0;f<nf;f++) { 1218135c375SStefano Zampini ierr = VecGetBlockSize((Vec)oXf[f],&bss[f]);CHKERRQ(ierr); 1228135c375SStefano Zampini ierr = VecGetArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr); 1238135c375SStefano Zampini } 1248135c375SStefano Zampini for (ke = kst, ii = 0; ke < kst + ken; ke++) { 1258135c375SStefano Zampini for (je = jst; je < jst + jen; je++) { 1268135c375SStefano Zampini for (ie = ist; ie < ist + ien; ie++) { 1278135c375SStefano Zampini PetscInt cf,b; 1288135c375SStefano Zampini i = ke * gm * gn + je * gm + ie; 1298135c375SStefano Zampini for (f=0,cf=0;f<nf;f++) 1308135c375SStefano Zampini for (b=0;b<bss[f];b++) 1318135c375SStefano Zampini arrayf[f][bss[f]*ii+b] = array[i*bs+cf++]; 1328135c375SStefano Zampini ii++; 1338135c375SStefano Zampini } 1348135c375SStefano Zampini } 1358135c375SStefano Zampini } 1368135c375SStefano Zampini for (f=0;f<nf;f++) { ierr = VecRestoreArray((Vec)oXf[f],&arrayf[f]);CHKERRQ(ierr); } 1378135c375SStefano Zampini ierr = VecRestoreArrayRead(ctx->xlocal,&array);CHKERRQ(ierr); 1388135c375SStefano Zampini ierr = PetscFree2(arrayf,bss);CHKERRQ(ierr); 1398135c375SStefano Zampini PetscFunctionReturn(0); 1408135c375SStefano Zampini } 1418135c375SStefano Zampini 1428135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer) 1438135c375SStefano Zampini { 144*46900a7fSStefano Zampini DM da = (DM)oda,daview; 1458135c375SStefano Zampini PetscErrorCode ierr; 1468135c375SStefano Zampini 1478135c375SStefano Zampini PetscFunctionBegin; 148*46900a7fSStefano Zampini ierr = PetscObjectQuery(oda,"GLVisGraphicsDMDAGhosted",(PetscObject*)&daview);CHKERRQ(ierr); 149*46900a7fSStefano Zampini if (!daview) { 150*46900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 151*46900a7fSStefano Zampini DM dacoord = NULL; 152*46900a7fSStefano Zampini Vec xcoor,xcoorl; 153*46900a7fSStefano Zampini PetscBool hashocoord = PETSC_FALSE; 154*46900a7fSStefano Zampini const PetscInt *lx,*ly,*lz; 155*46900a7fSStefano Zampini PetscInt dim,M,N,P,m,n,p,dof,s,i; 156*46900a7fSStefano Zampini 157*46900a7fSStefano Zampini ierr = PetscNew(&dactx);CHKERRQ(ierr); 158*46900a7fSStefano Zampini ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");CHKERRQ(ierr); 159*46900a7fSStefano Zampini ierr = PetscOptionsBool("-viewer_glvis_dm_da_ll","Left-looking subdomain view",NULL,dactx->ll,&dactx->ll,NULL);CHKERRQ(ierr); 160*46900a7fSStefano Zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 1618135c375SStefano Zampini /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */ 1628135c375SStefano Zampini ierr = DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 1638135c375SStefano Zampini ierr = DMDAGetOwnershipRanges(da,&lx,&ly,&lz);CHKERRQ(ierr); 16477eacf09SStefano Zampini ierr = PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor);CHKERRQ(ierr); 16577eacf09SStefano Zampini if (!xcoor) { 1668135c375SStefano Zampini ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 16777eacf09SStefano Zampini } else { 1688f07809aSStefano Zampini hashocoord = PETSC_TRUE; 16977eacf09SStefano Zampini } 1708135c375SStefano Zampini ierr = PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n");CHKERRQ(ierr); 1718135c375SStefano Zampini switch (dim) { 1728135c375SStefano Zampini case 1: 1738135c375SStefano Zampini ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview);CHKERRQ(ierr); 1748f07809aSStefano Zampini if (!hashocoord) { 1758f07809aSStefano Zampini ierr = DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord);CHKERRQ(ierr); 1768f07809aSStefano Zampini } 1778135c375SStefano Zampini break; 1788135c375SStefano Zampini case 2: 1798135c375SStefano 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); 1808f07809aSStefano Zampini if (!hashocoord) { 1818f07809aSStefano 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); 1828f07809aSStefano Zampini } 1838135c375SStefano Zampini break; 1848135c375SStefano Zampini case 3: 1858135c375SStefano 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); 1868f07809aSStefano Zampini if (!hashocoord) { 1878f07809aSStefano 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); 1888f07809aSStefano Zampini } 1898135c375SStefano Zampini break; 1908135c375SStefano Zampini default: 1918135c375SStefano Zampini SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim); 1928135c375SStefano Zampini break; 1938135c375SStefano Zampini } 194*46900a7fSStefano Zampini ierr = DMSetApplicationContext(daview,dactx);CHKERRQ(ierr); 195*46900a7fSStefano Zampini ierr = DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 1968135c375SStefano Zampini ierr = DMSetUp(daview);CHKERRQ(ierr); 1978135c375SStefano Zampini if (!xcoor) { 1988135c375SStefano Zampini ierr = DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); 1998135c375SStefano Zampini ierr = DMGetCoordinates(daview,&xcoor);CHKERRQ(ierr); 2008135c375SStefano Zampini } 2018f07809aSStefano Zampini if (dacoord) { 2028f07809aSStefano Zampini ierr = DMSetUp(dacoord);CHKERRQ(ierr); 2038135c375SStefano Zampini ierr = DMCreateLocalVector(dacoord,&xcoorl);CHKERRQ(ierr); 2048135c375SStefano Zampini ierr = DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 2058135c375SStefano Zampini ierr = DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 206*46900a7fSStefano Zampini ierr = DMDestroy(&dacoord);CHKERRQ(ierr); 2078f07809aSStefano Zampini } else { 2088f07809aSStefano Zampini PetscInt ien,jen,ken,nc,nl,cdof,deg; 2098f07809aSStefano Zampini char fecmesh[32]; 2108f07809aSStefano Zampini 2118f07809aSStefano Zampini ierr = DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken);CHKERRQ(ierr); 2128f07809aSStefano Zampini nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1); 2138f07809aSStefano Zampini 2148f07809aSStefano Zampini ierr = VecGetLocalSize(xcoor,&nl);CHKERRQ(ierr); 2158f07809aSStefano Zampini if (nc && nl % nc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc); 2168f07809aSStefano Zampini ierr = VecDuplicate(xcoor,&xcoorl);CHKERRQ(ierr); 2178f07809aSStefano Zampini ierr = VecCopy(xcoor,xcoorl);CHKERRQ(ierr); 2188f07809aSStefano Zampini ierr = VecSetDM(xcoorl,NULL);CHKERRQ(ierr); 2198f07809aSStefano Zampini cdof = nl/(nc*dim); 2208f07809aSStefano Zampini deg = 1; 2218f07809aSStefano Zampini while (1) { 2228f07809aSStefano Zampini PetscInt degd = 1; 2238f07809aSStefano Zampini for (i=0;i<dim;i++) degd *= (deg+1); 2248f07809aSStefano Zampini if (degd > cdof) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof); 2258f07809aSStefano Zampini if (degd == cdof) break; 2268f07809aSStefano Zampini deg++; 2278f07809aSStefano Zampini } 2288f07809aSStefano Zampini ierr = PetscSNPrintf(fecmesh,sizeof(fecmesh),"L2_T1_%dD_P%d",dim,deg);CHKERRQ(ierr); 2298f07809aSStefano Zampini ierr = PetscObjectSetName((PetscObject)xcoorl,fecmesh);CHKERRQ(ierr); 2308f07809aSStefano Zampini } 2318135c375SStefano Zampini 232*46900a7fSStefano Zampini /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */ 233*46900a7fSStefano Zampini ierr = PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl);CHKERRQ(ierr); 2348135c375SStefano Zampini ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr); 2358135c375SStefano Zampini 23677eacf09SStefano Zampini /* daview is composed with the original DMDA */ 23777eacf09SStefano Zampini ierr = PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview);CHKERRQ(ierr); 23877eacf09SStefano Zampini ierr = PetscObjectDereference((PetscObject)daview);CHKERRQ(ierr); 239*46900a7fSStefano Zampini } 24077eacf09SStefano Zampini 241c9493c35SStefano Zampini /* customize the viewer if present */ 242c9493c35SStefano Zampini if (viewer) { 243*46900a7fSStefano Zampini DMDAFieldGLVisViewerCtx *ctx; 244*46900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 245*46900a7fSStefano Zampini char fec[64]; 2468f07809aSStefano Zampini Vec xlocal; 247*46900a7fSStefano Zampini const char **dafieldname; 248*46900a7fSStefano Zampini char **fec_type,**fieldname; 249*46900a7fSStefano Zampini PetscInt *nlocal,*bss,*dims; 250*46900a7fSStefano Zampini PetscInt dim,M,N,P,dof,s,i,nf; 251*46900a7fSStefano Zampini PetscBool bsset; 2528f07809aSStefano Zampini 253*46900a7fSStefano Zampini ierr = DMGetApplicationContext(daview,(void**)&dactx);CHKERRQ(ierr); 2548f07809aSStefano Zampini ierr = DMCreateLocalVector(daview,&xlocal);CHKERRQ(ierr); 2558135c375SStefano Zampini ierr = DMDAGetFieldNames(da,(const char * const **)&dafieldname);CHKERRQ(ierr); 2568135c375SStefano Zampini ierr = DMDAGetNumVerticesGhosted(daview,&M,&N,&P);CHKERRQ(ierr); 257*46900a7fSStefano Zampini ierr = DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 258*46900a7fSStefano Zampini ierr = PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%dD_P1",dim);CHKERRQ(ierr); 259bb77a09fSStefano Zampini ierr = PetscMalloc5(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname);CHKERRQ(ierr); 2608135c375SStefano Zampini for (i=0;i<dof;i++) bss[i] = 1; 2618135c375SStefano Zampini nf = dof; 2628135c375SStefano Zampini 263*46900a7fSStefano Zampini ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");CHKERRQ(ierr); 2648135c375SStefano Zampini ierr = PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset);CHKERRQ(ierr); 2658135c375SStefano Zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 2668135c375SStefano Zampini if (bsset) { 2678135c375SStefano Zampini PetscInt t; 2688135c375SStefano Zampini for (i=0,t=0;i<nf;i++) t += bss[i]; 2698135c375SStefano Zampini if (t != dof) SETERRQ2(PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof); 2708135c375SStefano Zampini } else nf = dof; 2718135c375SStefano Zampini 2728135c375SStefano Zampini for (i=0,s=0;i<nf;i++) { 2738135c375SStefano Zampini ierr = PetscStrallocpy(fec,&fec_type[i]);CHKERRQ(ierr); 2748135c375SStefano Zampini if (bss[i] == 1) { 2758135c375SStefano Zampini ierr = PetscStrallocpy(dafieldname[s],&fieldname[i]);CHKERRQ(ierr); 2768135c375SStefano Zampini } else { 2778135c375SStefano Zampini PetscInt b; 2788135c375SStefano Zampini size_t tlen = 9; /* "Vector-" + end */ 2798135c375SStefano Zampini for (b=0;b<bss[i];b++) { 2808135c375SStefano Zampini size_t len; 2818135c375SStefano Zampini ierr = PetscStrlen(dafieldname[s+b],&len);CHKERRQ(ierr); 2828135c375SStefano Zampini tlen += len + 1; /* field + "-" */ 2838135c375SStefano Zampini } 2848135c375SStefano Zampini ierr = PetscMalloc1(tlen,&fieldname[i]);CHKERRQ(ierr); 2858135c375SStefano Zampini ierr = PetscStrcpy(fieldname[i],"Vector-");CHKERRQ(ierr); 2868135c375SStefano Zampini for (b=0;b<bss[i]-1;b++) { 2878135c375SStefano Zampini ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr); 2888135c375SStefano Zampini ierr = PetscStrcat(fieldname[i],"-");CHKERRQ(ierr); 2898135c375SStefano Zampini } 2908135c375SStefano Zampini ierr = PetscStrcat(fieldname[i],dafieldname[s+b]);CHKERRQ(ierr); 2918135c375SStefano Zampini } 292bb77a09fSStefano Zampini dims[i] = dim; 2938135c375SStefano Zampini nlocal[i] = M*N*P*bss[i]; 2948135c375SStefano Zampini s += bss[i]; 2958135c375SStefano Zampini } 2968135c375SStefano Zampini 297*46900a7fSStefano Zampini /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */ 2988135c375SStefano Zampini ierr = PetscNew(&ctx);CHKERRQ(ierr); 2998135c375SStefano Zampini ctx->xlocal = xlocal; 3008135c375SStefano Zampini 301*46900a7fSStefano Zampini ierr = PetscViewerGLVisSetFields(viewer,nf,(const char**)fieldname,(const char**)fec_type,nlocal,bss,dims,DMDASampleGLVisFields_Private,ctx,DMDAFieldDestroyGLVisViewerCtx_Private);CHKERRQ(ierr); 3028135c375SStefano Zampini for (i=0;i<nf;i++) { 3038135c375SStefano Zampini ierr = PetscFree(fec_type[i]);CHKERRQ(ierr); 3048135c375SStefano Zampini ierr = PetscFree(fieldname[i]);CHKERRQ(ierr); 3058135c375SStefano Zampini } 306bb77a09fSStefano Zampini ierr = PetscFree5(fec_type,nlocal,bss,dims,fieldname);CHKERRQ(ierr); 307c9493c35SStefano Zampini } 3088135c375SStefano Zampini PetscFunctionReturn(0); 3098135c375SStefano Zampini } 3108135c375SStefano Zampini 3118135c375SStefano Zampini static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer) 3128135c375SStefano Zampini { 31377eacf09SStefano Zampini DM da,cda; 3148135c375SStefano Zampini Vec xcoorl; 3158135c375SStefano Zampini PetscMPIInt commsize; 3168135c375SStefano Zampini const PetscScalar *array; 3178135c375SStefano Zampini PetscContainer glvis_container; 31877eacf09SStefano Zampini PetscInt dim,sdim,i,vid[8],mid,cid,cdof; 3198135c375SStefano Zampini PetscInt sx,sy,sz,ie,je,ke,ien,jen,ken; 3208135c375SStefano Zampini PetscInt gsx,gsy,gsz,gm,gn,gp,kst,jst,ist; 3218135c375SStefano Zampini PetscBool enabled = PETSC_TRUE, isascii; 3228135c375SStefano Zampini PetscErrorCode ierr; 3238f07809aSStefano Zampini const char *fmt; 3248135c375SStefano Zampini 3258135c375SStefano Zampini PetscFunctionBegin; 3268135c375SStefano Zampini PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3278135c375SStefano Zampini PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3288135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 3298135c375SStefano Zampini if (!isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 3308135c375SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&commsize);CHKERRQ(ierr); 3318135c375SStefano Zampini if (commsize > 1) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 3328135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 3338135c375SStefano Zampini 3348135c375SStefano Zampini /* get container: determines if a process visualizes is portion of the data or not */ 3358135c375SStefano Zampini ierr = PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container);CHKERRQ(ierr); 3368135c375SStefano Zampini if (!glvis_container) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 3378135c375SStefano Zampini { 3388135c375SStefano Zampini PetscViewerGLVisInfo glvis_info; 3398135c375SStefano Zampini ierr = PetscContainerGetPointer(glvis_container,(void**)&glvis_info);CHKERRQ(ierr); 3408135c375SStefano Zampini enabled = glvis_info->enabled; 34177eacf09SStefano Zampini fmt = glvis_info->fmt; 3428135c375SStefano Zampini } 343c9493c35SStefano Zampini /* this can happen if we are calling DMView outside of VecView_GLVis */ 344*46900a7fSStefano Zampini ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);CHKERRQ(ierr); 345*46900a7fSStefano Zampini if (!da) {ierr = DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL);CHKERRQ(ierr);} 34677eacf09SStefano Zampini ierr = PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da);CHKERRQ(ierr); 3478135c375SStefano Zampini if (!da) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA"); 3488135c375SStefano Zampini ierr = DMGetCoordinateDim(da,&sdim);CHKERRQ(ierr); 3498135c375SStefano Zampini 3508135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.1\n");CHKERRQ(ierr); 3518135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\ndimension\n");CHKERRQ(ierr); 3528135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",dim);CHKERRQ(ierr); 3538135c375SStefano Zampini 3548135c375SStefano Zampini if (!enabled) { 3558135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 3568135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 3578135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 3588135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 3598135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 3608135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 3618135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 3628135c375SStefano Zampini PetscFunctionReturn(0); 3638135c375SStefano Zampini } 3648135c375SStefano Zampini 3658135c375SStefano Zampini ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr); 3668135c375SStefano Zampini i = ien; 3678135c375SStefano Zampini if (dim > 1) i *= jen; 3688135c375SStefano Zampini if (dim > 2) i *= ken; 3698135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nelements\n");CHKERRQ(ierr); 3708135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",i);CHKERRQ(ierr); 3718135c375SStefano Zampini switch (dim) { 3728135c375SStefano Zampini case 1: 3738135c375SStefano Zampini for (ie = 0; ie < ien; ie++) { 3748135c375SStefano Zampini vid[0] = ie; 3758135c375SStefano Zampini vid[1] = ie+1; 3768135c375SStefano Zampini mid = 1; /* material id */ 3778135c375SStefano Zampini cid = 1; /* segment */ 3788135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1]);CHKERRQ(ierr); 3798135c375SStefano Zampini } 3808135c375SStefano Zampini break; 3818135c375SStefano Zampini case 2: 3828135c375SStefano Zampini for (je = 0; je < jen; je++) { 3838135c375SStefano Zampini for (ie = 0; ie < ien; ie++) { 3848135c375SStefano Zampini vid[0] = je*(ien+1) + ie; 3858135c375SStefano Zampini vid[1] = je*(ien+1) + ie+1; 3868135c375SStefano Zampini vid[2] = (je+1)*(ien+1) + ie+1; 3878135c375SStefano Zampini vid[3] = (je+1)*(ien+1) + ie; 3888135c375SStefano Zampini mid = 1; /* material id */ 3898135c375SStefano Zampini cid = 3; /* quad */ 3908135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3]);CHKERRQ(ierr); 3918135c375SStefano Zampini } 3928135c375SStefano Zampini } 3938135c375SStefano Zampini break; 3948135c375SStefano Zampini case 3: 3958135c375SStefano Zampini for (ke = 0; ke < ken; ke++) { 3968135c375SStefano Zampini for (je = 0; je < jen; je++) { 3978135c375SStefano Zampini for (ie = 0; ie < ien; ie++) { 3988135c375SStefano Zampini vid[0] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie; 3998135c375SStefano Zampini vid[1] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie+1; 4008135c375SStefano Zampini vid[2] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1; 4018135c375SStefano Zampini vid[3] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie; 4028135c375SStefano Zampini vid[4] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie; 4038135c375SStefano Zampini vid[5] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie+1; 4048135c375SStefano Zampini vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1; 4058135c375SStefano Zampini vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie; 4068135c375SStefano Zampini mid = 1; /* material id */ 4078135c375SStefano Zampini cid = 5; /* hex */ 4088135c375SStefano 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); 4098135c375SStefano Zampini } 4108135c375SStefano Zampini } 4118135c375SStefano Zampini } 4128135c375SStefano Zampini break; 4138135c375SStefano Zampini default: 4148135c375SStefano Zampini SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim); 4158135c375SStefano Zampini break; 4168135c375SStefano Zampini } 4178135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nboundary\n");CHKERRQ(ierr); 4188135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",0);CHKERRQ(ierr); 4198135c375SStefano Zampini 42077eacf09SStefano Zampini /* vertex coordinates */ 421*46900a7fSStefano Zampini ierr = PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr); 42277eacf09SStefano Zampini if (!xcoorl) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords"); 4238135c375SStefano Zampini ierr = DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr); 42477eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\nvertices\n");CHKERRQ(ierr); 42577eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken);CHKERRQ(ierr); 426*46900a7fSStefano Zampini ierr = VecGetDM(xcoorl,&cda);CHKERRQ(ierr); 42777eacf09SStefano Zampini ierr = VecGetArrayRead(xcoorl,&array);CHKERRQ(ierr); 4288f07809aSStefano Zampini if (!cda) { /* HO viz */ 4298f07809aSStefano Zampini const char *fecname; 4308f07809aSStefano Zampini PetscInt nc,nl; 4318f07809aSStefano Zampini 4328f07809aSStefano Zampini ierr = PetscObjectGetName((PetscObject)xcoorl,&fecname);CHKERRQ(ierr); 43377eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"nodes\n");CHKERRQ(ierr); 43477eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n");CHKERRQ(ierr); 43577eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"FiniteElementCollection: %s\n",fecname);CHKERRQ(ierr); 43677eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim);CHKERRQ(ierr); 43777eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n");CHKERRQ(ierr); /*Ordering::byVDIM*/ 4388f07809aSStefano Zampini /* L2 coordinates */ 4398f07809aSStefano Zampini ierr = DMDAGetNumElementsGhosted(da,&ien,&jen,&ken);CHKERRQ(ierr); 4408f07809aSStefano Zampini ierr = VecGetLocalSize(xcoorl,&nl);CHKERRQ(ierr); 4418f07809aSStefano Zampini nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1); 4428f07809aSStefano Zampini cdof = nl/nc; 4438f07809aSStefano Zampini if (!ien) ien++; 4448f07809aSStefano Zampini if (!jen) jen++; 4458f07809aSStefano Zampini if (!ken) ken++; 4468f07809aSStefano Zampini ist = jst = kst = 0; 4478f07809aSStefano Zampini gm = ien; 4488f07809aSStefano Zampini gn = jen; 4498f07809aSStefano Zampini gp = ken; 45077eacf09SStefano Zampini } else { 451*46900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 452*46900a7fSStefano Zampini 453*46900a7fSStefano Zampini ierr = DMGetApplicationContext(da,(void**)&dactx);CHKERRQ(ierr); 45477eacf09SStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"%D\n",sdim);CHKERRQ(ierr); 4558f07809aSStefano Zampini cdof = sdim; 4568135c375SStefano Zampini ierr = DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL);CHKERRQ(ierr); 4578135c375SStefano Zampini ierr = DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp);CHKERRQ(ierr); 458*46900a7fSStefano Zampini if (dactx->ll) { 459*46900a7fSStefano Zampini kst = jst = ist = 0; 460*46900a7fSStefano Zampini } else { 4618135c375SStefano Zampini kst = gsz != sz ? 1 : 0; 4628135c375SStefano Zampini jst = gsy != sy ? 1 : 0; 4638135c375SStefano Zampini ist = gsx != sx ? 1 : 0; 4648f07809aSStefano Zampini } 465*46900a7fSStefano Zampini } 4668135c375SStefano Zampini for (ke = kst; ke < kst + ken; ke++) { 4678135c375SStefano Zampini for (je = jst; je < jst + jen; je++) { 4688135c375SStefano Zampini for (ie = ist; ie < ist + ien; ie++) { 46977eacf09SStefano Zampini PetscInt c; 4708135c375SStefano Zampini 4718135c375SStefano Zampini i = ke * gm * gn + je * gm + ie; 47277eacf09SStefano Zampini for (c=0;c<cdof/sdim;c++) { 47377eacf09SStefano Zampini PetscInt d; 47477eacf09SStefano Zampini for (d=0;d<sdim;d++) { 4758f07809aSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]));CHKERRQ(ierr); 4768135c375SStefano Zampini } 4778f07809aSStefano Zampini ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 47877eacf09SStefano Zampini } 4798135c375SStefano Zampini } 4808135c375SStefano Zampini } 4818135c375SStefano Zampini } 4828135c375SStefano Zampini ierr = VecRestoreArrayRead(xcoorl,&array);CHKERRQ(ierr); 4838135c375SStefano Zampini PetscFunctionReturn(0); 4848135c375SStefano Zampini } 4858135c375SStefano Zampini 4868135c375SStefano 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 */ 4878135c375SStefano Zampini PETSC_INTERN PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer) 4888135c375SStefano Zampini { 4898135c375SStefano Zampini PetscErrorCode ierr; 4908135c375SStefano Zampini PetscBool isglvis,isascii; 4918135c375SStefano Zampini 4928135c375SStefano Zampini PetscFunctionBegin; 4938135c375SStefano Zampini PetscValidHeaderSpecific(dm,DM_CLASSID,1); 4948135c375SStefano Zampini PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4958135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 4968135c375SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 4978135c375SStefano Zampini if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII"); 4988135c375SStefano Zampini if (isglvis) { 4998135c375SStefano Zampini PetscViewer view; 5008135c375SStefano Zampini PetscViewerGLVisType type; 5018135c375SStefano Zampini 5028135c375SStefano Zampini ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr); 5038135c375SStefano Zampini ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr); 5048135c375SStefano Zampini if (view) { /* in the socket case, it may happen that the connection failed */ 5058135c375SStefano Zampini if (type == PETSC_VIEWER_GLVIS_SOCKET) { 5068135c375SStefano Zampini PetscMPIInt size,rank; 5078135c375SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr); 5088135c375SStefano Zampini ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 5098135c375SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"parallel %D %D\nmesh\n",size,rank);CHKERRQ(ierr); 5108135c375SStefano Zampini } 5118135c375SStefano Zampini ierr = DMDAView_GLVis_ASCII(dm,view);CHKERRQ(ierr); 5128135c375SStefano Zampini ierr = PetscViewerFlush(view);CHKERRQ(ierr); 5138135c375SStefano Zampini if (type == PETSC_VIEWER_GLVIS_SOCKET) { 5148135c375SStefano Zampini PetscInt dim; 5158135c375SStefano Zampini const char* name; 5168135c375SStefano Zampini 5178135c375SStefano Zampini ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr); 5188135c375SStefano Zampini ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 5198135c375SStefano Zampini ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,dim,name);CHKERRQ(ierr); 5208135c375SStefano Zampini ierr = PetscBarrier((PetscObject)dm);CHKERRQ(ierr); 5218135c375SStefano Zampini } 5228135c375SStefano Zampini } 5238135c375SStefano Zampini } else { 5248135c375SStefano Zampini ierr = DMDAView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr); 5258135c375SStefano Zampini } 5268135c375SStefano Zampini PetscFunctionReturn(0); 5278135c375SStefano Zampini } 528