SCALE-RM
Data Types | Functions/Subroutines
scale_atmos_hydrostatic Module Reference

module atmosphere / hydrostatic barance More...

Functions/Subroutines

subroutine, public atmos_hydrostatic_setup
 Setup. More...
 
subroutine atmos_hydrostatic_buildrho_1d (KA, KS, KE, pott, qv, qc, pres_sfc, pott_sfc, qv_sfc, qc_sfc, cz, fz, dens, temp, pres, temp_sfc)
 Build up density from surface (1D) More...
 
subroutine atmos_hydrostatic_buildrho_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, pott, qv, qc, pres_sfc, pott_sfc, qv_sfc, qc_sfc, cz, fz, area, dens, temp, pres, temp_sfc)
 Build up density from surface (3D) More...
 
subroutine atmos_hydrostatic_buildrho_real_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, pott, qv, qc, cz, pres, dens, temp)
 Build up density from surface (3D), not to reverse from TOA. More...
 
subroutine atmos_hydrostatic_buildrho_atmos_0d (pott_L2, qv_L2, qc_L2, dens_L1, pott_L1, qv_L1, qc_L1, dz, k, dens_L2, temp_L2, pres_L2)
 Build up density (0D) More...
 
subroutine atmos_hydrostatic_buildrho_atmos_1d (KA, KS, KE, pott, qv, qc, dz, dens, temp, pres, kref)
 Build up density from lowermost atmosphere (1D) More...
 
subroutine atmos_hydrostatic_buildrho_atmos_rev_1d (KA, KS, KE, pott, qv, qc, dz, dens, temp, pres, kref)
 Build up density from upermost atmosphere (1D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, temp, qv, qc, pres_sfc, temp_sfc, qv_sfc, qc_sfc, cz, fz, dens, pott, pres, pott_sfc)
 Build up density from surface (3D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_rev_1d (KA, KS, KE, temp, qv, qc, dz, dens, pott, pres)
 Build up density from upermost atmosphere (1D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, temp, qv, qc, dz, dens, pott, pres)
 Build up density from lowermost atmosphere (3D) More...
 

Detailed Description

module atmosphere / hydrostatic barance

Description
make hydrostatic profile in the model
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_HYDROSTATIC
    nametypedefault valuecomment
    HYDROSTATIC_USELAPSERATE logical .false. use lapse rate?
    HYDROSTATIC_BUILDRHO_REAL_KREF integer 1
    HYDROSTATIC_BAROMETRIC_LAW_MSLP_KREF integer 1 reference layer for MSLP calculation

History Output
No history output

Function/Subroutine Documentation

◆ atmos_hydrostatic_setup()

subroutine, public scale_atmos_hydrostatic::atmos_hydrostatic_setup ( )

Setup.

Definition at line 116 of file scale_atmos_hydrostatic.F90.

References scale_const::const_eps, scale_io::io_fid_conf, and scale_prc::prc_abort().

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

116  use scale_prc, only: &
117  prc_abort
118  use scale_const, only: &
119  const_eps
120  implicit none
121 
122  namelist / param_atmos_hydrostatic / &
123  hydrostatic_uselapserate, &
124  hydrostatic_buildrho_real_kref, &
125  hydrostatic_barometric_law_mslp_kref
126 
127  integer :: ierr
128  !---------------------------------------------------------------------------
129 
130  log_newline
131  log_info("ATMOS_HYDROSTATIC_setup",*) 'Setup'
132 
133  !--- read namelist
134  rewind(io_fid_conf)
135  read(io_fid_conf,nml=param_atmos_hydrostatic,iostat=ierr)
136  if( ierr < 0 ) then !--- missing
137  log_info("ATMOS_HYDROSTATIC_setup",*) 'Not found namelist. Default used.'
138  elseif( ierr > 0 ) then !--- fatal error
139  log_error("ATMOS_HYDROSTATIC_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_HYDROSTATIC. Check!'
140  call prc_abort
141  endif
142  log_nml(param_atmos_hydrostatic)
143 
144  criteria = sqrt( const_eps )
145 
146  log_newline
147  log_info("ATMOS_HYDROSTATIC_setup",*) 'Use lapse rate for estimation of surface temperature? : ', hydrostatic_uselapserate
148  log_info("ATMOS_HYDROSTATIC_setup",*) 'Buildrho conversion criteria : ', criteria
149 
150  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_1d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_1d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), intent(in)  pott,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc,
real(rp), intent(in)  pres_sfc,
real(rp), intent(in)  pott_sfc,
real(rp), intent(in)  qv_sfc,
real(rp), intent(in)  qc_sfc,
real(rp), dimension (ka), intent(in)  cz,
real(rp), dimension (0:ka), intent(in)  fz,
real(rp), dimension(ka), intent(out)  dens,
real(rp), dimension(ka), intent(out)  temp,
real(rp), dimension(ka), intent(out)  pres,
real(rp), intent(out)  temp_sfc 
)

Build up density from surface (1D)

Parameters
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in]pres_sfcsurface pressure [Pa]
[in]pott_sfcsurface potential temperature [K]
[in]qv_sfcsurface water vapor [kg/kg]
[in]qc_sfcsurface liquid water [kg/kg]
[out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]
[out]temp_sfcsurface temperature [K]

Definition at line 163 of file scale_atmos_hydrostatic.F90.

References atmos_hydrostatic_buildrho_atmos_0d(), atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_prc::prc_abort().

Referenced by atmos_hydrostatic_buildrho_3d().

163  use scale_prc, only: &
164  prc_abort
165  use scale_atmos_hydrometeor, only: &
166  cv_vapor, &
167  cv_water, &
168  cp_vapor, &
169  cp_water
170  implicit none
171 
172  integer, intent(in) :: ka, ks, ke
173  real(RP), intent(in) :: pott(ka)
174  real(RP), intent(in) :: qv (ka)
175  real(RP), intent(in) :: qc (ka)
176  real(RP), intent(in) :: pres_sfc
177  real(RP), intent(in) :: pott_sfc
178  real(RP), intent(in) :: qv_sfc
179  real(RP), intent(in) :: qc_sfc
180  real(RP), intent(in) :: cz (ka)
181  real(RP), intent(in) :: fz (0:ka)
182  real(RP), intent(out) :: dens(ka)
183  real(RP), intent(out) :: temp(ka)
184  real(RP), intent(out) :: pres(ka)
185  real(RP), intent(out) :: temp_sfc
186 
187  real(RP) :: dens_sfc
188  real(RP) :: rtot_sfc
189  real(RP) :: cvtot_sfc
190  real(RP) :: cptot_sfc
191  real(RP) :: cpovcv_sfc
192  real(RP) :: rtot
193  real(RP) :: cvtot
194  real(RP) :: cptot
195  real(RP) :: cpovcv
196  real(RP) :: dz(ka)
197 
198  integer :: k
199  !---------------------------------------------------------------------------
200 
201  !--- from surface to lowermost atmosphere
202 
203  rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
204  + rvap * qv_sfc
205  cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
206  + cv_vapor * qv_sfc &
207  + cv_water * qc_sfc
208  cptot_sfc = cpdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
209  + cp_vapor * qv_sfc &
210  + cp_water * qc_sfc
211  cpovcv_sfc = cptot_sfc / cvtot_sfc
212 
213  rtot = rdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
214  + rvap * qv(ks)
215  cvtot = cvdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
216  + cv_vapor * qv(ks) &
217  + cv_water * qc(ks)
218  cptot = cpdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
219  + cp_vapor * qv(ks) &
220  + cp_water * qc(ks)
221  cpovcv = cptot / cvtot
222 
223  ! density at surface
224  dens_sfc = p00 / rtot_sfc / pott_sfc * ( pres_sfc/p00 )**(cvtot_sfc/cptot_sfc)
225  temp_sfc = pres_sfc / ( dens_sfc * rtot_sfc )
226 
227  dz(ks-1) = cz(ks) - fz(ks-1)
228  do k = ks, ke-1
229  dz(k) = cz(k+1) - cz(k)
230  end do
231 
232  ! make density at lowermost cell center
233  if ( hydrostatic_uselapserate ) then
234 
235  temp(ks) = pott_sfc - lapsdry * dz(ks-1) ! use dry lapse rate
236  pres(ks) = p00 * ( temp(ks)/pott(ks) )**(cptot/rtot)
237  dens(ks) = p00 / rtot / pott(ks) * ( pres(ks)/p00 )**(cvtot/cptot)
238 
239  else ! use itelation
240 
241  call atmos_hydrostatic_buildrho_atmos_0d( pott(ks), qv(ks), qc(ks), & ! [IN]
242  dens_sfc, pott_sfc, qv_sfc, qc_sfc, & ! [IN]
243  dz(ks-1), ks-1, & ! [IN]
244  dens(ks), temp(ks), pres(ks) ) ! [OUT]
245 
246  endif
247 
248  !--- from lowermost atmosphere to top of atmosphere
249  call atmos_hydrostatic_buildrho_atmos_1d( ka, ks, ke, &
250  pott(:), qv(:), qc(:), & ! [IN]
251  dz(:), & ! [IN]
252  dens(:), & ! [INOUT]
253  temp(:), pres(:) ) ! [OUT]
254  ! fill dummy
255  dens( 1:ks-1) = 0.0_rp
256  dens(ke+1:ka ) = 0.0_rp
257 
258  return
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public ks
start point of inner domain: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public ka
of whole cells: z, local, with HALO
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_3d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(in)  pott,
real(rp), dimension (ka,ia,ja), intent(in)  qv,
real(rp), dimension (ka,ia,ja), intent(in)  qc,
real(rp), dimension(ia,ja), intent(in)  pres_sfc,
real(rp), dimension(ia,ja), intent(in)  pott_sfc,
real(rp), dimension (ia,ja), intent(in)  qv_sfc,
real(rp), dimension (ia,ja), intent(in)  qc_sfc,
real(rp), dimension( ka,ia,ja), intent(in)  cz,
real(rp), dimension(0:ka,ia,ja), intent(in)  fz,
real(rp), dimension(ia,ja), intent(in)  area,
real(rp), dimension(ka,ia,ja), intent(out)  dens,
real(rp), dimension(ka,ia,ja), intent(out)  temp,
real(rp), dimension(ka,ia,ja), intent(out)  pres,
real(rp), dimension(ia,ja), intent(out)  temp_sfc 
)

Build up density from surface (3D)

Parameters
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in]pres_sfcsurface pressure [Pa]
[in]pott_sfcsurface potential temperature [K]
[in]qv_sfcsurface water vapor [kg/kg]
[in]qc_sfcsurface liquid water [kg/kg]
[out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]
[out]temp_sfcsurface temperature [K]

Definition at line 271 of file scale_atmos_hydrostatic.F90.

References atmos_hydrostatic_buildrho_1d().

271  use scale_statistics, only: &
272  statistics_horizontal_mean
273  implicit none
274 
275  integer, intent(in) :: ka, ks, ke
276  integer, intent(in) :: ia, is, ie
277  integer, intent(in) :: ja, js, je
278 
279  real(RP), intent(in) :: pott(ka,ia,ja)
280  real(RP), intent(in) :: qv (ka,ia,ja)
281  real(RP), intent(in) :: qc (ka,ia,ja)
282  real(RP), intent(in) :: pres_sfc(ia,ja)
283  real(RP), intent(in) :: pott_sfc(ia,ja)
284  real(RP), intent(in) :: qv_sfc (ia,ja)
285  real(RP), intent(in) :: qc_sfc (ia,ja)
286  real(RP), intent(in) :: cz( ka,ia,ja)
287  real(RP), intent(in) :: fz(0:ka,ia,ja)
288  real(RP), intent(in) :: area(ia,ja)
289 
290  real(RP), intent(out) :: dens(ka,ia,ja)
291  real(RP), intent(out) :: temp(ka,ia,ja)
292  real(RP), intent(out) :: pres(ka,ia,ja)
293  real(RP), intent(out) :: temp_sfc(ia,ja)
294 
295  real(RP) :: dz(ka,ia,ja), dz_top(ia,ja)
296 
297  ! TOA
298  real(RP) :: pott_toa(ia,ja)
299  real(RP) :: qv_toa (ia,ja)
300  real(RP) :: qc_toa (ia,ja)
301  real(RP) :: dens_toa(ia,ja)
302  real(RP) :: temp_toa(ia,ja)
303  real(RP) :: pres_toa(ia,ja)
304 
305  real(RP) :: dens_mean
306 
307  integer :: k, i, j
308  !---------------------------------------------------------------------------
309 
310  !--- from surface to lowermost atmosphere
311  !$omp parallel do OMP_SCHEDULE_ collapse(2)
312  do j = js, je
313  do i = is, ie
314  call atmos_hydrostatic_buildrho_1d( ka, ks, ke, &
315  pott(:,i,j), qv(:,i,j), qc(:,i,j), & ! [IN]
316  pres_sfc(i,j), pott_sfc(i,j), qv_sfc(i,j), qc_sfc(i,j), & ! [IN]
317  cz(:,i,j), fz(:,i,j), & ! [IN]
318  dens(:,i,j), temp(:,i,j), pres(:,i,j), & ! [OUT]
319  temp_sfc(i,j) ) ! [OUT]
320  end do
321  end do
322 
323  !$omp parallel do OMP_SCHEDULE_ collapse(2)
324  do j = js, je
325  do i = is, ie
326  do k = ks, ke-1
327  dz(k,i,j) = cz(k+1,i,j) - cz(k,i,j)
328  end do
329  dz_top(i,j) = fz(ke,i,j) - cz(ke,i,j) ! distance from cell center to TOA
330 
331  ! value at TOA
332  pott_toa(i,j) = pott(ke,i,j)
333  qv_toa(i,j) = qv(ke,i,j)
334  qc_toa(i,j) = qc(ke,i,j)
335  end do
336  end do
337 
338  call atmos_hydrostatic_buildrho_atmos_2d( ia, is, ie, ja, js, je, &
339  pott_toa(:,:), qv_toa(:,:), qc_toa(:,:), & ! [IN]
340  dens(ke,:,:), pott(ke,:,:), qv(ke,:,:), qc(ke,:,:), & ! [IN]
341  dz_top(:,:), ke+1, & ! [IN]
342  dens_toa(:,:), temp_toa(:,:), pres_toa(:,:) ) ! [OUT]
343 
344  call statistics_horizontal_mean( ia, is, ie, ja, js, je, &
345  dens_toa(:,:), area(:,:), & ! [IN]
346  dens_mean ) ! [OUT]
347 
348  !$omp parallel do OMP_SCHEDULE_ collapse(2)
349  do j = js, je
350  do i = is, ie
351  dens_toa(i,j) = dens_mean
352  enddo
353  enddo
354 
355  call atmos_hydrostatic_buildrho_atmos_rev_2d( ia, is, ie, ja, js, je, &
356  pott(ke,:,:), qv(ke,:,:), qc(ke,:,:), & ! [IN]
357  dens_toa(:,:), pott_toa(:,:), qv_toa(:,:), qc_toa(:,:), & ! [IN]
358  dz_top(:,:), ke+1, & ! [IN]
359  dens(ke,:,:), temp(ke,:,:), pres(ke,:,:) ) ! [OUT]
360 
361  !--- from top of atmosphere to lowermost atmosphere
362  call atmos_hydrostatic_buildrho_atmos_rev_3d( ka, ks, ke, ia, is, ie, ja, js, je, &
363  pott(:,:,:), qv(:,:,:), qc(:,:,:), & ! [IN]
364  dz(:,:,:), & ! [IN]
365  dens(:,:,:), & ! [INOUT]
366  temp(:,:,:), pres(:,:,:) ) ! [OUT]
367 
368  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
module Statistics
Here is the call graph for this function:

◆ atmos_hydrostatic_buildrho_real_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_real_3d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(in)  pott,
real(rp), dimension (ka,ia,ja), intent(in)  qv,
real(rp), dimension (ka,ia,ja), intent(in)  qc,
real(rp), dimension (ka,ia,ja), intent(in)  cz,
real(rp), dimension(ka,ia,ja), intent(inout)  pres,
real(rp), dimension(ka,ia,ja), intent(out)  dens,
real(rp), dimension(ka,ia,ja), intent(out)  temp 
)

Build up density from surface (3D), not to reverse from TOA.

Parameters
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in,out]prespressure [Pa]
[out]densdensity [kg/m3]
[out]temptemperature [K]

