計算ルーチン: DGTTRFで LU 分解済みの結果を用いて実三重対角連立方程式を解く

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

ホーム > LAPACKサンプルプログラム目次 > 計算ルーチン > DGTTRFで LU 分解済みの結果を用いて実三重対角連立方程式を解く

概要

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

入力データ

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

このデータをダウンロード
DGTTRS Example Program Data
  5     2                      :Values of N and NRHS
        2.1  -1.0   1.9   8.0
  3.0   2.3  -5.0  -0.9   7.1
  3.4   3.6   7.0  -6.0        :End of matrix A
  2.7   6.6
 -0.5  10.8
  2.6  -3.2
  0.6 -11.2
  2.7  19.1                    :End of matrix B

出力結果

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

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

 Solution(s)
             1          2
 1     -4.0000     5.0000
 2      7.0000    -4.0000
 3      3.0000    -3.0000
 4     -4.0000    -2.0000
 5     -3.0000     1.0000

ソースコード

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

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


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

!     DGTTRS 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: dgttrf, dgttrs
      Use lapack_precision, Only: dp
!     .. Implicit None Statement ..
      Implicit None
!     .. Parameters ..
      Integer, Parameter :: nin = 5, nout = 6
!     .. Local Scalars ..
      Integer :: i, ifail, info, ldb, n, nrhs
!     .. Local Arrays ..
      Real (Kind=dp), Allocatable :: b(:, :), d(:), dl(:), du(:), du2(:)
      Integer, Allocatable :: ipiv(:)
!     .. Executable Statements ..
      Write (nout, *) 'DGTTRS Example Program Results'
      Write (nout, *)
      Flush (nout)
!     Skip heading in data file
      Read (nin, *)
      Read (nin, *) n, nrhs
      ldb = n
      Allocate (b(ldb,nrhs), d(n), dl(n-1), du(n-1), du2(n-2), ipiv(n))

!     Read the tridiagonal matrix A from data file

      Read (nin, *) du(1:n-1)
      Read (nin, *) d(1:n)
      Read (nin, *) dl(1:n-1)

!     Read the right hand matrix B
      Read (nin, *)(b(i,1:nrhs), i=1, n)

!     Factorize the tridiagonal matrix A
      Call dgttrf(n, dl, d, du, du2, ipiv, info)

      If (info==0) Then

!       Solve the equations AX = B
        Call dgttrs('No transpose', n, nrhs, dl, d, du, du2, ipiv, b, ldb, &
          info)

!       Print the solution

!       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, nrhs, b, ldb, &
          'Solution(s)', ifail)

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

100   Format (1X, A, I3, A, I3, A, A)
    End Program


ご案内
関連情報
Privacy Policy  /  Trademarks