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 options that can be set via the command-line
3
      or procedurally.
4
 
5
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2116 eromero 7
   Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
1887 jroman 8
 
9
   This file is part of SLEPc.
10
 
11
   SLEPc is free software: you can redistribute it and/or modify it under  the
12
   terms of version 3 of the GNU Lesser General Public License as published by
13
   the Free Software Foundation.
14
 
15
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
16
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
17
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
18
   more details.
19
 
20
   You  should have received a copy of the GNU Lesser General  Public  License
21
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23
*/
24
 
25
#include "private/qepimpl.h"   /*I "slepcqep.h" I*/
26
 
27
#undef __FUNCT__  
28
#define __FUNCT__ "QEPSetFromOptions"
29
/*@
30
   QEPSetFromOptions - Sets QEP options from the options database.
31
   This routine must be called before QEPSetUp() if the user is to be
32
   allowed to set the solver type.
33
 
34
   Collective on QEP
35
 
36
   Input Parameters:
37
.  qep - the quadratic eigensolver context
38
 
39
   Notes:  
40
   To see all options, run your program with the -help option.
41
 
42
   Level: beginner
43
@*/
44
PetscErrorCode QEPSetFromOptions(QEP qep)
45
{
46
  PetscErrorCode ierr;
47
  char           type[256],monfilename[PETSC_MAX_PATH_LEN];
2216 jroman 48
  PetscBool      flg,val;
1887 jroman 49
  PetscReal      r;
50
  PetscInt       i,j,k;
51
  PetscViewerASCIIMonitor monviewer;
1951 jroman 52
  QEPMONITOR_CONV *ctx;
1887 jroman 53
 
54
  PetscFunctionBegin;
2213 jroman 55
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 56
  ierr = PetscOptionsBegin(((PetscObject)qep)->comm,((PetscObject)qep)->prefix,"Quadratic Eigenvalue Problem (QEP) Solver Options","QEP");CHKERRQ(ierr);
57
    ierr = PetscOptionsList("-qep_type","Quadratic Eigenvalue Problem method","QEPSetType",QEPList,(char*)(((PetscObject)qep)->type_name?((PetscObject)qep)->type_name:QEPLINEAR),type,256,&flg);CHKERRQ(ierr);
58
    if (flg) {
59
      ierr = QEPSetType(qep,type);CHKERRQ(ierr);
60
    } else if (!((PetscObject)qep)->type_name) {
61
      ierr = QEPSetType(qep,QEPLINEAR);CHKERRQ(ierr);
62
    }
63
 
2216 jroman 64
    ierr = PetscOptionsBoolGroupBegin("-qep_general","general quadratic eigenvalue problem","QEPSetProblemType",&flg);CHKERRQ(ierr);
1907 jroman 65
    if (flg) {ierr = QEPSetProblemType(qep,QEP_GENERAL);CHKERRQ(ierr);}
2216 jroman 66
    ierr = PetscOptionsBoolGroup("-qep_hermitian","hermitian quadratic eigenvalue problem","QEPSetProblemType",&flg);CHKERRQ(ierr);
1907 jroman 67
    if (flg) {ierr = QEPSetProblemType(qep,QEP_HERMITIAN);CHKERRQ(ierr);}
2216 jroman 68
    ierr = PetscOptionsBoolGroupEnd("-qep_gyroscopic","gyroscopic quadratic eigenvalue problem","QEPSetProblemType",&flg);CHKERRQ(ierr);
1907 jroman 69
    if (flg) {ierr = QEPSetProblemType(qep,QEP_GYROSCOPIC);CHKERRQ(ierr);}
70
 
1918 jroman 71
    r = PETSC_IGNORE;
72
    ierr = PetscOptionsReal("-qep_scale","Scale factor","QEPSetScaleFactor",qep->sfactor,&r,PETSC_NULL);CHKERRQ(ierr);
73
    ierr = QEPSetScaleFactor(qep,r);CHKERRQ(ierr);
74
 
1887 jroman 75
    r = i = PETSC_IGNORE;
76
    ierr = PetscOptionsInt("-qep_max_it","Maximum number of iterations","QEPSetTolerances",qep->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
77
    ierr = PetscOptionsReal("-qep_tol","Tolerance","QEPSetTolerances",qep->tol,&r,PETSC_NULL);CHKERRQ(ierr);
78
    ierr = QEPSetTolerances(qep,r,i);CHKERRQ(ierr);
2216 jroman 79
    ierr = PetscOptionsBoolGroupBegin("-qep_convergence_default","Default (relative error) convergence test","QEPSetConvergenceTest",&flg);CHKERRQ(ierr);
1887 jroman 80
    if (flg) {ierr = QEPSetConvergenceTest(qep,QEPDefaultConverged,PETSC_NULL);CHKERRQ(ierr);}
2216 jroman 81
    ierr = PetscOptionsBoolGroupEnd("-qep_convergence_absolute","Absolute error convergence test","QEPSetConvergenceTest",&flg);CHKERRQ(ierr);
1887 jroman 82
    if (flg) {ierr = QEPSetConvergenceTest(qep,QEPAbsoluteConverged,PETSC_NULL);CHKERRQ(ierr);}
83
 
84
    i = j = k = PETSC_IGNORE;
85
    ierr = PetscOptionsInt("-qep_nev","Number of eigenvalues to compute","QEPSetDimensions",qep->nev,&i,PETSC_NULL);CHKERRQ(ierr);
86
    ierr = PetscOptionsInt("-qep_ncv","Number of basis vectors","QEPSetDimensions",qep->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
87
    ierr = PetscOptionsInt("-qep_mpd","Maximum dimension of projected problem","QEPSetDimensions",qep->mpd,&k,PETSC_NULL);CHKERRQ(ierr);
88
    ierr = QEPSetDimensions(qep,i,j,k);CHKERRQ(ierr);
89
 
90
    /* -----------------------------------------------------------------------*/
91
    /*
92
      Cancels all monitors hardwired into code before call to QEPSetFromOptions()
93
    */
94
    flg  = PETSC_FALSE;
2216 jroman 95
    ierr = PetscOptionsBool("-qep_monitor_cancel","Remove any hardwired monitor routines","QEPMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1887 jroman 96
    if (flg) {
97
      ierr = QEPMonitorCancel(qep); CHKERRQ(ierr);
98
    }
99
    /*
100
      Prints approximate eigenvalues and error estimates at each iteration
101
    */
2054 eromero 102
    ierr = PetscOptionsString("-qep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","QEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1887 jroman 103
    if (flg) {
104
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)qep)->comm,monfilename,((PetscObject)qep)->tablevel,&monviewer);CHKERRQ(ierr);
2054 eromero 105
      ierr = QEPMonitorSet(qep,QEPMonitorFirst,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
1887 jroman 106
    }
107
    ierr = PetscOptionsString("-qep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","QEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
108
    if (flg) {
1951 jroman 109
      ierr = PetscNew(QEPMONITOR_CONV,&ctx);CHKERRQ(ierr);
110
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)qep)->comm,monfilename,((PetscObject)qep)->tablevel,&ctx->viewer);CHKERRQ(ierr);
111
      ierr = QEPMonitorSet(qep,QEPMonitorConverged,ctx,(PetscErrorCode (*)(void*))QEPMonitorDestroy_Converged);CHKERRQ(ierr);
1887 jroman 112
    }
2054 eromero 113
    ierr = PetscOptionsString("-qep_monitor_all","Monitor approximate eigenvalues and error estimates","QEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1887 jroman 114
    if (flg) {
115
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)qep)->comm,monfilename,((PetscObject)qep)->tablevel,&monviewer);CHKERRQ(ierr);
2054 eromero 116
      ierr = QEPMonitorSet(qep,QEPMonitorAll,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
117
      ierr = QEPSetTrackAll(qep,PETSC_TRUE);CHKERRQ(ierr);
1887 jroman 118
    }
