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
1887 jroman 1
/*
2
      QEP routines related to monitors.
3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 6
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
1887 jroman 7
 
8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22
*/
23
 
2284 jroman 24
#include <private/qepimpl.h>      /*I "slepcqep.h" I*/
1887 jroman 25
 
26
#undef __FUNCT__  
2313 jroman 27
#define __FUNCT__ "QEPMonitor"
28
/*
29
   Runs the user provided monitor routines, if any.
30
*/
31
PetscErrorCode QEPMonitor(QEP qep,PetscInt it,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest)
32
{
33
  PetscErrorCode ierr;
34
  PetscInt       i,n = qep->numbermonitors;
35
 
36
  PetscFunctionBegin;
37
  for (i=0;i<n;i++) {
38
    ierr = (*qep->monitor[i])(qep,it,nconv,eigr,eigi,errest,nest,qep->monitorcontext[i]);CHKERRQ(ierr);
39
  }
40
  PetscFunctionReturn(0);
41
}
42
 
43
#undef __FUNCT__  
1887 jroman 44
#define __FUNCT__ "QEPMonitorSet"
45
/*@C
46
   QEPMonitorSet - Sets an ADDITIONAL function to be called at every
47
   iteration to monitor the error estimates for each requested eigenpair.
48
 
2328 jroman 49
   Logically Collective on QEP
1887 jroman 50
 
51
   Input Parameters:
52
+  qep     - eigensolver context obtained from QEPCreate()
53
.  monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring)
54
.  mctx    - [optional] context for private data for the
55
             monitor routine (use PETSC_NULL if no context is desired)
56
-  monitordestroy - [optional] routine that frees monitor context
57
          (may be PETSC_NULL)
58
 
59
   Calling Sequence of monitor:
60
$     monitor (QEP qep, int its, int nconv, PetscScalar *eigr, PetscScalar *eigi, PetscReal* errest, int nest, void *mctx)
61
 
62
+  qep    - quadratic eigensolver context obtained from QEPCreate()
63
.  its    - iteration number
64
.  nconv  - number of converged eigenpairs
65
.  eigr   - real part of the eigenvalues
66
.  eigi   - imaginary part of the eigenvalues
67
.  errest - relative error estimates for each eigenpair
68
.  nest   - number of error estimates
69
-  mctx   - optional monitoring context, as set by QEPMonitorSet()
70
 
71
   Options Database Keys:
2054 eromero 72
+    -qep_monitor          - print only the first error estimate
73
.    -qep_monitor_all      - print error estimates at each iteration
74
.    -qep_monitor_conv     - print the eigenvalue approximations only when
1887 jroman 75
      convergence has been reached
2054 eromero 76
.    -qep_monitor_draw     - sets line graph monitor for the first unconverged
77
      approximate eigenvalue
78
.    -qep_monitor_draw_all - sets line graph monitor for all unconverged
79
      approximate eigenvalue
80
-    -qep_monitor_cancel   - cancels all monitors that have been hardwired into
1887 jroman 81
      a code by calls to QEPMonitorSet(), but does not cancel those set via
82
      the options database.
83
 
84
   Notes:  
85
   Several different monitoring routines may be set by calling
86
   QEPMonitorSet() multiple times; all will be called in the
87
   order in which they were set.
88
 
89
   Level: intermediate
90
 
2654 jroman 91
.seealso: QEPMonitorFirst(), QEPMonitorAll(), QEPMonitorCancel()
1887 jroman 92
@*/
2351 jroman 93
PetscErrorCode QEPMonitorSet(QEP qep,PetscErrorCode (*monitor)(QEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1887 jroman 94
{
95
  PetscFunctionBegin;
2213 jroman 96
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 97
  if (qep->numbermonitors >= MAXQEPMONITORS) {
2214 jroman 98
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Too many QEP monitors set");
1887 jroman 99
  }
100
  qep->monitor[qep->numbermonitors]           = monitor;
101
  qep->monitorcontext[qep->numbermonitors]    = (void*)mctx;
102
  qep->monitordestroy[qep->numbermonitors++]  = monitordestroy;
103
  PetscFunctionReturn(0);
104
}
105
 
106
#undef __FUNCT__  
107
#define __FUNCT__ "QEPMonitorCancel"
108
/*@
109
   QEPMonitorCancel - Clears all monitors for a QEP object.
110
 
2328 jroman 111
   Logically Collective on QEP
1887 jroman 112
 
113
   Input Parameters:
114
.  qep - eigensolver context obtained from QEPCreate()
115
 
116
   Options Database Key:
117
.    -qep_monitor_cancel - Cancels all monitors that have been hardwired
118
      into a code by calls to QEPMonitorSet(),
119
      but does not cancel those set via the options database.
120
 
121
   Level: intermediate
122
 
123
.seealso: QEPMonitorSet()
124
@*/
125
PetscErrorCode QEPMonitorCancel(QEP qep)
126
{
127
  PetscErrorCode ierr;
128
  PetscInt       i;
129
 
130
  PetscFunctionBegin;
2213 jroman 131
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 132
  for (i=0; i<qep->numbermonitors; i++) {
133
    if (qep->monitordestroy[i]) {
2311 jroman 134
      ierr = (*qep->monitordestroy[i])(&qep->monitorcontext[i]);CHKERRQ(ierr);
1887 jroman 135
    }
136
  }
137
  qep->numbermonitors = 0;
138
  PetscFunctionReturn(0);
139
}
140
 
141
#undef __FUNCT__  
142
#define __FUNCT__ "QEPGetMonitorContext"
143
/*@C
144
   QEPGetMonitorContext - Gets the monitor context, as set by
2242 jroman 145
   QEPMonitorSet() for the FIRST monitor only.
1887 jroman 146
 
147
   Not Collective
148
 
149
   Input Parameter:
150
.  qep - eigensolver context obtained from QEPCreate()
151
 
152
   Output Parameter:
153
.  ctx - monitor context
154
 
155
   Level: intermediate
156
 
2242 jroman 157
.seealso: QEPMonitorSet(), QEPDefaultMonitor()
1887 jroman 158
@*/
2331 jroman 159
PetscErrorCode QEPGetMonitorContext(QEP qep,void **ctx)
1887 jroman 160
{
161
  PetscFunctionBegin;
2213 jroman 162
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 163
  *ctx =      (qep->monitorcontext[0]);
164
  PetscFunctionReturn(0);
165
}
166
 
167
#undef __FUNCT__  
2054 eromero 168
#define __FUNCT__ "QEPMonitorAll"
1887 jroman 169
/*@C
2054 eromero 170
   QEPMonitorAll - Print the current approximate values and
1887 jroman 171
   error estimates at each iteration of the quadratic eigensolver.
172
 
173
   Collective on QEP
174
 
175
   Input Parameters:
176
+  qep    - quadratic eigensolver context
177
.  its    - iteration number
178
.  nconv  - number of converged eigenpairs so far
179
.  eigr   - real part of the eigenvalues
180
.  eigi   - imaginary part of the eigenvalues
181
.  errest - error estimates
182
.  nest   - number of error estimates to display
183
-  dummy  - unused monitor context
184
 
185
   Level: intermediate
186
 
2328 jroman 187
.seealso: QEPMonitorSet(), QEPMonitorFirst(), QEPMonitorConverged()
1887 jroman 188
@*/
2054 eromero 189
PetscErrorCode QEPMonitorAll(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *dummy)
1887 jroman 190
{
2419 jroman 191
  PetscErrorCode ierr;
192
  PetscInt       i;
193
  PetscViewer    viewer = dummy? (PetscViewer)dummy: PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);
1887 jroman 194
 
195
  PetscFunctionBegin;
196
  if (its) {
2419 jroman 197
    ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)qep)->tablevel);CHKERRQ(ierr);
