xref: /petsc/src/vec/vec/utils/vecglvis.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1 #include <petsc/private/glvisviewerimpl.h>
2 #include <petsc/private/glvisvecimpl.h>
3 
4 static PetscErrorCode PetscViewerGLVisVecInfoDestroy_Private(void *ptr) {
5   PetscViewerGLVisVecInfo info = (PetscViewerGLVisVecInfo)ptr;
6 
7   PetscFunctionBeginUser;
8   PetscCall(PetscFree(info->fec_type));
9   PetscCall(PetscFree(info));
10   PetscFunctionReturn(0);
11 }
12 
13 /* the main function to visualize vectors using GLVis */
14 PetscErrorCode VecView_GLVis(Vec U, PetscViewer viewer) {
15   PetscErrorCode (*g2lfields)(PetscObject, PetscInt, PetscObject[], void *);
16   Vec                   *Ufield;
17   const char           **fec_type;
18   PetscViewerGLVisStatus sockstatus;
19   PetscViewerGLVisType   socktype;
20   void                  *userctx;
21   PetscInt               i, nfields, *spacedim;
22   PetscBool              pause = PETSC_FALSE;
23 
24   PetscFunctionBegin;
25   PetscCall(PetscViewerGLVisGetStatus_Private(viewer, &sockstatus));
26   if (sockstatus == PETSCVIEWERGLVIS_DISABLED) PetscFunctionReturn(0);
27   /* if the user did not customize the viewer through the API, we need extra data that can be attached to the Vec */
28   PetscCall(PetscViewerGLVisGetFields_Private(viewer, &nfields, NULL, NULL, NULL, NULL, NULL));
29   if (!nfields) {
30     PetscObject dm;
31 
32     PetscCall(PetscObjectQuery((PetscObject)U, "__PETSc_dm", &dm));
33     if (dm) {
34       PetscCall(PetscViewerGLVisSetDM_Private(viewer, dm));
35     } else SETERRQ(PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "You need to provide a DM or use PetscViewerGLVisSetFields()");
36   }
37   PetscCall(PetscViewerGLVisGetFields_Private(viewer, &nfields, &fec_type, &spacedim, &g2lfields, (PetscObject **)&Ufield, &userctx));
38   if (!nfields) PetscFunctionReturn(0);
39 
40   PetscCall(PetscViewerGLVisGetType_Private(viewer, &socktype));
41 
42   for (i = 0; i < nfields; i++) {
43     PetscObject    fdm;
44     PetscContainer container;
45 
46     /* attach visualization info to the vector */
47     PetscCall(PetscObjectQuery((PetscObject)Ufield[i], "_glvis_info_container", (PetscObject *)&container));
48     if (!container) {
49       PetscViewerGLVisVecInfo info;
50 
51       PetscCall(PetscNew(&info));
52       PetscCall(PetscStrallocpy(fec_type[i], &info->fec_type));
53       PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)U), &container));
54       PetscCall(PetscContainerSetPointer(container, (void *)info));
55       PetscCall(PetscContainerSetUserDestroy(container, PetscViewerGLVisVecInfoDestroy_Private));
56       PetscCall(PetscObjectCompose((PetscObject)Ufield[i], "_glvis_info_container", (PetscObject)container));
57       PetscCall(PetscContainerDestroy(&container));
58     }
59     /* attach the mesh to the viz vectors */
60     PetscCall(PetscObjectQuery((PetscObject)Ufield[i], "__PETSc_dm", &fdm));
61     if (!fdm) {
62       PetscObject dm;
63 
64       PetscCall(PetscViewerGLVisGetDM_Private(viewer, &dm));
65       if (!dm) PetscCall(PetscObjectQuery((PetscObject)U, "__PETSc_dm", &dm));
66       PetscCheck(dm, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Mesh not present");
67       PetscCall(PetscObjectCompose((PetscObject)Ufield[i], "__PETSc_dm", dm));
68     }
69   }
70 
71   /* user-provided sampling */
72   if (g2lfields) {
73     PetscCall((*g2lfields)((PetscObject)U, nfields, (PetscObject *)Ufield, userctx));
74   } else {
75     PetscCheck(nfields <= 1, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Don't know how to sample %" PetscInt_FMT " fields", nfields);
76     PetscCall(VecCopy(U, Ufield[0]));
77   }
78 
79   /* TODO callback to user routine to disable/enable subdomains */
80   for (i = 0; i < nfields; i++) {
81     PetscObject dm;
82     PetscViewer view;
83 
84     PetscCall(PetscObjectQuery((PetscObject)Ufield[i], "__PETSc_dm", &dm));
85     PetscCall(PetscViewerGLVisGetWindow_Private(viewer, i, &view));
86     if (!view) continue; /* socket window has been closed */
87     if (socktype == PETSC_VIEWER_GLVIS_SOCKET) {
88       PetscMPIInt size, rank;
89       const char *name;
90 
91       PetscCallMPI(MPI_Comm_size(PetscObjectComm(dm), &size));
92       PetscCallMPI(MPI_Comm_rank(PetscObjectComm(dm), &rank));
93       PetscCall(PetscObjectGetName((PetscObject)Ufield[i], &name));
94 
95       PetscCall(PetscGLVisCollectiveBegin(PetscObjectComm(dm), &view));
96       PetscCall(PetscViewerASCIIPrintf(view, "parallel %d %d\nsolution\n", size, rank));
97       PetscCall(PetscObjectView(dm, view));
98       PetscCall(VecView(Ufield[i], view));
99       PetscCall(PetscViewerGLVisInitWindow_Private(view, PETSC_FALSE, spacedim[i], name));
100       PetscCall(PetscGLVisCollectiveEnd(PetscObjectComm(dm), &view));
101       if (view) pause = PETSC_TRUE; /* at least one window is connected */
102     } else {
103       PetscCall(PetscObjectView(dm, view));
104       PetscCall(VecView(Ufield[i], view));
105     }
106     PetscCall(PetscViewerGLVisRestoreWindow_Private(viewer, i, &view));
107   }
108   if (pause) PetscCall(PetscViewerGLVisPause_Private(viewer));
109   PetscFunctionReturn(0);
110 }
111