Class Surfaceflux_bulk
In: ../src/physics/surfaceflux_bulk.f90

下部境界でのフラックスの計算モジュール

Methods

Included Modules

dc_types dc_iounit dc_message gtool_historyauto mpi_wrapper gridset axesset basicset constants composition chemcalc namelist_util timeset DExnerDt

Public Instance methods

Subroutine :
pyz_VelX(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)
: 水平風速
xqz_VelY(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)
: 水平風速
xyz_PTemp(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)
: 温位の擾乱成分
xyz_Exner(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(in)
: 温位の擾乱成分
xyzf_QMix(imin:imax,jmin:jmax,kmin:kmax, ncmax) :real(DP), intent(in)
: 温位の擾乱成分
pyz_DVelXDt(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(inout)
xqz_DVelYDt(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(inout)
xyz_DPTempDt(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(inout)
xyz_DExnerDt(imin:imax,jmin:jmax,kmin:kmax) :real(DP), intent(inout)
xyzf_DQMixDt(imin:imax,jmin:jmax,kmin:kmax, ncmax) :real(DP), intent(inout)

下部境界からのフラックスによる温度の変化率を, バルク方法に基づいて計算する.

[Source]

  subroutine Surfaceflux_Bulk_forcing( pyz_VelX, xqz_VelY, xyz_PTemp, xyz_Exner, xyzf_QMix, pyz_DVelXDt, xqz_DVelYDt, xyz_DPTempDt, xyz_DExnerDt, xyzf_DQMixDt )
    ! 
    ! 下部境界からのフラックスによる温度の変化率を,
    ! バルク方法に基づいて計算する.
    !

    !暗黙の型宣言禁止
    implicit none

    !変数定義
    real(DP), intent(in)   :: pyz_VelX(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速
    real(DP), intent(in)   :: xqz_VelY(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速
    real(DP), intent(in)   :: xyz_PTemp(imin:imax,jmin:jmax,kmin:kmax)
                                           !温位の擾乱成分    
    real(DP), intent(in)   :: xyz_Exner(imin:imax,jmin:jmax,kmin:kmax)
                                           !温位の擾乱成分    
    real(DP), intent(in)   :: xyzf_QMix(imin:imax,jmin:jmax,kmin:kmax, ncmax)
                                           !温位の擾乱成分    
    real(DP), intent(inout):: pyz_DVelXDt(imin:imax,jmin:jmax,kmin:kmax)
    real(DP), intent(inout):: xqz_DVelYDt(imin:imax,jmin:jmax,kmin:kmax)
    real(DP), intent(inout):: xyz_DPTempDt(imin:imax,jmin:jmax,kmin:kmax)
    real(DP), intent(inout):: xyz_DExnerDt(imin:imax,jmin:jmax,kmin:kmax)
    real(DP), intent(inout):: xyzf_DQMixDt(imin:imax,jmin:jmax,kmin:kmax, ncmax)
    real(DP)               :: py_VelXflux (imin:imax,jmin:jmax)
                                           !運動量フラックス
                                           !(strictly speaking, this value is not 
                                           !momentum flux, but is it divided by density)
    real(DP)               :: xq_VelYflux (imin:imax,jmin:jmax)
                                           !運動量フラックス
                                           !(strictly speaking, this value is not 
                                           !momentum flux, but is it divided by density)
    real(DP)               :: xy_PTempFlux(imin:imax,jmin:jmax)
                                           !地表面熱フラックス
                                           !(strictly speaking, this value is not 
                                           !heat flux, but is it divided by density and 
                                           !specific heat)
    real(DP)               :: xy_ExnerFlux(imin:imax,jmin:jmax)
                                           !地表面フラックス
                                           !(strictly speaking, this value is not 
                                           !heat flux, but is it divided by density and 
                                           !specific heat)
    real(DP)               :: xyf_QMixFlux(imin:imax,jmin:jmax,ncmax)
                                           !物質的フラックス
    real(DP)               :: xyz_VelX(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速 (xyz 格子)
    real(DP)               :: xyz_VelY(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速 (xyz 格子)
    real(DP)               :: xyz_AbsVel(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速 (xyz 格子)
    real(DP)               :: pyz_AbsVel(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速 (pyz 格子)
    real(DP)               :: xqz_AbsVel(imin:imax,jmin:jmax,kmin:kmax)
                                           !水平風速 (xqz 格子)
    real(DP)               :: xyz_PTempAll(imin:imax,jmin:jmax,kmin:kmax)
                                           !Total value of potential temperature
    real(DP)               :: xyz_TempAll (imin:imax,jmin:jmax,kmin:kmax)
                                           !Total value of temperature
    real(DP)               :: xyzf_QMixAll(imin:imax,jmin:jmax,kmin:kmax, ncmax)
                                           !Total value of mixing ratios
    real(DP)               :: xy_DPTempDtBulk(imin:imax,jmin:jmax)
                                        !potential temperature tendency by surface flux
    real(DP)               :: xy_DExnerDtBulk(imin:imax,jmin:jmax)
                                        !exner function tendency by surface flux
    real(DP)               :: xyf_DQMixDtBulk(imin:imax,jmin:jmax, ncmax)
                                        !mixing ratio tendency by surface flux
    real(DP)               :: py_DVelXDtBulk (imin:imax,jmin:jmax)
                                        !x-component velocity tendency by surface flux
    real(DP)               :: xq_DVelYDtBulk (imin:imax,jmin:jmax)
                                        !y-component velocity tendency by surface flux

    real(DP)               :: xyz_DPTempDtBulk(imin:imax,jmin:jmax,kmin:kmax)
                                        ! variable for output
    real(DP)               :: xyz_DExnerDtBulk(imin:imax,jmin:jmax,kmin:kmax)
                                        ! variable for output
    real(DP)               :: xyzf_DQMixDtBulk(imin:imax,jmin:jmax,kmin:kmax, ncmax)
                                        ! variable for output
    real(DP)               :: pyz_DVelXDtBulk (imin:imax,jmin:jmax,kmin:kmax)
                                        ! variable for output
    real(DP)               :: xqz_DVelYDtBulk (imin:imax,jmin:jmax,kmin:kmax)
                                        ! variable for output

    real(DP)               :: ExnerBZSfc
                                        ! Basic state Exner function at the surface
    real(DP)               :: xy_PressSfc(imin:imax,jmin:jmax)
                                        ! Total pressure at the surface

    integer                :: kz            !配列添字
    integer                :: s             !ループ変数

    ! 初期化
    !
    kz = 1

    xyz_PTempAll  = xyz_PTemp + xyz_PTempBZ
    xyz_TempAll   = (xyz_Exner + xyz_ExnerBZ) * xyz_PTempAll
    xyzf_QMixAll  = xyzf_QMix + xyzf_QMixBZ

    ExnerBZSfc    = (PressSfc / PressBasis) ** (GasRDry / CpDry)
    xy_PressSfc   = PressBasis * ( ExnerBZSfc + xyz_Exner(:,:,kz) )**(CpDry / GasRDry)
                    ! Perturbation component of Exner function at the surface is assumed 
                    ! to be same as that at the lowest layer. (YOT, 2011/09/03)


    ! Velocities at xyz grid points are calculated.
    xyz_VelX = xyz_avr_pyz(pyz_VelX)
    xyz_VelY = xyz_avr_xqz(xqz_VelY)

    xyz_AbsVel = SQRT( xyz_VelX**2 + xyz_VelY**2 + Vel0**2 )
    pyz_AbsVel = pyz_avr_xyz(xyz_AbsVel)
    xqz_AbsVel = xqz_avr_xyz(xyz_AbsVel)


    ! Something like heat, mass, and momentum fluxes are calculated.
    ! The values below are not heat, mass, and momentum fluxes, but are those divided by
    ! by density and specific heat, density, and density, respectively.
    !
    xy_PTempFlux = - Bulk * xyz_AbsVel(:,:,kz) * ( xyz_PTempAll(:,:,kz) - TempSfc / ( ExnerBZSfc + xyz_Exner(:,:,kz) ) )
    !
    xyf_QMixFlux = 0.0d0
    do s = 1, CondNum
      xyf_QMixFlux(:,:,IdxCG(s)) = - Bulk * xyz_AbsVel(:,:,kz) * ( xyzf_QMixAll(:,:,kz,s) - SvapPress( SpcWetID(IdxCC(s)), TempSfc ) / ( xy_PressSfc ) * (MolWtWet(IdxCG(s)) / MolWtDry) )
    end do
    !
    xy_ExnerFlux = xy_DExnerDt_xy_xyf(xy_PTempFlux, xyf_QMixFlux, kz)
    !
    py_VelXFlux = - Bulk * pyz_AbsVel(:,:,kz) * pyz_VelX(:,:,kz)
    !
    xq_VelYFlux = - Bulk * xqz_AbsVel(:,:,kz) * xqz_VelY(:,:,kz)


    ! Something like heat flux and mass flux are restricted.
    ! This would be arbitrary treatment. 
    xy_PTempFlux = max( 0.0d0, xy_PTempFlux )
    xy_ExnerFlux = max( 0.0d0, xy_ExnerFlux )
    xyf_QMixFlux = max( 0.0d0, xyf_QMixFlux )


    ! Tendencies by surface fluxes (convergences of fluxes) are calculated.
    !
    xy_DPTempDtBulk = - ( 0.0d0 - xy_PTempFlux ) / z_dz(kz)
    !
    xy_DExnerDtBulk = - ( 0.0d0 - xy_ExnerFlux ) / z_dz(kz)
    !
    xyf_DQMixDtBulk = - ( 0.0d0 - xyf_QMixFlux ) / z_dz(kz)
    !
    py_DVelXDtBulk = - ( 0.0d0 - py_VelXFlux ) / z_dz(kz)
    !
    xq_DVelYDtBulk = - ( 0.0d0 - xq_VelYFlux ) / z_dz(kz)


    ! Add tendency by surface flux convergence
    !
    xyz_DPTempDt(:,:,kz) = xyz_DPTempDt(:,:,kz) + xy_DPTempDtBulk
    xyz_DExnerDt(:,:,kz) = xyz_DExnerDt(:,:,kz) + xy_DExnerDtBulk
    do s = 1, ncmax
      xyzf_DQMixDt(:,:,kz,s) = xyzf_DQMixDt(:,:,kz,s) + xyf_DQMixDtBulk(:,:,s)
    end do
    pyz_DVelXDt (:,:,kz) = pyz_DVelXDt (:,:,kz) + py_DVelXDtBulk
    xqz_DVelYDt (:,:,kz) = xqz_DVelYDt (:,:,kz) + xq_DVelYDtBulk

    ! Output
    !
    xyz_DPTempDtBulk = 0.0d0
    xyz_DExnerDtBulk = 0.0d0
    xyzf_DQMixDtBulk = 0.0d0
    pyz_DVelXDtBulk  = 0.0d0
    xqz_DVelYDtBulk  = 0.0d0

    xyz_DPTempDtBulk(:,:,kz) = xy_DPTempDtBulk
    xyz_DExnerDtBulk(:,:,kz) = xy_DExnerDtBulk
    do s = 1, ncmax
      xyzf_DQMixDtBulk(:,:,kz,s) = xyf_DQMixDtBulk(:,:,s)
    end do
    pyz_DVelXDtBulk (:,:,kz) = py_DVelXDtBulk
    xqz_DVelYDtBulk (:,:,kz) = xq_DVelYDtBulk
    !
    call HistoryAutoPut(TimeN, 'PTempSfc', xyz_DPTempDtBulk(1:nx,1:ny,1:nz))
    call HistoryAutoPut(TimeN, 'ExnerSfc', xyz_DExnerDtBulk(1:nx,1:ny,1:nz))
    call HistoryAutoPut(TimeN, 'VelXSfc',  pyz_DVelXDtBulk (1:nx,1:ny,1:nz))
    call HistoryAutoPut(TimeN, 'VelYSfc',  xqz_DVelYDtBulk (1:nx,1:ny,1:nz))
    do s = 1, ncmax
      call HistoryAutoPut(TimeN, trim(SpcWetSymbol(s))//'_Sfc', xyzf_DQMixDtBulk(1:nx,1:ny,1:nz,s))
    end do


    call HistoryAutoPut(TimeN, 'PTempSfcFlux', xy_PTempFlux(1:nx,1:ny))
    call HistoryAutoPut(TimeN, 'ExnerSfcFlux', xy_ExnerFlux(1:nx,1:ny))
    call HistoryAutoPut(TimeN, 'VelXSfcFlux',  py_VelXFlux (1:nx,1:ny))
    call HistoryAutoPut(TimeN, 'VelYSfcFlux',  xq_VelYFlux (1:nx,1:ny))
    do s = 1, ncmax
      call HistoryAutoPut(TimeN, trim(SpcWetSymbol(s))//'_SfcFlux', xyf_QMixFlux(1:nx,1:ny,s))
    end do

    call HistoryAutoPut(TimeN, 'SfcHeatFlux', CpDry * xyz_DensBZ(1:nx,1:ny,1) * xy_PTempFlux(1:nx,1:ny) * ( ExnerBZSfc + xyz_Exner(1:nx,1:ny,1) ) )
    call HistoryAutoPut(TimeN, 'SfcXMomFlux', xyz_DensBZ(1:nx,1:ny,1) * py_VelXFlux (1:nx,1:ny))
    call HistoryAutoPut(TimeN, 'SfcYMomFlux', xyz_DensBZ(1:nx,1:ny,1) * xq_VelYFlux (1:nx,1:ny))
    do s = 1, ncmax
      call HistoryAutoPut(TimeN, trim(SpcWetSymbol(s))//'_SfcMassFlux', xyz_DensBZ(1:nx,1:ny,1) * xyf_QMixFlux(1:nx,1:ny,s))
    end do

    ! Set Margin
    !
!    call SetMargin_xyz( xyz_DPTempDt )
!    call SetMargin_pyz( pyz_DVelXDt )
!    call SetMargin_xqz( xqz_DVelYDt )
!    call SetMargin_xyzf(xyzf_DQMixDt )

  end subroutine Surfaceflux_Bulk_forcing
Subroutine :

NAMELIST から必要な情報を読み取り, 時間関連の変数の設定を行う.

This procedure input/output NAMELIST#surfaceflux_bulk_nml .

[Source]

  subroutine Surfaceflux_Bulk_init
    !
    !NAMELIST から必要な情報を読み取り, 時間関連の変数の設定を行う. 
    !

    !暗黙の型宣言禁止
    implicit none

    !内部変数
    integer    :: l, unit

    !---------------------------------------------------------------
    ! NAMELIST から情報を取得
    !
    NAMELIST /surfaceflux_bulk_nml/ Bulk, Vel0

    call FileOpen(unit, file=namelist_filename, mode='r')
    read(unit, NML=surfaceflux_bulk_nml)
    close(unit)  

    if (myrank == 0) then 
      call MessageNotify( "M", module_name, "Bulk = %f", d=(/Bulk/) )
      call MessageNotify( "M", module_name, "Vel0 = %f", d=(/Vel0/))
    end if


    call HistoryAutoAddVariable( varname='PTempSfcFlux', dims=(/'x','y','t'/), longname='surface potential temperature flux (heat flux divided by density and specific heat)', units='K.m.s-1', xtype='float')

    call HistoryAutoAddVariable( varname='ExnerSfcFlux', dims=(/'x','y','t'/), longname='surface exner function flux (heat flux divided by density and specific heat)', units='s-1', xtype='float')

    call HistoryAutoAddVariable( varname='VelXSfcFlux', dims=(/'x','y','t'/), longname='surface flux of x-component of velocity (momentum flux divided by density)', units='m2.s-2', xtype='float')

    call HistoryAutoAddVariable( varname='VelYSfcFlux', dims=(/'x','y','t'/), longname='surface flux of y-component of velocity (momentum flux divided by density)', units='m2.s-2', xtype='float')

    do l = 1, ncmax
      call HistoryAutoAddVariable( varname=trim(SpcWetSymbol(l))//'_SfcFlux', dims=(/'x','y','t'/), longname='surface flux of ' //trim(SpcWetSymbol(l))//' mixing ratio (mass flux divided by density)', units='m.s-1', xtype='float')
    end do


    call HistoryAutoAddVariable( varname='PTempSfc', dims=(/'x','y','z','t'/), longname='potential temperature tendency by surface flux', units='K.s-1', xtype='float')

    call HistoryAutoAddVariable( varname='ExnerSfc', dims=(/'x','y','z','t'/), longname='exner function tendency by surface flux', units='K.s-1', xtype='float')

    call HistoryAutoAddVariable( varname='VelXSfc', dims=(/'x','y','z','t'/), longname='x-component velocity tendency by surface flux', units='m.s-2', xtype='float')

    call HistoryAutoAddVariable( varname='VelYSfc', dims=(/'x','y','z','t'/), longname='y-component velocity tendency by surface flux', units='m.s-2', xtype='float')

    do l = 1, ncmax
      call HistoryAutoAddVariable( varname=trim(SpcWetSymbol(l))//'_Sfc', dims=(/'x','y','z','t'/), longname=trim(SpcWetSymbol(l))//' mixing ratio tendency by surface flux', units='s-1', xtype='float')
    end do


    call HistoryAutoAddVariable( varname='SfcHeatFlux', dims=(/'x','y','t'/), longname='surface heat flux', units='W.m-2', xtype='float')

    call HistoryAutoAddVariable( varname='SfcXMomFlux', dims=(/'x','y','t'/), longname='surface x-component momentum flux', units='kg.m-2.s-1', xtype='float')

    call HistoryAutoAddVariable( varname='SfcYMomFlux', dims=(/'x','y','t'/), longname='surface y-component momentum flux', units='kg.m-2.s-1', xtype='float')

    do l = 1, ncmax
      call HistoryAutoAddVariable( varname=trim(SpcWetSymbol(l))//'_SfcMassFlux', dims=(/'x','y','t'/), longname=trim(SpcWetSymbol(l))//' surface mass flux', units='kg.m-2.s-1', xtype='float')
    end do

  end subroutine Surfaceflux_Bulk_init