SCALE-RM
scale_land_sfc_const.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: land_sfc_const_setup
27  public :: land_sfc_const
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public parameters & variables
32  !
33  !-----------------------------------------------------------------------------
34  !
35  !++ Private procedure
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private parameters & variables
40  !
41  !-----------------------------------------------------------------------------
42 contains
43  !-----------------------------------------------------------------------------
45  subroutine land_sfc_const_setup( LAND_TYPE )
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
55  end subroutine land_sfc_const_setup
56 
57  !-----------------------------------------------------------------------------
58  subroutine land_sfc_const( &
59  LST_t, &
60  ZMFLX, &
61  XMFLX, &
62  YMFLX, &
63  SHFLX, &
64  LHFLX, &
65  GHFLX, &
66  U10, &
67  V10, &
68  T2, &
69  Q2, &
70  TMPA, &
71  PRSA, &
72  WA, &
73  UA, &
74  VA, &
75  RHOA, &
76  QVA, &
77  Z1, &
78  PBL, &
79  PRSS, &
80  LWD, &
81  SWD, &
82  TG, &
83  LST, &
84  QVEF, &
85  ALB_LW, &
86  ALB_SW, &
87  DZG, &
88  TCS, &
89  Z0M, &
90  Z0H, &
91  Z0E, &
92  dt_DP )
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
266  end subroutine land_sfc_const
267 
268 end module scale_land_sfc_const
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
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.
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:51
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module STDIO
Definition: scale_stdio.F90:12
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
integer, public ia
of whole cells: x, local, with HALO
module LANDUSE
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
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)
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module PRECISION
subroutine, public land_sfc_const_setup(LAND_TYPE)
Setup.
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
module LAND / Surface fluxes with constant land model
integer, public ja
of whole cells: y, local, with HALO