SCALE-RM
Functions/Subroutines
scale_atmos_phy_sf_const Module Reference

module atmosphere / physics / surface / const More...

Functions/Subroutines

subroutine, public atmos_phy_sf_const_setup
 Setup. More...
 
subroutine, public atmos_phy_sf_const_flux (IA, IS, IE, JA, JS, JE, ATM_W, ATM_U, ATM_V, ATM_TEMP, ATM_Z1, SFC_DENS, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_QV, U10, V10)
 Constant flux. More...
 

Detailed Description

module atmosphere / physics / surface / const

Description
Flux from/to bottom wall of atmosphere (surface) Constant flux, domain-uniform
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_SF_CONST
    nametypedefault valuecomment
    ATMOS_PHY_SF_FLG_MOM_FLUX integer 0 application type for momentum flux
    ATMOS_PHY_SF_U_MINM real(RP) 0.0_RP minimum limit of absolute velocity for momentum [m/s]
    ATMOS_PHY_SF_CM_MIN real(RP) 1.0E-5_RP minimum limit of bulk coefficient for momentum [NIL]
    ATMOS_PHY_SF_CONST_USTAR real(RP) 0.25_RP constant friction velocity [m/s]
    ATMOS_PHY_SF_CONST_CM real(RP) 0.0011_RP constant bulk coefficient for momentum [NIL]
    ATMOS_PHY_SF_CONST_SH real(RP) 15.0_RP constant surface sensible heat flux [W/m2]
    ATMOS_PHY_SF_CONST_LH real(RP) 115.0_RP constant surface latent heat flux [W/m2]
    ATMOS_PHY_SF_FLG_SH_DIURNAL logical .false. diurnal modulation for sensible heat flux?
    ATMOS_PHY_SF_CONST_FREQ real(RP) 24.0_RP frequency of sensible heat flux modulation [hour]

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_sf_const_setup()

subroutine, public scale_atmos_phy_sf_const::atmos_phy_sf_const_setup ( )

Setup.

Definition at line 65 of file scale_atmos_phy_sf_const.F90.

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

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup().

65  use scale_prc, only: &
66  prc_abort
67  implicit none
68 
69  namelist / param_atmos_phy_sf_const / &
70  atmos_phy_sf_flg_mom_flux, &
71  atmos_phy_sf_u_minm, &
72  atmos_phy_sf_cm_min, &
73  atmos_phy_sf_const_ustar, &
74  atmos_phy_sf_const_cm, &
75  atmos_phy_sf_const_sh, &
76  atmos_phy_sf_const_lh, &
77  atmos_phy_sf_flg_sh_diurnal, &
78  atmos_phy_sf_const_freq
79 
80  integer :: ierr
81  !---------------------------------------------------------------------------
82 
83  log_newline
84  log_info("ATMOS_PHY_SF_const_setup",*) 'Setup'
85  log_info("ATMOS_PHY_SF_const_setup",*) 'Constant flux'
86 
87  !--- read namelist
88  rewind(io_fid_conf)
89  read(io_fid_conf,nml=param_atmos_phy_sf_const,iostat=ierr)
90  if( ierr < 0 ) then !--- missing
91  log_info("ATMOS_PHY_SF_const_setup",*) 'Not found namelist. Default used.'
92  elseif( ierr > 0 ) then !--- fatal error
93  log_error("ATMOS_PHY_SF_const_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_SF_CONST. Check!'
94  call prc_abort
95  endif
96  log_nml(param_atmos_phy_sf_const)
97 
98  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_sf_const_flux()

subroutine, public scale_atmos_phy_sf_const::atmos_phy_sf_const_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_Z1,
real(rp), dimension(ia,ja), intent(in)  SFC_DENS,
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)  U10,
real(rp), dimension (ia,ja), intent(out)  V10 
)

Constant flux.

Definition at line 110 of file scale_atmos_phy_sf_const.F90.

References scale_const::const_pi, and scale_time::time_nowsec.

Referenced by mod_atmos_phy_sf_driver::atmos_phy_sf_driver_calc_tendency().

110  use scale_const, only: &
111  pi => const_pi
112  use scale_atmos_hydrometeor, only: &
113  hydrometeor_lhv => atmos_hydrometeor_lhv
114  use scale_time, only: &
116  implicit none
117  integer, intent(in) :: ia, is, ie
118  integer, intent(in) :: ja, js, je
119 
120  real(RP), intent(in) :: atm_w (ia,ja) ! velocity w at the lowermost layer (cell center) [m/s]
121  real(RP), intent(in) :: atm_u (ia,ja) ! velocity u at the lowermost layer (cell center) [m/s]
122  real(RP), intent(in) :: atm_v (ia,ja) ! velocity v at the lowermost layer (cell center) [m/s]
123  real(RP), intent(in) :: atm_temp(ia,ja) ! temperature at the lowermost layer (cell center) [K]
124  real(RP), intent(in) :: atm_z1 (ia,ja) ! height of the lowermost grid from surface (cell center) [m]
125  real(RP), intent(in) :: sfc_dens(ia,ja) ! density at the surface atmosphere [kg/m3]
126 
127  real(RP), intent(out) :: sflx_mw(ia,ja) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
128  real(RP), intent(out) :: sflx_mu(ia,ja) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
129  real(RP), intent(out) :: sflx_mv(ia,ja) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
130  real(RP), intent(out) :: sflx_sh(ia,ja) ! surface flux for sensible heat (area center) [J/m2/s]
131  real(RP), intent(out) :: sflx_lh(ia,ja) ! surface flux for latent heat (area center) [J/m2/s]
132  real(RP), intent(out) :: sflx_qv(ia,ja) ! surface flux for qv (area center) [kg/m2/s]
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 
136  real(RP) :: atm_uabs(ia,ja) ! absolute velocity at z1 [m/s]
137 
138  real(RP) :: cm(ia,ja) ! bulk coefficient for momentum
139  real(RP) :: r10
140 
141  real(RP) :: modulation
142  real(RP) :: lhv(ia,ja)
143 
144  integer :: i, j
145  !---------------------------------------------------------------------------
146 
147  log_progress(*) 'atmosphere / physics / surface flux / const'
148 
149  !$omp parallel do
150  do j = js, je
151  do i = is, ie
152  atm_uabs(i,j) = min( atmos_phy_sf_u_maxm, max( atmos_phy_sf_u_minm, &
153  sqrt( atm_w(i,j)**2 + atm_u(i,j)**2 + atm_v(i,j)**2 ) ) ) ! at cell center
154  enddo
155  enddo
156 
157  if ( atmos_phy_sf_flg_mom_flux == 0 ) then ! Bulk coefficient is constant
158  !$omp parallel do
159  do j = js, je
160  do i = is, ie
161  cm(i,j) = atmos_phy_sf_const_cm
162  enddo
163  enddo
164  elseif( atmos_phy_sf_flg_mom_flux == 1 ) then ! Friction velocity is constant
165  !$omp parallel do
166  do j = js, je
167  do i = is, ie
168  cm(i,j) = ( atmos_phy_sf_const_ustar / atm_uabs(i,j) )**2
169  cm(i,j) = min( max( cm(i,j), atmos_phy_sf_cm_min ), atmos_phy_sf_cm_max )
170  enddo
171  enddo
172  endif
173 
174  !-----< momentum >-----
175 
176  !$omp parallel do
177  do j = js, je
178  do i = is, ie
179  sflx_mw(i,j) = -cm(i,j) * atm_uabs(i,j) * sfc_dens(i,j) * atm_w(i,j)
180  sflx_mu(i,j) = -cm(i,j) * atm_uabs(i,j) * sfc_dens(i,j) * atm_u(i,j)
181  sflx_mv(i,j) = -cm(i,j) * atm_uabs(i,j) * sfc_dens(i,j) * atm_v(i,j)
182  enddo
183  enddo
184 
185  !-----< heat flux >-----
186 
187  if ( atmos_phy_sf_flg_sh_diurnal ) then
188  modulation = sin( 2.0_rp * pi * time_nowsec / 3600.0_rp / atmos_phy_sf_const_freq )
189  else
190  modulation = 1.0_rp
191  endif
192 
193  !$omp parallel do
194  do j = js, je
195  do i = is, ie
196  sflx_sh(i,j) = atmos_phy_sf_const_sh * modulation
197  sflx_lh(i,j) = atmos_phy_sf_const_lh
198  enddo
199  enddo
200 
201  !-----< mass flux >-----
202  call hydrometeor_lhv( &
203  ia, is, ie, ja, js, je, &
204  atm_temp(:,:), & ! [IN]
205  lhv(:,:) ) ! [OUT]
206 
207  !$omp parallel do
208  do j = js, je
209  do i = is, ie
210  sflx_qv(i,j) = sflx_lh(i,j) / lhv(i,j)
211  enddo
212  enddo
213 
214  !-----< U10, V10 >-----
215 
216  !$omp parallel do
217  do j = js, je
218  do i = is, ie
219  r10 = 10.0_rp / atm_z1(i,j)
220 
221  u10(i,j) = r10 * atm_u(i,j)
222  v10(i,j) = r10 * atm_v(i,j)
223  enddo
224  enddo
225 
226  return
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:72
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module atmosphere / hydrometeor
integer, public je
end point of inner domain: y, local
module TIME
Definition: scale_time.F90:16
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
real(rp), public const_pi
pi
Definition: scale_const.F90:31
Here is the caller graph for this function: