SCALE-RM
Functions/Subroutines | Variables
scale_ocean_phy_ice_simple Module Reference

module ocean / physics / ice / simple More...

Functions/Subroutines

subroutine, public ocean_phy_ice_setup
 
subroutine, public ocean_phy_ice_fraction (OIA, OIS, OIE, OJA, OJS, OJE, ICE_MASS, ICE_FRAC)
 
subroutine, public ocean_phy_ice_adjustment (OIA, OIS, OIE, OJA, OJS, OJE, calc_flag, OCEAN_DEPTH, OCEAN_TEMP, ICE_TEMP, ICE_MASS, MASS_FLUX, ENGI_FLUX, MASS_SUPL, ENGI_SUPL)
 
subroutine, public ocean_phy_ice_simple (OIA, OIS, OIE, OJA, OJS, OJE, iflx_water, iflx_hbalance, subsfc_temp, TC_dz, ICE_TEMP, ICE_MASS, ICE_FRAC, calc_flag, dt, ICE_TEMP_t, ICE_MASS_t, sflx_G, sflx_water, sflx_RHOE)
 Slab ocean model. More...
 

Variables

real(rp), public ocean_phy_ice_freezetemp = 271.35_RP
 
real(rp), public ocean_phy_ice_density = 1000.0_RP
 

Detailed Description

module ocean / physics / ice / simple

Description
ocean common module
Author
Team SCALE
NAMELIST
  • PARAM_OCEAN_PHY_ICE
    nametypedefault valuecomment
    OCEAN_PHY_ICE_DENSITY real(RP) 1000.0_RP density of sea ice [kg/m3]
    OCEAN_PHY_ICE_MASS_CRITICAL real(RP) 1600.0_RP ice amount for fraction = 1 [kg/m2]
    OCEAN_PHY_ICE_MASS_LIMIT real(RP) 50000.0_RP maximum ice amount [kg/m2]
    OCEAN_PHY_ICE_FRACTION_LIMIT real(RP) 1.0_RP maximum ice fraction [1]
    OCEAN_PHY_ICE_DT_MAX real(RP) 5.E-2_RP maximum delta ice temperature [K/s]

History Output
No history output

Function/Subroutine Documentation

◆ ocean_phy_ice_setup()

subroutine, public scale_ocean_phy_ice_simple::ocean_phy_ice_setup

Definition at line 60 of file scale_ocean_phy_ice_simple.F90.

60  use scale_prc, only: &
61  prc_abort
62  use scale_const, only: &
63  undef => const_undef
64  use scale_calendar, only: &
66  use scale_file_external_input, only: &
67  file_external_input_regist
68  implicit none
69 
70  real(DP) :: OCEAN_PHY_ICE_nudging_tau = 0.0_dp ! Relaxation time
71  character(len=H_SHORT) :: OCEAN_PHY_ICE_nudging_tau_unit = "SEC"
72  character(len=H_LONG) :: OCEAN_PHY_ICE_nudging_basename = ''
73  logical :: OCEAN_PHY_ICE_nudging_basename_add_num = .false.
74  integer :: OCEAN_PHY_ICE_nudging_number_of_files = 1
75  logical :: OCEAN_PHY_ICE_nudging_enable_periodic_year = .false.
76  logical :: OCEAN_PHY_ICE_nudging_enable_periodic_month = .false.
77  logical :: OCEAN_PHY_ICE_nudging_enable_periodic_day = .false.
78  integer :: OCEAN_PHY_ICE_nudging_step_fixed = 0
79  real(RP) :: OCEAN_PHY_ICE_nudging_defval != UNDEF
80  logical :: OCEAN_PHY_ICE_nudging_check_coordinates = .true.
81  integer :: OCEAN_PHY_ICE_nudging_step_limit = 0
82 
83  namelist / param_ocean_phy_ice / &
84  ocean_phy_ice_density, &
85  ocean_phy_ice_mass_critical, &
86  ocean_phy_ice_mass_limit, &
87  ocean_phy_ice_fraction_limit, &
88  ocean_phy_ice_dt_max!, &
89 ! OCEAN_PHY_ICE_nudging, &
90 ! OCEAN_PHY_ICE_nudging_tau, &
91 ! OCEAN_PHY_ICE_nudging_tau_unit, &
92 ! OCEAN_PHY_ICE_nudging_basename, &
93 ! OCEAN_PHY_ICE_nudging_basename_add_num, &
94 ! OCEAN_PHY_ICE_nudging_number_of_files, &
95 ! OCEAN_PHY_ICE_nudging_enable_periodic_year, &
96 ! OCEAN_PHY_ICE_nudging_enable_periodic_month, &
97 ! OCEAN_PHY_ICE_nudging_enable_periodic_day, &
98 ! OCEAN_PHY_ICE_nudging_step_fixed, &
99 ! OCEAN_PHY_ICE_nudging_defval, &
100 ! OCEAN_PHY_ICE_nudging_check_coordinates, &
101 ! OCEAN_PHY_ICE_nudging_step_limit
102 
103  integer :: ierr
104  !---------------------------------------------------------------------------
105 
106  log_newline
107  log_info("OCEAN_PHY_ICE_setup",*) 'Setup'
108 
109  ocean_phy_ice_nudging_defval = undef
110 
111  !--- read namelist
112  rewind(io_fid_conf)
113  read(io_fid_conf,nml=param_ocean_phy_ice,iostat=ierr)
114  if( ierr < 0 ) then !--- missing
115  log_info("OCEAN_PHY_ICE_setup",*) 'Not found namelist. Default used.'
116  elseif( ierr > 0 ) then !--- fatal error
117  log_error("OCEAN_PHY_ICE_setup",*) 'Not appropriate names in namelist PARAM_OCEAN_PHY_ICE. Check!'
118  call prc_abort
119  endif
120  log_nml(param_ocean_phy_ice)
121 
122  log_newline
123  log_info("OCEAN_PHY_ICE_setup",*) 'Ice amount for frac. = 1 [kg/m2] : ', ocean_phy_ice_mass_critical
124  if ( ocean_phy_ice_nudging ) then
125  call calendar_unit2sec( ocean_phy_ice_nudging_tausec, ocean_phy_ice_nudging_tau, ocean_phy_ice_nudging_tau_unit )
126 
127  log_info("OCEAN_PHY_ICE_setup",*) 'Use nudging for sea ice fraction : ON'
128  log_info("OCEAN_PHY_ICE_setup",*) 'Relaxation time Tau [sec] : ', ocean_phy_ice_nudging_tausec
129 
130  if ( ocean_phy_ice_nudging_tausec == 0.0_rp ) then
131  ocean_phy_ice_offline_mode = .true.
132  log_info("OCEAN_PHY_ICE_setup",*) 'Tau=0 means that sea ice is completely replaced by the external data.'
133  endif
134 
135  if ( ocean_phy_ice_nudging_basename == '' ) then
136  log_error("OCEAN_PHY_ICE_setup",*) 'OCEAN_PHY_ICE_nudging_basename is necessary. STOP'
137  call prc_abort
138  endif
139  else
140  log_info("OCEAN_PHY_ICE_setup",*) 'Use nudging for sea ice fraction : OFF'
141  endif
142 
143  if ( ocean_phy_ice_nudging ) then
144  call file_external_input_regist( ocean_phy_ice_nudging_basename, & ! [IN]
145  ocean_phy_ice_nudging_basename_add_num, & ! [IN]
146  ocean_phy_ice_nudging_number_of_files, & ! [IN]
147  'OCEAN_ICE_FRAC', & ! [IN]
148  'XY', & ! [IN]
149  ocean_phy_ice_nudging_enable_periodic_year, & ! [IN]
150  ocean_phy_ice_nudging_enable_periodic_month, & ! [IN]
151  ocean_phy_ice_nudging_enable_periodic_day, & ! [IN]
152  ocean_phy_ice_nudging_step_fixed, & ! [IN]
153  ocean_phy_ice_nudging_defval, & ! [IN]
154  check_coordinates = ocean_phy_ice_nudging_check_coordinates, & ! [IN]
155  step_limit = ocean_phy_ice_nudging_step_limit ) ! [IN]
156  endif
157 
158  return

References scale_calendar::calendar_unit2sec(), scale_const::const_undef, scale_io::io_fid_conf, ocean_phy_ice_density, and scale_prc::prc_abort().

Referenced by mod_ocean_driver::ocean_driver_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ocean_phy_ice_fraction()

subroutine, public scale_ocean_phy_ice_simple::ocean_phy_ice_fraction ( integer, intent(in)  OIA,
integer, intent(in)  OIS,
integer, intent(in)  OIE,
integer, intent(in)  OJA,
integer, intent(in)  OJS,
integer, intent(in)  OJE,
real(rp), dimension(oia,oja), intent(in)  ICE_MASS,
real(rp), dimension(oia,oja), intent(out)  ICE_FRAC 
)

Definition at line 167 of file scale_ocean_phy_ice_simple.F90.

167  implicit none
168 
169  integer, intent(in) :: OIA, OIS, OIE
170  integer, intent(in) :: OJA, OJS, OJE
171  real(RP), intent(in) :: ICE_MASS(OIA,OJA) ! sea ice amount [kg/m2]
172  real(RP), intent(out) :: ICE_FRAC(OIA,OJA) ! sea ice area fraction [1]
173 
174  integer :: i, j
175  !---------------------------------------------------------------------------
176  !$acc data copyin(ICE_MASS) copyout(ICE_FRAC)
177 
178  !$acc kernels
179  do j = ojs, oje
180  do i = ois, oie
181  ice_frac(i,j) = ice_mass(i,j) / ocean_phy_ice_mass_critical
182 
183  ice_frac(i,j) = min( sqrt( max( ice_frac(i,j), 0.0_rp ) ), ocean_phy_ice_fraction_limit )
184  enddo
185  enddo
186  !$acc end kernels
187 
188  !$acc end data
189  return

Referenced by mod_ocean_driver::ocean_driver_update(), and mod_ocean_vars::ocean_vars_restart_read().

Here is the caller graph for this function:

◆ ocean_phy_ice_adjustment()

subroutine, public scale_ocean_phy_ice_simple::ocean_phy_ice_adjustment ( integer, intent(in)  OIA,
integer, intent(in)  OIS,
integer, intent(in)  OIE,
integer, intent(in)  OJA,
integer, intent(in)  OJS,
integer, intent(in)  OJE,
logical, dimension (oia,oja), intent(in)  calc_flag,
real(rp), intent(in)  OCEAN_DEPTH,
real(rp), dimension(oia,oja), intent(inout)  OCEAN_TEMP,
real(rp), dimension (oia,oja), intent(inout)  ICE_TEMP,
real(rp), dimension (oia,oja), intent(inout)  ICE_MASS,
real(rp), dimension(oia,oja), intent(out)  MASS_FLUX,
real(rp), dimension(oia,oja), intent(out)  ENGI_FLUX,
real(rp), dimension (oia,oja), intent(out)  MASS_SUPL,
real(rp), dimension (oia,oja), intent(out)  ENGI_SUPL 
)

Definition at line 205 of file scale_ocean_phy_ice_simple.F90.

205  use scale_const, only: &
206  dwatr => const_dwatr, &
207  eps => const_eps
208  use scale_atmos_hydrometeor, only: &
209  cv_water, &
210  cv_ice, &
211  lhf
212  implicit none
213 
214  integer, intent(in) :: OIA, OIS, OIE
215  integer, intent(in) :: OJA, OJS, OJE
216  logical, intent(in) :: calc_flag (OIA,OJA) ! to decide calculate or not
217  real(RP), intent(in) :: OCEAN_DEPTH ! depth of the first layer of the ocean
218  real(RP), intent(inout) :: OCEAN_TEMP(OIA,OJA) ! ocean temperature [K]
219  real(RP), intent(inout) :: ICE_TEMP (OIA,OJA) ! sea ice temperature [K]
220  real(RP), intent(inout) :: ICE_MASS (OIA,OJA) ! sea ice amount [kg/m2]
221  real(RP), intent(out) :: MASS_FLUX(OIA,OJA)
222  real(RP), intent(out) :: ENGI_FLUX(OIA,OJA)
223  real(RP), intent(out) :: MASS_SUPL (OIA,OJA)
224  real(RP), intent(out) :: ENGI_SUPL (OIA,OJA)
225 
226  real(RP) :: C_w
227  real(RP) :: ICE_MASS_frz
228  real(RP) :: ICE_MASS_prev
229 
230  integer :: i, j
231  !---------------------------------------------------------------------------
232  !$acc data copy(OCEAN_TEMP,ICE_TEMP,ICE_MASS) &
233  !$acc copyin(calc_flag) &
234  !$acc copyout(MASS_FLUX,ENGI_FLUX,MASS_SUPL,ENGI_SUPL)
235 
236  c_w = cv_water * dwatr * ocean_depth
237 
238  !$acc kernels
239  !$omp parallel do &
240  !$omp private(ICE_MASS_frz,ICE_MASS_prev)
241  do j = ojs, oje
242  do i = ois, oie
243  if ( calc_flag(i,j) &
244  .and. ocean_temp(i,j) < ocean_phy_ice_freezetemp .and. ice_mass(i,j) < ocean_phy_ice_mass_limit ) then
245  ice_mass_frz = c_w * ( ocean_phy_ice_freezetemp - ocean_temp(i,j) ) &
247  ice_mass_frz = min( ice_mass_frz, dwatr * ocean_depth )
248 
249  ! update ice mass
250  ice_mass_prev = ice_mass(i,j)
251  ice_mass(i,j) = ice_mass(i,j) + ice_mass_frz
252  ice_mass(i,j) = min( ice_mass(i,j), ocean_phy_ice_mass_limit ) ! apply limiter
253  ice_mass_frz = ice_mass(i,j) - ice_mass_prev
254 
255  ! update ice temperature
256  if ( ice_mass(i,j) > eps ) then
257  ice_temp(i,j) = ice_temp(i,j) &
258  + ( ocean_phy_ice_freezetemp - ice_temp(i,j) ) * ice_mass_frz / ice_mass(i,j)
259  else ! all ice melt
260  ice_temp(i,j) = ocean_phy_ice_freezetemp ! dummy
261  end if
262 
263  ! update ocean temperature
264  if ( c_w - cv_water * ice_mass_frz > eps ) then
265  ocean_temp(i,j) = ocean_temp(i,j) &
266  + ( cv_water * ocean_temp(i,j) - cv_ice * ocean_phy_ice_freezetemp + lhf ) * ice_mass_frz &
267  / ( c_w - cv_water * ice_mass_frz )
268  else ! all water freeze
269  ocean_temp(i,j) = ocean_phy_ice_freezetemp
270  end if
271 
272  mass_flux(i,j) = ice_mass_frz
273  engi_flux(i,j) = ( cv_ice * ocean_phy_ice_freezetemp - lhf ) * ice_mass_frz
274  mass_supl(i,j) = ice_mass_frz
275  engi_supl(i,j) = ice_mass_frz * cv_water * ocean_temp(i,j)
276  else
277  mass_flux(i,j) = 0.0_rp
278  engi_flux(i,j) = 0.0_rp
279  mass_supl(i,j) = 0.0_rp
280  engi_supl(i,j) = 0.0_rp
281  endif
282  enddo
283  enddo
284  !$acc end kernels
285 
286  !$acc end data
287  return

References scale_const::const_dwatr, scale_const::const_eps, scale_atmos_hydrometeor::cv_ice, scale_atmos_hydrometeor::cv_water, scale_atmos_hydrometeor::lhf, and ocean_phy_ice_freezetemp.

Referenced by mod_ocean_driver::ocean_driver_update().

Here is the caller graph for this function:

◆ ocean_phy_ice_simple()

subroutine, public scale_ocean_phy_ice_simple::ocean_phy_ice_simple ( integer, intent(in)  OIA,
integer, intent(in)  OIS,
integer, intent(in)  OIE,
integer, intent(in)  OJA,
integer, intent(in)  OJS,
integer, intent(in)  OJE,
real(rp), dimension (oia,oja), intent(in)  iflx_water,
real(rp), dimension(oia,oja), intent(in)  iflx_hbalance,
real(rp), dimension (oia,oja), intent(in)  subsfc_temp,
real(rp), dimension (oia,oja), intent(in)  TC_dz,
real(rp), dimension (oia,oja), intent(in)  ICE_TEMP,
real(rp), dimension (oia,oja), intent(in)  ICE_MASS,
real(rp), dimension (oia,oja), intent(in)  ICE_FRAC,
logical, dimension (oia,oja), intent(in)  calc_flag,
real(dp), intent(in)  dt,
real(rp), dimension (oia,oja), intent(out)  ICE_TEMP_t,
real(rp), dimension (oia,oja), intent(out)  ICE_MASS_t,
real(rp), dimension (oia,oja), intent(out)  sflx_G,
real(rp), dimension (oia,oja), intent(out)  sflx_water,
real(rp), dimension (oia,oja), intent(out)  sflx_RHOE 
)

Slab ocean model.

Definition at line 309 of file scale_ocean_phy_ice_simple.F90.

309  use scale_prc, only: &
310  prc_abort
311  use scale_atmos_hydrometeor, only: &
312  cv_water, &
313  cv_ice, &
314  lhf
315  use scale_file_external_input, only: &
316  file_external_input_update
317  implicit none
318 
319  integer, intent(in) :: OIA, OIS, OIE
320  integer, intent(in) :: OJA, OJS, OJE
321  real(RP), intent(in) :: iflx_water (OIA,OJA) ! input mass flux [kg/m2/s] (downward)
322  real(RP), intent(in) :: iflx_hbalance(OIA,OJA) ! input heat flux [J/m2/s] (downward)
323  real(RP), intent(in) :: subsfc_temp (OIA,OJA) ! subsurface temperature [K]
324  real(RP), intent(in) :: TC_dz (OIA,OJA) ! Thermal conductance [K/m]
325  real(RP), intent(in) :: ICE_TEMP (OIA,OJA) ! sea ice temperature [K]
326  real(RP), intent(in) :: ICE_MASS (OIA,OJA) ! sea ice amount [kg/m2]
327  real(RP), intent(in) :: ICE_FRAC (OIA,OJA) ! sea ice fraction [0-1]
328  logical, intent(in) :: calc_flag (OIA,OJA) ! to decide calculate or not
329  real(DP), intent(in) :: dt
330  real(RP), intent(out) :: ICE_TEMP_t (OIA,OJA) ! tendency of sea ice temperature [K/s]
331  real(RP), intent(out) :: ICE_MASS_t (OIA,OJA) ! tendency of sea ice amount [kg/m2/s]
332  real(RP), intent(out) :: SFLX_G (OIA,OJA) ! heat flux from sea ice to subsurface
333  real(RP), intent(out) :: SFLX_water (OIA,OJA) ! mass flux from sea ice to subsurface
334  real(RP), intent(out) :: SFLX_RHOE (OIA,OJA) ! internal energy flux from sea ice to subsurface
335 
336  real(RP) :: ICE_MASS_new ! [kg/m2]
337  real(RP) :: ICE_TEMP_new ! [K]
338  real(RP) :: mass_budget ! [kg/m2/s]
339  real(RP) :: heat_budget ! [J/m2/s]
340  real(RP) :: G ! [J/m2/s]
341  real(RP) :: dM ! [kg/m2]
342  real(RP) :: dE ! [J/m2]
343  real(RP) :: M_mlt ! [kg/m2]
344  real(RP) :: dt_RP
345 
346  integer :: i, j
347  !---------------------------------------------------------------------------
348  !$acc data copyin (iflx_water,iflx_hbalance,subsfc_temp, &
349  !$acc TC_dz,ICE_TEMP,ICE_MASS,ICE_FRAC,calc_flag) &
350  !$acc copyout(ICE_TEMP_t,ICE_MASS_t,SFLX_G,SFLX_water,SFLX_RHOE)
351 
352  log_progress(*) 'ocean / physics / seaice'
353 
354  dt_rp = real(dt,kind=rp)
355 
356  !$acc kernels
357  !$omp parallel do &
358  !$omp private(mass_budget,heat_budget,dM,dE,G,M_mlt, &
359  !$omp ICE_TEMP_new,ICE_MASS_new)
360  do j = ojs, oje
361  do i = ois, oie
362  if ( calc_flag(i,j) ) then
363 
364  ! mass change
365  dm = iflx_water(i,j) * ice_frac(i,j) * dt_rp
366  ice_mass_new = ice_mass(i,j) + dm
367 
368  if ( ice_mass_new > 0.0_rp ) then
369  ! internal energy change
370  g = ( subsfc_temp(i,j) - ice_temp(i,j) ) * tc_dz(i,j) ! heat flux from ocean
371  de = ( iflx_hbalance(i,j) + g ) * ice_frac(i,j) * dt_rp
372  ice_temp_new = ice_temp(i,j) &
373  + ( de - ( cv_ice * ice_temp(i,j) - lhf ) * dm ) / ( cv_ice * ice_mass_new )
374 
375  ! melting ice
376  m_mlt = cv_ice * ( ice_temp_new - ocean_phy_ice_freezetemp ) * ice_mass_new &
378  m_mlt = min( max( m_mlt, 0.0_rp ), ice_mass_new )
379 
380  ice_mass_new = ice_mass_new - m_mlt
381  ice_temp_new = ice_temp_new &
382  + ( cv_ice * ice_temp_new - lhf - cv_water * ocean_phy_ice_freezetemp ) * m_mlt &
383  / ( cv_ice * ice_mass_new )
384 
385  ! ice to ocean flux
386  mass_budget = m_mlt / dt_rp
387  sflx_rhoe(i,j) = cv_water * ocean_phy_ice_freezetemp * mass_budget
388  sflx_g(i,j) = - g * ice_frac(i,j)
389  sflx_water(i,j) = mass_budget
390 
391  else
392 
393  ice_mass_new = 0.0_rp
394  ice_temp_new = ocean_phy_ice_freezetemp ! dummy
395 
396  sflx_rhoe(i,j) = cv_water * subsfc_temp(i,j) * ice_mass_new / dt_rp
397  sflx_g(i,j) = ( cv_ice * ice_temp(i,j) - lhf ) * ice_mass(i,j) &
398  - sflx_rhoe(i,j) + iflx_hbalance(i,j) * ice_frac(i,j)
399  sflx_water(i,j) = ice_mass_new ! (negative)
400 
401  endif
402 
403  ice_mass_t(i,j) = ( ice_mass_new - ice_mass(i,j) ) / dt_rp
404  ice_temp_t(i,j) = ( ice_temp_new - ice_temp(i,j) ) / dt_rp
405 
406  else
407  sflx_g(i,j) = 0.0_rp
408  sflx_water(i,j) = 0.0_rp
409  sflx_rhoe(i,j) = 0.0_rp
410  ice_mass_t(i,j) = 0.0_rp
411  ice_temp_t(i,j) = 0.0_rp
412  endif
413  enddo
414  enddo
415  !$acc end kernels
416 
417  !$acc end data
418  return

References scale_atmos_hydrometeor::cv_ice, scale_atmos_hydrometeor::cv_water, scale_atmos_hydrometeor::lhf, ocean_phy_ice_freezetemp, scale_prc::prc_abort(), and scale_precision::rp.

Referenced by mod_ocean_driver::ocean_driver_calc_tendency().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ ocean_phy_ice_freezetemp

real(rp), public scale_ocean_phy_ice_simple::ocean_phy_ice_freezetemp = 271.35_RP

Definition at line 44 of file scale_ocean_phy_ice_simple.F90.

44  real(RP), public :: OCEAN_PHY_ICE_freezetemp = 271.35_rp ! freezing temperature of sea ice [K]

Referenced by ocean_phy_ice_adjustment(), ocean_phy_ice_simple(), mod_ocean_vars::ocean_vars_restart_read(), and mod_realinput::realinput_surface().

◆ ocean_phy_ice_density

real(rp), public scale_ocean_phy_ice_simple::ocean_phy_ice_density = 1000.0_RP

Definition at line 45 of file scale_ocean_phy_ice_simple.F90.

45  real(RP), public :: OCEAN_PHY_ICE_density = 1000.0_rp ! density of sea ice [kg/m3]

Referenced by ocean_phy_ice_setup(), and scale_ocean_phy_tc::ocean_phy_tc_seaice().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_calendar
module CALENDAR
Definition: scale_calendar.F90:13
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_ocean_phy_ice_simple::ocean_phy_ice_freezetemp
real(rp), public ocean_phy_ice_freezetemp
Definition: scale_ocean_phy_ice_simple.F90:44
scale_const::const_dwatr
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:89
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:146
scale_file_external_input
module file / external_input
Definition: scale_file_external_input.F90:12
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:151
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:153
scale_calendar::calendar_unit2sec
subroutine, public calendar_unit2sec(second, value, unit)
Convert several units to second.
Definition: scale_calendar.F90:486