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
545 dsic.upv.es!jroman 1
/*
2
      EPS routines related to options that can be set via the command-line
3
      or procedurally.
1376 slepc 4
 
5
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 6
   SLEPc - Scalable Library for Eigenvalue Problem Computations
7
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
1376 slepc 8
 
1672 slepc 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/>.
1376 slepc 22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
545 dsic.upv.es!jroman 23
*/
1376 slepc 24
 
1521 slepc 25
#include "private/epsimpl.h"   /*I "slepceps.h" I*/
526 dsic.upv.es!antodo 26
 
27
#undef __FUNCT__  
28
#define __FUNCT__ "EPSSetFromOptions"
29
/*@
30
   EPSSetFromOptions - Sets EPS options from the options database.
31
   This routine must be called before EPSSetUp() if the user is to be
32
   allowed to set the solver type.
33
 
34
   Collective on EPS
35
 
36
   Input Parameters:
37
.  eps - the eigensolver context
38
 
39
   Notes:  
40
   To see all options, run your program with the -help option.
41
 
42
   Level: beginner
43
@*/
44
PetscErrorCode EPSSetFromOptions(EPS eps)
45
{
46
  PetscErrorCode ierr;
1076 slepc 47
  char           type[256],monfilename[PETSC_MAX_PATH_LEN];
1944 jroman 48
  PetscTruth     flg,val;
1957 jroman 49
  PetscReal      r,nrma,nrmb;
1425 slepc 50
  PetscScalar    s;
1575 slepc 51
  PetscInt       i,j,k;
1804 jroman 52
  const char     *bal_list[4] = { "none", "oneside", "twoside","user" };
1532 slepc 53
  PetscViewerASCIIMonitor monviewer;
1950 jroman 54
  EPSMONITOR_CONV *ctx;
526 dsic.upv.es!antodo 55
 
56
  PetscFunctionBegin;
57
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1422 slepc 58
  ierr = PetscOptionsBegin(((PetscObject)eps)->comm,((PetscObject)eps)->prefix,"Eigenproblem Solver (EPS) Options","EPS");CHKERRQ(ierr);
59
    ierr = PetscOptionsList("-eps_type","Eigenproblem Solver method","EPSSetType",EPSList,(char*)(((PetscObject)eps)->type_name?((PetscObject)eps)->type_name:EPSKRYLOVSCHUR),type,256,&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 60
    if (flg) {
61
      ierr = EPSSetType(eps,type);CHKERRQ(ierr);
62
    }
63
 
830 dsic.upv.es!antodo 64
    ierr = PetscOptionsTruthGroupBegin("-eps_hermitian","hermitian eigenvalue problem","EPSSetProblemType",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 65
    if (flg) {ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);}
830 dsic.upv.es!antodo 66
    ierr = PetscOptionsTruthGroup("-eps_gen_hermitian","generalized hermitian eigenvalue problem","EPSSetProblemType",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 67
    if (flg) {ierr = EPSSetProblemType(eps,EPS_GHEP);CHKERRQ(ierr);}
830 dsic.upv.es!antodo 68
    ierr = PetscOptionsTruthGroup("-eps_non_hermitian","non-hermitian eigenvalue problem","EPSSetProblemType",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 69
    if (flg) {ierr = EPSSetProblemType(eps,EPS_NHEP);CHKERRQ(ierr);}
1358 slepc 70
    ierr = PetscOptionsTruthGroup("-eps_gen_non_hermitian","generalized non-hermitian eigenvalue problem","EPSSetProblemType",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 71
    if (flg) {ierr = EPSSetProblemType(eps,EPS_GNHEP);CHKERRQ(ierr);}
1915 jroman 72
    ierr = PetscOptionsTruthGroup("-eps_pos_gen_non_hermitian","generalized non-hermitian eigenvalue problem with positive semi-definite B","EPSSetProblemType",&flg);CHKERRQ(ierr);
1358 slepc 73
    if (flg) {ierr = EPSSetProblemType(eps,EPS_PGNHEP);CHKERRQ(ierr);}
1915 jroman 74
    ierr = PetscOptionsTruthGroupEnd("-eps_gen_indefinite","generalized hermitian-indefinite eigenvalue problem","EPSSetProblemType",&flg);CHKERRQ(ierr);
75
    if (flg) {ierr = EPSSetProblemType(eps,EPS_GHIEP);CHKERRQ(ierr);}
526 dsic.upv.es!antodo 76
 
657 dsic.upv.es!antodo 77
    /*
78
      Set the type if it was never set.
79
    */
1422 slepc 80
    if (!((PetscObject)eps)->type_name) {
1198 slepc 81
      ierr = EPSSetType(eps,EPSKRYLOVSCHUR);CHKERRQ(ierr);
657 dsic.upv.es!antodo 82
    }
83
 
1560 slepc 84
    ierr = PetscOptionsTruthGroupBegin("-eps_ritz","Rayleigh-Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
85
    if (flg) {ierr = EPSSetExtraction(eps,EPS_RITZ);CHKERRQ(ierr);}
86
    ierr = PetscOptionsTruthGroup("-eps_harmonic","harmonic Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
87
    if (flg) {ierr = EPSSetExtraction(eps,EPS_HARMONIC);CHKERRQ(ierr);}
1987 eromero 88
    ierr = PetscOptionsTruthGroup("-eps_harmonic_relative","relative harmonic Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
89
    if (flg) {ierr = EPSSetExtraction(eps,EPS_HARMONIC_RELATIVE);CHKERRQ(ierr);}
90
    ierr = PetscOptionsTruthGroup("-eps_harmonic_right","right harmonic Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
91
    if (flg) {ierr = EPSSetExtraction(eps,EPS_HARMONIC_RIGHT);CHKERRQ(ierr);}
92
    ierr = PetscOptionsTruthGroup("-eps_harmonic_largest","largest harmonic Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
93
    if (flg) {ierr = EPSSetExtraction(eps,EPS_HARMONIC_LARGEST);CHKERRQ(ierr);}
1560 slepc 94
    ierr = PetscOptionsTruthGroup("-eps_refined","refined Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
95
    if (flg) {ierr = EPSSetExtraction(eps,EPS_REFINED);CHKERRQ(ierr);}
96
    ierr = PetscOptionsTruthGroupEnd("-eps_refined_harmonic","refined harmonic Ritz extraction","EPSSetExtraction",&flg);CHKERRQ(ierr);
97
    if (flg) {ierr = EPSSetExtraction(eps,EPS_REFINED_HARMONIC);CHKERRQ(ierr);}
1426 slepc 98
 
1940 jroman 99
    if (!eps->balance) eps->balance = EPS_BALANCE_NONE;
100
    ierr = PetscOptionsEList("-eps_balance", "Balancing method","EPSSetBalance",bal_list,4,bal_list[eps->balance-EPS_BALANCE_NONE],&i,&flg);CHKERRQ(ierr);
101
    if (flg) { eps->balance = (EPSBalance)(i+EPS_BALANCE_NONE); }
1799 jroman 102
    r = j = PETSC_IGNORE;
103
    ierr = PetscOptionsInt("-eps_balance_its","Number of iterations in balancing","EPSSetBalance",eps->balance_its,&j,PETSC_NULL);CHKERRQ(ierr);
104
    ierr = PetscOptionsReal("-eps_balance_cutoff","Cutoff value in balancing","EPSSetBalance",eps->balance_cutoff,&r,PETSC_NULL);CHKERRQ(ierr);
1824 antodo 105
    ierr = EPSSetBalance(eps,(EPSBalance)PETSC_IGNORE,j,r);CHKERRQ(ierr);
1799 jroman 106
 
1282 slepc 107
    r = i = PETSC_IGNORE;
108
    ierr = PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&i,PETSC_NULL);CHKERRQ(ierr);
109
    ierr = PetscOptionsReal("-eps_tol","Tolerance","EPSSetTolerances",eps->tol,&r,PETSC_NULL);CHKERRQ(ierr);
1284 slepc 110
    ierr = EPSSetTolerances(eps,r,i);CHKERRQ(ierr);
1785 antodo 111
    ierr = PetscOptionsTruthGroupBegin("-eps_convergence_default","Default (relative error) convergence test","EPSSetConvergenceTest",&flg);CHKERRQ(ierr);
112
    if (flg) {ierr = EPSSetConvergenceTest(eps,EPSDefaultConverged,PETSC_NULL);CHKERRQ(ierr);}
1886 jroman 113
    ierr = PetscOptionsTruthGroup("-eps_convergence_absolute","Absolute error convergence test","EPSSetConvergenceTest",&flg);CHKERRQ(ierr);
1785 antodo 114
    if (flg) {ierr = EPSSetConvergenceTest(eps,EPSAbsoluteConverged,PETSC_NULL);CHKERRQ(ierr);}
1794 antodo 115
    ierr = PetscOptionsTruthGroupEnd("-eps_convergence_residual","Residual convergence test","EPSSetConvergenceTest",&flg);CHKERRQ(ierr);
116
    if (flg) {ierr = EPSSetConvergenceTest(eps,EPSResidualConverged,PETSC_NULL);CHKERRQ(ierr);}
526 dsic.upv.es!antodo 117
 
1599 slepc 118
    i = j = k = PETSC_IGNORE;
1282 slepc 119
    ierr = PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,PETSC_NULL);CHKERRQ(ierr);
120
    ierr = PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,PETSC_NULL);CHKERRQ(ierr);
1575 slepc 121
    ierr = PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,PETSC_NULL);CHKERRQ(ierr);
122
    ierr = EPSSetDimensions(eps,i,j,k);CHKERRQ(ierr);
1282 slepc 123
 
526 dsic.upv.es!antodo 124
    /* -----------------------------------------------------------------------*/
125
    /*
126
      Cancels all monitors hardwired into code before call to EPSSetFromOptions()
127
    */
1713 antodo 128
    flg  = PETSC_FALSE;
129
    ierr = PetscOptionsTruth("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
526 dsic.upv.es!antodo 130
    if (flg) {
1331 slepc 131
      ierr = EPSMonitorCancel(eps); CHKERRQ(ierr);
526 dsic.upv.es!antodo 132
    }
133
    /*
134
      Prints approximate eigenvalues and error estimates at each iteration
135
    */
1331 slepc 136
    ierr = PetscOptionsString("-eps_monitor","Monitor approximate eigenvalues and error estimates","EPSMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 137
    if (flg) {
1532 slepc 138
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)eps)->comm,monfilename,((PetscObject)eps)->tablevel,&monviewer);CHKERRQ(ierr);
139
      ierr = EPSMonitorSet(eps,EPSMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
526 dsic.upv.es!antodo 140
    }
1721 antodo 141
    ierr = PetscOptionsString("-eps_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","EPSMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
142
    if (flg) {
1950 jroman 143
      ierr = PetscNew(EPSMONITOR_CONV,&ctx);CHKERRQ(ierr);
144
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)eps)->comm,monfilename,((PetscObject)eps)->tablevel,&ctx->viewer);CHKERRQ(ierr);
145
      ierr = EPSMonitorSet(eps,EPSMonitorConverged,ctx,(PetscErrorCode (*)(void*))EPSMonitorDestroy_Converged);CHKERRQ(ierr);
1721 antodo 146
    }
147
    ierr = PetscOptionsString("-eps_monitor_first","Monitor first unconverged approximate eigenvalue and error estimate","EPSMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
148
    if (flg) {
149
      ierr = PetscViewerASCIIMonitorCreate(((PetscObject)eps)->comm,monfilename,((PetscObject)eps)->tablevel,&monviewer);CHKERRQ(ierr);
150
      ierr = EPSMonitorSet(eps,EPSMonitorFirst,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
151
    }
1713 antodo 152
    flg = PETSC_FALSE;
153
    ierr = PetscOptionsTruth("-eps_monitor_draw","Monitor error estimates graphically","EPSMonitorSet",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
623 dsic.upv.es!antodo 154
    if (flg) {
1331 slepc 155
      ierr = EPSMonitorSet(eps,EPSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
623 dsic.upv.es!antodo 156
    }
526 dsic.upv.es!antodo 157
  /* -----------------------------------------------------------------------*/
1782 antodo 158
    ierr = PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg);CHKERRQ(ierr);
159
    if (flg) {
160
      ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);CHKERRQ(ierr);
161
      ierr = EPSSetTarget(eps,s);CHKERRQ(ierr);
162
    }
163
 
830 dsic.upv.es!antodo 164
    ierr = PetscOptionsTruthGroupBegin("-eps_largest_magnitude","compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 165
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE);CHKERRQ(ierr);}
830 dsic.upv.es!antodo 166
    ierr = PetscOptionsTruthGroup("-eps_smallest_magnitude","compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 167
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE);CHKERRQ(ierr);}
830 dsic.upv.es!antodo 168
    ierr = PetscOptionsTruthGroup("-eps_largest_real","compute largest real parts","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 169
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);CHKERRQ(ierr);}
830 dsic.upv.es!antodo 170
    ierr = PetscOptionsTruthGroup("-eps_smallest_real","compute smallest real parts","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 171
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);CHKERRQ(ierr);}
830 dsic.upv.es!antodo 172
    ierr = PetscOptionsTruthGroup("-eps_largest_imaginary","compute largest imaginary parts","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 173
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY);CHKERRQ(ierr);}
1782 antodo 174
    ierr = PetscOptionsTruthGroup("-eps_smallest_imaginary","compute smallest imaginary parts","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
526 dsic.upv.es!antodo 175
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY);CHKERRQ(ierr);}
1782 antodo 176
    ierr = PetscOptionsTruthGroup("-eps_target_magnitude","compute nearest eigenvalues to target","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
177
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);CHKERRQ(ierr);}
178
    ierr = PetscOptionsTruthGroup("-eps_target_real","compute eigenvalues with real parts close to target","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
179
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL);CHKERRQ(ierr);}
180
    ierr = PetscOptionsTruthGroupEnd("-eps_target_imaginary","compute eigenvalues with imaginary parts close to target","EPSSetWhichEigenpairs",&flg);CHKERRQ(ierr);
