SCALE-RM
Functions/Subroutines | Variables
scale_const Module Reference

module CONSTANT More...

Functions/Subroutines

subroutine, public const_setup
 Setup. More...
 
subroutine, public const_finalize
 Finalize. More...
 

Variables

real(rp), parameter, public const_pi = 4.0_RP * atan( 1.0_RP )
 pi More...
 
real(rp), public const_d2r
 degree to radian More...
 
real(rp), public const_r2d
 radian to degree More...
 
real(rp), public const_eps = 1.E-16_RP
 small number More...
 
real(rp), public const_eps1 = 0.99999999999999_RP
 small number More...
 
real(rp), public const_huge = 1.E+30_RP
 huge number More...
 
integer, parameter, public const_undef2 = -32768
 undefined value (INT2) More...
 
real(sp), parameter, public const_undef4 = -9.9999E30
 undefined value (REAL4) More...
 
real(dp), parameter, public const_undef8 = -9.9999D30
 undefined value (REAL8) More...
 
real(rp), public const_undef
 
real(rp), public const_radius
 radius of the planet [m] More...
 
real(rp), public const_ohm
 angular velocity of the planet [1/s] More...
 
real(rp), public const_grav
 standard acceleration of gravity [m/s2] More...
 
real(rp), parameter, public const_stb = 5.670373E-8_RP
 Stefan-Boltzman constant [W/m2/K4]. More...
 
real(rp), parameter, public const_karman = 0.4_RP
 von Karman constant More...
 
real(rp), parameter, public const_r = 8.31436_RP
 universal gas constant [J/mol/K] More...
 
real(rp), public const_mdry = 28.966_RP
 mass weight (dry air) [g/mol] More...
 
real(rp), public const_rdry
 specific gas constant (dry air) [J/kg/K] More...
 
real(rp), public const_cpdry
 specific heat (dry air,constant pressure) [J/kg/K] More...
 
real(rp), public const_cvdry
 specific heat (dry air,constant volume) [J/kg/K] More...
 
real(rp), public const_laps
 lapse rate of ISA [K/m] More...
 
real(rp), public const_lapsdry
 dry adiabatic lapse rate [K/m] More...
 
real(rp), public const_mvap = 18.016_RP
 mass weight (water vapor) [g/mol] More...
 
real(rp), parameter, public const_rvap = 461.50_RP
 specific gas constant (water vapor) [J/kg/K] More...
 
real(rp), parameter, public const_cpvap = 1846.00_RP
 specific heat (water vapor, constant pressure) [J/kg/K] More...
 
real(rp), public const_cvvap
 specific heat (water vapor, constant volume) [J/kg/K] More...
 
real(rp), parameter, public const_cl = 4218.0_RP
 specific heat (liquid water) [J/kg/K] More...
 
real(rp), parameter, public const_ci = 2106.0_RP
 specific heat (ice) [J/kg/K] More...
 
real(rp), public const_epsvap
 Rdry / Rvap. More...
 
real(rp), public const_epstvap
 1 / epsilon - 1 More...
 
real(rp), parameter, public const_emelt = 3.4E5_RP
 
real(rp), parameter, public const_tmelt = 273.15_RP
 
real(rp), parameter, public const_lhv0 = 2.501E+6_RP
 latent heat of vaporizaion at 0C [J/kg] More...
 
real(rp), public const_lhv00
 latent heat of vaporizaion at 0K [J/kg] More...
 
real(rp), parameter, public const_lhs0 = 2.834E+6_RP
 latent heat of sublimation at 0C [J/kg] More...
 
real(rp), public const_lhs00
 latent heat of sublimation at 0K [J/kg] More...
 
real(rp), public const_lhf0
 latent heat of fusion at 0C [J/kg] More...
 
real(rp), public const_lhf00
 latent heat of fusion at 0K [J/kg] More...
 
real(rp), parameter, public const_psat0 = 610.78_RP
 saturate pressure of water vapor at 0C [Pa] More...
 
real(rp), parameter, public const_dwatr = 1000.0_RP
 density of water [kg/m3] More...
 
real(rp), parameter, public const_dice = 916.8_RP
 density of ice [kg/m3] More...
 
real(rp), public const_sound
 speed of sound (dry air at 0C) [m/s] More...
 
real(rp), public const_pstd
 standard pressure [Pa] More...
 
real(rp), public const_pre00
 pressure reference [Pa] More...
 
real(rp), public const_tstd
 standard temperature (15C) [K] More...
 
real(rp), parameter, public const_tem00 = 273.15_RP
 temperature reference (0C) [K] More...
 
real(rp), parameter, public const_ppm = 1.E-6_RP
 parts par million More...
 
real(rp), parameter, public const_epsvac = 8.854187817E-12_RP
 parts par million More...
 
real(rp), public const_epsair = 1.00059_RP
 parts par million More...
 
integer, public const_i_lw = 1
 long-wave radiation index More...
 
integer, public const_i_sw = 2
 short-wave radiation index More...
 
character(len=h_short), public const_thermodyn_type = 'EXACT'
 internal energy type More...
 

Detailed Description

module CONSTANT

Description
Physical constants module
Author
Team SCALE
NAMELIST
  • PARAM_CONST
    nametypedefault valuecomment
    CONST_RADIUS real(RP) radius of the planet [m]
    CONST_OHM real(RP) angular velocity of the planet [1/s]
    CONST_GRAV real(RP) standard acceleration of gravity [m/s2]
    CONST_RDRY real(RP) specific gas constant (dry air) [J/kg/K]
    CONST_CPDRY real(RP) specific heat (dry air,constant pressure) [J/kg/K]
    CONST_LAPS real(RP) lapse rate of ISA [K/m]
    CONST_PSTD real(RP) standard pressure [Pa]
    CONST_PRE00 real(RP) pressure reference [Pa]
    CONST_TSTD real(RP) standard temperature (15C) [K]
    CONST_THERMODYN_TYPE character(len=H_SHORT) 'EXACT' internal energy type
    CONST_SMALLPLANETFACTOR real(RP) factor for small planet

