Class MoistAdjust_3d
In: moist/moistadjust_3d.f90

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

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

Methods

Included Modules

dc_types dc_message gridset_3d basicset_3d moistset ChemCalc_3d MoistFunc_3d

Public Instance methods

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

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

[Source]

  subroutine MoistAdjustNH4SH(xyz_Exner, xyz_PotTemp, xyza_MixRt ) 
    !
    ! NH3 + H2S --> NH4SH の生成反応の圧平衡定数 Kp を用いた飽和湿潤調節法
    !
    
    !暗黙の型宣言禁止
    implicit none
    
    !変数定義
    real(8),intent(in)   :: xyz_Exner(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
                                                          !エクスナー関数
    real(8),intent(inout):: xyz_PotTemp(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
                                                          !温位
    real(8),intent(inout):: xyza_MixRt(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax, SpcNum)
                                                          !凝縮成分の混合比
    real(DP):: xyz_PotTemp_pre(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
    real(DP):: xyz_PotTemp_nxt(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtNH3_pre(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtNH3_nxt(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtH2S_pre(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtH2S_nxt(DimXMin:DimXMax,DimXMin:DimXMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtNH4SH_pre(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtNH4SH_nxt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    
    real(DP):: xyz_ExnerAll(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_PressAll(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_TempAll(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_Gamma(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_DelMixRt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_Cond(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    real(DP)::xyz_Evap(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    
    integer            :: i
    integer, parameter :: ItrNum = 2
            
    !---------------------------------------------------------------------
    ! 初期化
    !---------------------------------------------------------------------
    xyz_PotTemp_nxt = 0.0d0
    xyz_MixRtNH3_nxt   = 0.0d0
    xyz_MixRtH2S_nxt   = 0.0d0
    xyz_MixRtNH4SH_nxt = 0.0d0
    
    xyz_Gamma     = 0.0d0
    xyz_DelMixRt  = 0.0d0

    !湿潤飽和法では圧力は変化させない. 
    xyz_ExnerAll = xyz_Exner + xyz_ExnerBasicZ
    xyz_PressAll = PressBasis * (xyz_ExnerAll ** (CpDry / GasRDry))
        
    !今までに得られた擾乱成分の値を暫定量とみなす. 添え字を追加
    xyz_MixRtNH3_pre(:,:,:)   = xyza_MixRt(:,:,:,IdxNH3) + xyza_MixRtBasicZ(:,:,:,IdxNH3)
    xyz_MixRtH2S_pre(:,:,:)   = xyza_MixRt(:,:,:,IdxH2S) + xyza_MixRtBasicZ(:,:,:,IdxH2S)
    xyz_MixRtNH4SH_pre(:,:,:) = xyza_MixRt(:,:,:,IdxNH4SHc) + xyza_MixRtBasicZ(:,:,:,IdxNH4SHc)
    xyz_PotTemp_pre          = xyz_PotTemp
    
    AdjustNH4SH: do i = 1, ItrNum
      !---------------------------------------------------------------
      ! 変数の初期化
      !---------------------------------------------------------------
      !温度
      xyz_TempAll = ( xyz_PotTemp_pre + xyz_PotTempBasicZ ) * xyz_ExnerAll
      
      !規格化された反応熱 (NH4SH 1kg に対する熱量)
      xyz_Gamma = ReactHeatNH4SH * xyz_EffMolWtBasicZ / ( xyz_ExnerAll * CpDry )

      !NH4SH の生成量
      xyz_DelMixRt = xyz_DelMixRtNH4SH( xyz_TempAll, xyz_PressAll, xyz_MixRtNH3_pre, xyz_MixRtH2S_pre, MolWtWet(IdxNH3), MolWtWet(IdxH2S) )

      xyz_Cond = max( 0.0d0, xyz_DelMixRt )
      xyz_Evap = max( 0.0d0, min( - xyz_DelMixRt, xyz_MixRtNH4SH_pre ) )

      !---------------------------------------------------------------
      ! より真に近い値を求める飽和蒸気圧から飽和混合比を求める
      !---------------------------------------------------------------
      ! NH4SH の混合比を修正
      xyz_MixRtNH4SH_nxt  = xyz_MixRtNH4SH_pre + xyz_Cond - xyz_Evap
      
      ! DelPress を元に, NH3 と H2S の混合比を修正
      xyz_MixRtNH3_nxt = xyz_MixRtNH3_pre - ( xyz_Cond - xyz_Evap ) * MolWtWet(IdxNH3) / MolWtWet(IdxNH4SHc)
      xyz_MixRtH2S_nxt = xyz_MixRtH2S_pre - ( xyz_Cond - xyz_Evap ) * MolWtWet(IdxH2S) / MolWtWet(IdxNH4SHc)
          
      !温位を修正
      xyz_PotTemp_nxt = xyz_PotTemp_pre + xyz_Gamma * ( xyz_Cond - xyz_Evap )
      
      !ループを回すための変数変化
      xyz_PotTemp_pre    = xyz_PotTemp_nxt
      xyz_MixRtNH3_pre   = xyz_MixRtNH3_nxt 
      xyz_MixRtH2S_pre   = xyz_MixRtH2S_nxt 
      xyz_MixRtNH4SH_pre = xyz_MixRtNH4SH_nxt 

    end do AdjustNH4SH

    xyz_PotTemp              = xyz_PotTemp_nxt                 
    xyza_MixRt(:,:,:,IdxNH3)   = xyz_MixRtNH3_nxt   - xyza_MixRtBasicZ(:,:,:,IdxNH3)
    xyza_MixRt(:,:,:,IdxH2S)   = xyz_MixRtH2S_nxt   - xyza_MixRtBasicZ(:,:,:,IdxH2S)
    xyza_MixRt(:,:,:,IdxNH4SHc) = xyz_MixRtNH4SH_nxt - xyza_MixRtBasicZ(:,:,:,IdxNH4SHc)
    
  end subroutine MoistAdjustNH4SH
Subroutine :
xyz_Exner(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(8),intent(in)
: エクスナー関数
xyz_PotTemp(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) :real(8),intent(inout)
: 温位
xyza_MixRt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax, SpcNum) :real(8),intent(inout)
: 混合比

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

[Source]

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

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

    real(DP):: xyz_MixRtV_pre(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtV_nxt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtC_pre(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtC_nxt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_PotTemp_pre(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_PotTemp_nxt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_ExnerAll(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_TempAll(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_DelMixRt(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_MixRtSat(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    real(DP):: xyz_Cond(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    real(DP):: xyz_Evap(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax)  
    real(DP):: xyz_Gamma(DimXMin:DimXMax,DimYMin:DimYMax,DimZMin:DimZMax) 
    integer :: i, s                                       ! 添え字  

    !---------------------------------------------------------------------
    ! 初期化
    !---------------------------------------------------------------------
    xyz_MixRtV_pre  = 0.0d0
    xyz_MixRtV_nxt  = 0.0d0
    xyz_MixRtC_pre  = 0.0d0
    xyz_MixRtC_nxt  = 0.0d0
    xyz_PotTemp_pre = 0.0d0
    xyz_PotTemp_nxt = 0.0d0
    
    xyz_ExnerAll    = 0.0d0
    xyz_TempAll     = 0.0d0
    xyz_DelMixRt    = 0.0d0
    xyz_MixRtSat    = 0.0d0 
    xyz_Gamma       = 0.0d0
    
    !---------------------------------------------------------------------
    ! 湿潤飽和調節法の実行
    !   ループを回すのは, 雲についてだけ.  
    !---------------------------------------------------------------------
    LoopSvapPress: do s = 1, CondNum
      
      !湿潤飽和法では圧力は変化させない. 
      xyz_ExnerAll = xyz_Exner + xyz_ExnerBasicZ
          
      !今までに得られた擾乱成分の値を暫定量とみなす. 添え字を追加
      xyz_MixRtV_pre  = xyza_MixRt(:,:,:,IdxCG(s)) + xyza_MixRtBasicZ(:,:,:,IdxCG(s))
      xyz_MixRtC_pre  = xyza_MixRt(:,:,:,IdxCC(s)) + xyza_MixRtBasicZ(:,:,:,IdxCC(s))
      xyz_PotTemp_pre = xyz_PotTemp
          
      Adjusting: do i = 1, ItrNum
        !---------------------------------------------------------------
        ! 飽和蒸気圧から飽和混合比を求める
        !---------------------------------------------------------------
        !温度
        xyz_TempAll = ( xyz_PotTemp_pre + xyz_PotTempBasicZ ) * xyz_ExnerAll

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

        !規格化された潜熱
        xyz_Gamma = xyz_LatentHeat(SpcWetID(IdxCC(s)), xyz_TempAll) * xyz_EffMolWtBasicZ / (xyz_ExnerAll * CpDry)
        
        !凝結量を求める. 
        !  凝結が生じる場合には, xz_MixRtV_pre - xz_MixRtSat は必ず正の値となる.
        !  蒸発が生じる場合には, 蒸発量は - MixRtC を超えることはない. 
        xyz_DelMixRt = ( xyz_MixRtV_pre - xyz_MixRtSat ) / (1.0d0 + xyz_Gamma * xyz_DMixRtSatDPotTemp( SpcWetID(IdxCC(s)), MolWtWet(IdxCC(s)), xyz_TempAll, xyz_ExnerAll ) )

        xyz_Cond = max( 0.0d0, min( xyz_MixRtV_pre,   xyz_DelMixRt ) )
        xyz_Evap = max( 0.0d0, min( xyz_MixRtC_pre, - xyz_DelMixRt ) )

        !より真に近い値を計算
        xyz_PotTemp_nxt = xyz_PotTemp_pre + xyz_Gamma * ( xyz_Cond - xyz_Evap )
        xyz_MixRtV_nxt  = xyz_MixRtV_pre  - xyz_Cond + xyz_Evap
        xyz_MixRtC_nxt  = xyz_MixRtC_pre  + xyz_Cond - xyz_Evap
        
        !繰り返しのための変数定義
        xyz_PotTemp_pre = xyz_PotTemp_nxt
        xyz_MixRtV_pre  = xyz_MixRtV_nxt
        xyz_MixRtC_pre  = xyz_MixRtC_nxt

      end do Adjusting
      
      xyz_PotTemp              = xyz_PotTemp_nxt                 
      xyza_MixRt(:,:,:,IdxCG(s)) = xyz_MixRtV_nxt - xyza_MixRtBasicZ(:,:,:,IdxCG(s))
      xyza_MixRt(:,:,:,IdxCC(s)) = xyz_MixRtC_nxt - xyza_MixRtBasicZ(:,:,:,IdxCC(s))
      
    end do LoopSvapPress
    
  end subroutine MoistAdjustSvapPress

[Validate]