198
    ierr = PetscViewerASCIIPrintf(viewer,"%3D QEP nconv=%D Values (Errors)",its,nconv);CHKERRQ(ierr);
1887 jroman 199
    for (i=0;i<nest;i++) {
200
#if defined(PETSC_USE_COMPLEX)
2419 jroman 201
      ierr = PetscViewerASCIIPrintf(viewer," %G%+Gi",PetscRealPart(eigr[i]),PetscImaginaryPart(eigr[i]));CHKERRQ(ierr);
1887 jroman 202
#else
2419 jroman 203
      ierr = PetscViewerASCIIPrintf(viewer," %G",eigr[i]);CHKERRQ(ierr);
204
      if (eigi[i]!=0.0) { ierr = PetscViewerASCIIPrintf(viewer,"%+Gi",eigi[i]);CHKERRQ(ierr); }
1887 jroman 205
#endif
2419 jroman 206
      ierr = PetscViewerASCIIPrintf(viewer," (%10.8e)",(double)errest[i]);CHKERRQ(ierr);
1887 jroman 207
    }
2419 jroman 208
    ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
209
    ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)qep)->tablevel);CHKERRQ(ierr);
1887 jroman 210
  }
211
  PetscFunctionReturn(0);
212
}
213
 
214
#undef __FUNCT__  
215
#define __FUNCT__ "QEPMonitorFirst"
216
/*@C
217
   QEPMonitorFirst - Print the first unconverged approximate value and
218
   error estimate at each iteration of the quadratic eigensolver.
219
 
220
   Collective on QEP
221
 
222
   Input Parameters:
223
+  qep    - quadratic eigensolver context
224
.  its    - iteration number
225
.  nconv  - number of converged eigenpairs so far
226
.  eigr   - real part of the eigenvalues
227
.  eigi   - imaginary part of the eigenvalues
228
.  errest - error estimates
229
.  nest   - number of error estimates to display
230
-  dummy  - unused monitor context
231
 
232
   Level: intermediate
233
 
2328 jroman 234
.seealso: QEPMonitorSet(), QEPMonitorAll(), QEPMonitorConverged()
1887 jroman 235
@*/
236
PetscErrorCode QEPMonitorFirst(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *dummy)
237
{
2419 jroman 238
  PetscErrorCode ierr;
239
  PetscViewer    viewer = dummy? (PetscViewer)dummy: PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);
