SCALE-RM
scale_atmos_phy_sf_bulk.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 #include "scalelib.h"
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_io
20  use scale_prof
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: atmos_phy_sf_bulk_setup
29  public :: atmos_phy_sf_bulk_flux
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  !-----------------------------------------------------------------------------
36  !
37  !++ Private procedure
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private parameters & variables
42  !
43  real(RP), private :: ATMOS_PHY_SF_BULK_beta = 1.0_rp ! evaporation efficiency (0-1)
44 
45  !-----------------------------------------------------------------------------
46 contains
47  !-----------------------------------------------------------------------------
49  subroutine atmos_phy_sf_bulk_setup
50  use scale_prc, only: &
51  prc_abort
52  implicit none
53 
54  namelist / param_atmos_phy_sf_bulk / &
55  atmos_phy_sf_bulk_beta
56 
57  integer :: ierr
58  !---------------------------------------------------------------------------
59 
60  log_newline
61  log_info("ATMOS_PHY_SF_bulk_setup",*) 'Setup'
62  log_info("ATMOS_PHY_SF_bulk_setup",*) 'Bulk scheme'
63 
64  !--- read namelist
65  rewind(io_fid_conf)
66  read(io_fid_conf,nml=param_atmos_phy_sf_bulk,iostat=ierr)
67  if( ierr < 0 ) then !--- missing
68  log_info("ATMOS_PHY_SF_bulk_setup",*) 'Not found namelist. Default used.'
69  elseif( ierr > 0 ) then !--- fatal error
70  log_error("ATMOS_PHY_SF_bulk_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_SF_BULK. Check!'
71  call prc_abort
72  endif
73  log_nml(param_atmos_phy_sf_bulk)
74 
75  return
76  end subroutine atmos_phy_sf_bulk_setup
77 
78  !-----------------------------------------------------------------------------
80  subroutine atmos_phy_sf_bulk_flux( &
81  IA, IS, IE, JA, JS, JE, &
82  ATM_W, ATM_U, ATM_V, &
83  ATM_TEMP, ATM_PRES, ATM_QV, &
84  SFC_DENS, SFC_TEMP, SFC_PRES, &
85  SFC_Z0M, SFC_Z0H, SFC_Z0E, &
86  PBL, ATM_Z1, &
87  SFLX_MW, SFLX_MU, SFLX_MV, &
88  SFLX_SH, SFLX_LH, SFLX_QV, &
89  Ustar, Tstar, Qstar, Wstar, &
90  RLmo, &
91  U10, V10, T2, Q2 )
92  use scale_const, only: &
93  eps => const_eps, &
94  cpdry => const_cpdry, &
95  epsvap => const_epsvap
96  use scale_atmos_hydrometeor, only: &
97  hydrometeor_lhv => atmos_hydrometeor_lhv
98  use scale_atmos_saturation, only: &
99  saturation_psat_all => atmos_saturation_psat_all
100  use scale_bulkflux, only: &
101  bulkflux, &
102  bulkflux_diagnose_surface
103  implicit none
104  integer, intent(in) :: ia, is, ie
105  integer, intent(in) :: ja, js, je
106 
107  real(rp), intent(in) :: atm_w (ia,ja) ! velocity w at the lowermost layer (cell center) [m/s]
108  real(rp), intent(in) :: atm_u (ia,ja) ! velocity u at the lowermost layer (cell center) [m/s]
109  real(rp), intent(in) :: atm_v (ia,ja) ! velocity v at the lowermost layer (cell center) [m/s]
110  real(rp), intent(in) :: atm_temp(ia,ja) ! temperature at the lowermost layer (cell center) [K]
111  real(rp), intent(in) :: atm_pres(ia,ja) ! pressure at the lowermost layer (cell center) [Pa]
112  real(rp), intent(in) :: atm_qv (ia,ja) ! qv at the lowermost layer (cell center) [kg/kg]
113  real(rp), intent(in) :: sfc_dens(ia,ja) ! density at the surface atmosphere [kg/m3]
114  real(rp), intent(in) :: sfc_temp(ia,ja) ! temperature at the surface skin [K]
115  real(rp), intent(in) :: sfc_pres(ia,ja) ! pressure at the surface atmosphere [Pa]
116  real(rp), intent(in) :: sfc_z0m (ia,ja) ! surface roughness length (momentum) [m]
117  real(rp), intent(in) :: sfc_z0h (ia,ja) ! surface roughness length (heat) [m]
118  real(rp), intent(in) :: sfc_z0e (ia,ja) ! surface roughness length (vapor) [m]
119  real(rp), intent(in) :: pbl (ia,ja) ! depth of the PBL [m]
120  real(rp), intent(in) :: atm_z1 (ia,ja) ! height of the lowermost grid from surface (cell center) [m]
121 
122  real(rp), intent(out) :: sflx_mw(ia,ja) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
123  real(rp), intent(out) :: sflx_mu(ia,ja) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
124  real(rp), intent(out) :: sflx_mv(ia,ja) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
125  real(rp), intent(out) :: sflx_sh(ia,ja) ! surface flux for sensible heat (area center) [J/m2/s]
126  real(rp), intent(out) :: sflx_lh(ia,ja) ! surface flux for latent heat (area center) [J/m2/s]
127  real(rp), intent(out) :: sflx_qv(ia,ja) ! surface flux for qv (area center) [kg/m2/s]
128  real(rp), intent(out) :: ustar (ia,ja) ! friction velocity
129  real(rp), intent(out) :: tstar (ia,ja) ! temperatuer scale
130  real(rp), intent(out) :: qstar (ia,ja) ! moisture scale
131  real(rp), intent(out) :: wstar (ia,ja) ! convective veolocity scale
132  real(rp), intent(out) :: rlmo (ia,ja) ! inverse of Obukhov length
133  real(rp), intent(out) :: u10 (ia,ja) ! velocity u at 10m height
134  real(rp), intent(out) :: v10 (ia,ja) ! velocity v at 10m height
135  real(rp), intent(out) :: t2 (ia,ja) ! temperature t at 2m height
136  real(rp), intent(out) :: q2 (ia,ja) ! water vapor q at 2m height
137 
138  real(rp) :: sfc_psat (ia,ja) ! saturatad water vapor pressure [Pa]
139  real(rp) :: lhv (ia,ja)
140 
141  real(rp) :: uabs ! modified absolute velocity [m/s]
142  real(rp) :: ra ! Aerodynamic resistance (=1/Ce) [1/s]
143  real(rp) :: sfc_qsat ! saturatad water vapor mixing ratio [kg/kg]
144  real(rp) :: sfc_qv(ia,ja) ! water vapor mixing ratio [kg/kg]
145 
146  real(rp) :: fracu10(ia,ja) ! calculation parameter for U10 [-]
147  real(rp) :: fract2 (ia,ja) ! calculation parameter for T2 [-]
148  real(rp) :: fracq2 (ia,ja) ! calculation parameter for Q2 [-]
149 
150  real(rp) :: mflux
151 
152  integer :: i, j
153  !---------------------------------------------------------------------------
154 
155  log_progress(*) 'atmosphere / physics / surface flux / bulk'
156 
157  call hydrometeor_lhv( ia, is, ie, ja, js, je, &
158  sfc_temp(:,:), & ! [IN]
159  lhv(:,:) ) ! [OUT]
160 
161  call saturation_psat_all( ia, is, ie, ja, js, je, &
162  sfc_temp(:,:), & ! [IN]
163  sfc_psat(:,:) ) ! [OUT]
164 
165  !$omp parallel do &
166 #ifndef __GFORTRAN__
167  !$omp default(none) &
168  !$omp shared (IS,IE,JS,JE,EPSvap,ATMOS_PHY_SF_BULK_beta,EPS,CPdry,LHV,bulkflux,&
169  !$omp ATM_TEMP,ATM_PRES,ATM_QV,ATM_W,ATM_U,ATM_V,ATM_Z1,PBL, &
170  !$omp SFC_DENS,SFC_TEMP,SFC_PRES,SFC_QV,SFC_PSAT,SFC_Z0M,SFC_Z0H,SFC_Z0E, &
171  !$omp SFLX_MW,SFLX_MU,SFLX_MV,SFLX_SH,SFLX_LH,SFLX_QV, &
172  !$omp FracU10,FracT2,FracQ2, &
173  !$omp Ustar,Tstar,Qstar,Wstar,RLmo,U10,V10,T2,Q2) &
174 #else
175  !$omp default(shared) &
176 #endif
177  !$omp private(SFC_QSAT,Uabs,Ra,MFLUX)
178  do j = js, je
179  do i = is, ie
180  ! qdry = 1 - psat
181  sfc_qsat = epsvap * sfc_psat(i,j) / ( sfc_pres(i,j) - ( 1.0_rp-epsvap ) * sfc_psat(i,j) )
182 
183  sfc_qv(i,j) = ( 1.0_rp - atmos_phy_sf_bulk_beta ) * atm_qv(i,j) + atmos_phy_sf_bulk_beta * sfc_qsat
184  uabs = sqrt( atm_w(i,j)**2 + atm_u(i,j)**2 + atm_v(i,j)**2 )
185 
186  call bulkflux( atm_temp(i,j), sfc_temp(i,j), & ! [IN]
187  atm_pres(i,j), sfc_pres(i,j), & ! [IN]
188  atm_qv(i,j), sfc_qv(i,j), & ! [IN]
189  uabs, atm_z1(i,j), pbl(i,j), & ! [IN]
190  sfc_z0m(i,j), sfc_z0h(i,j), sfc_z0e(i,j), & ! [IN]
191  ustar(i,j), tstar(i,j), qstar(i,j), & ! [OUT]
192  wstar(i,j), rlmo(i,j), ra, & ! [OUT]
193  fracu10(i,j), fract2(i,j), fracq2(i,j) ) ! [OUT]
194 
195  !-----< momentum >-----
196  if ( uabs < eps ) then
197  sflx_mw(i,j) = 0.0_rp
198  sflx_mu(i,j) = 0.0_rp
199  sflx_mv(i,j) = 0.0_rp
200  else
201  mflux = -sfc_dens(i,j) * ustar(i,j)**2
202  sflx_mw(i,j) = mflux * atm_w(i,j) / uabs
203  sflx_mu(i,j) = mflux * atm_u(i,j) / uabs
204  sflx_mv(i,j) = mflux * atm_v(i,j) / uabs
205  end if
206 
207  !-----< heat flux >-----
208  sflx_sh(i,j) = -sfc_dens(i,j) * ustar(i,j) * tstar(i,j) * cpdry
209  sflx_lh(i,j) = -sfc_dens(i,j) * ustar(i,j) * qstar(i,j) * lhv(i,j)
210 
211  !-----< mass flux >-----
212  sflx_qv(i,j) = sflx_lh(i,j) / lhv(i,j)
213 
214  enddo
215  enddo
216 
217  call bulkflux_diagnose_surface( ia, is, ie, ja, js, je, &
218  atm_u(:,:), atm_v(:,:), &
219  atm_temp(:,:), atm_qv(:,:), &
220  sfc_temp(:,:), sfc_qv(:,:), &
221  atm_z1(:,:), &
222  sfc_z0m(:,:), sfc_z0h(:,:), sfc_z0e(:,:), &
223  u10(:,:), v10(:,:), t2(:,:), q2(:,:), &
224  fracu10 = fracu10(:,:), &
225  fract2 = fract2(:,:), fracq2 = fracq2(:,:) )
226 
227  return
228  end subroutine atmos_phy_sf_bulk_flux
229 
230 end module scale_atmos_phy_sf_bulk
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_flux
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.
Definition: scale_atmos_phy_sf_bulk.F90:92
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:69
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_phy_sf_bulk
module atmosphere / physics / surface / bulk
Definition: scale_atmos_phy_sf_bulk.F90:13
scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_setup
subroutine, public atmos_phy_sf_bulk_setup
Setup.
Definition: scale_atmos_phy_sf_bulk.F90:50
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_io
module STDIO
Definition: scale_io.F90:10
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:56
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:56
scale_bulkflux::bulkflux
procedure(bc), pointer, public bulkflux
Definition: scale_bulkflux.F90:78