Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
525 dsic.upv.es!antodo 1
/*
545 dsic.upv.es!jroman 2
      EPS routines related to monitors.
525 dsic.upv.es!antodo 3
*/
4
#include "src/eps/epsimpl.h"   /*I "slepceps.h" I*/
5
 
6
#undef __FUNCT__  
7
#define __FUNCT__ "EPSSetMonitor"
8
/*@C
9
   EPSSetMonitor - Sets an ADDITIONAL function to be called at every
10
   iteration to monitor the error estimates for each requested eigenpair.
11
 
12
   Collective on EPS
13
 
14
   Input Parameters:
15
+  eps     - eigensolver context obtained from EPSCreate()
1021 slepc 16
.  monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring)
525 dsic.upv.es!antodo 17
-  mctx    - [optional] context for private data for the
18
             monitor routine (use PETSC_NULL if no context is desired)
19
 
20
   Calling Sequence of monitor:
617 dsic.upv.es!antodo 21
$     monitor (EPS eps, int its, int nconv, PetscScalar *eigr, PetscScalar *eigi, PetscReal* errest, int nest, void *mctx)
525 dsic.upv.es!antodo 22
 
23
+  eps    - eigensolver context obtained from EPSCreate()
24
.  its    - iteration number
25
.  nconv  - number of converged eigenpairs
617 dsic.upv.es!antodo 26
.  eigr   - real part of the eigenvalues
27
.  eigi   - imaginary part of the eigenvalues
663 dsic.upv.es!antodo 28
.  errest - relative error estimates for each eigenpair
525 dsic.upv.es!antodo 29
.  nest   - number of error estimates
30
-  mctx   - optional monitoring context, as set by EPSSetMonitor()
31
 
32
   Options Database Keys:
33
+    -eps_monitor        - print error estimates at each iteration
34
-    -eps_cancelmonitors - cancels all monitors that have been hardwired into
35
      a code by calls to EPSetMonitor(), but does not cancel those set via
36
      the options database.
37
 
38
   Notes:  
39
   Several different monitoring routines may be set by calling
40
   EPSSetMonitor() multiple times; all will be called in the
41
   order in which they were set.
42
 
43
   Level: intermediate
44
 
1021 slepc 45
.seealso: EPSDefaultMonitor(), EPSClearMonitor()
525 dsic.upv.es!antodo 46
@*/
47
PetscErrorCode EPSSetMonitor(EPS eps, int (*monitor)(EPS,int,int,PetscScalar*,PetscScalar*,PetscReal*,int,void*), void *mctx)
48
{
49
  PetscFunctionBegin;
50
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
51
  if (eps->numbermonitors >= MAXEPSMONITORS) {
52
    SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
53
  }
54
  eps->monitor[eps->numbermonitors]           = monitor;
55
  eps->monitorcontext[eps->numbermonitors++]  = (void*)mctx;
56
  PetscFunctionReturn(0);
57
}
58
 
59
#undef __FUNCT__  
60
#define __FUNCT__ "EPSClearMonitor"
1021 slepc 61
/*@
525 dsic.upv.es!antodo 62
   EPSClearMonitor - Clears all monitors for an EPS object.
63
 
64
   Collective on EPS
65
 
66
   Input Parameters:
67
.  eps - eigensolver context obtained from EPSCreate()
68
 
69
   Options Database Key:
70
.    -eps_cancelmonitors - Cancels all monitors that have been hardwired
1021 slepc 71
      into a code by calls to EPSSetMonitor(),
525 dsic.upv.es!antodo 72
      but does not cancel those set via the options database.
73
 
74
   Level: intermediate
75
 
623 dsic.upv.es!antodo 76
.seealso: EPSSetMonitor()
525 dsic.upv.es!antodo 77
@*/
78
PetscErrorCode EPSClearMonitor(EPS eps)
79
{
80
  PetscFunctionBegin;
81
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
82
  eps->numbermonitors = 0;
83
  PetscFunctionReturn(0);
84
}
85
 
