xref: /petsc/src/ksp/ksp/interface/xmon.c (revision b41ce5d507ea9a58bfa83cf403107a702e77a67d)
1 
2 #include <petsc/private/kspimpl.h>              /*I  "petscksp.h"   I*/
3 #include <petscdraw.h>
4 
5 /*@C
6    KSPMonitorLGResidualNormCreate - Creates a line graph context for use with
7    KSP to monitor convergence of preconditioned residual norms.
8 
9    Collective on KSP
10 
11    Input Parameters:
12 +  comm - communicator context
13 .  host - the X display to open, or null for the local machine
14 .  label - the title to put in the title bar
15 .  x, y - the screen coordinates of the upper left coordinate of
16           the window
17 -  m, n - the screen width and height in pixels
18 
19    Output Parameter:
20 .  lgctx - the drawing context
21 
22    Options Database Key:
23 .  -ksp_monitor_lg_residualnorm - Sets line graph monitor
24 
25    Notes:
26    Use PetscDrawLGDestroy() to destroy this line graph.
27 
28    Level: intermediate
29 
30 .keywords: KSP, monitor, line graph, residual, create
31 
32 .seealso: KSPMonitorSet(), KSPMonitorLGTrueResidualCreate()
33 @*/
34 PetscErrorCode  KSPMonitorLGResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
35 {
36   PetscDraw      draw;
37   PetscErrorCode ierr;
38   PetscDrawAxis  axis;
39   PetscDrawLG    lg;
40 
41   PetscFunctionBegin;
42   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&draw);CHKERRQ(ierr);
43   ierr = PetscDrawSetFromOptions(draw);CHKERRQ(ierr);
44   ierr = PetscDrawLGCreate(draw,1,&lg);CHKERRQ(ierr);
45   ierr = PetscDrawLGSetFromOptions(lg);CHKERRQ(ierr);
46   ierr = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr);
47   ierr = PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norm");CHKERRQ(ierr);
48   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
49   *lgctx = lg;
50   PetscFunctionReturn(0);
51 }
52 
53 PetscErrorCode  KSPMonitorLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx)
54 {
55   PetscDrawLG    lg = (PetscDrawLG) ctx;
56   PetscReal      x,y;
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,4);
61   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
62   x = (PetscReal) n;
63   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
64   else y = -15.0;
65   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
66   if (n <= 20 || !(n % 5) || ksp->reason) {
67     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
68     ierr = PetscDrawLGSave(lg);CHKERRQ(ierr);
69   }
70   PetscFunctionReturn(0);
71 }
72 
73 extern PetscErrorCode  KSPMonitorRange_Private(KSP,PetscInt,PetscReal*);
74 PetscErrorCode  KSPMonitorLGRange(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx)
75 {
76   PetscDrawLG      lg;
77   PetscErrorCode   ierr;
78   PetscReal        x,y,per;
79   PetscViewer      v = (PetscViewer)monctx;
80   static PetscReal prev; /* should be in the context */
81   PetscDraw        draw;
82 
83   PetscFunctionBegin;
84   PetscValidHeaderSpecific(v,PETSC_VIEWER_CLASSID,4);
85 
86   ierr = KSPMonitorRange_Private(ksp,n,&per);CHKERRQ(ierr);
87   if (!n) prev = rnorm;
88 
89   ierr = PetscViewerDrawGetDrawLG(v,0,&lg);CHKERRQ(ierr);
90   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
91   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
92   ierr = PetscDrawSetTitle(draw,"Residual norm");CHKERRQ(ierr);
93   x    = (PetscReal) n;
94   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
95   else y = -15.0;
96   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
97   if (n < 20 || !(n % 5) || ksp->reason) {
98     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
99     ierr = PetscDrawLGSave(lg);CHKERRQ(ierr);
100   }
101 
102   ierr = PetscViewerDrawGetDrawLG(v,1,&lg);CHKERRQ(ierr);
103   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
104   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
105   ierr = PetscDrawSetTitle(draw,"% elemts > .2*max elemt");CHKERRQ(ierr);
106   x    = (PetscReal) n;
107   y    = 100.0*per;
108   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
109   if (n < 20 || !(n % 5) || ksp->reason) {
110     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
111     ierr = PetscDrawLGSave(lg);CHKERRQ(ierr);
112   }
113 
114   ierr = PetscViewerDrawGetDrawLG(v,2,&lg);CHKERRQ(ierr);
115   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
116   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
117   ierr = PetscDrawSetTitle(draw,"(norm-oldnorm)/oldnorm");CHKERRQ(ierr);
118   x    = (PetscReal) n;
119   y    = (prev - rnorm)/prev;
120   ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
121   if (n < 20 || !(n % 5) || ksp->reason) {
122     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
123     ierr = PetscDrawLGSave(lg);CHKERRQ(ierr);
124   }
125 
126   ierr = PetscViewerDrawGetDrawLG(v,3,&lg);CHKERRQ(ierr);
127   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
128   ierr = PetscDrawLGGetDraw(lg,&draw);CHKERRQ(ierr);
129   ierr = PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");CHKERRQ(ierr);
130   x    = (PetscReal) n;
131   y    = (prev - rnorm)/(prev*per);
132   if (n > 5) {
133     ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
134   }
135   if (n < 20 || !(n % 5) || ksp->reason) {
136     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
137     ierr = PetscDrawLGSave(lg);CHKERRQ(ierr);
138   }
139 
140   prev = rnorm;
141   PetscFunctionReturn(0);
142 }
143 
144 /*@C
145    KSPMonitorLGTrueResidualNormCreate - Creates a line graph context for use with
146    KSP to monitor convergence of true residual norms (as opposed to
147    preconditioned residual norms).
148 
149    Collective on KSP
150 
151    Input Parameters:
152 +  comm - communicator context
153 .  host - the X display to open, or null for the local machine
154 .  label - the title to put in the title bar
155 .  x, y - the screen coordinates of the upper left coordinate of
156           the window
157 -  m, n - the screen width and height in pixels
158 
159    Output Parameter:
160 .  lgctx - the drawing context
161 
162    Options Database Key:
163 .  -ksp_monitor_lg_true_residualnorm - Sets true line graph monitor
164 
165    Notes:
166    Use PetscDrawLGDestroy() to destroy this line graph.
167 
168    Level: intermediate
169 
170 .keywords: KSP, monitor, line graph, residual, create, true
171 
172 .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
173 @*/
174 PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
175 {
176   PetscDraw      draw;
177   PetscErrorCode ierr;
178   PetscDrawAxis  axis;
179   PetscDrawLG    lg;
180   const char     *names[] = {"Preconditioned","True"};
181 
182   PetscFunctionBegin;
183   ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&draw);CHKERRQ(ierr);
184   ierr = PetscDrawSetFromOptions(draw);CHKERRQ(ierr);
185   ierr = PetscDrawLGCreate(draw,2,&lg);CHKERRQ(ierr);
186   ierr = PetscDrawLGSetLegend(lg,names);CHKERRQ(ierr);
187   ierr = PetscDrawLGSetFromOptions(lg);CHKERRQ(ierr);
188   ierr = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr);
189   ierr = PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norm");CHKERRQ(ierr);
190   ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr);
191   *lgctx = lg;
192   PetscFunctionReturn(0);
193 }
194 
195 PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx)
196 {
197   PetscDrawLG    lg = (PetscDrawLG) ctx;
198   PetscReal      x[2],y[2],scnorm;
199   Vec            resid,work;
200   PetscErrorCode ierr;
201 
202   PetscFunctionBegin;
203   PetscValidHeaderSpecific(ksp,KSP_CLASSID,1);
204   PetscValidHeaderSpecific(lg,PETSC_DRAWLG_CLASSID,4);
205   if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);}
206   x[0] = x[1] = (PetscReal) n;
207   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
208   else y[0] = -15.0;
209   ierr = VecDuplicate(ksp->vec_rhs,&work);CHKERRQ(ierr);
210   ierr = KSPBuildResidual(ksp,NULL,work,&resid);CHKERRQ(ierr);
211   ierr = VecNorm(resid,NORM_2,&scnorm);CHKERRQ(ierr);
212   ierr = VecDestroy(&work);CHKERRQ(ierr);
213   if (scnorm > 0.0) y[1] = PetscLog10Real(scnorm);
214   else y[1] = -15.0;
215   ierr = PetscDrawLGAddPoint(lg,x,y);CHKERRQ(ierr);
216   if (n <= 20 || !(n % 5)) {
217     ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
218     ierr = PetscDrawLGSave(lg);CHKERRQ(ierr);
219   }
220   PetscFunctionReturn(0);
221 }
222