xref: /petsc/src/ksp/ksp/interface/xmon.c (revision 6ac5842e34eedc6428162d8d42bedaaf46eae34c)
1 
2 #include <petsc-private/kspimpl.h>              /*I  "petscksp.h"   I*/
3 #include <petscdraw.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "KSPMonitorLGResidualNormCreate"
7 /*@C
8    KSPMonitorLGResidualNormCreate - Creates a line graph context for use with
9    KSP to monitor convergence of preconditioned residual norms.
10 
11    Collective on KSP
12 
13    Input Parameters:
14 +  host - the X display to open, or null for the local machine
15 .  label - the title to put in the title bar
16 .  x, y - the screen coordinates of the upper left coordinate of
17           the window
18 -  m, n - the screen width and height in pixels
19 
20    Output Parameter:
21 .  draw - the drawing context
22 
23    Options Database Key:
24 .  -ksp_monitor_lg_residualnorm - Sets line graph monitor
25 
26    Notes:
27    Use KSPMonitorLGResidualNormDestroy() to destroy this line graph; do not use PetscDrawLGDestroy().
28 
29    Level: intermediate
30 
31 .keywords: KSP, monitor, line graph, residual, create
32 
33 .seealso: KSPMonitorLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorLGTrueResidualCreate()
34 @*/
35 PetscErrorCode  KSPMonitorLGResidualNormCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
36 {
37   PetscDraw      win;
38   PetscErrorCode ierr;
39   PetscDrawAxis  axis;
40 
41   PetscFunctionBegin;
42   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
43   ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr);
44   ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr);
45   ierr = PetscDrawLGGetAxis(*draw,&axis);CHKERRQ(ierr);
46   ierr = PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norm");CHKERRQ(ierr);
47   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 #undef __FUNCT__
52 #define __FUNCT__ "KSPMonitorLGResidualNorm"
53 PetscErrorCode  KSPMonitorLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx)
54 {
55   PetscDrawLG    lg = (PetscDrawLG)monctx;
56   PetscErrorCode ierr;
57   PetscReal      x,y;
58 
59   PetscFunctionBegin;
60   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
61   x = (PetscReal) n;
62   if (rnorm > 0.0) y = log10(rnorm);
63   else y = -15.0;
64   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
65   if (n < 20 || !(n % 5)) {
66     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
67   }
68   PetscFunctionReturn(0);
69 }
70 
71 #undef __FUNCT__
72 #define __FUNCT__ "KSPMonitorLGResidualNormDestroy"
73 /*@
74    KSPMonitorLGResidualNormDestroy - Destroys a line graph context that was created
75    with KSPMonitorLGResidualNormCreate().
76 
77    Collective on KSP
78 
79    Input Parameter:
80 .  draw - the drawing context
81 
82    Level: intermediate
83 
84 .keywords: KSP, monitor, line graph, destroy
85 
86 .seealso: KSPMonitorLGResidualNormCreate(), KSPMonitorLGTrueResidualDestroy(), KSPMonitorSet()
87 @*/
88 PetscErrorCode  KSPMonitorLGResidualNormDestroy(PetscDrawLG *drawlg)
89 {
90   PetscDraw      draw;
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
95   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
96   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99 
100 extern PetscErrorCode  KSPMonitorRange_Private(KSP,PetscInt,PetscReal*);
101 #undef __FUNCT__
102 #define __FUNCT__ "KSPMonitorLGRange"
103 PetscErrorCode  KSPMonitorLGRange(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx)
104 {
105   PetscDrawLG      lg;
106   PetscErrorCode   ierr;
107   PetscReal        x,y,per;
108   PetscViewer      v = (PetscViewer)monctx;
109   static PetscReal prev; /* should be in the context */
110   PetscDraw        draw;
111 
112   PetscFunctionBegin;
113   ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
114   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
115   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
116   ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
117   x    = (PetscReal) n;
118   if (rnorm > 0.0) y = log10(rnorm);
119   else y = -15.0;
120   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
121   if (n < 20 || !(n % 5)) {
122     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
123   }
124 
125   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
126   ierr =  KSPMonitorRange_Private(ksp,n,&per);CHKERRQ(ierr);
127   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
128   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
129   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
130   x    = (PetscReal) n;
131   y    = 100.0*per;
132   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
133   if (n < 20 || !(n % 5)) {
134     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
135   }
136 
137   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
138   if (!n) {prev = rnorm;ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
139   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
140   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
141   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
142   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm");CHKERRQ(ierr);
143   x    = (PetscReal) n;
144   y    = (prev - rnorm)/prev;
145   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
146   if (n < 20 || !(n % 5)) {
147     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
148   }
149 
150   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
151   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
152   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
153   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
154   x    = (PetscReal) n;
155   y    = (prev - rnorm)/(prev*per);
156   if (n > 5) {
157     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
158   }
159   if (n < 20 || !(n % 5)) {
160     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
161   }
162   prev = rnorm;
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "KSPMonitorLGTrueResidualNormCreate"
168 /*@C
169    KSPMonitorLGTrueResidualNormCreate - Creates a line graph context for use with
170    KSP to monitor convergence of true residual norms (as opposed to
171    preconditioned residual norms).
172 
173    Collective on KSP
174 
175    Input Parameters:
176 +  host - the X display to open, or null for the local machine
177 .  label - the title to put in the title bar
178 .  x, y - the screen coordinates of the upper left coordinate of
179           the window
180 -  m, n - the screen width and height in pixels
181 
182    Output Parameter:
183 .  draw - the drawing context
184 
185    Options Database Key:
186 .  -ksp_monitor_lg_true_residualnorm - Sets true line graph monitor
187 
188    Notes:
189    Use KSPMonitorLGTrueResidualNormDestroy() to destroy this line graph, not
190    PetscDrawLGDestroy().
191 
192    Level: intermediate
193 
194 .keywords: KSP, monitor, line graph, residual, create, true
195 
196 .seealso: KSPMonitorLGResidualNormDestroy(), KSPMonitorSet(), KSPMonitorDefault()
197 @*/
198 PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
199 {
200   PetscDraw      win;
201   PetscErrorCode ierr;
202   PetscMPIInt    rank;
203   PetscDrawAxis  axis;
204 
205   PetscFunctionBegin;
206   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
207   if (rank) { *draw = 0; PetscFunctionReturn(0);}
208 
209   ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr);
210   ierr = PetscDrawSetFromOptions(win);CHKERRQ(ierr);
211   ierr = PetscDrawLGCreate(win,2,draw);CHKERRQ(ierr);
212   ierr = PetscDrawLGGetAxis(*draw,&axis);CHKERRQ(ierr);
213   ierr = PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norms");CHKERRQ(ierr);
214   ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "KSPMonitorLGTrueResidualNorm"
220 PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx)
221 {
222   PetscDrawLG    lg = (PetscDrawLG) monctx;
223   PetscReal      x[2],y[2],scnorm;
224   PetscErrorCode ierr;
225   PetscMPIInt    rank;
226   Vec            resid,work;
227 
228   PetscFunctionBegin;
229   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)ksp),&rank);CHKERRQ(ierr);
230   if (!rank) {
231     if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
232     x[0] = x[1] = (PetscReal) n;
233     if (rnorm > 0.0) y[0] = log10(rnorm);
234     else y[0] = -15.0;
235   }
236 
237   ierr = VecDuplicate(ksp->vec_rhs,&work);CHKERRQ(ierr);
238   ierr = KSPBuildResidual(ksp,0,work,&resid);CHKERRQ(ierr);
239   ierr = VecNorm(resid,NORM_2,&scnorm);CHKERRQ(ierr);
240   ierr = VecDestroy(&work);CHKERRQ(ierr);
241 
242   if (!rank) {
243     if (scnorm > 0.0) y[1] = log10(scnorm);
244     else y[1] = -15.0;
245     ierr = PetscDrawLGAddPoint(lg,x,y);CHKERRQ(ierr);
246     if (n <= 20 || (n % 3)) {
247       ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
248     }
249   }
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "KSPMonitorLGTrueResidualNormDestroy"
255 /*@C
256    KSPMonitorLGTrueResidualNormDestroy - Destroys a line graph context that was created
257    with KSPMonitorLGTrueResidualNormCreate().
258 
259    Collective on KSP
260 
261    Input Parameter:
262 .  draw - the drawing context
263 
264    Level: intermediate
265 
266 .keywords: KSP, monitor, line graph, destroy, true
267 
268 .seealso: KSPMonitorLGTrueResidualNormCreate(), KSPMonitorSet()
269 @*/
270 PetscErrorCode  KSPMonitorLGTrueResidualNormDestroy(PetscDrawLG *drawlg)
271 {
272   PetscErrorCode ierr;
273   PetscDraw      draw;
274 
275   PetscFunctionBegin;
276   ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr);
277   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
278   ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 
283