Actual source code: dupl.c
petsc-3.8.4 2018-03-24
2: #include <petsc/private/viewerimpl.h>
4: /*@
5: PetscViewerGetSubViewer - Creates a new PetscViewer (same type as the old)
6: that lives on a subcommunicator
8: Collective on PetscViewer
10: Input Parameter:
11: . viewer - the PetscViewer to be reproduced
13: Output Parameter:
14: . outviewer - new PetscViewer
16: Level: advanced
18: Notes: Call PetscViewerRestoreSubViewer() to return this PetscViewer, NOT PetscViewerDestroy()
20: This is most commonly used to view a sequential object that is part of a
21: parallel object. For example block Jacobi PC view could use this to obtain a
22: PetscViewer that is used with the sequential KSP on one block of the preconditioner.
24: Concepts: PetscViewer^sequential version
26: .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerRestoreSubViewer()
27: @*/
28: PetscErrorCode PetscViewerGetSubViewer(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer)
29: {
31: PetscMPIInt size;
36: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
37: if (size == 1) {
38: PetscObjectReference((PetscObject)viewer);
39: *outviewer = viewer;
40: } else if (viewer->ops->getsubviewer) {
41: (*viewer->ops->getsubviewer)(viewer,comm,outviewer);
42: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get SubViewer PetscViewer for type %s",((PetscObject)viewer)->type_name);
43: PetscViewerASCIIPushSynchronized(viewer);
44: return(0);
45: }
47: /*@
48: PetscViewerRestoreSubViewer - Restores a new PetscViewer obtained with PetscViewerGetSubViewer().
50: Collective on PetscViewer
52: Input Parameters:
53: + viewer - the PetscViewer that was reproduced
54: - outviewer - new PetscViewer
56: Level: advanced
58: Notes: Call PetscViewerGetSubViewer() to get this PetscViewer, NOT PetscViewerCreate()
60: .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerGetSubViewer()
61: @*/
62: PetscErrorCode PetscViewerRestoreSubViewer(PetscViewer viewer,MPI_Comm comm,PetscViewer *outviewer)
63: {
65: PetscMPIInt size;
70: MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);
71: if (size == 1) {
72: PetscObjectDereference((PetscObject)viewer);
73: if (outviewer) *outviewer = NULL;
74: } else if (viewer->ops->restoresubviewer) {
75: (*viewer->ops->restoresubviewer)(viewer,comm,outviewer);
76: }
77: PetscViewerASCIIPopSynchronized(viewer);
78: return(0);
79: }