SCALE-RM
Functions/Subroutines
scale_land_sfc_thick_slab Module Reference

module LAND / Surface fluxes with thick-slab land model More...

Functions/Subroutines

subroutine, public land_sfc_thick_slab_setup (LAND_TYPE)
 Setup. More...
 
subroutine, public land_sfc_thick_slab (LST_t, ZMFLX, XMFLX, YMFLX, SHFLX, LHFLX, GHFLX, U10, V10, T2, Q2, TMPA, PRSA, WA, UA, VA, RHOA, QVA, Z1, PBL, PRSS, LWD, SWD, TG, LST, QVEF, ALB_LW, ALB_SW, DZG, TCS, Z0M, Z0H, Z0E, dt)
 

Detailed Description

module LAND / Surface fluxes with thick-slab land model

Description
Surface flux with thick-slab land model
Author
Team SCALE
NAMELIST
  • PARAM_LAND_SFC_THICK_SLAB
    nametypedefault valuecomment
    DUMMY logical

History Output
No history output

Function/Subroutine Documentation

◆ land_sfc_thick_slab_setup()

subroutine, public scale_land_sfc_thick_slab::land_sfc_thick_slab_setup ( character(len=*), intent(in)  LAND_TYPE)

Setup.

Definition at line 46 of file scale_land_sfc_thick-slab.F90.

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

Referenced by scale_land_sfc::land_sfc_setup().

46  use scale_process, only: &
48  implicit none
49 
50  character(len=*), intent(in) :: LAND_TYPE
51 
52  logical :: dummy
53 
54  namelist / param_land_sfc_thick_slab / &
55  dummy
56 
57  integer :: ierr
58  !---------------------------------------------------------------------------
59 
60  if( io_l ) write(io_fid_log,*)
61  if( io_l ) write(io_fid_log,*) '++++++ Module[THICK-SLAB] / Categ[LAND SFC] / Origin[SCALElib]'
62 
63  !--- read namelist
64  rewind(io_fid_conf)
65  read(io_fid_conf,nml=param_land_sfc_thick_slab,iostat=ierr)
66  if( ierr < 0 ) then !--- missing
67  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
68  elseif( ierr > 0 ) then !--- fatal error
69  write(*,*) 'xxx Not appropriate names in namelist PARAM_LAND_SFC_THICK_SLAB. Check!'
70  call prc_mpistop
71  endif
72  if( io_nml ) write(io_fid_nml,nml=param_land_sfc_thick_slab)
73 
74  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_sfc_thick_slab()

subroutine, public scale_land_sfc_thick_slab::land_sfc_thick_slab ( real(rp), dimension(ia,ja), intent(out)  LST_t,
real(rp), dimension(ia,ja), intent(out)  ZMFLX,
real(rp), dimension(ia,ja), intent(out)  XMFLX,
real(rp), dimension(ia,ja), intent(out)  YMFLX,
real(rp), dimension(ia,ja), intent(out)  SHFLX,
real(rp), dimension(ia,ja), intent(out)  LHFLX,
real(rp), dimension(ia,ja), intent(out)  GHFLX,
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,
real(rp), dimension(ia,ja), intent(in)  TMPA,
real(rp), dimension(ia,ja), intent(in)  PRSA,
real(rp), dimension (ia,ja), intent(in)  WA,
real(rp), dimension (ia,ja), intent(in)  UA,
real(rp), dimension (ia,ja), intent(in)  VA,
real(rp), dimension(ia,ja), intent(in)  RHOA,
real(rp), dimension (ia,ja), intent(in)  QVA,
real(rp), dimension (ia,ja), intent(in)  Z1,
real(rp), dimension (ia,ja), intent(in)  PBL,
real(rp), dimension(ia,ja), intent(in)  PRSS,
real(rp), dimension (ia,ja), intent(in)  LWD,
real(rp), dimension (ia,ja), intent(in)  SWD,
real(rp), dimension (ia,ja), intent(in)  TG,
real(rp), dimension (ia,ja), intent(in)  LST,
real(rp), dimension (ia,ja), intent(in)  QVEF,
real(rp), dimension(ia,ja), intent(in)  ALB_LW,
real(rp), dimension(ia,ja), intent(in)  ALB_SW,
real(rp), dimension (ia,ja), intent(in)  DZG,
real(rp), dimension (ia,ja), intent(in)  TCS,
real(rp), dimension (ia,ja), intent(in)  Z0M,
real(rp), dimension (ia,ja), intent(in)  Z0H,
real(rp), dimension (ia,ja), intent(in)  Z0E,
real(dp), intent(in)  dt 
)

