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 - the polynomial 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)=K+\lambda C+\lambda^2M\), with symmetric (or Hermitian) coefficient matrices \(K\), \(C\), \(M\).
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 by Niendorf and Voss [2010]. Furthermore, if
hyperbolic=1 then only xi is computed.
References#
V. Niendorf and H. Voss. Detecting hyperbolic and definite matrix polynomials. Linear Algebra Appl., 432(4):1017–1035, 2010. doi:10.1016/j.laa.2009.10.014.
See Also#
Level#
advanced
Location#
Examples#
Index of all PEP routines Table of Contents for all manual pages Index of all manual pages