119
    flg = PETSC_FALSE;
2216 jroman 120
    ierr = PetscOptionsBool("-qep_monitor_draw","Monitor first unconverged approximate error estimate graphically","QEPMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1887 jroman 121
    if (flg) {
122
      ierr = QEPMonitorSet(qep,QEPMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
123
    }
2054 eromero 124
    flg = PETSC_FALSE;
2216 jroman 125
    ierr = PetscOptionsBool("-qep_monitor_draw_all","Monitor error estimates graphically","QEPMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
2054 eromero 126
    if (flg) {
127
      ierr = QEPMonitorSet(qep,QEPMonitorLGAll,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
128
      ierr = QEPSetTrackAll(qep,PETSC_TRUE);CHKERRQ(ierr);
129
    }
1887 jroman 130
  /* -----------------------------------------------------------------------*/
131
 
2216 jroman 132
    ierr = PetscOptionsBoolGroupBegin("-qep_largest_magnitude","compute largest eigenvalues in magnitude","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
1887 jroman 133
    if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_LARGEST_MAGNITUDE);CHKERRQ(ierr);}
2216 jroman 134
    ierr = PetscOptionsBoolGroup("-qep_smallest_magnitude","compute smallest eigenvalues in magnitude","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
1887 jroman 135
    if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_SMALLEST_MAGNITUDE);CHKERRQ(ierr);}
2216 jroman 136
    ierr = PetscOptionsBoolGroup("-qep_largest_real","compute largest real parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
1887 jroman 137
    if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_LARGEST_REAL);CHKERRQ(ierr);}
2216 jroman 138
    ierr = PetscOptionsBoolGroup("-qep_smallest_real","compute smallest real parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
1887 jroman 139
    if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_SMALLEST_REAL);CHKERRQ(ierr);}