Definition at line 379 of file scale_atmos_hydrostatic.F90.

References scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, and scale_atmos_hydrometeor::cv_water.

379  use scale_atmos_hydrometeor, only: &
380  cv_vapor, &
381  cv_water, &
382  cp_vapor, &
383  cp_water
384  implicit none
385  integer, intent(in) :: ka, ks, ke
386  integer, intent(in) :: ia, is, ie
387  integer, intent(in) :: ja, js, je
388 
389  real(RP), intent(in) :: pott(ka,ia,ja)
390  real(RP), intent(in) :: qv (ka,ia,ja)
391  real(RP), intent(in) :: qc (ka,ia,ja)
392  real(RP), intent(in) :: cz (ka,ia,ja)
393 
394  real(RP), intent(inout) :: pres(ka,ia,ja)
395 
396  real(RP), intent(out) :: dens(ka,ia,ja)
397  real(RP), intent(out) :: temp(ka,ia,ja)
398 
399  real(RP) :: dz(ka,ia,ja)
400 
401  real(RP) :: pott_toa(ia,ja)
402  real(RP) :: qv_toa (ia,ja)
403  real(RP) :: qc_toa (ia,ja)
404 
405  real(RP) :: rtot
406  real(RP) :: cvtot
407  real(RP) :: cptot
408 
409  integer :: kref
410  integer :: k, i, j
411  !---------------------------------------------------------------------------
412 
413  !--- from surface to lowermost atmosphere
414 
415  !$omp parallel do OMP_SCHEDULE_ collapse(2)
416  do j = js, je
417  do i = is, ie
418  do k = ks, ke-1
419  dz(k,i,j) = cz(k+1,i,j) - cz(k,i,j) ! distance from cell center to cell center
420  enddo
421  enddo
422  enddo
423 
424  !$omp parallel do OMP_SCHEDULE_ collapse(2)
425  do j = js, je
426  do i = is, ie
427  pott_toa(i,j) = pott(ke,i,j)
428  qv_toa(i,j) = qv(ke,i,j)
429  qc_toa(i,j) = qc(ke,i,j)
430  enddo
431  enddo
432 
433  kref = hydrostatic_buildrho_real_kref + ks - 1
434 
435  ! calc density at reference level
436  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
437  !$omp private(Rtot,CVtot,CPtot)
438  do j = js, je
439  do i = is, ie
440  rtot = rdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
441  + rvap * qv(kref,i,j)
442  cvtot = cvdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
443  + cv_vapor * qv(kref,i,j) &
444  + cv_water * qc(kref,i,j)
445  cptot = cpdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
446  + cp_vapor * qv(kref,i,j) &
447  + cp_water * qc(kref,i,j)
448  dens(kref,i,j) = p00 / ( rtot * pott(kref,i,j) ) * ( pres(kref,i,j)/p00 )**(cvtot/cptot)
449  enddo
450  enddo
451 
452  !--- from lowermost atmosphere to top of atmosphere
453  call atmos_hydrostatic_buildrho_atmos_3d( ka, ks, ke, ia, is, ie, ja, js, je, &
454  pott(:,:,:), qv(:,:,:), qc(:,:,:), & ! [IN]
455  dz(:,:,:), & ! [IN]
456  dens(:,:,:), & ! [INOUT]
457  temp(:,:,:), pres(:,:,:), & ! [OUT]
458  kref = kref ) ! [IN]
459  call atmos_hydrostatic_buildrho_atmos_rev_3d( ka, ks, ke, ia, is, ie, ja, js, je, &
460  pott(:,:,:), qv(:,:,:), qc(:,:,:), & ! [IN]
461  dz(:,:,:), & ! [IN]
462  dens(:,:,:), & ! [INOUT]
463  temp(:,:,:), pres(:,:,:), & ! [OUT]
464  kref = kref ) ! [IN]
465 
466 !!$ call ATMOS_HYDROSTATIC_buildrho_atmos_2D( IA, IS, IE, JA, JS, JE, &
467 !!$ pott_toa(:,:), qv_toa(:,:), qc_toa(:,:), & ! [IN]
468 !!$ dens(KE,:,:), & ! [IN]
469 !!$ pott(KE,:,:), qv(KE,:,:), qc(KE,:,:), & ! [IN]
470 !!$ dz(KE+1,:,:), KE+1, & ! [IN]
471 !!$ dens_toa(:,:), temp_toa(:,:), pres_toa(:,:) ) ! [OUT]
472 
473  ! density at TOA
474  dens( 1:ks-1,:,:) = 0.0_rp ! fill dummy
475 !!$ dens(KE+2:KA ,:,:) = 0.0_RP ! fill dummy
476  dens(ke+1:ka ,:,:) = 0.0_rp ! fill dummy
477 
478  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].

◆ atmos_hydrostatic_buildrho_atmos_0d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_0d ( real(rp), intent(in)  pott_L2,
real(rp), intent(in)  qv_L2,
real(rp), intent(in)  qc_L2,
real(rp), intent(in)  dens_L1,
real(rp), intent(in)  pott_L1,
real(rp), intent(in)  qv_L1,
real(rp), intent(in)  qc_L1,
real(rp), intent(in)  dz,
integer, intent(in)  k,
real(rp), intent(out)  dens_L2,
real(rp), intent(out)  temp_L2,
real(rp), intent(out)  pres_L2 
)

Build up density (0D)

Parameters
[in]pott_l2potential temperature at layer 2 [K]
[in]qv_l2water vapor at layer 2 [kg/kg]
[in]qc_l2liquid water at layer 2 [kg/kg]
[in]dens_l1density at layer 1 [Pa]
[in]pott_l1potential temperature at layer 1 [K]
[in]qv_l1water vapor at layer 1 [kg/kg]
[in]qc_l1liquid water at layer 1 [kg/kg]
[in]dzdistance from layer 1 to layer 2 [m]
[in]kfor monitor
[out]dens_l2density at layer 2 [kg/m3]
[out]temp_l2temperature at layer 2 [K]
[out]pres_l2pressure at layer 2 [Pa]

Definition at line 490 of file scale_atmos_hydrostatic.F90.

References scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_prc::prc_abort().

Referenced by atmos_hydrostatic_buildrho_1d(), and atmos_hydrostatic_buildrho_atmos_rev_1d().

490  use scale_prc, only: &
491  prc_abort
492  use scale_atmos_hydrometeor, only: &
493  cv_vapor, &
494  cv_water, &
495  cp_vapor, &
496  cp_water
497  implicit none
498  real(RP), intent(in) :: pott_l2
499  real(RP), intent(in) :: qv_l2
500  real(RP), intent(in) :: qc_l2
501  real(RP), intent(in) :: dens_l1
502  real(RP), intent(in) :: pott_l1
503  real(RP), intent(in) :: qv_l1
504  real(RP), intent(in) :: qc_l1
505  real(RP), intent(in) :: dz
506  integer, intent(in) :: k
507 
508  real(RP), intent(out) :: dens_l2
509  real(RP), intent(out) :: temp_l2
510  real(RP), intent(out) :: pres_l2
511 
512  real(RP) :: rtot_l1 , rtot_l2
513  real(RP) :: cvtot_l1 , cvtot_l2
514  real(RP) :: cptot_l1 , cptot_l2
515  real(RP) :: cpovcv_l1, cpovcv_l2
516 
517  real(RP) :: pres_l1
518  real(RP) :: dens_s, dhyd, dgrd
519  integer :: ite
520  logical :: converged
521  !---------------------------------------------------------------------------
522 
523  rtot_l1 = rdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
524  + rvap * qv_l1
525  cvtot_l1 = cvdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
526  + cv_vapor * qv_l1 &
527  + cv_water * qc_l1
528  cptot_l1 = cpdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
529  + cp_vapor * qv_l1 &
530  + cp_water * qc_l1
531  cpovcv_l1 = cptot_l1 / cvtot_l1
532 
533  rtot_l2 = rdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
534  + rvap * qv_l2
535  cvtot_l2 = cvdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
536  + cv_vapor * qv_l2 &
537  + cv_water * qc_l2
538  cptot_l2 = cpdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
539  + cp_vapor * qv_l2 &
540  + cp_water * qc_l2
541  cpovcv_l2 = cptot_l2 / cvtot_l2
542 
543  dens_s = 0.0_rp
544  dens_l2 = dens_l1 ! first guess
545 
546  pres_l1 = p00 * ( dens_l1 * rtot_l1 * pott_l1 / p00 )**cpovcv_l1
547 
548  converged = .false.
549  do ite = 1, itelim
550  if ( abs(dens_l2-dens_s) <= criteria ) then
551  converged = .true.
552  exit
553  endif
554 
555  dens_s = dens_l2
556  pres_l2 = p00 * ( dens_s * rtot_l2 * pott_l2 / p00 )**cpovcv_l2
557 
558  dhyd = + ( pres_l1 - pres_l2 ) / dz & ! dp/dz
559  - grav * 0.5_rp * ( dens_l1 + dens_s ) ! rho*g
560 
561  dgrd = - pres_l2 * cpovcv_l2 / dens_s / dz &
562  - grav * 0.5_rp
563 
564  dens_l2 = dens_s - dhyd/dgrd
565 
566  if ( dens_l2*0.0_rp /= 0.0_rp ) exit
567  enddo
568 
569  if ( .NOT. converged ) then
570  log_error("ATMOS_HYDROSTATIC_buildrho_atmos_0D",*) 'iteration not converged!', &
571  k,dens_l2,ite,dens_s,dhyd,dgrd
572  call prc_abort
573  endif
574 
575  pres_l2 = p00 * ( dens_l2 * rtot_l2 * pott_l2 / p00 )**cpovcv_l2
576  temp_l2 = pres_l2 / ( dens_l2 * rtot_l2 )
577 
578  return
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_atmos_1d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), intent(in)  pott,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc,
real(rp), dimension (ka), intent(in)  dz,
real(rp), dimension(ka), intent(inout)  dens,
real(rp), dimension(ka), intent(out)  temp,
real(rp), dimension(ka), intent(out)  pres,
integer, intent(in), optional  kref 
)

Build up density from lowermost atmosphere (1D)

Parameters
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in,out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]

Definition at line 591 of file scale_atmos_hydrostatic.F90.

References scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_prc::prc_abort().

Referenced by atmos_hydrostatic_buildrho_1d(), and atmos_hydrostatic_buildrho_atmos_rev_1d().

591  use scale_prc, only: &
592  prc_abort
593  use scale_atmos_hydrometeor, only: &
594  cv_vapor, &
595  cv_water, &
596  cp_vapor, &
597  cp_water
598  implicit none
599  integer, intent(in) :: ka, ks, ke
600 
601  real(RP), intent(in) :: pott(ka)
602  real(RP), intent(in) :: qv (ka)
603  real(RP), intent(in) :: qc (ka)
604  real(RP), intent(in) :: dz (ka)
605 
606  real(RP), intent(inout) :: dens(ka)
607 
608  real(RP), intent(out) :: temp(ka)
609  real(RP), intent(out) :: pres(ka)
610 
611  integer, intent(in), optional :: kref
612 
613  real(RP) :: rtot (ka)
614  real(RP) :: cpovcv(ka)
615  real(RP) :: cvtot
616  real(RP) :: cptot
617 
618  real(RP) :: dens_s, dhyd, dgrd
619  integer :: ite
620  logical :: converged
621 
622  integer :: kref_
623  integer :: k
624  !---------------------------------------------------------------------------
625 
626  if ( present(kref) ) then
627  kref_ = kref
628  else
629  kref_ = ks
630  end if
631 
632  do k = kref_, ke
633  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
634  + rvap * qv(k)
635  cvtot = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
636  + cv_vapor * qv(k) &
637  + cv_water * qc(k)
638  cptot = cpdry * ( 1.0_rp - qv(k) - qc(k) ) &
639  + cp_vapor * qv(k) &
640  + cp_water * qc(k)
641  cpovcv(k) = cptot / cvtot
642  enddo
643 
644  pres(kref_) = p00 * ( dens(kref_) * rtot(kref_) * pott(kref_) / p00 )**cpovcv(kref_)
645 
646  do k = kref_+1, ke
647 
648  dens_s = 0.0_rp
649  dens(k) = dens(k-1) ! first guess
650 
651  converged = .false.
652  do ite = 1, itelim
653  if ( abs(dens(k)-dens_s) <= criteria ) then
654  converged = .true.
655  exit
656  endif
657 
658  dens_s = dens(k)
659  pres(k) = p00 * ( dens_s * rtot(k) * pott(k) / p00 )**cpovcv(k)
660 
661  dhyd = + ( pres(k-1) - pres(k) ) / dz(k-1) & ! dpdz
662  - grav * 0.5_rp * ( dens(k-1) + dens_s ) ! rho*g
663 
664  dgrd = - cpovcv(k) * pres(k) / dens_s / dz(k-1) &
665  - grav * 0.5_rp
666 
667  dens(k) = dens_s - dhyd/dgrd
668 
669  if ( dens(k)*0.0_rp /= 0.0_rp ) exit
670  enddo
671 
672  if ( .NOT. converged ) then
673  log_error("ATMOS_HYDROSTATIC_buildrho_atmos_1D",*) 'iteration not converged!', &
674  k,dens(k),ite,dens_s,dhyd,dgrd
675  call prc_abort
676  endif
677 
678  pres(k) = p00 * ( dens(k) * rtot(k) * pott(k) / p00 )**cpovcv(k)
679 
680  enddo
681 
682  do k = kref_, ke
683  temp(k) = pres(k) / ( dens(k) * rtot(k) )
684  enddo
685 
686  dens( 1:ks-1) = dens(ks)
687  dens(ke+1:ka ) = dens(ke)
688  pres( 1:ks-1) = pres(ks)
689  pres(ke+1:ka ) = pres(ke)
690  temp( 1:ks-1) = temp(ks)
691  temp(ke+1:ka ) = temp(ke)
692 
693  return
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public ks
start point of inner domain: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public ka
of whole cells: z, local, with HALO
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_atmos_rev_1d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_1d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), intent(in)  pott,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc,
real(rp), dimension (ka), intent(in)  dz,
real(rp), dimension(ka), intent(inout)  dens,
real(rp), dimension(ka), intent(out)  temp,
real(rp), dimension(ka), intent(out)  pres,
integer, intent(in), optional  kref 
)

Build up density from upermost atmosphere (1D)

Parameters
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in,out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]

Definition at line 706 of file scale_atmos_hydrostatic.F90.

References atmos_hydrostatic_buildrho_atmos_0d(), atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_prc::prc_abort().

706  use scale_prc, only: &
707  prc_abort
708  use scale_atmos_hydrometeor, only: &
709  cv_vapor, &
710  cv_water, &
711  cp_vapor, &
712  cp_water
713  implicit none
714  integer, intent(in) :: ka, ks, ke
715 
716  real(RP), intent(in) :: pott(ka)
717  real(RP), intent(in) :: qv (ka)
718  real(RP), intent(in) :: qc (ka)
719  real(RP), intent(in) :: dz (ka)
720 
721  real(RP), intent(inout) :: dens(ka)
722 
723  real(RP), intent(out) :: temp(ka)
724  real(RP), intent(out) :: pres(ka)
725 
726  integer, intent(in), optional :: kref
727 
728  real(RP) :: rtot (ka)
729  real(RP) :: cpovcv(ka)
730  real(RP) :: cvtot
731  real(RP) :: cptot
732 
733  real(RP) :: dens_s, dhyd, dgrd
734  integer :: ite
735  logical :: converged
736 
737  integer :: kref_
738  integer :: k
739  !---------------------------------------------------------------------------
740 
741  if ( present(kref) ) then
742  kref_ = kref
743  else
744  kref_ = ke
745  end if
746 
747  do k = ks, kref_
748  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
749  + rvap * qv(k)
750  cvtot = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
751  + cv_vapor * qv(k) &
752  + cv_water * qc(k)
753  cptot = cpdry * ( 1.0_rp - qv(k) - qc(k) ) &
754  + cp_vapor * qv(k) &
755  + cp_water * qc(k)
756  cpovcv(k) = cptot / cvtot
757  enddo
758 
759  pres(kref_) = p00 * ( dens(kref_) * rtot(kref_) * pott(kref_) / p00 )**cpovcv(kref_)
760 
761  do k = kref_-1, ks, -1
762 
763  dens_s = 0.0_rp
764  dens(k) = dens(k+1) ! first guess
765 
766  converged = .false.
767  do ite = 1, itelim
768  if ( abs(dens(k)-dens_s) <= criteria ) then
769  converged = .true.
770  exit
771  endif
772 
773  dens_s = dens(k)
774  pres(k) = p00 * ( dens_s * rtot(k) * pott(k) / p00 )**cpovcv(k)
775 
776  dhyd = - ( pres(k+1) - pres(k) ) / dz(k) & ! dpdz
777  - grav * 0.5_rp * ( dens(k+1) + dens_s ) ! rho*g
778 
779  dgrd = + cpovcv(k) * pres(k) / dens_s / dz(k) &
780  - grav * 0.5_rp
781 
782  dens(k) = dens_s - dhyd/dgrd
783 
784  if ( dens(k)*0.0_rp /= 0.0_rp ) exit
785  enddo
786 
787  if ( .NOT. converged ) then
788  log_error("ATMOS_HYDROSTATIC_buildrho_atmos_rev_1D",*) 'iteration not converged!', &
789  k,dens(k),ite,dens_s,dhyd,dgrd
790  call prc_abort
791  endif
792 
793  pres(k) = p00 * ( dens(k) * rtot(k) * pott(k) / p00 )**cpovcv(k)
794 
795  enddo
796 
797  do k = ks, kref_
798  temp(k) = pres(k) / ( dens(k) * rtot(k) )
799  enddo
800 
801  dens( 1:ks-1) = dens(ks)
802  dens(ke+1:ka ) = dens(ke)
803  pres( 1:ks-1) = pres(ks)
804  pres(ke+1:ka ) = pres(ke)
805  temp( 1:ks-1) = temp(ks)
806  temp(ke+1:ka ) = temp(ke)
807 
808  return
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public ks
start point of inner domain: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public ka
of whole cells: z, local, with HALO
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
Here is the call graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_3d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(in)  temp,
real(rp), dimension (ka,ia,ja), intent(in)  qv,
real(rp), dimension (ka,ia,ja), intent(in)  qc,
real(rp), dimension(ia,ja), intent(in)  pres_sfc,
real(rp), dimension(ia,ja), intent(in)  temp_sfc,
real(rp), dimension (ia,ja), intent(in)  qv_sfc,
real(rp), dimension (ia,ja), intent(in)  qc_sfc,
real(rp), dimension(ka,ia,ja), intent(in)  cz,
real(rp), dimension(ka,ia,ja), intent(in)  fz,
real(rp), dimension(ka,ia,ja), intent(out)  dens,
real(rp), dimension(ka,ia,ja), intent(out)  pott,
real(rp), dimension(ka,ia,ja), intent(out)  pres,
real(rp), dimension(ia,ja), intent(out)  pott_sfc 
)

Build up density from surface (3D)

Parameters
[in]temptemperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in]pres_sfcsurface pressure [Pa]
[in]temp_sfcsurface temperature [K]
[in]qv_sfcsurface water vapor [kg/kg]
[in]qc_sfcsurface liquid water [kg/kg]
[out]densdensity [kg/m3]
[out]pottpotential temperature [K]
[out]prespressure [Pa]
[out]pott_sfcsurface potential temperature [K]

Definition at line 1114 of file scale_atmos_hydrostatic.F90.

References scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_prc::prc_abort().

1114  implicit none
1115  integer, intent(in) :: ka, ks, ke
1116  integer, intent(in) :: ia, is, ie
1117  integer, intent(in) :: ja, js, je
1118 
1119  real(RP), intent(in) :: temp(ka,ia,ja)
1120  real(RP), intent(in) :: qv (ka,ia,ja)
1121  real(RP), intent(in) :: qc (ka,ia,ja)
1122  real(RP), intent(in) :: pres_sfc(ia,ja)
1123  real(RP), intent(in) :: temp_sfc(ia,ja)
1124  real(RP), intent(in) :: qv_sfc (ia,ja)
1125  real(RP), intent(in) :: qc_sfc (ia,ja)
1126  real(RP), intent(in) :: cz(ka,ia,ja)
1127  real(RP), intent(in) :: fz(ka,ia,ja)
1128 
1129  real(RP), intent(out) :: dens(ka,ia,ja)
1130  real(RP), intent(out) :: pott(ka,ia,ja)
1131  real(RP), intent(out) :: pres(ka,ia,ja)
1132  real(RP), intent(out) :: pott_sfc(ia,ja)
1133 
1134  integer :: i, j
1135  !---------------------------------------------------------------------------
1136 
1137  !--- from surface to lowermost atmosphere
1138  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1139  do j = js, je
1140  do i = is, ie
1141  call atmos_hydrostatic_buildrho_bytemp_1d( ka, ks, ke, &
1142  temp(:,i,j), qv(:,i,j), qc(:,i,j), & ! [IN]
1143  pres_sfc(i,j), temp_sfc(i,j), qv_sfc(i,j), qc_sfc(i,j), & ! [IN]
1144  cz(:,i,j), fz(:,i,j), & ! [IN]
1145  dens(:,i,j), pott(:,i,j), pres(:,i,j), & ! [OUT]
1146  pott_sfc(i,j) ) ! [OUT]
1147  enddo
1148  enddo
1149 
1150  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
Here is the call graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_atmos_rev_1d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_rev_1d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), intent(in)  temp,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc,
real(rp), dimension (ka), intent(in)  dz,
real(rp), dimension(ka), intent(inout)  dens,
real(rp), dimension(ka), intent(out)  pott,
real(rp), dimension(ka), intent(out)  pres 
)

Build up density from upermost atmosphere (1D)

Parameters
[in]temptemperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in,out]densdensity [kg/m3]
[out]pottpotential temperature [K]
[out]prespressure [Pa]

Definition at line 1259 of file scale_atmos_hydrostatic.F90.

References scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_prc::prc_abort().

1259  use scale_prc, only: &
1260  prc_abort
1261  use scale_atmos_hydrometeor, only: &
1262  cv_vapor, &
1263  cv_water, &
1264  cp_vapor, &
1265  cp_water
1266  implicit none
1267  integer, intent(in) :: ka, ks, ke
1268 
1269  real(RP), intent(in) :: temp(ka)
1270  real(RP), intent(in) :: qv (ka)
1271  real(RP), intent(in) :: qc (ka)
1272  real(RP), intent(in) :: dz (ka)
1273 
1274  real(RP), intent(inout) :: dens(ka)
1275 
1276  real(RP), intent(out) :: pott(ka)
1277  real(RP), intent(out) :: pres(ka)
1278 
1279  real(RP) :: rtot (ka)
1280  real(RP) :: cvtot (ka)
1281  real(RP) :: cptot (ka)
1282 
1283  real(RP) :: rovcp
1284  real(RP) :: dens_s, dhyd, dgrd
1285  integer :: ite
1286  logical :: converged
1287 
1288  integer :: k
1289  !---------------------------------------------------------------------------
1290 
1291  do k = ks, ke
1292  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
1293  + rvap * qv(k)
1294  cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
1295  + cv_vapor * qv(k) &
1296  + cv_water * qc(k)
1297  cptot(k) = cpdry * ( 1.0_rp - qv(k) - qc(k) ) &
1298  + cp_vapor * qv(k) &
1299  + cp_water * qc(k)
1300  enddo
1301 
1302  pres(ke) = dens(ke) * rtot(ke) * temp(ke)
1303 
1304  do k = ke-1, ke, -1
1305 
1306  dens_s = 0.0_rp
1307  dens(k) = dens(k+1) ! first guess
1308 
1309  converged = .false.
1310  do ite = 1, itelim
1311  if ( abs(dens(k)-dens_s) <= criteria ) then
1312  converged = .true.
1313  exit
1314  endif
1315 
1316  dens_s = dens(k)
1317 
1318  dhyd = - ( pres(k+1) - dens_s * rtot(k) * temp(k) ) / dz(k) & ! dp/dz
1319  - grav * 0.5_rp * ( dens(k+1) + dens_s ) ! rho*g
1320 
1321  dgrd = + rtot(k) * temp(k) / dz(k) &
1322  - 0.5_rp * grav
1323 
1324  dens(k) = dens_s - dhyd/dgrd
1325 
1326  if ( dens(k)*0.0_rp /= 0.0_rp ) exit
1327  enddo
1328 
1329  if ( .NOT. converged ) then
1330  log_error("ATMOS_HYDROSTATIC_buildrho_bytemp_atmos_rev_1D",*) 'iteration not converged!', &
1331  k,dens(k),ite,dens_s,dhyd,dgrd
1332  call prc_abort
1333  endif
1334 
1335  pres(k) = dens(k) * rtot(k) * temp(k)
1336 
1337  enddo
1338 
1339  do k = ks, ke
1340  rovcp = rtot(k) / cptot(k)
1341  pott(k) = temp(k) * ( p00 / pres(k) )**rovcp
1342  enddo
1343 
1344  return
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
module PROCESS
Definition: scale_prc.F90:11
integer, public ks
start point of inner domain: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public ka
of whole cells: z, local, with HALO
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
Here is the call graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_atmos_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_3d ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(in)  temp,
real(rp), dimension (ka,ia,ja), intent(in)  qv,
real(rp), dimension (ka,ia,ja), intent(in)  qc,
real(rp), dimension (ka,ia,ja), intent(in)  dz,
real(rp), dimension(ka,ia,ja), intent(inout)  dens,
real(rp), dimension(ka,ia,ja), intent(out)  pott,
real(rp), dimension(ka,ia,ja), intent(out)  pres 
)

Build up density from lowermost atmosphere (3D)

Parameters
[in]temptemperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in,out]densdensity [kg/m3]
[out]pottpotential temperature [K]
[out]prespressure [Pa]

Definition at line 1355 of file scale_atmos_hydrostatic.F90.

1355  implicit none
1356  integer, intent(in) :: ka, ks, ke
1357  integer, intent(in) :: ia, is, ie
1358  integer, intent(in) :: ja, js, je
1359 
1360  real(RP), intent(in) :: temp(ka,ia,ja)
1361  real(RP), intent(in) :: qv (ka,ia,ja)
1362  real(RP), intent(in) :: qc (ka,ia,ja)
1363  real(RP), intent(in) :: dz (ka,ia,ja)
1364 
1365  real(RP), intent(inout) :: dens(ka,ia,ja)
1366  real(RP), intent(out) :: pott(ka,ia,ja)
1367  real(RP), intent(out) :: pres(ka,ia,ja)
1368 
1369  integer :: i, j
1370  !---------------------------------------------------------------------------
1371 
1372  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1373  do j = js, je
1374  do i = is, ie
1375  call atmos_hydrostatic_buildrho_bytemp_atmos_1d( ka, ks, ke, &
1376  temp(:,i,j), qv(:,i,j), qc(:,i,j), & ! [IN]
1377  dz(:,i,j), & ! [IN]
1378  dens(:,i,j), & ! [INOUT]
1379  pott(:,i,j), pres(:,i,j) ) ! [OUT]
1380  enddo
1381  enddo
1382 
1383  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO