Actual source code: nepopts.c

slepc-3.17.2 2022-08-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    NEP routines related to options that can be set via the command-line
 12:    or procedurally
 13: */

 15: #include <slepc/private/nepimpl.h>
 16: #include <petscdraw.h>

 18: /*@C
 19:    NEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 20:    indicated by the user.

 22:    Collective on nep

 24:    Input Parameters:
 25: +  nep      - the nonlinear eigensolver context
 26: .  opt      - the command line option for this monitor
 27: .  name     - the monitor type one is seeking
 28: .  ctx      - an optional user context for the monitor, or NULL
 29: -  trackall - whether this monitor tracks all eigenvalues or not

 31:    Level: developer

 33: .seealso: NEPMonitorSet(), NEPSetTrackAll()
 34: @*/
 35: PetscErrorCode NEPMonitorSetFromOptions(NEP nep,const char opt[],const char name[],void *ctx,PetscBool trackall)
 36: {
 37:   PetscErrorCode       (*mfunc)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 38:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
 39:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
 40:   PetscViewerAndFormat *vf;
 41:   PetscViewer          viewer;
 42:   PetscViewerFormat    format;
 43:   PetscViewerType      vtype;
 44:   char                 key[PETSC_MAX_PATH_LEN];
 45:   PetscBool            flg;

 47:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,opt,&viewer,&format,&flg);
 48:   if (!flg) PetscFunctionReturn(0);

 50:   PetscViewerGetType(viewer,&vtype);
 51:   SlepcMonitorMakeKey_Internal(name,vtype,format,key);
 52:   PetscFunctionListFind(NEPMonitorList,key,&mfunc);
 54:   PetscFunctionListFind(NEPMonitorCreateList,key,&cfunc);
 55:   PetscFunctionListFind(NEPMonitorDestroyList,key,&dfunc);
 56:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 57:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

 59:   (*cfunc)(viewer,format,ctx,&vf);
 60:   PetscObjectDereference((PetscObject)viewer);
 61:   NEPMonitorSet(nep,mfunc,vf,(PetscErrorCode(*)(void **))dfunc);
 62:   if (trackall) NEPSetTrackAll(nep,PETSC_TRUE);
 63:   PetscFunctionReturn(0);
 64: }

 66: /*@
 67:    NEPSetFromOptions - Sets NEP options from the options database.
 68:    This routine must be called before NEPSetUp() if the user is to be
 69:    allowed to set the solver type.

 71:    Collective on nep

 73:    Input Parameters:
 74: .  nep - the nonlinear eigensolver context

 76:    Notes:
 77:    To see all options, run your program with the -help option.

 79:    Level: beginner

 81: .seealso: NEPSetOptionsPrefix()
 82: @*/
 83: PetscErrorCode NEPSetFromOptions(NEP nep)
 84: {
 85:   PetscErrorCode  ierr;
 86:   char            type[256];
 87:   PetscBool       set,flg,flg1,flg2,flg3,flg4,flg5,bval;
 88:   PetscReal       r;
 89:   PetscScalar     s;
 90:   PetscInt        i,j,k;
 91:   NEPRefine       refine;
 92:   NEPRefineScheme scheme;

 95:   NEPRegisterAll();
 96:   ierr = PetscObjectOptionsBegin((PetscObject)nep);
 97:     PetscOptionsFList("-nep_type","Nonlinear eigensolver method","NEPSetType",NEPList,(char*)(((PetscObject)nep)->type_name?((PetscObject)nep)->type_name:NEPRII),type,sizeof(type),&flg);
 98:     if (flg) NEPSetType(nep,type);
 99:     else if (!((PetscObject)nep)->type_name) NEPSetType(nep,NEPRII);

101:     PetscOptionsBoolGroupBegin("-nep_general","General nonlinear eigenvalue problem","NEPSetProblemType",&flg);
102:     if (flg) NEPSetProblemType(nep,NEP_GENERAL);
103:     PetscOptionsBoolGroupEnd("-nep_rational","Rational eigenvalue problem","NEPSetProblemType",&flg);
104:     if (flg) NEPSetProblemType(nep,NEP_RATIONAL);

106:     refine = nep->refine;
107:     PetscOptionsEnum("-nep_refine","Iterative refinement method","NEPSetRefine",NEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1);
108:     i = nep->npart;
109:     PetscOptionsInt("-nep_refine_partitions","Number of partitions of the communicator for iterative refinement","NEPSetRefine",nep->npart,&i,&flg2);
110:     r = nep->rtol;
111:     PetscOptionsReal("-nep_refine_tol","Tolerance for iterative refinement","NEPSetRefine",nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:nep->rtol,&r,&flg3);
112:     j = nep->rits;
113:     PetscOptionsInt("-nep_refine_its","Maximum number of iterations for iterative refinement","NEPSetRefine",nep->rits,&j,&flg4);
114:     scheme = nep->scheme;
115:     PetscOptionsEnum("-nep_refine_scheme","Scheme used for linear systems within iterative refinement","NEPSetRefine",NEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5);
116:     if (flg1 || flg2 || flg3 || flg4 || flg5) NEPSetRefine(nep,refine,i,r,j,scheme);

118:     i = nep->max_it;
119:     PetscOptionsInt("-nep_max_it","Maximum number of iterations","NEPSetTolerances",nep->max_it,&i,&flg1);
120:     r = nep->tol;
121:     PetscOptionsReal("-nep_tol","Tolerance","NEPSetTolerances",SlepcDefaultTol(nep->tol),&r,&flg2);
122:     if (flg1 || flg2) NEPSetTolerances(nep,r,i);

124:     PetscOptionsBoolGroupBegin("-nep_conv_rel","Relative error convergence test","NEPSetConvergenceTest",&flg);
125:     if (flg) NEPSetConvergenceTest(nep,NEP_CONV_REL);
126:     PetscOptionsBoolGroup("-nep_conv_norm","Convergence test relative to the matrix norms","NEPSetConvergenceTest",&flg);
127:     if (flg) NEPSetConvergenceTest(nep,NEP_CONV_NORM);
128:     PetscOptionsBoolGroup("-nep_conv_abs","Absolute error convergence test","NEPSetConvergenceTest",&flg);
129:     if (flg) NEPSetConvergenceTest(nep,NEP_CONV_ABS);
130:     PetscOptionsBoolGroupEnd("-nep_conv_user","User-defined convergence test","NEPSetConvergenceTest",&flg);
131:     if (flg) NEPSetConvergenceTest(nep,NEP_CONV_USER);

133:     PetscOptionsBoolGroupBegin("-nep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","NEPSetStoppingTest",&flg);
134:     if (flg) NEPSetStoppingTest(nep,NEP_STOP_BASIC);
135:     PetscOptionsBoolGroupEnd("-nep_stop_user","User-defined stopping test","NEPSetStoppingTest",&flg);
136:     if (flg) NEPSetStoppingTest(nep,NEP_STOP_USER);

138:     i = nep->nev;
139:     PetscOptionsInt("-nep_nev","Number of eigenvalues to compute","NEPSetDimensions",nep->nev,&i,&flg1);
140:     j = nep->ncv;
141:     PetscOptionsInt("-nep_ncv","Number of basis vectors","NEPSetDimensions",nep->ncv,&j,&flg2);
142:     k = nep->mpd;
143:     PetscOptionsInt("-nep_mpd","Maximum dimension of projected problem","NEPSetDimensions",nep->mpd,&k,&flg3);
144:     if (flg1 || flg2 || flg3) NEPSetDimensions(nep,i,j,k);

146:     PetscOptionsBoolGroupBegin("-nep_largest_magnitude","Compute largest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
147:     if (flg) NEPSetWhichEigenpairs(nep,NEP_LARGEST_MAGNITUDE);
148:     PetscOptionsBoolGroup("-nep_smallest_magnitude","Compute smallest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
149:     if (flg) NEPSetWhichEigenpairs(nep,NEP_SMALLEST_MAGNITUDE);
150:     PetscOptionsBoolGroup("-nep_largest_real","Compute eigenvalues with largest real parts","NEPSetWhichEigenpairs",&flg);
151:     if (flg) NEPSetWhichEigenpairs(nep,NEP_LARGEST_REAL);
152:     PetscOptionsBoolGroup("-nep_smallest_real","Compute eigenvalues with smallest real parts","NEPSetWhichEigenpairs",&flg);
153:     if (flg) NEPSetWhichEigenpairs(nep,NEP_SMALLEST_REAL);
154:     PetscOptionsBoolGroup("-nep_largest_imaginary","Compute eigenvalues with largest imaginary parts","NEPSetWhichEigenpairs",&flg);
155:     if (flg) NEPSetWhichEigenpairs(nep,NEP_LARGEST_IMAGINARY);
156:     PetscOptionsBoolGroup("-nep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","NEPSetWhichEigenpairs",&flg);
157:     if (flg) NEPSetWhichEigenpairs(nep,NEP_SMALLEST_IMAGINARY);
158:     PetscOptionsBoolGroup("-nep_target_magnitude","Compute eigenvalues closest to target","NEPSetWhichEigenpairs",&flg);
159:     if (flg) NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
160:     PetscOptionsBoolGroup("-nep_target_real","Compute eigenvalues with real parts closest to target","NEPSetWhichEigenpairs",&flg);
161:     if (flg) NEPSetWhichEigenpairs(nep,NEP_TARGET_REAL);
162:     PetscOptionsBoolGroup("-nep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","NEPSetWhichEigenpairs",&flg);
163:     if (flg) NEPSetWhichEigenpairs(nep,NEP_TARGET_IMAGINARY);
164:     PetscOptionsBoolGroupEnd("-nep_all","Compute all eigenvalues in a region","NEPSetWhichEigenpairs",&flg);
165:     if (flg) NEPSetWhichEigenpairs(nep,NEP_ALL);

167:     PetscOptionsScalar("-nep_target","Value of the target","NEPSetTarget",nep->target,&s,&flg);
168:     if (flg) {
169:       if (nep->which!=NEP_TARGET_REAL && nep->which!=NEP_TARGET_IMAGINARY) NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
170:       NEPSetTarget(nep,s);
171:     }

173:     PetscOptionsBool("-nep_two_sided","Use two-sided variant (to compute left eigenvectors)","NEPSetTwoSided",nep->twosided,&bval,&flg);
174:     if (flg) NEPSetTwoSided(nep,bval);

176:     /* -----------------------------------------------------------------------*/
177:     /*
178:       Cancels all monitors hardwired into code before call to NEPSetFromOptions()
179:     */
180:     PetscOptionsBool("-nep_monitor_cancel","Remove any hardwired monitor routines","NEPMonitorCancel",PETSC_FALSE,&flg,&set);
181:     if (set && flg) NEPMonitorCancel(nep);
182:     NEPMonitorSetFromOptions(nep,"-nep_monitor","first_approximation",NULL,PETSC_FALSE);
183:     NEPMonitorSetFromOptions(nep,"-nep_monitor_all","all_approximations",NULL,PETSC_TRUE);
184:     NEPMonitorSetFromOptions(nep,"-nep_monitor_conv","convergence_history",NULL,PETSC_FALSE);

186:     /* -----------------------------------------------------------------------*/
187:     PetscOptionsName("-nep_view","Print detailed information on solver used","NEPView",NULL);
188:     PetscOptionsName("-nep_view_vectors","View computed eigenvectors","NEPVectorsView",NULL);
189:     PetscOptionsName("-nep_view_values","View computed eigenvalues","NEPValuesView",NULL);
190:     PetscOptionsName("-nep_converged_reason","Print reason for convergence, and number of iterations","NEPConvergedReasonView",NULL);
191:     PetscOptionsName("-nep_error_absolute","Print absolute errors of each eigenpair","NEPErrorView",NULL);
192:     PetscOptionsName("-nep_error_relative","Print relative errors of each eigenpair","NEPErrorView",NULL);

194:     if (nep->ops->setfromoptions) (*nep->ops->setfromoptions)(PetscOptionsObject,nep);
195:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)nep);
196:   ierr = PetscOptionsEnd();

198:   if (!nep->V) NEPGetBV(nep,&nep->V);
199:   BVSetFromOptions(nep->V);
200:   if (!nep->rg) NEPGetRG(nep,&nep->rg);
201:   RGSetFromOptions(nep->rg);
202:   if (nep->useds) {
203:     if (!nep->ds) NEPGetDS(nep,&nep->ds);
204:     DSSetFromOptions(nep->ds);
205:   }
206:   if (!nep->refineksp) NEPRefineGetKSP(nep,&nep->refineksp);
207:   KSPSetFromOptions(nep->refineksp);
208:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) for (i=0;i<nep->nt;i++) FNSetFromOptions(nep->f[i]);
209:   PetscFunctionReturn(0);
210: }

212: /*@C
213:    NEPGetTolerances - Gets the tolerance and maximum iteration count used
214:    by the NEP convergence tests.

216:    Not Collective

218:    Input Parameter:
219: .  nep - the nonlinear eigensolver context

221:    Output Parameters:
222: +  tol - the convergence tolerance
223: -  maxits - maximum number of iterations

225:    Notes:
226:    The user can specify NULL for any parameter that is not needed.

228:    Level: intermediate

230: .seealso: NEPSetTolerances()
231: @*/
232: PetscErrorCode NEPGetTolerances(NEP nep,PetscReal *tol,PetscInt *maxits)
233: {
235:   if (tol)    *tol    = nep->tol;
236:   if (maxits) *maxits = nep->max_it;
237:   PetscFunctionReturn(0);
238: }

240: /*@
241:    NEPSetTolerances - Sets the tolerance and maximum iteration count used
242:    by the NEP convergence tests.

244:    Logically Collective on nep

246:    Input Parameters:
247: +  nep    - the nonlinear eigensolver context
248: .  tol    - the convergence tolerance
249: -  maxits - maximum number of iterations to use

251:    Options Database Keys:
252: +  -nep_tol <tol> - Sets the convergence tolerance
253: -  -nep_max_it <maxits> - Sets the maximum number of iterations allowed

255:    Notes:
256:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

258:    Level: intermediate

260: .seealso: NEPGetTolerances()
261: @*/
262: PetscErrorCode NEPSetTolerances(NEP nep,PetscReal tol,PetscInt maxits)
263: {
267:   if (tol == PETSC_DEFAULT) {
268:     nep->tol   = PETSC_DEFAULT;
269:     nep->state = NEP_STATE_INITIAL;
270:   } else {
272:     nep->tol = tol;
273:   }
274:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
275:     nep->max_it = PETSC_DEFAULT;
276:     nep->state  = NEP_STATE_INITIAL;
277:   } else {
279:     nep->max_it = maxits;
280:   }
281:   PetscFunctionReturn(0);
282: }

284: /*@C
285:    NEPGetDimensions - Gets the number of eigenvalues to compute
286:    and the dimension of the subspace.

288:    Not Collective

290:    Input Parameter:
291: .  nep - the nonlinear eigensolver context

293:    Output Parameters:
294: +  nev - number of eigenvalues to compute
295: .  ncv - the maximum dimension of the subspace to be used by the solver
296: -  mpd - the maximum dimension allowed for the projected problem

298:    Notes:
299:    The user can specify NULL for any parameter that is not needed.

301:    Level: intermediate

303: .seealso: NEPSetDimensions()
304: @*/
305: PetscErrorCode NEPGetDimensions(NEP nep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
306: {
308:   if (nev) *nev = nep->nev;
309:   if (ncv) *ncv = nep->ncv;
310:   if (mpd) *mpd = nep->mpd;
311:   PetscFunctionReturn(0);
312: }

314: /*@
315:    NEPSetDimensions - Sets the number of eigenvalues to compute
316:    and the dimension of the subspace.

318:    Logically Collective on nep

320:    Input Parameters:
321: +  nep - the nonlinear eigensolver context
322: .  nev - number of eigenvalues to compute
323: .  ncv - the maximum dimension of the subspace to be used by the solver
324: -  mpd - the maximum dimension allowed for the projected problem

326:    Options Database Keys:
327: +  -nep_nev <nev> - Sets the number of eigenvalues
328: .  -nep_ncv <ncv> - Sets the dimension of the subspace
329: -  -nep_mpd <mpd> - Sets the maximum projected dimension

331:    Notes:
332:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
333:    dependent on the solution method.

335:    The parameters ncv and mpd are intimately related, so that the user is advised
336:    to set one of them at most. Normal usage is that
337:    (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
338:    (b) in cases where nev is large, the user sets mpd.

340:    The value of ncv should always be between nev and (nev+mpd), typically
341:    ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
342:    a smaller value should be used.

344:    Level: intermediate

346: .seealso: NEPGetDimensions()
347: @*/
348: PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd)
349: {
355:   nep->nev = nev;
356:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
357:     nep->ncv = PETSC_DEFAULT;
358:   } else {
360:     nep->ncv = ncv;
361:   }
362:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
363:     nep->mpd = PETSC_DEFAULT;
364:   } else {
366:     nep->mpd = mpd;
367:   }
368:   nep->state = NEP_STATE_INITIAL;
369:   PetscFunctionReturn(0);
370: }

372: /*@
373:     NEPSetWhichEigenpairs - Specifies which portion of the spectrum is
374:     to be sought.

376:     Logically Collective on nep

378:     Input Parameters:
379: +   nep   - eigensolver context obtained from NEPCreate()
380: -   which - the portion of the spectrum to be sought

382:     Possible values:
383:     The parameter 'which' can have one of these values

385: +     NEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
386: .     NEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
387: .     NEP_LARGEST_REAL - largest real parts
388: .     NEP_SMALLEST_REAL - smallest real parts
389: .     NEP_LARGEST_IMAGINARY - largest imaginary parts
390: .     NEP_SMALLEST_IMAGINARY - smallest imaginary parts
391: .     NEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
392: .     NEP_TARGET_REAL - eigenvalues with real part closest to target
393: .     NEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
394: .     NEP_ALL - all eigenvalues contained in a given region
395: -     NEP_WHICH_USER - user defined ordering set with NEPSetEigenvalueComparison()

397:     Options Database Keys:
398: +   -nep_largest_magnitude - Sets largest eigenvalues in magnitude
399: .   -nep_smallest_magnitude - Sets smallest eigenvalues in magnitude
400: .   -nep_largest_real - Sets largest real parts
401: .   -nep_smallest_real - Sets smallest real parts
402: .   -nep_largest_imaginary - Sets largest imaginary parts
403: .   -nep_smallest_imaginary - Sets smallest imaginary parts
404: .   -nep_target_magnitude - Sets eigenvalues closest to target
405: .   -nep_target_real - Sets real parts closest to target
406: .   -nep_target_imaginary - Sets imaginary parts closest to target
407: -   -nep_all - Sets all eigenvalues in a region

409:     Notes:
410:     Not all eigensolvers implemented in NEP account for all the possible values
411:     stated above. If SLEPc is compiled for real numbers NEP_LARGEST_IMAGINARY
412:     and NEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
413:     for eigenvalue selection.

415:     The target is a scalar value provided with NEPSetTarget().

417:     NEP_ALL is intended for use in the context of the CISS solver for
418:     computing all eigenvalues in a region.

420:     Level: intermediate

422: .seealso: NEPGetWhichEigenpairs(), NEPSetTarget(), NEPSetEigenvalueComparison(), NEPWhich
423: @*/
424: PetscErrorCode NEPSetWhichEigenpairs(NEP nep,NEPWhich which)
425: {
428:   switch (which) {
429:     case NEP_LARGEST_MAGNITUDE:
430:     case NEP_SMALLEST_MAGNITUDE:
431:     case NEP_LARGEST_REAL:
432:     case NEP_SMALLEST_REAL:
433:     case NEP_LARGEST_IMAGINARY:
434:     case NEP_SMALLEST_IMAGINARY:
435:     case NEP_TARGET_MAGNITUDE:
436:     case NEP_TARGET_REAL:
437: #if defined(PETSC_USE_COMPLEX)
438:     case NEP_TARGET_IMAGINARY:
439: #endif
440:     case NEP_ALL:
441:     case NEP_WHICH_USER:
442:       if (nep->which != which) {
443:         nep->state = NEP_STATE_INITIAL;
444:         nep->which = which;
445:       }
446:       break;
447: #if !defined(PETSC_USE_COMPLEX)
448:     case NEP_TARGET_IMAGINARY:
449:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEP_TARGET_IMAGINARY can be used only with complex scalars");
450: #endif
451:     default:
452:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
453:   }
454:   PetscFunctionReturn(0);
455: }

457: /*@
458:     NEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
459:     sought.

461:     Not Collective

463:     Input Parameter:
464: .   nep - eigensolver context obtained from NEPCreate()

466:     Output Parameter:
467: .   which - the portion of the spectrum to be sought

469:     Notes:
470:     See NEPSetWhichEigenpairs() for possible values of 'which'.

472:     Level: intermediate

474: .seealso: NEPSetWhichEigenpairs(), NEPWhich
475: @*/
476: PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which)
477: {
480:   *which = nep->which;
481:   PetscFunctionReturn(0);
482: }

484: /*@C
485:    NEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
486:    when NEPSetWhichEigenpairs() is set to NEP_WHICH_USER.

488:    Logically Collective on nep

490:    Input Parameters:
491: +  nep  - eigensolver context obtained from NEPCreate()
492: .  func - a pointer to the comparison function
493: -  ctx  - a context pointer (the last parameter to the comparison function)

495:    Calling Sequence of func:
496: $   func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)

498: +   ar     - real part of the 1st eigenvalue
499: .   ai     - imaginary part of the 1st eigenvalue
500: .   br     - real part of the 2nd eigenvalue
501: .   bi     - imaginary part of the 2nd eigenvalue
502: .   res    - result of comparison
503: -   ctx    - optional context, as set by NEPSetEigenvalueComparison()

505:    Note:
506:    The returning parameter 'res' can be
507: +  negative - if the 1st eigenvalue is preferred to the 2st one
508: .  zero     - if both eigenvalues are equally preferred
509: -  positive - if the 2st eigenvalue is preferred to the 1st one

511:    Level: advanced

513: .seealso: NEPSetWhichEigenpairs(), NEPWhich
514: @*/
515: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
516: {
518:   nep->sc->comparison    = func;
519:   nep->sc->comparisonctx = ctx;
520:   nep->which             = NEP_WHICH_USER;
521:   PetscFunctionReturn(0);
522: }

524: /*@
525:    NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.

527:    Logically Collective on nep

529:    Input Parameters:
530: +  nep  - the nonlinear eigensolver context
531: -  type - a known type of nonlinear eigenvalue problem

533:    Options Database Keys:
534: +  -nep_general - general problem with no particular structure
535: -  -nep_rational - a rational eigenvalue problem defined in split form with all f_i rational

537:    Notes:
538:    Allowed values for the problem type are general (NEP_GENERAL), and rational
539:    (NEP_RATIONAL).

541:    This function is used to provide a hint to the NEP solver to exploit certain
542:    properties of the nonlinear eigenproblem. This hint may be used or not,
543:    depending on the solver. By default, no particular structure is assumed.

545:    Level: intermediate

547: .seealso: NEPSetType(), NEPGetProblemType(), NEPProblemType
548: @*/
549: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)
550: {
554:   if (type != nep->problem_type) {
555:     nep->problem_type = type;
556:     nep->state = NEP_STATE_INITIAL;
557:   }
558:   PetscFunctionReturn(0);
559: }

561: /*@
562:    NEPGetProblemType - Gets the problem type from the NEP object.

564:    Not Collective

566:    Input Parameter:
567: .  nep - the nonlinear eigensolver context

569:    Output Parameter:
570: .  type - the problem type

572:    Level: intermediate

574: .seealso: NEPSetProblemType(), NEPProblemType
575: @*/
576: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)
577: {
580:   *type = nep->problem_type;
581:   PetscFunctionReturn(0);
582: }

584: /*@
585:    NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
586:    eigenvectors are also computed.

588:    Logically Collective on nep

590:    Input Parameters:
591: +  nep      - the eigensolver context
592: -  twosided - whether the two-sided variant is to be used or not

594:    Options Database Keys:
595: .  -nep_two_sided <boolean> - Sets/resets the twosided flag

597:    Notes:
598:    If the user sets twosided=PETSC_TRUE then the solver uses a variant of
599:    the algorithm that computes both right and left eigenvectors. This is
600:    usually much more costly. This option is not available in all solvers.

602:    When using two-sided solvers, the problem matrices must have both the
603:    MatMult and MatMultTranspose operations defined.

605:    Level: advanced

607: .seealso: NEPGetTwoSided(), NEPGetLeftEigenvector()
608: @*/
609: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)
610: {
613:   if (twosided!=nep->twosided) {
614:     nep->twosided = twosided;
615:     nep->state    = NEP_STATE_INITIAL;
616:   }
617:   PetscFunctionReturn(0);
618: }

620: /*@
621:    NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
622:    of the algorithm is being used or not.

624:    Not Collective

626:    Input Parameter:
627: .  nep - the eigensolver context

629:    Output Parameter:
630: .  twosided - the returned flag

632:    Level: advanced

634: .seealso: NEPSetTwoSided()
635: @*/
636: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)
637: {
640:   *twosided = nep->twosided;
641:   PetscFunctionReturn(0);
642: }

644: /*@C
645:    NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
646:    used in the convergence test.

648:    Logically Collective on nep

650:    Input Parameters:
651: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
652: .  func    - a pointer to the convergence test function
653: .  ctx     - context for private data for the convergence routine (may be null)
654: -  destroy - a routine for destroying the context (may be null)

656:    Calling Sequence of func:
657: $   func(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)

659: +   nep    - nonlinear eigensolver context obtained from NEPCreate()
660: .   eigr   - real part of the eigenvalue
661: .   eigi   - imaginary part of the eigenvalue
662: .   res    - residual norm associated to the eigenpair
663: .   errest - (output) computed error estimate
664: -   ctx    - optional context, as set by NEPSetConvergenceTestFunction()

666:    Note:
667:    If the error estimate returned by the convergence test function is less than
668:    the tolerance, then the eigenvalue is accepted as converged.

670:    Level: advanced

672: .seealso: NEPSetConvergenceTest(), NEPSetTolerances()
673: @*/
674: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
675: {
677:   if (nep->convergeddestroy) (*nep->convergeddestroy)(nep->convergedctx);
678:   nep->convergeduser    = func;
679:   nep->convergeddestroy = destroy;
680:   nep->convergedctx     = ctx;
681:   if (func == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
682:   else if (func == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
683:   else if (func == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
684:   else {
685:     nep->conv      = NEP_CONV_USER;
686:     nep->converged = nep->convergeduser;
687:   }
688:   PetscFunctionReturn(0);
689: }

691: /*@
692:    NEPSetConvergenceTest - Specifies how to compute the error estimate
693:    used in the convergence test.

695:    Logically Collective on nep

697:    Input Parameters:
698: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
699: -  conv - the type of convergence test

701:    Options Database Keys:
702: +  -nep_conv_abs  - Sets the absolute convergence test
703: .  -nep_conv_rel  - Sets the convergence test relative to the eigenvalue
704: -  -nep_conv_user - Selects the user-defined convergence test

706:    Note:
707:    The parameter 'conv' can have one of these values
708: +     NEP_CONV_ABS  - absolute error ||r||
709: .     NEP_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
710: .     NEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(|f_i(l)|*||A_i||)
711: -     NEP_CONV_USER - function set by NEPSetConvergenceTestFunction()

713:    Level: intermediate

715: .seealso: NEPGetConvergenceTest(), NEPSetConvergenceTestFunction(), NEPSetStoppingTest(), NEPConv
716: @*/
717: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)
718: {
721:   switch (conv) {
722:     case NEP_CONV_ABS:  nep->converged = NEPConvergedAbsolute; break;
723:     case NEP_CONV_REL:  nep->converged = NEPConvergedRelative; break;
724:     case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
725:     case NEP_CONV_USER:
727:       nep->converged = nep->convergeduser;
728:       break;
729:     default:
730:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
731:   }
732:   nep->conv = conv;
733:   PetscFunctionReturn(0);
734: }

736: /*@
737:    NEPGetConvergenceTest - Gets the method used to compute the error estimate
738:    used in the convergence test.

740:    Not Collective

742:    Input Parameters:
743: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

745:    Output Parameters:
746: .  conv  - the type of convergence test

748:    Level: intermediate

750: .seealso: NEPSetConvergenceTest(), NEPConv
751: @*/
752: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)
753: {
756:   *conv = nep->conv;
757:   PetscFunctionReturn(0);
758: }

760: /*@C
761:    NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
762:    iteration of the eigensolver.

764:    Logically Collective on nep

766:    Input Parameters:
767: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
768: .  func    - pointer to the stopping test function
769: .  ctx     - context for private data for the stopping routine (may be null)
770: -  destroy - a routine for destroying the context (may be null)

772:    Calling Sequence of func:
773: $   func(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx)

775: +   nep    - nonlinear eigensolver context obtained from NEPCreate()
776: .   its    - current number of iterations
777: .   max_it - maximum number of iterations
778: .   nconv  - number of currently converged eigenpairs
779: .   nev    - number of requested eigenpairs
780: .   reason - (output) result of the stopping test
781: -   ctx    - optional context, as set by NEPSetStoppingTestFunction()

783:    Note:
784:    Normal usage is to first call the default routine NEPStoppingBasic() and then
785:    set reason to NEP_CONVERGED_USER if some user-defined conditions have been
786:    met. To let the eigensolver continue iterating, the result must be left as
787:    NEP_CONVERGED_ITERATING.

789:    Level: advanced

791: .seealso: NEPSetStoppingTest(), NEPStoppingBasic()
792: @*/
793: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
794: {
796:   if (nep->stoppingdestroy) (*nep->stoppingdestroy)(nep->stoppingctx);
797:   nep->stoppinguser    = func;
798:   nep->stoppingdestroy = destroy;
799:   nep->stoppingctx     = ctx;
800:   if (func == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
801:   else {
802:     nep->stop     = NEP_STOP_USER;
803:     nep->stopping = nep->stoppinguser;
804:   }
805:   PetscFunctionReturn(0);
806: }

808: /*@
809:    NEPSetStoppingTest - Specifies how to decide the termination of the outer
810:    loop of the eigensolver.

812:    Logically Collective on nep

814:    Input Parameters:
815: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
816: -  stop - the type of stopping test

818:    Options Database Keys:
819: +  -nep_stop_basic - Sets the default stopping test
820: -  -nep_stop_user  - Selects the user-defined stopping test

822:    Note:
823:    The parameter 'stop' can have one of these values
824: +     NEP_STOP_BASIC - default stopping test
825: -     NEP_STOP_USER  - function set by NEPSetStoppingTestFunction()

827:    Level: advanced

829: .seealso: NEPGetStoppingTest(), NEPSetStoppingTestFunction(), NEPSetConvergenceTest(), NEPStop
830: @*/
831: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)
832: {
835:   switch (stop) {
836:     case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
837:     case NEP_STOP_USER:
839:       nep->stopping = nep->stoppinguser;
840:       break;
841:     default:
842:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
843:   }
844:   nep->stop = stop;
845:   PetscFunctionReturn(0);
846: }

848: /*@
849:    NEPGetStoppingTest - Gets the method used to decide the termination of the outer
850:    loop of the eigensolver.

852:    Not Collective

854:    Input Parameters:
855: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

857:    Output Parameters:
858: .  stop  - the type of stopping test

860:    Level: advanced

862: .seealso: NEPSetStoppingTest(), NEPStop
863: @*/
864: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)
865: {
868:   *stop = nep->stop;
869:   PetscFunctionReturn(0);
870: }

