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, Ustar, Tstar, Qstar, Wstar, RLmo, 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)  Ustar,
real(rp), dimension(lia,lja), intent(out)  Tstar,
real(rp), dimension(lia,lja), intent(out)  Qstar,
real(rp), dimension(lia,lja), intent(out)  Wstar,
real(rp), dimension (lia,lja), intent(out)  RLmo,
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 59 of file scale_land_phy_snow_diagnos.F90.

59  use scale_const, only: &
60  eps => const_eps, &
61  undef => const_undef, &
62  pre00 => const_pre00, &
63  rdry => const_rdry, &
64  cpdry => const_cpdry, &
65  rvap => const_rvap, &
66  stb => const_stb
67  use scale_atmos_saturation, only: &
68  qsat => atmos_saturation_dens2qsat_all
69 ! qsat => ATMOS_SATURATION_pres2qsat_all
70  use scale_bulkflux, only: &
71  bulkflux, &
72  bulkflux_diagnose_surface
73  implicit none
74  integer, intent(in) :: LIA, LIS, LIE
75  integer, intent(in) :: LJA, LJS, LJE
76 
77  real(RP), intent(in) :: SNOW_frac(LIA,LJA) ! snow fraction
78  real(RP), intent(in) :: TMPA(LIA,LJA) ! temperature at the lowest atmospheric layer [K]
79  real(RP), intent(in) :: PRSA(LIA,LJA) ! pressure at the lowest atmospheric layer [Pa]
80  real(RP), intent(in) :: WA (LIA,LJA) ! velocity w at the lowest atmospheric layer [m/s]
81  real(RP), intent(in) :: UA (LIA,LJA) ! velocity u at the lowest atmospheric layer [m/s]
82  real(RP), intent(in) :: VA (LIA,LJA) ! velocity v at the lowest atmospheric layer [m/s]
83  real(RP), intent(in) :: RHOA(LIA,LJA) ! density at the lowest atmospheric layer [kg/m3]
84  real(RP), intent(in) :: QVA (LIA,LJA) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
85  real(RP), intent(in) :: Z1 (LIA,LJA) ! cell center height at the lowest atmospheric layer [m]
86  real(RP), intent(in) :: PBL (LIA,LJA) ! the top of atmospheric mixing layer [m]
87  real(RP), intent(in) :: RHOS(LIA,LJA) ! density at the surface [kg/m3]
88  real(RP), intent(in) :: PRSS(LIA,LJA) ! pressure at the surface [Pa]
89 
90  real(RP), intent(in) :: LST1 (LIA,LJA) ! land surface temperature [K]
91  real(RP), intent(in) :: QVEF (LIA,LJA) ! efficiency of evaporation (0-1)
92 
93  real(RP), intent(in) :: Z0M (LIA,LJA) ! roughness length for momemtum [m]
94  real(RP), intent(in) :: Z0H (LIA,LJA) ! roughness length for heat [m]
95  real(RP), intent(in) :: Z0E (LIA,LJA) ! roughness length for vapor [m]
96 
97  real(RP), intent(out) :: ZMFLX(LIA,LJA) ! z-momentum flux at the surface [kg/m/s2]
98  real(RP), intent(out) :: XMFLX(LIA,LJA) ! x-momentum flux at the surface [kg/m/s2]
99  real(RP), intent(out) :: YMFLX(LIA,LJA) ! y-momentum flux at the surface [kg/m/s2]
100  real(RP), intent(out) :: Ustar(LIA,LJA) ! friction velocity [m/s]
101  real(RP), intent(out) :: Tstar(LIA,LJA) ! temperature scale [K]
102  real(RP), intent(out) :: Qstar(LIA,LJA) ! moisture scale [kg/kg]
103  real(RP), intent(out) :: Wstar(LIA,LJA) ! convective velocity scale [m/s]
104  real(RP), intent(out) :: RLmo (LIA,LJA) ! inversed Obukhov length [1/m]
105  real(RP), intent(out) :: U10 (LIA,LJA) ! velocity u at 10m [m/s]
106  real(RP), intent(out) :: V10 (LIA,LJA) ! velocity v at 10m [m/s]
107  real(RP), intent(out) :: T2 (LIA,LJA) ! temperature at 2m [K]
108  real(RP), intent(out) :: Q2 (LIA,LJA) ! water vapor at 2m [kg/kg]
109 
110 
111  ! works
112  real(RP) :: Uabs ! modified absolute velocity [m/s]
113  real(RP) :: Ra ! Aerodynamic resistance (=1/Ce) [1/s]
114 
115  real(RP) :: QVsat ! saturation water vapor mixing ratio at surface [kg/kg]
116  real(RP) :: QVS(LIA,LJA) ! water vapor mixing ratio at surface [kg/kg]
117 ! real(RP) :: Rtot ! total gas constant
118 ! real(RP) :: qdry ! dry air mass ratio [kg/kg]
119 
120  real(RP) :: FracU10(LIA,LJA) ! calculation parameter for U10 [-]
121  real(RP) :: FracT2 (LIA,LJA) ! calculation parameter for T2 [-]
122  real(RP) :: FracQ2 (LIA,LJA) ! calculation parameter for Q2 [-]
123 
124  real(RP) :: MFLUX
125 
126  logical :: calc_flag(LIA,LJA)
127 
128  integer :: i, j
129  !---------------------------------------------------------------------------
130 
131  log_info("LAND_PHY_SNOW_DIAGS",*) 'Snow surface diagnostic'
132 
133  ! calculate surface flux
134  !$omp parallel do &
135  !$omp private(QVsat,Uabs,MFLUX)
136  do j = ljs, lje
137  do i = lis, lie
138 
139  if( snow_frac(i,j) > 0.0_rp ) then
140 
141  calc_flag(i,j) = .true.
142 
143 ! qdry = 1.0_RP - QVA(i,j)
144 ! Rtot = qdry * Rdry + QVA(i,j) * Rvap
145 
146 ! call qsat( LST1(i,j), PRSS(i,j), qdry, & ! [IN]
147 ! QVsat ) ! [OUT]
148  call qsat( lst1(i,j), rhos(i,j), & ! [IN]
149  qvsat ) ! [OUT]
150 
151  qvs(i,j) = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * qvsat
152 
153  uabs = sqrt( wa(i,j)**2 + ua(i,j)**2 + va(i,j)**2 )
154 
155  call bulkflux( tmpa(i,j), lst1(i,j), & ! [IN]
156  prsa(i,j), prss(i,j), & ! [IN]
157  qva(i,j), qvs(i,j), & ! [IN]
158  uabs, z1(i,j), pbl(i,j), & ! [IN]
159  z0m(i,j), z0h(i,j), z0e(i,j), & ! [IN]
160  ustar(i,j), tstar(i,j), qstar(i,j), & ! [OUT]
161  wstar(i,j), rlmo(i,j), ra, & ! [OUT]
162  fracu10(i,j), fract2(i,j), fracq2(i,j) ) ! [OUT]
163 
164  if ( uabs < eps ) then
165  zmflx(i,j) = 0.0_rp
166  xmflx(i,j) = 0.0_rp
167  ymflx(i,j) = 0.0_rp
168  else
169  mflux = - rhos(i,j) * ustar(i,j)**2
170  zmflx(i,j) = mflux * wa(i,j) / uabs
171  xmflx(i,j) = mflux * ua(i,j) / uabs
172  ymflx(i,j) = mflux * va(i,j) / uabs
173  end if
174 
175  else
176 
177  ! not calculate surface flux
178 
179  calc_flag(i,j) = .false.
180 
181  zmflx(i,j) = undef
182  xmflx(i,j) = undef
183  ymflx(i,j) = undef
184  ustar(i,j) = undef
185  tstar(i,j) = undef
186  qstar(i,j) = undef
187  wstar(i,j) = undef
188  rlmo(i,j) = undef
189  u10(i,j) = undef
190  v10(i,j) = undef
191  t2(i,j) = undef
192  q2(i,j) = undef
193 
194  end if
195 
196  end do
197  end do
198 
199  call bulkflux_diagnose_surface( lia, lis, lie, lja, ljs, lje, &
200  ua(:,:), va(:,:), & ! (in)
201  tmpa(:,:), qva(:,:), & ! (in)
202  lst1(:,:), qvs(:,:), & ! (in)
203  z1(:,:), z0m(:,:), z0h(:,:), z0e(:,:), & ! (in)
204  u10(:,:), v10(:,:), t2(:,:), q2(:,:), & ! (out)
205  mask = calc_flag(:,:), & ! (in)
206  fracu10 = fracu10(:,:), & ! (in)
207  fract2 = fract2(:,:), & ! (in)
208  fracq2 = fracq2(:,:) ) ! (in)
209 
210  return

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

Referenced by mod_land_driver::land_driver_calc_tendency().

Here is the caller graph for this function:
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:68
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
scale_const::const_stb
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:53
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_bulkflux::bulkflux
procedure(bc), pointer, public bulkflux
Definition: scale_bulkflux.F90:84