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_water, SFLX_ENGI, PRSA, TA, QA, WA, UA, VA, DENS, SFLX_RAD_dn, exists_land, dt, TSNOW, SWE, SDepth, SDzero, nosnowsec, Salbedo, SFLX_SH, SFLX_LH, SFLX_QV, SFLX_QV_ENGI, 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, Evap_ENGI, QCC, QFUSION, MELT, SWEMELT, Gflux2land, SFLX_SNOW, SFLX_ENGI, 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
LAND_SNOW_MELT Heat used for snow melt J/m2 MELT
LAND_SNOW_QCC Heat used for changing temperature profile J/m2 QCC
LAND_SNOW_QFUSION Heat used for phase change of snow J/m2 QFUSION
LAND_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.

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

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

Referenced by mod_land_driver::land_driver_setup().

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_water,
real(rp), dimension (lia,lja), intent(in)  SFLX_ENGI,
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,
logical, dimension(lia,lja), intent(in)  exists_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_QV,
real(rp), dimension (lia,lja), intent(out)  SFLX_QV_ENGI,
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 151 of file scale_land_phy_snow_ky90.F90.

151  use scale_prc, only: &
152  prc_abort
153  use scale_file_history, only: &
154  file_history_in
155  use scale_atmos_saturation, only: &
156  qsatf => atmos_saturation_psat_all ! better to change name from qsatf to qsat
157  use scale_const, only: &
158  eps => const_eps, &
159  t0 => const_tem00, &
160  i_sw => const_i_sw, &
161  i_lw => const_i_lw, &
162  epsvap => const_epsvap
163  implicit none
164  integer, intent(in) :: LIA, LIS, LIE
165  integer, intent(in) :: LJA, LJS, LJE
166 
167  ! input data
168  real(RP), intent(in) :: SFLX_water(LIA,LJA)
169  real(RP), intent(in) :: SFLX_ENGI (LIA,LJA)
170  real(RP), intent(in) :: PRSA (LIA,LJA)
171  real(RP), intent(in) :: TA (LIA,LJA)
172  real(RP), intent(in) :: WA (LIA,LJA)
173  real(RP), intent(in) :: UA (LIA,LJA)
174  real(RP), intent(in) :: VA (LIA,LJA)
175  real(RP), intent(in) :: QA (LIA,LJA) ! specific humidity [kg/kg]
176  real(RP), intent(in) :: DENS (LIA,LJA)
177  real(RP), intent(in) :: SFLX_RAD_dn(LIA,LJA,N_RAD_DIR,N_RAD_RGN)
178  real(DP), intent(in) :: dt ! dt of land
179  logical, intent(in) :: exists_land(LIA,LJA)
180 
181  ! prognostic variables
182  real(RP), intent(inout) :: TSNOW (LIA,LJA) ! snow temperature [K]
183  real(RP), intent(inout) :: SWE (LIA,LJA) ! equivalent water [kg/m2]
184  real(RP), intent(inout) :: SDepth (LIA,LJA) ! depth of melting point [m]
185  real(RP), intent(inout) :: SDzero (LIA,LJA) ! total snow depth [m]
186  real(RP), intent(inout) :: nosnowsec (LIA,LJA) ! elapsed time of no snow condition [s]
187 
188  ! updated variables
189  real(RP), intent(out) :: Salbedo (LIA,LJA,2) ! snow albedo [-]
190  real(RP), intent(out) :: SFLX_SH (LIA,LJA) ! sensible heat flux between atmos and snow [W/m2]
191  real(RP), intent(out) :: SFLX_LH (LIA,LJA) ! latente heat flux between atmos and snow [W/m2]
192  real(RP), intent(out) :: SFLX_QV (LIA,LJA) ! evaporation due to LH [kg/m2/s]
193  real(RP), intent(out) :: SFLX_QV_ENGI (LIA,LJA) ! internal energy of evaporation [J/m2/s]
194  real(RP), intent(out) :: SFLX_GH (LIA,LJA) ! whole snowpack Ground flux [W/m2]
195  real(RP), intent(out) :: SNOW_LAND_GH (LIA,LJA) ! heat flux from snow to land [W/m2]
196  real(RP), intent(out) :: SNOW_LAND_water(LIA,LJA) ! water flux from snow to land [W/m2]
197  real(RP), intent(out) :: SNOW_frac (LIA,LJA) ! snow fraction, defined by time direction [-]
198 
199  real(RP) :: QCC (LIA,LJA)
200  real(RP) :: QFUSION (LIA,LJA)
201  real(RP) :: MELT (LIA,LJA)
202  real(RP) :: SWEMELT (LIA,LJA)
203 
204  ! works
205  real(RP) :: TSNOW1 ! updated snow surface temperature [K]
206  real(RP) :: ZNSNOW1 ! updated freezing depth [m]
207  real(RP) :: SWE1 ! updated snow water equivalence [kg/m2]
208  real(RP) :: DEPTH1 ! updated snow depth [m]
209 
210  real(RP), parameter :: Uabs_min = 0.1_rp
211  real(RP) :: Uabs
212  real(RP) :: QAsat
213  real(RP) :: RH
214  real(RP) :: qdry, psat
215  real(RP) :: SFLX_SW_dn, SFLX_LW_dn
216 
217  real(RP) :: w
218 
219  integer :: k, i, j
220  !---------------------------------------------------------------------------
221  log_progress(*) 'Land / physics / snow / KY90'
222 
223  do j = ljs, lje
224  do i = lis, lie
225 
226  if( ( exists_land(i,j) ) .and. &
227  ( swe(i,j)>0. .or. sflx_water(i,j)>0. ) )then
228 
229  uabs = sqrt( wa(i,j)**2 + ua(i,j)**2 + va(i,j)**2 )
230 
231  !qdry = 1.0_RP - QA(i,j)
232  !call qsatf( TA(i,j), PRSA(i,j), qdry, & ![IN]
233  ! QAsat ) ![OUT]
234  call qsatf( ta(i,j), psat )
235  qasat = epsvap * psat / ( prsa(i,j) - ( 1.0_rp-epsvap ) * psat )
236  rh = qa(i,j) / qasat
237  if ( debug ) then
238  log_info("LAND_PHY_SNOW_KY90",*) "RH, ",rh," DENS (1.289 org) ",dens(i,j)
239  end if
240 
241  sflx_lw_dn = sflx_rad_dn(i,j,i_r_diffuse,i_r_ir )
242  sflx_sw_dn = sflx_rad_dn(i,j,i_r_direct ,i_r_nir) &
243  + sflx_rad_dn(i,j,i_r_diffuse,i_r_nir) &
244  + sflx_rad_dn(i,j,i_r_direct ,i_r_vis) &
245  + sflx_rad_dn(i,j,i_r_diffuse,i_r_vis)
246 
247  tsnow1 = tsnow(i,j)
248  swe1 = swe(i,j)
249  depth1 = sdepth(i,j)
250  znsnow1 = sdzero(i,j)
251 
252  call snow_ky90_main( tsnow1, & ! [INOUT]
253  swe1, & ! [INOUT]
254  depth1, & ! [INOUT]
255  znsnow1, & ! [INOUT]
256  nosnowsec(i,j), & ! [INOUT]
257  salbedo(i,j,i_sw), & ! [OUT]
258  salbedo(i,j,i_lw), & ! [OUT]
259  sflx_sh(i,j), & ! [OUT]
260  sflx_lh(i,j), & ! [OUT]
261  sflx_gh(i,j), & ! [OUT]
262  sflx_qv(i,j), & ! [OUT]
263  sflx_qv_engi(i,j), & ! [OUT]
264  qcc(i,j), & ! [OUT]
265  qfusion(i,j), & ! [OUT]
266  melt(i,j), & ! [OUT]
267  swemelt(i,j), & ! [OUT]
268  snow_land_gh(i,j), & ! [OUT]
269  sflx_water(i,j), & ! [IN] ! [kg/m2/s]
270  sflx_engi(i,j), & ! [IN]
271  ta(i,j), & ! [IN]
272  uabs, & ! [IN]
273  rh, & ! [IN]
274  dens(i,j), & ! [IN]
275  sflx_sw_dn, & ! [IN]
276  sflx_lw_dn, & ! [IN]
277  dt )
278 
279 
280  snow_land_gh(i,j) = snow_land_gh(i,j) / dt ! [J/m2] -> [J/m2/s]
281  snow_land_water(i,j) = swemelt(i,j) / dt ! [kg/m2] -> [kg/m2/s]
282 
283  if ( swe1 <= 0. .and. swe(i,j) <= 0. ) then ! no accumulated snow during the time step
284  snow_frac(i,j) = 0.0_rp
285  else
286  snow_frac(i,j) = 1.0_rp
287  endif
288 
289  tsnow(i,j) = tsnow1
290  swe(i,j) = swe1
291  sdepth(i,j) = depth1
292  sdzero(i,j) = znsnow1
293 
294  if ( debug ) then
295  log_info("LAND_PHY_SNOW_KY90",*) "SNOW_frac, SWE, TSNOW", snow_frac(i,j), swe(i,j), tsnow(i,j)
296  end if
297 
298  else
299 
300  snow_land_water(i,j) = sflx_water(i,j)
301  snow_frac(i,j) = 0.0_rp
302 
303  tsnow(i,j) = t0 !!!
304  swe(i,j) = 0.0_rp !!!
305  sdepth(i,j) = 0.0_rp !!!
306  sdzero(i,j) = 0.0_rp !!!
307  salbedo(i,j,:) = 0.0_rp !!!
308  sflx_sh(i,j) = 0.0_rp
309  sflx_lh(i,j) = 0.0_rp
310  sflx_gh(i,j) = 0.0_rp
311  sflx_qv(i,j) = 0.0_rp
312  qcc(i,j) = 0.0_rp
313  qfusion(i,j) = 0.0_rp
314  melt(i,j) = 0.0_rp
315  swemelt(i,j) = 0.0_rp
316  snow_land_gh(i,j) = 0.0_rp
317 
318  endif
319 
320  call file_history_in( qcc(:,:), 'LAND_SNOW_QCC', 'Heat used for changing temperature profile', 'J/m2', dim_type='XY' )
321  call file_history_in( qfusion(:,:), 'LAND_SNOW_QFUSION', 'Heat used for phase change of snow', 'J/m2', dim_type='XY' )
322  call file_history_in( melt(:,:), 'LAND_SNOW_MELT', 'Heat used for snow melt', 'J/m2', dim_type='XY' )
323  call file_history_in( swemelt(:,:), 'LAND_SNOW_SWEMELT', 'Equivalent water of melt snow', 'kg/m2', dim_type='XY' )
324 
325  end do
326  end do
327 
328  return

References scale_const::const_eps, 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().

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)  Evap_ENGI,
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)  SFLX_ENGI,
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 361 of file scale_land_phy_snow_ky90.F90.

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

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

Referenced by land_phy_snow_ky90().

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 625 of file scale_land_phy_snow_ky90.F90.

625  implicit none
626 
627  real(RP), intent(in) :: TS, TA, UA, RH, rhoair, ALPHA, SW, LW
628  real(RP), intent(out) :: GFLUX, RFLUX, SFLUX, LINFLUX, LOUTFLUX, HFLUX, LATENTFLUX
629 
630  esat = 0.6112_rp * exp( (17.67_rp * (ta-273.15_rp)) / (ta-29.66_rp) )
631  qsat = 0.622_rp * esat / 101.325_rp
632  deltaqsat = 0.622_rp * lv * qsat /(287.04_rp * (ta**2))
633  cp = 1004.67_rp * (1.0_rp + (0.84_rp * qsat))
634 
635  ! Into snowpack is positive
636  sflux = (1.0_rp-alpha) * sw
637  rflux = epsilon * (lw-(sigma*(ts**4)))
638  linflux = epsilon * lw
639  loutflux = -1.0_rp * (epsilon * sigma * (ts**4))
640  hflux = cp * rhoair * ch * ua * (ts-ta)
641  latentflux = lv * rhoair * ce * ua * ((1-rh)*qsat + (deltaqsat*(ts-ta)))
642 
643  gflux = (sflux + rflux - hflux - latentflux)
644 
645  if ( debug ) then
646  log_info("LAND_PHY_SNOW_KY90_groundflux",*) "-------------- groundflux --------------"
647  log_info_cont(*) "GFLUX is: ", gflux
648  log_info_cont(*) "SFLUX: ", sflux
649  log_info_cont(*) "RFLUX: ", rflux, linflux+loutflux
650  log_info_cont(*) " (LONG in: ", linflux,")"
651  log_info_cont(*) " (LONG out: ", loutflux,")"
652  log_info_cont(*) "HFLUX is: ", hflux
653  log_info_cont(*) "LATENT FLUX: ", latentflux
654  end if
655 
656 return

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

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 693 of file scale_land_phy_snow_ky90.F90.

693  use scale_const, only: &
694  t0 => const_tem00
695 
696  implicit none
697  real(RP), intent(in) :: ZN1, TS1, TA, UA, RH, rhoair, LW
698  real(RP), intent(in) :: GFLUX
699  real(DP), intent(in) :: time
700 
701  c1 = csrhos*0.5_rp
702  c2 = (4.0_rp*epsilon*sigma*(ta**3)) + (cp*rhoair*ch*ua) + (lv*rhoair*ce*ua*deltaqsat)
703  c3 = w0*rhosnow*lf
704 
705  a0 = lambdas*(c3*zn1 + (c1*zn1*(t0-ts1)- gflux*time))
706  a1 = c2*(c3*zn1 + c1*zn1*(t0-ts1) - gflux*time) - c3*lambdas
707  a2 = c1* &
708  ( epsilon*(lw-(sigma*(ta**4))) + (c2*(ta-t0)) - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) ) &
709  - c2*c3
710 
711  if ( debug ) then
712  log_info("LAND_PHY_SNOW_KY90_cal_param",*) "-------------- snowdepth --------------"
713  log_info_cont(*) "C1",c1
714  log_info_cont(*) "C2",c2,(4.0_rp*epsilon*sigma*(ta**3)),(cp*rhoair*ch*ua), (lv*rhoair*ce*ua*deltaqsat)
715  log_info_cont(*) "C3",c3
716  log_info_cont(*) "A0",a0
717  log_info_cont(*) "A1",a1
718  log_info_cont(*) "A2",a2
719  end if
720 

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

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 725 of file scale_land_phy_snow_ky90.F90.

725  use scale_const, only: &
726  t0 => const_tem00
727 
728  implicit none
729  real(RP), intent(in) :: GFLUX, TS1, ZN1
730  real(RP), intent(in) :: TA, UA, RH, rhoair, LW
731  real(DP), intent(in) :: time
732  real(RP), intent(out) :: GFLUX_res, beta
733  real(RP) :: energy_in, energy_use_max
734 
735  energy_in = gflux * time
736  energy_use_max = 0.5_rp*csrhos*zn1*(t0-ts1) + w0*rhosnow*lf*zn1
737 
738  gflux_res = energy_in - energy_use_max ! residual energy after being used to melt all snow
739 
740  !
741  beta = (lambdas/c2)* &
742  ( t0 - ( ta*c2 - lv*rhoair*ce*ua*(1.0_rp-rh)*qsat + epsilon*(lw-sigma*(ta**4)))/c2 )
743 
744  return

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

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 749 of file scale_land_phy_snow_ky90.F90.

749 
750  implicit none
751 
752  real(RP), intent(in) :: GFLUX
753  real(RP), intent(in) :: ZN1
754  real(DP), intent(in) :: time
755  real(RP), intent(out) :: ZN2
756 
757 !calcultaing snowdepth
758 
759  !print *, -1.0_RP*A1
760  !print *, -1.0_RP*((A1**2.0_RP - 4.0_RP * A2 * A0)**0.5_RP)
761  !print *, 2.0_RP*A2
762  !print *, A1**2.0_RP, -1.0_RP * 4.0_RP * A2 * A0
763 
764  zn2 = ((-1.0_rp*a1) - ((a1**2.0_rp - 4.0_rp*a2*a0)**0.5_rp)) / (2.0_rp*a2)
765 
766  if ( debug ) then
767  log_info("LAND_PHY_SNOW_KY90_snowdepth",*) "ZN old = ",zn1, "ZN new = ",zn2
768  end if
769 
770 return

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

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 800 of file scale_land_phy_snow_ky90.F90.

800  use scale_const, only: &
801  t0 => const_tem00
802  implicit none
803 
804  real(RP), intent(in) :: ZN1
805  real(RP), intent(in) :: TS, GFLUX
806  real(DP), intent(in) :: time
807  real(RP), intent(out) :: ZN2
808 
809  zn2 = zn1 + ((c1*zn1*(t0-ts) - gflux*time) / c3)
810 
811  if ( debug ) then
812  log_info("LAND_PHY_SNOW_KY90_recalculateZ",*) "-------------- recalculate Z --------------"
813  end if
814 
815  return

References scale_const::const_tem00.

Referenced by snow_ky90_main().

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 821 of file scale_land_phy_snow_ky90.F90.

821  use scale_const, only: &
822  t0 => const_tem00
823  use scale_prc, only: &
824  prc_abort
825  implicit none
826  real(RP), intent(in) :: GFLUX, CSRHOS, ZN1, TS1, TS2, ZN2
827  real(DP), intent(in) :: time
828  real(RP), intent(out) :: MELT, QCC, QFUSION
829 
830  qcc = 0.5_rp*csrhos*(zn1*(t0-ts1)-zn2*(t0-ts2))
831  qfusion = w0*rhosnow*lf*(zn1-zn2)
832  melt = ( gflux*time - qcc - qfusion )
833 
834  if ( debug ) then
835  log_info("LAND_PHY_SNOW_KY90_calculationMO",*) "--------------------------MELT----------------"
836  log_info_cont(*) "GFLUX*time is: ", gflux*time
837  log_info_cont(*) "QCC is : ", qcc
838  log_info_cont(*) "QFUSION is : ", qfusion
839  log_info_cont(*) "QMELT is : ", melt
840 
841  log_info_cont(*) qcc+qfusion+melt
842  log_info_cont(*) "diff= ", qcc + qfusion + melt - (gflux*time)
843  end if
844 
845  if ( abs(qcc+qfusion+melt - (gflux*time)) > 10.) then
846  log_error("LAND_PHY_SNOW_KY90_calculationMO",*) "Calculation is fault. Model would include bugs. Please check! Melt"
847  call prc_abort
848  endif
849 
850 
851  return

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

Referenced by snow_ky90_main().

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 857 of file scale_land_phy_snow_ky90.F90.

857  use scale_const, only: &
858  t0 => const_tem00
859  implicit none
860  real(RP), intent(in) :: GFLUX, CSRHOS, ZN1, TS1, TS2, ZN2
861  real(DP), intent(in) :: time
862  real(RP), intent(out) :: MELT, QCC, QFUSION
863 
864  qcc = 0.5_rp*csrhos*(zn1*(t0-ts1)-zn2*(t0-ts2))
865  qfusion = w0*rhosnow*lf*(zn1-zn2)
866  melt = 0.0_rp
867 
868  if ( debug ) then
869  log_info("LAND_PHY_SNOW_KY90_calculationNoMO",*) "--------------------------NOMELT----------------"
870  log_info_cont(*) "GFLUX*time is: ", gflux*time
871  log_info_cont(*) "QCC is : ", qcc
872  log_info_cont(*) "QFUSION is : ", qfusion
873  log_info_cont(*) "QMELT is : ", melt
874 
875  log_info_cont(*) qcc+qfusion+melt
876  log_info_cont(*) "diff= ", qcc +qfusion - (gflux*time)
877  end if
878 
879  !if ( ABS(QCC+QFUSION - (GFLUX*time)) > 10.) then
880  ! LOG_ERROR("LAND_PHY_SNOW_KY90_calculationNoMO",*) "Calculation is fault. Model would include bugs. Please check! No Melt"
881  ! call PRC_abort
882  !endif
883 
884  return

References scale_const::const_tem00, rhosnow, and w0.

Referenced by snow_ky90_main().

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 889 of file scale_land_phy_snow_ky90.F90.

889  use scale_const, only: &
890  t0 => const_tem00
891  implicit none
892  real(RP), intent(in) :: ZN1, ZN2, TS1, TS2, GFLUX
893  real(RP), intent(in) :: TA, UA, RH, rhoair, LW
894  real(DP), intent(in) :: time
895  character(len=*) :: flag
896  real(RP) :: R1 ! Eq. (2)
897  real(RP) :: R2 ! Eq. (8)
898  real(RP) :: R3 ! Eq. (8)
899 
900  r3 = -999.
901  r1 = c1 * (zn1*(t0-ts1) - zn2*(t0-ts2)) + c3*(zn1-zn2)-gflux*time
902  r2 = zn2*( epsilon*lw-epsilon*sigma*(ta**4) - (c2*(ts2-ta)) - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) ) &
903  + lambdas*(t0-ts2)
904  r3 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
905  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))-ts2
906 
907 
908  if ( debug ) then
909  log_info("LAND_PHY_SNOW_KY90_check_res",*) "R1 is : ", r1, "flag = ", flag
910  log_info("LAND_PHY_SNOW_KY90_check_res",*) "R2 is : ", r2, r3,"flag = ", flag
911 
912  if(abs(r1)>10000)then
913  log_info("LAND_PHY_SNOW_KY90_check_res",*) c1 * (zn1*(t0-ts1) - zn2*(t0-ts2)), c3*(zn1-zn2), -gflux*time
914  endif
915  end if
916 
917  return

References scale_const::const_tem00.

Referenced by snow_ky90_main().

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 922 of file scale_land_phy_snow_ky90.F90.

922  use scale_const, only: &
923  t0 => const_tem00
924  implicit none
925  real(RP), intent(in) :: ZN1, TS1
926  real(RP), intent(in) :: GFLUX, TA, UA, RH, rhoair, LW
927  real(DP), intent(in) :: time
928  real(RP) :: ZN2, TS2
929  real(RP) :: R1 ! Eq. (2)
930  real(RP) :: R2 ! Eq. (8)
931  real(RP) :: R3 ! Eq. (8)
932 
933  real(RP) :: ts0,zn0
934  integer :: i,j
935  real(RP) :: ts_r1,ts_r2,zn_r1,zn_r2
936  character(len=3) :: ttt = ""
937 
938  real(RP) :: a,b,c,d,e,f,g
939 
940  ts0 = -50.0_rp + t0 ! -25 to + 25
941  zn0 = -10.0_rp ! -50 to + 50
942 
943  !write(ttt,'(i3.3)') t
944 
945  open(70,file='check_R1-R2_zn-base'//ttt//'.dat',status='unknown')
946  do j=1,10000
947  zn2 = zn0 + 9.0_rp*0.0001_rp*real((j-1),kind=rp)
948 
949  if (zn2/=0.0)then
950  ts_r1 = t0 - (c1*zn1*(t0-ts1) + c3*(zn1-zn2) - gflux*time)/(c1*zn2)
951 
952  ts_r2 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
953  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))
954  endif
955  if ( debug ) write(70,'(3f15.5)') zn2, ts_r1, ts_r2
956  enddo
957  do j=1,100000
958  zn2 = -1.0_rp + 2.0_rp*0.00001_rp*real((j-1),kind=rp)
959 
960  if (zn2/=0.0)then
961  ts_r1 = t0 - (c1*zn1*(t0-ts1) + c3*(zn1-zn2) - gflux*time)/(c1*zn2)
962 
963  ts_r2 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
964  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))
965  endif
966  if ( debug ) write(70,'(3f15.5)') zn2, ts_r1, ts_r2
967  enddo
968  do j=2,10000
969  zn2 = 1.0_rp + 9.0_rp*0.0001_rp*real((j-1),kind=rp)
970 
971  if (zn2/=0.0)then
972  ts_r1 = t0 - (c1*zn1*(t0-ts1) + c3*(zn1-zn2) - gflux*time)/(c1*zn2)
973 
974  ts_r2 = ((lambdas*t0) + (ta*c2*zn2) - (zn2*lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) &
975  + zn2*epsilon*( lw - (sigma*(ta**4))))/ (lambdas + (c2*zn2))
976  endif
977  if ( debug ) write(70,'(3f15.5)') zn2, ts_r1, ts_r2
978  enddo
979  close(70)
980 
981  a=t0+c3/c1
982  b=-1.0_rp * (c1*zn1*(t0-ts1) + c3*zn1 - gflux*time)
983  c=c1
984 
985  d=(lambdas*t0)
986  e=( ta*c2 - (lv*rhoair*ce*ua*(1.0_rp-rh)*qsat) + epsilon*( lw - (sigma*(ta**4)))) ! *ZN2
987  f=(lambdas + (c2*zn2))
988 
989 
990  if ( debug ) then
991  open(71,file='check_R1-R2-grad'//ttt//'.dat',status='unknown')
992  write(71,'(5f15.5)') b/c, d,e,lambdas,c2
993  close(71)
994 
995  open(70,file='check_R1-R2_ts-base'//ttt//'.dat',status='unknown')
996  do i=1,100000
997  ts2 = ts0 + 150.0_rp*0.00001_rp*real((i-1),kind=rp)
998 
999  zn_r1 = (c1*zn1*(t0-ts1) + c3*zn1 - gflux*time)/(c1*(t0-ts2)+c3)
1000  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))
1001 
1002  write(70,'(3f15.5)') ts2, zn_r1, zn_r2
1003  enddo
1004  close(70)
1005  end if
1006 
1007  !LOG_INFO("cal_R1R2",*) "aa",(LINFLUX-epsilon*sigma*(TA**4) + (C2*TA) - (LV*RHOAIR*CE*UA*(1.0_RP-RH)*QSAT))/C2
1008  return

References scale_const::const_tem00.

Variable Documentation

◆ w0

real(rp), public scale_land_phy_snow_ky90::w0

Definition at line 35 of file scale_land_phy_snow_ky90.F90.

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

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

◆ 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.

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

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

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:128
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:134