xref: /petsc/src/dm/interface/dmglvis.c (revision 67ea8cea998d0800565aa55bb307baf58d681386)
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 
DMView_GLVis(DM dm,PetscViewer viewer,PetscErrorCode (* DMView_GLVis_ASCII)(DM,PetscViewer))6d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_GLVis(DM dm, PetscViewer viewer, PetscErrorCode (*DMView_GLVis_ASCII)(DM, PetscViewer))
7d71ae5a4SJacob Faibussowitsch {
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);
139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
149566063dSJacob 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*34e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetType_Internal(viewer, &type));
21*34e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisGetDMWindow_Internal(viewer, &view));
223ba16761SJacob Faibussowitsch     if (!view) PetscFunctionReturn(PETSC_SUCCESS); /* 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 
289566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
299566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
309566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDim(dm, &sdim));
319566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)dm, &name));
320286d493SLisandro Dalcin 
339566063dSJacob Faibussowitsch       PetscCall(PetscGLVisCollectiveBegin(PetscObjectComm((PetscObject)dm), &view));
349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(view, "parallel %d %d\nmesh\n", size, rank));
359566063dSJacob Faibussowitsch       PetscCall(DMView_GLVis_ASCII(dm, view));
36*34e79e72SJacob Faibussowitsch       PetscCall(PetscViewerGLVisInitWindow_Internal(view, PETSC_TRUE, sdim, name));
379566063dSJacob Faibussowitsch       PetscCall(PetscGLVisCollectiveEnd(PetscObjectComm((PetscObject)dm), &view));
380286d493SLisandro Dalcin     } else {
399566063dSJacob Faibussowitsch       PetscCall(DMView_GLVis_ASCII(dm, view));
400286d493SLisandro Dalcin     }
41*34e79e72SJacob Faibussowitsch     PetscCall(PetscViewerGLVisRestoreDMWindow_Internal(viewer, &view));
420286d493SLisandro Dalcin   } else {
439566063dSJacob Faibussowitsch     PetscCall(DMView_GLVis_ASCII(dm, viewer));
440286d493SLisandro Dalcin   }
453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
460286d493SLisandro Dalcin }
47