SCALE-RM
Functions/Subroutines
scale_land_sfc_const Module Reference

module LAND / Surface fluxes with constant land model More...

Functions/Subroutines

subroutine, public land_sfc_const_setup (LAND_TYPE)
 Setup. More...
 
subroutine, public land_sfc_const (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_DP)
 

Detailed Description

module LAND / Surface fluxes with constant land model

Description
Surface flux with constant land model
Author
Team SCALE

Function/Subroutine Documentation

◆ land_sfc_const_setup()

subroutine, public scale_land_sfc_const::land_sfc_const_setup ( character(len=*), intent(in)  LAND_TYPE)

Setup.

Definition at line 46 of file scale_land_sfc_const.F90.

References scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by scale_land_sfc::land_sfc_setup().

46  implicit none
47 
48  character(len=*), intent(in) :: LAND_TYPE
49  !---------------------------------------------------------------------------
50 
51  if( io_l ) write(io_fid_log,*)
52  if( io_l ) write(io_fid_log,*) '++++++ Module[CONST] / Categ[LAND SFC] / Origin[SCALElib]'
53 
54  return
Here is the caller graph for this function:

◆ land_sfc_const()

subroutine, public scale_land_sfc_const::land_sfc_const ( 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_DP 
)

Definition at line 93 of file scale_land_sfc_const.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().

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