872: /*@
873:    NEPSetTrackAll - Specifies if the solver must compute the residual of all
874:    approximate eigenpairs or not.

876:    Logically Collective on nep

878:    Input Parameters:
879: +  nep      - the eigensolver context
880: -  trackall - whether compute all residuals or not

882:    Notes:
883:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
884:    the residual for each eigenpair approximation. Computing the residual is
885:    usually an expensive operation and solvers commonly compute the associated
886:    residual to the first unconverged eigenpair.

888:    The option '-nep_monitor_all' automatically activates this option.

890:    Level: developer

892: .seealso: NEPGetTrackAll()
893: @*/
894: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)
895: {
898:   nep->trackall = trackall;
899:   PetscFunctionReturn(0);
900: }

902: /*@
903:    NEPGetTrackAll - Returns the flag indicating whether all residual norms must
904:    be computed or not.

906:    Not Collective

908:    Input Parameter:
909: .  nep - the eigensolver context

911:    Output Parameter:
912: .  trackall - the returned flag

914:    Level: developer

916: .seealso: NEPSetTrackAll()
917: @*/
918: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)
919: {
922:   *trackall = nep->trackall;
923:   PetscFunctionReturn(0);
924: }

926: /*@
927:    NEPSetRefine - Specifies the refinement type (and options) to be used
928:    after the solve.

930:    Logically Collective on nep

932:    Input Parameters:
933: +  nep    - the nonlinear eigensolver context
934: .  refine - refinement type
935: .  npart  - number of partitions of the communicator
936: .  tol    - the convergence tolerance
937: .  its    - maximum number of refinement iterations
938: -  scheme - which scheme to be used for solving the involved linear systems

940:    Options Database Keys:
941: +  -nep_refine <type> - refinement type, one of <none,simple,multiple>
942: .  -nep_refine_partitions <n> - the number of partitions
943: .  -nep_refine_tol <tol> - the tolerance
944: .  -nep_refine_its <its> - number of iterations
945: -  -nep_refine_scheme - to set the scheme for the linear solves

947:    Notes:
948:    By default, iterative refinement is disabled, since it may be very
949:    costly. There are two possible refinement strategies, simple and multiple.
950:    The simple approach performs iterative refinement on each of the
951:    converged eigenpairs individually, whereas the multiple strategy works
952:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
953:    The latter may be required for the case of multiple eigenvalues.

955:    In some cases, especially when using direct solvers within the
956:    iterative refinement method, it may be helpful for improved scalability
957:    to split the communicator in several partitions. The npart parameter
958:    indicates how many partitions to use (defaults to 1).

960:    The tol and its parameters specify the stopping criterion. In the simple
961:    method, refinement continues until the residual of each eigenpair is
962:    below the tolerance (tol defaults to the NEP tol, but may be set to a
963:    different value). In contrast, the multiple method simply performs its
964:    refinement iterations (just one by default).

966:    The scheme argument is used to change the way in which linear systems are
967:    solved. Possible choices are explicit, mixed block elimination (MBE),
968:    and Schur complement.

970:    Level: intermediate

972: .seealso: NEPGetRefine()
973: @*/
974: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
975: {
976:   PetscMPIInt    size;

984:   nep->refine = refine;
985:   if (refine) {  /* process parameters only if not REFINE_NONE */
986:     if (npart!=nep->npart) {
987:       PetscSubcommDestroy(&nep->refinesubc);
988:       KSPDestroy(&nep->refineksp);
989:     }
990:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
991:       nep->npart = 1;
992:     } else {
993:       MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size);
995:       nep->npart = npart;
996:     }
997:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
998:       nep->rtol = PETSC_DEFAULT;
999:     } else {
1001:       nep->rtol = tol;
1002:     }
1003:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1004:       nep->rits = PETSC_DEFAULT;
1005:     } else {
1007:       nep->rits = its;
1008:     }
1009:     nep->scheme = scheme;
1010:   }
1011:   nep->state = NEP_STATE_INITIAL;
1012:   PetscFunctionReturn(0);
1013: }

1015: /*@C
1016:    NEPGetRefine - Gets the refinement strategy used by the NEP object, and the
1017:    associated parameters.

1019:    Not Collective

1021:    Input Parameter:
1022: .  nep - the nonlinear eigensolver context

1024:    Output Parameters:
1025: +  refine - refinement type
1026: .  npart  - number of partitions of the communicator
1027: .  tol    - the convergence tolerance
1028: .  its    - maximum number of refinement iterations
1029: -  scheme - the scheme used for solving linear systems

1031:    Level: intermediate

1033:    Note:
1034:    The user can specify NULL for any parameter that is not needed.

1036: .seealso: NEPSetRefine()
1037: @*/
1038: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)
1039: {
1041:   if (refine) *refine = nep->refine;
1042:   if (npart)  *npart  = nep->npart;
1043:   if (tol)    *tol    = nep->rtol;
1044:   if (its)    *its    = nep->rits;
1045:   if (scheme) *scheme = nep->scheme;
1046:   PetscFunctionReturn(0);
1047: }

1049: /*@C
1050:    NEPSetOptionsPrefix - Sets the prefix used for searching for all
1051:    NEP options in the database.

1053:    Logically Collective on nep

1055:    Input Parameters:
1056: +  nep - the nonlinear eigensolver context
1057: -  prefix - the prefix string to prepend to all NEP option requests

1059:    Notes:
1060:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1061:    The first character of all runtime options is AUTOMATICALLY the
1062:    hyphen.

1064:    For example, to distinguish between the runtime options for two
1065:    different NEP contexts, one could call
1066: .vb
1067:       NEPSetOptionsPrefix(nep1,"neig1_")
1068:       NEPSetOptionsPrefix(nep2,"neig2_")
1069: .ve

1071:    Level: advanced

1073: .seealso: NEPAppendOptionsPrefix(), NEPGetOptionsPrefix()
1074: @*/
1075: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char *prefix)
1076: {
1078:   if (!nep->V) NEPGetBV(nep,&nep->V);
1079:   BVSetOptionsPrefix(nep->V,prefix);
1080:   if (!nep->ds) NEPGetDS(nep,&nep->ds);
1081:   DSSetOptionsPrefix(nep->ds,prefix);
1082:   if (!nep->rg) NEPGetRG(nep,&nep->rg);
1083:   RGSetOptionsPrefix(nep->rg,prefix);
1084:   PetscObjectSetOptionsPrefix((PetscObject)nep,prefix);
1085:   PetscFunctionReturn(0);
1086: }

1088: /*@C
1089:    NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1090:    NEP options in the database.

1092:    Logically Collective on nep

1094:    Input Parameters:
1095: +  nep - the nonlinear eigensolver context
1096: -  prefix - the prefix string to prepend to all NEP option requests

1098:    Notes:
1099:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1100:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1102:    Level: advanced

1104: .seealso: NEPSetOptionsPrefix(), NEPGetOptionsPrefix()
1105: @*/
1106: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char *prefix)
1107: {
1109:   if (!nep->V) NEPGetBV(nep,&nep->V);
1110:   BVAppendOptionsPrefix(nep->V,prefix);
1111:   if (!nep->ds) NEPGetDS(nep,&nep->ds);
1112:   DSAppendOptionsPrefix(nep->ds,prefix);
1113:   if (!nep->rg) NEPGetRG(nep,&nep->rg);
1114:   RGAppendOptionsPrefix(nep->rg,prefix);
1115:   PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix);
1116:   PetscFunctionReturn(0);
1117: }

1119: /*@C
1120:    NEPGetOptionsPrefix - Gets the prefix used for searching for all
1121:    NEP options in the database.

1123:    Not Collective

1125:    Input Parameters:
1126: .  nep - the nonlinear eigensolver context

1128:    Output Parameters:
1129: .  prefix - pointer to the prefix string used is returned

1131:    Note:
1132:    On the Fortran side, the user should pass in a string 'prefix' of
1133:    sufficient length to hold the prefix.

1135:    Level: advanced

1137: .seealso: NEPSetOptionsPrefix(), NEPAppendOptionsPrefix()
1138: @*/
1139: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])
1140: {
1143:   PetscObjectGetOptionsPrefix((PetscObject)nep,prefix);
1144:   PetscFunctionReturn(0);
1145: }