1887 jroman 240
 
241
  PetscFunctionBegin;
242
  if (its && nconv<nest) {
2419 jroman 243
    ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)qep)->tablevel);CHKERRQ(ierr);
244
    ierr = PetscViewerASCIIPrintf(viewer,"%3D QEP nconv=%D first unconverged value (error)",its,nconv);CHKERRQ(ierr);
1887 jroman 245
#if defined(PETSC_USE_COMPLEX)
2419 jroman 246
    ierr = PetscViewerASCIIPrintf(viewer," %G%+Gi",PetscRealPart(eigr[nconv]),PetscImaginaryPart(eigr[nconv]));CHKERRQ(ierr);
1887 jroman 247
#else
2419 jroman 248
    ierr = PetscViewerASCIIPrintf(viewer," %G",eigr[nconv]);CHKERRQ(ierr);
249
    if (eigi[nconv]!=0.0) { ierr = PetscViewerASCIIPrintf(viewer,"%+Gi",eigi[nconv]);CHKERRQ(ierr); }
1887 jroman 250
#endif
2419 jroman 251
    ierr = PetscViewerASCIIPrintf(viewer," (%10.8e)\n",(double)errest[nconv]);CHKERRQ(ierr);
252
    ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)qep)->tablevel);CHKERRQ(ierr);
1887 jroman 253
  }
254
  PetscFunctionReturn(0);
255
}
256
 
257
#undef __FUNCT__  
258
#define __FUNCT__ "QEPMonitorConverged"
259
/*@C
260
   QEPMonitorConverged - Print the approximate values and
261
   error estimates as they converge.
262
 
263
   Collective on QEP
264
 
265
   Input Parameters:
266
+  qep    - quadratic eigensolver context
267
.  its    - iteration number
268
.  nconv  - number of converged eigenpairs so far
269
.  eigr   - real part of the eigenvalues
270
.  eigi   - imaginary part of the eigenvalues
271
.  errest - error estimates
272
.  nest   - number of error estimates to display
273
-  dummy  - unused monitor context
274
 
275
   Level: intermediate
276
 
2328 jroman 277
.seealso: QEPMonitorSet(), QEPMonitorFirst(), QEPMonitorAll()
1887 jroman 278
@*/
279
PetscErrorCode QEPMonitorConverged(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *dummy)
280
{
2311 jroman 281
  PetscErrorCode   ierr;
282
  PetscInt         i;
2419 jroman 283
  SlepcConvMonitor ctx = (SlepcConvMonitor)dummy;
1887 jroman 284
 
285
  PetscFunctionBegin;
286
  if (!its) {
1951 jroman 287
    ctx->oldnconv = 0;
1887 jroman 288
  } else {
1951 jroman 289
    for (i=ctx->oldnconv;i<nconv;i++) {
2419 jroman 290
      ierr = PetscViewerASCIIAddTab(ctx->viewer,((PetscObject)qep)->tablevel);CHKERRQ(ierr);
291
      ierr = PetscViewerASCIIPrintf(ctx->viewer,"%3D QEP converged value (error) #%D",its,i);CHKERRQ(ierr);
1887 jroman 292
#if defined(PETSC_USE_COMPLEX)
2419 jroman 293
      ierr = PetscViewerASCIIPrintf(ctx->viewer," %G%+Gi",PetscRealPart(eigr[i]),PetscImaginaryPart(eigr[i]));CHKERRQ(ierr);
1887 jroman 294
#else
2419 jroman 295
      ierr = PetscViewerASCIIPrintf(ctx->viewer," %G",eigr[i]);CHKERRQ(ierr);
296
      if (eigi[i]!=0.0) { ierr = PetscViewerASCIIPrintf(ctx->viewer,"%+Gi",eigi[i]);CHKERRQ(ierr); }
1887 jroman 297
#endif
2419 jroman 298
      ierr = PetscViewerASCIIPrintf(ctx->viewer," (%10.8e)\n",(double)errest[i]);CHKERRQ(ierr);
299
      ierr = PetscViewerASCIISubtractTab(ctx->viewer,((PetscObject)qep)->tablevel);CHKERRQ(ierr);
1887 jroman 300
    }
1951 jroman 301
    ctx->oldnconv = nconv;
1887 jroman 302
  }
303
  PetscFunctionReturn(0);
304
}
305
 