2216 jroman 140
    ierr = PetscOptionsBoolGroup("-qep_largest_imaginary","compute largest imaginary parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
1887 jroman 141
    if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_LARGEST_IMAGINARY);CHKERRQ(ierr);}
2216 jroman 142
    ierr = PetscOptionsBoolGroupEnd("-qep_smallest_imaginary","compute smallest imaginary parts","QEPSetWhichEigenpairs",&flg);CHKERRQ(ierr);
1887 jroman 143
    if (flg) {ierr = QEPSetWhichEigenpairs(qep,QEP_SMALLEST_IMAGINARY);CHKERRQ(ierr);}
144
 
2216 jroman 145
    ierr = PetscOptionsBool("-qep_left_vectors","Compute left eigenvectors also","QEPSetLeftVectorsWanted",qep->leftvecs,&val,&flg);CHKERRQ(ierr);
1944 jroman 146
    if (flg) {
147
      ierr = QEPSetLeftVectorsWanted(qep,val);CHKERRQ(ierr);
148
    }
149
 
1887 jroman 150
    ierr = PetscOptionsName("-qep_view","Print detailed information on solver used","QEPView",0);CHKERRQ(ierr);
151
    ierr = PetscOptionsName("-qep_view_binary","Save the matrices associated to the eigenproblem","QEPSetFromOptions",0);CHKERRQ(ierr);
152
    ierr = PetscOptionsName("-qep_plot_eigs","Make a plot of the computed eigenvalues","QEPSolve",0);CHKERRQ(ierr);
153
 
154
    if (qep->ops->setfromoptions) {
155
      ierr = (*qep->ops->setfromoptions)(qep);CHKERRQ(ierr);
156
    }
157
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
158
 
159
  ierr = IPSetFromOptions(qep->ip); CHKERRQ(ierr);
160
  PetscFunctionReturn(0);
161
}
162
 
163
#undef __FUNCT__  
164
#define __FUNCT__ "QEPGetTolerances"
165
/*@
166
   QEPGetTolerances - Gets the tolerance and maximum iteration count used
167
   by the QEP convergence tests.
168
 
169
   Not Collective
170
 
171
   Input Parameter:
172
.  qep - the quadratic eigensolver context
173
 
174
   Output Parameters:
175
+  tol - the convergence tolerance
176
-  maxits - maximum number of iterations
177
 
178
   Notes:
179
   The user can specify PETSC_NULL for any parameter that is not needed.
180
 
181
   Level: intermediate
182
 
183
.seealso: QEPSetTolerances()
184
@*/
185
PetscErrorCode QEPGetTolerances(QEP qep,PetscReal *tol,PetscInt *maxits)
186
{
187
  PetscFunctionBegin;
2213 jroman 188
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 189
  if (tol)    *tol    = qep->tol;
190
  if (maxits) *maxits = qep->max_it;
191
  PetscFunctionReturn(0);
192
}
193
 
