Class MoistAdjust
In: moist/moistadjust.f90

湿潤飽和調節法を行うためのパッケージ型モジュール

 * 飽和蒸気圧を使う化学種は, 一元的に実行する
 * 化学反応については, それぞれの化学反応毎に実行する.

Methods

Included Modules

dc_message gridset basicset moistset ChemCalc MoistFunc StorePotTemp

Public Instance methods

Subroutine :
xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax) :real(8),intent(in)
: エクスナー関数
xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) :real(8),intent(inout)
: 温位
xza_MixRt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8),intent(inout)
: 凝縮成分の混合比
DelTime :real(8),intent(in), optional

NH3 + H2S —> NH4SH の生成反応の圧平衡定数 Kp を用いた飽和湿潤調節法

[Source]

  subroutine MoistAdjustNH4SH(xz_Exner, xz_PotTemp, xza_MixRt, DelTime) 
    !
    ! NH3 + H2S --> NH4SH の生成反応の圧平衡定数 Kp を用いた飽和湿潤調節法
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8),intent(in)   :: xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax) 
                                                          !エクスナー関数
    real(8),intent(inout):: xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) 
                                                          !温位
    real(8),intent(inout):: xza_MixRt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
                                                          !凝縮成分の混合比
    real(8),intent(in), optional:: DelTime
    real(8) :: xz_PotTemp_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_PotTemp_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_PotTemp_org(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtNH3_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtNH3_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtH2S_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtH2S_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtNH4SH_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtNH4SH_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    
    real(8) :: xz_ExnerAll(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_PressAll(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_TempAll(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_Gamma(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_DelMixRt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_Cond(DimXMin:DimXMax, DimZMin:DimZMax)  
    real(8) :: xz_Evap(DimXMin:DimXMax, DimZMin:DimZMax)  
    
    integer            :: i
    integer, parameter :: ItrNum = 2
            
    !---------------------------------------------------------------------
    ! 初期化
    !---------------------------------------------------------------------
    xz_PotTemp_nxt    = 0.0d0
    xz_MixRtNH3_nxt   = 0.0d0
    xz_MixRtH2S_nxt   = 0.0d0
    xz_MixRtNH4SH_nxt = 0.0d0
    
    xz_Gamma     = 0.0d0
    xz_DelMixRt  = 0.0d0

    ! 温位の初期値を保管しておく
    xz_PotTemp_org = xz_PotTemp

    !湿潤飽和法では圧力は変化させない. 
    xz_ExnerAll = xz_Exner + xz_ExnerBasicZ
    xz_PressAll = PressBasis * (xz_ExnerAll ** (CpDry / GasRDry))
        
    !今までに得られた擾乱成分の値を暫定量とみなす. 添え字を追加
    xz_MixRtNH3_pre(:,:)    = xza_MixRt(:,:,IdxNH3)    + xza_MixRtBasicZ(:,:,IdxNH3)
    xz_MixRtH2S_pre(:,:)    = xza_MixRt(:,:,IdxH2S)    + xza_MixRtBasicZ(:,:,IdxH2S)
    xz_MixRtNH4SH_pre(:,:)  = xza_MixRt(:,:,IdxNH4SHc) + xza_MixRtBasicZ(:,:,IdxNH4SHc)
    xz_PotTemp_pre          = xz_PotTemp
    
    AdjustNH4SH: do i = 1, ItrNum
      !---------------------------------------------------------------
      ! 変数の初期化
      !---------------------------------------------------------------
      !温度
      xz_TempAll = ( xz_PotTemp_pre + xz_PotTempBasicZ ) * xz_ExnerAll
      
      !規格化された反応熱 (NH4SH 1kg に対する熱量)
      xz_Gamma = ReactHeatNH4SH / ( xz_ExnerAll * CpDry )

      !NH4SH の生成量
      xz_DelMixRt = xz_DelMixRtNH4SH( xz_TempAll, xz_PressAll, xz_MixRtNH3_pre, xz_MixRtH2S_pre, MolWtWet(IdxNH3), MolWtWet(IdxH2S) )
      
      xz_Cond = max( 0.0d0, xz_DelMixRt )
      xz_Evap = max( 0.0d0, min( - xz_DelMixRt, xz_MixRtNH4SH_pre ) )

      !---------------------------------------------------------------
      ! より真に近い値を求める飽和蒸気圧から飽和混合比を求める
      !---------------------------------------------------------------
      ! NH4SH の混合比を修正
      xz_MixRtNH4SH_nxt  = xz_MixRtNH4SH_pre + xz_Cond - xz_Evap
      
      ! DelPress を元に, NH3 と H2S の混合比を修正
      xz_MixRtNH3_nxt = xz_MixRtNH3_pre - ( xz_Cond - xz_Evap ) * MolWtWet(IdxNH3) / MolWtWet(IdxNH4SHc)
      xz_MixRtH2S_nxt = xz_MixRtH2S_pre - ( xz_Cond - xz_Evap ) * MolWtWet(IdxH2S) / MolWtWet(IdxNH4SHc)
          
      !温位を修正
      xz_PotTemp_nxt = xz_PotTemp_pre + xz_Gamma * ( xz_Cond - xz_Evap )
      
      !ループを回すための変数変化
      xz_PotTemp_pre    = xz_PotTemp_nxt
      xz_MixRtNH3_pre   = xz_MixRtNH3_nxt 
      xz_MixRtH2S_pre   = xz_MixRtH2S_nxt 
      xz_MixRtNH4SH_pre = xz_MixRtNH4SH_nxt 

    end do AdjustNH4SH

    xz_PotTemp              = xz_PotTemp_nxt                 
    xza_MixRt(:,:,IdxNH3)   = xz_MixRtNH3_nxt   - xza_MixRtBasicZ(:,:,IdxNH3)
    xza_MixRt(:,:,IdxH2S)   = xz_MixRtH2S_nxt   - xza_MixRtBasicZ(:,:,IdxH2S)
    xza_MixRt(:,:,IdxNH4SHc) = xz_MixRtNH4SH_nxt - xza_MixRtBasicZ(:,:,IdxNH4SHc)

    if (present(DelTime)) then 
      call StorePotTempCond3((xz_PotTemp - xz_PotTemp_org) / DelTime) 
    end if
    
  end subroutine MoistAdjustNH4SH
Subroutine :
xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax) :real(8),intent(in)
: エクスナー関数
xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) :real(8),intent(inout)
: 温位
xza_MixRt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum) :real(8),intent(inout)
: 混合比
DelTime :real(8),intent(in), optional

飽和蒸気圧を用いた湿潤飽和調整法の実行 この副プログラムでは, 予め決めた回数だけ反復改良を行う.

[Source]

  subroutine MoistAdjustSvapPress(xz_Exner, xz_PotTemp, xza_MixRt, DelTime)
    !
    ! 飽和蒸気圧を用いた湿潤飽和調整法の実行
    ! この副プログラムでは, 予め決めた回数だけ反復改良を行う. 
    !

    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8),intent(in)   :: xz_Exner(DimXMin:DimXMax, DimZMin:DimZMax) 
                                                          !エクスナー関数
    real(8),intent(inout):: xz_PotTemp(DimXMin:DimXMax, DimZMin:DimZMax) 
                                                          !温位
    real(8),intent(inout):: xza_MixRt(DimXMin:DimXMax, DimZMin:DimZMax, SpcNum)
                                                          !混合比
    real(8),intent(in), optional:: DelTime
    integer, parameter   :: ItrNum = 4                    !反復改良の回数

    real(8) :: xz_MixRtV_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtV_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtC_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtC_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_PotTemp_pre(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_PotTemp_nxt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_PotTemp_org(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_ExnerAll(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_TempAll(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_DelMixRt(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_MixRtSat(DimXMin:DimXMax, DimZMin:DimZMax) 
    real(8) :: xz_Cond(DimXMin:DimXMax, DimZMin:DimZMax)  
    real(8) :: xz_Evap(DimXMin:DimXMax, DimZMin:DimZMax)  
    real(8) :: xz_Gamma(DimXMin:DimXMax, DimZMin:DimZMax) 
    integer :: i, s                                       ! 添え字  

    !---------------------------------------------------------------------
    ! 初期化
    !---------------------------------------------------------------------
    xz_MixRtV_pre  = 0.0d0
    xz_MixRtV_nxt  = 0.0d0
    xz_MixRtC_pre  = 0.0d0
    xz_MixRtC_nxt  = 0.0d0
    xz_PotTemp_pre = 0.0d0
    xz_PotTemp_nxt = 0.0d0
    
    xz_ExnerAll    = 0.0d0
    xz_TempAll     = 0.0d0
    xz_DelMixRt    = 0.0d0
    xz_MixRtSat    = 0.0d0 
    xz_Gamma       = 0.0d0
    
    !---------------------------------------------------------------------
    ! 湿潤飽和調節法の実行
    !   ループを回すのは, 雲についてだけ.  
    !---------------------------------------------------------------------
    LoopSvapPress: do s = 1, CondNum
      !初期値を取っておく
      xz_PotTemp_org = xz_PotTemp 
      
      !湿潤飽和法では圧力は変化させない. 
      xz_ExnerAll = xz_Exner + xz_ExnerBasicZ
          
      !今までに得られた擾乱成分の値を暫定量とみなす. 添え字を追加
      xz_MixRtV_pre  = xza_MixRt(:,:,IdxCG(s)) + xza_MixRtBasicZ(:,:,IdxCG(s))
      xz_MixRtC_pre  = xza_MixRt(:,:,IdxCC(s)) + xza_MixRtBasicZ(:,:,IdxCC(s))
      xz_PotTemp_pre = xz_PotTemp
          
      Adjusting: do i = 1, ItrNum
        !---------------------------------------------------------------
        ! 飽和蒸気圧から飽和混合比を求める
        !---------------------------------------------------------------
        !温度
        xz_TempAll = ( xz_PotTemp_pre + xz_PotTempBasicZ ) * xz_ExnerAll

        !飽和蒸気圧から飽和混合比を計算(基本場からの差). 
        xz_MixRtSat = xz_SvapPress(SpcWetID(IdxCC(s)), xz_TempAll) * MolWtWet(IdxCC(s)) / (MolWtDry * PressBasis * (xz_ExnerAll ** (CpDry / GasRDry)) )

        !規格化された潜熱
        xz_Gamma = xz_LatentHeat(SpcWetID(IdxCC(s)), xz_TempAll) / (xz_ExnerAll * CpDry)
        
        !凝結量を求める. 
        !  凝結が生じる場合には, xz_MixRtV_pre - xz_MixRtSat は必ず正の値となる.
        !  蒸発が生じる場合には, 蒸発量は - MixRtC を超えることはない. 
        xz_DelMixRt = ( xz_MixRtV_pre - xz_MixRtSat ) / (1.0d0 + xz_Gamma * xz_DMixRtSatDPotTemp( SpcWetID(IdxCC(s)), MolWtWet(IdxCC(s)), xz_TempAll, xz_ExnerAll ) )
        
        xz_Cond = max( 0.0d0, min( xz_MixRtV_pre,   xz_DelMixRt ) )
        xz_Evap = max( 0.0d0, min( xz_MixRtC_pre, - xz_DelMixRt ) )

        !より真に近い値を計算
        xz_PotTemp_nxt = xz_PotTemp_pre + xz_Gamma * ( xz_Cond - xz_Evap )
        xz_MixRtV_nxt  = xz_MixRtV_pre  - xz_Cond + xz_Evap
        xz_MixRtC_nxt  = xz_MixRtC_pre  + xz_Cond - xz_Evap
        
        !繰り返しのための変数定義
        xz_PotTemp_pre = xz_PotTemp_nxt
        xz_MixRtV_pre  = xz_MixRtV_nxt
        xz_MixRtC_pre  = xz_MixRtC_nxt
        
      end do Adjusting
      
      xz_PotTemp              = xz_PotTemp_nxt                 
      xza_MixRt(:,:,IdxCG(s)) = xz_MixRtV_nxt - xza_MixRtBasicZ(:,:,IdxCG(s))
      xza_MixRt(:,:,IdxCC(s)) = xz_MixRtC_nxt - xza_MixRtBasicZ(:,:,IdxCC(s))
      
      if (present(DelTime)) then 
        call StorePotTempCond2((xz_PotTemp - xz_PotTemp_org) / DelTime, idxCG(s)) 
      end if

    end do LoopSvapPress 
  
  end subroutine MoistAdjustSvapPress