SCALE-RM
Functions/Subroutines | Variables
scale_land_phy_snow_ky90 Module Reference

module land / physics / snow / ky90 More...

Functions/Subroutines

subroutine, public land_phy_snow_ky90_setup
 Setup. More...
 
subroutine, public land_phy_snow_ky90 (LIA, LIS, LIE, LJA, LJS, LJE, SFLX_rain, SFLX_snow, PRSA, TA, QA, WA, UA, VA, DENS, SFLX_RAD_dn, LANDUSE_fact_land, dt, TSNOW, SWE, SDepth, SDzero, nosnowsec, Salbedo, SFLX_SH, SFLX_LH, SFLX_GH, SNOW_LAND_GH, SNOW_LAND_Water, SNOW_frac)
 Main routine for land submodel. More...
 
subroutine snow_ky90_main (TSNOW, SWE, DEPTH, ZNSNOW, nosnowsec, ALBEDO_out, Emiss, HFLUX, LATENTFLUX, GFLUX, EvapFLX, QCC, QFUSION, MELT, SWEMELT, Gflux2land, SFLX_SNOW, TA, UA, RH, DENS, SW, LW, time)
 snow model main routine More...
 
subroutine groundflux (TS, TA, UA, RH, rhoair, ALPHA, SW, LW, GFLUX, RFLUX, SFLUX, LINFLUX, LOUTFLUX, HFLUX, LATENTFLUX)
 
subroutine cal_param (ZN1, TS1, GFLUX, TA, UA, RH, rhoair, LW, time)
 
subroutine check_applicability (GFLUX, TS1, ZN1, TA, UA, RH, rhoair, LW, GFLUX_res, beta, time)
 
subroutine snowdepth (GFLUX, ZN1, ZN2, time)
 
subroutine recalculatez (ZN1, TS, GFLUX, ZN2, time)
 
subroutine calculationmo (GFLUX, CSRHOS, ZN1, TS1, ZN2, TS2, MELT, QCC, QFUSION, time)
 
subroutine calculationnomo (GFLUX, CSRHOS, ZN1, TS1, ZN2, TS2, MELT, QCC, QFUSION, time)
 
subroutine check_res (ZN1, ZN2, TS1, TS2, GFLUX, TA, UA, RH, rhoair, LW, flag, time)
 
subroutine cal_r1r2 (ZN1, TS1, GFLUX, TA, UA, RH, rhoair, LW, time)
 

Variables

real(rp), public w0
 
real(rp), public rhosnow = 400.0_RP
 

Detailed Description

module land / physics / snow / ky90

Description
Profile model for snowpack by Kondo and Yamazaki (1990) Kondo and Yamazaki (1990): Journal of Applied Meteorology, 29, pp375-384.
Author
Team SCALE
NAMELIST
  • PARAM_LAND_PHY_SNOW_KY90
    nametypedefault valuecomment
    ALBEDO_CONST logical .true.
    SNOW_CONDUCTIVITY real(RP) 0.42_RP
    WATER_CONTENT real(RP) 0.1_RP
    SNOW_HEAT_CAPACITYRHO real(RP) 8.4e+05_RP
    SNOW_RHO real(RP) 400.0_RP
    SNOWDEPTH_INITIAL real(RP) 0.0_RP
    ALBEDO_VALUE real(RP) 0.686_RP
    DEBUG logical .false.

History Output
namedescriptionunitvariable
SNOW_MELT Heat used for snow melt J/m2 MELT
SNOW_QCC Heat used for changing temperature profile J/m2 QCC
SNOW_QFUSION Heat used for phase change of snow J/m2 QFUSION
SNOW_SWEMELT Equivalent water of melt snow kg/m2 SWEMELT

Function/Subroutine Documentation

◆ land_phy_snow_ky90_setup()

subroutine, public scale_land_phy_snow_ky90::land_phy_snow_ky90_setup ( )

Setup.

Definition at line 82 of file scale_land_phy_snow_ky90.F90.

References scale_io::io_fid_conf, scale_prc::prc_abort(), rhosnow, and w0.

Referenced by mod_land_driver::land_driver_setup().

82  use scale_prc, only: &
83  prc_abort
84  implicit none
85  real(RP) :: snow_conductivity = 0.42_rp
86  real(RP) :: water_content = 0.1_rp
87  real(RP) :: snow_heat_capacityrho = 8.4e+05_rp
88  real(RP) :: snow_rho = 400.0_rp
89  real(RP) :: snowdepth_initial = 0.0_rp
90  real(RP) :: albedo_value = 0.686_rp
91 
92  namelist / param_land_phy_snow_ky90 / &
93  albedo_const, &
94  snow_conductivity, &
95  water_content, &
96  snow_heat_capacityrho, &
97  snow_rho, &
98  snowdepth_initial, &
99  albedo_value, &
100  debug
101 
102  integer :: ierr
103  !---------------------------------------------------------------------------
104 
105  log_newline
106  log_info("LAND_PHY_SNOW_KY90_setup",*) 'Setup'
107 
108  !--- read namelist
109  rewind(io_fid_conf)
110  read(io_fid_conf,nml=param_land_phy_snow_ky90,iostat=ierr)
111  if( ierr < 0 ) then !--- missing
112  log_info("LAND_PHY_SNOW_KY90_setup",*) 'Not found namelist. Default used.'
113  elseif( ierr > 0 ) then !--- fatal error
114  log_error("LAND_PHY_SNOW_KY90_setup",*) 'Not appropriate names in namelist PARAM_LAND_PHY_SNOW_KY90. Check!'
115  call prc_abort
116  endif
117  log_nml(param_land_phy_snow_ky90)
118 
119  lambdas = snow_conductivity
120  w0 = water_content
121  csrhos = snow_heat_capacityrho ! [J/m3/K]
122  rhosnow = snow_rho
123  cs = csrhos/rhosnow ! [J/kg/K]
124  !LOG_WARN("LAND_PHY_SNOW_KY90_setup",*) "Specific heat capacity of snow [J/kg/K]: ",CS
125  albedo = albedo_value
126 
127  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_phy_snow_ky90()

