EPSDenseHessenberg

Computes the Hessenberg form of a dense matrix.

Synopsis

#include "slepceps.h" 
PetscErrorCode EPSDenseHessenberg(int n,int k,PetscScalar *A,int lda,PetscScalar *Q)
Not Collective

Input Parameters

n - dimension of the matrix
k - first active column
lda - leading dimension of A

Input/Output Parameters

A - on entry, the full matrix; on exit, the upper Hessenberg matrix (H)
Q - on exit, orthogonal matrix of vectors A = Q*H*Q'

Notes

Only active columns (from k to n) are computed.

Both A and Q are overwritten.

This routine uses LAPACK routines xGEHRD and xORGHR/xUNGHR.

See Also

EPSDenseSchur(), EPSSortDenseSchur(), EPSDenseTridiagonal()

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