SCALE-RM
Functions/Subroutines
scale_atmos_phy_sf_bulk Module Reference

module ATMOSPHERE / Physics Surface fluxes More...

Functions/Subroutines

subroutine, public atmos_phy_sf_bulk_setup (ATMOS_PHY_SF_TYPE)
 Setup. More...
 
subroutine, public atmos_phy_sf_bulk (ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, ATM_DENS, ATM_QTRC, ATM_Z1, dt, SFC_DENS, SFC_PRES, SFLX_LW_dn, SFLX_SW_dn, SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QTRC, U10, V10, T2, Q2)
 Calculate surface flux. More...
 

Detailed Description

module ATMOSPHERE / Physics Surface fluxes

Description
Flux from/to bottom wall of atmosphere (surface) Bulk Method
Author
Team SCALE
History
  • 2011-12-03 (Y.Miyamoto) [new]
  • 2011-12-11 (H.Yashiro) [mod] integrate to SCALE-LES ver.3
  • 2012-03-23 (H.Yashiro) [mod] Explicit index parameter inclusion
  • 2012-04-10 (Y.Miyamoto) [mod] introduce coefficients for interpolation
  • 2012-09-11 (S.Nishizawa) [mod] bugfix based on the scale document
NAMELIST
  • PARAM_ATMOS_PHY_SF_BULK
    nametypedefault valuecomment
    ATMOS_PHY_SF_BETA real(RP) 1.0_RP evaporation efficiency (0-1)

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_sf_bulk_setup()

subroutine, public scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_setup ( character(len=*), intent(in)  ATMOS_PHY_SF_TYPE)

Setup.

Definition at line 58 of file scale_atmos_phy_sf_bulk.F90.

References scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, and scale_process::prc_mpistop().

Referenced by scale_atmos_phy_sf::atmos_phy_sf_setup().

58  use scale_process, only: &
60  implicit none
61 
62  character(len=*), intent(in) :: ATMOS_PHY_SF_TYPE
63 
64  namelist / param_atmos_phy_sf_bulk / &
65  atmos_phy_sf_beta
66 
67  integer :: ierr
68  !---------------------------------------------------------------------------
69 
70  if( io_l ) write(io_fid_log,*)
71  if( io_l ) write(io_fid_log,*) '++++++ Module[SURFACE FLUX] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
72  if( io_l ) write(io_fid_log,*) '*** Bulk scheme'
73 
74  if ( atmos_phy_sf_type /= 'BULK' ) then
75  write(*,*) 'xxx ATMOS_PHY_SF_TYPE is not BULK. Check!'
76  call prc_mpistop
77  endif
78 
79  !--- read namelist
80  rewind(io_fid_conf)
81  read(io_fid_conf,nml=param_atmos_phy_sf_bulk,iostat=ierr)
82  if( ierr < 0 ) then !--- missing
83  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
84  elseif( ierr > 0 ) then !--- fatal error
85  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_BULK. Check!'
86  call prc_mpistop
87  endif
88  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_sf_bulk)
89 
90  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_bulk()