2311 jroman 306
#undef __FUNCT__  
1887 jroman 307
#define __FUNCT__ "QEPMonitorLG"
308
PetscErrorCode QEPMonitorLG(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
309
{
310
  PetscViewer    viewer = (PetscViewer) monctx;
311
  PetscDraw      draw;
312
  PetscDrawLG    lg;
313
  PetscErrorCode ierr;
2054 eromero 314
  PetscReal      x,y;
315
 
316
  PetscFunctionBegin;
317
  if (!viewer) { viewer = PETSC_VIEWER_DRAW_(((PetscObject)qep)->comm); }
318
  ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
319
  ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
320
  if (!its) {
321
    ierr = PetscDrawSetTitle(draw,"Error estimates");CHKERRQ(ierr);
322
    ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr);
323
    ierr = PetscDrawLGSetDimension(lg,1);CHKERRQ(ierr);
324
    ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
325
    ierr = PetscDrawLGSetLimits(lg,0,1.0,log10(qep->tol)-2,0.0);CHKERRQ(ierr);
326
  }
327
 
328
  x = (PetscReal) its;
329
  if (errest[nconv] > 0.0) y = log10(errest[nconv]); else y = 0.0;
330
  ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr);
331
 
332
  ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
333
  PetscFunctionReturn(0);
334
}
335
 
336
#undef __FUNCT__  
337
#define __FUNCT__ "QEPMonitorLGAll"
338
PetscErrorCode QEPMonitorLGAll(QEP qep,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *monctx)
339
{
340
  PetscViewer    viewer = (PetscViewer) monctx;
341
  PetscDraw      draw;
342
  PetscDrawLG    lg;
343
  PetscErrorCode ierr;
1887 jroman 344
  PetscReal      *x,*y;
2317 jroman 345
  PetscInt       i,n = PetscMin(qep->nev,255);
1887 jroman 346
 
347
  PetscFunctionBegin;
348
  if (!viewer) { viewer = PETSC_VIEWER_DRAW_(((PetscObject)qep)->comm); }
349
  ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
350
  ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr);
351
  if (!its) {
352
    ierr = PetscDrawSetTitle(draw,"Error estimates");CHKERRQ(ierr);
353
    ierr = PetscDrawSetDoubleBuffer(draw);CHKERRQ(ierr);
354
    ierr = PetscDrawLGSetDimension(lg,n);CHKERRQ(ierr);
355
    ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);
356
    ierr = PetscDrawLGSetLimits(lg,0,1.0,log10(qep->tol)-2,0.0);CHKERRQ(ierr);
357
  }
358
 
359
  ierr = PetscMalloc(sizeof(PetscReal)*n,&x);CHKERRQ(ierr);
360
  ierr = PetscMalloc(sizeof(PetscReal)*n,&y);CHKERRQ(ierr);
361
  for (i=0;i<n;i++) {
362
    x[i] = (PetscReal) its;
2054 eromero 363
    if (i < nest && errest[i] > 0.0) y[i] = log10(errest[i]);
364
    else y[i] = 0.0;
1887 jroman 365
  }
366
  ierr = PetscDrawLGAddPoint(lg,x,y);CHKERRQ(ierr);
367
 
368
  ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr);
369
  ierr = PetscFree(x);CHKERRQ(ierr);
370
  ierr = PetscFree(y);CHKERRQ(ierr);  
371
  PetscFunctionReturn(0);
372
}
373