Actual source code: dupl.c

  1: #define PETSC_DLL

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

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

 11:     Collective on PetscViewer

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

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

 19:    Level: advanced

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

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

 27:    Concepts: PetscViewer^sequential version

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

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

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

 56:     Collective on PetscViewer

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

 62:    Level: advanced

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

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


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