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: /*@
  6:    PetscViewerGetSingleton - Creates a new PetscViewer (same type as the old)
  7:     that lives on a single processor (with MPI_comm PETSC_COMM_SELF)

  9:     Collective on PetscViewer

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

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

 17:    Level: advanced

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

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

 25:    Concepts: PetscViewer^sequential version

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

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

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

 51:     Collective on PetscViewer

 53:    Input Parameters:
 54: +  viewer - the PetscViewer to be duplicated
 55: -  outviewer - new PetscViewer

 57:    Level: advanced

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

 61: .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerGetSingleton()
 62: @*/
 63: int PetscViewerRestoreSingleton(PetscViewer viewer,PetscViewer *outviewer)
 64: {
 65:   int ierr,size;


 70:   MPI_Comm_size(viewer->comm,&size);
 71:   if (size == 1) {
 72:     PetscObjectDereference((PetscObject)viewer);
 73:     if (outviewer) *outviewer = 0;
 74:   } else if (viewer->ops->restoresingleton) {
 75:     (*viewer->ops->restoresingleton)(viewer,outviewer);
 76:   }
 77:   return(0);
 78: }