86
#undef __FUNCT__  
87
#define __FUNCT__ "EPSGetMonitorContext"
88
/*@C
1021 slepc 89
   EPSGetMonitorContext - Gets the monitor context, as set by
525 dsic.upv.es!antodo 90
   EPSSetMonitor() for the FIRST monitor only.
91
 
92
   Not Collective
93
 
94
   Input Parameter:
95
.  eps - eigensolver context obtained from EPSCreate()
96
 
97
   Output Parameter:
98
.  ctx - monitor context
99
 
100
   Level: intermediate
101
 
1021 slepc 102
.seealso: EPSSetMonitor(), EPSDefaultMonitor()
525 dsic.upv.es!antodo 103
@*/
104
PetscErrorCode EPSGetMonitorContext(EPS eps, void **ctx)
105
{
106
  PetscFunctionBegin;
107
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
108
  *ctx =      (eps->monitorcontext[0]);
109
  PetscFunctionReturn(0);
110
}
111
 
112
#undef __FUNCT__  
113
#define __FUNCT__ "EPSDefaultMonitor"
114
/*@C
1021 slepc 115
   EPSDefaultMonitor - Print the current approximate values and
525 dsic.upv.es!antodo 116
   error estimates at each iteration of the eigensolver.
117
 
118
   Collective on EPS
119
 
120
   Input Parameters:
121
+  eps    - eigensolver context
122
.  its    - iteration number
123
.  nconv  - number of converged eigenpairs so far
124
.  eigr   - real part of the eigenvalues
125
.  eigi   - imaginary part of the eigenvalues
126
.  errest - error estimates
127
.  nest   - number of error estimates to display
128
-  dummy  - unused monitor context
129
 
130
   Level: intermediate
131
 
132
.seealso: EPSSetMonitor()
133
@*/
134
PetscErrorCode EPSDefaultMonitor(EPS eps,int its,int nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,int nest,void *dummy)
135
{
136
  PetscErrorCode ierr;
137
  int            i;
138
  PetscViewer    viewer = (PetscViewer) dummy;
139
 
140
  PetscFunctionBegin;
141
  if (!viewer) viewer = PETSC_VIEWER_STDOUT_(eps->comm);
142
  ierr = PetscViewerASCIIPrintf(viewer,"%3d EPS nconv=%d Values (Errors)",its,nconv);CHKERRQ(ierr);
143
  for (i=0;i<nest;i++) {
144
#if defined(PETSC_USE_COMPLEX)
145
    ierr = PetscViewerASCIIPrintf(viewer," %g%+gi",PetscRealPart(eigr[i]),PetscImaginaryPart(eigr[i]));CHKERRQ(ierr);
146
#else
147
    ierr = PetscViewerASCIIPrintf(viewer," %g",eigr[i]);CHKERRQ(ierr);
148
    if (eigi[i]!=0.0) { ierr = PetscViewerASCIIPrintf(viewer,"%+gi",eigi[i]);CHKERRQ(ierr); }
149
#endif
150
    ierr = PetscViewerASCIIPrintf(viewer," (%10.8e)",errest[i]);CHKERRQ(ierr);
151
  }
152
  ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
153
  PetscFunctionReturn(0);
154
}
623 dsic.upv.es!antodo 155
 
