Subversion Repositories slepc-dev

Rev

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