EPSSortDenseSchurTarget

Reorders the Schur decomposition computed by EPSDenseSchur().

Synopsis

#include "slepceps.h" 
PetscErrorCode EPSSortDenseSchurTarget(PetscInt n_,PetscInt k,PetscScalar *T,PetscInt ldt_,PetscScalar *Z,PetscScalar *wr,PetscScalar *wi,PetscScalar target,EPSWhich which)
Not Collective

Input Parameters

n - dimension of the matrix
k - first active column
ldt - leading dimension of T
target - the target value
which - eigenvalue sort order

Input/Output Parameters

T - the upper (quasi-)triangular matrix
Z - the orthogonal matrix of Schur vectors
wr - pointer to the array to store the computed eigenvalues
wi - imaginary part of the eigenvalues (only when using real numbers)

Notes

This function reorders the eigenvalues in wr,wi located in positions k to n according to increasing distance to the target. The parameter which is used to determine if distance is relative to magnitude, real axis, or imaginary axis. The Schur decomposition Z*T*Z^T, is also reordered by means of rotations so that eigenvalues in the diagonal blocks of T follow the same order.

Both T and Z are overwritten.

This routine uses LAPACK routines xTREXC.

See Also

EPSDenseHessenberg(), EPSDenseSchur(), EPSDenseTridiagonal()

Location: src/eps/interface/dense.c
Index of all EPS routines
Table of Contents for all manual pages
Index of all manual pages