複素一般化エルミート固有値問題: エルミート定値一般化行列 : (分割統治法)

LAPACKサンプルソースコード : 使用ルーチン名:ZHEGVD

概要

本サンプルはFortran言語によりLAPACKルーチンZHEGVDを利用するサンプルプログラムです。

一般化エルミート固有値問題 $ A B x = \lambda x$のすべての固有値と固有ベクトルを求めます。

\begin{displaymath}
A = \left(
\begin{array}{cccc}
-7.36 & 0.77 - 0.43i & -0....
...97i & 1.90 - 3.73i & 2.88 + 3.17i & -2.54
\end{array} \right)
\end{displaymath}

及び

\begin{displaymath}
B = \left(
\begin{array}{cccc}
3.23 & 1.51 - 1.92i & 1.90...
...0i & -1.18 - 1.37i & 2.33 + 0.14i & 4.29
\end{array} \right),
\end{displaymath}

$ B$の条件数の推定値と計算された固有値と固有ベクトルの誤差限界の近似値も合わせて求めます。

ZHEGVの例題プログラムは一般化エルミート固有値問題$ A x = \lambda B x$の解き方を示します。

入力データ

(本ルーチンの詳細はZHEGVD のマニュアルページを参照)

このデータをダウンロード
ZHEGVD Example Program Data

   4                                                        :Value of N

 (-7.36, 0.00) ( 0.77, -0.43) (-0.64, -0.92) ( 3.01, -6.97)
               ( 3.49,  0.00) ( 2.19,  4.45) ( 1.90,  3.73)
                              ( 0.12,  0.00) ( 2.88, -3.17)
                                             (-2.54,  0.00) :End of matrix A

 ( 3.23, 0.00) ( 1.51, -1.92) ( 1.90,  0.84) ( 0.42,  2.50)
               ( 3.58,  0.00) (-0.23,  1.11) (-1.18,  1.37)
                              ( 4.09,  0.00) ( 2.33, -0.14)
                                             ( 4.29,  0.00) :End of matrix B

出力結果

(本ルーチンの詳細はZHEGVD のマニュアルページを参照)

この出力例をダウンロード
 ZHEGVD Example Program Results

 Eigenvalues
      -61.7321    -6.6195     0.0725    43.1883
 Eigenvectors
             1          2          3          4
 1      0.3903    -0.1560     2.2909    -0.1943
        0.0000    -0.0404     0.0000    -0.0690
 
 2     -0.1814    -0.1552    -0.5042     0.3884
        0.0114    -0.3651    -0.7120     0.0000
 
 3      0.0438     0.5364    -1.2701     0.0657
        0.0338     0.0000    -0.4547    -0.2095
 
 4     -0.2221    -0.1298     0.5706     0.2924
       -0.2272    -0.1880     1.3132    -0.0675

 Estimate of reciprocal condition number for B
        2.5E-03

 Error estimates (relative to machine precision)
 for the eigenvalues:
        2.4E+04    2.8E+03    2.3E+02    1.7E+04

 for the eigenvectors:
        4.7E+02    1.0E+03    1.0E+03    4.9E+02

ソースコード

(本ルーチンの詳細はZHEGVD のマニュアルページを参照)

※本サンプルソースコードのご利用手順は「サンプルのコンパイル及び実行方法」をご参照下さい。


このソースコードをダウンロード
    Program zhegvd_example

!     ZHEGVD Example Program Text

!     Copyright 2017, Numerical Algorithms Group Ltd. http://www.nag.com

!     .. Use Statements ..
      Use lapack_example_aux, Only: nagf_file_print_matrix_complex_gen
      Use lapack_interfaces, Only: ddisna, zhegvd, zlanhe, ztrcon
      Use lapack_precision, Only: dp
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Real (Kind=dp), Parameter :: one = 1.0E+0_dp
      Integer, Parameter :: nb = 64, nin = 5, nout = 6
!     .. Local Scalars ..
      Complex (Kind=dp) :: scal
      Real (Kind=dp) :: anorm, bnorm, eps, rcond, rcondb, t1, t2
      Integer :: i, ifail, info, k, lda, ldb, liwork, lrwork, lwork, n
!     .. Local Arrays ..
      Complex (Kind=dp), Allocatable :: a(:, :), b(:, :), work(:)
      Complex (Kind=dp) :: cdum(1)
      Real (Kind=dp), Allocatable :: eerbnd(:), rcondz(:), rwork(:), w(:), &
        zerbnd(:)
      Real (Kind=dp) :: rdum(1)
      Integer :: idum(1)
      Integer, Allocatable :: iwork(:)
!     .. Intrinsic Procedures ..
      Intrinsic :: abs, conjg, epsilon, max, maxloc, nint, real
!     .. Executable Statements ..
      Write (nout, *) 'ZHEGVD Example Program Results'
      Write (nout, *)
!     Skip heading in data file
      Read (nin, *)
      Read (nin, *) n
      lda = n
      ldb = n
      Allocate (a(lda,n), b(ldb,n), eerbnd(n), rcondz(n), w(n), zerbnd(n))

!     Use routine workspace query to get optimal workspace.
      lwork = -1
      liwork = -1
      lrwork = -1
      Call zhegvd(2, 'Vectors', 'Upper', n, a, lda, b, ldb, w, cdum, lwork, &
        rdum, lrwork, idum, liwork, info)

!     Make sure that there is enough workspace for block size nb.
      lwork = max((nb+2+n)*n, nint(real(cdum(1))))
      lrwork = max(1+(5+2*n)*n, nint(rdum(1)))
      liwork = max(3+5*n, idum(1))
      Allocate (work(lwork), rwork(lrwork), iwork(liwork))

!     Read the upper triangular parts of the matrices A and B

      Read (nin, *)(a(i,i:n), i=1, n)
      Read (nin, *)(b(i,i:n), i=1, n)

!     Compute the one-norms of the symmetric matrices A and B

      anorm = zlanhe('One norm', 'Upper', n, a, lda, rwork)
      bnorm = zlanhe('One norm', 'Upper', n, b, ldb, rwork)

!     Solve the generalized Hermitian eigenvalue problem
!     A*B*x = lambda*x (itype = 2)
      Call zhegvd(2, 'Vectors', 'Upper', n, a, lda, b, ldb, w, work, lwork, &
        rwork, lrwork, iwork, liwork, info)

      If (info==0) Then

!       Print solution

        Write (nout, *) 'Eigenvalues'
        Write (nout, 100) w(1:n)
        Flush (nout)

!       Normalize the eigenvectors, largest element real
        Do i = 1, n
          rwork(1:n) = abs(a(1:n,i))
          k = maxloc(rwork(1:n), 1)
          scal = conjg(a(k,i))/abs(a(k,i))
          a(1:n, i) = a(1:n, i)*scal
        End Do

!       ifail: behaviour on error exit
!              =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft
        ifail = 0
        Call nagf_file_print_matrix_complex_gen('General', ' ', n, n, a, lda, &
          'Eigenvectors', ifail)

!       Call ZTRCON to estimate the reciprocal condition
!       number of the Cholesky factor of B.  Note that:
!       cond(B) = 1/rcond**2
        Call ztrcon('One norm', 'Upper', 'Non-unit', n, b, ldb, rcond, work, &
          rwork, info)

!       Print the reciprocal condition number of B

        rcondb = rcond**2
        Write (nout, *)
        Write (nout, *) 'Estimate of reciprocal condition number for B'
        Write (nout, 110) rcondb
        Flush (nout)

!       Get the machine precision, eps, and if rcondb is not less
!       than eps**2, compute error estimates for the eigenvalues and
!       eigenvectors

        eps = epsilon(1.0E0_dp)
        If (rcond>=eps) Then

!         Call DDISNA to estimate reciprocal condition
!         numbers for the eigenvectors of (A*B - lambda*I)

          Call ddisna('Eigenvectors', n, n, w, rcondz, info)

!         Compute the error estimates for the eigenvalues and
!         eigenvectors

          t1 = one/rcond
          t2 = anorm*bnorm
          Do i = 1, n
            eerbnd(i) = (t2+abs(w(i))/rcondb)
            zerbnd(i) = t1*(t2/rcondz(i)+t1)
          End Do

!         Print the approximate error bounds for the eigenvalues
!         and vectors

          Write (nout, *)
          Write (nout, *) 'Error estimates (relative to machine precision)'
          Write (nout, *) 'for the eigenvalues:'
          Write (nout, 110) eerbnd(1:n)
          Write (nout, *)
          Write (nout, *) 'for the eigenvectors:'
          Write (nout, 110) zerbnd(1:n)
        Else
          Write (nout, *)
          Write (nout, *) 'B is very ill-conditioned, error ', &
            'estimates have not been computed'
        End If
      Else If (info>n) Then
        i = info - n
        Write (nout, 120) 'The leading minor of order ', i, &
          ' of B is not positive definite'
      Else
        Write (nout, 130) 'Failure in ZHEGVD. INFO =', info
      End If

100   Format (3X, (6F11.4))
110   Format (4X, 1P, 6E11.1)
120   Format (1X, A, I4, A)
130   Format (1X, A, I4)
    End Program


ご案内
関連情報
Privacy Policy  /  Trademarks