SCALE-RM
scale_atmos_phy_sf_bulk.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
17 !-------------------------------------------------------------------------------
18 #include "inc_openmp.h"
20  !-----------------------------------------------------------------------------
21  !
22  !++ used modules
23  !
24  use scale_precision
25  use scale_stdio
26  use scale_prof
28  use scale_tracer
29  !-----------------------------------------------------------------------------
30  implicit none
31  private
32  !-----------------------------------------------------------------------------
33  !
34  !++ Public procedure
35  !
36  public :: atmos_phy_sf_bulk_setup
37  public :: atmos_phy_sf_bulk
38 
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public parameters & variables
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private procedure
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private parameters & variables
50  !
51  real(RP), private, parameter :: atmos_phy_sf_u_maxm = 100.0_rp ! maximum limit of absolute velocity for momentum [m/s]
52  real(RP), private, parameter :: atmos_phy_sf_u_maxh = 100.0_rp ! maximum limit of absolute velocity for heat [m/s]
53  real(RP), private, parameter :: atmos_phy_sf_u_maxe = 100.0_rp ! maximum limit of absolute velocity for vapor [m/s]
54  real(RP), private :: atmos_phy_sf_u_minm = 0.0_rp ! minimum limit of absolute velocity for momentum [m/s]
55  real(RP), private :: atmos_phy_sf_u_minh = 0.0_rp ! minimum limit of absolute velocity for heat [m/s]
56  real(RP), private :: atmos_phy_sf_u_mine = 0.0_rp ! minimum limit of absolute velocity for vapor [m/s]
57 
58  !-----------------------------------------------------------------------------
59 contains
60  !-----------------------------------------------------------------------------
62  subroutine atmos_phy_sf_bulk_setup( ATMOS_PHY_SF_TYPE )
63  use scale_process, only: &
65  use scale_atmos_phy_sf_bulkcoef, only: &
66  sf_bulkcoef_setup => atmos_phy_sf_bulkcoef_setup
67  implicit none
68 
69  character(len=*), intent(in) :: ATMOS_PHY_SF_TYPE
70 
71  namelist / param_atmos_phy_sf_bulk / &
72  atmos_phy_sf_u_minm, &
73  atmos_phy_sf_u_minh, &
74  atmos_phy_sf_u_mine
75 
76  integer :: ierr
77  !---------------------------------------------------------------------------
78 
79  if( io_l ) write(io_fid_log,*)
80  if( io_l ) write(io_fid_log,*) '++++++ Module[SURFACE FLUX] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
81  if( io_l ) write(io_fid_log,*) '+++ Bulk scheme'
82 
83  if ( atmos_phy_sf_type /= 'BULK' ) then
84  write(*,*) 'xxx ATMOS_PHY_SF_TYPE is not BULK. Check!'
85  call prc_mpistop
86  endif
87 
88  !--- read namelist
89  rewind(io_fid_conf)
90  read(io_fid_conf,nml=param_atmos_phy_sf_bulk,iostat=ierr)
91  if( ierr < 0 ) then !--- missing
92  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
93  elseif( ierr > 0 ) then !--- fatal error
94  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_BULK. Check!'
95  call prc_mpistop
96  endif
97  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_sf_bulk)
98 
99  call sf_bulkcoef_setup
100 
101  return
102  end subroutine atmos_phy_sf_bulk_setup
103 
104  !-----------------------------------------------------------------------------
106  subroutine atmos_phy_sf_bulk( &
107  ATM_TEMP, ATM_PRES, ATM_W, ATM_U, ATM_V, &
108  ATM_DENS, &
109  ATM_QTRC, &
110  ATM_Z1, dt, &
111  SFC_DENS, SFC_PRES, &
112  SFLX_LW_dn, SFLX_SW_dn, &
113  SFC_TEMP, SFC_albedo, SFC_beta, &
114  SFC_Z0M, SFC_Z0H, SFC_Z0E, &
115  SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, &
116  SFLX_QTRC, &
117  U10, V10, T2, Q2 )
119  use scale_tracer
120  use scale_const, only: &
121  cpdry => const_cpdry, &
122  rdry => const_rdry
123  use scale_atmos_phy_sf_bulkcoef, only: &
124  sf_bulkcoef => atmos_phy_sf_bulkcoef
125  use scale_atmos_saturation, only: &
126  saturation_pres2qsat_all => atmos_saturation_pres2qsat_all
127  use scale_atmos_thermodyn, only: &
128  atmos_thermodyn_templhv
129  use scale_roughness, only: &
130  roughness
131  use scale_bulkflux, only: &
132  bulkflux
133  implicit none
134 
135  real(RP), intent(in) :: ATM_TEMP (ia,ja) ! temperature at the lowermost layer (cell center) [K]
136  real(RP), intent(in) :: ATM_PRES (ia,ja) ! pressure at the lowermost layer (cell center) [Pa]
137  real(RP), intent(in) :: ATM_W (ia,ja) ! velocity w at the lowermost layer (cell center) [m/s]
138  real(RP), intent(in) :: ATM_U (ia,ja) ! velocity u at the lowermost layer (cell center) [m/s]
139  real(RP), intent(in) :: ATM_V (ia,ja) ! velocity v at the lowermost layer (cell center) [m/s]
140  real(RP), intent(in) :: ATM_DENS (ia,ja) ! density at the lowermost layer (cell center) [kg/m3]
141  real(RP), intent(in) :: ATM_QTRC (ia,ja,qa) ! tracer at the lowermost layer (cell center) [kg/kg]
142  real(RP), intent(in) :: ATM_Z1 (ia,ja) ! height of the lowermost grid from surface (cell center) [m]
143  real(DP), intent(in) :: dt ! delta time
144  real(RP), intent(in) :: SFC_DENS (ia,ja) ! density at the surface atmosphere [kg/m3]
145  real(RP), intent(in) :: SFC_PRES (ia,ja) ! pressure at the surface atmosphere [Pa]
146  real(RP), intent(in) :: SFLX_LW_dn(ia,ja) ! downward longwave radiation flux at the surface [J/m2/s]
147  real(RP), intent(in) :: SFLX_SW_dn(ia,ja) ! downward shortwave radiation flux at the surface [J/m2/s]
148  real(RP), intent(in) :: SFC_TEMP (ia,ja) ! temperature at the surface skin [K]
149  real(RP), intent(in) :: SFC_albedo(ia,ja,2) ! surface albedo (LW/SW) [0-1]
150  real(RP), intent(in) :: SFC_beta (ia,ja) ! evaporation efficiency [0-1]
151  real(RP), intent(inout) :: SFC_Z0M (ia,ja) ! surface roughness length (momentum) [m]
152  real(RP), intent(inout) :: SFC_Z0H (ia,ja) ! surface roughness length (heat) [m]
153  real(RP), intent(inout) :: SFC_Z0E (ia,ja) ! surface roughness length (vapor) [m]
154  real(RP), intent(out) :: SFLX_MW (ia,ja) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
155  real(RP), intent(out) :: SFLX_MU (ia,ja) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
156  real(RP), intent(out) :: SFLX_MV (ia,ja) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
157  real(RP), intent(out) :: SFLX_SH (ia,ja) ! surface flux for sensible heat (area center) [J/m2/s]
158  real(RP), intent(out) :: SFLX_LH (ia,ja) ! surface flux for latent heat (area center) [J/m2/s]
159  real(RP), intent(out) :: SFLX_QTRC (ia,ja,qa) ! surface flux for tracer mass (area center) [kg/m2/s]
160  real(RP), intent(out) :: U10 (ia,ja) ! velocity u at 10m height
161  real(RP), intent(out) :: V10 (ia,ja) ! velocity v at 10m height
162  real(RP), intent(out) :: T2 (ia,ja) ! temperature t at 2m height
163  real(RP), intent(out) :: Q2 (ia,ja) ! water vapor q at 2m height
164 
165  real(RP) :: SFC_Z0M_t(ia,ja)
166  real(RP) :: SFC_Z0H_t(ia,ja)
167  real(RP) :: SFC_Z0E_t(ia,ja)
168 
169  real(RP) :: SFC_QSAT(ia,ja) ! saturatad water vapor mixing ratio [kg/kg]
170 
171  real(RP) :: PBL(ia,ja)
172  real(RP) :: LHV (ia,ja)
173  real(RP) :: Ustar
174  real(RP) :: Tstar
175  real(RP) :: Qstar
176  real(RP) :: Uabs
177 
178  integer :: i, j
179  !---------------------------------------------------------------------------
180 
181  if( io_l ) write(io_fid_log,*) '*** Physics step: Surface flux(bulk)'
182 
183  call roughness( sfc_z0m_t(:,:), & ! [OUT]
184  sfc_z0h_t(:,:), & ! [OUT]
185  sfc_z0e_t(:,:), & ! [OUT]
186  sfc_z0m(:,:), & ! [IN]
187  sfc_z0h(:,:), & ! [IN]
188  sfc_z0e(:,:), & ! [IN]
189  atm_u(:,:), & ! [IN]
190  atm_v(:,:), & ! [IN]
191  atm_z1(:,:), & ! [IN]
192  dt ) ! [IN]
193 
194  sfc_z0m(:,:) = sfc_z0m(:,:) + sfc_z0m_t(:,:) * dt
195  sfc_z0h(:,:) = sfc_z0h(:,:) + sfc_z0h_t(:,:) * dt
196  sfc_z0e(:,:) = sfc_z0e(:,:) + sfc_z0e_t(:,:) * dt
197 
198  call saturation_pres2qsat_all( sfc_qsat(:,:), & ! [OUT]
199  sfc_temp(:,:), & ! [IN]
200  sfc_pres(:,:) ) ! [IN]
201 
202  call atmos_thermodyn_templhv( lhv, atm_temp )
203 
204  sflx_qtrc(:,:,:) = 0.0_rp
205  pbl(:,:) = 100.0_rp ! tentative
206  do j = js, je
207  do i = is, ie
208 
209  call bulkflux( &
210  ustar, & ! [OUT]
211  tstar, & ! [OUT]
212  qstar, & ! [OUT]
213  uabs, & ! [OUT]
214  atm_temp(i,j), & ! [IN]
215  sfc_temp(i,j), & ! [IN]
216  atm_pres(i,j), & ! [IN]
217  sfc_pres(i,j), & ! [IN]
218  atm_qtrc(i,j,i_qv), & ! [IN]
219  sfc_qsat(i,j), & ! [IN]
220  atm_u(i,j), & ! [IN]
221  atm_v(i,j), & ! [IN]
222  atm_z1(i,j), & ! [IN]
223  pbl(i,j), & ! [IN]
224  sfc_z0m(i,j), & ! [IN]
225  sfc_z0h(i,j), & ! [IN]
226  sfc_z0e(i,j) ) ! [IN]
227 
228  !-----< momentum >-----
229  sflx_mw(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_w(i,j)
230  sflx_mu(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_u(i,j)
231  sflx_mv(i,j) = -atm_dens(i,j) * ustar**2 / uabs * atm_v(i,j)
232 
233  !-----< heat flux >-----
234  sflx_sh(i,j) = -cpdry * atm_dens(i,j) * ustar * tstar
235  sflx_lh(i,j) = -lhv(i,j) * atm_dens(i,j) * ustar * qstar * sfc_beta(i,j)
236 
237  !-----< mass flux >-----
238  sflx_qtrc(i,j,i_qv) = sflx_lh(i,j) / lhv(i,j)
239  enddo
240  enddo
241 
242  !-----< U10, T2, q2 >-----
243 
244  do j = js, je
245  do i = is, ie
246  u10(i,j) = atm_u(i,j) * log( 10.0_rp / sfc_z0m(i,j) ) / log( atm_z1(i,j) / sfc_z0m(i,j) )
247  v10(i,j) = atm_v(i,j) * log( 10.0_rp / sfc_z0m(i,j) ) / log( atm_z1(i,j) / sfc_z0m(i,j) )
248  t2(i,j) = sfc_temp(i,j) + ( atm_temp(i,j) - sfc_temp(i,j) ) &
249  * ( log( 2.0_rp / sfc_z0m(i,j) ) * log( 2.0_rp / sfc_z0h(i,j) ) ) &
250  / ( log( atm_z1(i,j) / sfc_z0m(i,j) ) * log( atm_z1(i,j) / sfc_z0h(i,j) ) )
251  q2(i,j) = sfc_qsat(i,j) + ( atm_qtrc(i,j,i_qv) - sfc_qsat(i,j) ) &
252  * ( log( 2.0_rp / sfc_z0m(i,j) ) * log( 2.0_rp / sfc_z0e(i,j) ) ) &
253  / ( log( atm_z1(i,j) / sfc_z0m(i,j) ) * log( atm_z1(i,j) / sfc_z0e(i,j) ) )
254  enddo
255  enddo
256 
257  return
258  end subroutine atmos_phy_sf_bulk
259 
260 end module scale_atmos_phy_sf_bulk
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public atmos_phy_sf_bulk_setup(ATMOS_PHY_SF_TYPE)
Setup.
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
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_beta, 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.
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public i_qv
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
module ATMOSPHERE / Physics Surface fluxes
procedure(rl), pointer, public roughness
procedure(bc), pointer, public atmos_phy_sf_bulkcoef
module PROCESS
subroutine, public log(type, message)
Definition: dc_log.f90:133
module CONSTANT
Definition: scale_const.F90:14
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module ATMOSPHERE / Thermodynamics
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
module PRECISION
module Surface roughness length
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
module ATMOSPHERE / Physics Surface bulk coefficient
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)