SCALE-RM
scale_land_phy_snow_diagnos.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: land_phy_snow_diags
27 
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public parameters & variables
31  !
32  !-----------------------------------------------------------------------------
33  !
34  !++ Private procedure
35  !
36  !-----------------------------------------------------------------------------
37  !
38  !++ Private parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41 contains
42  !-----------------------------------------------------------------------------
43  subroutine land_phy_snow_diags( &
44  LIA, LIS, LIE, LJA, LJS, LJE, &
45  SNOW_frac, &
46  TMPA, PRSA, &
47  WA, UA, VA, &
48  RHOA, QVA, &
49  Z1, PBL, &
50  RHOS, PRSS, LST1, &
51  QVEF, &
52  Z0M, Z0H, Z0E, &
53  ZMFLX, XMFLX, YMFLX, &
54  U10, V10, &
55  T2, Q2 )
56  use scale_const, only: &
57  pre00 => const_pre00, &
58  rdry => const_rdry, &
59  cpdry => const_cpdry, &
60  rvap => const_rvap, &
61  stb => const_stb
62  use scale_atmos_saturation, only: &
63  qsat => atmos_saturation_pres2qsat_all
64  use scale_bulkflux, only: &
65  bulkflux
66  implicit none
67  integer, intent(in) :: LIA, LIS, LIE
68  integer, intent(in) :: LJA, LJS, LJE
69 
70  real(RP), intent(in) :: SNOW_frac(lia,lja) ! snow fraction
71  real(RP), intent(in) :: TMPA(lia,lja) ! temperature at the lowest atmospheric layer [K]
72  real(RP), intent(in) :: PRSA(lia,lja) ! pressure at the lowest atmospheric layer [Pa]
73  real(RP), intent(in) :: WA (lia,lja) ! velocity w at the lowest atmospheric layer [m/s]
74  real(RP), intent(in) :: UA (lia,lja) ! velocity u at the lowest atmospheric layer [m/s]
75  real(RP), intent(in) :: VA (lia,lja) ! velocity v at the lowest atmospheric layer [m/s]
76  real(RP), intent(in) :: RHOA(lia,lja) ! density at the lowest atmospheric layer [kg/m3]
77  real(RP), intent(in) :: QVA (lia,lja) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
78  real(RP), intent(in) :: Z1 (lia,lja) ! cell center height at the lowest atmospheric layer [m]
79  real(RP), intent(in) :: PBL (lia,lja) ! the top of atmospheric mixing layer [m]
80  real(RP), intent(in) :: RHOS(lia,lja) ! density at the surface [kg/m3]
81  real(RP), intent(in) :: PRSS(lia,lja) ! pressure at the surface [Pa]
82 
83  real(RP), intent(in) :: LST1 (lia,lja) ! land surface temperature [K]
84  real(RP), intent(in) :: QVEF (lia,lja) ! efficiency of evaporation (0-1)
85 
86  real(RP), intent(in) :: Z0M (lia,lja) ! roughness length for momemtum [m]
87  real(RP), intent(in) :: Z0H (lia,lja) ! roughness length for heat [m]
88  real(RP), intent(in) :: Z0E (lia,lja) ! roughness length for vapor [m]
89 
90  real(RP), intent(out) :: ZMFLX(lia,lja) ! z-momentum flux at the surface [kg/m/s2]
91  real(RP), intent(out) :: XMFLX(lia,lja) ! x-momentum flux at the surface [kg/m/s2]
92  real(RP), intent(out) :: YMFLX(lia,lja) ! y-momentum flux at the surface [kg/m/s2]
93  real(RP), intent(out) :: U10 (lia,lja) ! velocity u at 10m [m/s]
94  real(RP), intent(out) :: V10 (lia,lja) ! velocity v at 10m [m/s]
95  real(RP), intent(out) :: T2 (lia,lja) ! temperature at 2m [K]
96  real(RP), intent(out) :: Q2 (lia,lja) ! water vapor at 2m [kg/kg]
97 
98 
99  ! works
100  real(RP) :: Ustar ! friction velocity [m]
101  real(RP) :: Tstar ! friction potential temperature [K]
102  real(RP) :: Qstar ! friction water vapor mass ratio [kg/kg]
103  real(RP) :: Uabs ! modified absolute velocity [m/s]
104  real(RP) :: Ra ! Aerodynamic resistance (=1/Ce) [1/s]
105 
106  real(RP) :: QVsat ! saturation water vapor mixing ratio at surface [kg/kg]
107  real(RP) :: QVS ! water vapor mixing ratio at surface [kg/kg]
108  real(RP) :: Rtot ! total gas constant
109  real(RP) :: qdry ! dry air mass ratio [kg/kg]
110 
111  real(RP) :: FracU10 ! calculation parameter for U10 [-]
112  real(RP) :: FracT2 ! calculation parameter for T2 [-]
113  real(RP) :: FracQ2 ! calculation parameter for Q2 [-]
114 
115  integer :: i, j
116  !---------------------------------------------------------------------------
117 
118  log_info("LAND_PHY_SNOW_DIAGS",*) 'Snow surface diagnostic'
119 
120  ! calculate surface flux
121  do j = ljs, lje
122  do i = lis, lie
123 
124  if( snow_frac(i,j) > 0.0_rp ) then
125 
126  qdry = 1.0_rp - qva(i,j)
127  rtot = qdry * rdry + qva(i,j) * rvap
128 
129  call qsat( lst1(i,j), prss(i,j), qdry, & ! [IN]
130  qvsat ) ! [OUT]
131 
132  qvs = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * qvsat
133 
134  call bulkflux( &
135  ustar, & ! [OUT]
136  tstar, & ! [OUT]
137  qstar, & ! [OUT]
138  uabs, & ! [OUT]
139  ra, & ! [OUT]
140  fracu10, & ! [OUT]
141  fract2, & ! [OUT]
142  fracq2, & ! [OUT]
143  tmpa(i,j), & ! [IN]
144  lst1(i,j), & ! [IN]
145  prsa(i,j), & ! [IN]
146  prss(i,j), & ! [IN]
147  qva(i,j), & ! [IN]
148  qvs, & ! [IN]
149  ua(i,j), & ! [IN]
150  va(i,j), & ! [IN]
151  z1(i,j), & ! [IN]
152  pbl(i,j), & ! [IN]
153  z0m(i,j), & ! [IN]
154  z0h(i,j), & ! [IN]
155  z0e(i,j) ) ! [IN]
156 
157  zmflx(i,j) = -rhos(i,j) * ustar * ustar / uabs * wa(i,j)
158  xmflx(i,j) = -rhos(i,j) * ustar * ustar / uabs * ua(i,j)
159  ymflx(i,j) = -rhos(i,j) * ustar * ustar / uabs * va(i,j)
160 
161  ! diagnostic variables for neutral state
162  u10(i,j) = ua(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
163  v10(i,j) = va(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
164  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) ) ) &
165  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0h(i,j) ) )
166  q2(i,j) = qvs + ( qva(i,j) - qvs ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0e(i,j) ) ) &
167  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0e(i,j) ) )
168 
169  else
170 
171  ! not calculate surface flux
172  zmflx(i,j) = 0.0_rp
173  xmflx(i,j) = 0.0_rp
174  ymflx(i,j) = 0.0_rp
175  u10(i,j) = 0.0_rp
176  v10(i,j) = 0.0_rp
177  t2(i,j) = 0.0_rp
178  q2(i,j) = 0.0_rp
179 
180  end if
181 
182  end do
183  end do
184 
185  return
186  end subroutine land_phy_snow_diags
187 
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
module atmosphere / saturation
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:49
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
subroutine, public land_phy_snow_diags(LIA, LIS, LIE, LJA, LJS, LJE, SNOW_frac, TMPA, PRSA, WA, UA, VA, RHOA, QVA, Z1, PBL, RHOS, PRSS, LST1, QVEF, Z0M, Z0H, Z0E, ZMFLX, XMFLX, YMFLX, U10, V10, T2, Q2)
module land / physics / snow / diagnostics
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
procedure(bc), pointer, public bulkflux
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:63
module CONSTANT
Definition: scale_const.F90:11
module profiler
Definition: scale_prof.F90:11
module Surface bulk flux
module PRECISION
module STDIO
Definition: scale_io.F90:10