10286d493SLisandro Dalcin /* Routines to visualize DMs through GLVis */ 20286d493SLisandro Dalcin 30286d493SLisandro Dalcin #include <petsc/private/dmimpl.h> 40286d493SLisandro Dalcin #include <petsc/private/glvisviewerimpl.h> 50286d493SLisandro Dalcin 60286d493SLisandro Dalcin PetscErrorCode DMView_GLVis(DM dm, PetscViewer viewer, PetscErrorCode (*DMView_GLVis_ASCII)(DM,PetscViewer)) 70286d493SLisandro Dalcin { 80286d493SLisandro Dalcin PetscErrorCode ierr; 90286d493SLisandro Dalcin PetscBool isglvis,isascii; 100286d493SLisandro Dalcin 110286d493SLisandro Dalcin PetscFunctionBegin; 120286d493SLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13*064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 140286d493SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr); 150286d493SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 160286d493SLisandro Dalcin if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII"); 170286d493SLisandro Dalcin if (isglvis) { 180286d493SLisandro Dalcin PetscViewerGLVisType type; 190286d493SLisandro Dalcin PetscViewer view; 200286d493SLisandro Dalcin 210286d493SLisandro Dalcin ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr); 220286d493SLisandro Dalcin ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr); 230286d493SLisandro Dalcin if (!view) PetscFunctionReturn(0); /* socket window has been closed */ 240286d493SLisandro Dalcin if (type == PETSC_VIEWER_GLVIS_SOCKET) { 250286d493SLisandro Dalcin PetscMPIInt size,rank; 260286d493SLisandro Dalcin PetscInt sdim; 270286d493SLisandro Dalcin const char* name; 280286d493SLisandro Dalcin 29ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRMPI(ierr); 30ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr); 310286d493SLisandro Dalcin ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr); 320286d493SLisandro Dalcin ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr); 330286d493SLisandro Dalcin 340286d493SLisandro Dalcin ierr = PetscGLVisCollectiveBegin(PetscObjectComm((PetscObject)dm),&view);CHKERRQ(ierr); 350286d493SLisandro Dalcin ierr = PetscViewerASCIIPrintf(view,"parallel %d %d\nmesh\n",size,rank);CHKERRQ(ierr); 360286d493SLisandro Dalcin ierr = DMView_GLVis_ASCII(dm,view);CHKERRQ(ierr); 370286d493SLisandro Dalcin ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,sdim,name);CHKERRQ(ierr); 380286d493SLisandro Dalcin ierr = PetscGLVisCollectiveEnd(PetscObjectComm((PetscObject)dm),&view);CHKERRQ(ierr); 390286d493SLisandro Dalcin } else { 400286d493SLisandro Dalcin ierr = DMView_GLVis_ASCII(dm,view);CHKERRQ(ierr); 410286d493SLisandro Dalcin } 420286d493SLisandro Dalcin ierr = PetscViewerGLVisRestoreDMWindow_Private(viewer,&view);CHKERRQ(ierr); 430286d493SLisandro Dalcin } else { 440286d493SLisandro Dalcin ierr = DMView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr); 450286d493SLisandro Dalcin } 460286d493SLisandro Dalcin PetscFunctionReturn(0); 470286d493SLisandro Dalcin } 48