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