subroutine, public scale_land_phy_snow_ky90::land_phy_snow_ky90 ( 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)  SFLX_rain,
real(rp), dimension (lia,lja), intent(in)  SFLX_snow,
real(rp), dimension (lia,lja), intent(in)  PRSA,
real(rp), dimension (lia,lja), intent(in)  TA,
real(rp), dimension (lia,lja), intent(in)  QA,
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)  DENS,
real(rp), dimension(lia,lja,n_rad_dir,n_rad_rgn), intent(in)  SFLX_RAD_dn,
real(rp), dimension(lia,lja), intent(in)  LANDUSE_fact_land,
real(dp), intent(in)  dt,
real(rp), dimension (lia,lja), intent(inout)  TSNOW,
real(rp), dimension (lia,lja), intent(inout)  SWE,
real(rp), dimension (lia,lja), intent(inout)  SDepth,
real(rp), dimension (lia,lja), intent(inout)  SDzero,
real(rp), dimension (lia,lja), intent(inout)  nosnowsec,
real(rp), dimension (lia,lja,2), intent(out)  Salbedo,
real(rp), dimension (lia,lja), intent(out)  SFLX_SH,
real(rp), dimension (lia,lja), intent(out)  SFLX_LH,
real(rp), dimension (lia,lja), intent(out)  SFLX_GH,
real(rp), dimension (lia,lja), intent(out)  SNOW_LAND_GH,
real(rp), dimension(lia,lja), intent(out)  SNOW_LAND_Water,
real(rp), dimension (lia,lja), intent(out)  SNOW_frac 
)

Main routine for land submodel.

Definition at line 149 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_epsvap, scale_const::const_i_lw, scale_const::const_i_sw, scale_const::const_tem00, scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_prc::prc_abort(), and snow_ky90_main().

Referenced by mod_land_driver::land_driver_calc_tendency().

149  use scale_prc, only: &
150  prc_abort
151  use scale_file_history, only: &
152  file_history_in
153  use scale_atmos_saturation, only: &
154  qsatf => atmos_saturation_psat_all ! better to change name from qsatf to qsat
155  use scale_const, only: &
156  t0 => const_tem00, &
157  i_sw => const_i_sw, &
158  i_lw => const_i_lw, &
159  epsvap => const_epsvap
160  implicit none
161  integer, intent(in) :: lia, lis, lie
162  integer, intent(in) :: lja, ljs, lje
163 
164  ! input data
165  real(RP), intent(in) :: sflx_rain (lia,lja)
166  real(RP), intent(in) :: sflx_snow (lia,lja)
167  real(RP), intent(in) :: prsa (lia,lja)
168  real(RP), intent(in) :: ta (lia,lja)
169  real(RP), intent(in) :: wa (lia,lja)
170  real(RP), intent(in) :: ua (lia,lja)
171  real(RP), intent(in) :: va (lia,lja)
172  real(RP), intent(in) :: qa (lia,lja) ! specific humidity [kg/kg]
173  real(RP), intent(in) :: dens (lia,lja)
174  real(RP), intent(in) :: sflx_rad_dn(lia,lja,n_rad_dir,n_rad_rgn)
175  real(DP), intent(in) :: dt ! dt of land
176  real(RP), intent(in) :: landuse_fact_land(lia,lja)
177 
178  ! prognostic variables
179  real(RP), intent(inout) :: tsnow (lia,lja) ! snow temperature [K]
180  real(RP), intent(inout) :: swe (lia,lja) ! equivalent water [kg/m2]
181  real(RP), intent(inout) :: sdepth (lia,lja) ! depth of melting point [m]
182  real(RP), intent(inout) :: sdzero (lia,lja) ! total snow depth [m]
183  real(RP), intent(inout) :: nosnowsec (lia,lja) ! elapsed time of no snow condition [s]
184 
185  ! updated variables
186  real(RP), intent(out) :: salbedo (lia,lja,2) ! snow albedo [-]
187  real(RP), intent(out) :: sflx_sh (lia,lja) ! sensible heat flux between atmos and snow [W/m2]
188  real(RP), intent(out) :: sflx_lh (lia,lja) ! latente heat flux between atmos and snow [W/m2]
189  real(RP), intent(out) :: sflx_gh (lia,lja) ! whole snowpack Ground flux [W/m2]
190  real(RP), intent(out) :: snow_land_gh (lia,lja) ! heat flux from snow to land [W/m2]
191  real(RP), intent(out) :: snow_land_water(lia,lja) ! water flux from snow to land [W/m2]
192  real(RP), intent(out) :: snow_frac (lia,lja) ! snow fraction, defined by time direction [-]
193 
194  real(RP) :: qcc (lia,lja)
195  real(RP) :: qfusion (lia,lja)
196  real(RP) :: melt (lia,lja)
197  real(RP) :: swemelt (lia,lja)
198  real(RP) :: sflx_evap (lia,lja) ! evaporation due to LH [kg/m2/s]
199 
200  ! works
201  real(RP) :: tsnow1 ! updated snow surface temperature [K]
202  real(RP) :: znsnow1 ! updated freezing depth [m]
203  real(RP) :: swe1 ! updated snow water equivalence [kg/m2]
204  real(RP) :: depth1 ! updated snow depth [m]
205 
206  real(RP), parameter :: uabs_min = 0.1_rp
207  real(RP) :: uabs
208  real(RP) :: qasat
209  real(RP) :: rh
210  real(RP) :: qdry, psat
211  real(RP) :: sflx_sw_dn, sflx_lw_dn
212 
213  integer :: k, i, j
214  !---------------------------------------------------------------------------
215  log_progress(*) 'Land / physics / snow / KY90'
216 
217  do j = ljs, lje
218  do i = lis, lie
219 
220  if( ( landuse_fact_land(i,j) > 0.0_rp ) .and. &
221  ( swe(i,j)>0. .or. sflx_snow(i,j)>0. ) )then
222 
223  uabs = max( sqrt( ua(i,j)**2 + va(i,j)**2 + wa(i,j)**2 ), uabs_min )
224  !Uabs = sqrt( UA(i,j)**2 + VA(i,j)**2 + WA(i,j)**2 )
225 
226  !qdry = 1.0_RP - QA(i,j)
227  !call qsatf( TA(i,j), PRSA(i,j), qdry, & ![IN]
228  ! QAsat ) ![OUT]
229  call qsatf( ta(i,j), psat )
230  qasat = epsvap * psat / ( prsa(i,j) - ( 1.0_rp-epsvap ) * psat )
231  rh = qa(i,j) / qasat
232  if ( debug ) then
233  log_info("LAND_PHY_SNOW_KY90",*) "RH, ",rh," DENS (1.289 org) ",dens(i,j)
234  end if
235 
236  sflx_lw_dn = sflx_rad_dn(i,j,i_r_diffuse,i_r_ir )
237  sflx_sw_dn = sflx_rad_dn(i,j,i_r_direct ,i_r_nir) &
238  + sflx_rad_dn(i,j,i_r_diffuse,i_r_nir) &
239  + sflx_rad_dn(i,j,i_r_direct ,i_r_vis) &
240  + sflx_rad_dn(i,j,i_r_diffuse,i_r_vis)
241 
242  tsnow1 = tsnow(i,j)
243  swe1 = swe(i,j)
244  depth1 = sdepth(i,j)
245  znsnow1 = sdzero(i,j)
246 
247  call snow_ky90_main( tsnow1, & ! [INOUT]
248  swe1, & ! [INOUT]
249  depth1, & ! [INOUT]
250  znsnow1, & ! [INOUT]
251  nosnowsec(i,j), & ! [INOUT]
252  salbedo(i,j,i_sw), & ! [OUT]
253  salbedo(i,j,i_lw), & ! [OUT]
254  sflx_sh(i,j), & ! [OUT]
255  sflx_lh(i,j), & ! [OUT]
256  sflx_gh(i,j), & ! [OUT]
257  sflx_evap(i,j), & ! [OUT]
258  qcc(i,j), & ! [OUT]
259  qfusion(i,j), & ! [OUT]
260  melt(i,j), & ! [OUT]
261  swemelt(i,j), & ! [OUT]
262  snow_land_gh(i,j), & ! [OUT]
263  sflx_snow(i,j), & ! [IN] ! [kg/m2/s]
264  ta(i,j), & ! [IN]
265  uabs, & ! [IN]
266  rh, & ! [IN]
267  dens(i,j), & ! [IN]
268  sflx_sw_dn, & ! [IN]
269  sflx_lw_dn, & ! [IN]
270  dt )
271 
272 
273  snow_land_gh(i,j) = snow_land_gh(i,j) / dt ! [J/m2] -> [J/m2/s]
274  snow_land_water(i,j) = sflx_rain(i,j) + swemelt(i,j) / dt ! [kg/m2] -> [kg/m2/s]
275 
276  if ( swe1 <= 0. .and. swe(i,j) <= 0. ) then ! no accumulated snow during the time step
277  snow_frac(i,j) = 0.0_rp
278  else
279  snow_frac(i,j) = 1.0_rp
280  endif
281 
282  tsnow(i,j) = tsnow1
283  swe(i,j) = swe1
284  sdepth(i,j) = depth1
285  sdzero(i,j) = znsnow1
286 
287  if ( debug ) then
288  log_info("LAND_PHY_SNOW_KY90",*) "SNOW_frac, SWE, TSNOW", snow_frac(i,j), swe(i,j), tsnow(i,j)
289  end if
290 
291  else
292 
293  snow_land_water(i,j) = sflx_rain(i,j)
294  snow_frac(i,j) = 0.0_rp
295 
296  tsnow(i,j) = t0 !!!
297  swe(i,j) = 0.0_rp !!!
298  sdepth(i,j) = 0.0_rp !!!
299  sdzero(i,j) = 0.0_rp !!!
300  salbedo(i,j,:) = 0.0_rp !!!
301  sflx_sh(i,j) = 0.0_rp
302  sflx_lh(i,j) = 0.0_rp
303  sflx_gh(i,j) = 0.0_rp
304  sflx_evap(i,j) = 0.0_rp
305  qcc(i,j) = 0.0_rp
306  qfusion(i,j) = 0.0_rp
307  melt(i,j) = 0.0_rp
308  swemelt(i,j) = 0.0_rp
309  snow_land_gh(i,j) = 0.0_rp
310 
311  endif
312 
313  call file_history_in( qcc(:,:), 'SNOW_QCC', 'Heat used for changing temperature profile', 'J/m2', dim_type='XY' )
314  call file_history_in( qfusion(:,:), 'SNOW_QFUSION', 'Heat used for phase change of snow', 'J/m2', dim_type='XY' )
315  call file_history_in( melt(:,:), 'SNOW_MELT', 'Heat used for snow melt', 'J/m2', dim_type='XY' )
316  call file_history_in( swemelt(:,:), 'SNOW_SWEMELT', 'Equivalent water of melt snow', 'kg/m2', dim_type='XY' )
317 
318  end do
319  end do
320 
321  return
module atmosphere / saturation
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:93
integer, public va
Definition: scale_index.F90:35
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
integer, public qa
integer, parameter, public i_lw
integer, parameter, public i_sw
module PROCESS
Definition: scale_prc.F90:11
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:69
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:94
module file_history
Here is the call graph for this function:
Here is the caller graph for this function:

◆ snow_ky90_main()

subroutine scale_land_phy_snow_ky90::snow_ky90_main ( real(rp), intent(inout)  TSNOW,
real(rp), intent(inout)  SWE,
real(rp), intent(inout)  DEPTH,
real(rp), intent(inout)  ZNSNOW,
real(rp), intent(inout)  nosnowsec,
real(rp), intent(out)  ALBEDO_out,
real(rp), intent(out)  Emiss,
real(rp), intent(out)  HFLUX,
real(rp), intent(out)  LATENTFLUX,
real(rp), intent(out)  GFLUX,
real(rp), intent(out)  EvapFLX,
real(rp), intent(out)  QCC,
real(rp), intent(out)  QFUSION,
real(rp), intent(out)  MELT,
real(rp), intent(out)  SWEMELT,
real(rp), intent(out)  Gflux2land,
real(rp), intent(in)  SFLX_SNOW,
real(rp), intent(in)  TA,
real(rp), intent(in)  UA,
real(rp), intent(in)  RH,
real(rp), intent(in)  DENS,
real(rp), intent(in)  SW,
real(rp), intent(in)  LW,
real(dp), intent(in)  time 
)

snow model main routine

Definition at line 352 of file scale_land_phy_snow_ky90.F90.

References cal_param(), calculationmo(), calculationnomo(), check_applicability(), check_res(), scale_const::const_tem00, groundflux(), recalculatez(), rhosnow, snowdepth(), and w0.

Referenced by land_phy_snow_ky90().

352  use scale_const, only: &
353  t0 => const_tem00
354 
355  implicit none
356  ! prognostic variables
357  real(RP), intent(inout) :: tsnow ! snow temperature [K]
358  real(RP), intent(inout) :: swe ! equivalent water [kg/m2]
359  real(RP), intent(inout) :: depth ! depth of melting point [m]
360  real(RP), intent(inout) :: znsnow ! total snow depth [m]
361  ! update variables
362  real(RP), intent(inout) :: nosnowsec
363  real(RP), intent(out) :: albedo_out
364  real(RP), intent(out) :: emiss
365 
366  ! output variables
367  real(RP), intent(out) :: hflux ! HFLUX = whole snow Sensible heat flux [W/m2]
368  real(RP), intent(out) :: latentflux ! LATENTFLUX = whole snow Latent heat flux [W/m2]
369  real(RP), intent(out) :: gflux ! GFLUX = whole snow Ground flux [W/m2]
370  real(RP), intent(out) :: evapflx ! Evapolation due to LATENTFLUX [kg/m2/s]
371 
372  real(RP), intent(out) :: qcc ! QCC = heat used for change snow condition to isothermal [J m^-2]
373  real(RP), intent(out) :: qfusion ! QFUSION = heat used for change snow condition to melt point [J m^-2]
374  real(RP), intent(out) :: melt ! MELT = heat used for snow run off production [J m^-2 s^-1]
375  real(RP), intent(out) :: swemelt ! snow water equivalent ! [kg/m2]
376  real(RP), intent(out) :: gflux2land ! Residual heat, goes to land model ! [J/m2]
377 
378 
379  ! input data
380  real(RP), intent(in) :: sflx_snow
381  real(RP), intent(in) :: ta
382  real(RP), intent(in) :: ua
383  real(RP), intent(in) :: rh
384  real(RP), intent(in) :: dens
385  real(RP), intent(in) :: sw
386  real(RP), intent(in) :: lw
387  real(DP), intent(in) :: time
388 
389  ! initial value
390  real(RP) :: tsnow0 ! Initial time snow surface temperature [K]
391  real(RP) :: znsnow0 ! ZNSNOW0 = initial freezing depth [m]
392  real(RP) :: swe0 ! SWE0 = snow depth initial value in snow water equivalen [kg/m2]
393  real(RP) :: depth0 ! DEPTH0 = initial snow depth [m]
394 
395  ! works
396  real(RP) :: snow ! per dt
397  real(RP) :: rflux ! RFLUX = whole snow net long wave radiation flux [W/m2]
398  real(RP) :: linflux ! LINFLUX = whole snow downward long wave radiation flux [W/m2]
399  real(RP) :: loutflux ! LOUTFLUX = whole snow upward long wave radiation flux [W/m2]
400  real(RP) :: sflux ! SFLUX = whole snow net short wave radiation flux [W/m2]
401  real(RP) :: deltadepth
402 
403  real(RP) :: beta
404 
405 
406 !---------------------------------------------- !
407 ! ALL START HERE !
408 !---------------------------------------------- !
409 
410  ! initialize
411  zn_flag = 0
412  ts_flag = 0
413 
414  qcc = 0.0_rp
415  qfusion = 0.0_rp
416  melt = 0.0_rp
417  swemelt = 0.0_rp
418 
419  ! snowfall during timestep
420  snow = sflx_snow * time
421 
422  ! save previous timestep
423  tsnow0 = tsnow
424  znsnow0 = znsnow
425  swe0 = swe
426  depth0 = depth
427 
428  ! update
429  swe0 = swe0 + snow ! update according to snowfall
430  depth0 = depth0 + (snow /rhosnow)
431  znsnow0 = znsnow0 + (snow /rhosnow)
432 
433  if ( debug ) then
434  log_info("SNOW_ky90_main",*) "UA, SNOW,SFLX_SNOW,time : ", ua, snow, sflx_snow, time
435  log_info("SNOW_ky90_main",*) "SWE , TSNOW, and TA : ", swe0, tsnow0, ta
436  log_info("SNOW_ky90_main",*) "DEPTH is: ", depth0
437  log_info("SNOW_ky90_main",*) "ZN beginning: ", znsnow0
438  end if
439 
440 !----- Calculating albedo -------------------------------------------!
441 
442  if (snow > 0.0_rp) then ! snowfall
443  nosnowsec = 0.0_rp
444  else
445  nosnowsec = nosnowsec + time
446  endif
447  if ( .not. albedo_const ) then
448  albedo = albedomin + (albedomax-albedomin)*exp(-1.0_rp*nosnowsec/real(4.0_rp*24.0_rp*3600.0_rp))
449  endif
450 
451  albedo_out = albedo
452  emiss = 1.0_rp - epsilon
453  if ( debug ) then
454  log_info("SNOW_ky90_main",*) "Albedo ",albedo
455  end if
456 
457 !----- Energy balance at snow surface -------------------------------!
458 
459  call groundflux (tsnow0, ta, ua, rh, dens, albedo, sw, lw, & ! [IN]
460  gflux, rflux, sflux, linflux, loutflux, hflux, latentflux) ! [OUT]
461 
462  ! SWE change due to latent heat flux
463  evapflx = latentflux / lv ! [kg/m2/s] positive => snow decrease
464 
465  swe0 = swe0 - evapflx * time ! [kg/m2]
466  deltadepth = (evapflx * time) /rhosnow
467  depth0 = depth0 - deltadepth ! [m]
468  znsnow0 = znsnow0 - deltadepth ! [m]
469 
470 !! Check whether GFLUX (energy into snowpack) is enough to melt all snow.
471 !! If GFLUX is enough, the model melts all snow and then go to next timestep.
472 
473  call check_allsnowmelt (gflux, tsnow0, znsnow0, depth0, sflag, time)
474  if(sflag .eq. 1)then
475  log_info("SNOW_ky90_main",*) 'LAND/snow: All snow melt'
476  qcc = 0.5_rp*csrhos*znsnow0*(t0-tsnow0)
477  qfusion = w0*rhosnow*lf*znsnow0
478  melt = (1.0_rp-w0)*rhosnow*lf*depth0 ! [J/m2]
479  tsnow = t0
480  znsnow = 0.0_rp
481  depth = 0.0_rp
482  swe = 0.0_rp
483  !SWEMELT = MELT /((1.0_RP-W0)*LF)
484  swemelt = rhosnow*depth0
485  gflux2land = gflux*time - (qcc + qfusion + melt) ! [J/m2]
486 
487  else
488 
489  call cal_param (znsnow0, tsnow0, gflux, ta, ua, rh, dens, lw, time)
490 
491  ! check whether the model has solution
492  call check_applicability (gflux, tsnow0, znsnow0, ta, ua, rh, dens, lw, gres, beta, time)
493  if ((gres >= 0.0_rp).and.(beta >= 0.0_rp)) then
494  log_info("SNOW_ky90_main",*) 'LAND/snow model is not appropriate',gres,beta
495  qcc = 0.5_rp*csrhos*znsnow0*(t0-tsnow0)
496  qfusion = w0*rhosnow*lf*znsnow0
497  melt = gres
498  tsnow = t0
499  znsnow = 0.0_rp
500  else
501 
502  !if (t==1)then
503  ! call cal_R1R2 (ZNSNOW0, TSNOW0, GFLUX, TA, UA, RH, rhoair, LW, time)
504  !endif
505 
506  call snowdepth ( gflux, znsnow0, znsnow, time)
507 
508  if ( debug ) then
509  log_info("SNOW_ky90_main",*) "ZN is: ", znsnow
510  end if
511  IF(znsnow < zmin) THEN
512  zn_flag = 1
513  if ( debug ) then
514  log_info("SNOW_ky90_main",*) "ZN is replaced to: ", znsnow ," to ", zmin
515  end if
516  znsnow = zmin
517  ELSE IF(znsnow > depth0) THEN
518  zn_flag = 2
519  if ( debug ) then
520  log_info("SNOW_ky90_main",*) "ZN is replaced to: ", znsnow ," to ", depth0
521  end if
522  znsnow = depth0
523  END IF
524 
525  ! This equation is to calculate TSN
526  call equation415(lambdas, c2, znsnow, rh, qsat, tsnow0, znsnow0, gflux, ta, ua, dens, lw, tsnow, time)
527 
528  if ( debug ) then
529  log_info("SNOW_ky90_main",*) 'TSNOW is: ', tsnow
530  end if
531 
532  if (tsnow > t0) then
533  ts_flag = 1
534  tsnow = t0
535  call recalculatez(znsnow0, tsnow0, gflux, znsnow, time)
536  call check_res(znsnow0, znsnow, tsnow0, tsnow, gflux, ta, ua, rh, dens, lw, "1", time)
537  IF (znsnow < zmin) THEN
538  zn_flag = 4
539  if ( debug ) then
540  log_info("SNOW_ky90_main",*) "ZN is updated/replaced to: ", znsnow ," to ", zmin
541  end if
542  znsnow = zmin
543  ELSE IF(znsnow > depth0) THEN
544  zn_flag = 5
545  if ( debug ) then
546  log_info("SNOW_ky90_main",*) "ZN is updated/replaced to: ", znsnow ," to ", depth0
547  end if
548  znsnow = depth0
549  ELSE
550  zn_flag = 3
551  END IF
552 
553  if ( debug ) then
554  log_info("SNOW_ky90_main",*) 'TSNOW is updated: ', tsnow
555  log_info("SNOW_ky90_main",*) 'ZNSNOW is updated: ', znsnow
556  end if
557  else
558  call check_res( znsnow0, znsnow, tsnow0, tsnow, gflux, ta, ua, rh, dens, lw, "0", time)
559  endif
560 
561  IF ( (znsnow-zmin)/zmin < 0.00001 ) THEN
562  call calculationmo( &
563  gflux, csrhos, znsnow0, tsnow0, znsnow, tsnow, &
564  melt, qcc, qfusion, time)
565  ELSE
566  call calculationnomo( &
567  gflux, csrhos, znsnow0, tsnow0, znsnow, tsnow, &
568  melt, qcc, qfusion, time)
569  END IF
570 
571  endif ! Gres & beta
572 
573  gflux2land = gflux*time - (qcc + qfusion + melt)
574  if ( debug ) then
575  log_info("SNOW_ky90_main",*) "### ZN_flag = ", zn_flag, "TS_flag = ", ts_flag
576  log_info("SNOW_ky90_main",*) '### Heat flux from snowpack to land surface: ', gflux2land
577  end if
578 
579  deltadepth = melt / ((1.0_rp-w0)*lf*rhosnow)
580  swemelt = rhosnow*deltadepth
581  swe = swe0 - swemelt
582  depth = depth0 - deltadepth
583  if (depth < znsnow) then
584  ! NOTICE: energy budget has not been considered thought this process yet.
585  if ( debug ) then
586  log_info("SNOW_ky90_main",*) "replace ZNSNOW <= DEPTH"
587  end if
588  znsnow = depth
589  endif
590 
591  if ( debug ) then
592  log_info("SNOW_ky90_main",*) 'MELT in water equivalent is: ', swemelt
593  log_info("SNOW_ky90_main",*) 'SWE0 is: ', swe
594  log_info("SNOW_ky90_main",*) 'DELTADEPTH is: ', deltadepth
595  log_info("SNOW_ky90_main",*) 'DEPTH0 is: ', depth
596  log_info("SNOW_ky90_main",*) 'ZNSNOW0 is: ', znsnow
597  end if
598 
599  endif
600 
601  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
Here is the call graph for this function:
Here is the caller graph for this function:

◆ groundflux()

subroutine scale_land_phy_snow_ky90::groundflux ( real(rp), intent(in)  TS,
real(rp), intent(in)  TA,
real(rp), intent(in)  UA,
real(rp), intent(in)  RH,
real(rp), intent(in)  rhoair,
real(rp), intent(in)  ALPHA,
real(rp), intent(in)  SW,
real(rp), intent(in)  LW,
real(rp), intent(out)  GFLUX,
real(rp), intent(out)  RFLUX,
real(rp), intent(out)  SFLUX,
real(rp), intent(out)  LINFLUX,
real(rp), intent(out)  LOUTFLUX,
real(rp), intent(out)  HFLUX,
real(rp), intent(out)  LATENTFLUX 
)

Definition at line 608 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

608  implicit none
609 
610  real(RP), intent(in) :: ts, ta, ua, rh, rhoair, alpha, sw, lw
611  real(RP), intent(out) :: gflux, rflux, sflux, linflux, loutflux, hflux, latentflux
612 
613  esat = 0.6112_rp * exp( (17.67_rp * (ta-273.15_rp)) / (ta-29.66_rp) )
614  qsat = 0.622_rp * esat / 101.325_rp
615  deltaqsat = 0.622_rp * lv * qsat /(287.04_rp * (ta**2))
616  cp = 1004.67_rp * (1.0_rp + (0.84_rp * qsat))
617 
618  ! Into snowpack is positive
619  sflux = (1.0_rp-alpha) * sw
620  rflux = epsilon * (lw-(sigma*(ts**4)))
621  linflux = epsilon * lw
622  loutflux = -1.0_rp * (epsilon * sigma * (ts**4))
623  hflux = cp * rhoair * ch * ua * (ts-ta)
624  latentflux = lv * rhoair * ce * ua * ((1-rh)*qsat + (deltaqsat*(ts-ta)))
625 
626  gflux = (sflux + rflux - hflux - latentflux)
627 
628  if ( debug ) then
629  log_info("LAND_PHY_SNOW_KY90_groundflux",*) "-------------- groundflux --------------"
630  log_info_cont(*) "GFLUX is: ", gflux
631  log_info_cont(*) "SFLUX: ", sflux
632  log_info_cont(*) "RFLUX: ", rflux, linflux+loutflux
633  log_info_cont(*) " (LONG in: ", linflux,")"
634  log_info_cont(*) " (LONG out: ", loutflux,")"
635  log_info_cont(*) "HFLUX is: ", hflux
636  log_info_cont(*) "LATENT FLUX: ", latentflux
637  end if
638 
639 return
Here is the caller graph for this function:

