スパース固有値問題

Fortranによるサンプルソースコード : 使用ルーチン名:f12acf

Keyword: スパース, 固有値問題

概要

本サンプルはスパース固有値問題を解くFortranによるサンプルプログラムです。 本サンプルは以下に示されるスパース固有値問題を解いて固有値を出力します。本サンプルプログラムでは、関数 f12acf のほか、f12aaf、f12abf、f12adf や f12aef を呼び出しています。

スパース固有値問題のデータ 

※本サンプルはnAG Fortranライブラリに含まれるルーチン f12acf() のExampleコードです。本サンプル及びルーチンの詳細情報は f12acf のマニュアルページをご参照ください。
ご相談やお問い合わせはこちらまで

入力データ

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

このデータをダウンロード
F12ACF Example Program Data
 10 4 20 10.0 : Values for NX NEV NCV RHO

  • 1行目はタイトル行で読み飛ばされます。
  • 2行目に正方行列の行(列)の数(nx)、固有値の数(nev)、計算時に使用するためのArnoldi 規定ベクトルの数(ncv)、ρ(ロー)の値(rho)を指定しています。

出力結果

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

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


 The    4 generalized Ritz values of largest magnitude are:

        1     (   20383.0384 ,       0.0000 )
        2     (   20338.7563 ,       0.0000 )
        3     (   20265.2844 ,       0.0000 )
        4     (   20163.1142 ,       0.0000 )

  • 4行目には収束した固有値が4つあることが示されています。
  • 6〜9行目に収束した近似固有値の実部と虚部が出力されています。

ソースコード

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

※本サンプルソースコードは科学技術・統計計算ライブラリである「nAG Fortranライブラリ」のルーチンを呼び出します。
サンプルのコンパイル及び実行方法


このソースコードをダウンロード
!   F12ACF Example Program Text
!   Mark 23 Release. nAG Copyright 2011.

    MODULE f12acfe_mod

!      F12ACF Example Program Module:
!             Parameters and User-defined Routines

!      .. Use Statements ..
       USE nag_library, ONLY : nag_wp
!      .. Implicit None Statement ..
       IMPLICIT NONE
!      .. Parameters ..
       REAL (KIND=nag_wp), PARAMETER       :: one = 1.0_nag_wp
       INTEGER, PARAMETER                  :: imon = 0, nin = 5, nout = 6
    CONTAINS
       SUBROUTINE av(nx,rho,v,w)

!         .. Implicit None Statement ..
          IMPLICIT NONE
!         .. Parameters ..
          REAL (KIND=nag_wp), PARAMETER       :: two = 2.0_nag_wp
!         .. Scalar Arguments ..
          REAL (KIND=nag_wp), INTENT (IN)     :: rho
          INTEGER, INTENT (IN)                :: nx
!         .. Array Arguments ..
          REAL (KIND=nag_wp), INTENT (IN)     :: v(nx*nx)
          REAL (KIND=nag_wp), INTENT (OUT)    :: w(nx*nx)
!         .. Local Scalars ..
          REAL (KIND=nag_wp)                  :: dd, dl, du, h, s
          INTEGER                             :: j, n
!         .. Intrinsic Functions ..
          INTRINSIC                              real
!         .. Executable Statements ..
          n = nx*nx
          h = one/real(n+1,kind=nag_wp)
          s = rho/two
          dd = two/h
          dl = -one/h - s
          du = -one/h + s
          w(1) = dd*v(1) + du*v(2)
          DO j = 2, n - 1
             w(j) = dl*v(j-1) + dd*v(j) + du*v(j+1)
          END DO
          w(n) = dl*v(n-1) + dd*v(n)
          RETURN
       END SUBROUTINE av

       SUBROUTINE mv(nx,v,w)

!         .. Use Statements ..
          USE nag_library, ONLY : dscal
!         .. Implicit None Statement ..
          IMPLICIT NONE
!         .. Parameters ..
          REAL (KIND=nag_wp), PARAMETER       :: four = 4.0_nag_wp
!         .. Scalar Arguments ..
          INTEGER, INTENT (IN)                :: nx
!         .. Array Arguments ..
          REAL (KIND=nag_wp), INTENT (IN)     :: v(nx*nx)
          REAL (KIND=nag_wp), INTENT (OUT)    :: w(nx*nx)
!         .. Local Scalars ..
          REAL (KIND=nag_wp)                  :: h
          INTEGER                             :: j, n
!         .. Intrinsic Functions ..
          INTRINSIC                              real
!         .. Executable Statements ..
          n = nx*nx
          w(1) = four*v(1) + one*v(2)
          DO j = 2, n - 1
             w(j) = one*v(j-1) + four*v(j) + one*v(j+1)
          END DO
          w(n) = one*v(n-1) + four*v(n)
          h = one/real(n+1,kind=nag_wp)
!         The nAG name equivalent of dscal is f06edf
          CALL dscal(n,h,w,1)
          RETURN
       END SUBROUTINE mv
    END MODULE f12acfe_mod

    PROGRAM f12acfe

!      F12ACF Example Main Program

!      .. Use Statements ..
       USE nag_library, ONLY : dnrm2, dpttrf, dpttrs, f12aaf, f12abf, f12acf,  &
                               f12adf, f12aef
       USE f12acfe_mod, ONLY : av, imon, mv, nag_wp, nin, nout, one
!      .. Implicit None Statement ..
       IMPLICIT NONE
!      .. Local Scalars ..
       REAL (KIND=nag_wp)                  :: h, rho, sigmai, sigmar
       INTEGER                             :: ifail, ifail1, info, irevcm, j,  &
                                              lcomm, ldv, licomm, n, nconv,    &
                                              ncv, nev, niter, nshift, nx
!      .. Local Arrays ..
       REAL (KIND=nag_wp), ALLOCATABLE     :: comm(:), d(:,:), md(:), me(:),   &
                                              mx(:), resid(:), v(:,:), x(:)
       INTEGER, ALLOCATABLE                :: icomm(:)
!      .. Intrinsic Functions ..
       INTRINSIC                              real
!      .. Executable Statements ..
       WRITE (nout,*) 'F12ACF Example Program Results'
       WRITE (nout,*)
!      Skip heading in data file
       READ (nin,*)
       READ (nin,*) nx, nev, ncv, rho
       n = nx*nx
       ldv = n
       licomm = 140
       lcomm = 3*n + 3*ncv*ncv + 6*ncv + 60
       ALLOCATE (comm(lcomm),d(ncv,3),md(n),me(n-1),mx(n),resid(n),v(ldv,ncv), &
          x(n),icomm(licomm))

!      ifail: behaviour on error exit
!             =0 for hard exit, =1 for quiet-soft, =-1 for noisy-soft
       ifail = 0
       CALL f12aaf(n,nev,ncv,icomm,licomm,comm,lcomm,ifail)

!      Set the mode.
       ifail = 0
       CALL f12adf('REGULAR INVERSE',icomm,comm,ifail)
!      Set problem type.
       CALL f12adf('GENERALIZED',icomm,comm,ifail)
!      Use pointers to Workspace in calculating matrix vector
!      products rather than interfacing through the array X
       CALL f12adf('POINTERS=YES',icomm,comm,ifail)

!      Construct M, and factorize using DPTTRF/F07JDF.
       h = one/real(n+1,kind=nag_wp)
       md(1:n-1) = 4.0_nag_wp*h
       me(1:n-1) = h
       md(n) = 4.0_nag_wp*h

!      The nAG name equivalent of dpttrf is f07jdf
       CALL dpttrf(n,md,me,info)

       irevcm = 0

       ifail = -1
LOOP:  DO
          CALL f12abf(irevcm,resid,v,ldv,x,mx,nshift,comm,icomm,ifail)

          IF (irevcm/=5) THEN
             SELECT CASE (irevcm)
             CASE (-1,1)
!               Perform  y <--- OP*x = inv[M]*A*x using DPTTRS/F07JEF.
                CALL av(nx,rho,comm(icomm(1)),comm(icomm(2)))
!               The nAG name equivalent of dpttrs is f07jef
                CALL dpttrs(n,1,md,me,comm(icomm(2)),n,info)
             CASE (2)
!               Perform  y <--- M*x.
                CALL mv(nx,comm(icomm(1)),comm(icomm(2)))
             CASE (4)
                IF (imon/=0) THEN
!                  Output monitoring information if required.
                   CALL f12aef(niter,nconv,d,d(1,2),d(1,3),icomm,comm)
!                  The nAG name equivalent of dnrm2 is f06ejf
                   WRITE (6,99999) niter, nconv, dnrm2(nev,d(1,3),1)
                END IF
             END SELECT
          ELSE
             EXIT LOOP
          END IF
       END DO LOOP
       IF (ifail==0) THEN
!         Post-Process using F12ACF to compute eigenvalues/vectors.
          ifail1 = 0
          CALL f12acf(nconv,d,d(1,2),v,ldv,sigmar,sigmai,resid,v,ldv,comm, &
             icomm,ifail1)
!         Print computed eigenvalues.
          WRITE (nout,99998) nconv
          DO j = 1, nconv
             WRITE (nout,99997) j, d(j,1), d(j,2)
          END DO
       END IF

99999  FORMAT (1X,'Iteration',1X,I3,', No. converged =',1X,I3,', norm o', &
          'f estimates =',E16.8)
99998  FORMAT (1X/' The ',I4,' generalized Ritz values of largest ', &
          'magnitude are:'/)
99997  FORMAT (1X,I8,5X,'( ',F12.4,' , ',F12.4,' )')
    END PROGRAM f12acfe


関連情報
Privacy Policy  /  Trademarks