EPSDenseGNHEP

Solves a dense Generalized non-Hermitian Eigenvalue Problem.

Synopsis

#include "slepceps.h" 
PetscErrorCode EPSDenseGNHEP(PetscInt n_,PetscScalar *A,PetscScalar *B,PetscScalar *w,PetscScalar *wi,PetscScalar *V,PetscScalar *W)
Not Collective

Input Parameters

n - dimension of the eigenproblem
A - pointer to the array containing the matrix values for A
B - pointer to the array containing the matrix values for B

Output Parameters

w - pointer to the array to store the computed eigenvalues
wi - imaginary part of the eigenvalues (only when using real numbers)
V - pointer to the array to store right eigenvectors
W - pointer to the array to store left eigenvectors

Notes

If either V or W are PETSC_NULL then the corresponding eigenvectors are not computed.

Matrices A and B are overwritten.

This routine uses LAPACK routines xGGEVX.

See Also

EPSDenseNHEP(), EPSDenseHEP(), EPSDenseGHEP()

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