xref: /petsc/src/dm/interface/dmglvis.c (revision 0286d493dd5c8e4cb45d4ea43801100c0d0b6571)
1*0286d493SLisandro Dalcin /* Routines to visualize DMs through GLVis */
2*0286d493SLisandro Dalcin 
3*0286d493SLisandro Dalcin #include <petsc/private/dmimpl.h>
4*0286d493SLisandro Dalcin #include <petsc/private/glvisviewerimpl.h>
5*0286d493SLisandro Dalcin 
6*0286d493SLisandro Dalcin PetscErrorCode DMView_GLVis(DM dm, PetscViewer viewer, PetscErrorCode (*DMView_GLVis_ASCII)(DM,PetscViewer))
7*0286d493SLisandro Dalcin {
8*0286d493SLisandro Dalcin   PetscErrorCode ierr;
9*0286d493SLisandro Dalcin   PetscBool      isglvis,isascii;
10*0286d493SLisandro Dalcin 
11*0286d493SLisandro Dalcin   PetscFunctionBegin;
12*0286d493SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13*0286d493SLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
14*0286d493SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis);CHKERRQ(ierr);
15*0286d493SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
16*0286d493SLisandro Dalcin   if (!isglvis && !isascii) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
17*0286d493SLisandro Dalcin   if (isglvis) {
18*0286d493SLisandro Dalcin     PetscViewerGLVisType type;
19*0286d493SLisandro Dalcin     PetscViewer          view;
20*0286d493SLisandro Dalcin 
21*0286d493SLisandro Dalcin     ierr = PetscViewerGLVisGetType_Private(viewer,&type);CHKERRQ(ierr);
22*0286d493SLisandro Dalcin     ierr = PetscViewerGLVisGetDMWindow_Private(viewer,&view);CHKERRQ(ierr);
23*0286d493SLisandro Dalcin     if (!view) PetscFunctionReturn(0); /* socket window has been closed */
24*0286d493SLisandro Dalcin     if (type == PETSC_VIEWER_GLVIS_SOCKET) {
25*0286d493SLisandro Dalcin       PetscMPIInt size,rank;
26*0286d493SLisandro Dalcin       PetscInt    sdim;
27*0286d493SLisandro Dalcin       const char* name;
28*0286d493SLisandro Dalcin 
29*0286d493SLisandro Dalcin       ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size);CHKERRQ(ierr);
30*0286d493SLisandro Dalcin       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
31*0286d493SLisandro Dalcin       ierr = DMGetCoordinateDim(dm,&sdim);CHKERRQ(ierr);
32*0286d493SLisandro Dalcin       ierr = PetscObjectGetName((PetscObject)dm,&name);CHKERRQ(ierr);
33*0286d493SLisandro Dalcin 
34*0286d493SLisandro Dalcin       ierr = PetscGLVisCollectiveBegin(PetscObjectComm((PetscObject)dm),&view);CHKERRQ(ierr);
35*0286d493SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(view,"parallel %d %d\nmesh\n",size,rank);CHKERRQ(ierr);
36*0286d493SLisandro Dalcin       ierr = DMView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
37*0286d493SLisandro Dalcin       ierr = PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,sdim,name);CHKERRQ(ierr);
38*0286d493SLisandro Dalcin       ierr = PetscGLVisCollectiveEnd(PetscObjectComm((PetscObject)dm),&view);CHKERRQ(ierr);
39*0286d493SLisandro Dalcin     } else {
40*0286d493SLisandro Dalcin       ierr = DMView_GLVis_ASCII(dm,view);CHKERRQ(ierr);
41*0286d493SLisandro Dalcin     }
42*0286d493SLisandro Dalcin     ierr = PetscViewerGLVisRestoreDMWindow_Private(viewer,&view);CHKERRQ(ierr);
43*0286d493SLisandro Dalcin   } else {
44*0286d493SLisandro Dalcin     ierr = DMView_GLVis_ASCII(dm,viewer);CHKERRQ(ierr);
45*0286d493SLisandro Dalcin   }
46*0286d493SLisandro Dalcin   PetscFunctionReturn(0);
47*0286d493SLisandro Dalcin }
48