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 (dens, temp, pres, pott, qv, qc, temp_sfc, pres_sfc, pott_sfc, qv_sfc, qc_sfc)
 Build up density from surface (1D) More...
 
subroutine atmos_hydrostatic_buildrho_3d (dens, temp, pres, pott, qv, qc, temp_sfc, pres_sfc, pott_sfc, qv_sfc, qc_sfc)
 Build up density from surface (3D) More...
 
subroutine atmos_hydrostatic_buildrho_real_3d (dens, temp, pres, pott, qv, qc)
 Build up density from surface (3D), not to reverse from TOA. More...
 
subroutine, public atmos_hydrostatic_buildrho_atmos_0d (dens_L2, temp_L2, pres_L2, pott_L2, qv_L2, qc_L2, dens_L1, pott_L1, qv_L1, qc_L1, dz, k)
 Build up density (0D) More...
 
subroutine atmos_hydrostatic_buildrho_atmos_1d (dens, temp, pres, pott, qv, qc)
 Build up density from lowermost atmosphere (1D) More...
 
subroutine atmos_hydrostatic_buildrho_atmos_3d (dens, temp, pres, pott, qv, qc, dz, kref_in)
 Build up density from lowermost atmosphere (3D) More...
 
subroutine, public atmos_hydrostatic_buildrho_atmos_rev_2d (dens_L1, temp_L1, pres_L1, pott_L1, qv_L1, qc_L1, dens_L2, pott_L2, qv_L2, qc_L2, dz, k)
 Build up density (2D) More...
 
subroutine, public atmos_hydrostatic_buildrho_atmos_rev_3d (dens, temp, pres, pott, qv, qc, dz, kref_in)
 Build up density from lowermost atmosphere (3D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_1d (dens, pott, pres, temp, qv, qc, pott_sfc, pres_sfc, temp_sfc, qv_sfc, qc_sfc)
 Build up density from surface (1D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_3d (dens, pott, pres, temp, qv, qc, pott_sfc, pres_sfc, temp_sfc, qv_sfc, qc_sfc)
 Build up density from surface (3D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_1d (dens, pott, pres, temp, qv, qc)
 Build up density from lowermost atmosphere (1D) More...
 
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_3d (dens, pott, pres, temp, qv, qc)
 Build up density from lowermost atmosphere (3D) More...
 
subroutine atmos_hydrostatic_barometric_law_mslp_0d (mslp, pres, temp, dz)
 Calculate mean sea-level pressure from barometric law (0D) 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

History Output
No history output

Function/Subroutine Documentation

◆ atmos_hydrostatic_setup()

subroutine, public scale_atmos_hydrostatic::atmos_hydrostatic_setup ( )

Setup.

Definition at line 122 of file scale_atmos_sub_hydrostatic.F90.

References scale_const::const_eps, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, and scale_process::prc_mpistop().

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

122  use scale_process, only: &
124  use scale_const, only: &
125  const_eps
126  implicit none
127 
128  namelist / param_atmos_hydrostatic / &
129  hydrostatic_uselapserate, &
130  hydrostatic_buildrho_real_kref
131 
132  integer :: ierr
133  !---------------------------------------------------------------------------
134 
135  if( io_l ) write(io_fid_log,*)
136  if( io_l ) write(io_fid_log,*) '++++++ Module[HYDROSTATIC] / Categ[ATMOS SHARE] / Origin[SCALElib]'
137 
138  !--- read namelist
139  rewind(io_fid_conf)
140  read(io_fid_conf,nml=param_atmos_hydrostatic,iostat=ierr)
141  if( ierr < 0 ) then !--- missing
142  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
143  elseif( ierr > 0 ) then !--- fatal error
144  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_HYDROSTATIC. Check!'
145  call prc_mpistop
146  endif
147  if( io_lnml ) write(io_fid_log,nml=param_atmos_hydrostatic)
148 
149  criteria = sqrt( const_eps )
150 
151  if( io_l ) write(io_fid_log,*)
152  if( io_l ) write(io_fid_log,*) '*** use lapse rate for estimation of surface temperature? : ', hydrostatic_uselapserate
153  if( io_l ) write(io_fid_log,*) '*** buildrho conversion criteria : ', criteria
154 
155  if ( thermodyn_type == 'EXACT' ) then
156  cv_qv = cvvap
157  cv_qc = cl
158  elseif( thermodyn_type == 'SIMPLE' ) then
159  cv_qv = cvdry
160  cv_qc = cvdry
161  endif
162 
163  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 ( real(rp), dimension(ka), intent(out)  dens,
real(rp), dimension(ka), intent(out)  temp,
real(rp), dimension(ka), intent(out)  pres,
real(rp), dimension(ka), intent(in)  pott,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc,
real(rp), intent(out)  temp_sfc,
real(rp), intent(in)  pres_sfc,
real(rp), intent(in)  pott_sfc,
real(rp), intent(in)  qv_sfc,
real(rp), intent(in)  qc_sfc 
)

Build up density from surface (1D)

Parameters
[out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[out]temp_sfcsurface temperature [K]
[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]

Definition at line 180 of file scale_atmos_sub_hydrostatic.F90.

References atmos_hydrostatic_buildrho_atmos_1d(), scale_grid::grid_cz, scale_grid_index::ks, and scale_process::prc_mpistop().

180  use scale_process, only: &
182  implicit none
183 
184  real(RP), intent(out) :: dens(ka)
185  real(RP), intent(out) :: temp(ka)
186  real(RP), intent(out) :: pres(ka)
187  real(RP), intent(in) :: pott(ka)
188  real(RP), intent(in) :: qv (ka)
189  real(RP), intent(in) :: qc (ka)
190 
191  real(RP), intent(out) :: temp_sfc
192  real(RP), intent(in) :: pres_sfc
193  real(RP), intent(in) :: pott_sfc
194  real(RP), intent(in) :: qv_sfc
195  real(RP), intent(in) :: qc_sfc
196 
197  real(RP) :: dens_sfc
198 
199  real(RP) :: rtot_sfc
200  real(RP) :: cvtot_sfc
201  real(RP) :: cpovcv_sfc
202  real(RP) :: rtot
203  real(RP) :: cvtot
204  real(RP) :: cpovcv
205 
206  real(RP) :: cvovcp_sfc, cpovr, cvovcp, rovcv
207  real(RP) :: dens_s, dhyd, dgrd
208  integer :: ite
209  logical :: converged
210  !---------------------------------------------------------------------------
211 
212  !--- from surface to lowermost atmosphere
213 
214  rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
215  + rvap * qv_sfc
216  cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
217  + cv_qv * qv_sfc &
218  + cv_qc * qc_sfc
219  cpovcv_sfc = ( cvtot_sfc + rtot_sfc ) / cvtot_sfc
220 
221  rtot = rdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
222  + rvap * qv(ks)
223  cvtot = cvdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
224  + cv_qv * qv(ks) &
225  + cv_qc * qc(ks)
226  cpovcv = ( cvtot + rtot ) / cvtot
227 
228  ! density at surface
229  cvovcp_sfc = 1.0_rp / cpovcv_sfc
230  dens_sfc = p00 / rtot_sfc / pott_sfc * ( pres_sfc/p00 )**cvovcp_sfc
231  temp_sfc = pres_sfc / ( dens_sfc * rtot_sfc )
232 
233  ! make density at lowermost cell center
234  if ( hydrostatic_uselapserate ) then
235 
236  cpovr = ( cvtot + rtot ) / rtot
237  cvovcp = 1.0_rp / cpovcv
238 
239  temp(ks) = pott_sfc - lapsdry * grid_cz(ks) ! use dry lapse rate
240  pres(ks) = p00 * ( temp(ks)/pott(ks) )**cpovr
241  dens(ks) = p00 / rtot / pott(ks) * ( pres(ks)/p00 )**cvovcp
242 
243  else ! use itelation
244 
245  rovcv = rtot / cvtot
246 
247  dens_s = 0.0_rp
248  dens(ks) = dens_sfc ! first guess
249 
250  converged = .false.
251  do ite = 1, itelim
252  if ( abs(dens(ks)-dens_s) <= criteria ) then
253  converged = .true.
254  exit
255  endif
256 
257  dens_s = dens(ks)
258 
259  dhyd = + ( p00 * ( dens_sfc * rtot_sfc * pott_sfc / p00 )**cpovcv_sfc &
260  - p00 * ( dens_s * rtot * pott(ks) / p00 )**cpovcv ) / grid_cz(ks) & ! dp/dz
261  - grav * 0.5_rp * ( dens_sfc + dens_s ) ! rho*g
262 
263  dgrd = - p00 * ( rtot * pott(ks) / p00 )**cpovcv / grid_cz(ks) &
264  * cpovcv * dens_s**rovcv &
265  - 0.5_rp * grav
266 
267  dens(ks) = dens_s - dhyd/dgrd
268 
269  if( dens(ks)*0.0_rp /= 0.0_rp) exit
270  enddo
271 
272  if ( .NOT. converged ) then
273  write(*,*) 'xxx [buildrho 1D sfc] iteration not converged!', &
274  dens(ks),ite,dens_s,dhyd,dgrd
275  call prc_mpistop
276  endif
277 
278  endif
279 
280  !--- from lowermost atmosphere to top of atmosphere
281  call atmos_hydrostatic_buildrho_atmos_1d( dens(:), & ! [INOUT]
282  temp(:), & ! [OUT]
283  pres(:), & ! [OUT]
284  pott(:), & ! [IN]
285  qv(:), & ! [IN]
286  qc(:) ) ! [IN]
287 
288  return
subroutine, public prc_mpistop
Abort MPI.
integer, public ka
of z whole cells (local, with HALO)
module PROCESS
integer, public ks
start point of inner domain: z, local
Here is the call graph for this function:

◆ atmos_hydrostatic_buildrho_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_3d ( 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(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(1,ia,ja), intent(out)  temp_sfc,
real(rp), dimension(1,ia,ja), intent(in)  pres_sfc,
real(rp), dimension(1,ia,ja), intent(in)  pott_sfc,
real(rp), dimension (1,ia,ja), intent(in)  qv_sfc,
real(rp), dimension (1,ia,ja), intent(in)  qc_sfc 
)

Build up density from surface (3D)

Parameters
[out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[out]temp_sfcsurface temperature [K]
[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]

Definition at line 305 of file scale_atmos_sub_hydrostatic.F90.

References atmos_hydrostatic_buildrho_atmos_3d(), atmos_hydrostatic_buildrho_atmos_rev_2d(), atmos_hydrostatic_buildrho_atmos_rev_3d(), scale_comm::comm_horizontal_mean(), scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_process::prc_mpistop(), scale_grid_real::real_cz, and scale_grid_real::real_fz.

305  use scale_process, only: &
307  use scale_comm, only: &
309  implicit none
310 
311  real(RP), intent(out) :: dens(ka,ia,ja)
312  real(RP), intent(out) :: temp(ka,ia,ja)
313  real(RP), intent(out) :: pres(ka,ia,ja)
314  real(RP), intent(in) :: pott(ka,ia,ja)
315  real(RP), intent(in) :: qv (ka,ia,ja)
316  real(RP), intent(in) :: qc (ka,ia,ja)
317 
318  real(RP), intent(out) :: temp_sfc(1,ia,ja)
319  real(RP), intent(in) :: pres_sfc(1,ia,ja)
320  real(RP), intent(in) :: pott_sfc(1,ia,ja)
321  real(RP), intent(in) :: qv_sfc (1,ia,ja)
322  real(RP), intent(in) :: qc_sfc (1,ia,ja)
323 
324  real(RP) :: dz(ka,ia,ja)
325 
326  real(RP) :: dens_sfc (1,ia,ja)
327  real(RP) :: pott_toa (ia,ja)
328  real(RP) :: qv_toa (ia,ja)
329  real(RP) :: qc_toa (ia,ja)
330  real(RP) :: dens_1d (ka)
331 
332  real(RP) :: rtot_sfc (ia,ja)
333  real(RP) :: cvtot_sfc (ia,ja)
334  real(RP) :: cpovcv_sfc(ia,ja)
335  real(RP) :: rtot (ia,ja)
336  real(RP) :: cvtot (ia,ja)
337  real(RP) :: cpovcv (ia,ja)
338 
339  real(RP) :: cvovcp_sfc, cpovr, cvovcp
340 
341  integer :: k, i, j
342  !---------------------------------------------------------------------------
343 
344  !--- from surface to lowermost atmosphere
345 
346  do j = jsb, jeb
347  do i = isb, ieb
348  dz(ks,i,j) = real_cz(ks,i,j) - real_fz(ks-1,i,j) ! distance from surface to cell center
349  do k = ks+1, ke
350  dz(k,i,j) = real_cz(k,i,j) - real_cz(k-1,i,j) ! distance from cell center to cell center
351  enddo
352  dz(ke+1,i,j) = real_fz(ke,i,j) - real_cz(ke,i,j) ! distance from cell center to TOA
353  enddo
354  enddo
355 
356  do j = jsb, jeb
357  do i = isb, ieb
358  pott_toa(i,j) = pott(ke,i,j)
359  qv_toa(i,j) = qv(ke,i,j)
360  qc_toa(i,j) = qc(ke,i,j)
361  enddo
362  enddo
363 
364  ! density at surface
365  do j = jsb, jeb
366  do i = isb, ieb
367  rtot_sfc(i,j) = rdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
368  + rvap * qv_sfc(1,i,j)
369  cvtot_sfc(i,j) = cvdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
370  + cv_qv * qv_sfc(1,i,j) &
371  + cv_qc * qc_sfc(1,i,j)
372  cpovcv_sfc(i,j) = ( cvtot_sfc(i,j) + rtot_sfc(i,j) ) / cvtot_sfc(i,j)
373  enddo
374  enddo
375 
376  do j = jsb, jeb
377  do i = isb, ieb
378  rtot(i,j) = rdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
379  + rvap * qv(ks,i,j)
380  cvtot(i,j) = cvdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
381  + cv_qv * qv(ks,i,j) &
382  + cv_qc * qc(ks,i,j)
383  cpovcv(i,j) = ( cvtot(i,j) + rtot(i,j) ) / cvtot(i,j)
384  enddo
385  enddo
386 
387  ! density at surface
388  do j = jsb, jeb
389  do i = isb, ieb
390  cvovcp_sfc = 1.0_rp / cpovcv_sfc(i,j)
391  dens_sfc(1,i,j) = p00 / rtot_sfc(i,j) / pott_sfc(1,i,j) * ( pres_sfc(1,i,j)/p00 )**cvovcp_sfc
392  temp_sfc(1,i,j) = pres_sfc(1,i,j) / ( dens_sfc(1,i,j) * rtot_sfc(i,j) )
393  enddo
394  enddo
395 
396  ! make density at lowermost cell center
397  if ( hydrostatic_uselapserate ) then
398 
399  do j = jsb, jeb
400  do i = isb, ieb
401  cpovr = ( cvtot(i,j) + rtot(i,j) ) / rtot(i,j)
402  cvovcp = 1.0_rp / cpovcv(i,j)
403 
404  temp(ks,i,j) = pott_sfc(1,i,j) - lapsdry * ( real_cz(ks,i,j) - real_fz(ks-1,i,j) ) ! use dry lapse rate
405  pres(ks,i,j) = p00 * ( temp(ks,i,j)/pott(ks,i,j) )**cpovr
406  dens(ks,i,j) = p00 / rtot(i,j) / pott(ks,i,j) * ( pres(ks,i,j)/p00 )**cvovcp
407  enddo
408  enddo
409 
410  else ! use itelation
411 
412  call atmos_hydrostatic_buildrho_atmos_2d( dens(ks,:,:), & ! [OUT]
413  temp(ks,:,:), & ! [OUT]->not used
414  pres(ks,:,:), & ! [OUT]->not used
415  pott(ks,:,:), & ! [IN]
416  qv(ks,:,:), & ! [IN]
417  qc(ks,:,:), & ! [IN]
418  dens_sfc(1,:,:), & ! [IN]
419  pott_sfc(1,:,:), & ! [IN]
420  qv_sfc(1,:,:), & ! [IN]
421  qc_sfc(1,:,:), & ! [IN]
422  dz(ks,:,:), & ! [IN]
423  ks ) ! [IN]
424 
425  endif
426 
427  !--- from lowermost atmosphere to top of atmosphere
428  call atmos_hydrostatic_buildrho_atmos_3d( dens(:,:,:), & ! [INOUT]
429  temp(:,:,:), & ! [OUT]
430  pres(:,:,:), & ! [OUT]
431  pott(:,:,:), & ! [IN]
432  qv(:,:,:), & ! [IN]
433  qc(:,:,:), & ! [IN]
434  dz(:,:,:) ) ! [IN]
435 
436  call atmos_hydrostatic_buildrho_atmos_2d( dens(ke+1,:,:), & ! [OUT]
437  temp(ke+1,:,:), & ! [OUT]->not used
438  pres(ke+1,:,:), & ! [OUT]->not used
439  pott_toa(:,:), & ! [IN]
440  qv_toa(:,:), & ! [IN]
441  qc_toa(:,:), & ! [IN]
442  dens(ke ,:,:), & ! [IN]
443  pott(ke ,:,:), & ! [IN]
444  qv(ke ,:,:), & ! [IN]
445  qc(ke ,:,:), & ! [IN]
446  dz(ke+1,:,:), & ! [IN]
447  ke+1 ) ! [IN]
448 
449  ! density at TOA
450  dens( 1:ks-1,:,:) = 0.d0 ! fill dummy
451  dens(ke+2:ka ,:,:) = 0.d0 ! fill dummy
452 
453  call comm_horizontal_mean( dens_1d(:), dens(:,:,:) )
454  do j = jsb, jeb
455  do i = isb, ieb
456  dens(:,i,j) = dens_1d(:)
457  enddo
458  enddo
459 
460  call atmos_hydrostatic_buildrho_atmos_rev_2d( dens(ke ,:,:), & ! [OUT]
461  temp(ke ,:,:), & ! [OUT]->not used
462  pres(ke ,:,:), & ! [OUT]->not used
463  pott(ke ,:,:), & ! [IN]
464  qv(ke ,:,:), & ! [IN]
465  qc(ke ,:,:), & ! [IN]
466  dens(ke+1,:,:), & ! [IN]
467  pott_toa(:,:), & ! [IN]
468  qv_toa(:,:), & ! [IN]
469  qc_toa(:,:), & ! [IN]
470  dz(ke+1,:,:), & ! [IN]
471  ke+1 ) ! [IN]
472 
473  !--- from top of atmosphere to lowermost atmosphere
474  call atmos_hydrostatic_buildrho_atmos_rev_3d( dens(:,:,:), & ! [INOUT]
475  temp(:,:,:), & ! [OUT]
476  pres(:,:,:), & ! [OUT]
477  pott(:,:,:), & ! [IN]
478  qv(:,:,:), & ! [IN]
479  qc(:,:,:), & ! [IN]
480  dz(:,:,:) ) ! [IN]
481 
482  return
subroutine, public prc_mpistop
Abort MPI.
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)
subroutine, public comm_horizontal_mean(varmean, var)
calculate horizontal mean (global total with communication)
Definition: scale_comm.F90:516
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
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 call graph for this function:

◆ atmos_hydrostatic_buildrho_real_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_real_3d ( 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(inout)  pres,
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 
)

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

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

Definition at line 494 of file scale_atmos_sub_hydrostatic.F90.

References atmos_hydrostatic_buildrho_atmos_3d(), atmos_hydrostatic_buildrho_atmos_rev_3d(), scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_process::prc_mpistop(), scale_grid_real::real_cz, and scale_grid_real::real_fz.

494  use scale_process, only: &
496  implicit none
497 
498  real(RP), intent(out) :: dens(ka,ia,ja)
499  real(RP), intent(out) :: temp(ka,ia,ja)
500  real(RP), intent(inout) :: pres(ka,ia,ja)
501  real(RP), intent(in) :: pott(ka,ia,ja)
502  real(RP), intent(in) :: qv (ka,ia,ja)
503  real(RP), intent(in) :: qc (ka,ia,ja)
504 
505  real(RP) :: dz(ka,ia,ja)
506 
507  real(RP) :: pott_toa(ia,ja)
508  real(RP) :: qv_toa (ia,ja)
509  real(RP) :: qc_toa (ia,ja)
510 
511  real(RP) :: rtot
512  real(RP) :: cvtot
513  real(RP) :: cvovcp
514 
515  integer :: kref
516  integer :: k, i, j
517  !---------------------------------------------------------------------------
518 
519  !--- from surface to lowermost atmosphere
520 
521  do j = jsb, jeb
522  do i = isb, ieb
523  dz(ks,i,j) = real_cz(ks,i,j) - real_fz(ks-1,i,j) ! distance from surface to cell center
524  do k = ks+1, ke
525  dz(k,i,j) = real_cz(k,i,j) - real_cz(k-1,i,j) ! distance from cell center to cell center
526  enddo
527  dz(ke+1,i,j) = real_fz(ke,i,j) - real_cz(ke,i,j) ! distance from cell center to TOA
528  enddo
529  enddo
530 
531  do j = jsb, jeb
532  do i = isb, ieb
533  pott_toa(i,j) = pott(ke,i,j)
534  qv_toa(i,j) = qv(ke,i,j)
535  qc_toa(i,j) = qc(ke,i,j)
536  enddo
537  enddo
538 
539  kref = hydrostatic_buildrho_real_kref + ks - 1
540 
541  ! calc density at reference level
542  do j = jsb, jeb
543  do i = isb, ieb
544  rtot = rdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
545  + rvap * qv(kref,i,j)
546  cvtot = cvdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
547  + cv_qv * qv(kref,i,j) &
548  + cv_qc * qc(kref,i,j)
549  cvovcp = cvtot / ( cvtot + rtot )
550  dens(kref,i,j) = p00 / ( rtot * pott(kref,i,j) ) * ( pres(kref,i,j)/p00 )**cvovcp
551  enddo
552  enddo
553 
554  !--- from lowermost atmosphere to top of atmosphere
555  call atmos_hydrostatic_buildrho_atmos_3d( dens(:,:,:), & ! [INOUT]
556  temp(:,:,:), & ! [OUT]
557  pres(:,:,:), & ! [OUT]
558  pott(:,:,:), & ! [IN]
559  qv(:,:,:), & ! [IN]
560  qc(:,:,:), & ! [IN]
561  dz(:,:,:), & ! [IN]
562  kref ) ! [IN]
563  call atmos_hydrostatic_buildrho_atmos_rev_3d( dens(:,:,:), & ! [INOUT]
564  temp(:,:,:), & ! [OUT]
565  pres(:,:,:), & ! [OUT]
566  pott(:,:,:), & ! [IN]
567  qv(:,:,:), & ! [IN]
568  qc(:,:,:), & ! [IN]
569  dz(:,:,:), & ! [IN]
570  kref ) ! [IN]
571 
572  call atmos_hydrostatic_buildrho_atmos_2d( dens(ke+1,:,:), & ! [OUT]
573  temp(ke+1,:,:), & ! [OUT]->not used
574  pres(ke+1,:,:), & ! [OUT]->not used
575  pott_toa(:,:), & ! [IN]
576  qv_toa(:,:), & ! [IN]
577  qc_toa(:,:), & ! [IN]
578  dens(ke ,:,:), & ! [IN]
579  pott(ke ,:,:), & ! [IN]
580  qv(ke ,:,:), & ! [IN]
581  qc(ke ,:,:), & ! [IN]
582  dz(ke+1,:,:), & ! [IN]
583  ke+1 ) ! [IN]
584 
585  ! density at TOA
586  dens( 1:ks-1,:,:) = 0.d0 ! fill dummy
587  dens(ke+2:ka ,:,:) = 0.d0 ! fill dummy
588 
589  return
subroutine, public prc_mpistop
Abort MPI.
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)
module PROCESS
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 call graph for this function:

◆ atmos_hydrostatic_buildrho_atmos_0d()

subroutine, public scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_0d ( real(rp), intent(out)  dens_L2,
real(rp), intent(out)  temp_L2,
real(rp), intent(out)  pres_L2,
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 
)

Build up density (0D)

Parameters
[out]dens_l2density at layer 2 [kg/m3]
[out]temp_l2temperature at layer 2 [K]
[out]pres_l2pressure at layer 2 [Pa]
[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

Definition at line 607 of file scale_atmos_sub_hydrostatic.F90.

References scale_process::prc_mpistop().

Referenced by scale_atmos_refstate::atmos_refstate_calc3d().

607  use scale_process, only: &
609  implicit none
610 
611  real(RP), intent(out) :: dens_l2
612  real(RP), intent(out) :: temp_l2
613  real(RP), intent(out) :: pres_l2
614  real(RP), intent(in) :: pott_l2
615  real(RP), intent(in) :: qv_l2
616  real(RP), intent(in) :: qc_l2
617  real(RP), intent(in) :: dens_l1
618  real(RP), intent(in) :: pott_l1
619  real(RP), intent(in) :: qv_l1
620  real(RP), intent(in) :: qc_l1
621  real(RP), intent(in) :: dz
622  integer, intent(in) :: k
623 
624  real(RP) :: rtot_l1 , rtot_l2
625  real(RP) :: cvtot_l1 , cvtot_l2
626  real(RP) :: cpovcv_l1, cpovcv_l2
627 
628  real(RP) :: rovcv
629  real(RP) :: dens_s, dhyd, dgrd
630  integer :: ite
631  logical :: converged
632  !---------------------------------------------------------------------------
633 
634  rtot_l1 = rdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
635  + rvap * qv_l1
636  cvtot_l1 = cvdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
637  + cv_qv * qv_l1 &
638  + cv_qc * qc_l1
639  cpovcv_l1 = ( cvtot_l1 + rtot_l1 ) / cvtot_l1
640 
641  rtot_l2 = rdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
642  + rvap * qv_l2
643  cvtot_l2 = cvdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
644  + cv_qv * qv_l2 &
645  + cv_qc * qc_l2
646  cpovcv_l2 = ( cvtot_l2 + rtot_l2 ) / cvtot_l2
647 
648  rovcv = rtot_l2 / cvtot_l2
649 
650  dens_s = 0.0_rp
651  dens_l2 = dens_l1 ! first guess
652 
653  converged = .false.
654  do ite = 1, itelim
655  if ( abs(dens_l2-dens_s) <= criteria ) then
656  converged = .true.
657  exit
658  endif
659 
660  dens_s = dens_l2
661 
662  dhyd = + ( p00 * ( dens_l1 * rtot_l1 * pott_l1 / p00 )**cpovcv_l1 &
663  - p00 * ( dens_s * rtot_l2 * pott_l2 / p00 )**cpovcv_l2 ) / dz & ! dp/dz
664  - grav * 0.5_rp * ( dens_l1 + dens_s ) ! rho*g
665 
666  dgrd = - p00 * ( rtot_l2 * pott_l2 / p00 )**cpovcv_l2 / dz &
667  * cpovcv_l2 * dens_s**rovcv &
668  - 0.5_rp * grav
669 
670  dens_l2 = dens_s - dhyd/dgrd
671 
672  if( dens_l2*0.0_rp /= 0.0_rp) exit
673  enddo
674 
675  if ( .NOT. converged ) then
676  write(*,*) 'xxx [buildrho 0D atmos] iteration not converged!', &
677  k,dens_l2,ite,dens_s,dhyd,dgrd
678  call prc_mpistop
679  endif
680 
681  pres_l2 = p00 * ( dens_l2 * rtot_l2 * pott_l2 / p00 )**cpovcv_l2
682  temp_l2 = pres_l2 / ( dens_l2 * rtot_l2 )
683 
684  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
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 ( real(rp), dimension(ka), intent(inout)  dens,
real(rp), dimension(ka), intent(out)  temp,
real(rp), dimension(ka), intent(out)  pres,
real(rp), dimension(ka), intent(in)  pott,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc 
)

Build up density from lowermost atmosphere (1D)

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

Definition at line 696 of file scale_atmos_sub_hydrostatic.F90.

References scale_grid::grid_fdz, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, and scale_process::prc_mpistop().

Referenced by atmos_hydrostatic_buildrho_1d().

696  use scale_process, only: &
698  implicit none
699 
700  real(RP), intent(inout) :: dens(ka)
701  real(RP), intent(out) :: temp(ka)
702  real(RP), intent(out) :: pres(ka)
703  real(RP), intent(in) :: pott(ka)
704  real(RP), intent(in) :: qv (ka)
705  real(RP), intent(in) :: qc (ka)
706 
707  real(RP) :: rtot (ka)
708  real(RP) :: cvtot (ka)
709  real(RP) :: cpovcv(ka)
710 
711  real(RP) :: rovcv
712  real(RP) :: dens_s, dhyd, dgrd
713  integer :: ite
714  logical :: converged
715 
716  integer :: k
717  !---------------------------------------------------------------------------
718 
719  do k = ks, ke
720  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
721  + rvap * qv(k)
722  cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
723  + cv_qv * qv(k) &
724  + cv_qc * qc(k)
725  cpovcv(k) = ( cvtot(k) + rtot(k) ) / cvtot(k)
726  enddo
727 
728  do k = ks+1, ke
729  rovcv = rtot(k) / cvtot(k)
730 
731  dens_s = 0.0_rp
732  dens(k) = dens(k-1) ! first guess
733 
734  converged = .false.
735  do ite = 1, itelim
736  if ( abs(dens(k)-dens_s) <= criteria ) then
737  converged = .true.
738  exit
739  endif
740 
741  dens_s = dens(k)
742 
743  dhyd = + ( p00 * ( dens(k-1) * rtot(k-1) * pott(k-1) / p00 )**cpovcv(k-1) &
744  - p00 * ( dens_s * rtot(k ) * pott(k ) / p00 )**cpovcv(k ) ) / grid_fdz(k-1) & ! dpdz
745  - grav * 0.5_rp * ( dens(k-1) + dens_s ) ! rho*g
746 
747  dgrd = - p00 * ( rtot(k) * pott(k) / p00 )**cpovcv(k) / grid_fdz(k-1) &
748  * cpovcv(k) * dens_s**rovcv &
749  - 0.5_rp * grav
750 
751  dens(k) = dens_s - dhyd/dgrd
752 
753  if( dens(k)*0.0_rp /= 0.0_rp) exit
754  enddo
755 
756  if ( .NOT. converged ) then
757  write(*,*) 'xxx [buildrho 1D atmos] iteration not converged!', &
758  k,dens(k),ite,dens_s,dhyd,dgrd
759  call prc_mpistop
760  endif
761  enddo
762 
763  do k = ks, ke
764  pres(k) = p00 * ( dens(k) * rtot(k) * pott(k) / p00 )**cpovcv(k)
765  temp(k) = pres(k) / ( dens(k) * rtot(k) )
766  enddo
767 
768  dens( 1:ks-1) = dens(ks)
769  dens(ke+1:ka ) = dens(ke)
770  pres( 1:ks-1) = pres(ks)
771  pres(ke+1:ka ) = pres(ke)
772  temp( 1:ks-1) = temp(ks)
773  temp(ke+1:ka ) = temp(ke)
774 
775  return
subroutine, public prc_mpistop
Abort MPI.
integer, public ke
end point of inner domain: z, local
integer, public ka
of z whole cells (local, with HALO)
module PROCESS
integer, public ks
start point of inner domain: z, local
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_atmos_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_3d ( real(rp), dimension(ka,ia,ja), intent(inout)  dens,
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)  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)  dz,
integer, intent(in), optional  kref_in 
)

Build up density from lowermost atmosphere (3D)

Parameters
[in,out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in]dzdistance between the layer (center) [m]

Definition at line 898 of file scale_atmos_sub_hydrostatic.F90.

References 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 scale_process::prc_mpistop().

Referenced by atmos_hydrostatic_buildrho_3d(), and atmos_hydrostatic_buildrho_real_3d().

898  use scale_process, only: &
900  implicit none
901 
902  real(RP), intent(inout) :: dens(ka,ia,ja)
903  real(RP), intent(out) :: temp(ka,ia,ja)
904  real(RP), intent(out) :: pres(ka,ia,ja)
905  real(RP), intent(in) :: pott(ka,ia,ja)
906  real(RP), intent(in) :: qv (ka,ia,ja)
907  real(RP), intent(in) :: qc (ka,ia,ja)
908  real(RP), intent(in) :: dz (ka,ia,ja)
909  integer, intent(in), optional :: kref_in
910 
911  real(RP) :: rtot (ka,ia,ja)
912  real(RP) :: cvtot (ka,ia,ja)
913  real(RP) :: cpovcv(ka,ia,ja)
914 
915  real(RP) :: rovcv
916  real(RP) :: dens_s, dhyd, dgrd
917  integer :: ite
918  logical :: converged
919 
920  integer :: kref
921  integer :: k, i, j
922  !---------------------------------------------------------------------------
923 
924  if ( present(kref_in) ) then
925  kref = kref_in
926  else
927  kref = ks
928  end if
929 
930  do j = jsb, jeb
931  do i = isb, ieb
932  do k = kref, ke
933  rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
934  + rvap * qv(k,i,j)
935  cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
936  + cv_qv * qv(k,i,j) &
937  + cv_qc * qc(k,i,j)
938  cpovcv(k,i,j) = ( cvtot(k,i,j) + rtot(k,i,j) ) / cvtot(k,i,j)
939  enddo
940  enddo
941  enddo
942 
943  do j = jsb, jeb
944  do i = isb, ieb
945  do k = kref+1, ke
946  rovcv = rtot(k,i,j) / cvtot(k,i,j)
947 
948  dens_s = 0.0_rp
949  dens(k,i,j) = dens(k-1,i,j) ! first guess
950 
951  converged = .false.
952  do ite = 1, itelim
953  if ( abs(dens(k,i,j)-dens_s) <= criteria ) then
954  converged = .true.
955  exit
956  endif
957 
958  dens_s = dens(k,i,j)
959 
960  dhyd = + ( p00 * ( dens(k-1,i,j) * rtot(k-1,i,j) * pott(k-1,i,j) / p00 )**cpovcv(k-1,i,j) &
961  - p00 * ( dens_s * rtot(k ,i,j) * pott(k ,i,j) / p00 )**cpovcv(k ,i,j) ) / dz(k,i,j) & ! dpdz
962  - grav * 0.5_rp * ( dens(k-1,i,j) + dens_s ) ! rho*g
963 
964  dgrd = - p00 * ( rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j) / dz(k,i,j) &
965  * cpovcv(k,i,j) * dens_s**rovcv &
966  - 0.5_rp * grav
967 
968  dens(k,i,j) = dens_s - dhyd/dgrd
969 
970  if( dens(k,i,j)*0.0_rp /= 0.0_rp) exit
971  enddo
972 
973  if ( .NOT. converged ) then
974  write(*,*) 'xxx [buildrho 3D atmos] iteration not converged!', &
975  k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
976  call prc_mpistop
977  endif
978  enddo
979  enddo
980  enddo
981 
982  do j = jsb, jeb
983  do i = isb, ieb
984  do k = kref, ke
985  pres(k,i,j) = p00 * ( dens(k,i,j) * rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j)
986  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot(k,i,j) )
987  enddo
988  enddo
989  enddo
990 
991  return
subroutine, public prc_mpistop
Abort MPI.
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)
module PROCESS
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 call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_atmos_rev_2d()

subroutine, public scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_2d ( real(rp), dimension(ia,ja), intent(out)  dens_L1,
real(rp), dimension(ia,ja), intent(out)  temp_L1,
real(rp), dimension(ia,ja), intent(out)  pres_L1,
real(rp), dimension(ia,ja), intent(in)  pott_L1,
real(rp), dimension (ia,ja), intent(in)  qv_L1,
real(rp), dimension (ia,ja), intent(in)  qc_L1,
real(rp), dimension(ia,ja), intent(in)  dens_L2,
real(rp), dimension(ia,ja), intent(in)  pott_L2,
real(rp), dimension (ia,ja), intent(in)  qv_L2,
real(rp), dimension (ia,ja), intent(in)  qc_L2,
real(rp), dimension (ia,ja), intent(in)  dz,
integer, intent(in)  k 
)

Build up density (2D)

Parameters
[out]dens_l1density at layer 1 [kg/m3]
[out]temp_l1temperature at layer 1 [K]
[out]pres_l1pressure 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]dens_l2density at layer 2 [Pa]
[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]dzdistance from layer 1 to layer 2 [m]
[in]kfor monitor