Definition at line 113 of file scale_land_sfc_thick-slab.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_stb, 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, scale_landuse::landuse_fact_land, and scale_process::prc_mpistop().

Referenced by scale_land_sfc::land_sfc_setup().

113  use scale_process, only: &
115  use scale_const, only: &
116  pre00 => const_pre00, &
117  rdry => const_rdry, &
118  cpdry => const_cpdry, &
119  stb => const_stb
120  use scale_landuse, only: &
122  use scale_atmos_hydrometeor, only: &
123  hydrometeor_lhv => atmos_hydrometeor_lhv
124  use scale_atmos_saturation, only: &
125  qsat => atmos_saturation_pres2qsat_all
126  use scale_bulkflux, only: &
127  bulkflux
128  implicit none
129 
130  ! arguments
131  real(RP), intent(out) :: LST_t(IA,JA) ! tendency of LST
132  real(RP), intent(out) :: ZMFLX(IA,JA) ! z-momentum flux at the surface [kg/m/s2]
133  real(RP), intent(out) :: XMFLX(IA,JA) ! x-momentum flux at the surface [kg/m/s2]
134  real(RP), intent(out) :: YMFLX(IA,JA) ! y-momentum flux at the surface [kg/m/s2]
135  real(RP), intent(out) :: SHFLX(IA,JA) ! sensible heat flux at the surface [J/m2/s]
136  real(RP), intent(out) :: LHFLX(IA,JA) ! latent heat flux at the surface [J/m2/s]
137  real(RP), intent(out) :: GHFLX(IA,JA) ! ground heat flux at the surface [J/m2/s]
138  real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m [m/s]
139  real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m [m/s]
140  real(RP), intent(out) :: T2 (IA,JA) ! temperature at 2m [K]
141  real(RP), intent(out) :: Q2 (IA,JA) ! water vapor at 2m [kg/kg]
142 
143  real(RP), intent(in) :: TMPA(IA,JA) ! temperature at the lowest atmospheric layer [K]
144  real(RP), intent(in) :: PRSA(IA,JA) ! pressure at the lowest atmospheric layer [Pa]
145  real(RP), intent(in) :: WA (IA,JA) ! velocity w at the lowest atmospheric layer [m/s]
146  real(RP), intent(in) :: UA (IA,JA) ! velocity u at the lowest atmospheric layer [m/s]
147  real(RP), intent(in) :: VA (IA,JA) ! velocity v at the lowest atmospheric layer [m/s]
148  real(RP), intent(in) :: RHOA(IA,JA) ! density at the lowest atmospheric layer [kg/m3]
149  real(RP), intent(in) :: QVA (IA,JA) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
150  real(RP), intent(in) :: Z1 (IA,JA) ! cell center height at the lowest atmospheric layer [m]
151  real(RP), intent(in) :: PBL (IA,JA) ! the top of atmospheric mixing layer [m]
152  real(RP), intent(in) :: PRSS(IA,JA) ! pressure at the surface [Pa]
153  real(RP), intent(in) :: LWD (IA,JA) ! downward long-wave radiation flux at the surface [J/m2/s]
154  real(RP), intent(in) :: SWD (IA,JA) ! downward short-wave radiation flux at the surface [J/m2/s]
155 
156  real(RP), intent(in) :: TG (IA,JA) ! soil temperature [K]
157  real(RP), intent(in) :: LST (IA,JA) ! land surface temperature [K]
158  real(RP), intent(in) :: QVEF (IA,JA) ! efficiency of evaporation (0-1)
159  real(RP), intent(in) :: ALB_LW(IA,JA) ! surface albedo for LW (0-1)
160  real(RP), intent(in) :: ALB_SW(IA,JA) ! surface albedo for SW (0-1)
161  real(RP), intent(in) :: DZG (IA,JA) ! soil depth [m]
162  real(RP), intent(in) :: TCS (IA,JA) ! thermal conductivity for soil [J/m/K/s]
163  real(RP), intent(in) :: Z0M (IA,JA) ! roughness length for momemtum [m]
164  real(RP), intent(in) :: Z0H (IA,JA) ! roughness length for heat [m]
165  real(RP), intent(in) :: Z0E (IA,JA) ! roughness length for vapor [m]
166  real(DP), intent(in) :: dt ! delta time
167 
168  ! works
169  real(RP) :: LST1(IA,JA)
170 
171  real(RP) :: res ! residual
172 
173  real(RP) :: Ustar ! friction velocity [m]
174  real(RP) :: Tstar ! friction potential temperature [K]
175  real(RP) :: Qstar ! friction water vapor mass ratio [kg/kg]
176  real(RP) :: Uabs ! modified absolute velocity [m/s]
177  real(RP) :: QVsat ! saturation water vapor mixing ratio at surface [kg/kg]
178  real(RP) :: QVS ! water vapor mixing ratio at surface [kg/kg]
179 
180  real(RP) :: FracU10 ! calculation parameter for U10 [-]
181  real(RP) :: FracT2 ! calculation parameter for T2 [-]
182  real(RP) :: FracQ2 ! calculation parameter for Q2 [-]
183 
184  real(RP) :: LHV(IA,JA) ! latent heat of vaporization [J/kg]
185 
186  integer :: i, j
187  !---------------------------------------------------------------------------
188 
189  if( io_l ) write(io_fid_log,*) '*** Land surface step: Thick-Slab'
190 
191  call hydrometeor_lhv( lhv(:,:), tmpa(:,:) )
192 
193  ! update surface temperature
194  do j = js, je
195  do i = is, ie
196  lst1(i,j) = tg(i,j)
197  lst_t(i,j) = ( lst1(i,j) - lst(i,j) ) / dt
198  end do
199  end do
200 
201  ! calculate surface flux
202  do j = js, je
203  do i = is, ie
204 
205  if( landuse_fact_land(i,j) > 0.0_rp ) then
206 
207  call qsat( qvsat, & ! [OUT]
208  lst1(i,j), & ! [IN]
209  prss(i,j) ) ! [IN]
210 
211  qvs = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * qvsat
212 
213  call bulkflux( &
214  ustar, & ! [OUT]
215  tstar, & ! [OUT]
216  qstar, & ! [OUT]
217  uabs, & ! [OUT]
218  fracu10, & ! [OUT]
219  fract2, & ! [OUT]
220  fracq2, & ! [OUT]
221  tmpa(i,j), & ! [IN]
222  lst1(i,j), & ! [IN]
223  prsa(i,j), & ! [IN]
224  prss(i,j), & ! [IN]
225  qva(i,j), & ! [IN]
226  qvs, & ! [IN]
227  ua(i,j), & ! [IN]
228  va(i,j), & ! [IN]
229  z1(i,j), & ! [IN]
230  pbl(i,j), & ! [IN]
231  z0m(i,j), & ! [IN]
232  z0h(i,j), & ! [IN]
233  z0e(i,j) ) ! [IN]
234 
235  zmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * wa(i,j)
236  xmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * ua(i,j)
237  ymflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * va(i,j)
238 
239  shflx(i,j) = - rhoa(i,j) * ustar * tstar &
240  * cpdry * ( prss(i,j) / pre00 )**( rdry/cpdry )
241  lhflx(i,j) = - rhoa(i,j) * ustar * qstar * lhv(i,j)
242 
243  ghflx(i,j) = -2.0_rp * tcs(i,j) * ( lst1(i,j) - tg(i,j) ) / dzg(i,j)
244 
245  ! calculation for residual
246  res = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) &
247  + ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * lst1(i,j)**4 ) &
248  - shflx(i,j) - lhflx(i,j) + ghflx(i,j)
249 
250  ! put residual in ground heat flux
251  ghflx(i,j) = ghflx(i,j) - res
252 
253  ! diagnostic variables considering unstable/stable state
254  !U10(i,j) = FracU10 * UA(i,j)
255  !V10(i,j) = FracU10 * VA(i,j)
256  !T2 (i,j) = ( 1.0_RP - FracT2 ) * LST1(i,j) + FracT2 * TMPA(i,j)
257  !Q2 (i,j) = ( 1.0_RP - FracQ2 ) * QVS + FracQ2 * QVA (i,j)
258 
259  ! diagnostic variables for neutral state
260  u10(i,j) = ua(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
261  v10(i,j) = va(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
262  t2(i,j) = lst1(i,j) + ( tmpa(i,j) - lst1(i,j) ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0h(i,j) ) ) &
263  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0h(i,j) ) )
264  q2(i,j) = qvs + ( qva(i,j) - qvs ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0e(i,j) ) ) &
265  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0e(i,j) ) )
266 
267  else
268 
269  ! not calculate surface flux
270  zmflx(i,j) = 0.0_rp
271  xmflx(i,j) = 0.0_rp
272  ymflx(i,j) = 0.0_rp
273  shflx(i,j) = 0.0_rp
274  lhflx(i,j) = 0.0_rp
275  ghflx(i,j) = 0.0_rp
276  u10(i,j) = 0.0_rp
277  v10(i,j) = 0.0_rp
278  t2(i,j) = 0.0_rp
279  q2(i,j) = 0.0_rp
280 
281  end if
282 
283  end do
284  end do
285 
286  return
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
module LANDUSE
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module Surface bulk flux
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
Here is the call graph for this function:
Here is the caller graph for this function: