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 PetscErrorCode ierr; 1448135c375SStefano Zampini 1458135c375SStefano Zampini PetscFunctionBegin; 1469566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(oda,"GLVisGraphicsDMDAGhosted",(PetscObject*)&daview)); 14746900a7fSStefano Zampini if (!daview) { 14846900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 14946900a7fSStefano Zampini DM dacoord = NULL; 15046900a7fSStefano Zampini Vec xcoor,xcoorl; 15146900a7fSStefano Zampini PetscBool hashocoord = PETSC_FALSE; 15246900a7fSStefano Zampini const PetscInt *lx,*ly,*lz; 15346900a7fSStefano Zampini PetscInt dim,M,N,P,m,n,p,dof,s,i; 15446900a7fSStefano Zampini 1559566063dSJacob Faibussowitsch PetscCall(PetscNew(&dactx)); 156a7bda6aaSStefano Zampini dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */ 1579566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Options","PetscViewer");PetscCall(ierr); 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_glvis_dm_da_ll","Left-looking subdomain view",NULL,dactx->ll,&dactx->ll,NULL)); 1599566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 1608135c375SStefano Zampini /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */ 1619566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,&M,&N,&P,&m,&n,&p,&dof,&s,NULL,NULL,NULL,NULL)); 1629566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da,&lx,&ly,&lz)); 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da,"_glvis_mesh_coords",(PetscObject*)&xcoor)); 16477eacf09SStefano Zampini if (!xcoor) { 1659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&xcoor)); 16677eacf09SStefano Zampini } else { 1678f07809aSStefano Zampini hashocoord = PETSC_TRUE; 16877eacf09SStefano Zampini } 1699566063dSJacob Faibussowitsch PetscCall(PetscInfo(da,"Creating auxilary DMDA for managing GLVis graphics\n")); 1708135c375SStefano Zampini switch (dim) { 1718135c375SStefano Zampini case 1: 1729566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,dof,1,lx,&daview)); 1738f07809aSStefano Zampini if (!hashocoord) { 1749566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,M,1,1,lx,&dacoord)); 1758f07809aSStefano Zampini } 1768135c375SStefano Zampini break; 1778135c375SStefano Zampini case 2: 1789566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,lx,ly,&daview)); 1798f07809aSStefano Zampini if (!hashocoord) { 1809566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da),DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,2,1,lx,ly,&dacoord)); 1818f07809aSStefano Zampini } 1828135c375SStefano Zampini break; 1838135c375SStefano Zampini case 3: 1849566063dSJacob 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)); 1858f07809aSStefano Zampini if (!hashocoord) { 1869566063dSJacob 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)); 1878f07809aSStefano Zampini } 1888135c375SStefano Zampini break; 1898135c375SStefano Zampini default: 19098921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim); 1918135c375SStefano Zampini } 1929566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(daview,dactx)); 1939566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContextDestroy(daview,DMDAGhostedDestroyGLVisViewerCtx_Private)); 1949566063dSJacob Faibussowitsch PetscCall(DMSetUp(daview)); 1958135c375SStefano Zampini if (!xcoor) { 1969566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daview,0.0,1.0,0.0,1.0,0.0,1.0)); 1979566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daview,&xcoor)); 1988135c375SStefano Zampini } 1998f07809aSStefano Zampini if (dacoord) { 2009566063dSJacob Faibussowitsch PetscCall(DMSetUp(dacoord)); 2019566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dacoord,&xcoorl)); 2029566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dacoord,xcoor,INSERT_VALUES,xcoorl)); 2039566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dacoord,xcoor,INSERT_VALUES,xcoorl)); 2049566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dacoord)); 2058f07809aSStefano Zampini } else { 2068f07809aSStefano Zampini PetscInt ien,jen,ken,nc,nl,cdof,deg; 2073e7cab22SLisandro Dalcin char fecmesh[64]; 20821414b21SStefano Zampini const char *name; 20921414b21SStefano Zampini PetscBool flg; 2108135c375SStefano Zampini 2119566063dSJacob Faibussowitsch PetscCall(DMDAGetNumElementsGhosted(daview,&ien,&jen,&ken)); 2128f07809aSStefano Zampini nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1); 2138f07809aSStefano Zampini 2149566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xcoor,&nl)); 2152c71b3e2SJacob Faibussowitsch PetscCheckFalse(nc && nl % nc,PETSC_COMM_SELF,PETSC_ERR_SUP,"Incompatible local coordinate size %D and number of cells %D",nl,nc); 2169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xcoor,&xcoorl)); 2179566063dSJacob Faibussowitsch PetscCall(VecCopy(xcoor,xcoorl)); 2189566063dSJacob Faibussowitsch PetscCall(VecSetDM(xcoorl,NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xcoor,&name)); 2209566063dSJacob Faibussowitsch PetscCall(PetscStrbeginswith(name,"FiniteElementCollection:",&flg)); 22121414b21SStefano Zampini if (!flg) { 22221414b21SStefano Zampini deg = 0; 22321414b21SStefano Zampini if (nc && nl) { 2248f07809aSStefano Zampini cdof = nl/(nc*dim); 2258f07809aSStefano Zampini deg = 1; 2268f07809aSStefano Zampini while (1) { 2278f07809aSStefano Zampini PetscInt degd = 1; 2288f07809aSStefano Zampini for (i=0;i<dim;i++) degd *= (deg+1); 229*08401ef6SPierre Jolivet PetscCheck(degd <= cdof,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell dofs %D",cdof); 2308f07809aSStefano Zampini if (degd == cdof) break; 2318f07809aSStefano Zampini deg++; 2328f07809aSStefano Zampini } 23321414b21SStefano Zampini } 2349566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(fecmesh,sizeof(fecmesh),"FiniteElementCollection: L2_T1_%DD_P%D",dim,deg)); 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xcoorl,fecmesh)); 23621414b21SStefano Zampini } else { 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xcoorl,name)); 23821414b21SStefano Zampini } 2398f07809aSStefano Zampini } 2408135c375SStefano Zampini 24146900a7fSStefano Zampini /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */ 2429566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)daview,"GLVisGraphicsCoordsGhosted",(PetscObject)xcoorl)); 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)xcoorl)); 2448135c375SStefano Zampini 24577eacf09SStefano Zampini /* daview is composed with the original DMDA */ 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(oda,"GLVisGraphicsDMDAGhosted",(PetscObject)daview)); 2479566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)daview)); 24846900a7fSStefano Zampini } 24977eacf09SStefano Zampini 250c9493c35SStefano Zampini /* customize the viewer if present */ 251c9493c35SStefano Zampini if (viewer) { 25246900a7fSStefano Zampini DMDAFieldGLVisViewerCtx *ctx; 25346900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 25446900a7fSStefano Zampini char fec[64]; 255ab166b77SStefano Zampini Vec xlocal,*Ufield; 25646900a7fSStefano Zampini const char **dafieldname; 25746900a7fSStefano Zampini char **fec_type,**fieldname; 25846900a7fSStefano Zampini PetscInt *nlocal,*bss,*dims; 25946900a7fSStefano Zampini PetscInt dim,M,N,P,dof,s,i,nf; 26046900a7fSStefano Zampini PetscBool bsset; 2618f07809aSStefano Zampini 2629566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daview,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 2639566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(daview,&dactx)); 2649566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(daview,&xlocal)); 2659566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldNames(da,(const char * const **)&dafieldname)); 2669566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVerticesGhosted(daview,&M,&N,&P)); 2679566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec,sizeof(fec),"FiniteElementCollection: H1_%DD_P1",dim)); 2689566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(dof,&fec_type,dof,&nlocal,dof,&bss,dof,&dims,dof,&fieldname,dof,&Ufield)); 2698135c375SStefano Zampini for (i=0;i<dof;i++) bss[i] = 1; 2708135c375SStefano Zampini nf = dof; 2718135c375SStefano Zampini 2729566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PetscObjectComm(oda),oda->prefix,"GLVis PetscViewer DMDA Field options","PetscViewer");PetscCall(ierr); 2739566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-viewer_glvis_dm_da_bs","Block sizes for subfields; enable vector representation",NULL,bss,&nf,&bsset)); 2749566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 2758135c375SStefano Zampini if (bsset) { 2768135c375SStefano Zampini PetscInt t; 2778135c375SStefano Zampini for (i=0,t=0;i<nf;i++) t += bss[i]; 278*08401ef6SPierre Jolivet PetscCheck(t == dof,PetscObjectComm(oda),PETSC_ERR_USER,"Sum of block sizes %D should equal %D",t,dof); 2798135c375SStefano Zampini } else nf = dof; 2808135c375SStefano Zampini 2818135c375SStefano Zampini for (i=0,s=0;i<nf;i++) { 2829566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec,&fec_type[i])); 2838135c375SStefano Zampini if (bss[i] == 1) { 2849566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(dafieldname[s],&fieldname[i])); 2858135c375SStefano Zampini } else { 2868135c375SStefano Zampini PetscInt b; 2878135c375SStefano Zampini size_t tlen = 9; /* "Vector-" + end */ 2888135c375SStefano Zampini for (b=0;b<bss[i];b++) { 2898135c375SStefano Zampini size_t len; 2909566063dSJacob Faibussowitsch PetscCall(PetscStrlen(dafieldname[s+b],&len)); 2918135c375SStefano Zampini tlen += len + 1; /* field + "-" */ 2928135c375SStefano Zampini } 2939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tlen,&fieldname[i])); 2949566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(fieldname[i],"Vector-")); 2958135c375SStefano Zampini for (b=0;b<bss[i]-1;b++) { 2969566063dSJacob Faibussowitsch PetscCall(PetscStrcat(fieldname[i],dafieldname[s+b])); 2979566063dSJacob Faibussowitsch PetscCall(PetscStrcat(fieldname[i],"-")); 2988135c375SStefano Zampini } 2999566063dSJacob Faibussowitsch PetscCall(PetscStrcat(fieldname[i],dafieldname[s+b])); 3008135c375SStefano Zampini } 301bb77a09fSStefano Zampini dims[i] = dim; 3028135c375SStefano Zampini nlocal[i] = M*N*P*bss[i]; 3038135c375SStefano Zampini s += bss[i]; 3048135c375SStefano Zampini } 3058135c375SStefano Zampini 30646900a7fSStefano Zampini /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */ 3079566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 3088135c375SStefano Zampini ctx->xlocal = xlocal; 3098135c375SStefano Zampini 3104cac2994SStefano Zampini /* create work vectors */ 3114cac2994SStefano Zampini for (i=0;i<nf;i++) { 3129566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da),nlocal[i],PETSC_DECIDE,&Ufield[i])); 3139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Ufield[i],fieldname[i])); 3149566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(Ufield[i],bss[i])); 3159566063dSJacob Faibussowitsch PetscCall(VecSetDM(Ufield[i],da)); 3164cac2994SStefano Zampini } 3174cac2994SStefano Zampini 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer,nf,(const char**)fec_type,dims,DMDASampleGLVisFields_Private,(PetscObject*)Ufield,ctx,DMDAFieldDestroyGLVisViewerCtx_Private)); 3198135c375SStefano Zampini for (i=0;i<nf;i++) { 3209566063dSJacob Faibussowitsch PetscCall(PetscFree(fec_type[i])); 3219566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldname[i])); 3229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ufield[i])); 3238135c375SStefano Zampini } 3249566063dSJacob Faibussowitsch PetscCall(PetscFree6(fec_type,nlocal,bss,dims,fieldname,Ufield)); 325c9493c35SStefano Zampini } 3268135c375SStefano Zampini PetscFunctionReturn(0); 3278135c375SStefano Zampini } 3288135c375SStefano Zampini 3298135c375SStefano Zampini static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer) 3308135c375SStefano Zampini { 33177eacf09SStefano Zampini DM da,cda; 3328135c375SStefano Zampini Vec xcoorl; 3333924b612SStefano Zampini PetscMPIInt size; 3348135c375SStefano Zampini const PetscScalar *array; 3358135c375SStefano Zampini PetscContainer glvis_container; 33677eacf09SStefano Zampini PetscInt dim,sdim,i,vid[8],mid,cid,cdof; 33721414b21SStefano Zampini PetscInt sx,sy,sz,ie,je,ke,ien,jen,ken,nel; 3388135c375SStefano Zampini PetscInt gsx,gsy,gsz,gm,gn,gp,kst,jst,ist; 3398135c375SStefano Zampini PetscBool enabled = PETSC_TRUE, isascii; 3408f07809aSStefano Zampini const char *fmt; 3418135c375SStefano Zampini 3428135c375SStefano Zampini PetscFunctionBegin; 3438135c375SStefano Zampini PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3448135c375SStefano Zampini PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii)); 3467a8be351SBarry Smith PetscCheck(isascii,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERASCII"); 3479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size)); 348*08401ef6SPierre Jolivet PetscCheck(size <= 1,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Use single sequential viewers for parallel visualization"); 3499566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 3508135c375SStefano Zampini 3518135c375SStefano Zampini /* get container: determines if a process visualizes is portion of the data or not */ 3529566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)viewer,"_glvis_info_container",(PetscObject*)&glvis_container)); 3537a8be351SBarry Smith PetscCheck(glvis_container,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis container"); 3548135c375SStefano Zampini { 3558135c375SStefano Zampini PetscViewerGLVisInfo glvis_info; 3569566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(glvis_container,(void**)&glvis_info)); 3578135c375SStefano Zampini enabled = glvis_info->enabled; 35877eacf09SStefano Zampini fmt = glvis_info->fmt; 3598135c375SStefano Zampini } 360c9493c35SStefano Zampini /* this can happen if we are calling DMView outside of VecView_GLVis */ 3619566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da)); 3629566063dSJacob Faibussowitsch if (!da) PetscCall(DMSetUpGLVisViewer_DMDA((PetscObject)dm,NULL)); 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm,"GLVisGraphicsDMDAGhosted",(PetscObject*)&da)); 3647a8be351SBarry Smith PetscCheck(da,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted DMDA"); 3659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(da,&sdim)); 3668135c375SStefano Zampini 3679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"MFEM mesh v1.0\n")); 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\ndimension\n")); 3699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",dim)); 3708135c375SStefano Zampini 3718135c375SStefano Zampini if (!enabled) { 3729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\nelements\n")); 3739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",0)); 3749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\nboundary\n")); 3759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",0)); 3769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n")); 3779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",0)); 3789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",sdim)); 3798135c375SStefano Zampini PetscFunctionReturn(0); 3808135c375SStefano Zampini } 3818135c375SStefano Zampini 3829566063dSJacob Faibussowitsch PetscCall(DMDAGetNumElementsGhosted(da,&ien,&jen,&ken)); 38321414b21SStefano Zampini nel = ien; 38421414b21SStefano Zampini if (dim > 1) nel *= jen; 38521414b21SStefano Zampini if (dim > 2) nel *= ken; 3869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\nelements\n")); 3879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",nel)); 3888135c375SStefano Zampini switch (dim) { 3898135c375SStefano Zampini case 1: 3908135c375SStefano Zampini for (ie = 0; ie < ien; ie++) { 3918135c375SStefano Zampini vid[0] = ie; 3928135c375SStefano Zampini vid[1] = ie+1; 3938135c375SStefano Zampini mid = 1; /* material id */ 3948135c375SStefano Zampini cid = 1; /* segment */ 3959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D %D %D %D\n",mid,cid,vid[0],vid[1])); 3968135c375SStefano Zampini } 3978135c375SStefano Zampini break; 3988135c375SStefano Zampini case 2: 3998135c375SStefano Zampini for (je = 0; je < jen; je++) { 4008135c375SStefano Zampini for (ie = 0; ie < ien; ie++) { 4018135c375SStefano Zampini vid[0] = je*(ien+1) + ie; 4028135c375SStefano Zampini vid[1] = je*(ien+1) + ie+1; 4038135c375SStefano Zampini vid[2] = (je+1)*(ien+1) + ie+1; 4048135c375SStefano Zampini vid[3] = (je+1)*(ien+1) + ie; 4058135c375SStefano Zampini mid = 1; /* material id */ 4068135c375SStefano Zampini cid = 3; /* quad */ 4079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D %D %D %D %D %D\n",mid,cid,vid[0],vid[1],vid[2],vid[3])); 4088135c375SStefano Zampini } 4098135c375SStefano Zampini } 4108135c375SStefano Zampini break; 4118135c375SStefano Zampini case 3: 4128135c375SStefano Zampini for (ke = 0; ke < ken; ke++) { 4138135c375SStefano Zampini for (je = 0; je < jen; je++) { 4148135c375SStefano Zampini for (ie = 0; ie < ien; ie++) { 4158135c375SStefano Zampini vid[0] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie; 4168135c375SStefano Zampini vid[1] = ke*(jen+1)*(ien+1) + je*(ien+1) + ie+1; 4178135c375SStefano Zampini vid[2] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1; 4188135c375SStefano Zampini vid[3] = ke*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie; 4198135c375SStefano Zampini vid[4] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie; 4208135c375SStefano Zampini vid[5] = (ke+1)*(jen+1)*(ien+1) + je*(ien+1) + ie+1; 4218135c375SStefano Zampini vid[6] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie+1; 4228135c375SStefano Zampini vid[7] = (ke+1)*(jen+1)*(ien+1) + (je+1)*(ien+1) + ie; 4238135c375SStefano Zampini mid = 1; /* material id */ 4248135c375SStefano Zampini cid = 5; /* hex */ 4259566063dSJacob Faibussowitsch PetscCall(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])); 4268135c375SStefano Zampini } 4278135c375SStefano Zampini } 4288135c375SStefano Zampini } 4298135c375SStefano Zampini break; 4308135c375SStefano Zampini default: 43198921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Unsupported dimension %D",dim); 4328135c375SStefano Zampini } 4339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\nboundary\n")); 4349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",0)); 4358135c375SStefano Zampini 43677eacf09SStefano Zampini /* vertex coordinates */ 4379566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da,"GLVisGraphicsCoordsGhosted",(PetscObject*)&xcoorl)); 4387a8be351SBarry Smith PetscCheck(xcoorl,PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"Missing GLVis ghosted coords"); 4399566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVerticesGhosted(da,&ien,&jen,&ken)); 4409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\nvertices\n")); 4419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",ien*jen*ken)); 44221414b21SStefano Zampini if (nel) { 4439566063dSJacob Faibussowitsch PetscCall(VecGetDM(xcoorl,&cda)); 4449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoorl,&array)); 4458f07809aSStefano Zampini if (!cda) { /* HO viz */ 4468f07809aSStefano Zampini const char *fecname; 4478f07809aSStefano Zampini PetscInt nc,nl; 4488f07809aSStefano Zampini 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xcoorl,&fecname)); 4509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"nodes\n")); 4519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"FiniteElementSpace\n")); 4529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",fecname)); 4539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"VDim: %D\n",sdim)); 4549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Ordering: 1\n\n")); /*Ordering::byVDIM*/ 4558f07809aSStefano Zampini /* L2 coordinates */ 4569566063dSJacob Faibussowitsch PetscCall(DMDAGetNumElementsGhosted(da,&ien,&jen,&ken)); 4579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xcoorl,&nl)); 4588f07809aSStefano Zampini nc = ien*(jen>0 ? jen : 1)*(ken>0 ? ken : 1); 45921414b21SStefano Zampini cdof = nc ? nl/nc : 0; 4608f07809aSStefano Zampini if (!ien) ien++; 4618f07809aSStefano Zampini if (!jen) jen++; 4628f07809aSStefano Zampini if (!ken) ken++; 4638f07809aSStefano Zampini ist = jst = kst = 0; 4648f07809aSStefano Zampini gm = ien; 4658f07809aSStefano Zampini gn = jen; 4668f07809aSStefano Zampini gp = ken; 46777eacf09SStefano Zampini } else { 46846900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx; 46946900a7fSStefano Zampini 4709566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da,&dactx)); 4719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%D\n",sdim)); 4728f07809aSStefano Zampini cdof = sdim; 4739566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&sx,&sy,&sz,NULL,NULL,NULL)); 4749566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gsx,&gsy,&gsz,&gm,&gn,&gp)); 47546900a7fSStefano Zampini if (dactx->ll) { 47646900a7fSStefano Zampini kst = jst = ist = 0; 47746900a7fSStefano Zampini } else { 4788135c375SStefano Zampini kst = gsz != sz ? 1 : 0; 4798135c375SStefano Zampini jst = gsy != sy ? 1 : 0; 4808135c375SStefano Zampini ist = gsx != sx ? 1 : 0; 4818f07809aSStefano Zampini } 48246900a7fSStefano Zampini } 4838135c375SStefano Zampini for (ke = kst; ke < kst + ken; ke++) { 4848135c375SStefano Zampini for (je = jst; je < jst + jen; je++) { 4858135c375SStefano Zampini for (ie = ist; ie < ist + ien; ie++) { 48677eacf09SStefano Zampini PetscInt c; 4878135c375SStefano Zampini 4888135c375SStefano Zampini i = ke * gm * gn + je * gm + ie; 48977eacf09SStefano Zampini for (c=0;c<cdof/sdim;c++) { 49077eacf09SStefano Zampini PetscInt d; 49177eacf09SStefano Zampini for (d=0;d<sdim;d++) { 4929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,fmt,PetscRealPart(array[cdof*i+c*sdim+d]))); 4938135c375SStefano Zampini } 4949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 49577eacf09SStefano Zampini } 4968135c375SStefano Zampini } 4978135c375SStefano Zampini } 4988135c375SStefano Zampini } 4999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoorl,&array)); 50021414b21SStefano Zampini } 5018135c375SStefano Zampini PetscFunctionReturn(0); 5028135c375SStefano Zampini } 5038135c375SStefano Zampini 5040286d493SLisandro Dalcin PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer) 5058135c375SStefano Zampini { 5068135c375SStefano Zampini PetscFunctionBegin; 5079566063dSJacob Faibussowitsch PetscCall(DMView_GLVis(dm,viewer,DMDAView_GLVis_ASCII)); 5088135c375SStefano Zampini PetscFunctionReturn(0); 5098135c375SStefano Zampini } 510