Definition at line 1009 of file scale_atmos_sub_hydrostatic.F90.

References scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, and scale_process::prc_mpistop().

Referenced by atmos_hydrostatic_buildrho_3d(), and scale_atmos_refstate::atmos_refstate_calc3d().

1009  use scale_process, only: &
1010  prc_mpistop
1011  implicit none
1012 
1013  real(RP), intent(out) :: dens_l1(ia,ja)
1014  real(RP), intent(out) :: temp_l1(ia,ja)
1015  real(RP), intent(out) :: pres_l1(ia,ja)
1016  real(RP), intent(in) :: pott_l1(ia,ja)
1017  real(RP), intent(in) :: qv_l1 (ia,ja)
1018  real(RP), intent(in) :: qc_l1 (ia,ja)
1019  real(RP), intent(in) :: dens_l2(ia,ja)
1020  real(RP), intent(in) :: pott_l2(ia,ja)
1021  real(RP), intent(in) :: qv_l2 (ia,ja)
1022  real(RP), intent(in) :: qc_l2 (ia,ja)
1023  real(RP), intent(in) :: dz (ia,ja)
1024  integer, intent(in) :: k
1025 
1026  real(RP) :: rtot_l1 (ia,ja), rtot_l2 (ia,ja)
1027  real(RP) :: cvtot_l1 (ia,ja), cvtot_l2 (ia,ja)
1028  real(RP) :: cpovcv_l1(ia,ja), cpovcv_l2(ia,ja)
1029 
1030  real(RP) :: rovcv
1031  real(RP) :: dens_s, dhyd, dgrd
1032  integer :: ite
1033  logical :: converged
1034 
1035  integer :: i, j
1036  !---------------------------------------------------------------------------
1037 
1038  do j = jsb, jeb
1039  do i = isb, ieb
1040  rtot_l1(i,j) = rdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1041  + rvap * qv_l1(i,j)
1042  cvtot_l1(i,j) = cvdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1043  + cv_qv * qv_l1(i,j) &
1044  + cv_qc * qc_l1(i,j)
1045  cpovcv_l1(i,j) = ( cvtot_l1(i,j) + rtot_l1(i,j) ) / cvtot_l1(i,j)
1046 
1047  rtot_l2(i,j) = rdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1048  + rvap * qv_l2(i,j)
1049  cvtot_l2(i,j) = cvdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1050  + cv_qv * qv_l2(i,j) &
1051  + cv_qc * qc_l2(i,j)
1052  cpovcv_l2(i,j) = ( cvtot_l2(i,j) + rtot_l2(i,j) ) / cvtot_l2(i,j)
1053  enddo
1054  enddo
1055 
1056  do j = jsb, jeb
1057  do i = isb, ieb
1058  rovcv = rtot_l1(i,j) / cvtot_l1(i,j)
1059 
1060  dens_s = 0.0_rp
1061  dens_l1(i,j) = dens_l2(i,j) ! first guess
1062 
1063  converged = .false.
1064  do ite = 1, itelim
1065  if ( abs(dens_l1(i,j)-dens_s) <= criteria ) then
1066  converged = .true.
1067  exit
1068  endif
1069 
1070  dens_s = dens_l1(i,j)
1071 
1072  dhyd = + ( p00 * ( dens_s * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) &
1073  - p00 * ( dens_l2(i,j) * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) ) / dz(i,j) & ! dp/dz
1074  - grav * 0.5_rp * ( dens_s + dens_l2(i,j) ) ! rho*g
1075 
1076  dgrd = + p00 * ( rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) / dz(i,j) &
1077  * cpovcv_l1(i,j) * dens_s**rovcv &
1078  - 0.5_rp * grav
1079 
1080  dens_l1(i,j) = dens_s - dhyd/dgrd
1081 
1082  if( dens_l1(i,j)*0.0_rp /= 0.0_rp) exit
1083  enddo
1084 
1085  if ( .NOT. converged ) then
1086  write(*,*) 'xxx [buildrho 2D rev atmos] iteration not converged!', &
1087  i,j,k,dens_l1(i,j),ite,dens_s,dhyd,dgrd
1088  call prc_mpistop
1089  endif
1090  enddo
1091  enddo
1092 
1093  do j = jsb, jeb
1094  do i = isb, ieb
1095  pres_l1(i,j) = p00 * ( dens_l1(i,j) * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j)
1096  temp_l1(i,j) = pres_l1(i,j) / ( dens_l1(i,j) * rtot_l1(i,j) )
1097  enddo
1098  enddo
1099 
1100  return
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public ia
of x whole cells (local, with HALO)
module PROCESS
integer, public isb
integer, public jsb
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_atmos_rev_3d()

subroutine, public scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_3d ( real(rp), dimension(ka,ia,ja), intent(inout)  dens,
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)  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)  dz,
integer, intent(in), optional  kref_in 
)

Build up density from lowermost atmosphere (3D)

Parameters
[in,out]densdensity [kg/m3]
[out]temptemperature [K]
[out]prespressure [Pa]
[in]pottpotential temperature [K]
[in]qvwater vapor [kg/kg]
[in]qcliquid water [kg/kg]
[in]dzdistance between the layer (center) [m]

Definition at line 1114 of file scale_atmos_sub_hydrostatic.F90.

References 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 scale_process::prc_mpistop().

Referenced by atmos_hydrostatic_buildrho_3d(), atmos_hydrostatic_buildrho_real_3d(), and scale_atmos_refstate::atmos_refstate_calc3d().

1114  use scale_process, only: &
1115  prc_mpistop
1116  implicit none
1117 
1118  real(RP), intent(inout) :: dens(ka,ia,ja)
1119  real(RP), intent(out) :: temp(ka,ia,ja)
1120  real(RP), intent(out) :: pres(ka,ia,ja)
1121  real(RP), intent(in) :: pott(ka,ia,ja)
1122  real(RP), intent(in) :: qv (ka,ia,ja)
1123  real(RP), intent(in) :: qc (ka,ia,ja)
1124  real(RP), intent(in) :: dz (ka,ia,ja)
1125  integer, intent(in), optional :: kref_in
1126 
1127  real(RP) :: rtot (ka,ia,ja)
1128  real(RP) :: cvtot (ka,ia,ja)
1129  real(RP) :: cpovcv(ka,ia,ja)
1130 
1131  real(RP) :: rovcv
1132  real(RP) :: dens_s, dhyd, dgrd
1133  integer :: ite
1134  logical :: converged
1135 
1136  integer :: kref
1137  integer :: k, i, j
1138  !---------------------------------------------------------------------------
1139 
1140  if ( present(kref_in) ) then
1141  kref = kref_in
1142  else
1143  kref = ke
1144  end if
1145 
1146  do j = jsb, jeb
1147  do i = isb, ieb
1148  do k = ks, kref
1149  rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1150  + rvap * qv(k,i,j)
1151  cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1152  + cv_qv * qv(k,i,j) &
1153  + cv_qc * qc(k,i,j)
1154  cpovcv(k,i,j) = ( cvtot(k,i,j) + rtot(k,i,j) ) / cvtot(k,i,j)
1155  enddo
1156  enddo
1157  enddo
1158 
1159  do j = jsb, jeb
1160  do i = isb, ieb
1161  do k = kref-1, ks, -1
1162  rovcv = rtot(k,i,j) / cvtot(k,i,j)
1163 
1164  dens_s = 0.0_rp
1165  dens(k,i,j) = dens(k+1,i,j) ! first guess
1166 
1167  converged = .false.
1168  do ite = 1, itelim
1169  if ( abs(dens(k,i,j)-dens_s) <= criteria ) then
1170  converged = .true.
1171  exit
1172  endif
1173 
1174  dens_s = dens(k,i,j)
1175 
1176  dhyd = + ( p00 * ( dens_s * rtot(k ,i,j) * pott(k ,i,j) / p00 )**cpovcv(k ,i,j) &
1177  - p00 * ( dens(k+1,i,j) * rtot(k+1,i,j) * pott(k+1,i,j) / p00 )**cpovcv(k+1,i,j) ) / dz(k+1,i,j) & ! dpdz
1178  - grav * 0.5_rp * ( dens_s + dens(k+1,i,j) ) ! rho*g
1179 
1180  dgrd = + p00 * ( rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j) / dz(k+1,i,j) &
1181  * cpovcv(k,i,j) * dens_s**rovcv &
1182  - 0.5_rp * grav
1183 
1184  dens(k,i,j) = dens_s - dhyd/dgrd
1185 
1186  if( dens(k,i,j)*0.0_rp /= 0.0_rp) exit
1187  enddo
1188 
1189  if ( .NOT. converged ) then
1190  write(*,*) 'xxx [buildrho 3D rev atmos] iteration not converged!', &
1191  k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1192  call prc_mpistop
1193  endif
1194  enddo
1195  enddo
1196  enddo
1197 
1198  do j = jsb, jeb
1199  do i = isb, ieb
1200  do k = ks, kref
1201  pres(k,i,j) = p00 * ( dens(k,i,j) * rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j)
1202  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot(k,i,j) )
1203  enddo
1204  enddo
1205  enddo
1206 
1207  return
subroutine, public prc_mpistop
Abort MPI.
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)
module PROCESS
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 call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_1d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_1d ( real(rp), dimension(ka), intent(out)  dens,
real(rp), dimension(ka), intent(out)  pott,
real(rp), dimension(ka), intent(out)  pres,
real(rp), dimension(ka), intent(in)  temp,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc,
real(rp), intent(out)  pott_sfc,
real(rp), intent(in)  pres_sfc,
real(rp), intent(in)  temp_sfc,
real(rp), intent(in)  qv_sfc,
real(rp), intent(in)  qc_sfc 
)

Build up density from surface (1D)

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

Definition at line 1224 of file scale_atmos_sub_hydrostatic.F90.

References atmos_hydrostatic_buildrho_bytemp_atmos_1d(), scale_grid::grid_cz, scale_grid_index::ks, and scale_process::prc_mpistop().

1224  use scale_process, only: &
1225  prc_mpistop
1226  implicit none
1227 
1228  real(RP), intent(out) :: dens(ka)
1229  real(RP), intent(out) :: pott(ka)
1230  real(RP), intent(out) :: pres(ka)
1231  real(RP), intent(in) :: temp(ka)
1232  real(RP), intent(in) :: qv (ka)
1233  real(RP), intent(in) :: qc (ka)
1234 
1235  real(RP), intent(out) :: pott_sfc
1236  real(RP), intent(in) :: pres_sfc
1237  real(RP), intent(in) :: temp_sfc
1238  real(RP), intent(in) :: qv_sfc
1239  real(RP), intent(in) :: qc_sfc
1240 
1241  real(RP) :: dens_sfc
1242 
1243  real(RP) :: rtot_sfc
1244  real(RP) :: cvtot_sfc
1245  real(RP) :: rtot
1246  real(RP) :: cvtot
1247 
1248  real(RP) :: rovcp_sfc
1249  real(RP) :: dens_s, dhyd, dgrd
1250  integer :: ite
1251  logical :: converged
1252  !---------------------------------------------------------------------------
1253 
1254  !--- from surface to lowermost atmosphere
1255 
1256  rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1257  + rvap * qv_sfc
1258  cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1259  + cv_qv * qv_sfc &
1260  + cv_qc * qc_sfc
1261 
1262  rtot = rdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
1263  + rvap * qv(ks)
1264  cvtot = cvdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
1265  + cv_qv * qv(ks) &
1266  + cv_qc * qc(ks)
1267 
1268  ! density at surface
1269  rovcp_sfc = rtot_sfc / ( cvtot_sfc + rtot_sfc )
1270  dens_sfc = pres_sfc / ( rtot_sfc * temp_sfc )
1271  pott_sfc = temp_sfc * ( p00/pres_sfc )**rovcp_sfc
1272 
1273  ! make density at lowermost cell center
1274  dens_s = 0.0_rp
1275  dens(ks) = dens_sfc ! first guess
1276 
1277  converged = .false.
1278  do ite = 1, itelim
1279  if ( abs(dens(ks)-dens_s) <= criteria ) then
1280  converged = .true.
1281  exit
1282  endif
1283 
1284  dens_s = dens(ks)
1285 
1286  dhyd = + ( dens_sfc * rtot_sfc * temp_sfc &
1287  - dens_s * rtot * temp(ks) ) / grid_cz(ks) & ! dp/dz
1288  - grav * 0.5_rp * ( dens_sfc + dens_s ) ! rho*g
1289 
1290  dgrd = - rtot * temp(ks) / grid_cz(ks) &
1291  - 0.5_rp * grav
1292 
1293  dens(ks) = dens_s - dhyd/dgrd
1294 
1295  if( dens(ks)*0.0_rp /= 0.0_rp) exit
1296  enddo
1297 
1298  if ( .NOT. converged ) then
1299  write(*,*) 'xxx [buildrho bytemp 1D sfc] iteration not converged!', &
1300  dens(ks),ite,dens_s,dhyd,dgrd
1301  call prc_mpistop
1302  endif
1303 
1304  !--- from lowermost atmosphere to top of atmosphere
1305  call atmos_hydrostatic_buildrho_bytemp_atmos_1d( dens(:), & ! [INOUT]
1306  pott(:), & ! [OUT]
1307  pres(:), & ! [OUT]
1308  temp(:), & ! [IN]
1309  qv(:), & ! [IN]
1310  qc(:) ) ! [IN]
1311 
1312  return
subroutine, public prc_mpistop
Abort MPI.
integer, public ka
of z whole cells (local, with HALO)
module PROCESS
integer, public ks
start point of inner domain: z, local
Here is the call graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_3d ( 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(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(1,ia,ja), intent(out)  pott_sfc,
real(rp), dimension(1,ia,ja), intent(in)  pres_sfc,
real(rp), dimension(1,ia,ja), intent(in)  temp_sfc,
real(rp), dimension (1,ia,ja), intent(in)  qv_sfc,
real(rp), dimension (1,ia,ja), intent(in)  qc_sfc 
)

Build up density from surface (3D)

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

Definition at line 1329 of file scale_atmos_sub_hydrostatic.F90.

References atmos_hydrostatic_buildrho_bytemp_atmos_3d(), scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ks, scale_process::prc_mpistop(), scale_grid_real::real_cz, and scale_grid_real::real_fz.

1329  use scale_process, only: &
1330  prc_mpistop
1331  implicit none
1332 
1333  real(RP), intent(out) :: dens(ka,ia,ja)
1334  real(RP), intent(out) :: pott(ka,ia,ja)
1335  real(RP), intent(out) :: pres(ka,ia,ja)
1336  real(RP), intent(in) :: temp(ka,ia,ja)
1337  real(RP), intent(in) :: qv (ka,ia,ja)
1338  real(RP), intent(in) :: qc (ka,ia,ja)
1339 
1340  real(RP), intent(out) :: pott_sfc(1,ia,ja)
1341  real(RP), intent(in) :: pres_sfc(1,ia,ja)
1342  real(RP), intent(in) :: temp_sfc(1,ia,ja)
1343  real(RP), intent(in) :: qv_sfc (1,ia,ja)
1344  real(RP), intent(in) :: qc_sfc (1,ia,ja)
1345 
1346  real(RP) :: dens_sfc (1,ia,ja)
1347 
1348  real(RP) :: rtot_sfc (ia,ja)
1349  real(RP) :: cvtot_sfc (ia,ja)
1350  real(RP) :: rtot (ia,ja)
1351  real(RP) :: cvtot (ia,ja)
1352 
1353  real(RP) :: rovcp_sfc
1354  real(RP) :: dz
1355  real(RP) :: dens_s, dhyd, dgrd
1356  integer :: ite
1357  logical :: converged
1358 
1359  integer :: i, j
1360  !---------------------------------------------------------------------------
1361 
1362  !--- from surface to lowermost atmosphere
1363 
1364  do j = jsb, jeb
1365  do i = isb, ieb
1366  rtot_sfc(i,j) = rdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1367  + rvap * qv_sfc(1,i,j)
1368  cvtot_sfc(i,j) = cvdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1369  + cv_qv * qv_sfc(1,i,j) &
1370  + cv_qc * qc_sfc(1,i,j)
1371  enddo
1372  enddo
1373 
1374  do j = jsb, jeb
1375  do i = isb, ieb
1376  rtot(i,j) = rdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
1377  + rvap * qv(ks,i,j)
1378  cvtot(i,j) = cvdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
1379  + cv_qv * qv(ks,i,j) &
1380  + cv_qc * qc(ks,i,j)
1381  enddo
1382  enddo
1383 
1384  ! density at surface
1385  do j = jsb, jeb
1386  do i = isb, ieb
1387  rovcp_sfc = rtot_sfc(i,j) / ( cvtot_sfc(i,j) + rtot_sfc(i,j) )
1388  dens_sfc(1,i,j) = pres_sfc(1,i,j) / ( rtot_sfc(i,j) * temp_sfc(1,i,j) )
1389  pott_sfc(1,i,j) = temp_sfc(1,i,j) / ( p00/pres_sfc(1,i,j) )**rovcp_sfc
1390  enddo
1391  enddo
1392 
1393  ! make density at lowermost cell center
1394  do j = jsb, jeb
1395  do i = isb, ieb
1396  dz = real_cz(ks,i,j) - real_fz(ks-1,i,j)
1397 
1398  dens_s = 0.0_rp
1399  dens(ks,i,j) = dens_sfc(1,i,j) ! first guess
1400 
1401  converged = .false.
1402  do ite = 1, itelim
1403  if ( abs(dens(ks,i,j)-dens_s) <= criteria ) then
1404  converged = .true.
1405  exit
1406  endif
1407 
1408  dens_s = dens(ks,i,j)
1409 
1410  dhyd = + ( dens_sfc(1,i,j) * rtot_sfc(i,j) * temp_sfc(1,i,j) &
1411  - dens_s * rtot(i,j) * temp(ks,i,j) ) / dz & ! dp/dz
1412  - grav * 0.5_rp * ( dens_sfc(1,i,j) + dens_s ) ! rho*g
1413 
1414  dgrd = - rtot(i,j) * temp(ks,i,j) / dz &
1415  - 0.5_rp * grav
1416 
1417  dens(ks,i,j) = dens_s - dhyd/dgrd
1418 
1419  if( dens(ks,i,j)*0.0_rp /= 0.0_rp) exit
1420  enddo
1421 
1422  if ( .NOT. converged ) then
1423  write(*,*) 'xxx [buildrho bytemp 3D sfc] iteration not converged!', &
1424  i,j,dens(ks,i,j),ite,dens_s,dhyd,dgrd
1425  call prc_mpistop
1426  endif
1427  enddo
1428  enddo
1429 
1430  !--- from lowermost atmosphere to top of atmosphere
1431  call atmos_hydrostatic_buildrho_bytemp_atmos_3d( dens(:,:,:), & ! [INOUT]
1432  pott(:,:,:), & ! [OUT]
1433  pres(:,:,:), & ! [OUT]
1434  temp(:,:,:), & ! [IN]
1435  qv(:,:,:), & ! [IN]
1436  qc(:,:,:) ) ! [IN]
1437 
1438  return
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
integer, public ieb
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
module PROCESS
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 call graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_atmos_1d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_1d ( real(rp), dimension(ka), intent(inout)  dens,
real(rp), dimension(ka), intent(out)  pott,
real(rp), dimension(ka), intent(out)  pres,
real(rp), dimension(ka), intent(in)  temp,
real(rp), dimension (ka), intent(in)  qv,
real(rp), dimension (ka), intent(in)  qc 
)

Build up density from lowermost atmosphere (1D)

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

Definition at line 1450 of file scale_atmos_sub_hydrostatic.F90.

References scale_grid::grid_fdz, scale_grid_index::ke, scale_grid_index::ks, and scale_process::prc_mpistop().

Referenced by atmos_hydrostatic_buildrho_bytemp_1d().

1450  use scale_process, only: &
1451  prc_mpistop
1452  implicit none
1453 
1454  real(RP), intent(inout) :: dens(ka)
1455  real(RP), intent(out) :: pott(ka)
1456  real(RP), intent(out) :: pres(ka)
1457  real(RP), intent(in) :: temp(ka)
1458  real(RP), intent(in) :: qv (ka)
1459  real(RP), intent(in) :: qc (ka)
1460 
1461  real(RP) :: rtot (ka)
1462  real(RP) :: cvtot (ka)
1463 
1464  real(RP) :: rovcp
1465  real(RP) :: dens_s, dhyd, dgrd
1466  integer :: ite
1467  logical :: converged
1468 
1469  integer :: k
1470  !---------------------------------------------------------------------------
1471 
1472  do k = ks, ke
1473  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
1474  + rvap * qv(k)
1475  cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
1476  + cv_qv * qv(k) &
1477  + cv_qc * qc(k)
1478  enddo
1479 
1480  do k = ks+1, ke
1481 
1482  dens_s = 0.0_rp
1483  dens(k) = dens(k-1) ! first guess
1484 
1485  converged = .false.
1486  do ite = 1, itelim
1487  if ( abs(dens(k)-dens_s) <= criteria ) then
1488  converged = .true.
1489  exit
1490  endif
1491 
1492  dens_s = dens(k)
1493 
1494  dhyd = + ( dens(k-1) * rtot(k-1) * temp(k-1) &
1495  - dens_s * rtot(k ) * temp(k ) ) / grid_fdz(k-1) & ! dp/dz
1496  - grav * 0.5_rp * ( dens(k-1) + dens_s ) ! rho*g
1497 
1498  dgrd = - rtot(k) * temp(k) / grid_fdz(k-1) &
1499  - 0.5_rp * grav
1500 
1501  dens(k) = dens_s - dhyd/dgrd
1502 
1503  if( dens(k)*0.0_rp /= 0.0_rp) exit
1504  enddo
1505 
1506  if ( .NOT. converged ) then
1507  write(*,*) 'xxx [buildrho bytemp 1D atmos] iteration not converged!', &
1508  k,dens(k),ite,dens_s,dhyd,dgrd
1509  call prc_mpistop
1510  endif
1511  enddo
1512 
1513  do k = ks, ke
1514  rovcp = rtot(k) / ( cvtot(k) + rtot(k) )
1515  pres(k) = dens(k) * rtot(k) * temp(k)
1516  pott(k) = temp(k) * ( p00 / pres(k) )**rovcp
1517  enddo
1518 
1519  return
subroutine, public prc_mpistop
Abort MPI.
integer, public ke
end point of inner domain: z, local
integer, public ka
of z whole cells (local, with HALO)
module PROCESS
integer, public ks
start point of inner domain: z, local
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_buildrho_bytemp_atmos_3d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_3d ( 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,
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 
)

Build up density from lowermost atmosphere (3D)

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

Definition at line 1531 of file scale_atmos_sub_hydrostatic.F90.

References scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, scale_process::prc_mpistop(), and scale_grid_real::real_cz.

Referenced by atmos_hydrostatic_buildrho_bytemp_3d().

1531  use scale_process, only: &
1532  prc_mpistop
1533  implicit none
1534 
1535  real(RP), intent(inout) :: dens(ka,ia,ja)
1536  real(RP), intent(out) :: pott(ka,ia,ja)
1537  real(RP), intent(out) :: pres(ka,ia,ja)
1538  real(RP), intent(in) :: temp(ka,ia,ja)
1539  real(RP), intent(in) :: qv (ka,ia,ja)
1540  real(RP), intent(in) :: qc (ka,ia,ja)
1541 
1542  real(RP) :: rtot (ka,ia,ja)
1543  real(RP) :: cvtot (ka,ia,ja)
1544 
1545  real(RP) :: rovcp
1546  real(RP) :: dz
1547  real(RP) :: dens_s, dhyd, dgrd
1548  integer :: ite
1549  logical :: converged
1550 
1551  integer :: k, i, j
1552  !---------------------------------------------------------------------------
1553 
1554  do j = jsb, jeb
1555  do i = isb, ieb
1556  do k = ks, ke
1557  rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1558  + rvap * qv(k,i,j)
1559  cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1560  + cv_qv * qv(k,i,j) &
1561  + cv_qc * qc(k,i,j)
1562  enddo
1563  enddo
1564  enddo
1565 
1566  do j = jsb, jeb
1567  do i = isb, ieb
1568  do k = ks+1, ke
1569  dz = real_cz(k,i,j) - real_cz(k-1,i,j)
1570 
1571  dens_s = 0.0_rp
1572  dens(k,i,j) = dens(k-1,i,j) ! first guess
1573 
1574  converged = .false.
1575  do ite = 1, itelim
1576  if ( abs(dens(k,i,j)-dens_s) <= criteria ) then
1577  converged = .true.
1578  exit
1579  endif
1580 
1581  dens_s = dens(k,i,j)
1582 
1583  dhyd = + ( dens(k-1,i,j) * rtot(k-1,i,j) * temp(k-1,i,j) &
1584  - dens_s * rtot(k ,i,j) * temp(k ,i,j) ) / dz & ! dpdz
1585  - grav * 0.5_rp * ( dens(k-1,i,j) + dens_s ) ! rho*g
1586 
1587  dgrd = - rtot(k,i,j) * temp(k,i,j) / dz &
1588  - 0.5_rp * grav
1589 
1590  dens(k,i,j) = dens_s - dhyd/dgrd
1591 
1592  if( dens(k,i,j)*0.0_rp /= 0.0_rp) exit
1593  enddo
1594 
1595  if ( .NOT. converged ) then
1596  write(*,*) 'xxx [buildrho bytemp 3D atmos] iteration not converged!', &
1597  k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1598  call prc_mpistop
1599  endif
1600  enddo
1601  enddo
1602  enddo
1603 
1604  do j = jsb, jeb
1605  do i = isb, ieb
1606  do k = ks, ke
1607  rovcp = rtot(k,i,j) / ( cvtot(k,i,j) + rtot(k,i,j) )
1608  pres(k,i,j) = dens(k,i,j) * rtot(k,i,j) * temp(k,i,j)
1609  pott(k,i,j) = temp(k,i,j) * ( p00 / pres(k,i,j) )**rovcp
1610  enddo
1611  enddo
1612  enddo
1613 
1614  return
subroutine, public prc_mpistop
Abort MPI.
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)
module PROCESS
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 call graph for this function:
Here is the caller graph for this function:

◆ atmos_hydrostatic_barometric_law_mslp_0d()

subroutine scale_atmos_hydrostatic::atmos_hydrostatic_barometric_law_mslp_0d ( real(rp), intent(out)  mslp,
real(rp), intent(in)  pres,
real(rp), intent(in)  temp,
real(rp), intent(in)  dz 
)

Calculate mean sea-level pressure from barometric law (0D)

Parameters
[out]mslpmean sea-level pressure [Pa]
[in]pressurface pressure [Pa]
[in]tempsurface air temperature [K]
[in]dzsurface height from MSL [m]

Definition at line 1624 of file scale_atmos_sub_hydrostatic.F90.

References scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, and scale_grid_index::jsb.

1624  implicit none
1625 
1626  real(RP), intent(out) :: mslp
1627  real(RP), intent(in) :: pres
1628  real(RP), intent(in) :: temp
1629  real(RP), intent(in) :: dz
1630 
1631  ! work
1632  real(RP) :: tm
1633  !---------------------------------------------------------------------------
1634 
1635  tm = temp + laps * dz * 0.5_rp ! column-mean air temperature
1636 
1637  ! barometric law
1638  mslp = pres * exp( grav * dz / ( rdry * tm ) )
1639 
1640  return