xref: /petsc/src/dm/interface/dmglvis.c (revision 9566063d113dddea24716c546802770db7481bc0)
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   PetscBool      isglvis,isascii;
90286d493SLisandro Dalcin 
100286d493SLisandro Dalcin   PetscFunctionBegin;
110286d493SLisandro Dalcin   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
13*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
14*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
157a8be351SBarry Smith   PetscCheck(isglvis || isascii,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer must be of type VIEWERGLVIS or VIEWERASCII");
160286d493SLisandro Dalcin   if (isglvis) {
170286d493SLisandro Dalcin     PetscViewerGLVisType type;
180286d493SLisandro Dalcin     PetscViewer          view;
190286d493SLisandro Dalcin 
20*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetType_Private(viewer,&type));
21*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetDMWindow_Private(viewer,&view));
220286d493SLisandro Dalcin     if (!view) PetscFunctionReturn(0); /* socket window has been closed */
230286d493SLisandro Dalcin     if (type == PETSC_VIEWER_GLVIS_SOCKET) {
240286d493SLisandro Dalcin       PetscMPIInt size,rank;
250286d493SLisandro Dalcin       PetscInt    sdim;
260286d493SLisandro Dalcin       const char* name;
270286d493SLisandro Dalcin 
28*9566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm),&size));
29*9566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
30*9566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDim(dm,&sdim));
31*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)dm,&name));
320286d493SLisandro Dalcin 
33*9566063dSJacob Faibussowitsch       PetscCall(PetscGLVisCollectiveBegin(PetscObjectComm((PetscObject)dm),&view));
34*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(view,"parallel %d %d\nmesh\n",size,rank));
35*9566063dSJacob Faibussowitsch       PetscCall(DMView_GLVis_ASCII(dm,view));
36*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerGLVisInitWindow_Private(view,PETSC_TRUE,sdim,name));
37*9566063dSJacob Faibussowitsch       PetscCall(PetscGLVisCollectiveEnd(PetscObjectComm((PetscObject)dm),&view));
380286d493SLisandro Dalcin     } else {
39*9566063dSJacob Faibussowitsch       PetscCall(DMView_GLVis_ASCII(dm,view));
400286d493SLisandro Dalcin     }
41*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisRestoreDMWindow_Private(viewer,&view));
420286d493SLisandro Dalcin   } else {
43*9566063dSJacob Faibussowitsch     PetscCall(DMView_GLVis_ASCII(dm,viewer));
440286d493SLisandro Dalcin   }
450286d493SLisandro Dalcin   PetscFunctionReturn(0);
460286d493SLisandro Dalcin }
47