SCALE-RM
Functions/Subroutines
scale_land_phy_snow_diagnos Module Reference

module land / physics / snow / diagnostics More...

Functions/Subroutines

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)
 

Detailed Description

module land / physics / snow / diagnostics

Description
Momentum fluxes and Diagnostic variables
Author
Team SCALE

Function/Subroutine Documentation

◆ land_phy_snow_diags()

subroutine, public scale_land_phy_snow_diagnos::land_phy_snow_diags ( integer, intent(in)  LIA,
integer, intent(in)  LIS,
integer, intent(in)  LIE,
integer, intent(in)  LJA,
integer, intent(in)  LJS,
integer, intent(in)  LJE,
real(rp), dimension(lia,lja), intent(in)  SNOW_frac,
real(rp), dimension(lia,lja), intent(in)  TMPA,
real(rp), dimension(lia,lja), intent(in)  PRSA,
real(rp), dimension (lia,lja), intent(in)  WA,
real(rp), dimension (lia,lja), intent(in)  UA,
real(rp), dimension (lia,lja), intent(in)  VA,
real(rp), dimension(lia,lja), intent(in)  RHOA,
real(rp), dimension (lia,lja), intent(in)  QVA,
real(rp), dimension (lia,lja), intent(in)  Z1,
real(rp), dimension (lia,lja), intent(in)  PBL,
real(rp), dimension(lia,lja), intent(in)  RHOS,
real(rp), dimension(lia,lja), intent(in)  PRSS,
real(rp), dimension (lia,lja), intent(in)  LST1,
real(rp), dimension (lia,lja), intent(in)  QVEF,
real(rp), dimension (lia,lja), intent(in)  Z0M,
real(rp), dimension (lia,lja), intent(in)  Z0H,
real(rp), dimension (lia,lja), intent(in)  Z0E,
real(rp), dimension(lia,lja), intent(out)  ZMFLX,
real(rp), dimension(lia,lja), intent(out)  XMFLX,
real(rp), dimension(lia,lja), intent(out)  YMFLX,
real(rp), dimension (lia,lja), intent(out)  U10,
real(rp), dimension (lia,lja), intent(out)  V10,
real(rp), dimension (lia,lja), intent(out)  T2,
real(rp), dimension (lia,lja), intent(out)  Q2 
)

Definition at line 56 of file scale_land_phy_snow_diagnos.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_rvap, and scale_const::const_stb.

Referenced by mod_land_driver::land_driver_calc_tendency().

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
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
module atmosphere / saturation
integer, public va
Definition: scale_index.F90:35
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
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 Surface bulk flux
Here is the caller graph for this function: