Actual source code: dupl.c

  1: /*$Id: dupl.c,v 1.14 2001/04/10 19:34:10 bsmith Exp $*/

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

  5: #undef __FUNCT__  
  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 SLES on one block of the preconditioner.

 27:    Concepts: PetscViewer^sequential version

 29: .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerRestoreSingleton()
 30: @*/
 31: int PetscViewerGetSingleton(PetscViewer viewer,PetscViewer *outviewer)
 32: {
 33:   int ierr,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(1,"Cannot get singleton PetscViewer for type %s",viewer->type_name);
 46:   }
 47:   return(0);
 48: }

 50: #undef __FUNCT__  
 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: int PetscViewerRestoreSingleton(PetscViewer viewer,PetscViewer *outviewer)
 68: {
 69:   int ierr,size;


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