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