Actual source code: xmon.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <petsc/private/kspimpl.h>              /*I  "petscksp.h"   I*/
  3: #include <petscdraw.h>

  7: /*@C
  8:    KSPMonitorLGResidualNormCreate - Creates a line graph context for use with
  9:    KSP to monitor convergence of preconditioned residual norms.

 11:    Collective on KSP

 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

 21:    Output Parameter:
 22: .  lgctx - the drawing context

 24:    Options Database Key:
 25: .  -ksp_monitor_lg_residualnorm - Sets line graph monitor

 27:    Notes:
 28:    Use PetscDrawLGDestroy() to destroy this line graph.

 30:    Level: intermediate

 32: .keywords: KSP, monitor, line graph, residual, create

 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;
 40:   PetscDrawAxis  axis;
 41:   PetscDrawLG    lg;

 44:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
 45:   PetscDrawSetFromOptions(draw);
 46:   PetscDrawLGCreate(draw,1,&lg);
 47:   PetscDrawLGSetFromOptions(lg);
 48:   PetscDrawLGGetAxis(lg,&axis);
 49:   PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norm");
 50:   PetscDrawDestroy(&draw);
 51:   *lgctx = lg;
 52:   return(0);
 53: }

 57: PetscErrorCode  KSPMonitorLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx)
 58: {
 59:   PetscDrawLG    lg = (PetscDrawLG) ctx;
 60:   PetscReal      x,y;

 65:   if (!n) {PetscDrawLGReset(lg);}
 66:   x = (PetscReal) n;
 67:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
 68:   else y = -15.0;
 69:   PetscDrawLGAddPoint(lg,&x,&y);
 70:   if (n <= 20 || !(n % 5) || ksp->reason) {
 71:     PetscDrawLGDraw(lg);
 72:     PetscDrawLGSave(lg);
 73:   }
 74:   return(0);
 75: }

 77: extern PetscErrorCode  KSPMonitorRange_Private(KSP,PetscInt,PetscReal*);
 80: PetscErrorCode  KSPMonitorLGRange(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx)
 81: {
 82:   PetscDrawLG      lg;
 83:   PetscErrorCode   ierr;
 84:   PetscReal        x,y,per;
 85:   PetscViewer      v = (PetscViewer)monctx;
 86:   static PetscReal prev; /* should be in the context */
 87:   PetscDraw        draw;


 92:   KSPMonitorRange_Private(ksp,n,&per);
 93:   if (!n) prev = rnorm;

 95:   PetscViewerDrawGetDrawLG(v,0,&lg);
 96:   if (!n) {PetscDrawLGReset(lg);}
 97:   PetscDrawLGGetDraw(lg,&draw);
 98:   PetscDrawSetTitle(draw,"Residual norm");
 99:   x    = (PetscReal) n;
100:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
101:   else y = -15.0;
102:   PetscDrawLGAddPoint(lg,&x,&y);
103:   if (n < 20 || !(n % 5) || ksp->reason) {
104:     PetscDrawLGDraw(lg);
105:     PetscDrawLGSave(lg);
106:   }

108:   PetscViewerDrawGetDrawLG(v,1,&lg);
109:   if (!n) {PetscDrawLGReset(lg);}
110:   PetscDrawLGGetDraw(lg,&draw);
111:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
112:   x    = (PetscReal) n;
113:   y    = 100.0*per;
114:   PetscDrawLGAddPoint(lg,&x,&y);
115:   if (n < 20 || !(n % 5) || ksp->reason) {
116:     PetscDrawLGDraw(lg);
117:     PetscDrawLGSave(lg);
118:   }

120:   PetscViewerDrawGetDrawLG(v,2,&lg);
121:   if (!n) {PetscDrawLGReset(lg);}
122:   PetscDrawLGGetDraw(lg,&draw);
123:   PetscDrawSetTitle(draw,"(norm-oldnorm)/oldnorm");
124:   x    = (PetscReal) n;
125:   y    = (prev - rnorm)/prev;
126:   PetscDrawLGAddPoint(lg,&x,&y);
127:   if (n < 20 || !(n % 5) || ksp->reason) {
128:     PetscDrawLGDraw(lg);
129:     PetscDrawLGSave(lg);
130:   }

132:   PetscViewerDrawGetDrawLG(v,3,&lg);
133:   if (!n) {PetscDrawLGReset(lg);}
134:   PetscDrawLGGetDraw(lg,&draw);
135:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
136:   x    = (PetscReal) n;
137:   y    = (prev - rnorm)/(prev*per);
138:   if (n > 5) {
139:     PetscDrawLGAddPoint(lg,&x,&y);
140:   }
141:   if (n < 20 || !(n % 5) || ksp->reason) {
142:     PetscDrawLGDraw(lg);
143:     PetscDrawLGSave(lg);
144:   }

146:   prev = rnorm;
147:   return(0);
148: }

152: /*@C
153:    KSPMonitorLGTrueResidualNormCreate - Creates a line graph context for use with
154:    KSP to monitor convergence of true residual norms (as opposed to
155:    preconditioned residual norms).

157:    Collective on KSP

159:    Input Parameters:
160: +  comm - communicator context
161: .  host - the X display to open, or null for the local machine
162: .  label - the title to put in the title bar
163: .  x, y - the screen coordinates of the upper left coordinate of
164:           the window
165: -  m, n - the screen width and height in pixels

167:    Output Parameter:
168: .  lgctx - the drawing context

170:    Options Database Key:
171: .  -ksp_monitor_lg_true_residualnorm - Sets true line graph monitor

173:    Notes:
174:    Use PetscDrawLGDestroy() to destroy this line graph.

176:    Level: intermediate

178: .keywords: KSP, monitor, line graph, residual, create, true

180: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
181: @*/
182: PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
183: {
184:   PetscDraw      draw;
186:   PetscDrawAxis  axis;
187:   PetscDrawLG    lg;
188:   const char     *names[] = {"Preconditioned","True"};

191:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
192:   PetscDrawSetFromOptions(draw);
193:   PetscDrawLGCreate(draw,2,&lg);
194:   PetscDrawLGSetLegend(lg,names);
195:   PetscDrawLGSetFromOptions(lg);
196:   PetscDrawLGGetAxis(lg,&axis);
197:   PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norm");
198:   PetscDrawDestroy(&draw);
199:   *lgctx = lg;
200:   return(0);
201: }

205: PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx)
206: {
207:   PetscDrawLG    lg = (PetscDrawLG) ctx;
208:   PetscReal      x[2],y[2],scnorm;
209:   Vec            resid,work;

215:   if (!n) {PetscDrawLGReset(lg);}
216:   x[0] = x[1] = (PetscReal) n;
217:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
218:   else y[0] = -15.0;
219:   VecDuplicate(ksp->vec_rhs,&work);
220:   KSPBuildResidual(ksp,NULL,work,&resid);
221:   VecNorm(resid,NORM_2,&scnorm);
222:   VecDestroy(&work);
223:   if (scnorm > 0.0) y[1] = PetscLog10Real(scnorm);
224:   else y[1] = -15.0;
225:   PetscDrawLGAddPoint(lg,x,y);
226:   if (n <= 20 || !(n % 5)) {
227:     PetscDrawLGDraw(lg);
228:     PetscDrawLGSave(lg);
229:   }
230:   return(0);
231: }