181
    if (flg) {ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY);CHKERRQ(ierr);}
526 dsic.upv.es!antodo 182
 
1944 jroman 183
    ierr = PetscOptionsTruth("-eps_left_vectors","Compute left eigenvectors also","EPSSetLeftVectorsWanted",eps->leftvecs,&val,&flg);CHKERRQ(ierr);
184
    if (flg) {
185
      ierr = EPSSetLeftVectorsWanted(eps,val);CHKERRQ(ierr);
186
    }
2031 jroman 187
    ierr = PetscOptionsTruth("-eps_true_residual","Compute true residuals explicitly","EPSSetTrueResidual",eps->trueres,&val,&flg);CHKERRQ(ierr);
188
    if (flg) {
189
      ierr = EPSSetTrueResidual(eps,val);CHKERRQ(ierr);
190
    }
1944 jroman 191
 
1957 jroman 192
    nrma = nrmb = PETSC_IGNORE;
193
    ierr = PetscOptionsReal("-eps_norm_a","Norm of matrix A","EPSSetMatrixNorms",eps->nrma,&nrma,PETSC_NULL);CHKERRQ(ierr);
194
    ierr = PetscOptionsReal("-eps_norm_b","Norm of matrix B","EPSSetMatrixNorms",eps->nrmb,&nrmb,PETSC_NULL);CHKERRQ(ierr);
195
    ierr = EPSSetMatrixNorms(eps,nrma,nrmb,eps->adaptive);CHKERRQ(ierr);
196
    ierr = PetscOptionsTruth("-eps_norms_adaptive","Update the value of matrix norms adaptively","EPSSetMatrixNorms",eps->adaptive,&val,&flg);CHKERRQ(ierr);
197
    if (flg) {
198
      ierr = EPSSetMatrixNorms(eps,PETSC_IGNORE,PETSC_IGNORE,val);CHKERRQ(ierr);
199
    }
200
 
526 dsic.upv.es!antodo 201
    ierr = PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",0);CHKERRQ(ierr);
202
    ierr = PetscOptionsName("-eps_view_binary","Save the matrices associated to the eigenproblem","EPSSetFromOptions",0);CHKERRQ(ierr);
203
    ierr = PetscOptionsName("-eps_plot_eigs","Make a plot of the computed eigenvalues","EPSSolve",0);CHKERRQ(ierr);
1209 slepc 204
 
526 dsic.upv.es!antodo 205
    if (eps->ops->setfromoptions) {
206
      ierr = (*eps->ops->setfromoptions)(eps);CHKERRQ(ierr);
207
    }
208
  ierr = PetscOptionsEnd();CHKERRQ(ierr);
209
 
1345 slepc 210
  ierr = IPSetFromOptions(eps->ip); CHKERRQ(ierr);
526 dsic.upv.es!antodo 211
  ierr = STSetFromOptions(eps->OP); CHKERRQ(ierr);
212
  PetscFunctionReturn(0);
213
}
214
 
215
#undef __FUNCT__  
216
#define __FUNCT__ "EPSGetTolerances"
217
/*@
1811 jroman 218
   EPSGetTolerances - Gets the tolerance and maximum iteration count used
219
   by the EPS convergence tests.
526 dsic.upv.es!antodo 220
 
221
   Not Collective
222
 
223
   Input Parameter:
224
.  eps - the eigensolver context
225
 
226
   Output Parameters:
227
+  tol - the convergence tolerance
228
-  maxits - maximum number of iterations
229
 
230
   Notes:
231
   The user can specify PETSC_NULL for any parameter that is not needed.
232
 
233
   Level: intermediate
234
 
235
.seealso: EPSSetTolerances()
236
@*/
1509 slepc 237
PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)
526 dsic.upv.es!antodo 238
{
239
  PetscFunctionBegin;
240
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
241
  if (tol)    *tol    = eps->tol;
242
  if (maxits) *maxits = eps->max_it;
243
  PetscFunctionReturn(0);
244
}
245
 