194
#undef __FUNCT__  
1918 jroman 195
#define __FUNCT__ "QEPSetTolerances"
1887 jroman 196
/*@
197
   QEPSetTolerances - Sets the tolerance and maximum iteration count used
198
   by the QEP convergence tests.
199
 
200
   Collective on QEP
201
 
202
   Input Parameters:
203
+  qep - the quadratic eigensolver context
204
.  tol - the convergence tolerance
205
-  maxits - maximum number of iterations to use
206
 
207
   Options Database Keys:
208
+  -qep_tol <tol> - Sets the convergence tolerance
209
-  -qep_max_it <maxits> - Sets the maximum number of iterations allowed
210
 
211
   Notes:
212
   Use PETSC_IGNORE for an argument that need not be changed.
213
 
214
   Use PETSC_DECIDE for maxits to assign a reasonably good value, which is
215
   dependent on the solution method.
216
 
217
   Level: intermediate
218
 
219
.seealso: QEPGetTolerances()
220
@*/
221
PetscErrorCode QEPSetTolerances(QEP qep,PetscReal tol,PetscInt maxits)
222
{
223
  PetscFunctionBegin;
2213 jroman 224
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 225
  if (tol != PETSC_IGNORE) {
226
    if (tol == PETSC_DEFAULT) {
227
      qep->tol = 1e-7;
228
    } else {
2214 jroman 229
      if (tol < 0.0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1887 jroman 230
      qep->tol = tol;
231
    }
232
  }
233
  if (maxits != PETSC_IGNORE) {
234
    if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
235
      qep->max_it = 0;
236
      qep->setupcalled = 0;
237
    } else {
2214 jroman 238
      if (maxits < 0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
1887 jroman 239
      qep->max_it = maxits;
240
    }
241
  }
242
  PetscFunctionReturn(0);
243
}
244
 
245
#undef __FUNCT__  
246
#define __FUNCT__ "QEPGetDimensions"
247
/*@
248
   QEPGetDimensions - Gets the number of eigenvalues to compute
249
   and the dimension of the subspace.
250
 
251
   Not Collective
252
 
253
   Input Parameter:
254
.  qep - the quadratic eigensolver context
255
 
256
   Output Parameters:
257
+  nev - number of eigenvalues to compute
258
.  ncv - the maximum dimension of the subspace to be used by the solver
259
-  mpd - the maximum dimension allowed for the projected problem
260
 
261
   Notes:
262
   The user can specify PETSC_NULL for any parameter that is not needed.
263
 
264
   Level: intermediate
265
 
266
.seealso: QEPSetDimensions()
267
@*/
268
PetscErrorCode QEPGetDimensions(QEP qep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
269
{
270
  PetscFunctionBegin;
2213 jroman 271
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 272
  if (nev) *nev = qep->nev;
273
  if (ncv) *ncv = qep->ncv;
274
  if (mpd) *mpd = qep->mpd;
275
  PetscFunctionReturn(0);
276
}
277
 
278
#undef __FUNCT__  
279
#define __FUNCT__ "QEPSetDimensions"
280
/*@
281
   QEPSetDimensions - Sets the number of eigenvalues to compute
282
   and the dimension of the subspace.
283
 
284
   Collective on QEP
285
 
286
   Input Parameters:
287
+  qep - the quadratic eigensolver context
288
.  nev - number of eigenvalues to compute
289
.  ncv - the maximum dimension of the subspace to be used by the solver
290
-  mpd - the maximum dimension allowed for the projected problem
291
 
292
   Options Database Keys:
293
+  -qep_nev <nev> - Sets the number of eigenvalues
294
.  -qep_ncv <ncv> - Sets the dimension of the subspace
295
-  -qep_mpd <mpd> - Sets the maximum projected dimension
296
 
297
   Notes:
298
   Use PETSC_IGNORE to retain the previous value of any parameter.
299
 
300
   Use PETSC_DECIDE for ncv and mpd to assign a reasonably good value, which is
301
   dependent on the solution method.
302
 
303
   The parameters ncv and mpd are intimately related, so that the user is advised
304
   to set one of them at most. Normal usage is the following
305
+  - In cases where nev is small, the user sets ncv (a reasonable default is 2*nev).
306
-  - In cases where nev is large, the user sets mpd.
307
 
308
   The value of ncv should always be between nev and (nev+mpd), typically
309
   ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
310
   a smaller value should be used.
311
 
312
   Level: intermediate
313
 
314
.seealso: QEPGetDimensions()
315
@*/
316
PetscErrorCode QEPSetDimensions(QEP qep,PetscInt nev,PetscInt ncv,PetscInt mpd)
317
{
318
  PetscFunctionBegin;
2213 jroman 319
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 320
 
321
  if( nev != PETSC_IGNORE ) {
2214 jroman 322
    if (nev<1) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
1887 jroman 323
    qep->nev = nev;
324
    qep->setupcalled = 0;
325
  }
326
  if( ncv != PETSC_IGNORE ) {
327
    if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
328
      qep->ncv = 0;
329
    } else {
2214 jroman 330
      if (ncv<1) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
1887 jroman 331
      qep->ncv = ncv;
332
    }
333
    qep->setupcalled = 0;
334
  }
335
  if( mpd != PETSC_IGNORE ) {
336
    if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
337
      qep->mpd = 0;
338
    } else {
2214 jroman 339
      if (mpd<1) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
1887 jroman 340
      qep->mpd = mpd;
341
    }
342
  }
343
  PetscFunctionReturn(0);
344
}
345
 