subroutine, public scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk ( real(rp), dimension (ia,ja), intent(in)  ATM_TEMP,
real(rp), dimension (ia,ja), intent(in)  ATM_PRES,
real(rp), dimension (ia,ja), intent(in)  ATM_W,
real(rp), dimension (ia,ja), intent(in)  ATM_U,
real(rp), dimension (ia,ja), intent(in)  ATM_V,
real(rp), dimension (ia,ja), intent(in)  ATM_DENS,
real(rp), dimension (ia,ja,qa), intent(in)  ATM_QTRC,
real(rp), dimension (ia,ja), intent(in)  ATM_Z1,
real(dp), intent(in)  dt,
real(rp), dimension (ia,ja), intent(in)  SFC_DENS,
real(rp), dimension (ia,ja), intent(in)  SFC_PRES,
real(rp), dimension(ia,ja), intent(in)  SFLX_LW_dn,
real(rp), dimension(ia,ja), intent(in)  SFLX_SW_dn,
real(rp), dimension (ia,ja), intent(in)  SFC_TEMP,
real(rp), dimension(ia,ja,2), intent(in)  SFC_albedo,
real(rp), dimension (ia,ja), intent(inout)  SFC_Z0M,
real(rp), dimension (ia,ja), intent(inout)  SFC_Z0H,
real(rp), dimension (ia,ja), intent(inout)  SFC_Z0E,
real(rp), dimension (ia,ja), intent(out)  SFLX_MW,
real(rp), dimension (ia,ja), intent(out)  SFLX_MU,
real(rp), dimension (ia,ja), intent(out)  SFLX_MV,
real(rp), dimension (ia,ja), intent(out)  SFLX_SH,
real(rp), dimension (ia,ja), intent(out)  SFLX_LH,
real(rp), dimension (ia,ja,qa), intent(out)  SFLX_QTRC,
real(rp), dimension (ia,ja), intent(out)  U10,
real(rp), dimension (ia,ja), intent(out)  V10,
real(rp), dimension (ia,ja), intent(out)  T2,
real(rp), dimension (ia,ja), intent(out)  Q2 
)

Calculate surface flux.

Definition at line 107 of file scale_atmos_phy_sf_bulk.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_pre00, scale_const::const_rdry, scale_atmos_hydrometeor::i_qv, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, and scale_roughness::roughness.

Referenced by scale_atmos_phy_sf::atmos_phy_sf_setup().

