xref: /petsc/src/sys/classes/viewer/interface/dupl.c (revision 9895aa37ac365bac650f6bd8bf977519f7222510)
1 
2 #include <petsc-private/viewerimpl.h>  /*I "petscviewer.h" I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "PetscViewerGetSingleton"
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)
9 
10     Collective on PetscViewer
11 
12    Input Parameter:
13 .  viewer - the PetscViewer to be duplicated
14 
15    Output Parameter:
16 .  outviewer - new PetscViewer
17 
18    Level: advanced
19 
20    Notes: Call PetscViewerRestoreSingleton() to return this PetscViewer, NOT PetscViewerDestroy()
21 
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.
25 
26    Concepts: PetscViewer^sequential version
27 
28 .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerRestoreSingleton()
29 @*/
30 PetscErrorCode  PetscViewerGetSingleton(PetscViewer viewer,PetscViewer *outviewer)
31 {
32   PetscErrorCode ierr;
33   PetscMPIInt    size;
34 
35   PetscFunctionBegin;
36   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
37   PetscValidPointer(outviewer,2);
38   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
39   if (size == 1) {
40     ierr       = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
41     *outviewer = viewer;
42   } else if (viewer->ops->getsingleton) {
43     ierr = (*viewer->ops->getsingleton)(viewer,outviewer);CHKERRQ(ierr);
44   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get singleton PetscViewer for type %s",((PetscObject)viewer)->type_name);
45   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "PetscViewerRestoreSingleton"
51 /*@
52    PetscViewerRestoreSingleton - Restores a new PetscViewer obtained with PetscViewerGetSingleton().
53 
54     Collective on PetscViewer
55 
56    Input Parameters:
57 +  viewer - the PetscViewer to be duplicated
58 -  outviewer - new PetscViewer
59 
60    Level: advanced
61 
62    Notes: Call PetscViewerGetSingleton() to get this PetscViewer, NOT PetscViewerCreate()
63 
64 .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerGetSingleton()
65 @*/
66 PetscErrorCode  PetscViewerRestoreSingleton(PetscViewer viewer,PetscViewer *outviewer)
67 {
68   PetscErrorCode ierr;
69   PetscMPIInt    size;
70 
71   PetscFunctionBegin;
72   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
73 
74   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
75   if (size == 1) {
76     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
77     if (outviewer) *outviewer = 0;
78   } else if (viewer->ops->restoresingleton) {
79     ierr = (*viewer->ops->restoresingleton)(viewer,outviewer);CHKERRQ(ierr);
80   }
81   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
82   PetscFunctionReturn(0);
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "PetscViewerGetSubcomm"
87 /*@
88    PetscViewerGetSubcomm - Creates a new PetscViewer (same type as the old)
89     that lives on a subgroup of processors
90 
91     Collective on PetscViewer
92 
93    Input Parameter:
94 +  viewer - the PetscViewer to be duplicated
95 -  subcomm - MPI communicator
96 
97    Output Parameter:
98 .  outviewer - new PetscViewer
99 
100    Level: advanced
101 
102    Notes: Call PetscViewerRestoreSubcomm() to return this PetscViewer, NOT PetscViewerDestroy()
103 
104      This is used to view a sequential or a parallel object that is part of a larger
105     parallel object. For example redundant PC view could use this to obtain a
106     PetscViewer that is used within a subcommunicator on one duplicated preconditioner.
107 
108    Concepts: PetscViewer^sequential version
109 
110 .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerRestoreSubcomm()
111 @*/
112 PetscErrorCode  PetscViewerGetSubcomm(PetscViewer viewer,MPI_Comm subcomm,PetscViewer *outviewer)
113 {
114   PetscErrorCode ierr;
115   PetscMPIInt    size;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
119   PetscValidPointer(outviewer,3);
120   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
121   if (size == 1) {
122     ierr       = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr);
123     *outviewer = viewer;
124   } else if (viewer->ops->getsubcomm) {
125     ierr = (*viewer->ops->getsubcomm)(viewer,subcomm,outviewer);CHKERRQ(ierr);
126   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get subcommunicator PetscViewer for type %s",((PetscObject)viewer)->type_name);
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "PetscViewerRestoreSubcomm"
132 /*@
133    PetscViewerRestoreSubcomm - Restores a new PetscViewer obtained with PetscViewerGetSubcomm().
134 
135     Collective on PetscViewer
136 
137    Input Parameters:
138 +  viewer - the PetscViewer to be duplicated
139 .  subcomm - MPI communicator
140 -  outviewer - new PetscViewer
141 
142    Level: advanced
143 
144    Notes: Call PetscViewerGetSubcomm() to get this PetscViewer, NOT PetscViewerCreate()
145 
146 .seealso: PetscViewerSocketOpen(), PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscViewerGetSubcomm()
147 @*/
148 PetscErrorCode  PetscViewerRestoreSubcomm(PetscViewer viewer,MPI_Comm subcomm,PetscViewer *outviewer)
149 {
150   PetscErrorCode ierr;
151   PetscMPIInt    size;
152 
153   PetscFunctionBegin;
154   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
155 
156   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)viewer),&size);CHKERRQ(ierr);
157   if (size == 1 || (outviewer && viewer == *outviewer)) {
158     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
159     if (outviewer) *outviewer = 0;
160   } else if (viewer->ops->restoresubcomm) {
161     ierr = (*viewer->ops->restoresubcomm)(viewer,subcomm,outviewer);CHKERRQ(ierr);
162   }
163   PetscFunctionReturn(0);
164 }
165 
166