고유 공간
이 단원에서는 행렬의 고유 값과 고유 벡터를 계산하는 함수들에 대해 다룹니다. 실수 위 대칭, 비대칭, 일반회된 준대칭, 일반화된 비대칭 복소수 위 에르미트, 일반화된 준 에르미트, 그리고 행렬들의 풀이 기능들을 제공합니다. 고유 값들은 고유 벡터들과 함께 계산되거나 단독으로 계산될 수도 있습니다. 에르미트와 실수 대칭 행렬을 계산하는 알고리즘들은 대칭 bidigonalization 과 QR 감소을 이용합니다. 비대칭 알고리즘은 Francis-QR 이중 전이 방법을, 일반화된 비대칭 알고리즘은 Moler와 Strwart의 QZ 방법을 사용합니다.
실수 위 대칭 행렬
실수 위 대칭 행렬의 고유 값, 벡터 문제는 bidiagonalization과 QR-감소 방법을 사용합니다. 이 방법들은 Golub & van Loan, section 8.3 에 기술되어 있습니다. 고유 값은 \(\epsilon ||A||_2\) 의 절대 정확도를 가집니다. \(\epsilon\) 은 계산기 정밀도입니다.
-
type gsl_eigen_symm_workspace
실수 위 대칭 행렬의 고유 값을 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_symm_workspace *gsl_eigen_symm_alloc(const size_t n)
n
-by-n
크기의 실수 위 대칭 행렬의 고유값 문제를 풀기위한 작업 공간을 할당합니다. 작업 공간의 크기는 \(O(2n)\) 입니다.
-
void gsl_eigen_symm_free(gsl_eigen_symm_workspace *w)
할당된 작업 공간
w
를 해제합니다.
-
int gsl_eigen_symm(gsl_matrix *A, gsl_vector *eval, gsl_eigen_symm_workspace *w)
대칭 행렬
A
의 고유값을 계산합니다. 함수를 호출과정에서 반드시 적절한 크기의 작업 공간을w
에 전해주어야 합니다. 행렬A
의 대각, 하삼각 성분들은 계산과정에서 수정됩니다. 반면, 상삼각 성분들은 참조되지 않습니다. 고유 값들은 벡터eval
에 무작위 순서로 저장됩니다.
-
type gsl_eigen_symmv_workspace
실수 대칭 행렬의 고유 값과 고유 벡터를 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_symmv_workspace *gsl_eigen_symmv_alloc(const size_t n)
크기
n
-by-n
의 실수 대칭 행렬의 고유 값, 벡터 문제를 풀기 위한 작업공간을 할당합니다. 공간의 크기는 \(O(4n)\) 입니다.
-
void gsl_eigen_symmv_free(gsl_eigen_symmv_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_symmv(gsl_matrix *A, gsl_vector *eval, gsl_matrix *evec, gsl_eigen_symmv_workspace *w)
실수 대칭 행렬
A
의 고유값, 고유 벡터를 계산합니다. 반드시, 함수 호출과정에서 적절한 크기의 작업 공간을w
에 제공해 주어야합니다. 행렬A
의 대각, 하삼각 성분들은 계산과정에서 수정됩니다. 반면, 상삼각 성분들은 참조되지 않습니다. 고유값들은 벡터eval
에 무작위 순서로 저장됩니다. 대응 되는 고유 벡터들은 행렬evec
의 열에 고유값과 동일 순서로 저장됩니다. 예를 들어, 행렬의 첫번째 열은 첫번째 고유값에 해당하는 고유 벡터입니다. 계산된 고유 벡터들은 서로 직교하고 크기 1로 정규화 되어있습니다.
복소수 위 에르미트 행렬
에르미트 행렬의 계산은 대칭 bidiagonalization과 QR-감소 방법을 사용합니다.
-
type gsl_eigen_herm_workspace
에르미트 행렬의 고유 값을 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_herm_workspace *gsl_eigen_herm_alloc(const size_t n)
n
\(\times\)n
크기의 복소수 에르미트 행렬의 고유 값을 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(3n)\) 입니다.
-
void gsl_eigen_herm_free(gsl_eigen_herm_workspace *w)
T:data:w 가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_herm(gsl_matrix_complex *A, gsl_vector *eval, gsl_eigen_herm_workspace *w)
복소수 에르미트 행렬
A
의 고유 값을 계산합니다. 반드시, 함수 호출과정에서 적절한 크기의 작업 공간을w
에 제공해 주어야합니다. 행렬A
의 대각 성분과 하 삼각 성분은 계산 과정에서 수정됩니다. 하지만, 상 삼각 성분은 계산 과정에서 수정되지 않습니다. 대각 성분의 허수부는 0으로 간주되고 사용되지 않습니다. 고유 값들은 벡터eval
에 무작위로 저장됩니다.
-
type gsl_eigen_hermv_workspace
에르미트 행렬의 고유 값과 고유 벡터를 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_hermv_workspace *gsl_eigen_hermv_alloc(const size_t n)
n
\(\times\)n
크기의 복소수 에르미트 행렬의 고유 값과 고유 벡터를 찾기 위한 작업 공간을 메모리에 할당합니다. 작업 공간의 크기는 \(O(5n)\) 입니다.
-
void gsl_eigen_hermv_free(gsl_eigen_hermv_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_hermv(gsl_matrix_complex *A, gsl_vector *eval, gsl_matrix_complex *evec, gsl_eigen_hermv_workspace *w)
복소수 에르미트 행렬
A
의 고유 벡터와 고유 값을 계산합니다. 반드시, 함수 호출과정에서 적절한 크기의 작업 공간을w
에 제공해 주어야합니다. 행렬A
의 대각 성분과 하 삼각 성분은 계산 과정에서 수정됩니다. 대각 성분의 허수부는 0으로 간주되며 계산 과정에서 사용되지 않습니다. 고유 값들은eval
이 가르키는 벡터에 무자구이 순서로 저장됩니다. 이 고유 값들에 대응되는 고유 벡터들은 행렬evec
에 열 벡터로 저장됩니다. 예로 행렬evec
의 첫번째 열 벡터는eval
의 첫번째 값에 대응되는 고유 벡터입니다. 이 고유 벡터들은 서로 직교하고 크기 1로 정규화 되어 있습니다.
실수 위 비대칭 행렬
실수 비대칭 행렬 \(A\) 의 고유 값과 고유 벡터들을 찾는 문제는 슈어 분해(Schur decomposition)를 이용해 계산할 수 있습니다.
\(Z\) 는 슈어 벡터들의 직교 행렬, \(T\) 은 Schur 형태로 Quasi-상 삼각 행렬입니다. Quasi-상 삼각 행렬은 블럭 상 삼각 행렬로 대각 성분이 \(1 \times 1\) 이나 \(2 \times 2\) 크기의 블럭 행렬인 행렬을 의미합니다. 이 때, \(2 \times 2\) 크기 행렬의 고유 값들은 \(A\) 행렬의 고유 값들의 켤례 복소수입니다. 이 과정에서 사용되는 알고리즘에는 이중 전이 Francis 방법(double-shiftted Francis method)이 사용됩니다.
-
type gsl_eigen_nonsymm_workspace
비대칭 행렬의 고유 값 문제를 풀기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_nonsymm_workspace *gsl_eigen_nonsymm_alloc(const size_t n)
n
\(\times\)n
크기의 실수 비대칭 행렬의 고유 값을 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(2n)\) 입니다.
-
void gsl_eigen_nonsymm_free(gsl_eigen_nonsymm_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
void gsl_eigen_nonsymm_params(const int compute_t, const int balance, gsl_eigen_nonsymm_workspace *w)
gsl_eigen_nonsymm()
함수를 사용하기 위한 기본 계수들을 설정합니다.compute_t
가 1로 정해지면 슈어 형태 \(T\) 가gsl_eigen_nonsymm()
를 이용해 계산됩니다. 0으로 정해지면 \(T\) 가 0인 상황입니다. 완전한 슈어 형태 \(T\) 의 계산은 근사적으로 1.5-2 Flops의 시간을 필요로 합니다.balance
가 1로 정해지면 고유 값을 계산하기 전 군형 변환이 행렬에 적용됩니다. 이 변환은 행렬의 행과 열들이 comparable 노름을 가지도록 만들어 더 정확한 고유 값들을 얻을 수 있습니다. 더 자세한 내용은 Balancing 을 참고할 수 있습니다. 유의할 점은 균형 변환이 슈어 벡터들의 직교성을 보존하지 않는다는 점입니다. 따라서gsl_eigen_nonsymm_Z()
를 이용해 슈어 벡터를 계산하고자 한다면, 원본 행렬의 슈어 벡터가 아니라 균형 변환된 행렬의 슈어 벡터를 얻게됩니다. 이 두 벡터들 사이의 관계는 다음과 같습니다.\[T = Q^T D^{-1} A D Q\]Q
는 행렬로 균형 변환된 행렬에 대한 슈어 벡터들을 나타냅니다.D
는 균형 변환에 대응되는 행렬입니다.gsl_eigen_nonsymm_Z()
함수는 다음을 만족하는 행렬Z
를 계산합니다.\[T = Z^{-1} A Z\]\(Z = D Q\) 를 의미합니다. \(Z\) 는 직교하지 않음을 유의해야 합니다. 이러한 이유로 균형 변환을 포함하면 더 정확한 고유 값들을 계산할 수 있음에도 불구하고 기본 설정에서 제거하고 선택 사항으로 제공합니다.
-
int gsl_eigen_nonsymm(gsl_matrix *A, gsl_vector_complex *eval, gsl_eigen_nonsymm_workspace *w)
실수 위 비대칭 행렬
A
의 고유 값을 계산합니다. 계산된 고유 값들은 벡터eval
에 저장됩니다. \(T\) 값이 필요한 경우, 해당 값들은 꼐산이 끝난 후 행렬A
의 대각 성분을 포함한 상 삼각 부분에 저장됩니다. 계산 후A
의 대각 성분은 \(1 \times 1\) 크기의 고유 값과 고유 값이A
고유 값의 켤레 복소수 값을 가지는 크기 :math:` 2 times 2` 행렬들을 가지고 있습니다.A
의 나머지 부분들의 값은 계산 과정에서 훼손됩니다. 드믄 일이지만 이 함수가 모든 고유 값을 찾지 못하는 경우가 있습니다. 이 경우 오류 값이 반환되고 수렴한 고유 값들의 숫자가w -> n_evals
에 저장됩니다. 수렴된 고유 값들은eval
에 앞 자리부터 채워져 저장됩니다.
-
int gsl_eigen_nonsymm_Z(gsl_matrix *A, gsl_vector_complex *eval, gsl_matrix *Z, gsl_eigen_nonsymm_workspace *w)
gsl_eigen_nonsymm()
와 동일합니다. 다만 추가로 슈어 벡터를 계산하고Z
에 저장한다는 차이가 있습니다.
-
type gsl_eigen_nonsymmv_workspace
비대칭 행렬의 고유 값과 고유 벡터를 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_nonsymmv_workspace *gsl_eigen_nonsymmv_alloc(const size_t n)
n
\(\times\)n
크기의 실수 비대칭 행렬의 고유 값과 고유 벡터를 계산하기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(5n)\) 입니다.
-
void gsl_eigen_nonsymmv_free(gsl_eigen_nonsymmv_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
void gsl_eigen_nonsymmv_params(const int balance, gsl_eigen_nonsymm_workspace *w)
gsl_eigen_nonsymmv()
을 사용하기 위한 계수들을 설정합니다.balance
가 1 로 정해지면 행렬에 균형 변환이 적용됩니다. 자세한 내용은gsl_eigen_nonsymm_params()
를 참고하길 바랍니다. 균형 변환이 슈어 벡터들의 직교성을 보존하지 않기 때문에 기본적으로 비활성화 되어 있습니다.
-
int gsl_eigen_nonsymmv(gsl_matrix *A, gsl_vector_complex *eval, gsl_matrix_complex *evec, gsl_eigen_nonsymmv_workspace *w)
\(n \times n\) 크기의 실수 위 비대칭 행렬
A
의 고유 값과 고유 벡터들을 계산합니다. 먼저gsl_eigen_nonsymm()
를 호출해 고유 값, 슈어 형태 \(T\) 그리고 슈어 벡터를 계산합니다. 그 다음 \(T\) 의 고유 값을 계산하고 슈어 벡터들을 이용해 이들을 역변환합니다. 계산 과정에서 슈어 벡터들의 값은 유지되지 않습니다. 하지만, 함수gsl_eigen_nonsymmv_Z()
를 이용해 저장할 수 있습니다. 계산이 끝난 후A
의 상 삼각부는 슈어 형태 \(T\) 와 같습니다.gsl_eigen_nonsymm()
에서 계산이 실패하면 고유 벡터들이 계산되지 않고 오류 값이 반환됩니다.
-
int gsl_eigen_nonsymmv_Z(gsl_matrix *A, gsl_vector_complex *eval, gsl_matrix_complex *evec, gsl_matrix *Z, gsl_eigen_nonsymmv_workspace *w)
gsl_eigen_nonsymmv()
와 동일합니다. 다만 슈어 벡터를Z
에 저장합니다.
실수 위 일반화된 대칭-정부호 고유계
실수 위에서 일반화된 대칭-정부호 행렬의 고유 문제는 다음과 같은 고유 값과 고유 벡터들을 계산하는 문제입니다.
\(A\) 와 \(B\) 는 대칭 행렬이고, \(B\) 는 양의 정부호 행렬입니다. 이 문제는 \(B\) 에 춈스키 분해를 적용해 표준적인 대칭 행렬의 고유 값 문제로 바꿀 수 있습니다.
따라서, 원래 문제는 대칭 행렬 \(C = L^{-1} A L^{-T}\) 와 \(y = L^T x\) 에 대해 \(C y = \lambda y\) 형태로 표현할 수 있습니다. 이 경우 표준적인 대칭 행렬의 고유 문제 풀이 기능을 행렬 \(C\) 에 적용할 수 있습니다. 해당 계산 결과로 나온 고유 값들을 역변환해 원래 행렬의 고유 값들을 찾을 수 있습니다. 일반화된 대칭-정부호 행렬의 고유 값과 고유 벡터들은 항상 실수 있습니다.
-
type gsl_eigen_gensymm_workspace
일반화된 대칭 행렬의 고유 값 문제를 풀기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_gensymm_workspace *gsl_eigen_gensymm_alloc(const size_t n)
n
\(\times\)n
크기의 일반화된 대칭-정부호 행렬의 고유 값을 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(2n)\) 입니다.
-
void gsl_eigen_gensymm_free(gsl_eigen_gensymm_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_gensymm(gsl_matrix *A, gsl_matrix *B, gsl_vector *eval, gsl_eigen_gensymm_workspace *w)
일반화된 대칭-정부호 행렬쌍 (
A
,B
)의 고유 값을 계산해eval
에 저장합니다. 이 함수는 위에 기술된 일반 대칭 행렬 문제로 변환하는 방법을 사용합니다. 계산후B
에는B
의 춈스키 분해가 저장되고A
는 본래 값을 잃습니다.
-
type gsl_eigen_gensymmv_workspace
일반화된 대칭 행렬의 고유 값과 고유 벡터를 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_gensymmv_workspace *gsl_eigen_gensymmv_alloc(const size_t n)
n
\(\times\)n
크기의 실수 일반화된 대칭-정부호 행렬의 고유 값과 고유 벡터를 계산하기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(4n)\) 입니다.
-
void gsl_eigen_gensymmv_free(gsl_eigen_gensymmv_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_gensymmv(gsl_matrix *A, gsl_matrix *B, gsl_vector *eval, gsl_matrix *evec, gsl_eigen_gensymmv_workspace *w)
실수 위 일반화된 대칭-정부호 행렬 쌍(
A
,B
)의 고유 값과 고유 벡터들을 계산합니다. 계산된 고유 값과 벡터들은 각각eval
과evec
에 저장됩니다. 계산된 고유 벡터들은 서로 직교하고 크기 1로 정규화 되어있습니다. 계산 후B
에는B
의 춈스키 분해가 저장되어 있고A
는 계산 과정에서 훼손됩니다.
복소수 위 일반화된 에르미트-정부호 고유계
복소수 위 일반화딘 에르미트-정부호 행렬의 고유 값 문제는 다음을 만족하는 고유 값 \(\lambda\) 와 고유 벡터 \(x\) 를 찾는 문제입니다.
\(A\) 와 \(B\) 는 에르미트 행렬이고, \(B\) 는 추가로 정부호 행렬이기도 합니다. 실수 위와 비슷하게 이 문제는 에르미트 행렬 \(C = L^{-1} A L^{-\dagger}\) 와 \(y = L^{\dagger} x\) 로 \(C y = \lambda y\) 형태로 변환할 수 있습니다. 에르미트 고유 문제의 표준 해결 방법을 \(C\) 에 적용해 이 문제를 풀수 있고, 풀이로 나온 고유 벡터들에 역변환을 적용해 원래 문제의 해답을 얻을 수 있습니다. 일반화된 에르미트-정부호 고유 문제의 답은 언제나 실수입니다.
-
type gsl_eigen_genherm_workspace
일반화된 에르미트-정부호 행렬의 고유 값 문제를 풀기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_genherm_workspace *gsl_eigen_genherm_alloc(const size_t n)
복소수 위
n
\(\times\)n
크기의 일반화된 에르미트-정부호 행렬의 고유 값을 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(3n)\) 입니다.
-
void gsl_eigen_genherm_free(gsl_eigen_genherm_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_genherm(gsl_matrix_complex *A, gsl_matrix_complex *B, gsl_vector *eval, gsl_eigen_genherm_workspace *w)
복소수 위에서 일반화된 에르미트-정부호 행렬 쌍 (
A
,B
)의 고유 값을 계산해eval
에 저장합니다. 사용하는 방법은 위에 기술되어 있습니다. 계산 후B
에는B
의 춈스키 분해가 저장되고A
는 계산과정에서 훼손됩니다.
-
type gsl_eigen_genhermv_workspace
일반화된 에르미트-정부호 행렬의 고유 값가 고유 벡터 문제를 풀기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_genhermv_workspace *gsl_eigen_genhermv_alloc(const size_t n)
복소수 위
n
\(\times\)n
크기의 일반화된 에르미트-정부호 행렬의 고유 값과 고유 벡터를 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(5n)\) 입니다.
-
void gsl_eigen_genhermv_free(gsl_eigen_genhermv_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_genhermv(gsl_matrix_complex *A, gsl_matrix_complex *B, gsl_vector *eval, gsl_matrix_complex *evec, gsl_eigen_genhermv_workspace *w)
복소수 위 일반화된 에르미트-정부호 행렬 쌍 (
A
,B
)의 고유 값과 고유 벡터를 계산하고 이들을 각각eval
와evec
에 저장합니다. 계산된 고유 벡터들은 서로 직교하고 크기 1로 정규화 되어있습니다. 계산 후B
에는B
의 춈스키 분해가 저장되고A
는 계산과정에서 훼손됩니다.
실수 위 일반화된 비대칭 고유계
주어진 두 정사각 행렬 ( \(A\), \(B\) ) 에 대해, 일반화된 비대칭 계의 풀이는 다음을 만족하는 고유 값 \(\lambda\) 와 고유 벡터 \(x\) 를 찾는 문제입니다.
고유 값 \(mu\) 와 고유 벡터 \(y\) 에 대해 다음과 같은 형식도 가능합니다.
이 두 문제는 \(\lambda = 1/\mu\) 식을 이용하면 동일한 식임을 알 수 있습니다. 이때, \(lambda\) 와 \(mu\) 는 모두 0이 아니여야 합니다. \(lambda\) 가 0이여도 첫번째 식은 여전히 잘 정의된 고유계를 형성하지만 두번째 식은 \(mu\) 가 정의되지 않습니다. 고유 값으로 0과 무한대를 사용하기 위해 라이브러리에서 제공하는 함수들은 다음과 같은 형태의 고유 계를 풀게 됩니다.
풀이 함수들은 \(\alpha\) 와 \(\beta\) 두 값을 반환해 사용자가 \(\lambda\) 나 \(\mu\) 를 계산해 사용할 수 있게 합니다. 각각은 나눗셈을 이용해서 \(\lambda = \alpha / \beta\) 와 \(\mu = \beta / \alpha\) 로 계산할 수 있습니다.
선형 행렬 pencil \(A - \lambda B\) 의 행렬식이 모든 \(\lambda\) 에 대해 0이면 해당 고유 계는 특이성을 가집니다. 특이성을 가지는 고유 계는 \(\alpha = \beta = 0\) 를 가지게 됩니다. 이는 해당 고유 계가 제대로 정의된 고유 값과 벡터를 가지지 않음을 의미합니다. 이 단원의 함수들은 특이성을 가지지 않는 일반적인 문제들에 대한 풀이를 제공합니다. 이러한 특이 pencil 값을 가지는 문제에 함수들을 적용할 경우 예기치 못한 결과가 나올 수 있습니다.
행렬 쌍 \((A, B)\) 의 실수 위의 일반화된 비 대칭 고유 계 문제는 일반화된 슈어 분해와 관련되어 있습니다.
\(Q\) 와 \(Z\) 는 각각 좌, 우 슈어 벡터들의 직교 행렬입니다. \((S, T)\) 는 일반화된 슈어 형태로 이들의 대각 성분은 각각 \(\alpha\) 와 \(\beta\) 값들을 제공해줍니다. 이 알고리즘들은 Moler & StewartQZ 의 QZ 방법에 기반해 있습니다. (참고 문헌을 확인할 수 있습니다.)
-
type gsl_eigen_gen_workspace
일반화된 비대칭 행렬의 고유 값 문제를 풀기 위한 내부 인자들을 가지고 있는 작업 공간입니다. problems.
-
gsl_eigen_gen_workspace *gsl_eigen_gen_alloc(const size_t n)
n
\(\times\)n
크기의 실수 일반화된 비대칭 행렬의 고유 값을 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(n)\) 입니다.
-
void gsl_eigen_gen_free(gsl_eigen_gen_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
void gsl_eigen_gen_params(const int compute_s, const int compute_t, const int balance, gsl_eigen_gen_workspace *w)
gsl_eigen_gen()
을 사용하기 위한 계수들을 설정합니다.compute_s
가 1로 정해지면 완전한 슈어 형태 \(S\) 함수gsl_eigen_gen()
에 의해 계산됩니다. 0이 주어지면 \(S\) 는 계산되지 않습니다. 0이 기본 설정입니다. \(S\) 는 Schur 형태로 Quasi-상 삼각 행렬입니다. Quasi-상 삼각 행렬은 블럭 상 삼각 행렬로 대각 성분이 \(1 \times 1\) 이나 \(2 \times 2\) 크기의 블럭 행렬인 행렬을 의미합니다. 이 때, \(2 \times 2\) 크기 행렬의 고유 값들은 \(A\) 행렬의 고유 값들의 켤례 복소수입니다. \(1 \times 1\) 행렬은 실수 고유 값과 대응됩니다.compute_t
가 1로 정해지면 완전한 슈어 형태 \(T\) 가 함수gsl_eigen_gen()
에 의해 계산됩니다. 0이 주어지면 \(T\) 는 계산되지 않습니다. 0이 기본 설정입니다. \(T\) 는 상 삼각 행렬로 대각선에 음수 값이 없는 행렬입니다. \(S\) 의 \(2 \times 2\) 크기 블럭 행렬은 \(T\) 의 \(2 \times 2\) 크기 블럭 행렬에 대응됩니다.balance
값은 현재 사용되지 않습니다. 잃반화된 균형 행렬기능은 아직 구현되지 않았습니다.
-
int gsl_eigen_gen(gsl_matrix *A, gsl_matrix *B, gsl_vector_complex *alpha, gsl_vector *beta, gsl_eigen_gen_workspace *w)
- 실수 위 일반화된 비대칭 행렬 쌍 (
A
,B
) 의 고유 값을 계산해 (
alpha
,beta
)에 저장합니다.alpha
는 복소수,beta
는 실수입니다. \(\beta_i\) 가 0이 아니라면 \(\lambda = \alpha_i / \beta_i\) 이 고유 값이 됩니다. 비슷하게 \(\alpha_i\) 가 0이 아니라면 \(\mu = \beta_i / \alpha_i\) 이 \(\mu A y = B y\) 의 고유 값이 됩니다.beta
의 원소들은 모두 음수가 되지 않도록 정규화되어 있습니다.\(S\) 가 필요하다면 이 행렬은 계산이 끝난 후
A
에 저장되어 있습니다. \(T\) 도 계산이 끝난 후B
에 저장됩니다. (alpha
,beta
) 의 고유 값 순서는 슈어 형태 \(S\) 와 \(T\) 의 대각선 블럭들 순서를 따릅니다. 드문 일이지만 함수에서 모든 고유 값들을 찾지 못하는 경우가 있습니다. 이 경우 오류 값이 반환됩니다.
- 실수 위 일반화된 비대칭 행렬 쌍 (
-
int gsl_eigen_gen_QZ(gsl_matrix *A, gsl_matrix *B, gsl_vector_complex *alpha, gsl_vector *beta, gsl_matrix *Q, gsl_matrix *Z, gsl_eigen_gen_workspace *w)
gsl_eigen_gen()
함수와 동일합니다. 다만 좌, 우 슈어 벡터들을 계산해 각각Q
와Z
에 저장합니다.
-
type gsl_eigen_genv_workspace
일반화된 고유 값과 고유 벡터를 계산하기 위한 내부 인자들을 가지고 있는 작업 공간입니다.
-
gsl_eigen_genv_workspace *gsl_eigen_genv_alloc(const size_t n)
실수 위 크기 \(n \times n\) 의 일반화된 비 대칭 고유 계의 고유 값과 고유 벡터를 찾기 위한 작업 공간을 메모리에 할당합니다. 공간의 크기는 \(O(7n)\) 입니다.
-
void gsl_eigen_genv_free(gsl_eigen_genv_workspace *w)
w
가 가르키는 작업 공간 영역을 해제합니다.
-
int gsl_eigen_genv(gsl_matrix *A, gsl_matrix *B, gsl_vector_complex *alpha, gsl_vector *beta, gsl_matrix_complex *evec, gsl_eigen_genv_workspace *w)
실수 위 크기 \(n \times n\) 의 일반화된 비 대칭 행렬 쌍(
A
,B
)의 고유 값과 우 고유 벡터를 찾습니다. 고유 벡터들은 (alpha
,beta
)에 저장되고 고유 벡터들은evec
에 저장됩니다. 이 함수는 먼저 함수gsl_eigen_gen()
를 호출해 고유 값들과 슈어 형태 그리고 슈어 벡터를 계산합니다. 그 후 슈어 형태에 대한 고유 벡터를 찾고 이들을 슈어 벡터를 이용해 역변환 시켜 원래 계의 고유 벡터들을 얻습니다. 슈어 벡터들은 계산과정에서 사용되 훼손됩니다. 함수gsl_eigen_genv_QZ()
를 사용하면 별도로 저장해둘 수 있습니다. 계산된 고유 벡터는 크기 1을 가지도록 정규화 되어 있습니다. 계산 후 (A
,B
) 는 일반화된 슈어 형태 ( \(S\) , \(T\) )가 됩니다. 함수gsl_eigen_gen()
의 구동이 실패하면 고유 벡터들은 계산되지 않고 오류 코드가 반환됩니다.
-
int gsl_eigen_genv_QZ(gsl_matrix *A, gsl_matrix *B, gsl_vector_complex *alpha, gsl_vector *beta, gsl_matrix_complex *evec, gsl_matrix *Q, gsl_matrix *Z, gsl_eigen_genv_workspace *w)
gsl_eigen_genv()
와 동일합니다. 다만, 좌우 슈어 벡터들을 계산해 각각Q
와Z
에 저장합니다.
고유 값과 벡터의 정렬
-
int gsl_eigen_symmv_sort(gsl_vector *eval, gsl_matrix *evec, gsl_eigen_sort_t sort_type)
eval
벡터에 저장된 고유 값과, 대응되는 실수 고유 벡터들을 동시에 정렬합니다. 해당 벡터들은 행렬evec
에 열 벡터로 저장되어 있습니다.sort_type
의 값에 따라 내림, 오름차순 정렬이 결정됩니다.-
type gsl_eigen_sort_t
GSL_EIGEN_SORT_VAL_ASC
수치 값에 따른 오름차순 정렬
GSL_EIGEN_SORT_VAL_DESC
수치 값에 따른 내림차순 정렬
GSL_EIGEN_SORT_ABS_ASC
크기에 따른 오름차순 정렬
GSL_EIGEN_SORT_ABS_DESC
크기에 따른 내림차순 정렬
-
type gsl_eigen_sort_t
-
int gsl_eigen_hermv_sort(gsl_vector *eval, gsl_matrix_complex *evec, gsl_eigen_sort_t sort_type)
eval
벡터에 저장된 고유 값과, 행렬evec
에 열 벡터로 저장된 복소수 고유 벡터들을 동시에 정렬합니다.sort_type
의 값에 따라 내림, 오름차순 정렬이 결정됩니다. 해당 값들은 위 표에 있습니다.
-
int gsl_eigen_nonsymmv_sort(gsl_vector_complex *eval, gsl_matrix_complex *evec, gsl_eigen_sort_t sort_type)
eval
벡터에 저장된 고유 값과, 대응되는 복소수 고유 벡터들을 동시에 정렬합니다. 해당 벡터들은 행렬evec
에 열 벡터로 저장되어 있습니다.sort_type
의 값에 따라 내림, 오름차순 정렬이 결정됩니다. 해당 값들은 위 표에 있습니다. 고유 값들이 복소수이기 때문에GSL_EIGEN_SORT_ABS_ASC
와GSL_EIGEN_SORT_ABS_DESC
매크로만 사용 가능합니다.
-
int gsl_eigen_gensymmv_sort(gsl_vector *eval, gsl_matrix *evec, gsl_eigen_sort_t sort_type)
eval
벡터에 저장된 고유 값과, 대응되는 실수 고유 벡터들을 동시에 정렬합니다. 해당 벡터들은 행렬evec
에 열 벡터로 저장되어 있습니다.sort_type
의 값에 따라 내림, 오름차순 정렬이 결정됩니다.sort_type
값은 이전 표에 기술되어 있습니다.
-
int gsl_eigen_genhermv_sort(gsl_vector *eval, gsl_matrix_complex *evec, gsl_eigen_sort_t sort_type)
eval
벡터에 저장된 고유 값과, 대응되는 복소수 고유 벡터들을 동시에 정렬합니다. 해당 벡터들은 행렬evec
에 열 벡터로 저장되어 있습니다.sort_type
의 값에 따라 내림, 오름차순 정렬이 결정됩니다.sort_type
값은 이전 표에 기술되어 있습니다.
-
int gsl_eigen_genv_sort(gsl_vector_complex *alpha, gsl_vector *beta, gsl_matrix_complex *evec, gsl_eigen_sort_t sort_type)
두 벡터
alpha
,beta
에 저장된 고유 값들과 대응 되는 복소수 고유 벡터들을 동시에 정렬합니다. 대응되는 고유 벡터들은 행렬evec
에 열 벡터로 저장되어 있습니다.sort_type
의 값에 따라 내림, 오름차순 정렬이 결정됩니다.sort_type
값은 이전 표에 기술되어 있습니다. 고유 값들이 복소수이기 때문에GSL_EIGEN_SORT_ABS_ASC
와GSL_EIGEN_SORT_ABS_DESC
매크로만 사용 가능합니다.
예제
다음 프로그램은 4차 힐베르트 행렬 \(H(i,j) = 1/(i+j+1)\) 의 고유 값과 고유 벡터를 계산합니다.
#include <stdio.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_eigen.h>
int
main (void)
{
double data[] = { 1.0 , 1/2.0, 1/3.0, 1/4.0,
1/2.0, 1/3.0, 1/4.0, 1/5.0,
1/3.0, 1/4.0, 1/5.0, 1/6.0,
1/4.0, 1/5.0, 1/6.0, 1/7.0 };
gsl_matrix_view m
= gsl_matrix_view_array (data, 4, 4);
gsl_vector *eval = gsl_vector_alloc (4);
gsl_matrix *evec = gsl_matrix_alloc (4, 4);
gsl_eigen_symmv_workspace * w =
gsl_eigen_symmv_alloc (4);
gsl_eigen_symmv (&m.matrix, eval, evec, w);
gsl_eigen_symmv_free (w);
gsl_eigen_symmv_sort (eval, evec,
GSL_EIGEN_SORT_ABS_ASC);
{
int i;
for (i = 0; i < 4; i++)
{
double eval_i
= gsl_vector_get (eval, i);
gsl_vector_view evec_i
= gsl_matrix_column (evec, i);
printf ("eigenvalue = %g\n", eval_i);
printf ("eigenvector = \n");
gsl_vector_fprintf (stdout,
&evec_i.vector, "%g");
}
}
gsl_vector_free (eval);
gsl_matrix_free (evec);
return 0;
}
다음은 프로그램의 실행 결과입니다.
$ ./a.out
eigenvalue = 9.67023e-05
eigenvector =
-0.0291933
0.328712
-0.791411
0.514553
...
GNU octave 의 결과 값과 비교해 볼 수 있습니다.
octave> [v,d] = eig(hilb(4));
octave> diag(d)
ans =
9.6702e-05
6.7383e-03
1.6914e-01
1.5002e+00
octave> v
v =
0.029193 0.179186 -0.582076 0.792608
-0.328712 -0.741918 0.370502 0.451923
0.791411 0.100228 0.509579 0.322416
-0.514553 0.638283 0.514048 0.252161
참고
고유 벡터들은 부호에 의해 달라질 수 있습니다. 고유 벡터들의 부호는 무작위로 정해집니다.
다음 프로그램은 비대칭 행렬의 고유문제를 계산합니다. 이 프로그램은 \(x = (-1,-2,3,4)\) 에 대해, Vandermonde 행렬 \(V(x;i,j) = x_i^{n - j}\) 을 사용했습니다.
#include <stdio.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_eigen.h>
int
main (void)
{
double data[] = { -1.0, 1.0, -1.0, 1.0,
-8.0, 4.0, -2.0, 1.0,
27.0, 9.0, 3.0, 1.0,
64.0, 16.0, 4.0, 1.0 };
gsl_matrix_view m
= gsl_matrix_view_array (data, 4, 4);
gsl_vector_complex *eval = gsl_vector_complex_alloc (4);
gsl_matrix_complex *evec = gsl_matrix_complex_alloc (4, 4);
gsl_eigen_nonsymmv_workspace * w =
gsl_eigen_nonsymmv_alloc (4);
gsl_eigen_nonsymmv (&m.matrix, eval, evec, w);
gsl_eigen_nonsymmv_free (w);
gsl_eigen_nonsymmv_sort (eval, evec,
GSL_EIGEN_SORT_ABS_DESC);
{
int i, j;
for (i = 0; i < 4; i++)
{
gsl_complex eval_i
= gsl_vector_complex_get (eval, i);
gsl_vector_complex_view evec_i
= gsl_matrix_complex_column (evec, i);
printf ("eigenvalue = %g + %gi\n",
GSL_REAL(eval_i), GSL_IMAG(eval_i));
printf ("eigenvector = \n");
for (j = 0; j < 4; ++j)
{
gsl_complex z =
gsl_vector_complex_get(&evec_i.vector, j);
printf("%g + %gi\n", GSL_REAL(z), GSL_IMAG(z));
}
}
}
gsl_vector_complex_free(eval);
gsl_matrix_complex_free(evec);
return 0;
}
다음은 프로그램의 실행 결과입니다.
$ ./a.out
eigenvalue = -6.41391 + 0i
eigenvector =
-0.0998822 + 0i
-0.111251 + 0i
0.292501 + 0i
0.944505 + 0i
eigenvalue = 5.54555 + 3.08545i
eigenvector =
-0.043487 + -0.0076308i
0.0642377 + -0.142127i
-0.515253 + 0.0405118i
-0.840592 + -0.00148565i
...
GNU octave 의 결과 값과 비교해 볼 수 있습니다.
octave> [v,d] = eig(vander([-1 -2 3 4]));
octave> diag(d)
ans =
-6.4139 + 0.0000i
5.5456 + 3.0854i
5.5456 - 3.0854i
2.3228 + 0.0000i
octave> v
v =
Columns 1 through 3:
-0.09988 + 0.00000i -0.04350 - 0.00755i -0.04350 + 0.00755i
-0.11125 + 0.00000i 0.06399 - 0.14224i 0.06399 + 0.14224i
0.29250 + 0.00000i -0.51518 + 0.04142i -0.51518 - 0.04142i
0.94451 + 0.00000i -0.84059 + 0.00000i -0.84059 - 0.00000i
Column 4:
-0.14493 + 0.00000i
0.35660 + 0.00000i
0.91937 + 0.00000i
0.08118 + 0.00000i
참고
고유 값 \(5.54555 + 3.08545i\) 의 대응 고유 벡터는 곱 계수 \(0.9999984 + 0.0017674i\) 에 의해 다릅니다. 이 계수는 크기 1의 임의의 위상 계수입니다.
참고 자료와 추가 문헌
이 단원에 기술된 알고리즘들은 다음 책에 자세히 기술되어 있습니다.
G. H. Golub, C. F. Van Loan, “Matrix Computations” (3rd Ed, 1996), Johns Hopkins University Press, ISBN 0-8018-5414-8.
일반화된 고유계의 QZ 알고리즘은 다음 논문을 참고할 수 있습니다.
C. Moler, G. Stewart, “An Algorithm for Generalized Matrix Eigenvalue Problems”, SIAM J. Numer. Anal., Vol 10, No 2, 1973.
큰 규모의 행렬들에 대한 고유 문제 풀이 기능들은 포트란 라이브러리 LAPACK 에서 찾을 수 있습니다. LAPACK 은 다음 문헌에 기술되어 있습니다.
LAPACK Users’ Guide (Third Edition, 1999), Published by SIAM, ISBN 0-89871-447-8.
LAPACK 의 소스코드는 사용자 설명서와 함께 http://www.netlib.org/lapack 에서 찾을 수 있습니다.