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