246
#undef __FUNCT__  
247
#define __FUNCT__ "EPSSetTolerances"
248
/*@
1811 jroman 249
   EPSSetTolerances - Sets the tolerance and maximum iteration count used
250
   by the EPS convergence tests.
526 dsic.upv.es!antodo 251
 
252
   Collective on EPS
253
 
254
   Input Parameters:
255
+  eps - the eigensolver context
256
.  tol - the convergence tolerance
257
-  maxits - maximum number of iterations to use
258
 
259
   Options Database Keys:
260
+  -eps_tol <tol> - Sets the convergence tolerance
261
-  -eps_max_it <maxits> - Sets the maximum number of iterations allowed
262
 
263
   Notes:
1284 slepc 264
   Use PETSC_IGNORE for an argument that need not be changed.
526 dsic.upv.es!antodo 265
 
1282 slepc 266
   Use PETSC_DECIDE for maxits to assign a reasonably good value, which is
267
   dependent on the solution method.
268
 
526 dsic.upv.es!antodo 269
   Level: intermediate
270
 
271
.seealso: EPSGetTolerances()
272
@*/
1509 slepc 273
PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)
526 dsic.upv.es!antodo 274
{
275
  PetscFunctionBegin;
276
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1282 slepc 277
  if (tol != PETSC_IGNORE) {
278
    if (tol == PETSC_DEFAULT) {
279
      eps->tol = 1e-7;
280
    } else {
281
      if (tol < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
282
      eps->tol = tol;
283
    }
1276 slepc 284
  }
1282 slepc 285
  if (maxits != PETSC_IGNORE) {
286
    if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
287
      eps->max_it = 0;
1886 jroman 288
      eps->setupcalled = 0;
1282 slepc 289
    } else {
290
      if (maxits < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
291
      eps->max_it = maxits;
292
    }
1276 slepc 293
  }
526 dsic.upv.es!antodo 294
  PetscFunctionReturn(0);
295
}
296
 
297
#undef __FUNCT__  
298
#define __FUNCT__ "EPSGetDimensions"
299
/*@
300
   EPSGetDimensions - Gets the number of eigenvalues to compute
301
   and the dimension of the subspace.
302
 
303
   Not Collective
304
 
305
   Input Parameter:
306
.  eps - the eigensolver context
307
 
308
   Output Parameters:
309
+  nev - number of eigenvalues to compute
1575 slepc 310
.  ncv - the maximum dimension of the subspace to be used by the solver
311
-  mpd - the maximum dimension allowed for the projected problem
526 dsic.upv.es!antodo 312
 
313
   Notes:
314
   The user can specify PETSC_NULL for any parameter that is not needed.
315
 
316
   Level: intermediate
317
 
318
.seealso: EPSSetDimensions()
319
@*/
1575 slepc 320
PetscErrorCode EPSGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
526 dsic.upv.es!antodo 321
{
322
  PetscFunctionBegin;
323
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1575 slepc 324
  if (nev) *nev = eps->nev;
325
  if (ncv) *ncv = eps->ncv;
326
  if (mpd) *mpd = eps->mpd;
526 dsic.upv.es!antodo 327
  PetscFunctionReturn(0);
328
}
329
 
330
#undef __FUNCT__  
331
#define __FUNCT__ "EPSSetDimensions"
332
/*@
333
   EPSSetDimensions - Sets the number of eigenvalues to compute
334
   and the dimension of the subspace.
335
 
336
   Collective on EPS
337
 
338
   Input Parameters:
339
+  eps - the eigensolver context
340
.  nev - number of eigenvalues to compute
1575 slepc 341
.  ncv - the maximum dimension of the subspace to be used by the solver
342
-  mpd - the maximum dimension allowed for the projected problem
526 dsic.upv.es!antodo 343
 
344
   Options Database Keys:
345
+  -eps_nev <nev> - Sets the number of eigenvalues
1575 slepc 346
.  -eps_ncv <ncv> - Sets the dimension of the subspace
347
-  -eps_mpd <mpd> - Sets the maximum projected dimension
526 dsic.upv.es!antodo 348
 
349
   Notes:
1282 slepc 350
   Use PETSC_IGNORE to retain the previous value of any parameter.
526 dsic.upv.es!antodo 351
 
1575 slepc 352
   Use PETSC_DECIDE for ncv and mpd to assign a reasonably good value, which is
526 dsic.upv.es!antodo 353
   dependent on the solution method.
354
 
1575 slepc 355
   The parameters ncv and mpd are intimately related, so that the user is advised
1799 jroman 356
   to set one of them at most. Normal usage is the following
1575 slepc 357
+  - In cases where nev is small, the user sets ncv (a reasonable default is 2*nev).
358
-  - In cases where nev is large, the user sets mpd.
359
 
360
   The value of ncv should always be between nev and (nev+mpd), typically
361
   ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
362
   a smaller value should be used.
363
 
526 dsic.upv.es!antodo 364
   Level: intermediate
365
 
366
.seealso: EPSGetDimensions()
367
@*/
1575 slepc 368
PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
526 dsic.upv.es!antodo 369
{
370
  PetscFunctionBegin;
371
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
372
 
1282 slepc 373
  if( nev != PETSC_IGNORE ) {
374
    if (nev<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
526 dsic.upv.es!antodo 375
    eps->nev = nev;
376
    eps->setupcalled = 0;
377
  }
1282 slepc 378
  if( ncv != PETSC_IGNORE ) {
379
    if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
380
      eps->ncv = 0;
381
    } else {
382
      if (ncv<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
526 dsic.upv.es!antodo 383
      eps->ncv = ncv;
384
    }
385
    eps->setupcalled = 0;
386
  }
1575 slepc 387
  if( mpd != PETSC_IGNORE ) {
388
    if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
389
      eps->mpd = 0;
390
    } else {
391
      if (mpd<1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
392
      eps->mpd = mpd;
393
    }
394
  }
526 dsic.upv.es!antodo 395
  PetscFunctionReturn(0);
396
}
397
 
398
#undef __FUNCT__  
399
#define __FUNCT__ "EPSSetWhichEigenpairs"
400
/*@
401
    EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
402
    to be sought.
403
 
404
    Collective on EPS
405
 
1782 antodo 406
    Input Parameters:
1942 jroman 407
+   eps   - eigensolver context obtained from EPSCreate()
1782 antodo 408
-   which - the portion of the spectrum to be sought
526 dsic.upv.es!antodo 409
 
410
    Possible values:
1799 jroman 411
    The parameter 'which' can have one of these values
526 dsic.upv.es!antodo 412
 
413
+     EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
414
.     EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
415
.     EPS_LARGEST_REAL - largest real parts
416
.     EPS_SMALLEST_REAL - smallest real parts
417
.     EPS_LARGEST_IMAGINARY - largest imaginary parts
1782 antodo 418
.     EPS_SMALLEST_IMAGINARY - smallest imaginary parts
1811 jroman 419
.     EPS_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
420
.     EPS_TARGET_REAL - eigenvalues with real part closest to target
421
.     EPS_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
1945 jroman 422
-     EPS_WHICH_USER - user defined ordering set with EPSSetEigenvalueComparison()
526 dsic.upv.es!antodo 423
 
424
    Options Database Keys:
425
+   -eps_largest_magnitude - Sets largest eigenvalues in magnitude
426
.   -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
427
.   -eps_largest_real - Sets largest real parts
428
.   -eps_smallest_real - Sets smallest real parts
1811 jroman 429
.   -eps_largest_imaginary - Sets largest imaginary parts
430
.   -eps_smallest_imaginary - Sets smallest imaginary parts
431
.   -eps_target_magnitude - Sets eigenvalues closest to target
432
.   -eps_target_real - Sets real parts closest to target
433
-   -eps_target_imaginary - Sets imaginary parts closest to target
526 dsic.upv.es!antodo 434
 
435
    Notes:
436
    Not all eigensolvers implemented in EPS account for all the possible values
437
    stated above. Also, some values make sense only for certain types of
438
    problems. If SLEPc is compiled for real numbers EPS_LARGEST_IMAGINARY
439
    and EPS_SMALLEST_IMAGINARY use the absolute value of the imaginary part
1811 jroman 440
    for eigenvalue selection.
526 dsic.upv.es!antodo 441
 
1811 jroman 442
    The target is a scalar value provided with EPSSetTarget().
443
 
1815 jroman 444
    The criterion EPS_TARGET_IMAGINARY is available only in case PETSc and
445
    SLEPc have been built with complex scalars.
446
 
526 dsic.upv.es!antodo 447
    Level: intermediate
448
 
1811 jroman 449
.seealso: EPSGetWhichEigenpairs(), EPSSetTarget(), EPSSetEigenvalueComparison(), EPSSortEigenvalues(), EPSWhich
526 dsic.upv.es!antodo 450
@*/
451
PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)
452
{
453
  PetscFunctionBegin;
454
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1942 jroman 455
  if (which!=PETSC_IGNORE) {
456
    if (which==PETSC_DECIDE || which==PETSC_DEFAULT) eps->which = (EPSWhich)0;
457
    else switch (which) {
458
      case EPS_LARGEST_MAGNITUDE:
459
      case EPS_SMALLEST_MAGNITUDE:
460
      case EPS_LARGEST_REAL:
461
      case EPS_SMALLEST_REAL:
462
      case EPS_LARGEST_IMAGINARY:
463
      case EPS_SMALLEST_IMAGINARY:
464
      case EPS_TARGET_MAGNITUDE:
465
      case EPS_TARGET_REAL:
1815 jroman 466
#if defined(PETSC_USE_COMPLEX)
1942 jroman 467
      case EPS_TARGET_IMAGINARY:
1815 jroman 468
#endif
1945 jroman 469
      case EPS_WHICH_USER:
1942 jroman 470
        if (eps->which != which) {
471
          eps->setupcalled = 0;
472
          eps->which = which;
473
        }
474
        break;
475
      default:
476
        SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
477
    }
1810 jroman 478
  }
526 dsic.upv.es!antodo 479
  PetscFunctionReturn(0);
480
}
481
 
482
#undef __FUNCT__  
483
#define __FUNCT__ "EPSGetWhichEigenpairs"
707 dsic.upv.es!antodo 484
/*@C
526 dsic.upv.es!antodo 485
    EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
486
    sought.
487
 
488
    Not Collective
489
 
490
    Input Parameter:
491
.   eps - eigensolver context obtained from EPSCreate()
492
 
493
    Output Parameter:
494
.   which - the portion of the spectrum to be sought
495
 
496
    Notes:
1811 jroman 497
    See EPSSetWhichEigenpairs() for possible values of 'which'.
526 dsic.upv.es!antodo 498
 
499
    Level: intermediate
500
 
1364 slepc 501
.seealso: EPSSetWhichEigenpairs(), EPSWhich
526 dsic.upv.es!antodo 502
@*/
503
PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
504
{
505
  PetscFunctionBegin;
506
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1273 slepc 507
  PetscValidPointer(which,2);
526 dsic.upv.es!antodo 508
  *which = eps->which;
509
  PetscFunctionReturn(0);
510
}
511
 
512
#undef __FUNCT__  
1944 jroman 513
#define __FUNCT__ "EPSSetLeftVectorsWanted"
514
/*@
515
    EPSSetLeftVectorsWanted - Specifies which eigenvectors are required.
516
 
517
    Collective on EPS
518
 
519
    Input Parameters:
520
+   eps      - the eigensolver context
521
-   leftvecs - whether left eigenvectors are required or not
522
 
523
    Options Database Keys:
524
.   -eps_left_vectors <boolean> - Sets/resets the boolean flag 'leftvecs'
525
 
526
    Notes:
527
    If the user sets leftvecs=PETSC_TRUE then the solver uses a variant of
528
    the algorithm that computes both right and left eigenvectors. This is
529
    usually much more costly. This option is not available in all solvers.
530
 
531
    Level: intermediate
532
 
533
.seealso: EPSGetLeftVectorsWanted(), EPSGetEigenvectorLeft()
534
@*/
535
PetscErrorCode EPSSetLeftVectorsWanted(EPS eps,PetscTruth leftvecs)
536
{
537
  PetscFunctionBegin;
538
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1948 jroman 539
  if (eps->leftvecs != leftvecs) {
540
    eps->leftvecs = leftvecs;
541
    eps->setupcalled = 0;
542
  }
1944 jroman 543
  PetscFunctionReturn(0);
544
}
545
 
546
#undef __FUNCT__  
547
#define __FUNCT__ "EPSGetLeftVectorsWanted"
548
/*@C
549
    EPSGetLeftVectorsWanted - Returns the flag indicating whether left
550
    eigenvectors are required or not.
551
 
552
    Not Collective
553
 
554
    Input Parameter:
555
.   eps - the eigensolver context
556
 
557
    Output Parameter:
558
.   leftvecs - the returned flag
559
 
560
    Level: intermediate
561
 
562
.seealso: EPSSetLeftVectorsWanted(), EPSWhich
563
@*/
564
PetscErrorCode EPSGetLeftVectorsWanted(EPS eps,PetscTruth *leftvecs)
565
{
566
  PetscFunctionBegin;
567
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
568
  PetscValidPointer(leftvecs,2);
569
  *leftvecs = eps->leftvecs;
570
  PetscFunctionReturn(0);
571
}
572
 
573
#undef __FUNCT__  
1957 jroman 574
#define __FUNCT__ "EPSSetMatrixNorms"
575
/*@
576
    EPSSetMatrixNorms - Gives the reference values of the matrix norms
577
    and specifies whether these values should be improved adaptively.
578
 
579
    Collective on EPS
580
 
581
    Input Parameters:
582
+   eps      - the eigensolver context
583
.   nrma     - a reference value for the norm of matrix A
584
.   nrmb     - a reference value for the norm of matrix B
585
-   adaptive - whether matrix norms are improved adaptively
586
 
587
    Options Database Keys:
588
+   -eps_norm_a <nrma> - norm of A
589
.   -eps_norm_b <nrma> - norm of B
590
-   -eps_norms_adaptive <boolean> - Sets/resets the boolean flag 'adaptive'
591
 
592
    Notes:
593
    If the user sets adaptive=PETSC_FALSE then the solver uses the values
594
    of nrma and nrmb for the matrix norms, and these values do not change
595
    throughout the iteration.
596
 
597
    If the user sets adaptive=PETSC_TRUE then the solver tries to adaptively
598
    improve the supplied values, with the numerical information generated
599
    during the iteration. This option is not available in all solvers.
600
 
601
    If a passed value is PETSC_DEFAULT, the corresponding norm will be set to 1.
602
    If a passed value is PETSC_DETERMINE, the corresponding norm will be computed
603
    as the NORM_INFINITY with MatNorm().
604
 
605
    Level: intermediate
606
 
607
.seealso: EPSGetMatrixNorms()
608
@*/
609
PetscErrorCode EPSSetMatrixNorms(EPS eps,PetscReal nrma,PetscReal nrmb,PetscTruth adaptive)
610
{
611
  PetscFunctionBegin;
612
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
613
  if (nrma != PETSC_IGNORE) {
614
    if (nrma == PETSC_DEFAULT) eps->nrma = 1.0;
615
    else if (nrma == PETSC_DETERMINE) {
616
      eps->nrma = nrma;
617
      eps->setupcalled = 0;
618
    } else {
619
      if (nrma < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nrma. Must be > 0");
620
      eps->nrma = nrma;
621
    }
622
  }
623
  if (nrmb != PETSC_IGNORE) {
624
    if (!eps->isgeneralized) SETERRQ(PETSC_ERR_ARG_WRONG,"Norm of B only allowed in generalized problems");
625
    if (nrmb == PETSC_DEFAULT) eps->nrmb = 1.0;
626
    else if (nrmb == PETSC_DETERMINE) {
627
      eps->nrmb = nrmb;
628
      eps->setupcalled = 0;
629
    } else {
630
      if (nrmb < 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nrmb. Must be > 0");
631
      eps->nrmb = nrmb;
632
    }
633
  }
634
  if (eps->adaptive != adaptive) {
635
    eps->adaptive = adaptive;
636
    eps->setupcalled = 0;
637
  }
638
  PetscFunctionReturn(0);
639
}
640
 
641
#undef __FUNCT__  
642
#define __FUNCT__ "EPSGetMatrixNorms"
643
/*@C
644
    EPSGetMatrixNorms - Returns the value of the matrix norms (either set
645
    by the user or estimated by the solver) and the flag indicating whether
646
    the norms are being adaptively improved.
647
 
648
    Not Collective
649
 
650
    Input Parameter:
651
.   eps - the eigensolver context
652
 
653
    Output Parameters:
654
+   nrma     - the norm of matrix A
655
.   nrmb     - the norm of matrix B
656
-   adaptive - whether matrix norms are improved adaptively
657
 
658
    Level: intermediate
659
 
660
.seealso: EPSSetMatrixNorms()
661
@*/
662
PetscErrorCode EPSGetMatrixNorms(EPS eps,PetscReal *nrma,PetscReal *nrmb,PetscTruth *adaptive)
663
{
664
  PetscFunctionBegin;
665
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
666
  PetscValidPointer(nrma,2);
667
  PetscValidPointer(nrmb,3);
668
  PetscValidPointer(adaptive,4);
669
  if (nrma) *nrma = eps->nrma;
670
  if (nrmb) *nrmb = eps->nrmb;
671
  if (adaptive) *adaptive = eps->adaptive;
672
  PetscFunctionReturn(0);
673
}
674
 
675
#undef __FUNCT__  
1782 antodo 676
#define __FUNCT__ "EPSSetEigenvalueComparison"
677
/*@C
1811 jroman 678
    EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
1945 jroman 679
    when EPSSetWhichEigenpairs() is set to EPS_WHICH_USER.
1811 jroman 680
 
1782 antodo 681
    Collective on EPS
682
 
683
    Input Parameters:
684
+   eps  - eigensolver context obtained from EPSCreate()
685
.   func - a pointer to the comparison function
686
-   ctx  - a context pointer (the last parameter to the comparison function)
687
 
1811 jroman 688
    Calling Sequence of func:
689
$   func(EPS eps,PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
690
 
691
+   eps    - eigensolver context obtained from EPSCreate()
692
.   ar     - real part of the 1st eigenvalue
693
.   ai     - imaginary part of the 1st eigenvalue
694
.   br     - real part of the 2nd eigenvalue
695
.   bi     - imaginary part of the 2nd eigenvalue
696
.   res    - result of comparison
697
-   ctx    - optional context, as set by EPSSetEigenvalueComparison()
698
 
699
    Note:
700
    The comparison function must return an integer less than, equal to, or
701
    greater than zero if the first eigenvalue is considered to be respectively
702
    less than, equal to, or greater than the second one.
1782 antodo 703
 
704
    Level: advanced
705
 
706
.seealso: EPSSetWhichEigenpairs(), EPSSortEigenvalues(), EPSWhich
707
@*/
708
PetscErrorCode EPSSetEigenvalueComparison(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
709
{
710
  PetscFunctionBegin;
711
  eps->which_func = func;
712
  eps->which_ctx = ctx;
713
  PetscFunctionReturn(0);
714
}
715
 
716
#undef __FUNCT__  
1785 antodo 717
#define __FUNCT__ "EPSSetConvergenceTest"
718
/*@C
719
    EPSSetConvergenceTest - Specifies the convergence test.
1811 jroman 720
 
1785 antodo 721
    Collective on EPS
722
 
723
    Input Parameters:
724
+   eps  - eigensolver context obtained from EPSCreate()
725
.   func - a pointer to the convergence test function
726
-   ctx  - a context pointer (the last parameter to the convergence test function)
727
 
1811 jroman 728
    Calling Sequence of func:
2030 jroman 729
$   func(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscTruth *conv,void *ctx)
1811 jroman 730
 
731
+   eps    - eigensolver context obtained from EPSCreate()
2030 jroman 732
.   eigr   - real part of the eigenvalue
733
.   eigi   - imaginary part of the eigenvalue
734
.   res    - computed or estimated residual (on output: error estimate)
735
.   conv   - (output) boolean value, true if the convergence criterion is satisfied
1811 jroman 736
-   ctx    - optional context, as set by EPSSetConvergenceTest()
737
 
1785 antodo 738
    Level: advanced
739
 
740
.seealso: EPSSetTolerances()
741
@*/
2030 jroman 742
EXTERN PetscErrorCode EPSSetConvergenceTest(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscReal*,PetscTruth*,void*),void* ctx)
1785 antodo 743
{
744
  PetscFunctionBegin;
745
  eps->conv_func = func;
746
  eps->conv_ctx = ctx;
747
  PetscFunctionReturn(0);
748
}
749
 
750
#undef __FUNCT__  
526 dsic.upv.es!antodo 751
#define __FUNCT__ "EPSSetProblemType"
752
/*@
753
   EPSSetProblemType - Specifies the type of the eigenvalue problem.
754
 
755
   Collective on EPS
756
 
757
   Input Parameters:
758
+  eps      - the eigensolver context
759
-  type     - a known type of eigenvalue problem
760
 
761
   Options Database Keys:
762
+  -eps_hermitian - Hermitian eigenvalue problem
763
.  -eps_gen_hermitian - generalized Hermitian eigenvalue problem
764
.  -eps_non_hermitian - non-Hermitian eigenvalue problem
1426 slepc 765
.  -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
766
-  -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
767
   with positive semi-definite B
526 dsic.upv.es!antodo 768
 
1426 slepc 769
   Notes:  
526 dsic.upv.es!antodo 770
   Allowed values for the problem type are: Hermitian (EPS_HEP), non-Hermitian
1697 slepc 771
   (EPS_NHEP), generalized Hermitian (EPS_GHEP), generalized non-Hermitian
1915 jroman 772
   (EPS_GNHEP), generalized non-Hermitian with positive semi-definite B
773
   (EPS_PGNHEP), and generalized Hermitian-indefinite (EPS_GHIEP).
526 dsic.upv.es!antodo 774
 
775
   This function must be used to instruct SLEPc to exploit symmetry. If no
776
   problem type is specified, by default a non-Hermitian problem is assumed
777
   (either standard or generalized). If the user knows that the problem is
1697 slepc 778
   Hermitian (i.e. A=A^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and
526 dsic.upv.es!antodo 779
   B positive definite) then it is recommended to set the problem type so
780
   that eigensolver can exploit these properties.
781
 
782
   Level: beginner
783
 
1364 slepc 784
.seealso: EPSSetOperators(), EPSSetType(), EPSGetProblemType(), EPSProblemType
526 dsic.upv.es!antodo 785
@*/
786
PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)
787
{
788
  PetscFunctionBegin;
789
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
790
 
791
  switch (type) {
792
    case EPS_HEP:
793
      eps->isgeneralized = PETSC_FALSE;
794
      eps->ishermitian = PETSC_TRUE;
1358 slepc 795
      eps->ispositive = PETSC_FALSE;
526 dsic.upv.es!antodo 796
      break;      
797
    case EPS_NHEP:
798
      eps->isgeneralized = PETSC_FALSE;
799
      eps->ishermitian = PETSC_FALSE;
1358 slepc 800
      eps->ispositive = PETSC_FALSE;
526 dsic.upv.es!antodo 801
      break;
802
    case EPS_GHEP:
803
      eps->isgeneralized = PETSC_TRUE;
804
      eps->ishermitian = PETSC_TRUE;
1358 slepc 805
      eps->ispositive = PETSC_TRUE;
526 dsic.upv.es!antodo 806
      break;
807
    case EPS_GNHEP:
808
      eps->isgeneralized = PETSC_TRUE;
809
      eps->ishermitian = PETSC_FALSE;
1358 slepc 810
      eps->ispositive = PETSC_FALSE;
526 dsic.upv.es!antodo 811
      break;
1358 slepc 812
    case EPS_PGNHEP:
813
      eps->isgeneralized = PETSC_TRUE;
814
      eps->ishermitian = PETSC_FALSE;
815
      eps->ispositive = PETSC_TRUE;
816
      break;
1915 jroman 817
    case EPS_GHIEP:
818
      eps->isgeneralized = PETSC_TRUE;
819
      eps->ishermitian = PETSC_TRUE;
820
      eps->ispositive = PETSC_FALSE;
821
      break;
526 dsic.upv.es!antodo 822
/*
823
    case EPS_CSEP:
824
      eps->isgeneralized = PETSC_FALSE;
825
      eps->ishermitian = PETSC_FALSE;
826
      ierr = STSetBilinearForm(eps->OP,STINNER_SYMMETRIC);CHKERRQ(ierr);
827
      break;
828
    case EPS_GCSEP:
829
      eps->isgeneralized = PETSC_TRUE;
830
      eps->ishermitian = PETSC_FALSE;
831
      ierr = STSetBilinearForm(eps->OP,STINNER_B_SYMMETRIC);CHKERRQ(ierr);
832
      break;
833
*/
834
    default:
835
      SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
836
  }
837
  eps->problem_type = type;
838
 
839
  PetscFunctionReturn(0);
840
}
841
 
842
#undef __FUNCT__  
843
#define __FUNCT__ "EPSGetProblemType"
707 dsic.upv.es!antodo 844
/*@C
526 dsic.upv.es!antodo 845
   EPSGetProblemType - Gets the problem type from the EPS object.
846
 
847
   Not Collective
848
 
849
   Input Parameter:
850
.  eps - the eigensolver context
851
 
852
   Output Parameter:
853
.  type - name of EPS problem type
854
 
855
   Level: intermediate
856
 
1364 slepc 857
.seealso: EPSSetProblemType(), EPSProblemType
526 dsic.upv.es!antodo 858
@*/
859
PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)
860
{
861
  PetscFunctionBegin;
862
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1273 slepc 863
  PetscValidPointer(type,2);
526 dsic.upv.es!antodo 864
  *type = eps->problem_type;
865
  PetscFunctionReturn(0);
866
}
867
 
868
#undef __FUNCT__  
1560 slepc 869
#define __FUNCT__ "EPSSetExtraction"
1426 slepc 870
/*@
1560 slepc 871
   EPSSetExtraction - Specifies the type of extraction technique to be employed
1426 slepc 872
   by the eigensolver.
873
 
874
   Collective on EPS
875
 
876
   Input Parameters:
877
+  eps  - the eigensolver context
1560 slepc 878
-  extr - a known type of extraction
1426 slepc 879
 
880
   Options Database Keys:
1560 slepc 881
+  -eps_ritz - Rayleigh-Ritz extraction
882
.  -eps_harmonic - hamonic Ritz extraction
883
.  -eps_refined - refined Ritz extraction
884
-  -eps_refined_harmonic - refined harmonic Ritz extraction
1426 slepc 885
 
886
   Notes:  
1560 slepc 887
   Not all eigensolvers support all types of extraction. See the SLEPc
1426 slepc 888
   Users Manual for details.
889
 
1560 slepc 890
   By default, a standard Rayleigh-Ritz extraction is used. Other extractions
1426 slepc 891
   may be useful when computing interior eigenvalues.
892
 
1560 slepc 893
   Harmonic-type extractions are used in combination with a 'target'.
1426 slepc 894
 
895
   Level: beginner
896
 
1560 slepc 897
.seealso: EPSSetTarget(), EPSGetExtraction(), EPSExtraction
1426 slepc 898
@*/
1560 slepc 899
PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)
1426 slepc 900
{
901
  PetscFunctionBegin;
902
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1560 slepc 903
  eps->extraction = extr;
1426 slepc 904
  PetscFunctionReturn(0);
905
}
906
 
907
#undef __FUNCT__  
1560 slepc 908
#define __FUNCT__ "EPSGetExtraction"
1426 slepc 909
/*@C
1560 slepc 910
   EPSGetExtraction - Gets the extraction type used by the EPS object.
1426 slepc 911
 
912
   Not Collective
913
 
914
   Input Parameter:
915
.  eps - the eigensolver context
916
 
917
   Output Parameter:
1560 slepc 918
.  extr - name of extraction type
1426 slepc 919
 
920
   Level: intermediate
921
 
1560 slepc 922
.seealso: EPSSetExtraction(), EPSExtraction
1426 slepc 923
@*/
1560 slepc 924
PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)
1426 slepc 925
{
926
  PetscFunctionBegin;
927
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1560 slepc 928
  PetscValidPointer(extr,2);
929
  *extr = eps->extraction;
1426 slepc 930
  PetscFunctionReturn(0);
931
}
932
 
933
#undef __FUNCT__  
1799 jroman 934
#define __FUNCT__ "EPSSetBalance"
935
/*@
936
   EPSSetBalance - Specifies the balancing technique to be employed by the
937
   eigensolver, and some parameters associated to it.
938
 
939
   Collective on EPS
940
 
941
   Input Parameters:
942
+  eps    - the eigensolver context
1940 jroman 943
.  bal    - the balancing method, one of EPS_BALANCE_NONE, EPS_BALANCE_ONESIDE,
944
            EPS_BALANCE_TWOSIDE, or EPS_BALANCE_USER
1799 jroman 945
.  its    - number of iterations of the balancing algorithm
946
-  cutoff - cutoff value
947
 
948
   Options Database Keys:
949
+  -eps_balance <method> - the balancing method, where <method> is one of
1804 jroman 950
                           'none', 'oneside', 'twoside', or 'user'
1799 jroman 951
.  -eps_balance_its <its> - number of iterations
952
-  -eps_balance_cutoff <cutoff> - cutoff value
953
 
954
   Notes:
955
   When balancing is enabled, the solver works implicitly with matrix DAD^-1,
956
   where D is an appropriate diagonal matrix. This improves the accuracy of
957
   the computed results in some cases. See the SLEPc Users Manual for details.
958
 
959
   Balancing makes sense only for non-Hermitian problems when the required
960
   precision is high (i.e. a small tolerance such as 1e-15).
961
 
962
   By default, balancing is disabled. The two-sided method is much more
963
   effective than the one-sided counterpart, but it requires the system
964
   matrices to have the MatMultTranspose operation defined.
965
 
966
   The parameter 'its' is the number of iterations performed by the method. The
967
   cutoff value is used only in the two-side variant. Use PETSC_IGNORE for an
968
   argument that need not be changed. Use PETSC_DECIDE to assign a reasonably
969
   good value.
970
 
1804 jroman 971
   User-defined balancing is allowed provided that the corresponding matrix
972
   is set via STSetBalanceMatrix.
973
 
1799 jroman 974
   Level: intermediate
975
 
1804 jroman 976
.seealso: EPSGetBalance(), EPSBalance, STSetBalanceMatrix()
1799 jroman 977
@*/
978
PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)
979
{
980
  PetscFunctionBegin;
981
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
982
  if (bal!=PETSC_IGNORE) {
1940 jroman 983
    if (bal==PETSC_DECIDE || bal==PETSC_DEFAULT) eps->balance = EPS_BALANCE_TWOSIDE;
1810 jroman 984
    else switch (bal) {
1940 jroman 985
      case EPS_BALANCE_NONE:
986
      case EPS_BALANCE_ONESIDE:
987
      case EPS_BALANCE_TWOSIDE:
988
      case EPS_BALANCE_USER:
1810 jroman 989
        eps->balance = bal;
990
        break;
991
      default:
992
        SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
993
    }
1799 jroman 994
  }
995
  if (its!=PETSC_IGNORE) {
996
    if (its==PETSC_DECIDE || its==PETSC_DEFAULT) eps->balance_its = 5;
997
    eps->balance_its = its;
998
  }
999
  if (cutoff!=PETSC_IGNORE) {
1000
    if (cutoff==PETSC_DECIDE || cutoff==PETSC_DEFAULT) eps->balance_cutoff = 1e-8;
1001
    eps->balance_cutoff = cutoff;
1002
  }
1003
  PetscFunctionReturn(0);
1004
}
1005
 
