Class surface_flux_bulk
In: surface_flux/surface_flux_bulk.f90

地表面フラックス

Surface flux

Note that Japanese and English are described in parallel.

地表面フラックスを計算します.

Surface fluxes are calculated.

Procedures List

SurfaceFlux :地表面フラックスの計算
———— :————
SurfaceFlux :Calculate surface fluxes

NAMELIST

NAMELIST#surface_flux_bulk_nml

Methods

Included Modules

gridset dc_types dc_message constants saturate_nha1992 timeset gtool_historyauto dc_trace namelist_util dc_iounit dc_string axesset

Public Instance methods

Subroutine :
xyz_U(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ u $ . 東西風速. Eastward wind
xyz_V(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ v $ . 南北風速. Northward wind
xyz_Temp(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ T $ . 温度 (整数レベル). Temperature (full level)
xyr_Temp(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ T $ . 温度 (半整数レベル). Temperature (half level)
xyz_QVap(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ q $ . 比湿. Specific humidity
xyz_Press(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ p_s $ . 地表面気圧 (整数レベル). Surface pressure (full level)
xyr_Press(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(in)
: $ p_s $ . 地表面気圧 (半整数レベル). Surface pressure (half level)
xyz_GeoPot(0:imax-1, 1:jmax, 1:kmax) :real(DP), intent(in)
: $ phi $ . ジオポテンシャル (整数レベル). Geo-potential (full level)
xy_SurfTemp(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表面温度. Surface temperature
xy_SurfHumidCoeff(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表湿潤度. Surface humidity coefficient
xy_SurfRoughLength(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表粗度長. Surface rough length
xy_SurfCondition(0:imax-1, 1:jmax) :integer, intent(in)
: 地表状態. Surface condition
xyr_UFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(inout)
: 東西風速フラックス. Eastward wind flux
xyr_VFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(inout)
: 南北風速フラックス. Northward wind flux
xyr_TempFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(inout)
: 温度フラックス. Temperature flux
xyr_QVapFlux(0:imax-1, 1:jmax, 0:kmax) :real(DP), intent(inout)
: 比湿フラックス. Specific humidity flux
xy_SurfUVMtx(0:imax-1, 1:jmax) :real(DP), intent(out)
: 速度陰解行列: 地表. Implicit matrix about velocity: surface
xyaa_SurfTempMtx(0:imax-1, 1:jmax, 0:1, -1:1) :real(DP), intent(out)
xyaa_SurfQVapMtx(0:imax-1, 1:jmax, 0:1, -1:1) :real(DP), intent(out)
: 比湿陰解行列: 地表. Implicit matrix about specific humidity: surface

温度, 比湿, 気圧から, 放射フラックスを計算します.

Calculate radiation flux from temperature, specific humidity, and air pressure.

[Source]

  subroutine SurfaceFlux( xyz_U, xyz_V, xyz_Temp, xyr_Temp, xyz_QVap, xyz_Press, xyr_Press, xyz_GeoPot, xy_SurfTemp, xy_SurfHumidCoeff, xy_SurfRoughLength, xy_SurfCondition, xyr_UFlux, xyr_VFlux, xyr_TempFlux, xyr_QVapFlux, xy_SurfUVMtx, xyaa_SurfTempMtx, xyaa_SurfQVapMtx )
    !
    ! 温度, 比湿, 気圧から, 放射フラックスを計算します. 
    !
    ! Calculate radiation flux from temperature, specific humidity, and 
    ! air pressure. 
    !

    ! モジュール引用 ; USE statements
    !

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: Grav, GasRDry, CpDry, LatentHeat
                              ! $ L $ [J kg-1] . 
                              ! 凝結の潜熱. 
                              ! Latent heat of condensation

    ! 飽和比湿計算 (Nakajima et al., 1992)
    ! Evaluate saturation specific humidity (Nakajima et al., 1992)
    !
    use saturate_nha1992, only: CalcQVapSat, CalcDQVapSatDTemp

    ! 時刻管理
    ! Time control
    !
    use timeset, only: TimeN, TimesetClockStart, TimesetClockStop

    ! ヒストリデータ出力
    ! History data output
    !
    use gtool_historyauto, only: HistoryAutoPut

    ! デバッグ用ユーティリティ
    ! Utilities for debug
    !
    use dc_trace, only: DbgMessage, BeginSub, EndSub

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(in):: xyz_U (0:imax-1, 1:jmax, 1:kmax)
                              ! $ u $ . 東西風速. Eastward wind
    real(DP), intent(in):: xyz_V (0:imax-1, 1:jmax, 1:kmax)
                              ! $ v $ . 南北風速. Northward wind

    real(DP), intent(in):: xyz_Temp (0:imax-1, 1:jmax, 1:kmax)
                              ! $ T $ . 温度 (整数レベル). 
                              ! Temperature (full level)
    real(DP), intent(in):: xyr_Temp (0:imax-1, 1:jmax, 0:kmax)
                              ! $ T $ . 温度 (半整数レベル). 
                              ! Temperature (half level)
    real(DP), intent(in):: xyz_QVap (0:imax-1, 1:jmax, 1:kmax)
                              ! $ q $ .     比湿. Specific humidity
    real(DP), intent(in):: xyz_Press (0:imax-1, 1:jmax, 1:kmax)
                              ! $ p_s $ . 地表面気圧 (整数レベル). 
                              ! Surface pressure (full level)
    real(DP), intent(in):: xyr_Press (0:imax-1, 1:jmax, 0:kmax)
                              ! $ p_s $ . 地表面気圧 (半整数レベル). 
                              ! Surface pressure (half level)
    real(DP), intent(in):: xyz_GeoPot (0:imax-1, 1:jmax, 1:kmax)
                              ! $ \phi $ . ジオポテンシャル (整数レベル). 
                              ! Geo-potential (full level)
    real(DP), intent(in):: xy_SurfTemp (0:imax-1, 1:jmax)
                              ! 地表面温度. 
                              ! Surface temperature
    real(DP), intent(in):: xy_SurfHumidCoeff (0:imax-1, 1:jmax)
                              ! 地表湿潤度. 
                              ! Surface humidity coefficient
    real(DP), intent(in):: xy_SurfRoughLength (0:imax-1, 1:jmax)
                              ! 地表粗度長. 
                              ! Surface rough length
    integer, intent(in):: xy_SurfCondition (0:imax-1, 1:jmax)
                              ! 地表状態. 
                              ! Surface condition
    real(DP), intent(inout):: xyr_UFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 東西風速フラックス. 
                              ! Eastward wind flux
    real(DP), intent(inout):: xyr_VFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 南北風速フラックス. 
                              ! Northward wind flux
    real(DP), intent(inout):: xyr_TempFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 温度フラックス. 
                              ! Temperature flux
    real(DP), intent(inout):: xyr_QVapFlux (0:imax-1, 1:jmax, 0:kmax)
                              ! 比湿フラックス. 
                              ! Specific humidity flux
    real(DP), intent(out):: xy_SurfUVMtx (0:imax-1, 1:jmax)
                              !  速度陰解行列: 地表. 
                              ! Implicit matrix about velocity: surface
    real(DP), intent(out):: xyaa_SurfTempMtx (0:imax-1, 1:jmax, 0:1, -1:1)

                              ! 温度陰解行列: 地表. 
                              ! Implicit matrix about temperature: surface
    real(DP), intent(out):: xyaa_SurfQVapMtx (0:imax-1, 1:jmax, 0:1, -1:1)
                              ! 比湿陰解行列: 地表. 
                              ! Implicit matrix about specific humidity: surface


    ! 作業変数
    ! Work variables
    !
    real(DP):: xy_SurfBulkRiNum (0:imax-1, 1:jmax)
                              ! バルク $ R_i $ 数. 
                              ! Bulk $ R_i $ number
    real(DP):: xy_SurfTempTransCoeff (0:imax-1, 1:jmax)
                              ! 輸送係数:温度. 
                              ! Transfer coefficient: temperature
    real(DP):: xy_SurfQVapTransCoeff (0:imax-1, 1:jmax)
                              ! 輸送係数:比湿. 
                              ! Transfer coefficient: specific humidity
    real(DP):: xy_SurfVelTransCoeff (0:imax-1, 1:jmax)
                              ! 輸送係数:運動量. 
                              ! Diffusion coefficient: velocity
    real(DP):: xy_SurfTempBulkCoeff (0:imax-1, 1:jmax)
                              ! バルク係数:温度. 
                              ! Bulk coefficient: temperature
    real(DP):: xy_SurfQVapBulkCoeff (0:imax-1, 1:jmax)
                              ! バルク係数:比湿. 
                              ! Bulk coefficient: specific humidity
    real(DP):: xy_SurfVelBulkCoeff (0:imax-1, 1:jmax)
                              ! バルク係数:運動量. 
                              ! Bulk coefficient: temperature
    real(DP):: xy_SurfExner (0:imax-1, 1:jmax)
                              ! Exner 関数. 
                              ! Exner function
    real(DP):: xy_SurfVelAbs (0:imax-1, 1:jmax)
                              ! 風速絶対値. 
                              ! Absolute velocity
    real(DP):: xy_SurfQVapSat (0:imax-1, 1:jmax)
                              ! 地表飽和比湿. 
                              ! Saturated specific humidity on surface
    real(DP):: xy_SurfDQVapSatDTemp (0:imax-1, 1:jmax)
                              ! 地表飽和比湿変化. 
                              ! Saturated specific humidity tendency on surface

    real(DP):: xy_UFluxSurf (0:imax-1, 1:jmax)
                              ! 地表面の東西風速フラックス. 
                              ! Eastward wind flux on surface
    real(DP):: xy_VFluxSurf (0:imax-1, 1:jmax)
                              ! 地表面の南北風速フラックス. 
                              ! Northward wind flux on surface
    real(DP):: xy_TempFluxSurf (0:imax-1, 1:jmax)
                              ! 地表面の温度フラックス. 
                              ! Temperature flux on surface
    real(DP):: xy_QVapFluxSurf (0:imax-1, 1:jmax)
                              ! 地表面の比湿フラックス. 
                              ! Specific humidity flux on surface

    integer:: i               ! 経度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in longitude
    integer:: j               ! 緯度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in latitude

    ! 実行文 ; Executable statement
    !

    ! 計算時間計測開始
    ! Start measurement of computation time
    !
    call TimesetClockStart( module_name )

    ! 初期化
    ! Initialization
    !
    if ( .not. surface_flux_bulk_inited ) call SurfFluxInit

    ! Exner 関数算出
    ! Calculate Exner functions
    !
    xy_SurfExner = ( xyz_Press(:,:,1) / xyr_Press(:,:,0)  )**( GasRDry / CpDry )

    ! バルク $ R_i $ 数算出
    ! Calculate bulk $ R_i $
    !
    do i = 0, imax-1
      do j = 1, jmax
        xy_SurfVelAbs(i,j) = sqrt ( xyz_U(i,j,1)**2 + xyz_V(i,j,1)**2 )

        xy_SurfBulkRiNum(i,j) = Grav / BasePotTemp * ( xyz_Temp(i,j,1) / xy_SurfExner(i,j) - xy_SurfTemp(i,j)  ) / max( xy_SurfVelAbs(i,j), VelMinForRi )**2 * xyz_GeoPot(i,j,1)
      end do
    end do
    
    ! バルク係数算出
    ! Calculate bulk coefficients
    !
    call BulkCoeff( xy_SurfBulkRiNum, xy_SurfRoughLength, xy_SurfRoughLength, xyz_GeoPot(:,:,1), xy_SurfVelBulkCoeff, xy_SurfTempBulkCoeff, xy_SurfQVapBulkCoeff )  ! (out)

    ! 輸送係数の計算
    ! Calculate transfer coefficient
    !
    do i = 0, imax-1
      do j = 1, jmax
        
        xy_SurfVelTransCoeff(i,j) = xy_SurfVelBulkCoeff(i,j) * xyr_Press(i,j,0) / ( GasRDry * xyr_Temp(i,j,0) ) * min( max( xy_SurfVelAbs(i,j), VelMinForVel ), VelMaxForVel )
        
        xy_SurfTempTransCoeff(i,j) = xy_SurfTempBulkCoeff(i,j) * xyr_Press(i,j,0) / ( GasRDry * xyr_Temp(i,j,0) ) * min( max( xy_SurfVelAbs(i,j), VelMinForTemp ), VelMaxForTemp )
        
        xy_SurfQVapTransCoeff(i,j) = xy_SurfQVapBulkCoeff(i,j) * xyr_Press(i,j,0) / ( GasRDry * xyr_Temp(i,j,0) ) * min( max( xy_SurfVelAbs(i,j), VelMinForQVap ), VelMaxForQVap )
        
      end do
    end do
    
    ! 飽和比湿の計算
    ! Calculate saturated specific humidity
    !
    call CalcQVapSat( xy_SurfTemp, xyr_Press(:,:,0), xy_SurfQVapSat )                 ! (out)

    call CalcDQVapSatDTemp( xy_SurfTemp, xyr_Press(:,:,0), xy_SurfDQVapSatDTemp )           ! (out)

!!$    xy_SurfQVapSat = &
!!$      &   EpsV * ES0  &
!!$      &   * exp( LatentHeat / RVap * ( 1.0_DP / 273.0_DP - 1.0_DP / xy_SurfTemp ) ) &
!!$      &   / xyr_Press(:,:,0)
!!$
!!$    xy_SurfDQVapSatDTemp = &
!!$      &   LatentHeat * xy_SurfQVapSat &
!!$      &   / ( RVap * xy_SurfTemp * xy_SurfTemp )

    ! 地表面フラックスの計算
    ! Calculate fluxes on flux
    !
    xy_UFluxSurf    = - xy_SurfVelTransCoeff * xyz_U(:,:,1)
    xy_VFluxSurf    = - xy_SurfVelTransCoeff * xyz_V(:,:,1)
    xy_TempFluxSurf =   CpDry * xy_SurfTempTransCoeff * (   xy_SurfTemp - xyz_Temp(:,:,1) / xy_SurfExner )
    xy_QVapFluxSurf =   LatentHeat * xy_SurfQVapTransCoeff * xy_SurfHumidCoeff * ( xy_SurfQVapSat - xyz_QVap(:,:,1) )

    ! フラックスの計算
    ! Calculate fluxes
    !
    xyr_UFlux(:,:,0)    = xyr_UFlux(:,:,0)    + xy_UFluxSurf
    xyr_VFlux(:,:,0)    = xyr_VFlux(:,:,0)    + xy_VFluxSurf
    xyr_TempFlux(:,:,0) = xyr_TempFlux(:,:,0) + xy_TempFluxSurf
    xyr_QVapFlux(:,:,0) = xyr_QVapFlux(:,:,0) + xy_QVapFluxSurf

!!$    xyr_UFlux(:,:,0) =   xyr_UFlux(:,:,0) &
!!$      &                - xy_SurfVelTransCoeff * xyz_U(:,:,1)
!!$
!!$    xyr_VFlux(:,:,0) =   xyr_VFlux(:,:,0) &
!!$      &                - xy_SurfVelTransCoeff * xyz_V(:,:,1)
!!$
!!$    xyr_TempFlux(:,:,0) =   xyr_TempFlux(:,:,0) &
!!$      &                   + CpDry * xy_SurfTempTransCoeff &
!!$      &                     * (   xy_SurfTemp           &
!!$      &                         - xyz_Temp(:,:,1) / xy_SurfExner )
!!$
!!$    xyr_QVapFlux(:,:,0) =   xyr_QVapFlux(:,:,0) &
!!$      &                   + LatentHeat * xy_SurfQVapTransCoeff * xy_SurfHumidCoeff &
!!$      &                     * ( xy_SurfQVapSat - xyz_QVap(:,:,1) )
    
    ! 陰解行列の計算
    ! Calculate implicit matrices
    !
    xyaa_SurfTempMtx = 0.0_DP
    xyaa_SurfQVapMtx = 0.0_DP
    
    xy_SurfUVMtx = xy_SurfVelTransCoeff
    
    xyaa_SurfTempMtx(:,:,1,0) =   CpDry * xy_SurfTempTransCoeff / xy_SurfExner
    xyaa_SurfTempMtx(:,:,0,1) = - CpDry * xy_SurfTempTransCoeff / xy_SurfExner
    
    xyaa_SurfQVapMtx(:,:,1,0) =   CpDry * xy_SurfQVapTransCoeff * xy_SurfHumidCoeff
    xyaa_SurfQVapMtx(:,:,0,1) = - CpDry * xy_SurfQVapTransCoeff * xy_SurfHumidCoeff
    
    do i = 0, imax-1
      do j = 1, jmax
        if ( xy_SurfCondition(i,j) >= 1 ) then
          
          xyaa_SurfTempMtx(i,j,1,-1) = - CpDry * xy_SurfTempTransCoeff(i,j)
          xyaa_SurfTempMtx(i,j,0,0)  =   CpDry * xy_SurfTempTransCoeff(i,j)
          
          xyaa_SurfQVapMtx(i,j,1,-1) = - LatentHeat * xy_SurfQVapTransCoeff(i,j) * xy_SurfHumidCoeff(i,j) * xy_SurfDQVapSatDTemp(i,j)

          xyaa_SurfQVapMtx(i,j,0,0)  = LatentHeat * xy_SurfQVapTransCoeff(i,j) * xy_SurfHumidCoeff(i,j) * xy_SurfDQVapSatDTemp(i,j)
          
        end if
      end do
    end do

    ! ヒストリデータ出力
    ! History data output
    !
    call HistoryAutoPut( TimeN, 'UFluxSurf', xy_UFluxSurf )
    call HistoryAutoPut( TimeN, 'VFluxSurf', xy_VFluxSurf )
    call HistoryAutoPut( TimeN, 'TempFluxSurf', xy_TempFluxSurf )
    call HistoryAutoPut( TimeN, 'QVapFluxSurf', xy_QVapFluxSurf )

    ! 計算時間計測一時停止
    ! Pause measurement of computation time
    !
    call TimesetClockStop( module_name )

  end subroutine SurfaceFlux
surface_flux_bulk_inited
Variable :
surface_flux_bulk_inited = .false. :logical, save, public
: 初期設定フラグ. Initialization flag

Private Instance methods

BasePotTemp
Variable :
BasePotTemp :real(DP), save
: 基本温位. Basic potential temperature
Subroutine :
xy_SurfBulkRiNum(0:imax-1, 1:jmax) :real(DP), intent(in)
: バルク $ R_i $ 数. Bulk $ R_i $ number
xy_SurfVelRoughLength(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表粗度長 (運動量). Surface rough length (momentum)
xy_SurfTempRoughLength(0:imax-1, 1:jmax) :real(DP), intent(in)
: 地表粗度長 (熱). Surface rough length (thermal)
xy_SurfGeoPot(0:imax-1, 1:jmax) :real(DP), intent(in)
: $ phi $ . ジオポテンシャル (地表面). Geo-potential (on surface)
xy_SurfVelBulkCoeff(0:imax-1, 1:jmax) :real(DP), intent(out)
: バルク係数:運動量. Bulk coefficient: temperature
xy_SurfTempBulkCoeff(0:imax-1, 1:jmax) :real(DP), intent(out)
: バルク係数:温度. Bulk coefficient: temperature
xy_SurfQVapBulkCoeff(0:imax-1, 1:jmax) :real(DP), intent(out)
: バルク係数:比湿. Bulk coefficient: specific humidity

バルク係数を算出します.

Bulk coefficients are calculated.

[Source]

  subroutine BulkCoeff( xy_SurfBulkRiNum, xy_SurfVelRoughLength, xy_SurfTempRoughLength, xy_SurfGeoPot, xy_SurfVelBulkCoeff, xy_SurfTempBulkCoeff, xy_SurfQVapBulkCoeff )
    !
    ! バルク係数を算出します.
    !
    ! Bulk coefficients are calculated.
    !

    ! モジュール引用 ; USE statements
    !

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: FKarm                 ! $ k $ .
                              ! カルマン定数. 
                              ! Karman constant

    ! 宣言文 ; Declaration statements
    !
    implicit none
    real(DP), intent(in):: xy_SurfBulkRiNum (0:imax-1, 1:jmax)
                              ! バルク $ R_i $ 数. 
                              ! Bulk $ R_i $ number

!!$    real(DP), intent(in):: xy_SurfVelAbs (0:imax-1, 1:jmax)
!!$                              ! 風速絶対値. 
!!$                              ! Absolute velocity
    real(DP), intent(in):: xy_SurfVelRoughLength (0:imax-1, 1:jmax)
                              ! 地表粗度長 (運動量). 
                              ! Surface rough length (momentum)
    real(DP), intent(in):: xy_SurfTempRoughLength (0:imax-1, 1:jmax)
                              ! 地表粗度長 (熱). 
                              ! Surface rough length (thermal)
    real(DP), intent(in):: xy_SurfGeoPot (0:imax-1, 1:jmax)
                              ! $ \phi $ . ジオポテンシャル (地表面). 
                              ! Geo-potential (on surface)
    real(DP), intent(out):: xy_SurfVelBulkCoeff (0:imax-1, 1:jmax)
                              ! バルク係数:運動量. 
                              ! Bulk coefficient: temperature
    real(DP), intent(out):: xy_SurfTempBulkCoeff (0:imax-1, 1:jmax)
                              ! バルク係数:温度. 
                              ! Bulk coefficient: temperature
    real(DP), intent(out):: xy_SurfQVapBulkCoeff (0:imax-1, 1:jmax)
                              ! バルク係数:比湿. 
                              ! Bulk coefficient: specific humidity

    ! 作業変数
    ! Work variables
    !
    integer:: i               ! 経度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in longitude
    integer:: j               ! 緯度方向に回る DO ループ用作業変数
                              ! Work variables for DO loop in latitude

    ! 実行文 ; Executable statement
    !

    ! 中立バルク係数の計算
    ! Calculate neutral bulk coefficients
    !
    if ( ConstBulkCoeff < 0.0_DP ) then
      
      xy_SurfVelBulkCoeff  = ( FKarm / log ( xy_SurfGeoPot / xy_SurfVelRoughLength ) )**2

      xy_SurfTempBulkCoeff = ( FKarm / log ( xy_SurfGeoPot / xy_SurfTempRoughLength ) )**2

      xy_SurfQVapBulkCoeff = xy_SurfTempBulkCoeff
      
    else
      xy_SurfVelBulkCoeff  = ConstBulkCoeff
      xy_SurfTempBulkCoeff = ConstBulkCoeff
      xy_SurfQVapBulkCoeff = ConstBulkCoeff
    end if
    
    ! 非中立バルク係数の計算
    ! Calculate non-neutral bulk coefficients
    !
    if ( .not. Neutral ) then
      
      do i = 0, imax-1
        do j = 1, jmax

          if ( xy_SurfBulkRiNum(i,j) > 0.0_DP ) then 

            xy_SurfVelBulkCoeff(i,j) = xy_SurfVelBulkCoeff(i,j) / (   1.0_DP + 10.0_DP * xy_SurfBulkRiNum(i,j) / sqrt( 1.0_DP + 5.0_DP * xy_SurfBulkRiNum(i,j) ) )

            xy_SurfTempBulkCoeff(i,j) = xy_SurfTempBulkCoeff(i,j) / (   1.0_DP + 15.0_DP * xy_SurfBulkRiNum(i,j) / sqrt( 1.0_DP + 5.0_DP * xy_SurfBulkRiNum(i,j) ) )

            xy_SurfQVapBulkCoeff(i,j) = xy_SurfTempBulkCoeff(i,j)

          else

            xy_SurfVelBulkCoeff(i,j) = xy_SurfVelBulkCoeff(i,j) * (   1.0_DP - 10.0_DP * xy_SurfBulkRiNum(i,j) / (   1.0_DP + 75.0_DP * xy_SurfVelBulkCoeff(i,j) * sqrt( - xy_SurfGeoPot(i,j) / xy_SurfVelRoughLength(i,j) * xy_SurfBulkRiNum(i,j) ) ) )
            
            xy_SurfTempBulkCoeff(i,j) = xy_SurfTempBulkCoeff(i,j) * (   1.0_DP - 15.0_DP * xy_SurfBulkRiNum(i,j) / (   1.0_DP + 75.0_DP * xy_SurfTempBulkCoeff(i,j) * sqrt( - xy_SurfGeoPot(i,j) / xy_SurfTempRoughLength(i,j) * xy_SurfBulkRiNum(i,j) ) ) )

            xy_SurfQVapBulkCoeff(i,j) = xy_SurfTempBulkCoeff(i,j)

          end if
        end do
      end do
      
    end if
    
    ! 最大/最小 判定
    ! Measure maximum/minimum
    !
    do i = 0, imax-1
      do j = 1, jmax

        xy_SurfVelBulkCoeff(i,j)  = max( min( xy_SurfVelBulkCoeff(i,j), VelBulkCoeffMax ), VelBulkCoeffMin )

        xy_SurfTempBulkCoeff(i,j) = max( min( xy_SurfTempBulkCoeff(i,j), TempBulkCoeffMax ), TempBulkCoeffMin )

        xy_SurfQVapBulkCoeff(i,j) = max( min( xy_SurfQVapBulkCoeff(i,j), QVapBulkCoeffMax ), QVapBulkCoeffMin )

      end do
    end do

  end subroutine BulkCoeff
ConstBulkCoeff
Variable :
ConstBulkCoeff :real(DP), save
: バルク係数一定値. Steady value of bulk coefficient
Subroutine :

依存モジュールの初期化チェック

Check initialization of dependency modules

[Source]

  subroutine InitCheck
    !
    ! 依存モジュールの初期化チェック
    !
    ! Check initialization of dependency modules

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_util_inited

    ! 格子点設定
    ! Grid points settings
    !
    use gridset, only: gridset_inited

    ! 物理定数設定
    ! Physical constants settings
    !
    use constants, only: constants_inited

    ! 座標データ設定
    ! Axes data settings
    !
    use axesset, only: axesset_inited

    ! 時刻管理
    ! Time control
    !
    use timeset, only: timeset_inited

    ! 実行文 ; Executable statement
    !

    if ( .not. namelist_util_inited ) call MessageNotify( 'E', module_name, '"namelist_util" module is not initialized.' )

    if ( .not. gridset_inited ) call MessageNotify( 'E', module_name, '"gridset" module is not initialized.' )

    if ( .not. constants_inited ) call MessageNotify( 'E', module_name, '"constants" module is not initialized.' )

    if ( .not. axesset_inited ) call MessageNotify( 'E', module_name, '"axesset" module is not initialized.' )

    if ( .not. timeset_inited ) call MessageNotify( 'E', module_name, '"timeset" module is not initialized.' )

  end subroutine InitCheck
Neutral
Variable :
Neutral :logical, save
: 中立であるか否か. Whether it was neutral or not?
QVapBulkCoeffMax
Variable :
QVapBulkCoeffMax :real(DP), save
: $ q $ バルク係数最大値. Maximum value of $ q $ bulk coefficient
QVapBulkCoeffMin
Variable :
QVapBulkCoeffMin :real(DP), save
: $ q $ バルク係数最小値. Minimum value of $ q $ bulk coefficient
Subroutine :

surface_flux_bulk モジュールの初期化を行います. NAMELIST#surface_flux_bulk_nml の読み込みはこの手続きで行われます.

"surface_flux_bulk" module is initialized. "NAMELIST#surface_flux_bulk_nml" is loaded in this procedure.

This procedure input/output NAMELIST#surface_flux_bulk_nml .

[Source]

  subroutine SurfFluxInit
    !
    ! surface_flux_bulk モジュールの初期化を行います. 
    ! NAMELIST#surface_flux_bulk_nml の読み込みはこの手続きで行われます. 
    !
    ! "surface_flux_bulk" module is initialized. 
    ! "NAMELIST#surface_flux_bulk_nml" is loaded in this procedure. 
    !

    ! モジュール引用 ; USE statements
    !

    ! NAMELIST ファイル入力に関するユーティリティ
    ! Utilities for NAMELIST file input
    !
    use namelist_util, only: namelist_filename, NmlutilMsg, NmlutilAryValid

    ! ファイル入出力補助
    ! File I/O support
    !
    use dc_iounit, only: FileOpen

    ! 種別型パラメタ
    ! Kind type parameter
    !
    use dc_types, only: STDOUT ! 標準出力の装置番号. Unit number of standard output

    ! 文字列操作
    ! Character handling
    !
    use dc_string, only: StoA

    ! ヒストリデータ出力
    ! History data output
    !
    use gtool_historyauto, only: HistoryAutoAddVariable

    ! 宣言文 ; Declaration statements
    !
    implicit none

    integer:: unit_nml        ! NAMELIST ファイルオープン用装置番号. 
                              ! Unit number for NAMELIST file open
    integer:: iostat_nml      ! NAMELIST 読み込み時の IOSTAT. 
                              ! IOSTAT of NAMELIST read

    ! NAMELIST 変数群
    ! NAMELIST group name
    !
    namelist /surface_flux_bulk_nml/ VelMinForRi, VelMinForVel, VelMinForTemp, VelMinForQVap, VelMaxForVel, VelMaxForTemp, VelMaxForQVap, Neutral, ConstBulkCoeff, VelBulkCoeffMin, TempBulkCoeffMin, QVapBulkCoeffMin, VelBulkCoeffMax, TempBulkCoeffMax, QVapBulkCoeffMax
          !
          ! デフォルト値については初期化手続 "surface_flux_bulk#SurfFluxInit" 
          ! のソースコードを参照のこと. 
          !
          ! Refer to source codes in the initialization procedure
          ! "surface_flux_bulk#SurfFluxInit" for the default values. 
          !

    ! 実行文 ; Executable statement
    !

    if ( surface_flux_bulk_inited ) return
    call InitCheck

    ! デフォルト値の設定
    ! Default values settings
    !
    BasePotTemp   = 300.0_DP
    VelMinForRi   = 0.01_DP
    VelMinForVel  = 0.01_DP
    VelMinForTemp = 0.01_DP
    VelMinForQVap = 0.01_DP
    VelMaxForVel  = 1000.0_DP
    VelMaxForTemp = 1000.0_DP
    VelMaxForQVap = 1000.0_DP

    Neutral          = .false.
    ConstBulkCoeff   = -1.0_DP
    VelBulkCoeffMin  =  0.0_DP
    TempBulkCoeffMin =  0.0_DP
    QVapBulkCoeffMin =  0.0_DP
    VelBulkCoeffMax  =  1.0_DP
    TempBulkCoeffMax =  1.0_DP
    QVapBulkCoeffMax =  1.0_DP


    ! NAMELIST の読み込み
    ! NAMELIST is input
    !
    if ( trim(namelist_filename) /= '' ) then
      call FileOpen( unit_nml, namelist_filename, mode = 'r' ) ! (in)

      rewind( unit_nml )
      read( unit_nml, nml = surface_flux_bulk_nml, iostat = iostat_nml )        ! (out)
      close( unit_nml )

      call NmlutilMsg( iostat_nml, module_name ) ! (in)
    end if

    ! ヒストリデータ出力のためのへの変数登録
    ! Register of variables for history data output
    !
    call HistoryAutoAddVariable( 'UFluxSurf' , (/ 'lon ', 'lat ', 'time'/), 'eastward wind flux by surface process', 'N m-2' )

    call HistoryAutoAddVariable( 'VFluxSurf' , (/ 'lon ', 'lat ', 'time'/), 'northward wind flux by surface process', 'N m-2' )

    call HistoryAutoAddVariable( 'TempFluxSurf' , (/ 'lon ', 'lat ', 'time'/), 'temperature flux by surface process', 'W m-2' )

    call HistoryAutoAddVariable( 'QVapFluxSurf' , (/ 'lon ', 'lat ', 'time'/), 'specific humidity flux by surface process', 'W m-2' )

    ! 印字 ; Print
    !
    call MessageNotify( 'M', module_name, '----- Initialization Messages -----' )

    call MessageNotify( 'M', module_name, '  VelMinForRi   = %f', d = (/ VelMinForRi   /) )
    call MessageNotify( 'M', module_name, '  VelMinForVel  = %f', d = (/ VelMinForVel  /) )
    call MessageNotify( 'M', module_name, '  VelMinForTemp = %f', d = (/ VelMinForTemp /) )
    call MessageNotify( 'M', module_name, '  VelMinForQVap = %f', d = (/ VelMinForQVap /) )
    call MessageNotify( 'M', module_name, '  VelMaxForVel  = %f', d = (/ VelMaxForVel  /) )
    call MessageNotify( 'M', module_name, '  VelMaxForTemp = %f', d = (/ VelMaxForTemp /) )
    call MessageNotify( 'M', module_name, '  VelMaxForQVap = %f', d = (/ VelMaxForQVap /) )

    call MessageNotify( 'M', module_name, 'Bulk coefficients:' )
    call MessageNotify( 'M', module_name, '  Neutral          = %b', l = (/ Neutral          /) )
    call MessageNotify( 'M', module_name, '  ConstBulkCoeff   = %f', d = (/ ConstBulkCoeff   /) )
    call MessageNotify( 'M', module_name, '  VelBulkCoeffMin  = %f', d = (/ VelBulkCoeffMin  /) )
    call MessageNotify( 'M', module_name, '  TempBulkCoeffMin = %f', d = (/ TempBulkCoeffMin /) )
    call MessageNotify( 'M', module_name, '  QVapBulkCoeffMin = %f', d = (/ QVapBulkCoeffMin /) )
    call MessageNotify( 'M', module_name, '  VelBulkCoeffMax  = %f', d = (/ VelBulkCoeffMax  /) )
    call MessageNotify( 'M', module_name, '  TempBulkCoeffMax = %f', d = (/ TempBulkCoeffMax /) )
    call MessageNotify( 'M', module_name, '  QVapBulkCoeffMax = %f', d = (/ QVapBulkCoeffMax /) )
    call MessageNotify( 'M', module_name, '-- version = %c', c1 = trim(version) )

    surface_flux_bulk_inited = .true.
  end subroutine SurfFluxInit
TempBulkCoeffMax
Variable :
TempBulkCoeffMax :real(DP), save
: $ T $ バルク係数最大値. Maximum value of $ T $ bulk coefficient
TempBulkCoeffMin
Variable :
TempBulkCoeffMin :real(DP), save
: $ T $ バルク係数最小値. Minimum value of $ T $ bulk coefficient
VelBulkCoeffMax
Variable :
VelBulkCoeffMax :real(DP), save
: $ u $ バルク係数最大値. Maximum value of $ u $ bulk coefficient
VelBulkCoeffMin
Variable :
VelBulkCoeffMin :real(DP), save
: $ u $ バルク係数最小値. Minimum value of $ u $ bulk coefficient
VelMaxForQVap
Variable :
VelMaxForQVap :real(DP), save
: 水蒸気用風最大値. Maximum value of velocity for vapor
VelMaxForTemp
Variable :
VelMaxForTemp :real(DP), save
: 熱用風最大値. Maximum value of velocity for thermal
VelMaxForVel
Variable :
VelMaxForVel :real(DP), save
: 運動量用風最大値. Maximum value of velocity for momentum
VelMinForQVap
Variable :
VelMinForQVap :real(DP), save
: 水蒸気用風最小値. Minimum value of velocity for vapor
VelMinForRi
Variable :
VelMinForRi :real(DP), save
: $ R_i $ 数用風最小値. Minimum value of velocity for $ R_i $ number
VelMinForTemp
Variable :
VelMinForTemp :real(DP), save
: 熱用風最小値. Minimum value of velocity for thermal
VelMinForVel
Variable :
VelMinForVel :real(DP), save
: 運動量用風最小値. Minimum value of velocity for momentum
module_name
Constant :
module_name = ‘surface_flux_bulk :character(*), parameter
: モジュールの名称. Module name
version
Constant :
version = ’$Name: dcpam5-20081109-1 $’ // ’$Id: surface_flux_bulk.f90,v 1.8 2008-10-06 16:30:12 morikawa Exp $’ :character(*), parameter
: モジュールのバージョン Module version

[Validate]