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