107  use scale_grid_index
108  use scale_tracer
109  use scale_const, only: &
110  pre00 => const_pre00, &
111  cpdry => const_cpdry, &
112  rdry => const_rdry
113  use scale_atmos_hydrometeor, only: &
114  hydrometeor_lhv => atmos_hydrometeor_lhv, &
115  i_qv
116  use scale_atmos_saturation, only: &
117  saturation_pres2qsat_all => atmos_saturation_pres2qsat_all
118  use scale_roughness, only: &
119  roughness
120  use scale_bulkflux, only: &
121  bulkflux
122  implicit none
123 
124  real(RP), intent(in) :: ATM_TEMP (IA,JA) ! temperature at the lowermost layer (cell center) [K]
125  real(RP), intent(in) :: ATM_PRES (IA,JA) ! pressure at the lowermost layer (cell center) [Pa]
126  real(RP), intent(in) :: ATM_W (IA,JA) ! velocity w at the lowermost layer (cell center) [m/s]
127  real(RP), intent(in) :: ATM_U (IA,JA) ! velocity u at the lowermost layer (cell center) [m/s]
128  real(RP), intent(in) :: ATM_V (IA,JA) ! velocity v at the lowermost layer (cell center) [m/s]
129  real(RP), intent(in) :: ATM_DENS (IA,JA) ! density at the lowermost layer (cell center) [kg/m3]
130  real(RP), intent(in) :: ATM_QTRC (IA,JA,QA) ! tracer at the lowermost layer (cell center) [kg/kg]
131  real(RP), intent(in) :: ATM_Z1 (IA,JA) ! height of the lowermost grid from surface (cell center) [m]
132  real(DP), intent(in) :: dt ! delta time
133  real(RP), intent(in) :: SFC_DENS (IA,JA) ! density at the surface atmosphere [kg/m3]
134  real(RP), intent(in) :: SFC_PRES (IA,JA) ! pressure at the surface atmosphere [Pa]
135  real(RP), intent(in) :: SFLX_LW_dn(IA,JA) ! downward longwave radiation flux at the surface [J/m2/s]
136  real(RP), intent(in) :: SFLX_SW_dn(IA,JA) ! downward shortwave radiation flux at the surface [J/m2/s]
137  real(RP), intent(in) :: SFC_TEMP (IA,JA) ! temperature at the surface skin [K]
138  real(RP), intent(in) :: SFC_albedo(IA,JA,2) ! surface albedo (LW/SW) (0-1)
139  real(RP), intent(inout) :: SFC_Z0M (IA,JA) ! surface roughness length (momentum) [m]
140  real(RP), intent(inout) :: SFC_Z0H (IA,JA) ! surface roughness length (heat) [m]
141  real(RP), intent(inout) :: SFC_Z0E (IA,JA) ! surface roughness length (vapor) [m]
142  real(RP), intent(out) :: SFLX_MW (IA,JA) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
143  real(RP), intent(out) :: SFLX_MU (IA,JA) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
144  real(RP), intent(out) :: SFLX_MV (IA,JA) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
145  real(RP), intent(out) :: SFLX_SH (IA,JA) ! surface flux for sensible heat (area center) [J/m2/s]
146  real(RP), intent(out) :: SFLX_LH (IA,JA) ! surface flux for latent heat (area center) [J/m2/s]
147  real(RP), intent(out) :: SFLX_QTRC (IA,JA,QA) ! surface flux for tracer mass (area center) [kg/m2/s]
148  real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m height
149  real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m height
150  real(RP), intent(out) :: T2 (IA,JA) ! temperature t at 2m height
151  real(RP), intent(out) :: Q2 (IA,JA) ! water vapor q at 2m height
152 
153  real(RP) :: ATM_QV (IA,JA)
154  real(RP) :: SFC_Z0M_t(IA,JA)
155  real(RP) :: SFC_Z0H_t(IA,JA)
156  real(RP) :: SFC_Z0E_t(IA,JA)
157  real(RP) :: SFC_QSAT (IA,JA) ! saturatad water vapor mixing ratio [kg/kg]
158  real(RP) :: SFC_QV (IA,JA) ! water vapor mixing ratio [kg/kg]
159  real(RP) :: LHV (IA,JA)
160  real(RP) :: PBL (IA,JA)
161 
162  real(RP) :: Ustar ! friction velocity [m]
163  real(RP) :: Tstar ! friction temperature [K]
164  real(RP) :: Qstar ! friction mixing rate [kg/kg]
165  real(RP) :: Uabs ! modified absolute velocity [m/s]
166  real(RP) :: FracU10 ! calculation parameter for U10 [-]
167  real(RP) :: FracT2 ! calculation parameter for T2 [-]
168  real(RP) :: FracQ2 ! calculation parameter for Q2 [-]
169 
170  integer :: i, j
171  !---------------------------------------------------------------------------
172 
173  if( io_l ) write(io_fid_log,*) '*** Atmos physics step: Surface flux(bulk)'
174 
175  call roughness( sfc_z0m_t(:,:), & ! [OUT]
176  sfc_z0h_t(:,:), & ! [OUT]
177  sfc_z0e_t(:,:), & ! [OUT]
178  sfc_z0m(:,:), & ! [IN]
179  sfc_z0h(:,:), & ! [IN]
180  sfc_z0e(:,:), & ! [IN]
181  atm_u(:,:), & ! [IN]
182  atm_v(:,:), & ! [IN]
183  atm_z1(:,:), & ! [IN]
184  dt ) ! [IN]
185 
186  sfc_z0m(:,:) = sfc_z0m(:,:) + sfc_z0m_t(:,:) * dt
187  sfc_z0h(:,:) = sfc_z0h(:,:) + sfc_z0h_t(:,:) * dt
188  sfc_z0e(:,:) = sfc_z0e(:,:) + sfc_z0e_t(:,:) * dt
189 
190  call hydrometeor_lhv( lhv(:,:), atm_temp(:,:) )
191 
192  if ( i_qv > 0 ) then
193  atm_qv(:,:) = atm_qtrc(:,:,i_qv)
194 
195  call saturation_pres2qsat_all( sfc_qsat(:,:), & ! [OUT]
196  sfc_temp(:,:), & ! [IN]
197  sfc_pres(:,:) ) ! [IN]
198 
199  sfc_qv(:,:) = ( 1.0_rp - atmos_phy_sf_beta ) * atm_qv(:,:) + atmos_phy_sf_beta * sfc_qsat(:,:)
200  else
201  atm_qv(:,:) = 0.0_rp
202  sfc_qv(:,:) = 0.0_rp
203  end if
204 
205  sflx_qtrc(:,:,:) = 0.0_rp
206  pbl(:,:) = 100.0_rp ! tentative
207  do j = js, je
208  do i = is, ie
209  call bulkflux( ustar, & ! [OUT]
210  tstar, & ! [OUT]
211  qstar, & ! [OUT]
212  uabs, & ! [OUT]
213  fracu10, & ! [OUT]
214  fract2, & ! [OUT]
215  fracq2, & ! [OUT]
216  atm_temp(i,j), & ! [IN]
217  sfc_temp(i,j), & ! [IN]
218  atm_pres(i,j), & ! [IN]
219  sfc_pres(i,j), & ! [IN]
220  atm_qv(i,j), & ! [IN]
221  sfc_qv(i,j), & ! [IN]
222  atm_u(i,j), & ! [IN]
223  atm_v(i,j), & ! [IN]
224  atm_z1(i,j), & ! [IN]
225  pbl(i,j), & ! [IN]
226  sfc_z0m(i,j), & ! [IN]
227  sfc_z0h(i,j), & ! [IN]
228  sfc_z0e(i,j) ) ! [IN]
229 
230  !-----< momentum >-----
231  sflx_mw(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_w(i,j)
232  sflx_mu(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_u(i,j)
233  sflx_mv(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_v(i,j)
234 
235  !-----< heat flux >-----
236  sflx_sh(i,j) = - atm_dens(i,j) * ustar * tstar &
237  * cpdry * ( sfc_pres(i,j) / pre00 )**( rdry/cpdry )
238  sflx_lh(i,j) = - atm_dens(i,j) * ustar * qstar * lhv(i,j)
239 
240  !-----< mass flux >-----
241  if ( i_qv > 0 ) then
242  sflx_qtrc(i,j,i_qv) = sflx_lh(i,j) / lhv(i,j)
243  endif
244 
245  !-----< U10, T2, q2 >-----
246  !U10(i,j) = FracU10 * ATM_U(i,j)
247  !V10(i,j) = FracU10 * ATM_V(i,j)
248  !T2 (i,j) = ( 1.0_RP - FracT2 ) * SFC_TEMP(i,j) + FracT2 * ATM_TEMP(i,j)
249  !Q2 (i,j) = ( 1.0_RP - FracQ2 ) * SFC_QV (i,j) + FracQ2 * ATM_QV (i,j)
250 
251  u10(i,j) = atm_u(i,j) * log( 10.0_rp / sfc_z0m(i,j) ) / log( atm_z1(i,j) / sfc_z0m(i,j) )
252  v10(i,j) = atm_v(i,j) * log( 10.0_rp / sfc_z0m(i,j) ) / log( atm_z1(i,j) / sfc_z0m(i,j) )
253  t2(i,j) = sfc_temp(i,j) + ( atm_temp(i,j) - sfc_temp(i,j) ) &
254  * ( log( 2.0_rp / sfc_z0m(i,j) ) * log( 2.0_rp / sfc_z0h(i,j) ) ) &
255  / ( log( atm_z1(i,j) / sfc_z0m(i,j) ) * log( atm_z1(i,j) / sfc_z0h(i,j) ) )
256  q2(i,j) = sfc_qv(i,j) + ( atm_qv(i,j) - sfc_qv(i,j) ) &
257  * ( log( 2.0_rp / sfc_z0m(i,j) ) * log( 2.0_rp / sfc_z0e(i,j) ) ) &
258  / ( log( atm_z1(i,j) / sfc_z0m(i,j) ) * log( atm_z1(i,j) / sfc_z0e(i,j) ) )
259  enddo
260  enddo
261 
262  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
module TRACER
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
procedure(rl), pointer, public roughness
module CONSTANT
Definition: scale_const.F90:14
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module Surface roughness length
Here is the caller graph for this function: