Class DynImpFunc_3d
In: dynamic/dynimpfunc_3d.f90

陰解法を用いた力学過程の各項の計算モジュール. エクスナー関数を陰解法で解く際に必要となる, 係数行列の要素の決定,

本モジュールでは科学計算ライブラリとして, 以下をサポートする.

 * MATRIX/MPP (HITACH 最適化 Fortran)
 * LAPACK

Methods

Included Modules

dc_types dc_iounit dc_message gridset_3d timeset damping_3d basicset_3d xyz_module xyz_deriv_module

Public Instance methods

Function :
xyr_GradPi(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP)
: 圧力傾度力
xyz_ExnerA(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: エクスナー関数の擾乱
xyz_ExnerN(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: エクスナー関数の擾乱
pyz_VelX(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 水平速度
xqz_VelY(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 水平速度
xyr_VelZ(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 鉛直速度

x 方向に半格子ずれた点での圧力傾度力項の計算. クランク・ニコルソン法を用いるために, 時刻 τ と τ+Δτでの エクスナー関数の値を引数として取る. 音波減衰項を含めた形式で定式化してあることに注意.

[Source]

  function xyr_GradPi(xyz_ExnerA, xyz_ExnerN, pyz_VelX, xqz_VelY, xyr_VelZ)
    !
    ! x 方向に半格子ずれた点での圧力傾度力項の計算. 
    ! クランク・ニコルソン法を用いるために, 時刻 τ と τ+Δτでの
    ! エクスナー関数の値を引数として取る.
    ! 音波減衰項を含めた形式で定式化してあることに注意.
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(DP), intent(in)  :: xyz_ExnerA (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !エクスナー関数の擾乱
    real(DP), intent(in)  :: xyz_ExnerN (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !エクスナー関数の擾乱
    real(DP), intent(in)  :: pyz_VelX (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !水平速度
    real(DP), intent(in)  :: xqz_VelY (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !水平速度
    real(DP), intent(in)  :: xyr_VelZ (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !鉛直速度
    real(DP)              :: xyr_GradPi (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !圧力傾度力
    real(DP)              :: xyz_DivVel (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
                                               !速度の収束

    !速度の収束
    xyz_DivVel =  xyz_dx_pyz( pyz_VelX ) + xyz_dy_xqz( xqz_VelY ) + xyz_dz_xyr( xyr_VelZ )
    
    !速度 w の圧力勾配
!    xyr_GradPi = 0.0d0
    xyr_GradPi = xyr_avr_xyz(CpDry * xyz_PotTempBasicZ ) * ( beta * xyr_dz_xyz( xyz_ExnerA ) + (1.0d0 - beta) * xyr_dz_xyz( xyz_ExnerN ) - xyr_dz_xyz( DampSound * xyz_DivVel ) )                                              
    
  end function xyr_GradPi
Function :
xyz_Exner(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP)
: 無次元圧力[τ+Δτ]
xyr_FzNl(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: Z 方向の外力項[t]
pyz_VelXNs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 速度 u [τ]
pyz_VelXAs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 速度 u [τ+Δτ]
xqz_VelYNs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 速度 v [τ]
xqz_VelYAs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 速度 v [τ+Δτ]
xyr_VelZNs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 速度 w [τ]
xyz_ExnerNs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(DP), intent(in)
: 無次元圧力

陰解法を用いたエクスナー関数の計算.

[Source]

  function xyz_Exner(xyr_FzNl, pyz_VelXNs, pyz_VelXAs, xqz_VelYNs, xqz_VelYAs, xyr_VelZNs, xyz_ExnerNs)
    !
    !陰解法を用いたエクスナー関数の計算. 
    !

    !暗黙の型宣言禁止
    implicit none
    
   !入出力変数
    real(DP), intent(in)   :: pyz_VelXNs (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !速度 u [τ]
    real(DP), intent(in)   :: pyz_VelXAs (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !速度 u [τ+Δτ]
    real(DP), intent(in)   :: xqz_VelYNs (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !速度 v [τ]
    real(DP), intent(in)   :: xqz_VelYAs (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !速度 v [τ+Δτ]
    real(DP), intent(in)   :: xyr_VelZNs (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !速度 w [τ]
    real(DP), intent(in)   :: xyr_FzNl (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !Z 方向の外力項[t]
    real(DP), intent(in)   :: xyz_ExnerNs (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !無次元圧力
    real(DP)               :: xyz_Exner (DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
                                                           !無次元圧力[τ+Δτ]

    !変数定義
    real(DP)               :: D(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    real(DP)               :: E(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    real(DP)               :: F(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    real(DP)               :: xyz_DivVelNs(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)
    real(DP)               :: X(RegXMin:RegXMax, RegZMin:RegZMax)

    real(DP)               :: DTS ! 短い時間格子間隔
    integer                :: jy

    DTS = DelTimeShort


    !変数の初期化
    xyz_Exner = 0.0d0
    
    !速度の収束を計算
    xyz_DivVelNs =  xyz_dx_pyz( pyz_VelXNs ) + xyz_dy_xqz( xqz_VelYNs ) + xyz_dz_xyr( xyr_VelZNs )
  
    !行列計算のための係数
    E = xyr_dz_xyz( DampSound * xyz_DivVelNs ) - ( 1.0d0 - beta ) * xyr_dz_xyz( xyz_ExnerNs ) + xyr_FzNl / xyr_avr_xyz( CpDry * xyz_VPotTempBasicZ ) 
    
    F = - beta * xyz_F1BasicZ * DTS * xyz_dz_xyr( xyr_avr_xyz( xyz_DensBasicZ * xyz_VPotTempBasicZ) * ( xyr_VelZNs - xyr_avr_xyz(CpDry * xyz_VPotTempBasicZ) * DTS * ( (1.0d0 - beta) * xyr_dz_xyz( xyz_ExnerNs ) - xyr_dz_xyz( DampSound * xyz_DivVelNs ) ) + xyr_FzNl * DTS ) ) 
    
    D =   xyz_ExnerNs - (1.0d0 - beta) * xyz_F1BasicZ * DTS * xyz_dz_xyr( xyr_avr_xyz(xyz_DensBasicZ * xyz_VPotTempBasicZ) * xyr_VelZNs ) - (xyz_VelSoundBasicZ ** 2.0d0) *DTS / (CpDry * xyz_VPotTempBasicZ) * xyz_dx_pyz( pyz_VelXAs ) - (xyz_VelSoundBasicZ ** 2.0d0) *DTS / (CpDry * xyz_VPotTempBasicZ) * xyz_dy_xqz( xqz_VelYAs ) + F
    
    D(:,:,RegZMin) = D(:,:,RegZMin) - beta * xyz_F1BasicZ(:,:,RegZMin) * (DTS ** 2.0d0) * xyr_F2BasicZ(:,:, RegZMin-1) * E(:,:,RegZMin-1) / z_dz(RegZMin)
    
    D(:,:,RegZMax) = D(:,:,RegZMax) + beta * xyz_F1BasicZ(:,:,RegZMax) * (DTS ** 2.0d0) * xyr_F2BasicZ(:,:,RegZMax) * E(:,:,RegZMax) / z_dz(RegZMax)     

    !-----------------------------------------------------------
    !連立一次方程式の解を求める
    !------------------------------------------------------------
    !配列の準備. 
    !  * D を求めるためには, 微分平均演算で領域の「糊代の部分」を
    !    利用するので, この段階で必要な部分を切り出す

    do jy = RegYMin, RegYMax

      X = D(RegXMin:RegXMax,jy,RegZMin:RegZMax) 

      !解の計算
      !  LAPACK 利用
      call LinSolv_Lapack( X )

      !戻り値を出力
      xyz_Exner(RegXMin:RegXMax,jy,RegZMin:RegZMax) = X
    
    end do
   
   
  end function xyz_Exner
Subroutine :

エクスナー関数を陰解法で解く際に必要となる, 係数行列の要素を決め, LU 分解を行う.

[Source]

  subroutine xyz_Exner_init()
    !
    !エクスナー関数を陰解法で解く際に必要となる, 係数行列の要素を決め, 
    !LU 分解を行う. 
    !

    !暗黙の型宣言禁止
    implicit none

    real(DP)  :: DTS ! 短い時間格子

    DTS = DelTimeShort


    !配列の割り付け
    allocate( A(RegZMin:RegZMax), B(RegZMin+1:RegZMax), C(RegZMin:RegZMax-1), xyz_F1BasicZ(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax), xyr_F2BasicZ(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax), xyz_VPotTempBasicZ(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) )

    !----------------------------------------------------------------
    ! 係数行列と共通して利用される配列の値を決める
    !----------------------------------------------------------------
    !仮温位の定義
    xyz_VPotTempBasicZ = xyz_PotTempBasicZ 

    !係数行列の計算
    !  A, B, C を求める際, F1BasicZ と F2BasicZ は X 方向に一様なので. 
    !  RegXMax, RegYMax の値で代表させることとした. 

    xyz_F1BasicZ = ( xyz_VelSoundBasicZ ** 2.0d0 ) / (CpDry * xyz_DensBasicZ * (xyz_VPotTempBasicZ ** 2.0d0))

    xyr_F2BasicZ = xyr_avr_xyz( CpDry * xyz_DensBasicZ * ( xyz_VPotTempBasicZ ** 2.0d0 ) )
    
    A(RegZMin+1: RegZMax-1) = (beta ** 2.0d0) * xyz_F1BasicZ(RegXMax,RegYMax,RegZMin+1: RegZMax-1) * ( xyr_F2BasicZ(RegXMax,RegYMax,RegZMin+1: RegZMax-1) / r_dz(RegZMin+1: RegZMax-1) + xyr_F2BasicZ(RegXMax,RegYMax,RegZMin  : RegZMax-2) / r_dz(RegZMin: RegZMax-2) ) * (DTS ** 2.0d0) / z_dz(RegZMin+1: RegZMax-1) + 1.0d0

    A(RegZMin) = (beta ** 2.0d0) * xyz_F1BasicZ(RegXMax,RegYMax,RegZMin) * xyr_F2BasicZ(RegXMax,RegYMax,RegZMin) / r_dz(RegZMin) * (DTS ** 2.0d0) / z_dz(RegZMin) + 1.0d0                                         

    A(RegZMax) = (beta ** 2.0d0) * xyz_F1BasicZ(RegXMax,RegYMax,RegZMax) * xyr_F2BasicZ(RegXMax,RegYMax,RegZMax-1) / r_dz(RegZMax-1) * (DTS ** 2.0d0) / z_dz(RegZMax) + 1.0d0                                         
    
    B(RegZMin+1:RegZMax) = - (beta ** 2.0d0) * xyz_F1BasicZ(RegXMax,RegYMax,RegZMin:RegZMax-1) * xyr_F2BasicZ(RegXMax,RegYMax,RegZMin:RegZMax-1) / r_dz(RegZmin:RegZMax-1) * (DTS ** 2.0d0) / z_dz(RegZMin:RegZMax-1)
    
    C(RegZMin: RegZMax-1) = - ( beta ** 2.0d0 ) * xyz_F1BasicZ(RegXMax,RegYMax,RegZMin+1:RegZMax) * xyr_F2BasicZ(RegXMax,RegYMax,RegZMin  :RegZMax-1) / r_dz(RegZmin:RegZMax-1) * (DTS ** 2.0d0) / z_dz(RegZMin+1:RegZMax)


    !デバッグ出力ファイルオープン. 情報取得. 
!    call FileOpen(unit, file="dynimpfunc_log", mode='w')
!    write(unit,*) "kz, A_k"
!    write(unit,*) "------------------------------------"
!    do kz = RegZmin, RegZmax
!      write(unit,*) kz, A(kz)
!    end do
!    write(unit,*) "------------------------------------"
!    write(unit,*) "kz, B_k"
!    write(unit,*) "------------------------------------"
!    do kz = RegZmin+1, RegZmax
!      write(unit,*) kz, B(kz)
!    end do
!    write(unit,*) "------------------------------------"
!    write(unit,*) "kz, C_k"
!    write(unit,*) "------------------------------------"
!    do kz = RegZmin, RegZmax-1
!      write(unit,*) kz, C(kz)
!    end do
!    write(unit,*) "------------------------------------"
!    close(unit)


    !----------------------------------------------------------------
    ! 係数行列を LU 分解
    !----------------------------------------------------------------
    !配列の大きさを保管
    N   = RegZMax - RegZMin +1 !係数行列/改行列の次数, 整合寸法
    M   = RegXMax - RegXMin +1 !方程式の組数
    NUD = 1                    !係数行列の上三角部分の帯幅
    NLD = 1                    !係数行列の下三角部分の帯幅
    NAL = NLD                  !LU 分解の結果 L の整合寸法
    NA  = NUD + NLD + 1

    !配列の割り当て
    allocate( AL1(N), AL2(NAL, N), AU2(NA, N), IP(N) )

    !LU 分解の実行
    !  LAPACK の利用
    call ResolvLU_Lapack( )

   
  end subroutine xyz_Exner_init

[Validate]