SCALE-RM
Functions/Subroutines
scale_atmos_phy_sf_bulk Module Reference

module atmosphere / physics / surface / bulk More...

Functions/Subroutines

subroutine, public atmos_phy_sf_bulk_setup
 Setup. More...
 
subroutine, public atmos_phy_sf_bulk_flux (IA, IS, IE, JA, JS, JE, ATM_W, ATM_U, ATM_V, ATM_TEMP, ATM_PRES, ATM_QV, SFC_DENS, SFC_TEMP, SFC_PRES, SFC_Z0M, SFC_Z0H, SFC_Z0E, PBL, ATM_Z1, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QV, Ustar, Tstar, Qstar, Wstar, RLmo, U10, V10, T2, Q2)
 Calculate surface flux. More...
 

Detailed Description

module atmosphere / physics / surface / bulk

Description
Flux from/to bottom wall of atmosphere (surface) Bulk Method
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_SF_BULK
    nametypedefault valuecomment
    ATMOS_PHY_SF_BULK_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

Setup.

Definition at line 51 of file scale_atmos_phy_sf_bulk.F90.

51  use scale_prc, only: &
52  prc_abort
53  implicit none
54 
55  namelist / param_atmos_phy_sf_bulk / &
56  atmos_phy_sf_bulk_beta
57 
58  integer :: ierr
59  !---------------------------------------------------------------------------
60 
61  log_newline
62  log_info("ATMOS_PHY_SF_bulk_setup",*) 'Setup'
63  log_info("ATMOS_PHY_SF_bulk_setup",*) 'Bulk scheme'
64 
65  !--- read namelist
66  rewind(io_fid_conf)
67  read(io_fid_conf,nml=param_atmos_phy_sf_bulk,iostat=ierr)
68  if( ierr < 0 ) then !--- missing
69  log_info("ATMOS_PHY_SF_bulk_setup",*) 'Not found namelist. Default used.'
70  elseif( ierr > 0 ) then !--- fatal error
71  log_error("ATMOS_PHY_SF_bulk_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_SF_BULK. Check!'
72  call prc_abort
73  endif
74  log_nml(param_atmos_phy_sf_bulk)
75 
76  !$acc update device(ATMOS_PHY_SF_BULK_beta)
77 
78  return

References scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_bulk_flux()

subroutine, public scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_flux ( integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
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_TEMP,
real(rp), dimension(ia,ja), intent(in)  ATM_PRES,
real(rp), dimension (ia,ja), intent(in)  ATM_QV,
real(rp), dimension(ia,ja), intent(in)  SFC_DENS,
real(rp), dimension(ia,ja), intent(in)  SFC_TEMP,
real(rp), dimension(ia,ja), intent(in)  SFC_PRES,
real(rp), dimension (ia,ja), intent(in)  SFC_Z0M,
real(rp), dimension (ia,ja), intent(in)  SFC_Z0H,
real(rp), dimension (ia,ja), intent(in)  SFC_Z0E,
real(rp), dimension (ia,ja), intent(in)  PBL,
real(rp), dimension (ia,ja), intent(in)  ATM_Z1,
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), intent(out)  SFLX_QV,
real(rp), dimension (ia,ja), intent(out)  Ustar,
real(rp), dimension (ia,ja), intent(out)  Tstar,
real(rp), dimension (ia,ja), intent(out)  Qstar,
real(rp), dimension (ia,ja), intent(out)  Wstar,
real(rp), dimension (ia,ja), intent(out)  RLmo,
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 95 of file scale_atmos_phy_sf_bulk.F90.

95  use scale_const, only: &
96  eps => const_eps, &
97  cpdry => const_cpdry, &
98  epsvap => const_epsvap
99  use scale_atmos_hydrometeor, only: &
100  hydrometeor_lhv => atmos_hydrometeor_lhv
101  use scale_atmos_saturation, only: &
102  saturation_psat_all => atmos_saturation_psat_all
103  use scale_bulkflux, only: &
104  bulkflux, &
105  bulkflux_diagnose_surface
106  implicit none
107  integer, intent(in) :: IA, IS, IE
108  integer, intent(in) :: JA, JS, JE
109 
110  real(RP), intent(in) :: ATM_W (IA,JA) ! velocity w at the lowermost layer (cell center) [m/s]
111  real(RP), intent(in) :: ATM_U (IA,JA) ! velocity u at the lowermost layer (cell center) [m/s]
112  real(RP), intent(in) :: ATM_V (IA,JA) ! velocity v at the lowermost layer (cell center) [m/s]
113  real(RP), intent(in) :: ATM_TEMP(IA,JA) ! temperature at the lowermost layer (cell center) [K]
114  real(RP), intent(in) :: ATM_PRES(IA,JA) ! pressure at the lowermost layer (cell center) [Pa]
115  real(RP), intent(in) :: ATM_QV (IA,JA) ! qv at the lowermost layer (cell center) [kg/kg]
116  real(RP), intent(in) :: SFC_DENS(IA,JA) ! density at the surface atmosphere [kg/m3]
117  real(RP), intent(in) :: SFC_TEMP(IA,JA) ! temperature at the surface skin [K]
118  real(RP), intent(in) :: SFC_PRES(IA,JA) ! pressure at the surface atmosphere [Pa]
119  real(RP), intent(in) :: SFC_Z0M (IA,JA) ! surface roughness length (momentum) [m]
120  real(RP), intent(in) :: SFC_Z0H (IA,JA) ! surface roughness length (heat) [m]
121  real(RP), intent(in) :: SFC_Z0E (IA,JA) ! surface roughness length (vapor) [m]
122  real(RP), intent(in) :: PBL (IA,JA) ! depth of the PBL [m]
123  real(RP), intent(in) :: ATM_Z1 (IA,JA) ! height of the lowermost grid from surface (cell center) [m]
124 
125  real(RP), intent(out) :: SFLX_MW(IA,JA) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
126  real(RP), intent(out) :: SFLX_MU(IA,JA) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
127  real(RP), intent(out) :: SFLX_MV(IA,JA) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
128  real(RP), intent(out) :: SFLX_SH(IA,JA) ! surface flux for sensible heat (area center) [J/m2/s]
129  real(RP), intent(out) :: SFLX_LH(IA,JA) ! surface flux for latent heat (area center) [J/m2/s]
130  real(RP), intent(out) :: SFLX_QV(IA,JA) ! surface flux for qv (area center) [kg/m2/s]
131  real(RP), intent(out) :: Ustar (IA,JA) ! friction velocity
132  real(RP), intent(out) :: Tstar (IA,JA) ! temperatuer scale
133  real(RP), intent(out) :: Qstar (IA,JA) ! moisture scale
134  real(RP), intent(out) :: Wstar (IA,JA) ! convective veolocity scale
135  real(RP), intent(out) :: RLmo (IA,JA) ! inverse of Obukhov length
136  real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m height
137  real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m height
138  real(RP), intent(out) :: T2 (IA,JA) ! temperature t at 2m height
139  real(RP), intent(out) :: Q2 (IA,JA) ! water vapor q at 2m height
140 
141  real(RP) :: SFC_PSAT (IA,JA) ! saturatad water vapor pressure [Pa]
142  real(RP) :: LHV (IA,JA)
143 
144  real(RP) :: Uabs ! modified absolute velocity [m/s]
145  real(RP) :: Ra ! Aerodynamic resistance (=1/Ce) [1/s]
146  real(RP) :: SFC_QSAT ! saturatad water vapor mixing ratio [kg/kg]
147  real(RP) :: SFC_QV(IA,JA) ! water vapor mixing ratio [kg/kg]
148 
149  real(RP) :: FracU10(IA,JA) ! calculation parameter for U10 [-]
150  real(RP) :: FracT2 (IA,JA) ! calculation parameter for T2 [-]
151  real(RP) :: FracQ2 (IA,JA) ! calculation parameter for Q2 [-]
152 
153  real(RP) :: MFLUX
154 
155  integer :: i, j
156  !---------------------------------------------------------------------------
157 
158  log_progress(*) 'atmosphere / physics / surface flux / bulk'
159 
160  !$acc data copyin(ATM_W,ATM_U,ATM_V,ATM_TEMP,ATM_PRES,ATM_QV,SFC_DENS,SFC_TEMP,SFC_PRES,SFC_Z0M,SFC_Z0H,SFC_Z0E,PBL,ATM_Z1) &
161  !$acc copyout(SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_QV,Ustar,Tstar,Qstar,Wstar,RLmo,U10,V10,T2,Q2) &
162  !$acc create(SFC_PSAT,LHV,SFC_QV,FracU10,FracT2,FracQ2)
163 
164 
165  call hydrometeor_lhv( ia, is, ie, ja, js, je, &
166  sfc_temp(:,:), & ! [IN]
167  lhv(:,:) ) ! [OUT]
168 
169  call saturation_psat_all( ia, is, ie, ja, js, je, &
170  sfc_temp(:,:), & ! [IN]
171  sfc_psat(:,:) ) ! [OUT]
172 
173  !$omp parallel do &
174 #ifndef __GFORTRAN__
175  !$omp default(none) &
176  !$omp shared (IS,IE,JS,JE,EPSvap,ATMOS_PHY_SF_BULK_beta,EPS,CPdry,LHV,bulkflux,&
177  !$omp ATM_TEMP,ATM_PRES,ATM_QV,ATM_W,ATM_U,ATM_V,ATM_Z1,PBL, &
178  !$omp SFC_DENS,SFC_TEMP,SFC_PRES,SFC_QV,SFC_PSAT,SFC_Z0M,SFC_Z0H,SFC_Z0E, &
179  !$omp SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_QV, &
180  !$omp FracU10,FracT2,FracQ2, &
181  !$omp Ustar,Tstar,Qstar,Wstar,RLmo,U10,V10,T2,Q2) &
182 #else
183  !$omp default(shared) &
184 #endif
185  !$omp private(SFC_QSAT,Uabs,Ra,MFLUX)
186  !$acc kernels
187  !$acc loop independent collapse(2) &
188  !$acc private(SFC_QSAT,Uabs,Ra,MFLUX)
189  do j = js, je
190  do i = is, ie
191  ! qdry = 1 - psat
192  sfc_qsat = epsvap * sfc_psat(i,j) / ( sfc_pres(i,j) - ( 1.0_rp-epsvap ) * sfc_psat(i,j) )
193 
194  sfc_qv(i,j) = ( 1.0_rp - atmos_phy_sf_bulk_beta ) * atm_qv(i,j) + atmos_phy_sf_bulk_beta * sfc_qsat
195  uabs = sqrt( atm_w(i,j)**2 + atm_u(i,j)**2 + atm_v(i,j)**2 )
196 
197  call bulkflux( atm_temp(i,j), sfc_temp(i,j), & ! [IN]
198  atm_pres(i,j), sfc_pres(i,j), & ! [IN]
199  atm_qv(i,j), sfc_qv(i,j), & ! [IN]
200  uabs, atm_z1(i,j), pbl(i,j), & ! [IN]
201  sfc_z0m(i,j), sfc_z0h(i,j), sfc_z0e(i,j), & ! [IN]
202  ustar(i,j), tstar(i,j), qstar(i,j), & ! [OUT]
203  wstar(i,j), rlmo(i,j), ra, & ! [OUT]
204  fracu10(i,j), fract2(i,j), fracq2(i,j) ) ! [OUT]
205 
206  !-----< momentum >-----
207  if ( uabs < eps ) then
208  sflx_mw(i,j) = 0.0_rp
209  sflx_mu(i,j) = 0.0_rp
210  sflx_mv(i,j) = 0.0_rp
211  else
212  mflux = -sfc_dens(i,j) * ustar(i,j)**2
213  sflx_mw(i,j) = mflux * atm_w(i,j) / uabs
214  sflx_mu(i,j) = mflux * atm_u(i,j) / uabs
215  sflx_mv(i,j) = mflux * atm_v(i,j) / uabs
216  end if
217 
218  !-----< heat flux >-----
219  sflx_sh(i,j) = -sfc_dens(i,j) * ustar(i,j) * tstar(i,j) * cpdry
220  sflx_lh(i,j) = -sfc_dens(i,j) * ustar(i,j) * qstar(i,j) * lhv(i,j)
221 
222  !-----< mass flux >-----
223  sflx_qv(i,j) = sflx_lh(i,j) / lhv(i,j)
224 
225  enddo
226  enddo
227  !$acc end kernels
228 
229  call bulkflux_diagnose_surface( ia, is, ie, ja, js, je, &
230  atm_u(:,:), atm_v(:,:), &
231  atm_temp(:,:), atm_qv(:,:), &
232  sfc_temp(:,:), sfc_qv(:,:), &
233  atm_z1(:,:), &
234  sfc_z0m(:,:), sfc_z0h(:,:), sfc_z0e(:,:), &
235  u10(:,:), v10(:,:), t2(:,:), q2(:,:), &
236  fracu10 = fracu10(:,:), &
237  fract2 = fract2(:,:), fracq2 = fracq2(:,:) )
238 
239  !$acc end data
240 
241  return

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_eps, and scale_const::const_epsvap.

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_mkinit(), and mod_atmos_phy_sf_driver::atmos_phy_sf_driver_calc_tendency().

Here is the caller graph for this function:
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:75
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_bulkflux::bulkflux
procedure(bc), pointer, public bulkflux
Definition: scale_bulkflux.F90:84