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:     PetscCheck(dm, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "You need to provide a DM or use PetscViewerGLVisSetFields()");
 36:     PetscCall(PetscViewerGLVisSetDM_Internal(viewer, dm));
 37:   }
 38:   PetscCall(PetscViewerGLVisGetFields_Internal(viewer, &nfields, &fec_type, &spacedim, &g2lfields, (PetscObject **)&Ufield, &userctx));
 39:   if (!nfields) PetscFunctionReturn(PETSC_SUCCESS);

 41:   PetscCall(PetscViewerGLVisGetType_Internal(viewer, &socktype));

 43:   for (i = 0; i < nfields; i++) {
 44:     PetscObject    fdm;
 45:     PetscContainer container;

 47:     /* attach visualization info to the vector */
 48:     PetscCall(PetscObjectQuery((PetscObject)Ufield[i], "_glvis_info_container", (PetscObject *)&container));
 49:     if (!container) {
 50:       PetscViewerGLVisVecInfo info;

 52:       PetscCall(PetscNew(&info));
 53:       PetscCall(PetscStrallocpy(fec_type[i], &info->fec_type));
 54:       PetscCall(PetscObjectContainerCompose((PetscObject)Ufield[i], "_glvis_info_container", info, PetscViewerGLVisVecInfoDestroy_Private));
 55:     }
 56:     /* attach the mesh to the viz vectors */
 57:     PetscCall(PetscObjectQuery((PetscObject)Ufield[i], "__PETSc_dm", &fdm));
 58:     if (!fdm) {
 59:       PetscObject dm;

 61:       PetscCall(PetscViewerGLVisGetDM_Internal(viewer, &dm));
 62:       if (!dm) PetscCall(PetscObjectQuery((PetscObject)U, "__PETSc_dm", &dm));
 63:       PetscCheck(dm, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Mesh not present");
 64:       PetscCall(PetscObjectCompose((PetscObject)Ufield[i], "__PETSc_dm", dm));
 65:     }
 66:   }

 68:   /* user-provided sampling */
 69:   if (g2lfields) {
 70:     PetscCall((*g2lfields)((PetscObject)U, nfields, (PetscObject *)Ufield, userctx));
 71:   } else {
 72:     PetscCheck(nfields <= 1, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Don't know how to sample %" PetscInt_FMT " fields", nfields);
 73:     PetscCall(VecCopy(U, Ufield[0]));
 74:   }

 76:   /* TODO callback to user routine to disable/enable subdomains */
 77:   for (i = 0; i < nfields; i++) {
 78:     PetscObject dm;
 79:     PetscViewer view;

 81:     PetscCall(PetscObjectQuery((PetscObject)Ufield[i], "__PETSc_dm", &dm));
 82:     PetscCall(PetscViewerGLVisGetWindow_Internal(viewer, i, &view));
 83:     if (!view) continue; /* socket window has been closed */
 84:     if (socktype == PETSC_VIEWER_GLVIS_SOCKET) {
 85:       PetscMPIInt size, rank;
 86:       const char *name;

 88:       PetscCallMPI(MPI_Comm_size(PetscObjectComm(dm), &size));
 89:       PetscCallMPI(MPI_Comm_rank(PetscObjectComm(dm), &rank));
 90:       PetscCall(PetscObjectGetName((PetscObject)Ufield[i], &name));

 92:       PetscCall(PetscGLVisCollectiveBegin(PetscObjectComm(dm), &view));
 93:       PetscCall(PetscViewerASCIIPrintf(view, "parallel %d %d\nsolution\n", size, rank));
 94:       PetscCall(PetscObjectView(dm, view));
 95:       PetscCall(VecView(Ufield[i], view));
 96:       PetscCall(PetscViewerGLVisInitWindow_Internal(view, PETSC_FALSE, spacedim[i], name));
 97:       PetscCall(PetscGLVisCollectiveEnd(PetscObjectComm(dm), &view));
 98:       if (view) pause = PETSC_TRUE; /* at least one window is connected */
 99:     } else {
100:       PetscCall(PetscObjectView(dm, view));
101:       PetscCall(VecView(Ufield[i], view));
102:     }
103:     PetscCall(PetscViewerGLVisRestoreWindow_Internal(viewer, i, &view));
104:   }
105:   if (pause) PetscCall(PetscViewerGLVisPause_Internal(viewer));
106:   PetscFunctionReturn(PETSC_SUCCESS);
107: }