◆ cal_param()

subroutine scale_land_phy_snow_ky90::cal_param ( real(rp), intent(in)  ZN1,
real(rp), intent(in)  TS1,
real(rp), intent(in)  GFLUX,
real(rp), intent(in)  TA,
real(rp), intent(in)  UA,
real(rp), intent(in)  RH,
real(rp), intent(in)  rhoair,
real(rp), intent(in)  LW,
real(dp), intent(in)  time 
)

Definition at line 676 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

676  use scale_const, only: &
677  t0 => const_tem00
678 
679  implicit none
680  real(RP), intent(in) :: zn1, ts1, ta, ua, rh, rhoair, lw
681  real(RP), intent(in) :: gflux
682  real(DP), intent(in) :: time
683 
684  c1 = csrhos*0.5_rp
685  c2 = (4.0_rp*epsilon*sigma*(ta**3)) + (cp*rhoair*ch*ua) + (lv*rhoair*ce*ua*deltaqsat)
686  c3 = w0*rhosnow*lf
687 
688  a0 = lambdas*(c3*zn1 + (c1*zn1*(t0-ts1)- gflux*time))
689  a1 = c2*(c3*zn1 + c1*zn1*(t0-ts1) - gflux*time) - c3*lambdas
690  a2 = c1* &
691  ( epsilon*(lw-(sigma*(ta**4))) + (c2*(ta-t0)) - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) ) &
692  - c2*c3
693 
694  if ( debug ) then
695  log_info("LAND_PHY_SNOW_KY90_cal_param",*) "-------------- snowdepth --------------"
696  log_info_cont(*) "C1",c1
697  log_info_cont(*) "C2",c2,(4.0_rp*epsilon*sigma*(ta**3)),(cp*rhoair*ch*ua), (lv*rhoair*ce*ua*deltaqsat)
698  log_info_cont(*) "C3",c3
699  log_info_cont(*) "A0",a0
700  log_info_cont(*) "A1",a1
701  log_info_cont(*) "A2",a2
702  end if
703 
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
Here is the caller graph for this function:

◆ check_applicability()

subroutine scale_land_phy_snow_ky90::check_applicability ( real(rp), intent(in)  GFLUX,
real(rp), intent(in)  TS1,
real(rp), intent(in)  ZN1,
real(rp), intent(in)  TA,
real(rp), intent(in)  UA,
real(rp), intent(in)  RH,
real(rp), intent(in)  rhoair,
real(rp), intent(in)  LW,
real(rp), intent(out)  GFLUX_res,
real(rp), intent(out)  beta,
real(dp), intent(in)  time 
)

Definition at line 708 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

708  use scale_const, only: &
709  t0 => const_tem00
710 
711  implicit none
712  real(RP), intent(in) :: gflux, ts1, zn1
713  real(RP), intent(in) :: ta, ua, rh, rhoair, lw
714  real(DP), intent(in) :: time
715  real(RP), intent(out) :: gflux_res, beta
716  real(RP) :: energy_in, energy_use_max
717 
718  energy_in = gflux * time
719  energy_use_max = 0.5_rp*csrhos*zn1*(t0-ts1) + w0*rhosnow*lf*zn1
720 
721  gflux_res = energy_in - energy_use_max ! residual energy after being used to melt all snow
722 
723  !
724  beta = (lambdas/c2)* &
725  ( t0 - ( ta*c2 - lv*rhoair*ce*ua*(1.0_rp-rh)*qsat + epsilon*(lw-sigma*(ta**4)))/c2 )
726 
727  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
Here is the caller graph for this function:

◆ snowdepth()

subroutine scale_land_phy_snow_ky90::snowdepth ( real(rp), intent(in)  GFLUX,
real(rp), intent(in)  ZN1,
real(rp), intent(out)  ZN2,
real(dp), intent(in)  time 
)

Definition at line 732 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

732 
733  implicit none
734 
735  real(RP), intent(in) :: gflux
736  real(RP), intent(in) :: zn1
737  real(DP), intent(in) :: time
738  real(RP), intent(out) :: zn2
739 
740 !calcultaing snowdepth
741 
742  !print *, -1.0_RP*A1
743  !print *, -1.0_RP*((A1**2.0_RP - 4.0_RP * A2 * A0)**0.5_RP)
744  !print *, 2.0_RP*A2
745  !print *, A1**2.0_RP, -1.0_RP * 4.0_RP * A2 * A0
746 
747  zn2 = ((-1.0_rp*a1) - ((a1**2.0_rp - 4.0_rp*a2*a0)**0.5_rp)) / (2.0_rp*a2)
748 
749  if ( debug ) then
750  log_info("LAND_PHY_SNOW_KY90_snowdepth",*) "ZN old = ",zn1, "ZN new = ",zn2
751  end if
752 
753 return
Here is the caller graph for this function:

◆ recalculatez()

subroutine scale_land_phy_snow_ky90::recalculatez ( real(rp), intent(in)  ZN1,
real(rp), intent(in)  TS,
real(rp), intent(in)  GFLUX,
real(rp), intent(out)  ZN2,
real(dp), intent(in)  time 
)

Definition at line 783 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00.

Referenced by snow_ky90_main().

783  use scale_const, only: &
784  t0 => const_tem00
785  implicit none
786 
787  real(RP), intent(in) :: zn1
788  real(RP), intent(in) :: ts, gflux
789  real(DP), intent(in) :: time
790  real(RP), intent(out) :: zn2
791 
792  zn2 = zn1 + ((c1*zn1*(t0-ts) - gflux*time) / c3)
793 
794  if ( debug ) then
795  log_info("LAND_PHY_SNOW_KY90_recalculateZ",*) "-------------- recalculate Z --------------"
796  end if
797 
798  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
Here is the caller graph for this function:

◆ calculationmo()

subroutine scale_land_phy_snow_ky90::calculationmo ( real(rp), intent(in)  GFLUX,
real(rp), intent(in)  CSRHOS,
real(rp), intent(in)  ZN1,
real(rp), intent(in)  TS1,
real(rp), intent(in)  ZN2,
real(rp), intent(in)  TS2,
real(rp), intent(out)  MELT,
real(rp), intent(out)  QCC,
real(rp), intent(out)  QFUSION,
real(dp), intent(in)  time 
)

Definition at line 804 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00, scale_prc::prc_abort(), rhosnow, and w0.

Referenced by snow_ky90_main().

804  use scale_const, only: &
805  t0 => const_tem00
806  use scale_prc, only: &
807  prc_abort
808  implicit none
809  real(RP), intent(in) :: gflux, csrhos, zn1, ts1, ts2, zn2
810  real(DP), intent(in) :: time
811  real(RP), intent(out) :: melt, qcc, qfusion
812 
813  qcc = 0.5_rp*csrhos*(zn1*(t0-ts1)-zn2*(t0-ts2))
814  qfusion = w0*rhosnow*lf*(zn1-zn2)
815  melt = ( gflux*time - qcc - qfusion )
816 
817  if ( debug ) then
818  log_info("LAND_PHY_SNOW_KY90_calculationMO",*) "--------------------------MELT----------------"
819  log_info_cont(*) "GFLUX*time is: ", gflux*time
820  log_info_cont(*) "QCC is : ", qcc
821  log_info_cont(*) "QFUSION is : ", qfusion
822  log_info_cont(*) "QMELT is : ", melt
823 
824  log_info_cont(*) qcc+qfusion+melt
825  log_info_cont(*) "diff= ", qcc + qfusion + melt - (gflux*time)
826  end if
827 
828  if ( abs(qcc+qfusion+melt - (gflux*time)) > 10.) then
829  log_error("LAND_PHY_SNOW_KY90_calculationMO",*) "Calculation is fault. Model would include bugs. Please check! Melt"
830  call prc_abort
831  endif
832 
833 
834  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calculationnomo()

subroutine scale_land_phy_snow_ky90::calculationnomo ( real(rp), intent(in)  GFLUX,
real(rp), intent(in)  CSRHOS,
real(rp), intent(in)  ZN1,
real(rp), intent(in)  TS1,
real(rp), intent(in)  ZN2,
real(rp), intent(in)  TS2,
real(rp), intent(out)  MELT,
real(rp), intent(out)  QCC,
real(rp), intent(out)  QFUSION,
real(dp), intent(in)  time 
)

Definition at line 840 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

840  use scale_const, only: &
841  t0 => const_tem00
842  implicit none
843  real(RP), intent(in) :: gflux, csrhos, zn1, ts1, ts2, zn2
844  real(DP), intent(in) :: time
845  real(RP), intent(out) :: melt, qcc, qfusion
846 
847  qcc = 0.5_rp*csrhos*(zn1*(t0-ts1)-zn2*(t0-ts2))
848  qfusion = w0*rhosnow*lf*(zn1-zn2)
849  melt = 0.0_rp
850 
851  if ( debug ) then
852  log_info("LAND_PHY_SNOW_KY90_calculationNoMO",*) "--------------------------NOMELT----------------"
853  log_info_cont(*) "GFLUX*time is: ", gflux*time
854  log_info_cont(*) "QCC is : ", qcc
855  log_info_cont(*) "QFUSION is : ", qfusion
856  log_info_cont(*) "QMELT is : ", melt
857 
858  log_info_cont(*) qcc+qfusion+melt
859  log_info_cont(*) "diff= ", qcc +qfusion - (gflux*time)
860  end if
861 
862  !if ( ABS(QCC+QFUSION - (GFLUX*time)) > 10.) then
863  ! LOG_ERROR("LAND_PHY_SNOW_KY90_calculationNoMO",*) "Calculation is fault. Model would include bugs. Please check! No Melt"
864  ! call PRC_abort
865  !endif
866 
867  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
Here is the caller graph for this function:

◆ check_res()

subroutine scale_land_phy_snow_ky90::check_res ( real(rp), intent(in)  ZN1,
real(rp), intent(in)  ZN2,
real(rp), intent(in)  TS1,
real(rp), intent(in)  TS2,
real(rp), intent(in)  GFLUX,
real(rp), intent(in)  TA,
real(rp), intent(in)  UA,
real(rp), intent(in)  RH,
real(rp), intent(in)  rhoair,
real(rp), intent(in)  LW,
character(len=*)  flag,
real(dp), intent(in)  time 
)

Definition at line 872 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00.

Referenced by snow_ky90_main().

872  use scale_const, only: &
873  t0 => const_tem00
874  implicit none
875  real(RP), intent(in) :: zn1, zn2, ts1, ts2, gflux
876  real(RP), intent(in) :: ta, ua, rh, rhoair, lw
877  real(DP), intent(in) :: time
878  character(len=*) :: flag
879  real(RP) :: r1 ! Eq. (2)
880  real(RP) :: r2 ! Eq. (8)
881  real(RP) :: r3 ! Eq. (8)
882 
883  r3 = -999.
884  r1 = c1 * (zn1*(t0-ts1) - zn2*(t0-ts2)) + c3*(zn1-zn2)-gflux*time
885  r2 = zn2*( epsilon*lw-epsilon*sigma*(ta**4) - (c2*(ts2-ta)) - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) ) &
886  + lambdas*(t0-ts2)
887  r3 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
888  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))-ts2
889 
890 
891  if ( debug ) then
892  log_info("LAND_PHY_SNOW_KY90_check_res",*) "R1 is : ", r1, "flag = ", flag
893  log_info("LAND_PHY_SNOW_KY90_check_res",*) "R2 is : ", r2, r3,"flag = ", flag
894 
895  if(abs(r1)>10000)then
896  log_info("LAND_PHY_SNOW_KY90_check_res",*) c1 * (zn1*(t0-ts1) - zn2*(t0-ts2)), c3*(zn1-zn2), -gflux*time
897  endif
898  end if
899 
900  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
Here is the caller graph for this function:

◆ cal_r1r2()

subroutine scale_land_phy_snow_ky90::cal_r1r2 ( real(rp), intent(in)  ZN1,
real(rp), intent(in)  TS1,
real(rp), intent(in)  GFLUX,
real(rp), intent(in)  TA,
real(rp), intent(in)  UA,
real(rp), intent(in)  RH,
real(rp), intent(in)  rhoair,
real(rp), intent(in)  LW,
real(dp), intent(in)  time 
)

Definition at line 905 of file scale_land_phy_snow_ky90.F90.

References scale_const::const_tem00.

905  use scale_const, only: &
906  t0 => const_tem00
907  implicit none
908  real(RP), intent(in) :: zn1, ts1
909  real(RP), intent(in) :: gflux, ta, ua, rh, rhoair, lw
910  real(DP), intent(in) :: time
911  real(RP) :: zn2, ts2
912  real(RP) :: r1 ! Eq. (2)
913  real(RP) :: r2 ! Eq. (8)
914  real(RP) :: r3 ! Eq. (8)
915 
916  real(RP) :: ts0,zn0
917  integer :: i,j
918  real(RP) :: ts_r1,ts_r2,zn_r1,zn_r2
919  character(len=3) :: ttt = ""
920 
921  real(RP) :: a,b,c,d,e,f,g
922 
923  ts0 = -50.0_rp + t0 ! -25 to + 25
924  zn0 = -10.0_rp ! -50 to + 50
925 
926  !write(ttt,'(i3.3)') t
927 
928  open(70,file='check_R1-R2_zn-base'//ttt//'.dat',status='unknown')
929  do j=1,10000
930  zn2 = zn0 + 9.0_rp*0.0001_rp*real((j-1),kind=rp)
931 
932  if (zn2/=0.0)then
933  ts_r1 = t0 - (c1*zn1*(t0-ts1) + c3*(zn1-zn2) - gflux*time)/(c1*zn2)
934 
935  ts_r2 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
936  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))
937  endif
938  if ( debug ) write(70,'(3f15.5)') zn2, ts_r1, ts_r2
939  enddo
940  do j=1,100000
941  zn2 = -1.0_rp + 2.0_rp*0.00001_rp*real((j-1),kind=rp)
942 
943  if (zn2/=0.0)then
944  ts_r1 = t0 - (c1*zn1*(t0-ts1) + c3*(zn1-zn2) - gflux*time)/(c1*zn2)
945 
946  ts_r2 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
947  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))
948  endif
949  if ( debug ) write(70,'(3f15.5)') zn2, ts_r1, ts_r2
950  enddo
951  do j=2,10000
952  zn2 = 1.0_rp + 9.0_rp*0.0001_rp*real((j-1),kind=rp)
953 
954  if (zn2/=0.0)then
955  ts_r1 = t0 - (c1*zn1*(t0-ts1) + c3*(zn1-zn2) - gflux*time)/(c1*zn2)
956 
957  ts_r2 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
958  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))
959  endif
960  if ( debug ) write(70,'(3f15.5)') zn2, ts_r1, ts_r2
961  enddo
962  close(70)
963 
964  a=t0+c3/c1
965  b=-1.0_rp * (c1*zn1*(t0-ts1) + c3*zn1 - gflux*time)
966  c=c1
967 
968  d=(lambdas*t0)
969  e=( ta*c2 - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) + epsilon*( lw - (sigma*(ta**4)))) ! *ZN2
970  f=(lambdas + (c2*zn2))
971 
972 
973  if ( debug ) then
974  open(71,file='check_R1-R2-grad'//ttt//'.dat',status='unknown')
975  write(71,'(5f15.5)') b/c, d,e,lambdas,c2
976  close(71)
977 
978  open(70,file='check_R1-R2_ts-base'//ttt//'.dat',status='unknown')
979  do i=1,100000
980  ts2 = ts0 + 150.0_rp*0.00001_rp*real((i-1),kind=rp)
981 
982  zn_r1 = (c1*zn1*(t0-ts1) + c3*zn1 - gflux*time)/(c1*(t0-ts2)+c3)
983  zn_r2 = -1.0_rp*lambdas*(t0-ts2)/(lw*epsilon-epsilon*sigma*(ta**4) - (c2*(ts2-ta)) - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat))
984 
985  write(70,'(3f15.5)') ts2, zn_r1, zn_r2
986  enddo
987  close(70)
988  end if
989 
990  !LOG_INFO("cal_R1R2",*) "aa",(LINFLUX-epsilon*sigma*(TA**4) + (C2*TA) - (LV*RHOAIR*CE*UA*(1.0_RP-RH)*QSAT))/C2
991  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
module CONSTANT
Definition: scale_const.F90:11
integer, parameter, public rp

Variable Documentation

◆ w0

real(rp), public scale_land_phy_snow_ky90::w0

Definition at line 35 of file scale_land_phy_snow_ky90.F90.

Referenced by cal_param(), calculationmo(), calculationnomo(), check_applicability(), groundflux(), land_phy_snow_ky90_setup(), snow_ky90_main(), and snowdepth().

35  real(RP), public :: w0 ! Maximum water content [ratio:0-1]

◆ rhosnow

real(rp), public scale_land_phy_snow_ky90::rhosnow = 400.0_RP

Definition at line 36 of file scale_land_phy_snow_ky90.F90.

Referenced by cal_param(), calculationmo(), calculationnomo(), check_applicability(), groundflux(), land_phy_snow_ky90_setup(), snow_ky90_main(), and snowdepth().

36  real(RP), public :: rhosnow = 400.0_rp ! Snow density [kg/m3]