1006
#undef __FUNCT__  
1007
#define __FUNCT__ "EPSGetBalance"
1008
/*@C
1009
   EPSGetBalance - Gets the balancing type used by the EPS object, and the associated
1010
   parameters.
1011
 
1012
   Not Collective
1013
 
1014
   Input Parameter:
1015
.  eps - the eigensolver context
1016
 
1017
   Output Parameters:
1018
+  bal    - the balancing method
1019
.  its    - number of iterations of the balancing algorithm
1020
-  cutoff - cutoff value
1021
 
1022
   Level: intermediate
1023
 
1024
   Note:
1025
   The user can specify PETSC_NULL for any parameter that is not needed.
1026
 
1027
.seealso: EPSSetBalance(), EPSBalance
1028
@*/
1029
PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)
1030
{
1031
  PetscFunctionBegin;
1032
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1033
  PetscValidPointer(bal,2);
1034
  PetscValidPointer(its,3);
1035
  PetscValidPointer(cutoff,4);
1036
  if (bal)    *bal = eps->balance;
1037
  if (its)    *its = eps->balance_its;
1038
  if (cutoff) *cutoff = eps->balance_cutoff;
1039
  PetscFunctionReturn(0);
1040
}
1041
 
1042
#undef __FUNCT__  
2031 jroman 1043
#define __FUNCT__ "EPSSetTrueResidual"
1044
/*@
1045
    EPSSetTrueResidual - Specifies if the solver must compute de true residual
1046
    explicitly or not.
1047
 
1048
    Collective on EPS
1049
 
1050
    Input Parameters:
1051
+   eps     - the eigensolver context
1052
-   trueres - whether true residuals are required or not
1053
 
1054
    Options Database Keys:
1055
.   -eps_true_residual <boolean> - Sets/resets the boolean flag 'trueres'
1056
 
1057
    Notes:
1058
    If the user sets trueres=PETSC_TRUE then the solver explicitly computes
1059
    the true residual for each eigenpair approximation, and uses it for
1060
    convergence testing. Computing the residual is usually an expensive
1061
    operation. Some solvers (e.g., Krylov solvers) can avoid this computation
1062
    by using a cheap estimate of the residual norm, but this may sometimes
1063
    give inaccurate results (especially if a spectral transform is being
1064
    used). On the contrary, preconditioned eigensolvers (e.g., Davidson solvers)
1065
    do rely on computing the true residual, so this option is irrelevant for them.
1066
 
1067
    Level: intermediate
1068
 
1069
.seealso: EPSGetTrueResidual()
1070
@*/
1071
PetscErrorCode EPSSetTrueResidual(EPS eps,PetscTruth trueres)
1072
{
1073
  PetscFunctionBegin;
1074
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1075
  eps->trueres = trueres;
1076
  PetscFunctionReturn(0);
1077
}
1078
 
1079
#undef __FUNCT__  
1080
#define __FUNCT__ "EPSGetTrueResidual"
1081
/*@C
1082
    EPSGetTrueResidual - Returns the flag indicating whether true
1083
    residuals must be computed explicitly or not.
1084
 
1085
    Not Collective
1086
 
1087
    Input Parameter:
1088
.   eps - the eigensolver context
1089
 
1090
    Output Parameter:
1091
.   trueres - the returned flag
1092
 
1093
    Level: intermediate
1094
 
1095
.seealso: EPSSetTrueResidual()
1096
@*/
1097
PetscErrorCode EPSGetTrueResidual(EPS eps,PetscTruth *trueres)
1098
{
1099
  PetscFunctionBegin;
1100
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1101
  PetscValidPointer(trueres,2);
1102
  *trueres = eps->trueres;
1103
  PetscFunctionReturn(0);
1104
}
1105
 
1106
 
1107
#undef __FUNCT__  
526 dsic.upv.es!antodo 1108
#define __FUNCT__ "EPSSetOptionsPrefix"
1109
/*@C
1110
   EPSSetOptionsPrefix - Sets the prefix used for searching for all
1111
   EPS options in the database.
1112
 
1113
   Collective on EPS
1114
 
1115
   Input Parameters:
1116
+  eps - the eigensolver context
1117
-  prefix - the prefix string to prepend to all EPS option requests
1118
 
1119
   Notes:
1120
   A hyphen (-) must NOT be given at the beginning of the prefix name.
1121
   The first character of all runtime options is AUTOMATICALLY the
1122
   hyphen.
1123
 
1124
   For example, to distinguish between the runtime options for two
1125
   different EPS contexts, one could call
1126
.vb
1127
      EPSSetOptionsPrefix(eps1,"eig1_")
1128
      EPSSetOptionsPrefix(eps2,"eig2_")
1129
.ve
1130
 
1131
   Level: advanced
1132
 
1133
.seealso: EPSAppendOptionsPrefix(), EPSGetOptionsPrefix()
1134
@*/
1248 slepc 1135
PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char *prefix)
526 dsic.upv.es!antodo 1136
{
1137
  PetscErrorCode ierr;
1138
  PetscFunctionBegin;
1139
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1140
  ierr = PetscObjectSetOptionsPrefix((PetscObject)eps, prefix);CHKERRQ(ierr);
1141
  ierr = STSetOptionsPrefix(eps->OP,prefix);CHKERRQ(ierr);
1345 slepc 1142
  ierr = IPSetOptionsPrefix(eps->ip,prefix);CHKERRQ(ierr);
1143
  ierr = IPAppendOptionsPrefix(eps->ip,"eps_");CHKERRQ(ierr);
526 dsic.upv.es!antodo 1144
  PetscFunctionReturn(0);  
1145
}
1146
 
1147
#undef __FUNCT__  
1148
#define __FUNCT__ "EPSAppendOptionsPrefix"
1149
/*@C
1150
   EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
1151
   EPS options in the database.
1152
 
1153
   Collective on EPS
1154
 
1155
   Input Parameters:
1156
+  eps - the eigensolver context
1157
-  prefix - the prefix string to prepend to all EPS option requests
1158
 
1159
   Notes:
1160
   A hyphen (-) must NOT be given at the beginning of the prefix name.
1161
   The first character of all runtime options is AUTOMATICALLY the hyphen.
1162
 
1163
   Level: advanced
1164
 
1165
.seealso: EPSSetOptionsPrefix(), EPSGetOptionsPrefix()
1166
@*/
1248 slepc 1167
PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char *prefix)
526 dsic.upv.es!antodo 1168
{
1169
  PetscErrorCode ierr;
1170
  PetscFunctionBegin;
1171
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1172
  ierr = PetscObjectAppendOptionsPrefix((PetscObject)eps, prefix);CHKERRQ(ierr);
1173
  ierr = STAppendOptionsPrefix(eps->OP,prefix); CHKERRQ(ierr);
1345 slepc 1174
  ierr = IPSetOptionsPrefix(eps->ip,prefix);CHKERRQ(ierr);
1175
  ierr = IPAppendOptionsPrefix(eps->ip,"eps_");CHKERRQ(ierr);
526 dsic.upv.es!antodo 1176
  PetscFunctionReturn(0);
1177
}
1178
 
1179
#undef __FUNCT__  
1180
#define __FUNCT__ "EPSGetOptionsPrefix"
1181
/*@C
1182
   EPSGetOptionsPrefix - Gets the prefix used for searching for all
1183
   EPS options in the database.
1184
 
1185
   Not Collective
1186
 
1187
   Input Parameters:
1188
.  eps - the eigensolver context
1189
 
1190
   Output Parameters:
1191
.  prefix - pointer to the prefix string used is returned
1192
 
1193
   Notes: On the fortran side, the user should pass in a string 'prefix' of
1194
   sufficient length to hold the prefix.
1195
 
1196
   Level: advanced
1197
 
1198
.seealso: EPSSetOptionsPrefix(), EPSAppendOptionsPrefix()
1199
@*/
812 dsic.upv.es!antodo 1200
PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])
526 dsic.upv.es!antodo 1201
{
1202
  PetscErrorCode ierr;
1203
  PetscFunctionBegin;
1204
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
1273 slepc 1205
  PetscValidPointer(prefix,2);
526 dsic.upv.es!antodo 1206
  ierr = PetscObjectGetOptionsPrefix((PetscObject)eps, prefix);CHKERRQ(ierr);
1207
  PetscFunctionReturn(0);
1208
}