Actual source code: xmon.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <petsc/private/kspimpl.h>
  3:  #include <petscdraw.h>

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

  9:    Collective

 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

 19:    Output Parameter:
 20: .  lgctx - the drawing context

 22:    Options Database Key:
 23: .  -ksp_monitor_lg_residualnorm - Sets line graph monitor

 25:    Notes:
 26:    Use PetscDrawLGDestroy() to destroy this line graph.

 28:    Level: intermediate

 30: .seealso: KSPMonitorSet(), KSPMonitorLGTrueResidualCreate()
 31: @*/
 32: PetscErrorCode  KSPMonitorLGResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
 33: {
 34:   PetscDraw      draw;
 36:   PetscDrawAxis  axis;
 37:   PetscDrawLG    lg;

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

 51: PetscErrorCode  KSPMonitorLGResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx)
 52: {
 53:   PetscDrawLG    lg = (PetscDrawLG) ctx;
 54:   PetscReal      x,y;

 59:   if (!n) {PetscDrawLGReset(lg);}
 60:   x = (PetscReal) n;
 61:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
 62:   else y = -15.0;
 63:   PetscDrawLGAddPoint(lg,&x,&y);
 64:   if (n <= 20 || !(n % 5) || ksp->reason) {
 65:     PetscDrawLGDraw(lg);
 66:     PetscDrawLGSave(lg);
 67:   }
 68:   return(0);
 69: }

 71: extern PetscErrorCode  KSPMonitorRange_Private(KSP,PetscInt,PetscReal*);
 72: PetscErrorCode  KSPMonitorLGRange(KSP ksp,PetscInt n,PetscReal rnorm,void *monctx)
 73: {
 74:   PetscDrawLG      lg;
 75:   PetscErrorCode   ierr;
 76:   PetscReal        x,y,per;
 77:   PetscViewer      v = (PetscViewer)monctx;
 78:   static PetscReal prev; /* should be in the context */
 79:   PetscDraw        draw;


 84:   KSPMonitorRange_Private(ksp,n,&per);
 85:   if (!n) prev = rnorm;

 87:   PetscViewerDrawGetDrawLG(v,0,&lg);
 88:   if (!n) {PetscDrawLGReset(lg);}
 89:   PetscDrawLGGetDraw(lg,&draw);
 90:   PetscDrawSetTitle(draw,"Residual norm");
 91:   x    = (PetscReal) n;
 92:   if (rnorm > 0.0) y = PetscLog10Real(rnorm);
 93:   else y = -15.0;
 94:   PetscDrawLGAddPoint(lg,&x,&y);
 95:   if (n < 20 || !(n % 5) || ksp->reason) {
 96:     PetscDrawLGDraw(lg);
 97:     PetscDrawLGSave(lg);
 98:   }

100:   PetscViewerDrawGetDrawLG(v,1,&lg);
101:   if (!n) {PetscDrawLGReset(lg);}
102:   PetscDrawLGGetDraw(lg,&draw);
103:   PetscDrawSetTitle(draw,"% elemts > .2*max elemt");
104:   x    = (PetscReal) n;
105:   y    = 100.0*per;
106:   PetscDrawLGAddPoint(lg,&x,&y);
107:   if (n < 20 || !(n % 5) || ksp->reason) {
108:     PetscDrawLGDraw(lg);
109:     PetscDrawLGSave(lg);
110:   }

112:   PetscViewerDrawGetDrawLG(v,2,&lg);
113:   if (!n) {PetscDrawLGReset(lg);}
114:   PetscDrawLGGetDraw(lg,&draw);
115:   PetscDrawSetTitle(draw,"(norm-oldnorm)/oldnorm");
116:   x    = (PetscReal) n;
117:   y    = (prev - rnorm)/prev;
118:   PetscDrawLGAddPoint(lg,&x,&y);
119:   if (n < 20 || !(n % 5) || ksp->reason) {
120:     PetscDrawLGDraw(lg);
121:     PetscDrawLGSave(lg);
122:   }

124:   PetscViewerDrawGetDrawLG(v,3,&lg);
125:   if (!n) {PetscDrawLGReset(lg);}
126:   PetscDrawLGGetDraw(lg,&draw);
127:   PetscDrawSetTitle(draw,"(norm -oldnorm)/oldnorm*(% > .2 max)");
128:   x    = (PetscReal) n;
129:   y    = (prev - rnorm)/(prev*per);
130:   if (n > 5) {
131:     PetscDrawLGAddPoint(lg,&x,&y);
132:   }
133:   if (n < 20 || !(n % 5) || ksp->reason) {
134:     PetscDrawLGDraw(lg);
135:     PetscDrawLGSave(lg);
136:   }

138:   prev = rnorm;
139:   return(0);
140: }

142: /*@C
143:    KSPMonitorLGTrueResidualNormCreate - Creates a line graph context for use with
144:    KSP to monitor convergence of true residual norms (as opposed to
145:    preconditioned residual norms).

147:    Collective

149:    Input Parameters:
150: +  comm - communicator context
151: .  host - the X display to open, or null for the local machine
152: .  label - the title to put in the title bar
153: .  x, y - the screen coordinates of the upper left coordinate of
154:           the window
155: -  m, n - the screen width and height in pixels

157:    Output Parameter:
158: .  lgctx - the drawing context

160:    Options Database Key:
161: .  -ksp_monitor_lg_true_residualnorm - Sets true line graph monitor

163:    Notes:
164:    Use PetscDrawLGDestroy() to destroy this line graph.

166:    Level: intermediate

168: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
169: @*/
170: PetscErrorCode  KSPMonitorLGTrueResidualNormCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *lgctx)
171: {
172:   PetscDraw      draw;
174:   PetscDrawAxis  axis;
175:   PetscDrawLG    lg;
176:   const char     *names[] = {"Preconditioned","True"};

179:   PetscDrawCreate(comm,host,label,x,y,m,n,&draw);
180:   PetscDrawSetFromOptions(draw);
181:   PetscDrawLGCreate(draw,2,&lg);
182:   PetscDrawLGSetLegend(lg,names);
183:   PetscDrawLGSetFromOptions(lg);
184:   PetscDrawLGGetAxis(lg,&axis);
185:   PetscDrawAxisSetLabels(axis,"Convergence","Iteration","Residual Norm");
186:   PetscDrawDestroy(&draw);
187:   *lgctx = lg;
188:   return(0);
189: }

191: PetscErrorCode  KSPMonitorLGTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *ctx)
192: {
193:   PetscDrawLG    lg = (PetscDrawLG) ctx;
194:   PetscReal      x[2],y[2],scnorm;
195:   Vec            resid,work;

201:   if (!n) {PetscDrawLGReset(lg);}
202:   x[0] = x[1] = (PetscReal) n;
203:   if (rnorm > 0.0) y[0] = PetscLog10Real(rnorm);
204:   else y[0] = -15.0;
205:   VecDuplicate(ksp->vec_rhs,&work);
206:   KSPBuildResidual(ksp,NULL,work,&resid);
207:   VecNorm(resid,NORM_2,&scnorm);
208:   VecDestroy(&work);
209:   if (scnorm > 0.0) y[1] = PetscLog10Real(scnorm);
210:   else y[1] = -15.0;
211:   PetscDrawLGAddPoint(lg,x,y);
212:   if (n <= 20 || !(n % 5)) {
213:     PetscDrawLGDraw(lg);
214:     PetscDrawLGSave(lg);
215:   }
216:   return(0);
217: }