SCALE-RM
Functions/Subroutines
scale_atmos_phy_sf_const Module Reference

module ATMOSPHERE / Physics Surface fluxes More...

Functions/Subroutines

subroutine, public atmos_phy_sf_const_setup (ATMOS_PHY_SF_TYPE)
 Setup. More...
 
subroutine, public atmos_phy_sf_const (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)
 Constant flux. More...
 

Detailed Description

module ATMOSPHERE / Physics Surface fluxes

Description
Flux from/to bottom wall of atmosphere (surface) Constant flux
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
  • 2012-09-12 (Y.Sato) [renew] constant FLUX version
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 ( character(len=*), intent(in)  ATMOS_PHY_SF_TYPE)

Setup.

Definition at line 74 of file scale_atmos_phy_sf_const.F90.

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

Referenced by scale_atmos_phy_sf::atmos_phy_sf_setup().

74  use scale_process, only: &
76  implicit none
77 
78  character(len=*), intent(in) :: atmos_phy_sf_type
79 
80  namelist / param_atmos_phy_sf_const / &
81  atmos_phy_sf_flg_mom_flux, &
82  atmos_phy_sf_u_minm, &
83  atmos_phy_sf_cm_min, &
84  atmos_phy_sf_const_ustar, &
85  atmos_phy_sf_const_cm, &
86  atmos_phy_sf_const_sh, &
87  atmos_phy_sf_const_lh, &
88  atmos_phy_sf_flg_sh_diurnal, &
89  atmos_phy_sf_const_freq
90 
91  integer :: ierr
92  !---------------------------------------------------------------------------
93 
94  if( io_l ) write(io_fid_log,*)
95  if( io_l ) write(io_fid_log,*) '++++++ Module[SURFACE FLUX] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
96  if( io_l ) write(io_fid_log,*) '+++ Constant flux'
97 
98  if ( atmos_phy_sf_type /= 'CONST' ) then
99  write(*,*) 'xxx ATMOS_PHY_SF_TYPE is not CONST. Check!'
100  call prc_mpistop
101  endif
102 
103  !--- read namelist
104  rewind(io_fid_conf)
105  read(io_fid_conf,nml=param_atmos_phy_sf_const,iostat=ierr)
106  if( ierr < 0 ) then !--- missing
107  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
108  elseif( ierr > 0 ) then !--- fatal error
109  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_SF_CONST. Check!'
110  call prc_mpistop
111  endif
112  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_sf_const)
113 
114  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_const()

subroutine, public scale_atmos_phy_sf_const::atmos_phy_sf_const ( 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 
)

Constant flux.

Definition at line 131 of file scale_atmos_phy_sf_const.F90.

References scale_const::const_pi, scale_tracer::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_time::time_nowsec.

Referenced by scale_atmos_phy_sf::atmos_phy_sf_setup().

131  use scale_grid_index
132  use scale_tracer
133  use scale_const, only: &
134  pi => const_pi
135  use scale_atmos_thermodyn, only: &
136  atmos_thermodyn_templhv
137  use scale_time, only: &
139  implicit none
140 
141  real(RP), intent(in) :: atm_temp (ia,ja) ! temperature at the lowermost layer (cell center) [K]
142  real(RP), intent(in) :: atm_pres (ia,ja) ! pressure at the lowermost layer (cell center) [Pa]
143  real(RP), intent(in) :: atm_w (ia,ja) ! velocity w at the lowermost layer (cell center) [m/s]
144  real(RP), intent(in) :: atm_u (ia,ja) ! velocity u at the lowermost layer (cell center) [m/s]
145  real(RP), intent(in) :: atm_v (ia,ja) ! velocity v at the lowermost layer (cell center) [m/s]
146  real(RP), intent(in) :: atm_dens (ia,ja) ! density at the lowermost layer (cell center) [kg/m3]
147  real(RP), intent(in) :: atm_qtrc (ia,ja,qa) ! tracer at the lowermost layer (cell center) [kg/kg]
148  real(RP), intent(in) :: atm_z1 (ia,ja) ! height of the lowermost grid from surface (cell center) [m]
149  real(DP), intent(in) :: dt ! delta time
150  real(RP), intent(in) :: sfc_dens (ia,ja) ! density at the surface atmosphere [kg/m3]
151  real(RP), intent(in) :: sfc_pres (ia,ja) ! pressure at the surface atmosphere [Pa]
152  real(RP), intent(in) :: sflx_lw_dn(ia,ja) ! downward longwave radiation flux at the surface [J/m2/s]
153  real(RP), intent(in) :: sflx_sw_dn(ia,ja) ! downward shortwave radiation flux at the surface [J/m2/s]
154  real(RP), intent(in) :: sfc_temp (ia,ja) ! temperature at the surface skin [K]
155  real(RP), intent(in) :: sfc_albedo(ia,ja,2) ! surface albedo (LW/SW) [0-1]
156  real(RP), intent(inout) :: sfc_z0m (ia,ja) ! surface roughness length (momentum) [m]
157  real(RP), intent(inout) :: sfc_z0h (ia,ja) ! surface roughness length (heat) [m]
158  real(RP), intent(inout) :: sfc_z0e (ia,ja) ! surface roughness length (vapor) [m]
159  real(RP), intent(out) :: sflx_mw (ia,ja) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
160  real(RP), intent(out) :: sflx_mu (ia,ja) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
161  real(RP), intent(out) :: sflx_mv (ia,ja) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
162  real(RP), intent(out) :: sflx_sh (ia,ja) ! surface flux for sensible heat (area center) [J/m2/s]
163  real(RP), intent(out) :: sflx_lh (ia,ja) ! surface flux for latent heat (area center) [J/m2/s]
164  real(RP), intent(out) :: sflx_qtrc (ia,ja,qa) ! surface flux for tracer mass (area center) [kg/m2/s]
165  real(RP), intent(out) :: u10 (ia,ja) ! velocity u at 10m height
166  real(RP), intent(out) :: v10 (ia,ja) ! velocity v at 10m height
167  real(RP), intent(out) :: t2 (ia,ja) ! temperature t at 2m height
168  real(RP), intent(out) :: q2 (ia,ja) ! water vapor q at 2m height
169 
170  real(RP) :: atm_uabs(ia,ja) ! absolute velocity at z1 [m/s]
171 
172  real(RP) :: cm(ia,ja) ! bulk coefficient for momentum
173  real(RP) :: r10
174 
175  real(RP) :: modulation
176  real(RP) :: uabs_lim
177  real(RP) :: lhv(ia,ja)
178  integer :: i, j
179  !---------------------------------------------------------------------------
180 
181  if( io_l ) write(io_fid_log,*) '*** Physics step: Surface flux(const)'
182 
183  do j = js, je
184  do i = is, ie
185  atm_uabs(i,j) = sqrt( atm_w(i,j)*atm_w(i,j) &
186  + atm_u(i,j)*atm_u(i,j) &
187  + atm_v(i,j)*atm_v(i,j) ) ! at cell center
188  enddo
189  enddo
190 
191  if ( atmos_phy_sf_flg_mom_flux == 0 ) then ! Bulk coefficient is constant
192  do j = js, je
193  do i = is, ie
194  cm(i,j) = atmos_phy_sf_const_cm
195  enddo
196  enddo
197  elseif( atmos_phy_sf_flg_mom_flux == 1 ) then ! Friction velocity is constant
198  do j = js, je
199  do i = is, ie
200  cm(i,j) = ( atmos_phy_sf_const_ustar / atm_uabs(i,j) )**2
201  cm(i,j) = min( max( cm(i,j), atmos_phy_sf_cm_min ), atmos_phy_sf_cm_max )
202  enddo
203  enddo
204  endif
205 
206  !-----< momentum >-----
207 
208  do j = js, je
209  do i = is, ie
210  uabs_lim = min( max( atm_uabs(i,j), atmos_phy_sf_u_minm ), atmos_phy_sf_u_maxm )
211 
212  sflx_mw(i,j) = -cm(i,j) * uabs_lim * sfc_dens(i,j) * atm_w(i,j)
213  sflx_mu(i,j) = -cm(i,j) * uabs_lim * sfc_dens(i,j) * atm_u(i,j)
214  sflx_mv(i,j) = -cm(i,j) * uabs_lim * sfc_dens(i,j) * atm_v(i,j)
215  enddo
216  enddo
217 
218  !-----< heat flux >-----
219 
220  if ( atmos_phy_sf_flg_sh_diurnal ) then
221  modulation = sin( 2.0_rp * pi * time_nowsec / 3600.0_rp / atmos_phy_sf_const_freq )
222  else
223  modulation = 1.0_rp
224  endif
225 
226  do j = js, je
227  do i = is, ie
228  sflx_sh(i,j) = atmos_phy_sf_const_sh * modulation
229  sflx_lh(i,j) = atmos_phy_sf_const_lh
230  enddo
231  enddo
232 
233  !-----< mass flux >-----
234  call atmos_thermodyn_templhv( lhv, atm_temp )
235 
236  sflx_qtrc(:,:,:) = 0.0_rp
237  do j = js, je
238  do i = is, ie
239  sflx_qtrc(i,j,i_qv) = sflx_lh(i,j) / lhv(i,j)
240  enddo
241  enddo
242 
243  !-----< U10, T2, q2 >-----
244 
245  do j = js, je
246  do i = is, ie
247  r10 = 10.0_rp / atm_z1(i,j)
248 
249  u10(i,j) = r10 * atm_u(i,j)
250  v10(i,j) = r10 * atm_v(i,j)
251  enddo
252  enddo
253 
254  do j = js, je
255  do i = is, ie
256  t2(i,j) = atm_temp(i,j)
257  q2(i,j) = atm_qtrc(i,j,i_qv)
258  enddo
259  enddo
260 
261  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:68
integer, public qa
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public i_qv
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module CONSTANT
Definition: scale_const.F90:14
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
real(rp), public const_pi
pi
Definition: scale_const.F90:34
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function: