Class BasicSet
In: src/setup/basicset.f90

Methods

Included Modules

dc_trace gt4_history gridset dc_message physset

Public Instance methods

basicfile :character(*), intent(in)
: Input

[Source]

  subroutine BasicSet_file_init(basicfile)

    !=== Dependency


    !=== Input
    character(*), intent(in) :: basicfile
                                                                 !=end
    
    character(20)  :: name            !変数名

    call BeginSub("basicset_file_init",                  fmt="%c",                        c1="Get basic state variables from netCDF file.")
    
    !--- 配列の割り当て
    allocate(   ss_CpBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),  ss_ExnerBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),  ss_PotTempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),  ss_DensBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),  ss_VelSoundBasicZ(DimXMin:DimXMax, DimZMin:DimZMax) )

                                                                 !=begin
    !=== Get a Value from netCDF File
    name = "CpBasicZ"
    call HistoryGet( basicfile, name, ss_CpBasicZ )

    name = "ExnerBasicZ"
    call HistoryGet( basicfile, name, ss_ExnerBasicZ )

    name = "PotTempBasicZ"
    call HistoryGet( basicfile, name, ss_PotTempBasicZ )

    name = "DensBasicZ"
    call HistoryGet( basicfile, name, ss_DensBasicZ )

    name = "VelSoundBasicZ"
    call HistoryGet( basicfile, name, ss_VelSoundBasicZ ) 
                                                                 !=end   

    call EndSub("basicset_file_init")

  end subroutine BasicSet_file_init
cfgfile :character(*), intent(in)
: Input 設定ファイル(NAMELIST)

典型的な基本場のを計算するためのサブルーチン. タイプとしては以下が設定可能.

  * 温度一様 or 温位一様
  * 圧力一様 or 静水圧平衡

[Source]

  subroutine basicset_calc_init(cfgfile)
    !=== Dependency

    
    !=== Input
    character(*), intent(in) :: cfgfile  !設定ファイル(NAMELIST)
                                                                 !=end
    !=== Work
    character(80) :: TempType
    character(80) :: PressType
    real(8)  :: ss_TempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax) ! 平均温度場
    real(8)  :: ss_CvBasicZ(DimXMin:DimXMax, DimZMin:DimZMax)
    real(8)  :: ss_Z(DimXMin:DimXMax, DimZMin:DimZMax)          ! 2-D 座標
    real(8)  :: GasR
    integer  :: k
                                                                 !=begin
    !=== NAMELIST 
    NAMELIST /basicset/ TempType, PressType
                                                                 !=end

    call BeginSub("basicset_calc_init",                  fmt="%c",                        c1="calculate the basic state variables.")

    open (10, FILE=cfgfile)
    read(10, NML=basicset)
    close(10)

    !==== 確認
!    write(*,*) "TempType",  TempType
!    write(*,*) "PressType", PressType
    

    !=== 変数の初期化
    allocate( ss_CpBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),         ss_DensBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),       ss_ExnerBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),      ss_PotTempBasicZ(DimXMin:DimXMax, DimZMin:DimZMax),     ss_VelSoundBasicZ(DimXMin:DimXMax, DimZMin:DimZMax) )

    ss_CpBasicZ = 0.0d0
    ss_DensBasicZ = 0.0d0
    ss_ExnerBasicZ = 0.0d0
    ss_PotTempBasicZ = 0.0d0
    ss_VelSoundBasicZ = 0.0d0

    !--- 座標を初期化
    do k = DimZMin, DimZMax
       ss_Z(:,k) = s_Z(k)
    end do
    
    !--------------------------------------------
    !--- 気体定数を決める
    !--------------------------------------------
    GasR = GasRUniv / MolWtDry
    
    !--------------------------------------------
    !--- 比熱を決める
    !---   比熱は乾燥大気を念頭に一様とする. 
    !--------------------------------------------
    ss_CpBasicZ = CpDry / MolWtDry
    ss_CvBasicZ = ss_CpBasicZ - GasR
  
    !--------------------------------------------
    !--- 温度/温位, 無次元圧力の基本場
    !--------------------------------------------
    if (TempType == "uniform" .AND. PressType == "uniform") then 
       
       !--- 温度場は TempSfc で一様
       ss_TempBasicZ = TempSfc
       
       !--- 地表面圧力で一定. 定義より値は 1 となる. 
       ss_ExnerBasicZ = 1.0d0
       
       !--- 温位
       ss_PotTempBasicZ = ss_TempBasicZ / ss_ExnerBasicZ 
       
    elseif (TempType == "uniform" .AND. PressType == "hydrostatic") then 
       
       !--- 温度場は TempSfc で一様
       ss_TempBasicZ = TempSfc
       
       !--- 無次元圧力と温位を連立して解く. 
       ss_ExnerBasicZ =  dexp( - Grav * ss_Z / ( ss_CpBasicZ * ss_TempBasicZ ) )

       !--- 無次元圧力を利用して温位を計算
       ss_PotTempBasicZ = ss_TempBasicZ / ss_ExnerBasicZ
       
    elseif (TempType == "adiabat" .AND. PressType == "uniform") then 
       
       !--- 地表面での温位(= TempSfc)で一定
       ss_PotTempBasicZ = TempSfc
       
       !--- 地表面圧力で一定. 定義より値は 1 となる. 
       ss_ExnerBasicZ = 1.0d0
       
    elseif (TempType == "adiabat" .AND. PressType == "hydrostatic") then 
       
       !--- 地表面での温位(= TempSfc)で一定
       ss_PotTempBasicZ = TempSfc
       
       !--- 無次元圧力を計算. 温位一定
       ss_ExnerBasicZ =  - Grav * ss_Z / ( ss_CpBasicZ * ss_PotTempBasicZ ) + 1.0d0
       
    else
       call MessageNotify("Message", "basicset_calc_init",                          "unknown TempType, or PressType")
    end if
    
    
    !--------------------------------------------
    !---- 密度
    !--------------------------------------------
    ss_DensBasicZ =  PressSfc * ( ss_ExnerBasicZ ** (ss_CvBasicZ / GasR) )    / ( GasR * ss_PotTempBasicZ )
    
    
    !--------------------------------------------    
    !--- 平均場の音速
    !--------------------------------------------
    ss_VelSoundBasicZ =  ss_CpBasicZ * GasR * ss_ExnerBasicZ * ss_PotTempBasicZ    / ss_CvBasicZ  
    ss_VelSoundBasicZ = sqrt( ss_VelSoundBasicZ )


    !==== 確認
!    write(*,*) "ss_CpBasicZ ", maxval( ss_CpBasicZ )
!    write(*,*) "ss_DensBasicZ ", maxval( ss_DensBasicZ )
!    write(*,*) "ss_ExnerBasicZ ", maxval( ss_ExnerBasicZ )
!    write(*,*) "ss_PotTempBasicZ ", maxval( ss_PotTempBasicZ )
!    write(*,*) "ss_VelSoundBasicZ ", maxval( ss_VelSoundBasicZ )
    

    call EndSub("basicset_calc_init")

  end subroutine basicset_calc_init

[Validate]