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