History Output
No history output

Function/Subroutine Documentation

◆ const_setup()

subroutine, public scale_const::const_setup

Setup.

Definition at line 128 of file scale_const.F90.

128  use scale_prc, only: &
129  prc_abort
130  implicit none
131 
132  real(RP) :: CONST_SmallPlanetFactor
133 
134  namelist / param_const / &
135  const_radius, &
136  const_ohm, &
137  const_grav, &
138  const_rdry, &
139  const_cpdry, &
140  const_laps, &
141  const_pstd, &
142  const_pre00, &
143  const_tstd, &
144  const_thermodyn_type, &
145  const_smallplanetfactor
146 
147  integer :: ierr
148  !---------------------------------------------------------------------------
149 
150  if ( initialized ) return
151  initialized = .true.
152 
153  log_newline
154  log_info("CONST_setup",*) 'Setup'
155 
156  const_radius = 6.37122e+6_rp
157  const_ohm = 7.2920e-5_rp
158  const_grav = 9.80665_rp
159  const_rdry = 287.04_rp
160  const_cpdry = 1004.64_rp
161  const_laps = 6.5e-3_rp
162  const_pstd = 101325.0_rp
163  const_pre00 = 100000.0_rp
164  const_tstd = 288.15_rp
165 
166  const_thermodyn_type = 'EXACT'
167 
168  const_smallplanetfactor = 1.0_rp
169 
170  !--- read namelist
171  rewind(io_fid_conf)
172  read(io_fid_conf,nml=param_const,iostat=ierr)
173  if( ierr < 0 ) then !--- missing
174  log_info("CONST_setup",*) 'Not found namelist. Default used.'
175  elseif( ierr > 0 ) then !--- fatal error
176  log_error("CONST_setup",*) 'Not appropriate names in namelist PARAM_CONST. Check!'
177  call prc_abort
178  endif
179  log_nml(param_const)
180 
181  if ( rp == sp ) then
182  const_undef = real(const_undef4,kind=rp)
183  elseif( rp == dp ) then
184  const_undef = real(const_undef8,kind=rp)
185  else
186  log_error("CONST_setup",*) 'unsupported precision: ', rp
187  call prc_abort
188  endif
189 
190  const_d2r = const_pi / 180.0_rp
191  const_r2d = 180.0_rp / const_pi
192  const_eps = epsilon(0.0_rp)
193  const_eps1 = 1.0_rp - epsilon(0.0_rp)
194  const_huge = huge(0.0_rp)
195 
196  const_radius = const_radius / const_smallplanetfactor
197  const_ohm = const_ohm * const_smallplanetfactor
198 
199  const_cvdry = const_cpdry - const_rdry
200  const_lapsdry = const_grav / const_cpdry
201 
202  const_cvvap = const_cpvap - const_rvap
203  const_epsvap = const_rdry / const_rvap
204  const_epstvap = 1.0_rp / const_epsvap - 1.0_rp
205 
206  const_lhf0 = const_lhs0 - const_lhv0
207 
208  const_lhv00 = const_lhv0 - ( const_cpvap - const_cl ) * const_tem00
209  const_lhs00 = const_lhs0 - ( const_cpvap - const_ci ) * const_tem00
210  const_lhf00 = const_lhf0 - ( const_cl - const_ci ) * const_tem00
211 
212  const_sound = sqrt( const_cpdry * const_rdry / ( const_cpdry - const_rdry ) * const_tem00 )
213 
214  log_newline
215  log_info("CONST_setup",*) 'List of constants '
216  log_info_cont(*) 'PI : PI = ', const_pi
217  log_info_cont(*) 'Small number : EPS = ', const_eps
218  log_info_cont(*) 'Small number (1-EPS) : EPS1 = ', const_eps1
219  log_info_cont(*) 'Huge number : HUGE = ', const_huge
220  log_info_cont(*) 'undefined number (INT2) : UNDEF2 = ', const_undef2
221  log_info_cont(*) 'undefined number (REAL,general use) : UNDEF = ', const_undef
222  log_info_cont(*) 'undefined number (REAL4) : UNDEF4 = ', const_undef4
223  log_info_cont(*) 'undefined number (REAL8) : UNDEF8 = ', const_undef8
224 
225  log_info_cont(*) 'radius of the planet [m] : RADIUS = ', const_radius
226  log_info_cont(*) 'angular velocity of the planet [1/s] : OHM = ', const_ohm
227  log_info_cont(*) 'standard acceleration of gravity [m/s2] : GRAV = ', const_grav
228 
229  log_info_cont(*) 'Stefan-Boltzman constant [W/m2/K4] : STB = ', const_stb
230  log_info_cont(*) 'von Karman constant : KARMAN = ', const_karman
231  log_info_cont(*) 'universal gas constant [J/mol/K] : R = ', const_r
232 
233  log_info_cont(*) 'mass weight (dry air) [g/mol] : Mdry = ', const_mdry
234  log_info_cont(*) 'specific gas constant (dry air) [J/kg/K] : Rdry = ', const_rdry
235  log_info_cont(*) 'specific heat (dry air, const. pressure) [J/kg/K] : CPdry = ', const_cpdry
236  log_info_cont(*) 'specific heat (dry air, const. volume) [J/kg/K] : Cvdry = ', const_cvdry
237  log_info_cont(*) 'lapse rate of ISA [K/m] : LAPS = ', const_laps
238  log_info_cont(*) 'dry adiabatic lapse rate [K/m] : LAPSdry = ', const_lapsdry
239 
240  log_info_cont(*) 'mass weight (water vapor) [g/mol] : Mvap = ', const_mvap
241  log_info_cont(*) 'specific gas constant (water vapor) [J/kg/K] : Rvap = ', const_rvap
242  log_info_cont(*) 'specific heat (vapor, const. pressure) [J/kg/K] : CPvap = ', const_cpvap
243  log_info_cont(*) 'specific heat (vapor, const. volume) [J/kg/K] : CVvap = ', const_cvvap
244  log_info_cont(*) 'specific heat (liquid water) [J/kg/K] : CL = ', const_cl
245  log_info_cont(*) 'specific heat (ice) [J/kg/K] : CI = ', const_ci
246  log_info_cont(*) 'Rdry / Rvap : EPSvap = ', const_epsvap
247  log_info_cont(*) '1 / EPSvap - 1 : EPSTvap = ', const_epstvap
248 
249  log_info_cont(*) 'latent heat of vaporizaion at 0C [J/kg] : LHV0 = ', const_lhv0
250  log_info_cont(*) 'latent heat of sublimation at 0C [J/kg] : LHS0 = ', const_lhs0
251  log_info_cont(*) 'latent heat of fusion at 0C [J/kg] : LHF0 = ', const_lhf0
252  log_info_cont(*) 'latent heat of vaporizaion at 0K [J/kg] : LHV00 = ', const_lhv00
253  log_info_cont(*) 'latent heat of sublimation at 0K [J/kg] : LHS00 = ', const_lhs00
254  log_info_cont(*) 'latent heat of fusion at 0K [J/kg] : LHF00 = ', const_lhf00
255  log_info_cont(*) 'Thermodynamics calculation type : ', trim(const_thermodyn_type)
256  log_info_cont(*) 'saturate pressure of water vapor at 0C [Pa] : PSAT0 = ', const_psat0
257  log_info_cont(*) 'density of water [kg/m3] : DWATR = ', const_dwatr
258  log_info_cont(*) 'density of ice [kg/m3] : DICE = ', const_dice
259 
260  log_info_cont(*) 'speed of sound (dry air at 0C) [m/s] : SOUND = ', const_sound
261  log_info_cont(*) 'standard pressure [Pa] : Pstd = ', const_pstd
262  log_info_cont(*) 'pressure reference [Pa] : PRE00 = ', const_pre00
263  log_info_cont(*) 'standard temperature (15C) [K] : Tstd = ', const_tstd
264  log_info_cont(*) 'temperature reference (0C) [K] : TEM00 = ', const_tem00
265 
266  !$acc update device(CONST_D2R, CONST_EPS, CONST_EPS1, CONST_HUGE)
267  !$acc update device(CONST_UNDEF)
268  !$acc update device(CONST_RADIUS, CONST_OHM, CONST_GRAV)
269  !$acc update device(CONST_Mdry, CONST_Rdry, CONST_CPdry, CONST_CVdry, CONST_LAPS, CONST_LAPSdry)
270  !$acc update device(CONST_Mvap)
271  !$acc update device(CONST_EPSvap, CONST_EPSTvap)
272  !$acc update device(CONST_LHV00, CONST_LHS00, CONST_LHF0, CONST_LHF00)
273  !$acc update device(CONST_SOUND)
274  !$acc update device(CONST_Pstd, CONST_PRE00, CONST_Tstd)
275  !$acc update device(CONST_EPSair)
276  !$acc update device(CONST_I_LW, CONST_I_SW)
277  !$acc update device(CONST_THERMODYN_TYPE)
278 
279  return

References const_ci, const_cl, const_cpdry, const_cpvap, const_cvdry, const_cvvap, const_d2r, const_dice, const_dwatr, const_eps, const_eps1, const_epstvap, const_epsvap, const_grav, const_huge, const_karman, const_laps, const_lapsdry, const_lhf0, const_lhf00, const_lhs0, const_lhs00, const_lhv0, const_lhv00, const_mdry, const_mvap, const_ohm, const_pi, const_pre00, const_psat0, const_pstd, const_r, const_r2d, const_radius, const_rdry, const_rvap, const_sound, const_stb, const_tem00, const_thermodyn_type, const_tstd, const_undef, const_undef2, const_undef4, const_undef8, scale_precision::dp, scale_io::io_fid_conf, scale_prc::prc_abort(), scale_precision::rp, and scale_precision::sp.

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_setup(), mod_rm_driver::rm_driver(), mod_rm_prep::rm_prep(), and scale::scale_init().

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

◆ const_finalize()

subroutine, public scale_const::const_finalize

Finalize.

Definition at line 285 of file scale_const.F90.

285 
286  initialized = .false.
287 
288  return

Referenced by mod_rm_driver::rm_driver(), mod_rm_prep::rm_prep(), and scale::scale_finalize().

Here is the caller graph for this function:

Variable Documentation

◆ const_pi

real(rp), parameter, public scale_const::const_pi = 4.0_RP * atan( 1.0_RP )

pi

Definition at line 32 of file scale_const.F90.

32  real(RP), public, parameter :: CONST_PI = 4.0_rp * atan( 1.0_rp )

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_set(), scale_atmos_dyn_common::atmos_dyn_wdamp_setup(), scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_setup(), scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_mkinit(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_setup(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qhyd2qtrc(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_setup(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tendency(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_effective_radius(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_sf_const::atmos_phy_sf_const_flux(), scale_atmos_solarins::atmos_solarins_ecliptic_longitude(), scale_atmos_solarins::atmos_solarins_orbit(), scale_bulkflux::bulkflux_univfunc(), mod_cnv2d::cnv2d_exec(), mod_cnvuser::cnvuser(), com_gamma(), const_setup(), mod_copytopo::copytopo(), scale_atmos_phy_mp_sn14::dep_vapor_ice_wrk(), mod_mkinit::faero(), scale_file_tiledata::file_tiledata_get_data_int1(), scale_file_tiledata::file_tiledata_get_domain_info(), scale_file_tiledata::file_tiledata_get_info(), scale_file_tiledata::file_tiledata_get_latlon(), scale_file_tiledata::file_tiledata_get_tile_info(), mod_realinput::land_interporation(), scale_letkf::letkf_add_inflation_setup(), scale_mapprojection::mapprojection_setup(), mod_mktopo::mktopo(), scale_atmos_phy_mp_sn14::nucleation_ice_hom(), mod_realinput::ocean_interporation(), scale_letkf::radar_superobing(), scale_random::random_knuth_shuffle(), mod_realinput::realinput_surface(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01(), and scale_vector::vectr_triangle().

◆ const_d2r

real(rp), public scale_const::const_d2r

degree to radian

Definition at line 33 of file scale_const.F90.

33  real(RP), public :: CONST_D2R

Referenced by scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_setup(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_finalize(), scale_atmos_sfc_ch_rn222::atmos_sfc_ch_rn222_setup(), scale_atmos_solarins::atmos_solarins_orbit(), scale_atmos_solarins::atmos_solarins_setup(), mod_cnv2d::cnv2d_grads_init(), mod_cnv2d::cnv2d_tile_init(), mod_cnvlanduse::cnvlanduse(), mod_cnvtopo::cnvtopo(), mod_cnvtopo::cnvtopo_setup(), mod_cnvuser::cnvuser(), com_gamma(), scale_comm_cartesc_nest::comm_cartesc_nest_setup(), const_setup(), mod_copytopo::copytopo_get_data_grads(), mod_copytopo::copytopo_get_data_scale(), mod_copytopo::copytopo_get_data_wrfarw(), scale_file_cartesc::file_cartesc_set_areavol_atmos(), scale_file_cartesc::file_cartesc_set_coordinates_atmos(), scale_file_history_cartesc::file_history_cartesc_truncate_3d(), scale_file_tiledata::file_tiledata_get_data_int1(), scale_file_tiledata::file_tiledata_get_latlon(), scale_file_tiledata::file_tiledata_read_catalog_file(), scale_interp::interp_domain_compatibility(), scale_letkf::letkf_add_inflation_setup(), scale_mapprojection::mapprojection_setup(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_nicam::parentatmosinputnicam(), mod_realinput_nicam::parentatmosopennicam(), mod_realinput_grads::parentatmossetupgrads(), mod_realinput_netcdf::parentatmossetupnetcdf(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_netcdf::parentlandinputnetcdf(), mod_realinput_nicam::parentlandinputnicam(), mod_realinput_grads::parentlandsetupgrads(), mod_realinput_netcdf::parentlandsetupnetcdf(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_netcdf::parentoceaninputnetcdf(), mod_realinput_nicam::parentoceaninputnicam(), mod_realinput_nicam::parentoceanopennicam(), mod_realinput_grads::parentoceansetupgrads(), mod_realinput_netcdf::parentoceansetupnetcdf(), scale_letkf::radar_superobing(), scale_atmos_phy_rd_mm5sw::swpara(), and mod_urban_driver::urban_driver_calc_tendency().

◆ const_r2d

real(rp), public scale_const::const_r2d

radian to degree

Definition at line 34 of file scale_const.F90.

34  real(RP), public :: CONST_R2D

Referenced by com_gamma(), const_setup(), scale_letkf::letkf_add_inflation_setup(), and scale_letkf::radar_superobing().

◆ const_eps

real(rp), public scale_const::const_eps = 1.E-16_RP

small number

Definition at line 35 of file scale_const.F90.

35  real(RP), public :: CONST_EPS = 1.e-16_rp

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_set(), scale_atmos_diagnostic::atmos_diagnostic_get_phyd(), scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_fct(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxy_xvz_ud3koren1993(), scale_atmos_dyn_fvm_flux_cd2::atmos_dyn_fvm_fluxz_xyz_cd2(), scale_atmos_dyn_fvm_flux_cd4::atmos_dyn_fvm_fluxz_xyz_cd4(), scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xyz_cd6(), scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxz_xyz_cd8(), scale_atmos_dyn_fvm_flux_ud1::atmos_dyn_fvm_fluxz_xyz_ud1(), scale_atmos_dyn_fvm_flux_ud3::atmos_dyn_fvm_fluxz_xyz_ud3(), scale_atmos_dyn_fvm_flux_ud3koren1993::atmos_dyn_fvm_fluxz_xyz_ud3koren1993(), scale_atmos_dyn_fvm_flux_ud5::atmos_dyn_fvm_fluxz_xyz_ud5(), scale_atmos_dyn_fvm_flux_ud7::atmos_dyn_fvm_fluxz_xyz_ud7(), scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), scale_atmos_dyn_tstep_short_fvm_heve::atmos_dyn_tstep_short_fvm_heve(), scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi(), scale_atmos_dyn_common::atmos_dyn_wdamp_setup(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tendency(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_electric_field(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_f2013(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_mg2001(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_adjustment(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_setup(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_effective_radius(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver_calc_tendency(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_flux(), scale_atmos_phy_tb_smg::atmos_phy_tb_smg(), scale_atmos_refstate::atmos_refstate_smoothing(), scale_atmos_saturation::atmos_saturation_moist_conversion_pres_liq_0d(), scale_atmos_solarins::atmos_solarins_insolation_0d(), scale_atmos_solarins::atmos_solarins_insolation_2d(), scale_atmos_phy_mp_sn14::aut_acc_slc_brk(), scale_bulkflux::bulkflux_diagnose_scales_0d(), scale_bulkflux::bulkflux_diagnose_scales_2d(), mod_cnvlanduse::cnvlanduse(), mod_cnvtopo::cnvtopo(), const_setup(), mod_copytopo::copytopo(), mod_cpl_vars::cpl_getsfc_atm(), scale_cpl_phy_sfc_fixed_temp::cpl_phy_sfc_fixed_temp(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), scale_interp::cross(), scale_file_cartesc::file_cartesc_write_var_4d(), scale_file_external_input::file_external_input_put_ref_1d(), scale_file_external_input::file_external_input_put_ref_2d(), scale_file_external_input::file_external_input_put_ref_3d(), scale_file_external_input::file_external_input_update_1d(), scale_file_external_input::file_external_input_update_2d(), scale_file_external_input::file_external_input_update_3d(), scale_file_grads::file_grads_close(), scale_file_history::file_history_in_0d(), scale_file_history::file_history_in_1d(), scale_file_history::file_history_in_2d(), scale_file_history::file_history_in_3d(), scale_file_history::file_history_reg(), mod_realinput::get_ijrange(), scale_interp::interp_bilinear_inv(), scale_interp::interp_div_block(), scale_interp::interp_factor1d(), scale_interp::interp_factor2d_linear_xy(), scale_interp::interp_insert_2d(), scale_interp::interp_interp1d(), scale_interp::interp_interp2d(), scale_interp::interp_interp3d(), scale_interp::interp_search_horiz_struct(), mod_realinput::land_interporation(), scale_land_phy_snow_diagnos::land_phy_snow_diags(), scale_land_phy_snow_ky90::land_phy_snow_ky90(), mod_realinput::make_mask(), mod_ocean_driver::ocean_driver_calc_tendency(), scale_ocean_phy_ice_simple::ocean_phy_ice_adjustment(), scale_ocean_phy_tc::ocean_phy_tc_seaice(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_nicam::parentatmosinputnicam(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_nicam::parentlandinputnicam(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_nicam::parentoceaninputnicam(), mod_realinput::realinput_surface(), mod_realinput::replace_misval_const(), mod_realinput::replace_misval_map(), scale_interp::spline_exec(), scale_statistics::statistics_total_2d(), scale_statistics::statistics_total_3d(), mod_mkinit::tke_setup(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01(), scale_vector::vectr_anticlockwise(), scale_vector::vectr_intersec(), scale_vector::vectr_triangle(), and scale_vector::vectr_xyz2latlon().

◆ const_eps1

real(rp), public scale_const::const_eps1 = 0.99999999999999_RP

small number

Definition at line 36 of file scale_const.F90.

36  real(RP), public :: CONST_EPS1 = 0.99999999999999_rp

Referenced by scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), and const_setup().

◆ const_huge

real(rp), public scale_const::const_huge = 1.E+30_RP

◆ const_undef2

integer, parameter, public scale_const::const_undef2 = -32768

undefined value (INT2)

Definition at line 40 of file scale_const.F90.

40  integer, public, parameter :: CONST_UNDEF2 = -32768

Referenced by scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_fct(), mod_cnvlanduse::cnvlanduse(), const_setup(), and scale_file_tiledata::file_tiledata_get_data_int1().

◆ const_undef4

real(sp), parameter, public scale_const::const_undef4 = -9.9999E30

undefined value (REAL4)

Definition at line 41 of file scale_const.F90.

41  real(SP), public, parameter :: CONST_UNDEF4 = -9.9999e30

Referenced by const_setup().

◆ const_undef8

real(dp), parameter, public scale_const::const_undef8 = -9.9999D30

undefined value (REAL8)

Definition at line 42 of file scale_const.F90.

42  real(DP), public, parameter :: CONST_UNDEF8 = -9.9999d30

Referenced by mod_admin_time::admin_time_setup(), const_setup(), mod_mkinit::mkinit(), and mod_mkinit::rect_setup().

◆ const_undef

real(rp), public scale_const::const_undef

Definition at line 43 of file scale_const.F90.

43  real(RP), public :: CONST_UNDEF

Referenced by mod_atmos_vars::allocate_1d(), mod_atmos_vars::allocate_2d(), mod_atmos_vars::allocate_3d(), mod_atmos_bnd_driver::atmos_boundary_driver_setup(), scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_fct(), scale_atmos_dyn::atmos_dyn_setup(), scale_atmos_dyn_tinteg_large_euler::atmos_dyn_tinteg_large_euler_setup(), scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o_setup(), scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3_setup(), scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4_setup(), scale_atmos_dyn_tinteg_short_rk7s6o::atmos_dyn_tinteg_short_rk7s6o_setup(), scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3_setup(), scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve_setup(), mod_atmos_dyn_vars::atmos_dyn_vars_setup(), scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_calc_areavol(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_real_3d(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_mkinit(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_qtrc2qaero(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_setup(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_tendency(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_setup(), mod_atmos_phy_bl_driver::atmos_phy_bl_driver_mkinit(), mod_atmos_phy_bl_vars::atmos_phy_bl_vars_setup(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_setup(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_setup(), mod_atmos_phy_lt_vars::atmos_phy_lt_vars_setup(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qhyd2qtrc(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_terminal_velocity(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup(), scale_atmos_phy_rd_offline::atmos_phy_rd_offline_setup(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_setup(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_calc_tendency(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_setup(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_setup(), scale_atmos_refstate::atmos_refstate_fillhalo(), scale_atmos_refstate::atmos_refstate_setup(), scale_atmos_saturation::atmos_saturation_dqs_dtem_dpre_ice_1d(), scale_atmos_solarins::atmos_solarins_orbit(), scale_atmos_solarins::atmos_solarins_setup(), mod_atmos_vars::atmos_vars_setup(), scale_debug::check(), mod_cnv2d::cnv2d_tile_init(), mod_cnvtopo::cnvtopo(), com_gamma(), scale_comm_cartesc::comm_horizontal_mean_2d(), scale_comm_cartesc::comm_horizontal_mean_3d(), const_setup(), scale_cpl_phy_sfc_fixed_temp::cpl_phy_sfc_fixed_temp(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), mod_cpl_vars::cpl_vars_setup(), scale_interp::cross(), mod_da_driver::da_driver_update(), mod_da_param_estimation::da_param_estimation_update(), scale_file_cartesc::file_cartesc_def_axes(), scale_file_cartesc::file_cartesc_write_axes(), scale_file_external_input::file_external_input_put_ref_1d(), scale_file_external_input::file_external_input_put_ref_2d(), scale_file_external_input::file_external_input_put_ref_3d(), scale_file_external_input::file_external_input_regist_external_1d(), scale_file_external_input::file_external_input_regist_external_2d(), scale_file_external_input::file_external_input_regist_external_3d(), scale_file_external_input::file_external_input_setup(), scale_file_external_input::file_external_input_update_1d(), scale_file_external_input::file_external_input_update_2d(), scale_file_external_input::file_external_input_update_3d(), scale_file_grads::file_grads_close(), scale_file_grads::file_grads_open(), scale_file_history_cartesc::file_history_cartesc_truncate_3d(), scale_file_history::file_history_in_0d(), scale_file_history::file_history_in_1d(), scale_file_history::file_history_in_2d(), scale_file_history::file_history_in_3d(), scale_file_history::file_history_reg(), scale_file_tiledata::file_tiledata_get_data_int1(), scale_file_tiledata::file_tiledata_get_latlon(), mod_atmos_phy_sf_driver::history_output(), scale_interp::interp_div_block(), scale_interp::interp_factor1d(), scale_interp::interp_factor2d_linear_xy(), scale_interp::interp_insert_2d(), scale_interp::interp_interp1d(), scale_interp::interp_interp2d(), scale_interp::interp_interp3d(), scale_interp::interp_search_horiz_struct(), scale_land_dyn_bucket::land_dyn_bucket(), scale_land_dyn_bucket::land_dyn_bucket_setup(), mod_realinput::land_interporation(), scale_land_phy_snow_diagnos::land_phy_snow_diags(), mod_land_vars::land_vars_setup(), scale_letkf::letkf_obs_initialize(), scale_letkf::letkf_param_estimation_system(), scale_letkf::letkf_system(), mod_realinput::make_mask(), scale_mapprojection::mapprojection_lonlat2xy_equidistantcylindrical(), scale_mapprojection::mapprojection_lonlat2xy_lambertconformal(), scale_mapprojection::mapprojection_lonlat2xy_mercator(), scale_mapprojection::mapprojection_lonlat2xy_none(), scale_mapprojection::mapprojection_lonlat2xy_polarstereographic(), scale_mapprojection::mapprojection_mapfactor_none(), scale_mapprojection::mapprojection_setup(), scale_mapprojection::mapprojection_xy2lonlat_equidistantcylindrical(), scale_mapprojection::mapprojection_xy2lonlat_lambertconformal(), scale_mapprojection::mapprojection_xy2lonlat_mercator(), scale_mapprojection::mapprojection_xy2lonlat_none(), scale_mapprojection::mapprojection_xy2lonlat_polarstereographic(), mod_mkinit::mkinit(), scale_ocean_dyn_offline::ocean_dyn_offline_setup(), scale_ocean_dyn_slab::ocean_dyn_slab(), scale_ocean_dyn_slab::ocean_dyn_slab_setup(), mod_realinput::ocean_interporation(), scale_ocean_phy_ice_simple::ocean_phy_ice_setup(), scale_ocean_phy_roughness_moon07::ocean_phy_roughness_moon07(), mod_ocean_vars::ocean_vars_setup(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_netcdf::parentatmosinputnetcdf(), mod_realinput_netcdf::parentatmossetupnetcdf(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_netcdf::parentlandinputnetcdf(), mod_realinput_nicam::parentlandinputnicam(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_netcdf::parentoceaninputnetcdf(), mod_realinput_nicam::parentoceaninputnicam(), mod_realinput_nicam::parentoceanopennicam(), mod_realinput::realinput_surface(), mod_realinput::replace_misval_map(), scale_interp::spline_exec(), scale_statistics::statistics_horizontal_max_2d(), scale_statistics::statistics_horizontal_max_3d(), scale_statistics::statistics_horizontal_mean_2d(), scale_statistics::statistics_horizontal_mean_3d(), scale_statistics::statistics_horizontal_min_2d(), scale_statistics::statistics_horizontal_min_3d(), scale_statistics::statistics_total_2d(), scale_statistics::statistics_total_3d(), mod_mkinit::tke_setup(), mod_urban_driver::urban_driver_setup(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01_setup(), and mod_urban_vars::urban_vars_setup().

◆ const_radius

real(rp), public scale_const::const_radius

◆ const_ohm

real(rp), public scale_const::const_ohm

angular velocity of the planet [1/s]

Definition at line 48 of file scale_const.F90.

48  real(RP), public :: CONST_OHM

Referenced by scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve_setup(), const_setup(), scale_coriolis::coriolis_setup(), and mod_mkinit::read_sounding().

◆ const_grav

real(rp), public scale_const::const_grav

standard acceleration of gravity [m/s2]

Definition at line 49 of file scale_const.F90.

49  real(RP), public :: CONST_GRAV

Referenced by scale_atmos_adiabat::atmos_adiabat_cape_1d(), scale_atmos_bottom::atmos_bottom_estimate(), scale_atmos_diagnostic::atmos_diagnostic_get_n2(), scale_atmos_diagnostic::atmos_diagnostic_get_phyd(), scale_atmos_dyn_tstep_short_fvm_heve::atmos_dyn_tstep_short_fvm_heve(), scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi(), scale_atmos_dyn_tstep_short_fvm_hivi::atmos_dyn_tstep_short_fvm_hivi(), scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_calc_z(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tendency(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_semilag(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_upwind(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_putvar(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_read(), scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980(), scale_atmos_phy_tb_dns::atmos_phy_tb_dns(), scale_atmos_phy_tb_smg::atmos_phy_tb_smg(), mod_atmos_vars::atmos_vars_check(), mod_atmos_vars::atmos_vars_get_diagnostic_3d(), mod_atmos_vars::atmos_vars_restart_open(), scale_bulkflux::bulkflux_diagnose_scales_0d(), scale_bulkflux::bulkflux_diagnose_scales_2d(), scale_bulkflux::bulkflux_diagnose_surface_2d(), scale_bulkflux::calc_scales_b91w01(), com_gamma(), const_setup(), scale_atmos_phy_bl_mynn::get_phi(), scale_atmos_phy_mp_sn14::nucleation_ice_hom(), scale_ocean_phy_roughness_miller92::ocean_phy_roughness_miller92(), scale_ocean_phy_roughness_moon07::ocean_phy_roughness_moon07(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_netcdf::parentatmossetupnetcdf(), mod_realinput_grads::parentlandinputgrads(), scale_letkf::radar_superobing(), mod_mkinit::read_sounding(), mod_realinput::realinput_surface(), mod_rm_driver::rm_driver(), mod_rm_prep::rm_prep(), and scale_atmos_phy_rd_mm5sw::swrad().

◆ const_stb

real(rp), parameter, public scale_const::const_stb = 5.670373E-8_RP

Stefan-Boltzman constant [W/m2/K4].

Definition at line 53 of file scale_const.F90.

53  real(RP), public, parameter :: CONST_STB = 5.670373e-8_rp

Referenced by const_setup(), scale_cpl_phy_sfc_fixed_temp::cpl_phy_sfc_fixed_temp(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), scale_land_phy_snow_diagnos::land_phy_snow_diags(), mod_rm_driver::rm_driver(), and mod_rm_prep::rm_prep().

◆ const_karman

real(rp), parameter, public scale_const::const_karman = 0.4_RP

◆ const_r

real(rp), parameter, public scale_const::const_r = 8.31436_RP

universal gas constant [J/mol/K]

Definition at line 55 of file scale_const.F90.

55  real(RP), public, parameter :: CONST_R = 8.31436_rp

Referenced by const_setup().

◆ const_mdry

real(rp), public scale_const::const_mdry = 28.966_RP

mass weight (dry air) [g/mol]

Definition at line 58 of file scale_const.F90.

58  real(RP), public :: CONST_Mdry = 28.966_rp

Referenced by scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup_zgrid(), and const_setup().

◆ const_rdry

real(rp), public scale_const::const_rdry

◆ const_cpdry

real(rp), public scale_const::const_cpdry

◆ const_cvdry

real(rp), public scale_const::const_cvdry

◆ const_laps

real(rp), public scale_const::const_laps

◆ const_lapsdry

real(rp), public scale_const::const_lapsdry

dry adiabatic lapse rate [K/m]

Definition at line 63 of file scale_const.F90.

63  real(RP), public :: CONST_LAPSdry

Referenced by const_setup().

◆ const_mvap

real(rp), public scale_const::const_mvap = 18.016_RP

mass weight (water vapor) [g/mol]

Definition at line 67 of file scale_const.F90.

67  real(RP), public :: CONST_Mvap = 18.016_rp

Referenced by scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup_zgrid(), and const_setup().

◆ const_rvap

real(rp), parameter, public scale_const::const_rvap = 461.50_RP

◆ const_cpvap

real(rp), parameter, public scale_const::const_cpvap = 1846.00_RP

◆ const_cvvap

real(rp), public scale_const::const_cvvap

specific heat (water vapor, constant volume) [J/kg/K]

Definition at line 70 of file scale_const.F90.

70  real(RP), public :: CONST_CVvap

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_setup(), scale_atmos_saturation::atmos_saturation_setup(), const_setup(), scale_atmos_phy_mp_sn14::nucleation_ice_hom(), and mod_mkinit::read_sounding().

◆ const_cl

real(rp), parameter, public scale_const::const_cl = 4218.0_RP

◆ const_ci

real(rp), parameter, public scale_const::const_ci = 2106.0_RP

◆ const_epsvap

real(rp), public scale_const::const_epsvap

◆ const_epstvap

real(rp), public scale_const::const_epstvap

◆ const_emelt

real(rp), parameter, public scale_const::const_emelt = 3.4E5_RP

◆ const_tmelt

real(rp), parameter, public scale_const::const_tmelt = 273.15_RP

Definition at line 80 of file scale_const.F90.

80  real(RP), public, parameter :: CONST_TMELT = 273.15_rp

Referenced by scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd().

◆ const_lhv0

real(rp), parameter, public scale_const::const_lhv0 = 2.501E+6_RP

◆ const_lhv00

real(rp), public scale_const::const_lhv00

latent heat of vaporizaion at 0K [J/kg]

Definition at line 83 of file scale_const.F90.

83  real(RP), public :: CONST_LHV00

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_setup(), scale_atmos_saturation::atmos_saturation_setup(), const_setup(), and scale_atmos_phy_mp_sn14::nucleation_ice_hom().

◆ const_lhs0

real(rp), parameter, public scale_const::const_lhs0 = 2.834E+6_RP

◆ const_lhs00

real(rp), public scale_const::const_lhs00

latent heat of sublimation at 0K [J/kg]

Definition at line 85 of file scale_const.F90.

85  real(RP), public :: CONST_LHS00

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_setup(), scale_atmos_saturation::atmos_saturation_setup(), const_setup(), and scale_atmos_phy_mp_sn14::nucleation_ice_hom().

◆ const_lhf0

real(rp), public scale_const::const_lhf0

◆ const_lhf00

real(rp), public scale_const::const_lhf00

latent heat of fusion at 0K [J/kg]

Definition at line 87 of file scale_const.F90.

87  real(RP), public :: CONST_LHF00

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_setup(), const_setup(), scale_atmos_phy_mp_sn14::dep_vapor_ice_wrk(), and scale_atmos_phy_mp_sn14::nucleation_ice_hom().

◆ const_psat0

real(rp), parameter, public scale_const::const_psat0 = 610.78_RP

◆ const_dwatr

real(rp), parameter, public scale_const::const_dwatr = 1000.0_RP

◆ const_dice

real(rp), parameter, public scale_const::const_dice = 916.8_RP

◆ const_sound

real(rp), public scale_const::const_sound

speed of sound (dry air at 0C) [m/s]

Definition at line 93 of file scale_const.F90.

93  real(RP), public :: CONST_SOUND

Referenced by const_setup().

◆ const_pstd

real(rp), public scale_const::const_pstd

◆ const_pre00

real(rp), public scale_const::const_pre00

◆ const_tstd

real(rp), public scale_const::const_tstd

standard temperature (15C) [K]

Definition at line 98 of file scale_const.F90.

98  real(RP), public :: CONST_Tstd

Referenced by const_setup().

◆ const_tem00

real(rp), parameter, public scale_const::const_tem00 = 273.15_RP

temperature reference (0C) [K]

Definition at line 99 of file scale_const.F90.

99  real(RP), public, parameter :: CONST_TEM00 = 273.15_rp

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_entr_2d(), scale_atmos_hydrometeor::atmos_hydrometeor_lhf_1d(), scale_atmos_hydrometeor::atmos_hydrometeor_lhs_1d(), scale_atmos_hydrometeor::atmos_hydrometeor_lhv_0d(), scale_atmos_hydrometeor::atmos_hydrometeor_lhv_1d(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_select_dqcrg_from_lut(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_setup(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tendency(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_effective_radius(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_terminal_velocity(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_land_phy_snow_ky90::cal_param(), scale_land_phy_snow_ky90::cal_r1r2(), scale_land_phy_snow_ky90::calculationmo(), scale_land_phy_snow_ky90::calculationnomo(), scale_land_phy_snow_ky90::check_applicability(), scale_land_phy_snow_ky90::check_res(), const_setup(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), scale_atmos_phy_mp_sn14::dep_vapor_ice_wrk(), scale_land_phy_snow_ky90::groundflux(), mod_land_driver::land_driver_calc_tendency(), scale_land_dyn_bucket::land_dyn_bucket_setup(), scale_land_phy_snow_ky90::land_phy_snow_ky90(), scale_letkf::letkf_obs_initialize(), scale_atmos_phy_mp_sn14::nucleation_ice_hom(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_nicam::parentlandinputnicam(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_nicam::parentoceaninputnicam(), mod_realinput::realinput_surface(), scale_land_phy_snow_ky90::recalculatez(), mod_rm_driver::rm_driver(), mod_rm_prep::rm_prep(), scale_land_phy_snow_ky90::snow_ky90_main(), and scale_land_phy_snow_ky90::snowdepth().

◆ const_ppm

real(rp), parameter, public scale_const::const_ppm = 1.E-6_RP

parts par million

Definition at line 100 of file scale_const.F90.

100  real(RP), public, parameter :: CONST_PPM = 1.e-6_rp

Referenced by scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), and scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup_zgrid().

◆ const_epsvac

real(rp), parameter, public scale_const::const_epsvac = 8.854187817E-12_RP

parts par million

Definition at line 103 of file scale_const.F90.

103  real(RP), public, parameter :: CONST_EPSvac = 8.854187817e-12_rp

Referenced by scale_atmos_phy_lt_sato2019::atmos_phy_lt_electric_field(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_f2013(), and scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_mg2001().

◆ const_epsair

real(rp), public scale_const::const_epsair = 1.00059_RP

◆ const_i_lw

integer, public scale_const::const_i_lw = 1

long-wave radiation index

Definition at line 107 of file scale_const.F90.

107  integer, public :: CONST_I_LW = 1

Referenced by mod_realinput::land_interporation(), scale_land_phy_snow_ky90::land_phy_snow_ky90(), and mod_realinput::realinput_surface().

◆ const_i_sw

integer, public scale_const::const_i_sw = 2

short-wave radiation index

Definition at line 108 of file scale_const.F90.

108  integer, public :: CONST_I_SW = 2

Referenced by mod_realinput::land_interporation(), scale_land_phy_snow_ky90::land_phy_snow_ky90(), and mod_realinput::realinput_surface().

◆ const_thermodyn_type

character(len=h_short), public scale_const::const_thermodyn_type = 'EXACT'

internal energy type

Definition at line 111 of file scale_const.F90.

111  character(len=H_SHORT), public :: CONST_THERMODYN_TYPE = 'EXACT'

Referenced by scale_atmos_hydrometeor::atmos_hydrometeor_setup(), scale_atmos_saturation::atmos_saturation_setup(), and const_setup().

scale_precision::sp
integer, parameter, public sp
Definition: scale_precision.F90:31
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33