PEPCheckDefiniteQEP#

Determines if a symmetric/Hermitian quadratic eigenvalue problem is definite or not.

Synopsis#

#include "slepcpep.h" 
PetscErrorCode PEPCheckDefiniteQEP(PEP pep,PetscReal *xi,PetscReal *mu,PetscInt *definite,PetscInt *hyperbolic)

Collective

Input Parameter#

  • pep - eigensolver context

Output Parameters#

  • xi - first computed parameter

  • mu - second computed parameter

  • definite - flag indicating that the problem is definite

  • hyperbolic - flag indicating that the problem is hyperbolic

Notes#

This function is intended for quadratic eigenvalue problems, Q(lambda)=Alambda^2+Blambda+C, with symmetric (or Hermitian) coefficient matrices A,B,C.

On output, the flag ‘definite’ may have the values -1 (meaning that the QEP is not definite), 1 (if the problem is definite), or 0 if the algorithm was not able to determine whether the problem is definite or not.

If definite=1, the output flag ‘hyperbolic’ informs in a similar way about whether the problem is hyperbolic or not.

If definite=1, the computed values xi and mu satisfy Q(xi)<0 and Q(mu)>0, as obtained via the method proposed in [Niendorf and Voss, LAA 2010]. Furthermore, if hyperbolic=1 then only xi is computed.

See Also#

PEPSetProblemType()

Level#

advanced

Location#

src/pep/impls/krylov/stoar/qslice.c

Examples#

src/pep/tutorials/ex40.c


Index of all PEP routines Table of Contents for all manual pages Index of all manual pages