156
#undef __FUNCT__  
157
#define __FUNCT__ "EPSLGMonitor"
158
PetscErrorCode EPSLGMonitor(EPS eps,int its,int nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,int nest,void *monctx)
159
{
626 dsic.upv.es!antodo 160
  PetscViewer    viewer = (PetscViewer) monctx;
161
  PetscDraw      draw;
162
  PetscDrawLG    lg;
623 dsic.upv.es!antodo 163
  PetscErrorCode ierr;
164
  PetscReal      *x,*y;
165
  int            i,n = eps->nev;
626 dsic.upv.es!antodo 166
#if !defined(PETSC_USE_COMPLEX)
1004 slepc 167
  int            p;
626 dsic.upv.es!antodo 168
  PetscDraw      draw1;
169
  PetscDrawLG    lg1;
170
#endif
623 dsic.upv.es!antodo 171
 
172
  PetscFunctionBegin;
173
 
626 dsic.upv.es!antodo 174
  if (!viewer) { viewer = PETSC_VIEWER_DRAW_(eps->comm); }
175
 
176
  ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
177
  ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
178
  if (!its) {
628 dsic.upv.es!antodo 179
    ierr = PetscDrawSetTitle(draw,"Error estimates");CHKERRQ(ierr);
626 dsic.upv.es!antodo 180
    ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr);
181
    ierr = PetscDrawLGSetDimension(lg,n);CHKERRQ(ierr);
182
    ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
183
    ierr = PetscDrawLGSetLimits(lg,0,1.0,log10(eps->tol)-2,0.0);CHKERRQ(ierr);
184
  }
185
 
186
#if !defined(PETSC_USE_COMPLEX)
187
  if (eps->ishermitian) {
188
    ierr = PetscViewerDrawGetDraw(viewer,1,&draw1);CHKERRQ(ierr);
189
    ierr = PetscViewerDrawGetDrawLG(viewer,1,&lg1);CHKERRQ(ierr);
623 dsic.upv.es!antodo 190
    if (!its) {
628 dsic.upv.es!antodo 191
      ierr = PetscDrawSetTitle(draw1,"Approximate eigenvalues");CHKERRQ(ierr);
626 dsic.upv.es!antodo 192
      ierr = PetscDrawSetDoubleBuffer(draw1);CHKERRQ(ierr);
193
      ierr = PetscDrawLGSetDimension(lg1,n);CHKERRQ(ierr);
194
      ierr = PetscDrawLGReset(lg1);CHKERRQ(ierr);
195
      ierr = PetscDrawLGSetLimits(lg1,0,1.0,1.e20,-1.e20);CHKERRQ(ierr);
623 dsic.upv.es!antodo 196
    }
197
  }
626 dsic.upv.es!antodo 198
#endif
623 dsic.upv.es!antodo 199
 
200
  ierr = PetscMalloc(sizeof(PetscReal)*n,&x);CHKERRQ(ierr);
201
  ierr = PetscMalloc(sizeof(PetscReal)*n,&y);CHKERRQ(ierr);
202
  for (i=0;i<n;i++) {
203
    x[i] = (PetscReal) its;
204
    if (errest[i] > 0.0) y[i] = log10(errest[i]); else y[i] = 0.0;
205
  }
206
  ierr = PetscDrawLGAddPoint(lg,x,y);CHKERRQ(ierr);
626 dsic.upv.es!antodo 207
#if !defined(PETSC_USE_COMPLEX)
208
  if (eps->ishermitian) {
209
    ierr = PetscDrawLGAddPoint(lg1,x,eps->eigr);CHKERRQ(ierr);
1004 slepc 210
    ierr = PetscDrawGetPause(draw1,&p);CHKERRQ(ierr);
626 dsic.upv.es!antodo 211
    ierr = PetscDrawSetPause(draw1,0);CHKERRQ(ierr);    
212
    ierr = PetscDrawLGDraw(lg1);CHKERRQ(ierr);
1004 slepc 213
    ierr = PetscDrawSetPause(draw1,p);CHKERRQ(ierr);    
626 dsic.upv.es!antodo 214
  }
215
#endif  
705 dsic.upv.es!antodo 216
  ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
623 dsic.upv.es!antodo 217
  ierr = PetscFree(x);CHKERRQ(ierr);
218
  ierr = PetscFree(y);CHKERRQ(ierr);  
219
  PetscFunctionReturn(0);
220
}
1021 slepc 221