346
#undef __FUNCT__  
347
#define __FUNCT__ "QEPSetWhichEigenpairs"
348
/*@
349
    QEPSetWhichEigenpairs - Specifies which portion of the spectrum is
350
    to be sought.
351
 
352
    Collective on QEP
353
 
354
    Input Parameters:
1942 jroman 355
+   qep   - eigensolver context obtained from QEPCreate()
1887 jroman 356
-   which - the portion of the spectrum to be sought
357
 
358
    Possible values:
359
    The parameter 'which' can have one of these values
360
 
361
+     QEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
362
.     QEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
363
.     QEP_LARGEST_REAL - largest real parts
364
.     QEP_SMALLEST_REAL - smallest real parts
365
.     QEP_LARGEST_IMAGINARY - largest imaginary parts
366
-     QEP_SMALLEST_IMAGINARY - smallest imaginary parts
367
 
368
    Options Database Keys:
369
+   -qep_largest_magnitude - Sets largest eigenvalues in magnitude
370
.   -qep_smallest_magnitude - Sets smallest eigenvalues in magnitude
371
.   -qep_largest_real - Sets largest real parts
372
.   -qep_smallest_real - Sets smallest real parts
373
.   -qep_largest_imaginary - Sets largest imaginary parts
374
-   -qep_smallest_imaginary - Sets smallest imaginary parts
375
 
376
    Notes:
377
    Not all eigensolvers implemented in QEP account for all the possible values
378
    stated above. If SLEPc is compiled for real numbers QEP_LARGEST_IMAGINARY
379
    and QEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
380
    for eigenvalue selection.
381
 
382
    Level: intermediate
383
 
384
.seealso: QEPGetWhichEigenpairs(), QEPWhich
385
@*/
386
PetscErrorCode QEPSetWhichEigenpairs(QEP qep,QEPWhich which)
387
{
388
  PetscFunctionBegin;
2213 jroman 389
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1942 jroman 390
  if (which!=PETSC_IGNORE) {
391
    if (which==PETSC_DECIDE || which==PETSC_DEFAULT) qep->which = (QEPWhich)0;
392
    else switch (which) {
393
      case QEP_LARGEST_MAGNITUDE:
394
      case QEP_SMALLEST_MAGNITUDE:
395
      case QEP_LARGEST_REAL:
396
      case QEP_SMALLEST_REAL:
397
      case QEP_LARGEST_IMAGINARY:
398
      case QEP_SMALLEST_IMAGINARY:
399
        if (qep->which != which) {
400
          qep->setupcalled = 0;
401
          qep->which = which;
402
        }
403
        break;
404
      default:
2214 jroman 405
        SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
1942 jroman 406
    }
1887 jroman 407
  }
408
  PetscFunctionReturn(0);
409
}
410
 
411
#undef __FUNCT__  
412
#define __FUNCT__ "QEPGetWhichEigenpairs"
413
/*@C
414
    QEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
415
    sought.
416
 
417
    Not Collective
418
 
419
    Input Parameter:
420
.   qep - eigensolver context obtained from QEPCreate()
421
 
422
    Output Parameter:
423
.   which - the portion of the spectrum to be sought
424
 
425
    Notes:
426
    See QEPSetWhichEigenpairs() for possible values of 'which'.
427
 
428
    Level: intermediate
429
 
430
.seealso: QEPSetWhichEigenpairs(), QEPWhich
431
@*/
432
PetscErrorCode QEPGetWhichEigenpairs(QEP qep,QEPWhich *which)
433
{
434
  PetscFunctionBegin;
2213 jroman 435
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 436
  PetscValidPointer(which,2);
437
  *which = qep->which;
438
  PetscFunctionReturn(0);
439
}
440
 
441
#undef __FUNCT__  
1944 jroman 442
#define __FUNCT__ "QEPSetLeftVectorsWanted"
443
/*@
444
    QEPSetLeftVectorsWanted - Specifies which eigenvectors are required.
445
 
446
    Collective on QEP
447
 
448
    Input Parameters:
449
+   qep      - the quadratic eigensolver context
450
-   leftvecs - whether left eigenvectors are required or not
451
 
452
    Options Database Keys:
453
.   -qep_left_vectors <boolean> - Sets/resets the boolean flag 'leftvecs'
454
 
455
    Notes:
456
    If the user sets leftvecs=PETSC_TRUE then the solver uses a variant of
457
    the algorithm that computes both right and left eigenvectors. This is
458
    usually much more costly. This option is not available in all solvers.
459
 
460
    Level: intermediate
461
 
462
.seealso: QEPGetLeftVectorsWanted(), QEPGetEigenvectorLeft()
463
@*/
2216 jroman 464
PetscErrorCode QEPSetLeftVectorsWanted(QEP qep,PetscBool leftvecs)
1944 jroman 465
{
466
  PetscFunctionBegin;
2213 jroman 467
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1949 jroman 468
  if (qep->leftvecs != leftvecs) {
469
    qep->leftvecs = leftvecs;
470
    qep->setupcalled = 0;
471
  }
1944 jroman 472
  PetscFunctionReturn(0);
473
}
474
 
475
#undef __FUNCT__  
476
#define __FUNCT__ "QEPGetLeftVectorsWanted"
477
/*@C
478
    QEPGetLeftVectorsWanted - Returns the flag indicating whether left
479
    eigenvectors are required or not.
480
 
481
    Not Collective
482
 
483
    Input Parameter:
484
.   qep - the eigensolver context
485
 
486
    Output Parameter:
487
.   leftvecs - the returned flag
488
 
489
    Level: intermediate
490
 
491
.seealso: QEPSetLeftVectorsWanted()
492
@*/
2216 jroman 493
PetscErrorCode QEPGetLeftVectorsWanted(QEP qep,PetscBool *leftvecs)
1944 jroman 494
{
495
  PetscFunctionBegin;
2213 jroman 496
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1944 jroman 497
  PetscValidPointer(leftvecs,2);
498
  *leftvecs = qep->leftvecs;
499
  PetscFunctionReturn(0);
500
}
501
 
502
#undef __FUNCT__  
1918 jroman 503
#define __FUNCT__ "QEPGetScaleFactor"
504
/*@
505
   QEPGetScaleFactor - Gets the factor used for scaling the quadratic eigenproblem.
506
 
507
   Not Collective
508
 
509
   Input Parameter:
510
.  qep - the quadratic eigensolver context
511
 
512
   Output Parameters:
513
.  alpha - the scaling factor
514
 
515
   Notes:
516
   If the user did not specify a scaling factor, then after QEPSolve() the
517
   default value is returned.
518
 
519
   Level: intermediate
520
 
521
.seealso: QEPSetScaleFactor(), QEPSolve()
522
@*/
523
PetscErrorCode QEPGetScaleFactor(QEP qep,PetscReal *alpha)
524
{
525
  PetscFunctionBegin;
2213 jroman 526
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2114 eromero 527
  if (alpha) *alpha = qep->sfactor;
1918 jroman 528
  PetscFunctionReturn(0);
529
}
530
 
531
#undef __FUNCT__  
532
#define __FUNCT__ "QEPSetScaleFactor"
533
/*@
534
   QEPSetScaleFactor - Sets the scaling factor to be used for scaling the
535
   quadratic problem before attempting to solve.
536
 
537
   Collective on QEP
538
 
539
   Input Parameters:
540
+  qep   - the quadratic eigensolver context
541
-  alpha - the scaling factor
542
 
543
   Options Database Keys:
544
.  -qep_scale <alpha> - Sets the scaling factor
545
 
546
   Notes:
547
   For the problem (l^2*M + l*C + K)*x = 0, the effect of scaling is to work
548
   with matrices (alpha^2*M, alpha*C, K), then scale the computed eigenvalue.
549
 
550
   The default is to scale with alpha = norm(K)/norm(M).
551
 
552
   Level: intermediate
553
 
554
.seealso: QEPGetScaleFactor()
555
@*/
556
PetscErrorCode QEPSetScaleFactor(QEP qep,PetscReal alpha)
557
{
558
  PetscFunctionBegin;
2213 jroman 559
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1918 jroman 560
  if (alpha != PETSC_IGNORE) {
561
    if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
562
      qep->sfactor = 0.0;
563
    } else {
2214 jroman 564
      if (alpha < 0.0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1918 jroman 565
      qep->sfactor = alpha;
566
    }
567
  }
568
  PetscFunctionReturn(0);
569
}
570
 
571
#undef __FUNCT__  
1906 jroman 572
#define __FUNCT__ "QEPSetProblemType"
573
/*@
574
   QEPSetProblemType - Specifies the type of the quadratic eigenvalue problem.
575
 
576
   Collective on QEP
577
 
578
   Input Parameters:
579
+  qep      - the quadratic eigensolver context
580
-  type     - a known type of quadratic eigenvalue problem
581
 
582
   Options Database Keys:
583
+  -qep_general - general problem with no particular structure
584
.  -qep_hermitian - problem whose coefficient matrices are Hermitian
585
-  -qep_gyroscopic - problem with Hamiltonian structure
586
 
587
   Notes:  
588
   Allowed values for the problem type are: general (QEP_GENERAL), Hermitian
589
   (QEP_HERMITIAN), and gyroscopic (QEP_GYROSCOPIC).
590
 
591
   This function is used to instruct SLEPc to exploit certain structure in
592
   the quadratic eigenproblem. By default, no particular structure is assumed.
593
 
594
   If the problem matrices are Hermitian (symmetric in the real case) or
595
   Hermitian/skew-Hermitian then the solver can exploit this fact to perform
596
   less operations or provide better stability.
597
 
598
   Level: intermediate
599
 
600
.seealso: QEPSetOperators(), QEPSetType(), QEPGetProblemType(), QEPProblemType
601
@*/
602
PetscErrorCode QEPSetProblemType(QEP qep,QEPProblemType type)
603
{
604
  PetscFunctionBegin;
2213 jroman 605
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1906 jroman 606
  if (type!=QEP_GENERAL && type!=QEP_HERMITIAN && type!=QEP_GYROSCOPIC)
2214 jroman 607
    SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
1906 jroman 608
  qep->problem_type = type;
609
  PetscFunctionReturn(0);
610
}
611
 
612
#undef __FUNCT__  
613
#define __FUNCT__ "QEPGetProblemType"
614
/*@C
615
   QEPGetProblemType - Gets the problem type from the QEP object.
616
 
617
   Not Collective
618
 
619
   Input Parameter:
620
.  qep - the quadratic eigensolver context
621
 
622
   Output Parameter:
623
.  type - name of QEP problem type
624
 
625
   Level: intermediate
626
 
627
.seealso: QEPSetProblemType(), QEPProblemType
628
@*/
629
PetscErrorCode QEPGetProblemType(QEP qep,QEPProblemType *type)
630
{
631
  PetscFunctionBegin;
2213 jroman 632
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1906 jroman 633
  PetscValidPointer(type,2);
634
  *type = qep->problem_type;
635
  PetscFunctionReturn(0);
636
}
637
 
638
#undef __FUNCT__  
1887 jroman 639
#define __FUNCT__ "QEPSetConvergenceTest"
640
/*@C
2070 jroman 641
    QEPSetConvergenceTest - Sets a function to compute the error estimate used in
642
    the convergence test.
1887 jroman 643
 
644
    Collective on QEP
645
 
646
    Input Parameters:
647
+   qep  - eigensolver context obtained from QEPCreate()
648
.   func - a pointer to the convergence test function
649
-   ctx  - a context pointer (the last parameter to the convergence test function)
650
 
651
    Calling Sequence of func:
2070 jroman 652
$   func(QEP qep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal* errest,void *ctx)
1887 jroman 653
 
654
+   qep    - eigensolver context obtained from QEPCreate()
2070 jroman 655
.   eigr   - real part of the eigenvalue
656
.   eigi   - imaginary part of the eigenvalue
657
.   res    - residual norm associated to the eigenpair
658
.   errest - (output) computed error estimate
1887 jroman 659
-   ctx    - optional context, as set by QEPSetConvergenceTest()
660
 
661
    Note:
2070 jroman 662
    If the error estimate returned by the convergence test function is less than
663
    the tolerance, then the eigenvalue is accepted as converged.
664
 
1887 jroman 665
    Level: advanced
666
 
667
.seealso: QEPSetTolerances()
668
@*/
2070 jroman 669
EXTERN PetscErrorCode QEPSetConvergenceTest(QEP qep,PetscErrorCode (*func)(QEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx)
1887 jroman 670
{
671
  PetscFunctionBegin;
672
  qep->conv_func = func;
673
  qep->conv_ctx = ctx;
674
  PetscFunctionReturn(0);
675
}
676
 
677
#undef __FUNCT__  
2054 eromero 678
#define __FUNCT__ "QEPSetTrackAll"
679
/*@
680
    QEPSetTrackAll - Specifies if the solver must compute the residual of all
681
    approximate eigenpairs or not.
682
 
683
    Collective on QEP
684
 
685
    Input Parameters:
686
+   qep      - the eigensolver context
687
-   trackall - whether compute all residuals or not
688
 
689
    Notes:
690
    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
691
    the residual for each eigenpair approximation. Computing the residual is
692
    usually an expensive operation and solvers commonly compute the associated
693
    residual to the first unconverged eigenpair.
694
    The options '-qep_monitor_all' and '-qep_monitor_draw_all' automatically
695
    activates this option.
696
 
697
    Level: intermediate
698
 
699
.seealso: EPSGetTrackAll()
700
@*/
2216 jroman 701
PetscErrorCode QEPSetTrackAll(QEP qep,PetscBool trackall)
2054 eromero 702
{
703
  PetscFunctionBegin;
2213 jroman 704
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2054 eromero 705
  qep->trackall = trackall;
706
  PetscFunctionReturn(0);
707
}
708
 
709
#undef __FUNCT__  
710
#define __FUNCT__ "QEPGetTrackAll"
711
/*@
2185 jroman 712
    QEPGetTrackAll - Returns the flag indicating whether all residual norms must
713
    be computed or not.
2054 eromero 714
 
715
    Not Collective
716
 
717
    Input Parameter:
718
.   qep - the eigensolver context
719
 
720
    Output Parameter:
721
.   trackall - the returned flag
722
 
723
    Level: intermediate
724
 
725
.seealso: EPSSetTrackAll()
726
@*/
2216 jroman 727
PetscErrorCode QEPGetTrackAll(QEP qep,PetscBool *trackall)
2054 eromero 728
{
729
  PetscFunctionBegin;
2213 jroman 730
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
2054 eromero 731
  PetscValidPointer(trackall,2);
732
  *trackall = qep->trackall;
733
  PetscFunctionReturn(0);
734
}
735
 
736
#undef __FUNCT__  
1887 jroman 737
#define __FUNCT__ "QEPSetOptionsPrefix"
738
/*@C
739
   QEPSetOptionsPrefix - Sets the prefix used for searching for all
740
   QEP options in the database.
741
 
742
   Collective on QEP
743
 
744
   Input Parameters:
745
+  qep - the quadratic eigensolver context
746
-  prefix - the prefix string to prepend to all QEP option requests
747
 
748
   Notes:
749
   A hyphen (-) must NOT be given at the beginning of the prefix name.
750
   The first character of all runtime options is AUTOMATICALLY the
751
   hyphen.
752
 
753
   For example, to distinguish between the runtime options for two
754
   different QEP contexts, one could call
755
.vb
756
      QEPSetOptionsPrefix(qep1,"qeig1_")
757
      QEPSetOptionsPrefix(qep2,"qeig2_")
758
.ve
759
 
760
   Level: advanced
761
 
762
.seealso: QEPAppendOptionsPrefix(), QEPGetOptionsPrefix()
763
@*/
764
PetscErrorCode QEPSetOptionsPrefix(QEP qep,const char *prefix)
765
{
766
  PetscErrorCode ierr;
767
  PetscFunctionBegin;
2213 jroman 768
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 769
  ierr = PetscObjectSetOptionsPrefix((PetscObject)qep,prefix);CHKERRQ(ierr);
770
  ierr = IPSetOptionsPrefix(qep->ip,prefix);CHKERRQ(ierr);
771
  ierr = IPAppendOptionsPrefix(qep->ip,"qep_");CHKERRQ(ierr);
772
  PetscFunctionReturn(0);  
773
}
774
 
775
#undef __FUNCT__  
776
#define __FUNCT__ "QEPAppendOptionsPrefix"
777
/*@C
778
   QEPAppendOptionsPrefix - Appends to the prefix used for searching for all
779
   QEP options in the database.
780
 
781
   Collective on QEP
782
 
783
   Input Parameters:
784
+  qep - the quadratic eigensolver context
785
-  prefix - the prefix string to prepend to all QEP option requests
786
 
787
   Notes:
788
   A hyphen (-) must NOT be given at the beginning of the prefix name.
789
   The first character of all runtime options is AUTOMATICALLY the hyphen.
790
 
791
   Level: advanced
792
 
793
.seealso: QEPSetOptionsPrefix(), QEPGetOptionsPrefix()
794
@*/
795
PetscErrorCode QEPAppendOptionsPrefix(QEP qep,const char *prefix)
796
{
797
  PetscErrorCode ierr;
798
  PetscFunctionBegin;
2213 jroman 799
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 800
  ierr = PetscObjectAppendOptionsPrefix((PetscObject)qep, prefix);CHKERRQ(ierr);
801
  ierr = IPSetOptionsPrefix(qep->ip,prefix);CHKERRQ(ierr);
802
  ierr = IPAppendOptionsPrefix(qep->ip,"qep_");CHKERRQ(ierr);
803
  PetscFunctionReturn(0);
804
}
805
 
806
#undef __FUNCT__  
807
#define __FUNCT__ "QEPGetOptionsPrefix"
808
/*@C
809
   QEPGetOptionsPrefix - Gets the prefix used for searching for all
810
   QEP options in the database.
811
 
812
   Not Collective
813
 
814
   Input Parameters:
815
.  qep - the quadratic eigensolver context
816
 
817
   Output Parameters:
818
.  prefix - pointer to the prefix string used is returned
819
 
820
   Notes: On the fortran side, the user should pass in a string 'prefix' of
821
   sufficient length to hold the prefix.
822
 
823
   Level: advanced
824
 
825
.seealso: QEPSetOptionsPrefix(), QEPAppendOptionsPrefix()
826
@*/
827
PetscErrorCode QEPGetOptionsPrefix(QEP qep,const char *prefix[])
828
{
829
  PetscErrorCode ierr;
830
  PetscFunctionBegin;
2213 jroman 831
  PetscValidHeaderSpecific(qep,QEP_CLASSID,1);
1887 jroman 832
  PetscValidPointer(prefix,2);
833
  ierr = PetscObjectGetOptionsPrefix((PetscObject)qep,prefix);CHKERRQ(ierr);
834
  PetscFunctionReturn(0);
835
}