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

module ATMOSPHERE / Thermodynamics More...

Functions/Subroutines

subroutine, public atmos_thermodyn_setup
 Setup. More...
 
subroutine atmos_thermodyn_qd_0d (qdry, q)
 calc dry air mass (0D) More...
 
subroutine, public atmos_thermodyn_tempre (temp, pres, Ein, dens, qdry, q)
 
subroutine, public atmos_thermodyn_tempre2 (temp, pres, dens, pott, qdry, q)
 

Variables

real(rp), public thermodyn_emask = 0.0_RP
 =0: SIMPLE, 1: EXACT More...
 
real(rp), dimension(:), allocatable, public aq_cp
 CP for each hydrometeors [J/kg/K]. More...
 
real(rp), dimension(:), allocatable, public aq_cv
 CV for each hydrometeors [J/kg/K]. More...
 

Detailed Description

module ATMOSPHERE / Thermodynamics

Description
Thermodynamics module
Author
Team SCALE
History
  • 2011-10-24 (T.Seiki) [new] Import from NICAM
  • 2012-02-10 (H.Yashiro) [mod] Reconstruction
  • 2012-03-23 (H.Yashiro) [mod] Explicit index parameter inclusion
  • 2012-12-22 (S.Nishizawa) [mod] Use thermodyn macro set

Function/Subroutine Documentation

◆ atmos_thermodyn_setup()

subroutine, public scale_atmos_thermodyn::atmos_thermodyn_setup ( )

Setup.

Definition at line 158 of file scale_atmos_sub_thermodyn.F90.

References aq_cp, aq_cv, scale_const::const_ci, scale_const::const_cl, scale_const::const_cpvap, scale_const::const_cvvap, scale_const::const_thermodyn_type, scale_tracer::i_qv, scale_stdio::io_fid_log, scale_stdio::io_l, scale_process::prc_mpistop(), scale_tracer::qie, scale_tracer::qis, scale_tracer::qqa, scale_tracer::qwe, scale_tracer::qws, and thermodyn_emask.

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

158  use scale_process, only: &
160  use scale_const, only: &
161  cpvap => const_cpvap, &
162  cvvap => const_cvvap, &
163  cl => const_cl, &
164  ci => const_ci, &
165  thermodyn_type => const_thermodyn_type
166  implicit none
167 
168  integer :: n
169  !---------------------------------------------------------------------------
170 
171  if( io_l ) write(io_fid_log,*)
172  if( io_l ) write(io_fid_log,*) '++++++ Module[THERMODYN] / Categ[ATMOS SHARE] / Origin[SCALElib]'
173 
174  allocate( aq_cp(qqa) )
175  allocate( aq_cv(qqa) )
176 
177  if ( thermodyn_type == 'EXACT' ) then
178  thermodyn_emask = 1.0_rp
179 
180  aq_cp(i_qv) = cpvap
181  aq_cv(i_qv) = cvvap
182 
183  if ( qws /= 0 ) then
184  do n = qws, qwe
185  aq_cp(n) = cl
186  aq_cv(n) = cl
187  enddo
188  endif
189 
190  if ( qis /= 0 ) then
191  do n = qis, qie
192  aq_cp(n) = ci
193  aq_cv(n) = ci
194  enddo
195  endif
196  elseif( thermodyn_type == 'SIMPLE' ) then
197  thermodyn_emask = 0.0_rp
198 
199  aq_cp(i_qv) = cpdry
200  aq_cv(i_qv) = cvdry
201 
202  if ( qws /= 0 ) then
203  do n = qws, qwe
204  aq_cp(n) = cvdry
205  aq_cv(n) = cvdry
206  enddo
207  endif
208 
209  if ( qis /= 0 ) then
210  do n = qis, qie
211  aq_cp(n) = cvdry
212  aq_cv(n) = cvdry
213  enddo
214  endif
215  endif
216 
217  return
integer, public qie
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:69
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public qwe
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:68
integer, public qws
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:67
integer, public qis
integer, public i_qv
integer, public qqa
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:66
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_thermodyn_qd_0d()

subroutine scale_atmos_thermodyn::atmos_thermodyn_qd_0d ( real(rp), intent(out)  qdry,
real(rp), dimension(qa), intent(in)  q 
)

calc dry air mass (0D)

Parameters
[out]qdrydry mass concentration [kg/kg]
[in]qmass concentration [kg/kg]

Definition at line 225 of file scale_atmos_sub_thermodyn.F90.

References aq_cp, aq_cv, scale_tracer::i_qc, scale_tracer::i_qg, scale_tracer::i_qi, scale_tracer::i_qr, scale_tracer::i_qs, scale_tracer::i_qv, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, dc_log::log(), scale_tracer::qqe, scale_tracer::qqs, and thermodyn_emask.

225  implicit none
226 
227  real(RP), intent(out) :: qdry
228  real(RP), intent(in) :: q(qa)
229 
230  integer :: iqw
231  !-----------------------------------------------------------------------------
232 
233  qdry = 1.0_rp
234 #ifndef DRY
235  do iqw = qqs, qqe
236  qdry = qdry - q(iqw)
237  enddo
238 #endif
239 
240  return
integer, public qqe
integer, public qa
integer, public qqs
Here is the call graph for this function:

◆ atmos_thermodyn_tempre()

subroutine, public scale_atmos_thermodyn::atmos_thermodyn_tempre ( real(rp), dimension(ka,ia,ja), intent(out)  temp,
real(rp), dimension(ka,ia,ja), intent(out)  pres,
real(rp), dimension (ka,ia,ja), intent(in)  Ein,
real(rp), dimension(ka,ia,ja), intent(in)  dens,
real(rp), dimension(ka,ia,ja), intent(in)  qdry,
real(rp), dimension (ka,ia,ja,qa), intent(in)  q 
)

Definition at line 1239 of file scale_atmos_sub_thermodyn.F90.

References aq_cv, scale_tracer::i_qv, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, and scale_grid_index::ks.

1239  implicit none
1240 
1241  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
1242  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
1243  real(RP), intent(in) :: ein (ka,ia,ja) ! internal energy
1244  real(RP), intent(in) :: dens(ka,ia,ja) ! density
1245  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry concentration
1246  real(RP), intent(in) :: q (ka,ia,ja,qa) ! water concentration
1247 
1248  real(RP) :: cv, rmoist
1249 
1250  integer :: i, j, k, iqw
1251  !---------------------------------------------------------------------------
1252 
1253  !$omp parallel do private(i,j,k,iqw,cv,Rmoist) OMP_SCHEDULE_ collapse(2)
1254  do j = jsb, jeb
1255  do i = isb, ieb
1256  do k = ks, ke
1257 
1258  calc_cv(cv, qdry(k,i,j), q, k, i, j, iqw, cvdry, aq_cv)
1259  calc_r(rmoist, q(k,i,j,i_qv), qdry(k,i,j), rdry, rvap)
1260 
1261  temp(k,i,j) = ein(k,i,j) / cv
1262 
1263  pres(k,i,j) = dens(k,i,j) * rmoist * temp(k,i,j)
1264 
1265  enddo
1266  enddo
1267  enddo
1268 
1269  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
integer, public qa
integer, public ieb
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
integer, public i_qv
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)

◆ atmos_thermodyn_tempre2()

subroutine, public scale_atmos_thermodyn::atmos_thermodyn_tempre2 ( real(rp), dimension(ka,ia,ja), intent(out)  temp,
real(rp), dimension(ka,ia,ja), intent(out)  pres,
real(rp), dimension(ka,ia,ja), intent(in)  dens,
real(rp), dimension(ka,ia,ja), intent(in)  pott,
real(rp), dimension(ka,ia,ja), intent(in)  qdry,
real(rp), dimension (ka,ia,ja,qa), intent(in)  q 
)

Definition at line 1276 of file scale_atmos_sub_thermodyn.F90.

References aq_cp, scale_tracer::i_qv, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, and scale_grid_index::ks.

1276  implicit none
1277 
1278  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
1279  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
1280  real(RP), intent(in) :: dens(ka,ia,ja) ! density
1281  real(RP), intent(in) :: pott(ka,ia,ja) ! potential temperature
1282  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry concentration
1283  real(RP), intent(in) :: q (ka,ia,ja,qa) ! water concentration
1284 
1285  real(RP) :: rmoist, cp
1286 
1287  integer :: i, j, k, iqw
1288  !---------------------------------------------------------------------------
1289 
1290  !$omp parallel do private(i,j,k,iqw,cp,Rmoist) OMP_SCHEDULE_ collapse(2)
1291  do j = jsb, jeb
1292  do i = isb, ieb
1293  do k = ks, ke
1294 
1295  calc_cp(cp, qdry(k,i,j), q, k, i, j, iqw, cpdry, aq_cp)
1296  calc_r(rmoist, q(k,i,j,i_qv), qdry(k,i,j), rdry, rvap)
1297  calc_pre(pres(k,i,j), dens(k,i,j), pott(k,i,j), rmoist, cp, pre00)
1298 
1299  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rmoist )
1300 
1301  enddo
1302  enddo
1303  enddo
1304 
1305  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
integer, public qa
integer, public ieb
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
integer, public i_qv
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)

Variable Documentation

◆ thermodyn_emask

real(rp), public scale_atmos_thermodyn::thermodyn_emask = 0.0_RP

=0: SIMPLE, 1: EXACT

Definition at line 139 of file scale_atmos_sub_thermodyn.F90.

Referenced by atmos_thermodyn_qd_0d(), and atmos_thermodyn_setup().

139  real(RP), public :: thermodyn_emask = 0.0_rp

◆ aq_cp

real(rp), dimension(:), allocatable, public scale_atmos_thermodyn::aq_cp

◆ aq_cv

real(rp), dimension(:), allocatable, public scale_atmos_thermodyn::aq_cv