Actual source code: dupl.c

 2:  #include src/sys/src/viewer/viewerimpl.h

  6: /*@
  7:    PetscViewerGetSingleton - Creates a new PetscViewer (same type as the old)
  8:     that lives on a single processor (with MPI_comm PETSC_COMM_SELF)

 10:     Collective on PetscViewer

 12:    Input Parameter:
 13: .  viewer - the PetscViewer to be duplicated

 15:    Output Parameter:
 16: .  outviewer - new PetscViewer

 18:    Level: advanced

 20:    Notes: Call PetscViewerRestoreSingleton() to return this PetscViewer, NOT PetscViewerDestroy()

 22:      This is most commonly used to view a sequential object that is part of a 
 23:     parallel object. For example block Jacobi PC view could use this to obtain a
 24:     PetscViewer that is used with the sequential KSP on one block of the preconditioner.

 26:    Concepts: PetscViewer^sequential version

 28: .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerRestoreSingleton()
 29: @*/
 30: PetscErrorCode PetscViewerGetSingleton(PetscViewer viewer,PetscViewer *outviewer)
 31: {
 33:   PetscMPIInt    size;

 38:   MPI_Comm_size(viewer->comm,&size);
 39:   if (size == 1) {
 40:     *outviewer = viewer;
 41:     PetscObjectReference((PetscObject)viewer);
 42:   } else if (viewer->ops->getsingleton) {
 43:     (*viewer->ops->getsingleton)(viewer,outviewer);
 44:   } else {
 45:     SETERRQ1(PETSC_ERR_SUP,"Cannot get singleton PetscViewer for type %s",viewer->type_name);
 46:   }
 47:   return(0);
 48: }

 52: /*@
 53:    PetscViewerRestoreSingleton - Restores a new PetscViewer obtained with PetscViewerGetSingleton().

 55:     Collective on PetscViewer

 57:    Input Parameters:
 58: +  viewer - the PetscViewer to be duplicated
 59: -  outviewer - new PetscViewer

 61:    Level: advanced

 63:    Notes: Call PetscViewerGetSingleton() to get this PetscViewer, NOT PetscViewerCreate()

 65: .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerGetSingleton()
 66: @*/
 67: PetscErrorCode PetscViewerRestoreSingleton(PetscViewer viewer,PetscViewer *outviewer)
 68: {
 70:   PetscMPIInt    size;


 75:   MPI_Comm_size(viewer->comm,&size);
 76:   if (size == 1) {
 77:     PetscObjectDereference((PetscObject)viewer);
 78:     if (outviewer) *outviewer = 0;
 79:   } else if (viewer->ops->restoresingleton) {
 80:     (*viewer->ops->restoresingleton)(viewer,outviewer);
 81:   }
 82:   return(0);
 83: }