実線形方程式: 一般行列 : (右辺はベクトル)

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

概要

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

:6x4の行列の特異値と左および右特異ベクトルを求めます。

\begin{displaymath}
A = \left(
\begin{array}{rrrr}
2.27 & -1.54 & 1.15 & -1.9...
... & 2.30 \\
0.62 & -7.39 & 1.03 & -2.57
\end{array} \right),
\end{displaymath}

計算された特異値と特異ベクトルの誤差限界近似値も合わせて求めます。

DGESDDの例題プログラムは $ m \le n$の場合の特異値分解を示します。

入力データ

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

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

  4                           :Value of N

  1.80   2.88   2.05  -0.89
  5.25  -2.95  -0.95  -3.80
  1.58  -2.69  -2.90  -1.04
 -1.11  -0.66  -0.59   0.80   :End of matrix A

  9.52  24.35   0.77  -6.22   :End of vector b

出力結果

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

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

 Solution
        1.0000    -1.0000     3.0000    -5.0000

 Details of factorization
             1          2          3          4
 1      5.2500    -2.9500    -0.9500    -3.8000
 2      0.3429     3.8914     2.3757     0.4129
 3      0.3010    -0.4631    -1.5139     0.2948
 4     -0.2114    -0.3299     0.0047     0.1314

 Pivot indices
             2          2          3          4

ソースコード

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

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


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

!     DGESV Example Program Text

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

!     .. Use Statements ..
      Use lapack_example_aux, Only: nagf_file_print_matrix_real_gen
      Use lapack_interfaces, Only: dgesv
      Use lapack_precision, Only: dp
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Integer, Parameter :: nin = 5, nout = 6
!     .. Local Scalars ..
      Integer :: i, ifail, info, lda, ldb, n
!     .. Local Arrays ..
      Real (Kind=dp), Allocatable :: a(:, :), b(:)
      Integer, Allocatable :: ipiv(:)
!     .. Executable Statements ..
      Write (nout, *) 'DGESV 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), ipiv(n))

!     Read A and B from data file

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

!     Solve the equations Ax = b for x

      Call dgesv(n, 1, a, lda, ipiv, b, ldb, info)

      If (info==0) Then

!       Print solution

        Write (nout, *) 'Solution'
        Write (nout, 100) b(1:n)

!       Print details of factorization

        Write (nout, *)
        Flush (nout)

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

!       Print pivot indices

        Write (nout, *)
        Write (nout, *) 'Pivot indices'
        Write (nout, 110) ipiv(1:n)

      Else
        Write (nout, 120) 'The (', info, ',', info, ')', &
          ' element of the factor U is zero'
      End If

100   Format ((3X,7F11.4))
110   Format ((3X,7I11))
120   Format (1X, A, I3, A, I3, A, A)
    End Program


ご案内
関連情報
Privacy Policy  /  Trademarks