SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_atmos_saturation Module Reference

module ATMOSPHERE / Saturation adjustment More...

Functions/Subroutines

subroutine, public atmos_saturation_setup
 Setup. More...
 
subroutine atmos_saturation_alpha_0d (alpha, temp)
 calc liquid/ice separation factor (0D) More...
 
subroutine, public atmos_saturation_dqsw_dtem_rho (dqsdtem, temp, dens)
 
subroutine, public atmos_saturation_dqsi_dtem_rho (dqsdtem, temp, dens)
 
subroutine, public atmos_saturation_dqsw_dtem_dpre (dqsdtem, dqsdpre, temp, pres)
 
subroutine, public atmos_saturation_dqsi_dtem_dpre (dqsdtem, dqsdpre, temp, pres)
 

Variables

real(rp), public cpovr_liq
 
real(rp), public cpovr_ice
 
real(rp), public cvovr_liq
 
real(rp), public cvovr_ice
 
real(rp), public lovr_liq
 
real(rp), public lovr_ice
 

Detailed Description

module ATMOSPHERE / Saturation adjustment

Description
Saturation adjustment module
Author
Team SCALE
History
  • 2011-10-24 (T.Seiki) [new] Import from NICAM
  • 2012-02-10 (H.Yashiro) [mod] Reconstruction
NAMELIST
  • PARAM_ATMOS_SATURATION
    nametypedefault valuecomment
    ATMOS_SATURATION_ULIMIT_TEMP real(RP) 273.15_RP upper limit temperature
    ATMOS_SATURATION_LLIMIT_TEMP real(RP) 233.15_RP lower limit temperature

History Output
No history output

Function/Subroutine Documentation

◆ atmos_saturation_setup()

subroutine, public scale_atmos_saturation::atmos_saturation_setup ( )

Setup.

Definition at line 165 of file scale_atmos_sub_saturation.F90.

References scale_const::const_thermodyn_type, cpovr_ice, cpovr_liq, cvovr_ice, cvovr_liq, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, lovr_ice, lovr_liq, and scale_process::prc_mpistop().

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

165  use scale_process, only: &
167  use scale_const, only: &
169  implicit none
170 
171  namelist / param_atmos_saturation / &
172  atmos_saturation_ulimit_temp, &
173  atmos_saturation_llimit_temp
174 
175  integer :: ierr
176  !---------------------------------------------------------------------------
177 
178  if( io_l ) write(io_fid_log,*)
179  if( io_l ) write(io_fid_log,*) '++++++ Module[SATURATION] / Categ[ATMOS SHARE] / Origin[SCALElib]'
180 
181  !--- read namelist
182  rewind(io_fid_conf)
183  read(io_fid_conf,nml=param_atmos_saturation,iostat=ierr)
184  if( ierr < 0 ) then !--- missing
185  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
186  elseif( ierr > 0 ) then !--- fatal error
187  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_SATURATION. Check!'
188  call prc_mpistop
189  endif
190  if( io_lnml ) write(io_fid_log,nml=param_atmos_saturation)
191 
192  rtem00 = 1.0_rp / tem00
193 
194  if ( const_thermodyn_type == 'EXACT' ) then
195  cpovr_liq = ( cpvap - cl ) / rvap
196  cpovr_ice = ( cpvap - ci ) / rvap
197  cvovr_liq = ( cvvap - cl ) / rvap
198  cvovr_ice = ( cvvap - ci ) / rvap
199  elseif( const_thermodyn_type == 'SIMPLE' ) then
200  cpovr_liq = 0.0_rp
201  cpovr_ice = 0.0_rp
202  cvovr_liq = 0.0_rp
203  cvovr_ice = 0.0_rp
204  endif
205  lovr_liq = lhv / rvap
206  lovr_ice = lhs / rvap
207 
208 
209  dalphadt_const = 1.0_rp / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
210 
211  if( io_l ) write(io_fid_log,*)
212  if( io_l ) write(io_fid_log,'(1x,A,F7.2,A,F7.2)') '*** Temperature range for ice : ', &
213  atmos_saturation_llimit_temp, ' - ', &
214  atmos_saturation_ulimit_temp
215 
216  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
character(len=h_short), public const_thermodyn_type
internal energy type
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_saturation_alpha_0d()

subroutine scale_atmos_saturation::atmos_saturation_alpha_0d ( real(rp), intent(out)  alpha,
real(rp), intent(in)  temp 
)

calc liquid/ice separation factor (0D)

Parameters
[out]alphaliquid/ice separation factor [0-1]
[in]temptemperature [K]

Definition at line 224 of file scale_atmos_sub_saturation.F90.

References cpovr_ice, cpovr_liq, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, lovr_ice, and lovr_liq.

224  implicit none
225 
226  real(RP), intent(out) :: alpha
227  real(RP), intent(in) :: temp
228  !---------------------------------------------------------------------------
229 
230  alpha = ( temp - atmos_saturation_llimit_temp ) &
231  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
232 
233  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
234 
235  return

◆ atmos_saturation_dqsw_dtem_rho()

subroutine, public scale_atmos_saturation::atmos_saturation_dqsw_dtem_rho ( real(rp), dimension(ka,ia,ja), intent(out)  dqsdtem,
real(rp), dimension (ka,ia,ja), intent(in)  temp,
real(rp), dimension (ka,ia,ja), intent(in)  dens 
)

Definition at line 1242 of file scale_atmos_sub_saturation.F90.

References scale_const::const_lhv0, cpovr_liq, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and lovr_liq.

Referenced by scale_atmos_phy_mp_sn14::update_by_phase_change_kij().

1242  use scale_const, only: &
1243  lhv0 => const_lhv0
1244  implicit none
1245 
1246  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1247  real(RP), intent(in) :: temp (ka,ia,ja)
1248  real(RP), intent(in) :: dens (ka,ia,ja)
1249 
1250  real(RP) :: psat ! saturation vapor pressure
1251  real(RP) :: lhv ! latent heat for condensation
1252 
1253  real(RP) :: rtem00, tem
1254 
1255  integer :: k, i, j
1256  !---------------------------------------------------------------------------
1257 
1258  rtem00 = 1.0_rp / tem00
1259 
1260  !$omp parallel do private(i,j,k,TEM,psat,lhv) OMP_SCHEDULE_ collapse(2)
1261  do j = jsb, jeb
1262  do i = isb, ieb
1263  do k = ks, ke
1264  tem = max( temp(k,i,j), tem_min )
1265 
1266  psat = psat0 &
1267  * ( tem * rtem00 )**cpovr_liq &
1268  * exp( lovr_liq * ( rtem00 - 1.0_rp/tem ) )
1269 
1270  lhv = lhv0 + ( cpvap-cl ) * ( temp(k,i,j)-tem00 )
1271 
1272  dqsdtem(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) * temp(k,i,j) ) &
1273  * ( lhv / ( rvap * temp(k,i,j) ) - 1.0_rp )
1274  enddo
1275  enddo
1276  enddo
1277 
1278  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
integer, public ieb
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:80
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

◆ atmos_saturation_dqsi_dtem_rho()

subroutine, public scale_atmos_saturation::atmos_saturation_dqsi_dtem_rho ( real(rp), dimension(ka,ia,ja), intent(out)  dqsdtem,
real(rp), dimension (ka,ia,ja), intent(in)  temp,
real(rp), dimension (ka,ia,ja), intent(in)  dens 
)

Definition at line 1285 of file scale_atmos_sub_saturation.F90.

References scale_const::const_lhs0, cpovr_ice, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and lovr_ice.

Referenced by scale_atmos_phy_mp_sn14::nucleation_kij(), and scale_atmos_phy_mp_sn14::update_by_phase_change_kij().

1285  use scale_const, only: &
1286  lhs0 => const_lhs0
1287  implicit none
1288 
1289  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1290  real(RP), intent(in) :: temp (ka,ia,ja)
1291  real(RP), intent(in) :: dens (ka,ia,ja)
1292 
1293  real(RP) :: psat ! saturation vapor pressure
1294  real(RP) :: lhv ! latent heat for condensation
1295 
1296  real(RP) :: rtem00, tem
1297 
1298  integer :: k, i, j
1299  !---------------------------------------------------------------------------
1300 
1301  rtem00 = 1.0_rp / tem00
1302 
1303  !$omp parallel do private(i,j,k,TEM,psat,lhv) OMP_SCHEDULE_ collapse(2)
1304  do j = jsb, jeb
1305  do i = isb, ieb
1306  do k = ks, ke
1307  tem = max( temp(k,i,j), tem_min )
1308 
1309  psat = psat0 &
1310  * ( tem * rtem00 )**cpovr_ice &
1311  * exp( lovr_ice * ( rtem00 - 1.0_rp/tem ) )
1312  lhv = lhs0 + ( cpvap-ci ) * ( temp(k,i,j)-tem00 )
1313 
1314  dqsdtem(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) * temp(k,i,j) ) &
1315  * ( lhv / ( rvap * temp(k,i,j) ) - 1.0_rp )
1316  enddo
1317  enddo
1318  enddo
1319 
1320  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
integer, public ieb
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:82
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

◆ atmos_saturation_dqsw_dtem_dpre()

subroutine, public scale_atmos_saturation::atmos_saturation_dqsw_dtem_dpre ( real(rp), dimension(ka,ia,ja), intent(out)  dqsdtem,
real(rp), dimension(ka,ia,ja), intent(out)  dqsdpre,
real(rp), dimension (ka,ia,ja), intent(in)  temp,
real(rp), dimension (ka,ia,ja), intent(in)  pres 
)

Definition at line 1327 of file scale_atmos_sub_saturation.F90.

References scale_const::const_lhv0, cpovr_liq, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and lovr_liq.

Referenced by scale_atmos_phy_mp_sn14::update_by_phase_change_kij().

1327  use scale_const, only: &
1328  lhv0 => const_lhv0
1329  implicit none
1330 
1331  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1332  real(RP), intent(out) :: dqsdpre(ka,ia,ja)
1333  real(RP), intent(in) :: temp (ka,ia,ja)
1334  real(RP), intent(in) :: pres (ka,ia,ja)
1335 
1336  real(RP) :: psat ! saturation vapor pressure
1337  real(RP) :: lhv ! latent heat for condensation
1338 
1339  real(RP) :: den1, den2 ! denominator
1340  real(RP) :: rtem00, tem
1341 
1342  integer :: k, i, j
1343  !---------------------------------------------------------------------------
1344 
1345  rtem00 = 1.0_rp / tem00
1346 
1347  !$omp parallel do private(i,j,k,TEM,psat,den1,den2,lhv) OMP_SCHEDULE_ collapse(2)
1348  do j = jsb, jeb
1349  do i = isb, ieb
1350  do k = ks, ke
1351  tem = max( temp(k,i,j), tem_min )
1352 
1353  psat = psat0 &
1354  * ( tem * rtem00 )**cpovr_liq &
1355  * exp( lovr_liq * ( rtem00 - 1.0_rp/tem ) )
1356 
1357  den1 = ( pres(k,i,j) - (1.0_rp-epsvap) * psat ) &
1358  * ( pres(k,i,j) - (1.0_rp-epsvap) * psat )
1359  den2 = den1 * rvap * temp(k,i,j) * temp(k,i,j)
1360  lhv = lhv0 + ( cpvap-cl ) * ( temp(k,i,j)-tem00 )
1361 
1362  dqsdpre(k,i,j) = - epsvap * psat / den1
1363  dqsdtem(k,i,j) = epsvap * psat / den2 * lhv * pres(k,i,j)
1364  enddo
1365  enddo
1366  enddo
1367 
1368  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
integer, public ieb
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:80
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

◆ atmos_saturation_dqsi_dtem_dpre()

subroutine, public scale_atmos_saturation::atmos_saturation_dqsi_dtem_dpre ( real(rp), dimension(ka,ia,ja), intent(out)  dqsdtem,
real(rp), dimension(ka,ia,ja), intent(out)  dqsdpre,
real(rp), dimension (ka,ia,ja), intent(in)  temp,
real(rp), dimension (ka,ia,ja), intent(in)  pres 
)

Definition at line 1375 of file scale_atmos_sub_saturation.F90.

References scale_const::const_lhs0, cpovr_ice, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, and lovr_ice.

Referenced by scale_atmos_phy_mp_sn14::update_by_phase_change_kij().

1375  use scale_const, only: &
1376  lhs0 => const_lhs0
1377  implicit none
1378 
1379  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1380  real(RP), intent(out) :: dqsdpre(ka,ia,ja)
1381  real(RP), intent(in) :: temp (ka,ia,ja)
1382  real(RP), intent(in) :: pres (ka,ia,ja)
1383 
1384  real(RP) :: psat ! saturation vapor pressure
1385  real(RP) :: lhv ! latent heat for condensation
1386 
1387  real(RP) :: den1, den2 ! denominator
1388  real(RP) :: rtem00, tem
1389 
1390  integer :: k, i, j
1391  !---------------------------------------------------------------------------
1392 
1393  rtem00 = 1.0_rp / tem00
1394 
1395  !$omp parallel do private(i,j,k,TEM,psat,den1,den2,lhv) OMP_SCHEDULE_ collapse(2)
1396  do j = jsb, jeb
1397  do i = isb, ieb
1398  do k = ks, ke
1399  tem = max( temp(k,i,j), tem_min )
1400 
1401  psat = psat0 &
1402  * ( tem * rtem00 )**cpovr_ice &
1403  * exp( lovr_ice * ( rtem00 - 1.0_rp/tem ) )
1404 
1405  den1 = ( pres(k,i,j) - (1.0_rp-epsvap) * psat ) &
1406  * ( pres(k,i,j) - (1.0_rp-epsvap) * psat )
1407  den2 = den1 * rvap * temp(k,i,j) * temp(k,i,j)
1408  lhv = lhs0 + ( cpvap-ci ) * ( temp(k,i,j)-tem00 )
1409 
1410  dqsdpre(k,i,j) = - epsvap * psat / den1
1411  dqsdtem(k,i,j) = epsvap * psat / den2 * lhv * pres(k,i,j)
1412  enddo
1413  enddo
1414  enddo
1415 
1416  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
integer, public ieb
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:82
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

Variable Documentation

◆ cpovr_liq

real(rp), public scale_atmos_saturation::cpovr_liq

◆ cpovr_ice

real(rp), public scale_atmos_saturation::cpovr_ice

◆ cvovr_liq

real(rp), public scale_atmos_saturation::cvovr_liq

Definition at line 139 of file scale_atmos_sub_saturation.F90.

Referenced by scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment(), and atmos_saturation_setup().

139  real(RP), public :: cvovr_liq

◆ cvovr_ice

real(rp), public scale_atmos_saturation::cvovr_ice

Definition at line 140 of file scale_atmos_sub_saturation.F90.

Referenced by scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment(), and atmos_saturation_setup().

140  real(RP), public :: cvovr_ice

◆ lovr_liq

real(rp), public scale_atmos_saturation::lovr_liq

◆ lovr_ice

real(rp), public scale_atmos_saturation::lovr_ice