SCALE-RM
scale_atmos_sub_hydrostatic.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
18  use scale_prof
20 
21  use scale_const, only: &
22  grav => const_grav, &
23  rdry => const_rdry, &
24  rvap => const_rvap, &
25  cvdry => const_cvdry, &
26  cvvap => const_cvvap, &
27  cpdry => const_cpdry, &
28  cpvap => const_cpvap, &
29  cl => const_cl, &
30  laps => const_laps, &
31  lapsdry => const_lapsdry, &
32  p00 => const_pre00, &
33  thermodyn_type => const_thermodyn_type
34  use scale_grid, only: &
35  grid_cz, &
36  grid_fdz
37  use scale_grid_real, only: &
38  real_cz, &
39  real_fz
40  !-----------------------------------------------------------------------------
41  implicit none
42  private
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public procedure
46  !
47  public :: atmos_hydrostatic_setup
48  public :: atmos_hydrostatic_buildrho
49  public :: atmos_hydrostatic_buildrho_real
50  public :: atmos_hydrostatic_buildrho_atmos
51  public :: atmos_hydrostatic_buildrho_bytemp
52  public :: atmos_hydrostatic_buildrho_bytemp_atmos
53 
57 
58  public :: atmos_hydrostatic_barometric_law_mslp
59  public :: atmos_hydrostatic_barometric_law_pres
60 
61  interface atmos_hydrostatic_buildrho
62  module procedure atmos_hydrostatic_buildrho_1d
63  module procedure atmos_hydrostatic_buildrho_3d
64  end interface atmos_hydrostatic_buildrho
65 
66  interface atmos_hydrostatic_buildrho_real
67 ! module procedure ATMOS_HYDROSTATIC_buildrho_1D ! buildrho_real_1D is completely equal to buildrho_1D
68  module procedure atmos_hydrostatic_buildrho_real_3d
69  end interface atmos_hydrostatic_buildrho_real
70 
71  interface atmos_hydrostatic_buildrho_atmos
74  end interface atmos_hydrostatic_buildrho_atmos
75 
76  interface atmos_hydrostatic_buildrho_bytemp
79  end interface atmos_hydrostatic_buildrho_bytemp
80 
81  interface atmos_hydrostatic_buildrho_bytemp_atmos
84  end interface atmos_hydrostatic_buildrho_bytemp_atmos
85 
86  interface atmos_hydrostatic_barometric_law_mslp
88  module procedure atmos_hydrostatic_barometric_law_mslp_2d
89  end interface atmos_hydrostatic_barometric_law_mslp
90 
91  interface atmos_hydrostatic_barometric_law_pres
92  module procedure atmos_hydrostatic_barometric_law_pres_0d
93  module procedure atmos_hydrostatic_barometric_law_pres_2d
94  end interface atmos_hydrostatic_barometric_law_pres
95 
96  !-----------------------------------------------------------------------------
97  !
98  !++ Public parameters & variables
99  !
100  !-----------------------------------------------------------------------------
101  !
102  !++ Private procedure
103  !
104  private :: atmos_hydrostatic_buildrho_atmos_2d
105 
106  !-----------------------------------------------------------------------------
107  !
108  !++ Private parameters & variables
109  !
110  integer, private, parameter :: itelim = 100
111  real(RP), private :: criteria
112  logical, private :: HYDROSTATIC_uselapserate = .false.
113  integer, private :: HYDROSTATIC_buildrho_real_kref = 1
114 
115  real(RP), private :: CV_qv
116  real(RP), private :: CP_qv
117  real(RP), private :: CV_qc
118  real(RP), private :: CP_qc
119 
120  !-----------------------------------------------------------------------------
121 contains
122  !-----------------------------------------------------------------------------
124  subroutine atmos_hydrostatic_setup
125  use scale_process, only: &
127  use scale_const, only: &
128  const_eps
129  use scale_tracer, only: &
130  tracer_cv, &
131  tracer_cp
132  use scale_atmos_hydrometeor, only: &
133  i_qv, &
134  i_qc
135  implicit none
136 
137  namelist / param_atmos_hydrostatic / &
138  hydrostatic_uselapserate, &
139  hydrostatic_buildrho_real_kref
140 
141  integer :: ierr
142  !---------------------------------------------------------------------------
143 
144  if( io_l ) write(io_fid_log,*)
145  if( io_l ) write(io_fid_log,*) '++++++ Module[HYDROSTATIC] / Categ[ATMOS SHARE] / Origin[SCALElib]'
146 
147  !--- read namelist
148  rewind(io_fid_conf)
149  read(io_fid_conf,nml=param_atmos_hydrostatic,iostat=ierr)
150  if( ierr < 0 ) then !--- missing
151  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
152  elseif( ierr > 0 ) then !--- fatal error
153  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_HYDROSTATIC. Check!'
154  call prc_mpistop
155  endif
156  if( io_nml ) write(io_fid_nml,nml=param_atmos_hydrostatic)
157 
158  criteria = sqrt( const_eps )
159 
160  if( io_l ) write(io_fid_log,*)
161  if( io_l ) write(io_fid_log,*) '*** Use lapse rate for estimation of surface temperature? : ', hydrostatic_uselapserate
162  if( io_l ) write(io_fid_log,*) '*** Buildrho conversion criteria : ', criteria
163 
164  if ( i_qv > 0 ) then
165  cv_qv = tracer_cv(i_qv)
166  cp_qv = tracer_cp(i_qv)
167  else
168  cv_qv = 0.0_rp
169  cp_qv = 0.0_rp
170  end if
171  if ( i_qc > 0 ) then
172  cv_qc = tracer_cp(i_qc)
173  else
174  cv_qc = 0.0_rp
175  end if
176 
177  return
178  end subroutine atmos_hydrostatic_setup
179 
180  !-----------------------------------------------------------------------------
182  subroutine atmos_hydrostatic_buildrho_1d( &
183  dens, &
184  temp, &
185  pres, &
186  pott, &
187  qv, &
188  qc, &
189  temp_sfc, &
190  pres_sfc, &
191  pott_sfc, &
192  qv_sfc, &
193  qc_sfc )
194  use scale_process, only: &
196  implicit none
197 
198  real(RP), intent(out) :: dens(KA)
199  real(RP), intent(out) :: temp(KA)
200  real(RP), intent(out) :: pres(KA)
201  real(RP), intent(in) :: pott(KA)
202  real(RP), intent(in) :: qv (KA)
203  real(RP), intent(in) :: qc (KA)
204 
205  real(RP), intent(out) :: temp_sfc
206  real(RP), intent(in) :: pres_sfc
207  real(RP), intent(in) :: pott_sfc
208  real(RP), intent(in) :: qv_sfc
209  real(RP), intent(in) :: qc_sfc
210 
211  real(RP) :: dens_sfc
212 
213  real(RP) :: Rtot_sfc
214  real(RP) :: CVtot_sfc
215  real(RP) :: CPovCV_sfc
216  real(RP) :: Rtot
217  real(RP) :: CVtot
218  real(RP) :: CPtot
219  real(RP) :: CPovCV
220 
221  real(RP) :: CVovCP_sfc, CPovR, CVovCP, RovCV
222  real(RP) :: dens_s, dhyd, dgrd
223  integer :: ite
224  logical :: converged
225  !---------------------------------------------------------------------------
226 
227  !--- from surface to lowermost atmosphere
228 
229  rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
230  + rvap * qv_sfc
231  cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
232  + cv_qv * qv_sfc &
233  + cv_qc * qc_sfc
234  cpovcv_sfc = ( cvtot_sfc + rtot_sfc ) / cvtot_sfc
235 
236  rtot = rdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
237  + rvap * qv(ks)
238  cvtot = cvdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
239  + cv_qv * qv(ks) &
240  + cv_qc * qc(ks)
241  cptot = cpdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
242  + cp_qv * qv(ks) &
243  + cv_qc * qc(ks)
244  cpovcv = cptot / cvtot
245 
246  ! density at surface
247  cvovcp_sfc = 1.0_rp / cpovcv_sfc
248  dens_sfc = p00 / rtot_sfc / pott_sfc * ( pres_sfc/p00 )**cvovcp_sfc
249  temp_sfc = pres_sfc / ( dens_sfc * rtot_sfc )
250 
251  ! make density at lowermost cell center
252  if ( hydrostatic_uselapserate ) then
253 
254  cpovr = cptot / rtot
255  cvovcp = 1.0_rp / cpovcv
256 
257  temp(ks) = pott_sfc - lapsdry * grid_cz(ks) ! use dry lapse rate
258  pres(ks) = p00 * ( temp(ks)/pott(ks) )**cpovr
259  dens(ks) = p00 / rtot / pott(ks) * ( pres(ks)/p00 )**cvovcp
260 
261  else ! use itelation
262 
263  rovcv = rtot / cvtot
264 
265  dens_s = 0.0_rp
266  dens(ks) = dens_sfc ! first guess
267 
268  converged = .false.
269  do ite = 1, itelim
270  if ( abs(dens(ks)-dens_s) <= criteria ) then
271  converged = .true.
272  exit
273  endif
274 
275  dens_s = dens(ks)
276 
277  dhyd = + ( p00 * ( dens_sfc * rtot_sfc * pott_sfc / p00 )**cpovcv_sfc &
278  - p00 * ( dens_s * rtot * pott(ks) / p00 )**cpovcv ) / grid_cz(ks) & ! dp/dz
279  - grav * 0.5_rp * ( dens_sfc + dens_s ) ! rho*g
280 
281  dgrd = - p00 * ( rtot * pott(ks) / p00 )**cpovcv / grid_cz(ks) &
282  * cpovcv * dens_s**rovcv &
283  - 0.5_rp * grav
284 
285  dens(ks) = dens_s - dhyd/dgrd
286 
287  if( dens(ks)*0.0_rp /= 0.0_rp) exit
288  enddo
289 
290  if ( .NOT. converged ) then
291  write(*,*) 'xxx [buildrho 1D sfc] iteration not converged!', &
292  dens(ks),ite,dens_s,dhyd,dgrd
293  call prc_mpistop
294  endif
295 
296  endif
297 
298  !--- from lowermost atmosphere to top of atmosphere
299  call atmos_hydrostatic_buildrho_atmos_1d( dens(:), & ! [INOUT]
300  temp(:), & ! [OUT]
301  pres(:), & ! [OUT]
302  pott(:), & ! [IN]
303  qv(:), & ! [IN]
304  qc(:) ) ! [IN]
305 
306  return
307  end subroutine atmos_hydrostatic_buildrho_1d
308 
309  !-----------------------------------------------------------------------------
311  subroutine atmos_hydrostatic_buildrho_3d( &
312  dens, &
313  temp, &
314  pres, &
315  pott, &
316  qv, &
317  qc, &
318  temp_sfc, &
319  pres_sfc, &
320  pott_sfc, &
321  qv_sfc, &
322  qc_sfc )
323  use scale_process, only: &
325  use scale_comm, only: &
327  implicit none
328 
329  real(RP), intent(out) :: dens(KA,IA,JA)
330  real(RP), intent(out) :: temp(KA,IA,JA)
331  real(RP), intent(out) :: pres(KA,IA,JA)
332  real(RP), intent(in) :: pott(KA,IA,JA)
333  real(RP), intent(in) :: qv (KA,IA,JA)
334  real(RP), intent(in) :: qc (KA,IA,JA)
335 
336  real(RP), intent(out) :: temp_sfc(1,IA,JA)
337  real(RP), intent(in) :: pres_sfc(1,IA,JA)
338  real(RP), intent(in) :: pott_sfc(1,IA,JA)
339  real(RP), intent(in) :: qv_sfc (1,IA,JA)
340  real(RP), intent(in) :: qc_sfc (1,IA,JA)
341 
342  real(RP) :: dz(KA,IA,JA)
343 
344  real(RP) :: dens_sfc (1,IA,JA)
345  real(RP) :: pott_toa (IA,JA)
346  real(RP) :: qv_toa (IA,JA)
347  real(RP) :: qc_toa (IA,JA)
348  real(RP) :: dens_1D (KA)
349 
350  real(RP) :: Rtot_sfc (IA,JA)
351  real(RP) :: CVtot_sfc (IA,JA)
352  real(RP) :: CPtot_sfc (IA,JA)
353  real(RP) :: CPovCV_sfc(IA,JA)
354  real(RP) :: Rtot (IA,JA)
355  real(RP) :: CVtot (IA,JA)
356  real(RP) :: CPtot (IA,JA)
357  real(RP) :: CPovCV (IA,JA)
358 
359  real(RP) :: CVovCP_sfc, CPovR, CVovCP
360 
361  integer :: k, i, j
362  !---------------------------------------------------------------------------
363 
364  !--- from surface to lowermost atmosphere
365 
366  do j = jsb, jeb
367  do i = isb, ieb
368  dz(ks,i,j) = real_cz(ks,i,j) - real_fz(ks-1,i,j) ! distance from surface to cell center
369  do k = ks+1, ke
370  dz(k,i,j) = real_cz(k,i,j) - real_cz(k-1,i,j) ! distance from cell center to cell center
371  enddo
372  dz(ke+1,i,j) = real_fz(ke,i,j) - real_cz(ke,i,j) ! distance from cell center to TOA
373  enddo
374  enddo
375 
376  do j = jsb, jeb
377  do i = isb, ieb
378  pott_toa(i,j) = pott(ke,i,j)
379  qv_toa(i,j) = qv(ke,i,j)
380  qc_toa(i,j) = qc(ke,i,j)
381  enddo
382  enddo
383 
384  ! density at surface
385  do j = jsb, jeb
386  do i = isb, ieb
387  rtot_sfc(i,j) = rdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
388  + rvap * qv_sfc(1,i,j)
389  cvtot_sfc(i,j) = cvdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
390  + cv_qv * qv_sfc(1,i,j) &
391  + cv_qc * qc_sfc(1,i,j)
392  cptot_sfc(i,j) = cpdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
393  + cp_qv * qv_sfc(1,i,j) &
394  + cv_qc * qc_sfc(1,i,j)
395  cpovcv_sfc(i,j) = cptot_sfc(i,j) / cvtot_sfc(i,j)
396  enddo
397  enddo
398 
399  do j = jsb, jeb
400  do i = isb, ieb
401  rtot(i,j) = rdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
402  + rvap * qv(ks,i,j)
403  cvtot(i,j) = cvdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
404  + cv_qv * qv(ks,i,j) &
405  + cv_qc * qc(ks,i,j)
406  cptot(i,j) = cpdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
407  + cp_qv * qv(ks,i,j) &
408  + cv_qc * qc(ks,i,j)
409  cpovcv(i,j) = cptot(i,j) / cvtot(i,j)
410  enddo
411  enddo
412 
413  ! density at surface
414  do j = jsb, jeb
415  do i = isb, ieb
416  cvovcp_sfc = 1.0_rp / cpovcv_sfc(i,j)
417  dens_sfc(1,i,j) = p00 / rtot_sfc(i,j) / pott_sfc(1,i,j) * ( pres_sfc(1,i,j)/p00 )**cvovcp_sfc
418  temp_sfc(1,i,j) = pres_sfc(1,i,j) / ( dens_sfc(1,i,j) * rtot_sfc(i,j) )
419  enddo
420  enddo
421 
422  ! make density at lowermost cell center
423  if ( hydrostatic_uselapserate ) then
424 
425  do j = jsb, jeb
426  do i = isb, ieb
427  cpovr = cptot(i,j) / rtot(i,j)
428  cvovcp = 1.0_rp / cpovcv(i,j)
429 
430  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
431  pres(ks,i,j) = p00 * ( temp(ks,i,j)/pott(ks,i,j) )**cpovr
432  dens(ks,i,j) = p00 / rtot(i,j) / pott(ks,i,j) * ( pres(ks,i,j)/p00 )**cvovcp
433  enddo
434  enddo
435 
436  else ! use itelation
437 
438  call atmos_hydrostatic_buildrho_atmos_2d( dens(ks,:,:), & ! [OUT]
439  temp(ks,:,:), & ! [OUT]->not used
440  pres(ks,:,:), & ! [OUT]->not used
441  pott(ks,:,:), & ! [IN]
442  qv(ks,:,:), & ! [IN]
443  qc(ks,:,:), & ! [IN]
444  dens_sfc(1,:,:), & ! [IN]
445  pott_sfc(1,:,:), & ! [IN]
446  qv_sfc(1,:,:), & ! [IN]
447  qc_sfc(1,:,:), & ! [IN]
448  dz(ks,:,:), & ! [IN]
449  ks ) ! [IN]
450 
451  endif
452 
453  !--- from lowermost atmosphere to top of atmosphere
454  call atmos_hydrostatic_buildrho_atmos_3d( dens(:,:,:), & ! [INOUT]
455  temp(:,:,:), & ! [OUT]
456  pres(:,:,:), & ! [OUT]
457  pott(:,:,:), & ! [IN]
458  qv(:,:,:), & ! [IN]
459  qc(:,:,:), & ! [IN]
460  dz(:,:,:) ) ! [IN]
461 
462  call atmos_hydrostatic_buildrho_atmos_2d( dens(ke+1,:,:), & ! [OUT]
463  temp(ke+1,:,:), & ! [OUT]->not used
464  pres(ke+1,:,:), & ! [OUT]->not used
465  pott_toa(:,:), & ! [IN]
466  qv_toa(:,:), & ! [IN]
467  qc_toa(:,:), & ! [IN]
468  dens(ke ,:,:), & ! [IN]
469  pott(ke ,:,:), & ! [IN]
470  qv(ke ,:,:), & ! [IN]
471  qc(ke ,:,:), & ! [IN]
472  dz(ke+1,:,:), & ! [IN]
473  ke+1 ) ! [IN]
474 
475  ! density at TOA
476  dens( 1:ks-1,:,:) = 0.d0 ! fill dummy
477  dens(ke+2:ka ,:,:) = 0.d0 ! fill dummy
478 
479  call comm_horizontal_mean( dens_1d(:), dens(:,:,:) )
480  do j = jsb, jeb
481  do i = isb, ieb
482  dens(:,i,j) = dens_1d(:)
483  enddo
484  enddo
485 
486  call atmos_hydrostatic_buildrho_atmos_rev_2d( dens(ke ,:,:), & ! [OUT]
487  temp(ke ,:,:), & ! [OUT]->not used
488  pres(ke ,:,:), & ! [OUT]->not used
489  pott(ke ,:,:), & ! [IN]
490  qv(ke ,:,:), & ! [IN]
491  qc(ke ,:,:), & ! [IN]
492  dens(ke+1,:,:), & ! [IN]
493  pott_toa(:,:), & ! [IN]
494  qv_toa(:,:), & ! [IN]
495  qc_toa(:,:), & ! [IN]
496  dz(ke+1,:,:), & ! [IN]
497  ke+1 ) ! [IN]
498 
499  !--- from top of atmosphere to lowermost atmosphere
500  call atmos_hydrostatic_buildrho_atmos_rev_3d( dens(:,:,:), & ! [INOUT]
501  temp(:,:,:), & ! [OUT]
502  pres(:,:,:), & ! [OUT]
503  pott(:,:,:), & ! [IN]
504  qv(:,:,:), & ! [IN]
505  qc(:,:,:), & ! [IN]
506  dz(:,:,:) ) ! [IN]
507 
508  return
509  end subroutine atmos_hydrostatic_buildrho_3d
510 
511  !-----------------------------------------------------------------------------
514  dens, &
515  temp, &
516  pres, &
517  pott, &
518  qv, &
519  qc )
520  use scale_process, only: &
522  implicit none
523 
524  real(RP), intent(out) :: dens(KA,IA,JA)
525  real(RP), intent(out) :: temp(KA,IA,JA)
526  real(RP), intent(inout) :: pres(KA,IA,JA)
527  real(RP), intent(in) :: pott(KA,IA,JA)
528  real(RP), intent(in) :: qv (KA,IA,JA)
529  real(RP), intent(in) :: qc (KA,IA,JA)
530 
531  real(RP) :: dz(KA,IA,JA)
532 
533  real(RP) :: pott_toa(IA,JA)
534  real(RP) :: qv_toa (IA,JA)
535  real(RP) :: qc_toa (IA,JA)
536 
537  real(RP) :: Rtot
538  real(RP) :: CVtot
539  real(RP) :: CPtot
540  real(RP) :: CVovCP
541 
542  integer :: kref
543  integer :: k, i, j
544  !---------------------------------------------------------------------------
545 
546  !--- from surface to lowermost atmosphere
547 
548  do j = jsb, jeb
549  do i = isb, ieb
550  dz(ks,i,j) = real_cz(ks,i,j) - real_fz(ks-1,i,j) ! distance from surface to cell center
551  do k = ks+1, ke
552  dz(k,i,j) = real_cz(k,i,j) - real_cz(k-1,i,j) ! distance from cell center to cell center
553  enddo
554  dz(ke+1,i,j) = real_fz(ke,i,j) - real_cz(ke,i,j) ! distance from cell center to TOA
555  enddo
556  enddo
557 
558  do j = jsb, jeb
559  do i = isb, ieb
560  pott_toa(i,j) = pott(ke,i,j)
561  qv_toa(i,j) = qv(ke,i,j)
562  qc_toa(i,j) = qc(ke,i,j)
563  enddo
564  enddo
565 
566  kref = hydrostatic_buildrho_real_kref + ks - 1
567 
568  ! calc density at reference level
569  do j = jsb, jeb
570  do i = isb, ieb
571  rtot = rdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
572  + rvap * qv(kref,i,j)
573  cvtot = cvdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
574  + cv_qv * qv(kref,i,j) &
575  + cv_qc * qc(kref,i,j)
576  cptot = cpdry * ( 1.0_rp - qv(kref,i,j) - qc(kref,i,j) ) &
577  + cp_qv * qv(kref,i,j) &
578  + cv_qc * qc(kref,i,j)
579  cvovcp = cvtot / cptot
580  dens(kref,i,j) = p00 / ( rtot * pott(kref,i,j) ) * ( pres(kref,i,j)/p00 )**cvovcp
581  enddo
582  enddo
583 
584  !--- from lowermost atmosphere to top of atmosphere
585  call atmos_hydrostatic_buildrho_atmos_3d( dens(:,:,:), & ! [INOUT]
586  temp(:,:,:), & ! [OUT]
587  pres(:,:,:), & ! [OUT]
588  pott(:,:,:), & ! [IN]
589  qv(:,:,:), & ! [IN]
590  qc(:,:,:), & ! [IN]
591  dz(:,:,:), & ! [IN]
592  kref ) ! [IN]
593  call atmos_hydrostatic_buildrho_atmos_rev_3d( dens(:,:,:), & ! [INOUT]
594  temp(:,:,:), & ! [OUT]
595  pres(:,:,:), & ! [OUT]
596  pott(:,:,:), & ! [IN]
597  qv(:,:,:), & ! [IN]
598  qc(:,:,:), & ! [IN]
599  dz(:,:,:), & ! [IN]
600  kref ) ! [IN]
601 
602  call atmos_hydrostatic_buildrho_atmos_2d( dens(ke+1,:,:), & ! [OUT]
603  temp(ke+1,:,:), & ! [OUT]->not used
604  pres(ke+1,:,:), & ! [OUT]->not used
605  pott_toa(:,:), & ! [IN]
606  qv_toa(:,:), & ! [IN]
607  qc_toa(:,:), & ! [IN]
608  dens(ke ,:,:), & ! [IN]
609  pott(ke ,:,:), & ! [IN]
610  qv(ke ,:,:), & ! [IN]
611  qc(ke ,:,:), & ! [IN]
612  dz(ke+1,:,:), & ! [IN]
613  ke+1 ) ! [IN]
614 
615  ! density at TOA
616  dens( 1:ks-1,:,:) = 0.d0 ! fill dummy
617  dens(ke+2:ka ,:,:) = 0.d0 ! fill dummy
618 
619  return
621 
622  !-----------------------------------------------------------------------------
625  dens_L2, &
626  temp_L2, &
627  pres_L2, &
628  pott_L2, &
629  qv_L2, &
630  qc_L2, &
631  dens_L1, &
632  pott_L1, &
633  qv_L1, &
634  qc_L1, &
635  dz, &
636  k )
637  use scale_process, only: &
639  implicit none
640 
641  real(RP), intent(out) :: dens_l2
642  real(RP), intent(out) :: temp_l2
643  real(RP), intent(out) :: pres_l2
644  real(RP), intent(in) :: pott_l2
645  real(RP), intent(in) :: qv_l2
646  real(RP), intent(in) :: qc_l2
647  real(RP), intent(in) :: dens_l1
648  real(RP), intent(in) :: pott_l1
649  real(RP), intent(in) :: qv_l1
650  real(RP), intent(in) :: qc_l1
651  real(RP), intent(in) :: dz
652  integer, intent(in) :: k
653 
654  real(RP) :: rtot_l1 , rtot_l2
655  real(RP) :: cvtot_l1 , cvtot_l2
656  real(RP) :: cptot_l1 , cptot_l2
657  real(RP) :: cpovcv_l1, cpovcv_l2
658 
659  real(RP) :: rovcv
660  real(RP) :: dens_s, dhyd, dgrd
661  integer :: ite
662  logical :: converged
663  !---------------------------------------------------------------------------
664 
665  rtot_l1 = rdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
666  + rvap * qv_l1
667  cvtot_l1 = cvdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
668  + cv_qv * qv_l1 &
669  + cv_qc * qc_l1
670  cptot_l1 = cpdry * ( 1.0_rp - qv_l1 - qc_l1 ) &
671  + cp_qv * qv_l1 &
672  + cv_qc * qc_l1
673  cpovcv_l1 = cptot_l1 / cvtot_l1
674 
675  rtot_l2 = rdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
676  + rvap * qv_l2
677  cvtot_l2 = cvdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
678  + cv_qv * qv_l2 &
679  + cv_qc * qc_l2
680  cptot_l2 = cpdry * ( 1.0_rp - qv_l2 - qc_l2 ) &
681  + cp_qv * qv_l2 &
682  + cp_qc * qc_l2
683  cpovcv_l2 = cptot_l2 / cvtot_l2
684 
685  rovcv = rtot_l2 / cvtot_l2
686 
687  dens_s = 0.0_rp
688  dens_l2 = dens_l1 ! first guess
689 
690  converged = .false.
691  do ite = 1, itelim
692  if ( abs(dens_l2-dens_s) <= criteria ) then
693  converged = .true.
694  exit
695  endif
696 
697  dens_s = dens_l2
698 
699  dhyd = + ( p00 * ( dens_l1 * rtot_l1 * pott_l1 / p00 )**cpovcv_l1 &
700  - p00 * ( dens_s * rtot_l2 * pott_l2 / p00 )**cpovcv_l2 ) / dz & ! dp/dz
701  - grav * 0.5_rp * ( dens_l1 + dens_s ) ! rho*g
702 
703  dgrd = - p00 * ( rtot_l2 * pott_l2 / p00 )**cpovcv_l2 / dz &
704  * cpovcv_l2 * dens_s**rovcv &
705  - 0.5_rp * grav
706 
707  dens_l2 = dens_s - dhyd/dgrd
708 
709  if( dens_l2*0.0_rp /= 0.0_rp) exit
710  enddo
711 
712  if ( .NOT. converged ) then
713  write(*,*) 'xxx [buildrho 0D atmos] iteration not converged!', &
714  k,dens_l2,ite,dens_s,dhyd,dgrd
715  call prc_mpistop
716  endif
717 
718  pres_l2 = p00 * ( dens_l2 * rtot_l2 * pott_l2 / p00 )**cpovcv_l2
719  temp_l2 = pres_l2 / ( dens_l2 * rtot_l2 )
720 
721  return
723 
724  !-----------------------------------------------------------------------------
727  dens, &
728  temp, &
729  pres, &
730  pott, &
731  qv, &
732  qc )
733  use scale_process, only: &
735  implicit none
736 
737  real(RP), intent(inout) :: dens(KA)
738  real(RP), intent(out) :: temp(KA)
739  real(RP), intent(out) :: pres(KA)
740  real(RP), intent(in) :: pott(KA)
741  real(RP), intent(in) :: qv (KA)
742  real(RP), intent(in) :: qc (KA)
743 
744  real(RP) :: Rtot (KA)
745  real(RP) :: CVtot (KA)
746  real(RP) :: CPtot (KA)
747  real(RP) :: CPovCV(KA)
748 
749  real(RP) :: RovCV
750  real(RP) :: dens_s, dhyd, dgrd
751  integer :: ite
752  logical :: converged
753 
754  integer :: k
755  !---------------------------------------------------------------------------
756 
757  do k = ks, ke
758  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
759  + rvap * qv(k)
760  cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
761  + cv_qv * qv(k) &
762  + cv_qc * qc(k)
763  cptot(k) = cpdry * ( 1.0_rp - qv(k) - qc(k) ) &
764  + cp_qv * qv(k) &
765  + cv_qc * qc(k)
766  cpovcv(k) = cptot(k) / cvtot(k)
767  enddo
768 
769  do k = ks+1, ke
770  rovcv = rtot(k) / cvtot(k)
771 
772  dens_s = 0.0_rp
773  dens(k) = dens(k-1) ! first guess
774 
775  converged = .false.
776  do ite = 1, itelim
777  if ( abs(dens(k)-dens_s) <= criteria ) then
778  converged = .true.
779  exit
780  endif
781 
782  dens_s = dens(k)
783 
784  dhyd = + ( p00 * ( dens(k-1) * rtot(k-1) * pott(k-1) / p00 )**cpovcv(k-1) &
785  - p00 * ( dens_s * rtot(k ) * pott(k ) / p00 )**cpovcv(k ) ) / grid_fdz(k-1) & ! dpdz
786  - grav * 0.5_rp * ( dens(k-1) + dens_s ) ! rho*g
787 
788  dgrd = - p00 * ( rtot(k) * pott(k) / p00 )**cpovcv(k) / grid_fdz(k-1) &
789  * cpovcv(k) * dens_s**rovcv &
790  - 0.5_rp * grav
791 
792  dens(k) = dens_s - dhyd/dgrd
793 
794  if( dens(k)*0.0_rp /= 0.0_rp) exit
795  enddo
796 
797  if ( .NOT. converged ) then
798  write(*,*) 'xxx [buildrho 1D atmos] iteration not converged!', &
799  k,dens(k),ite,dens_s,dhyd,dgrd
800  call prc_mpistop
801  endif
802  enddo
803 
804  do k = ks, ke
805  pres(k) = p00 * ( dens(k) * rtot(k) * pott(k) / p00 )**cpovcv(k)
806  temp(k) = pres(k) / ( dens(k) * rtot(k) )
807  enddo
808 
809  dens( 1:ks-1) = dens(ks)
810  dens(ke+1:ka ) = dens(ke)
811  pres( 1:ks-1) = pres(ks)
812  pres(ke+1:ka ) = pres(ke)
813  temp( 1:ks-1) = temp(ks)
814  temp(ke+1:ka ) = temp(ke)
815 
816  return
818 
819  !-----------------------------------------------------------------------------
821  subroutine atmos_hydrostatic_buildrho_atmos_2d( &
822  dens_L2, &
823  temp_L2, &
824  pres_L2, &
825  pott_L2, &
826  qv_L2, &
827  qc_L2, &
828  dens_L1, &
829  pott_L1, &
830  qv_L1, &
831  qc_L1, &
832  dz, &
833  k )
834  use scale_process, only: &
836  implicit none
837 
838  real(RP), intent(out) :: dens_L2(IA,JA)
839  real(RP), intent(out) :: temp_L2(IA,JA)
840  real(RP), intent(out) :: pres_L2(IA,JA)
841  real(RP), intent(in) :: pott_L2(IA,JA)
842  real(RP), intent(in) :: qv_L2 (IA,JA)
843  real(RP), intent(in) :: qc_L2 (IA,JA)
844  real(RP), intent(in) :: dens_L1(IA,JA)
845  real(RP), intent(in) :: pott_L1(IA,JA)
846  real(RP), intent(in) :: qv_L1 (IA,JA)
847  real(RP), intent(in) :: qc_L1 (IA,JA)
848  real(RP), intent(in) :: dz (IA,JA)
849  integer, intent(in) :: k
850 
851  real(RP) :: Rtot_L1 (IA,JA), Rtot_L2 (IA,JA)
852  real(RP) :: CVtot_L1 (IA,JA), CVtot_L2 (IA,JA)
853  real(RP) :: CPtot_L1 (IA,JA), CPtot_L2 (IA,JA)
854  real(RP) :: CPovCV_L1(IA,JA), CPovCV_L2(IA,JA)
855 
856  real(RP) :: RovCV
857  real(RP) :: dens_s, dhyd, dgrd
858  integer :: ite
859  logical :: converged
860 
861  integer :: i, j
862  !---------------------------------------------------------------------------
863 
864  do j = jsb, jeb
865  do i = isb, ieb
866  rtot_l1(i,j) = rdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
867  + rvap * qv_l1(i,j)
868  cvtot_l1(i,j) = cvdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
869  + cv_qv * qv_l1(i,j) &
870  + cv_qc * qc_l1(i,j)
871  cptot_l1(i,j) = cpdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
872  + cp_qv * qv_l1(i,j) &
873  + cv_qc * qc_l1(i,j)
874  cpovcv_l1(i,j) = cptot_l1(i,j) / cvtot_l1(i,j)
875 
876  rtot_l2(i,j) = rdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
877  + rvap * qv_l2(i,j)
878  cvtot_l2(i,j) = cvdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
879  + cv_qv * qv_l2(i,j) &
880  + cv_qc * qc_l2(i,j)
881  cptot_l2(i,j) = cpdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
882  + cp_qv * qv_l2(i,j) &
883  + cv_qc * qc_l2(i,j)
884  cpovcv_l2(i,j) = cptot_l2(i,j) / cvtot_l2(i,j)
885  enddo
886  enddo
887 
888  do j = jsb, jeb
889  do i = isb, ieb
890  rovcv = rtot_l2(i,j) / cvtot_l2(i,j)
891 
892  dens_s = 0.0_rp
893  dens_l2(i,j) = dens_l1(i,j) ! first guess
894 
895  converged = .false.
896  do ite = 1, itelim
897  if ( abs(dens_l2(i,j)-dens_s) <= criteria ) then
898  converged = .true.
899  exit
900  endif
901 
902  dens_s = dens_l2(i,j)
903 
904  dhyd = + ( p00 * ( dens_l1(i,j) * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) &
905  - p00 * ( dens_s * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) ) / dz(i,j) & ! dp/dz
906  - grav * 0.5_rp * ( dens_l1(i,j) + dens_s ) ! rho*g
907 
908  dgrd = - p00 * ( rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) / dz(i,j) &
909  * cpovcv_l2(i,j) * dens_s**rovcv &
910  - 0.5_rp * grav
911 
912  dens_l2(i,j) = dens_s - dhyd/dgrd
913 
914  if( dens_l2(i,j)*0.0_rp /= 0.0_rp) exit
915  enddo
916 
917  if ( .NOT. converged ) then
918  write(*,*) 'xxx [buildrho 2D atmos] iteration not converged!', &
919  i,j,k,dens_l2(i,j),ite,dens_s,dhyd,dgrd
920  call prc_mpistop
921  endif
922  enddo
923  enddo
924 
925  do j = jsb, jeb
926  do i = isb, ieb
927  pres_l2(i,j) = p00 * ( dens_l2(i,j) * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j)
928  temp_l2(i,j) = pres_l2(i,j) / ( dens_l2(i,j) * rtot_l2(i,j) )
929  enddo
930  enddo
931 
932  return
933  end subroutine atmos_hydrostatic_buildrho_atmos_2d
934 
935  !-----------------------------------------------------------------------------
938  dens, &
939  temp, &
940  pres, &
941  pott, &
942  qv, &
943  qc, &
944  dz, &
945  kref_in )
946  use scale_process, only: &
948  implicit none
949 
950  real(RP), intent(inout) :: dens(KA,IA,JA)
951  real(RP), intent(out) :: temp(KA,IA,JA)
952  real(RP), intent(out) :: pres(KA,IA,JA)
953  real(RP), intent(in) :: pott(KA,IA,JA)
954  real(RP), intent(in) :: qv (KA,IA,JA)
955  real(RP), intent(in) :: qc (KA,IA,JA)
956  real(RP), intent(in) :: dz (KA,IA,JA)
957  integer, intent(in), optional :: kref_in
958 
959  real(RP) :: Rtot (KA,IA,JA)
960  real(RP) :: CVtot (KA,IA,JA)
961  real(RP) :: CPtot (KA,IA,JA)
962  real(RP) :: CPovCV(KA,IA,JA)
963 
964  real(RP) :: RovCV
965  real(RP) :: dens_s, dhyd, dgrd
966  integer :: ite
967  logical :: converged
968 
969  integer :: kref
970  integer :: k, i, j
971  !---------------------------------------------------------------------------
972 
973  if ( present(kref_in) ) then
974  kref = kref_in
975  else
976  kref = ks
977  end if
978 
979  do j = jsb, jeb
980  do i = isb, ieb
981  do k = kref, ke
982  rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
983  + rvap * qv(k,i,j)
984  cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
985  + cv_qv * qv(k,i,j) &
986  + cv_qc * qc(k,i,j)
987  cptot(k,i,j) = cpdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
988  + cp_qv * qv(k,i,j) &
989  + cv_qc * qc(k,i,j)
990  cpovcv(k,i,j) = cptot(k,i,j) / cvtot(k,i,j)
991  enddo
992  enddo
993  enddo
994 
995  do j = jsb, jeb
996  do i = isb, ieb
997  do k = kref+1, ke
998  rovcv = rtot(k,i,j) / cvtot(k,i,j)
999 
1000  dens_s = 0.0_rp
1001  dens(k,i,j) = dens(k-1,i,j) ! first guess
1002 
1003  converged = .false.
1004  do ite = 1, itelim
1005  if ( abs(dens(k,i,j)-dens_s) <= criteria ) then
1006  converged = .true.
1007  exit
1008  endif
1009 
1010  dens_s = dens(k,i,j)
1011 
1012  dhyd = + ( p00 * ( dens(k-1,i,j) * rtot(k-1,i,j) * pott(k-1,i,j) / p00 )**cpovcv(k-1,i,j) &
1013  - p00 * ( dens_s * rtot(k ,i,j) * pott(k ,i,j) / p00 )**cpovcv(k ,i,j) ) / dz(k,i,j) & ! dpdz
1014  - grav * 0.5_rp * ( dens(k-1,i,j) + dens_s ) ! rho*g
1015 
1016  dgrd = - p00 * ( rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j) / dz(k,i,j) &
1017  * cpovcv(k,i,j) * dens_s**rovcv &
1018  - 0.5_rp * grav
1019 
1020  dens(k,i,j) = dens_s - dhyd/dgrd
1021 
1022  if( dens(k,i,j)*0.0_rp /= 0.0_rp) exit
1023  enddo
1024 
1025  if ( .NOT. converged ) then
1026  write(*,*) 'xxx [buildrho 3D atmos] iteration not converged!', &
1027  k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1028  call prc_mpistop
1029  endif
1030  enddo
1031  enddo
1032  enddo
1033 
1034  do j = jsb, jeb
1035  do i = isb, ieb
1036  do k = kref, ke
1037  pres(k,i,j) = p00 * ( dens(k,i,j) * rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j)
1038  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot(k,i,j) )
1039  enddo
1040  enddo
1041  enddo
1042 
1043  return
1045 
1046  !-----------------------------------------------------------------------------
1049  dens_L1, &
1050  temp_L1, &
1051  pres_L1, &
1052  pott_L1, &
1053  qv_L1, &
1054  qc_L1, &
1055  dens_L2, &
1056  pott_L2, &
1057  qv_L2, &
1058  qc_L2, &
1059  dz, &
1060  k )
1061  use scale_process, only: &
1062  prc_mpistop
1063  implicit none
1064 
1065  real(RP), intent(out) :: dens_l1(ia,ja)
1066  real(RP), intent(out) :: temp_l1(ia,ja)
1067  real(RP), intent(out) :: pres_l1(ia,ja)
1068  real(RP), intent(in) :: pott_l1(ia,ja)
1069  real(RP), intent(in) :: qv_l1 (ia,ja)
1070  real(RP), intent(in) :: qc_l1 (ia,ja)
1071  real(RP), intent(in) :: dens_l2(ia,ja)
1072  real(RP), intent(in) :: pott_l2(ia,ja)
1073  real(RP), intent(in) :: qv_l2 (ia,ja)
1074  real(RP), intent(in) :: qc_l2 (ia,ja)
1075  real(RP), intent(in) :: dz (ia,ja)
1076  integer, intent(in) :: k
1077 
1078  real(RP) :: rtot_l1 (ia,ja), rtot_l2 (ia,ja)
1079  real(RP) :: cvtot_l1 (ia,ja), cvtot_l2 (ia,ja)
1080  real(RP) :: cptot_l1 (ia,ja), cptot_l2 (ia,ja)
1081  real(RP) :: cpovcv_l1(ia,ja), cpovcv_l2(ia,ja)
1082 
1083  real(RP) :: rovcv
1084  real(RP) :: dens_s, dhyd, dgrd
1085  integer :: ite
1086  logical :: converged
1087 
1088  integer :: i, j
1089  !---------------------------------------------------------------------------
1090 
1091  do j = jsb, jeb
1092  do i = isb, ieb
1093  rtot_l1(i,j) = rdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1094  + rvap * qv_l1(i,j)
1095  cvtot_l1(i,j) = cvdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1096  + cv_qv * qv_l1(i,j) &
1097  + cv_qc * qc_l1(i,j)
1098  cptot_l1(i,j) = cpdry * ( 1.0_rp - qv_l1(i,j) - qc_l1(i,j) ) &
1099  + cp_qv * qv_l1(i,j) &
1100  + cv_qc * qc_l1(i,j)
1101  cpovcv_l1(i,j) = cptot_l1(i,j) / cvtot_l1(i,j)
1102 
1103  rtot_l2(i,j) = rdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1104  + rvap * qv_l2(i,j)
1105  cvtot_l2(i,j) = cvdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1106  + cv_qv * qv_l2(i,j) &
1107  + cv_qc * qc_l2(i,j)
1108  cptot_l2(i,j) = cpdry * ( 1.0_rp - qv_l2(i,j) - qc_l2(i,j) ) &
1109  + cp_qv * qv_l2(i,j) &
1110  + cv_qc * qc_l2(i,j)
1111  cpovcv_l2(i,j) = cptot_l2(i,j) / cvtot_l2(i,j)
1112  enddo
1113  enddo
1114 
1115  do j = jsb, jeb
1116  do i = isb, ieb
1117  rovcv = rtot_l1(i,j) / cvtot_l1(i,j)
1118 
1119  dens_s = 0.0_rp
1120  dens_l1(i,j) = dens_l2(i,j) ! first guess
1121 
1122  converged = .false.
1123  do ite = 1, itelim
1124  if ( abs(dens_l1(i,j)-dens_s) <= criteria ) then
1125  converged = .true.
1126  exit
1127  endif
1128 
1129  dens_s = dens_l1(i,j)
1130 
1131  dhyd = + ( p00 * ( dens_s * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) &
1132  - p00 * ( dens_l2(i,j) * rtot_l2(i,j) * pott_l2(i,j) / p00 )**cpovcv_l2(i,j) ) / dz(i,j) & ! dp/dz
1133  - grav * 0.5_rp * ( dens_s + dens_l2(i,j) ) ! rho*g
1134 
1135  dgrd = + p00 * ( rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j) / dz(i,j) &
1136  * cpovcv_l1(i,j) * dens_s**rovcv &
1137  - 0.5_rp * grav
1138 
1139  dens_l1(i,j) = dens_s - dhyd/dgrd
1140 
1141  if( dens_l1(i,j)*0.0_rp /= 0.0_rp) exit
1142  enddo
1143 
1144  if ( .NOT. converged ) then
1145  write(*,*) 'xxx [buildrho 2D rev atmos] iteration not converged!', &
1146  i,j,k,dens_l1(i,j),ite,dens_s,dhyd,dgrd
1147  call prc_mpistop
1148  endif
1149  enddo
1150  enddo
1151 
1152  do j = jsb, jeb
1153  do i = isb, ieb
1154  pres_l1(i,j) = p00 * ( dens_l1(i,j) * rtot_l1(i,j) * pott_l1(i,j) / p00 )**cpovcv_l1(i,j)
1155  temp_l1(i,j) = pres_l1(i,j) / ( dens_l1(i,j) * rtot_l1(i,j) )
1156  enddo
1157  enddo
1158 
1159  return
1161 
1162  !-----------------------------------------------------------------------------
1165  dens, &
1166  temp, &
1167  pres, &
1168  pott, &
1169  qv, &
1170  qc, &
1171  dz, &
1172  kref_in )
1173  use scale_process, only: &
1174  prc_mpistop
1175  implicit none
1176 
1177  real(RP), intent(inout) :: dens(ka,ia,ja)
1178  real(RP), intent(out) :: temp(ka,ia,ja)
1179  real(RP), intent(out) :: pres(ka,ia,ja)
1180  real(RP), intent(in) :: pott(ka,ia,ja)
1181  real(RP), intent(in) :: qv (ka,ia,ja)
1182  real(RP), intent(in) :: qc (ka,ia,ja)
1183  real(RP), intent(in) :: dz (ka,ia,ja)
1184  integer, intent(in), optional :: kref_in
1185 
1186  real(RP) :: rtot (ka,ia,ja)
1187  real(RP) :: cvtot (ka,ia,ja)
1188  real(RP) :: cptot (ka,ia,ja)
1189  real(RP) :: cpovcv(ka,ia,ja)
1190 
1191  real(RP) :: rovcv
1192  real(RP) :: dens_s, dhyd, dgrd
1193  integer :: ite
1194  logical :: converged
1195 
1196  integer :: kref
1197  integer :: k, i, j
1198  !---------------------------------------------------------------------------
1199 
1200  if ( present(kref_in) ) then
1201  kref = kref_in
1202  else
1203  kref = ke
1204  end if
1205 
1206  do j = jsb, jeb
1207  do i = isb, ieb
1208  do k = ks, kref
1209  rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1210  + rvap * qv(k,i,j)
1211  cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1212  + cv_qv * qv(k,i,j) &
1213  + cv_qc * qc(k,i,j)
1214  cptot(k,i,j) = cpdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1215  + cp_qv * qv(k,i,j) &
1216  + cv_qc * qc(k,i,j)
1217  cpovcv(k,i,j) = cptot(k,i,j) / cvtot(k,i,j)
1218  enddo
1219  enddo
1220  enddo
1221 
1222  do j = jsb, jeb
1223  do i = isb, ieb
1224  do k = kref-1, ks, -1
1225  rovcv = rtot(k,i,j) / cvtot(k,i,j)
1226 
1227  dens_s = 0.0_rp
1228  dens(k,i,j) = dens(k+1,i,j) ! first guess
1229 
1230  converged = .false.
1231  do ite = 1, itelim
1232  if ( abs(dens(k,i,j)-dens_s) <= criteria ) then
1233  converged = .true.
1234  exit
1235  endif
1236 
1237  dens_s = dens(k,i,j)
1238 
1239  dhyd = + ( p00 * ( dens_s * rtot(k ,i,j) * pott(k ,i,j) / p00 )**cpovcv(k ,i,j) &
1240  - 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
1241  - grav * 0.5_rp * ( dens_s + dens(k+1,i,j) ) ! rho*g
1242 
1243  dgrd = + p00 * ( rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j) / dz(k+1,i,j) &
1244  * cpovcv(k,i,j) * dens_s**rovcv &
1245  - 0.5_rp * grav
1246 
1247  dens(k,i,j) = dens_s - dhyd/dgrd
1248 
1249  if( dens(k,i,j)*0.0_rp /= 0.0_rp) exit
1250  enddo
1251 
1252  if ( .NOT. converged ) then
1253  write(*,*) 'xxx [buildrho 3D rev atmos] iteration not converged!', &
1254  k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1255  call prc_mpistop
1256  endif
1257  enddo
1258  enddo
1259  enddo
1260 
1261  do j = jsb, jeb
1262  do i = isb, ieb
1263  do k = ks, kref
1264  pres(k,i,j) = p00 * ( dens(k,i,j) * rtot(k,i,j) * pott(k,i,j) / p00 )**cpovcv(k,i,j)
1265  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot(k,i,j) )
1266  enddo
1267  enddo
1268  enddo
1269 
1270  return
1272 
1273  !-----------------------------------------------------------------------------
1276  dens, &
1277  pott, &
1278  pres, &
1279  temp, &
1280  qv, &
1281  qc, &
1282  pott_sfc, &
1283  pres_sfc, &
1284  temp_sfc, &
1285  qv_sfc, &
1286  qc_sfc )
1287  use scale_process, only: &
1288  prc_mpistop
1289  implicit none
1290 
1291  real(RP), intent(out) :: dens(KA)
1292  real(RP), intent(out) :: pott(KA)
1293  real(RP), intent(out) :: pres(KA)
1294  real(RP), intent(in) :: temp(KA)
1295  real(RP), intent(in) :: qv (KA)
1296  real(RP), intent(in) :: qc (KA)
1297 
1298  real(RP), intent(out) :: pott_sfc
1299  real(RP), intent(in) :: pres_sfc
1300  real(RP), intent(in) :: temp_sfc
1301  real(RP), intent(in) :: qv_sfc
1302  real(RP), intent(in) :: qc_sfc
1303 
1304  real(RP) :: dens_sfc
1305 
1306  real(RP) :: Rtot_sfc
1307  real(RP) :: CVtot_sfc
1308  real(RP) :: CPtot_sfc
1309  real(RP) :: Rtot
1310  real(RP) :: CVtot
1311  real(RP) :: CPtot
1312 
1313  real(RP) :: RovCP_sfc
1314  real(RP) :: dens_s, dhyd, dgrd
1315  integer :: ite
1316  logical :: converged
1317  !---------------------------------------------------------------------------
1318 
1319  !--- from surface to lowermost atmosphere
1320 
1321  rtot_sfc = rdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1322  + rvap * qv_sfc
1323  cvtot_sfc = cvdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1324  + cv_qv * qv_sfc &
1325  + cv_qc * qc_sfc
1326  cptot_sfc = cpdry * ( 1.0_rp - qv_sfc - qc_sfc ) &
1327  + cp_qv * qv_sfc &
1328  + cv_qc * qc_sfc
1329 
1330  rtot = rdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
1331  + rvap * qv(ks)
1332  cvtot = cvdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
1333  + cv_qv * qv(ks) &
1334  + cv_qc * qc(ks)
1335  cptot = cpdry * ( 1.0_rp - qv(ks) - qc(ks) ) &
1336  + cp_qv * qv(ks) &
1337  + cv_qc * qc(ks)
1338 
1339  ! density at surface
1340  rovcp_sfc = rtot_sfc / cptot_sfc
1341  dens_sfc = pres_sfc / ( rtot_sfc * temp_sfc )
1342  pott_sfc = temp_sfc * ( p00/pres_sfc )**rovcp_sfc
1343 
1344  ! make density at lowermost cell center
1345  dens_s = 0.0_rp
1346  dens(ks) = dens_sfc ! first guess
1347 
1348  converged = .false.
1349  do ite = 1, itelim
1350  if ( abs(dens(ks)-dens_s) <= criteria ) then
1351  converged = .true.
1352  exit
1353  endif
1354 
1355  dens_s = dens(ks)
1356 
1357  dhyd = + ( dens_sfc * rtot_sfc * temp_sfc &
1358  - dens_s * rtot * temp(ks) ) / grid_cz(ks) & ! dp/dz
1359  - grav * 0.5_rp * ( dens_sfc + dens_s ) ! rho*g
1360 
1361  dgrd = - rtot * temp(ks) / grid_cz(ks) &
1362  - 0.5_rp * grav
1363 
1364  dens(ks) = dens_s - dhyd/dgrd
1365 
1366  if( dens(ks)*0.0_rp /= 0.0_rp) exit
1367  enddo
1368 
1369  if ( .NOT. converged ) then
1370  write(*,*) 'xxx [buildrho bytemp 1D sfc] iteration not converged!', &
1371  dens(ks),ite,dens_s,dhyd,dgrd
1372  call prc_mpistop
1373  endif
1374 
1375  !--- from lowermost atmosphere to top of atmosphere
1376  call atmos_hydrostatic_buildrho_bytemp_atmos_1d( dens(:), & ! [INOUT]
1377  pott(:), & ! [OUT]
1378  pres(:), & ! [OUT]
1379  temp(:), & ! [IN]
1380  qv(:), & ! [IN]
1381  qc(:) ) ! [IN]
1382 
1383  return
1385 
1386  !-----------------------------------------------------------------------------
1389  dens, &
1390  pott, &
1391  pres, &
1392  temp, &
1393  qv, &
1394  qc, &
1395  pott_sfc, &
1396  pres_sfc, &
1397  temp_sfc, &
1398  qv_sfc, &
1399  qc_sfc )
1400  use scale_process, only: &
1401  prc_mpistop
1402  implicit none
1403 
1404  real(RP), intent(out) :: dens(KA,IA,JA)
1405  real(RP), intent(out) :: pott(KA,IA,JA)
1406  real(RP), intent(out) :: pres(KA,IA,JA)
1407  real(RP), intent(in) :: temp(KA,IA,JA)
1408  real(RP), intent(in) :: qv (KA,IA,JA)
1409  real(RP), intent(in) :: qc (KA,IA,JA)
1410 
1411  real(RP), intent(out) :: pott_sfc(1,IA,JA)
1412  real(RP), intent(in) :: pres_sfc(1,IA,JA)
1413  real(RP), intent(in) :: temp_sfc(1,IA,JA)
1414  real(RP), intent(in) :: qv_sfc (1,IA,JA)
1415  real(RP), intent(in) :: qc_sfc (1,IA,JA)
1416 
1417  real(RP) :: dens_sfc (1,IA,JA)
1418 
1419  real(RP) :: Rtot_sfc (IA,JA)
1420  real(RP) :: CVtot_sfc (IA,JA)
1421  real(RP) :: CPtot_sfc (IA,JA)
1422  real(RP) :: Rtot (IA,JA)
1423  real(RP) :: CVtot (IA,JA)
1424  real(RP) :: CPtot (IA,JA)
1425 
1426  real(RP) :: RovCP_sfc
1427  real(RP) :: DZ
1428  real(RP) :: dens_s, dhyd, dgrd
1429  integer :: ite
1430  logical :: converged
1431 
1432  integer :: i, j
1433  !---------------------------------------------------------------------------
1434 
1435  !--- from surface to lowermost atmosphere
1436 
1437  do j = jsb, jeb
1438  do i = isb, ieb
1439  rtot_sfc(i,j) = rdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1440  + rvap * qv_sfc(1,i,j)
1441  cvtot_sfc(i,j) = cvdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1442  + cv_qv * qv_sfc(1,i,j) &
1443  + cv_qc * qc_sfc(1,i,j)
1444  cptot_sfc(i,j) = cpdry * ( 1.0_rp - qv_sfc(1,i,j) - qc_sfc(1,i,j) ) &
1445  + cp_qv * qv_sfc(1,i,j) &
1446  + cv_qc * qc_sfc(1,i,j)
1447  enddo
1448  enddo
1449 
1450  do j = jsb, jeb
1451  do i = isb, ieb
1452  rtot(i,j) = rdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
1453  + rvap * qv(ks,i,j)
1454  cvtot(i,j) = cvdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
1455  + cv_qv * qv(ks,i,j) &
1456  + cv_qc * qc(ks,i,j)
1457  cptot(i,j) = cpdry * ( 1.0_rp - qv(ks,i,j) - qc(ks,i,j) ) &
1458  + cp_qv * qv(ks,i,j) &
1459  + cv_qc * qc(ks,i,j)
1460  enddo
1461  enddo
1462 
1463  ! density at surface
1464  do j = jsb, jeb
1465  do i = isb, ieb
1466  rovcp_sfc = rtot_sfc(i,j) / cptot_sfc(i,j)
1467  dens_sfc(1,i,j) = pres_sfc(1,i,j) / ( rtot_sfc(i,j) * temp_sfc(1,i,j) )
1468  pott_sfc(1,i,j) = temp_sfc(1,i,j) / ( p00/pres_sfc(1,i,j) )**rovcp_sfc
1469  enddo
1470  enddo
1471 
1472  ! make density at lowermost cell center
1473  do j = jsb, jeb
1474  do i = isb, ieb
1475  dz = real_cz(ks,i,j) - real_fz(ks-1,i,j)
1476 
1477  dens_s = 0.0_rp
1478  dens(ks,i,j) = dens_sfc(1,i,j) ! first guess
1479 
1480  converged = .false.
1481  do ite = 1, itelim
1482  if ( abs(dens(ks,i,j)-dens_s) <= criteria ) then
1483  converged = .true.
1484  exit
1485  endif
1486 
1487  dens_s = dens(ks,i,j)
1488 
1489  dhyd = + ( dens_sfc(1,i,j) * rtot_sfc(i,j) * temp_sfc(1,i,j) &
1490  - dens_s * rtot(i,j) * temp(ks,i,j) ) / dz & ! dp/dz
1491  - grav * 0.5_rp * ( dens_sfc(1,i,j) + dens_s ) ! rho*g
1492 
1493  dgrd = - rtot(i,j) * temp(ks,i,j) / dz &
1494  - 0.5_rp * grav
1495 
1496  dens(ks,i,j) = dens_s - dhyd/dgrd
1497 
1498  if( dens(ks,i,j)*0.0_rp /= 0.0_rp) exit
1499  enddo
1500 
1501  if ( .NOT. converged ) then
1502  write(*,*) 'xxx [buildrho bytemp 3D sfc] iteration not converged!', &
1503  i,j,dens(ks,i,j),ite,dens_s,dhyd,dgrd
1504  call prc_mpistop
1505  endif
1506  enddo
1507  enddo
1508 
1509  !--- from lowermost atmosphere to top of atmosphere
1510  call atmos_hydrostatic_buildrho_bytemp_atmos_3d( dens(:,:,:), & ! [INOUT]
1511  pott(:,:,:), & ! [OUT]
1512  pres(:,:,:), & ! [OUT]
1513  temp(:,:,:), & ! [IN]
1514  qv(:,:,:), & ! [IN]
1515  qc(:,:,:) ) ! [IN]
1516 
1517  return
1519 
1520  !-----------------------------------------------------------------------------
1523  dens, &
1524  pott, &
1525  pres, &
1526  temp, &
1527  qv, &
1528  qc )
1529  use scale_process, only: &
1530  prc_mpistop
1531  implicit none
1532 
1533  real(RP), intent(inout) :: dens(KA)
1534  real(RP), intent(out) :: pott(KA)
1535  real(RP), intent(out) :: pres(KA)
1536  real(RP), intent(in) :: temp(KA)
1537  real(RP), intent(in) :: qv (KA)
1538  real(RP), intent(in) :: qc (KA)
1539 
1540  real(RP) :: Rtot (KA)
1541  real(RP) :: CVtot (KA)
1542  real(RP) :: CPtot (KA)
1543 
1544  real(RP) :: RovCP
1545  real(RP) :: dens_s, dhyd, dgrd
1546  integer :: ite
1547  logical :: converged
1548 
1549  integer :: k
1550  !---------------------------------------------------------------------------
1551 
1552  do k = ks, ke
1553  rtot(k) = rdry * ( 1.0_rp - qv(k) - qc(k) ) &
1554  + rvap * qv(k)
1555  cvtot(k) = cvdry * ( 1.0_rp - qv(k) - qc(k) ) &
1556  + cv_qv * qv(k) &
1557  + cv_qc * qc(k)
1558  cptot(k) = cpdry * ( 1.0_rp - qv(k) - qc(k) ) &
1559  + cp_qv * qv(k) &
1560  + cv_qc * qc(k)
1561  enddo
1562 
1563  do k = ks+1, ke
1564 
1565  dens_s = 0.0_rp
1566  dens(k) = dens(k-1) ! first guess
1567 
1568  converged = .false.
1569  do ite = 1, itelim
1570  if ( abs(dens(k)-dens_s) <= criteria ) then
1571  converged = .true.
1572  exit
1573  endif
1574 
1575  dens_s = dens(k)
1576 
1577  dhyd = + ( dens(k-1) * rtot(k-1) * temp(k-1) &
1578  - dens_s * rtot(k ) * temp(k ) ) / grid_fdz(k-1) & ! dp/dz
1579  - grav * 0.5_rp * ( dens(k-1) + dens_s ) ! rho*g
1580 
1581  dgrd = - rtot(k) * temp(k) / grid_fdz(k-1) &
1582  - 0.5_rp * grav
1583 
1584  dens(k) = dens_s - dhyd/dgrd
1585 
1586  if( dens(k)*0.0_rp /= 0.0_rp) exit
1587  enddo
1588 
1589  if ( .NOT. converged ) then
1590  write(*,*) 'xxx [buildrho bytemp 1D atmos] iteration not converged!', &
1591  k,dens(k),ite,dens_s,dhyd,dgrd
1592  call prc_mpistop
1593  endif
1594  enddo
1595 
1596  do k = ks, ke
1597  rovcp = rtot(k) / cptot(k)
1598  pres(k) = dens(k) * rtot(k) * temp(k)
1599  pott(k) = temp(k) * ( p00 / pres(k) )**rovcp
1600  enddo
1601 
1602  return
1604 
1605  !-----------------------------------------------------------------------------
1608  dens, &
1609  pott, &
1610  pres, &
1611  temp, &
1612  qv, &
1613  qc )
1614  use scale_process, only: &
1615  prc_mpistop
1616  implicit none
1617 
1618  real(RP), intent(inout) :: dens(KA,IA,JA)
1619  real(RP), intent(out) :: pott(KA,IA,JA)
1620  real(RP), intent(out) :: pres(KA,IA,JA)
1621  real(RP), intent(in) :: temp(KA,IA,JA)
1622  real(RP), intent(in) :: qv (KA,IA,JA)
1623  real(RP), intent(in) :: qc (KA,IA,JA)
1624 
1625  real(RP) :: Rtot (KA,IA,JA)
1626  real(RP) :: CVtot (KA,IA,JA)
1627  real(RP) :: CPtot (KA,IA,JA)
1628 
1629  real(RP) :: RovCP
1630  real(RP) :: DZ
1631  real(RP) :: dens_s, dhyd, dgrd
1632  integer :: ite
1633  logical :: converged
1634 
1635  integer :: k, i, j
1636  !---------------------------------------------------------------------------
1637 
1638  do j = jsb, jeb
1639  do i = isb, ieb
1640  do k = ks, ke
1641  rtot(k,i,j) = rdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1642  + rvap * qv(k,i,j)
1643  cvtot(k,i,j) = cvdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1644  + cv_qv * qv(k,i,j) &
1645  + cv_qc * qc(k,i,j)
1646  cptot(k,i,j) = cpdry * ( 1.0_rp - qv(k,i,j) - qc(k,i,j) ) &
1647  + cp_qv * qv(k,i,j) &
1648  + cv_qc * qc(k,i,j)
1649  enddo
1650  enddo
1651  enddo
1652 
1653  do j = jsb, jeb
1654  do i = isb, ieb
1655  do k = ks+1, ke
1656  dz = real_cz(k,i,j) - real_cz(k-1,i,j)
1657 
1658  dens_s = 0.0_rp
1659  dens(k,i,j) = dens(k-1,i,j) ! first guess
1660 
1661  converged = .false.
1662  do ite = 1, itelim
1663  if ( abs(dens(k,i,j)-dens_s) <= criteria ) then
1664  converged = .true.
1665  exit
1666  endif
1667 
1668  dens_s = dens(k,i,j)
1669 
1670  dhyd = + ( dens(k-1,i,j) * rtot(k-1,i,j) * temp(k-1,i,j) &
1671  - dens_s * rtot(k ,i,j) * temp(k ,i,j) ) / dz & ! dpdz
1672  - grav * 0.5_rp * ( dens(k-1,i,j) + dens_s ) ! rho*g
1673 
1674  dgrd = - rtot(k,i,j) * temp(k,i,j) / dz &
1675  - 0.5_rp * grav
1676 
1677  dens(k,i,j) = dens_s - dhyd/dgrd
1678 
1679  if( dens(k,i,j)*0.0_rp /= 0.0_rp) exit
1680  enddo
1681 
1682  if ( .NOT. converged ) then
1683  write(*,*) 'xxx [buildrho bytemp 3D atmos] iteration not converged!', &
1684  k,i,j,dens(k,i,j),ite,dens_s,dhyd,dgrd
1685  call prc_mpistop
1686  endif
1687  enddo
1688  enddo
1689  enddo
1690 
1691  do j = jsb, jeb
1692  do i = isb, ieb
1693  do k = ks, ke
1694  rovcp = rtot(k,i,j) / cptot(k,i,j)
1695  pres(k,i,j) = dens(k,i,j) * rtot(k,i,j) * temp(k,i,j)
1696  pott(k,i,j) = temp(k,i,j) * ( p00 / pres(k,i,j) )**rovcp
1697  enddo
1698  enddo
1699  enddo
1700 
1701  return
1703 
1704  !-----------------------------------------------------------------------------
1707  mslp, &
1708  pres, &
1709  temp, &
1710  dz )
1711  implicit none
1712 
1713  real(RP), intent(out) :: mslp
1714  real(RP), intent(in) :: pres
1715  real(RP), intent(in) :: temp
1716  real(RP), intent(in) :: dz
1717 
1718  ! work
1719  real(RP) :: TM
1720  !---------------------------------------------------------------------------
1721 
1722  tm = temp + laps * dz * 0.5_rp ! column-mean air temperature
1723 
1724  ! barometric law
1725  mslp = pres * exp( grav * dz / ( rdry * tm ) )
1726 
1727  return
1729 
1730  !-----------------------------------------------------------------------------
1732  subroutine atmos_hydrostatic_barometric_law_mslp_2d( &
1733  mslp, &
1734  pres, &
1735  temp, &
1736  dz )
1737  implicit none
1738 
1739  real(RP), intent(out) :: mslp(IA,JA)
1740  real(RP), intent(in) :: pres(IA,JA)
1741  real(RP), intent(in) :: temp(IA,JA)
1742  real(RP), intent(in) :: dz (IA,JA)
1743 
1744  ! work
1745  real(RP) :: TM
1746 
1747  integer :: i, j
1748  !---------------------------------------------------------------------------
1749 
1750  do j = jsb, jeb
1751  do i = isb, ieb
1752  tm = temp(i,j) + laps * dz(i,j) * 0.5_rp ! column-mean air temperature
1753 
1754  ! barometric law
1755  mslp(i,j) = pres(i,j) * exp( grav * dz(i,j) / ( rdry * tm ) )
1756  enddo
1757  enddo
1758 
1759  return
1760  end subroutine atmos_hydrostatic_barometric_law_mslp_2d
1761 
1762  !-----------------------------------------------------------------------------
1764  subroutine atmos_hydrostatic_barometric_law_pres_0d( &
1765  pres, &
1766  mslp, &
1767  temp, &
1768  dz )
1769  implicit none
1770 
1771  real(RP), intent(out) :: pres
1772  real(RP), intent(in) :: mslp
1773  real(RP), intent(in) :: temp
1774  real(RP), intent(in) :: dz
1775 
1776  ! work
1777  real(RP) :: TM
1778  !---------------------------------------------------------------------------
1779 
1780  tm = temp + laps * dz * 0.5_rp ! column-mean air temperature
1781 
1782  ! barometric law
1783  pres = mslp / exp( grav * dz / ( rdry * tm ) )
1784 
1785  return
1786  end subroutine atmos_hydrostatic_barometric_law_pres_0d
1787 
1788  !-----------------------------------------------------------------------------
1790  subroutine atmos_hydrostatic_barometric_law_pres_2d( &
1791  pres, &
1792  mslp, &
1793  temp, &
1794  dz )
1795  implicit none
1796 
1797  real(RP), intent(out) :: pres(IA,JA)
1798  real(RP), intent(in) :: mslp(IA,JA)
1799  real(RP), intent(in) :: temp(IA,JA)
1800  real(RP), intent(in) :: dz (IA,JA)
1801 
1802  ! work
1803  real(RP) :: TM
1804 
1805  integer :: i, j
1806  !---------------------------------------------------------------------------
1807 
1808  do j = jsb, jeb
1809  do i = isb, ieb
1810  tm = temp(i,j) + laps * dz(i,j) * 0.5_rp ! column-mean air temperature
1811 
1812  ! barometric law
1813  pres(i,j) = mslp(i,j) / exp( grav * dz(i,j) / ( rdry * tm ) )
1814  enddo
1815  enddo
1816 
1817  return
1818  end subroutine atmos_hydrostatic_barometric_law_pres_2d
1819 
1820 end module scale_atmos_hydrostatic
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)
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_1d(dens, pott, pres, temp, qv, qc)
Build up density from lowermost atmosphere (1D)
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
subroutine, public atmos_hydrostatic_buildrho_atmos_rev_3d(dens, temp, pres, pott, qv, qc, dz, kref_in)
Build up density from lowermost atmosphere (3D)
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:68
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:67
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
real(rp), dimension(qa_max), public tracer_cv
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
integer, public ieb
module grid index
real(rp), dimension(qa_max), public tracer_cp
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
module TRACER
integer, public ia
of whole cells: x, local, with HALO
subroutine, public comm_horizontal_mean(varmean, var)
calculate horizontal mean (global total with communication)
Definition: scale_comm.F90:482
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
subroutine atmos_hydrostatic_barometric_law_mslp_0d(mslp, pres, temp, dz)
Calculate mean sea-level pressure from barometric law (0D)
module GRID (real space)
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_lapsdry
dry adiabatic lapse rate [K/m]
Definition: scale_const.F90:61
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
module COMMUNICATION
Definition: scale_comm.F90:23
module ATMOSPHERE / Hydrostatic barance
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module PROCESS
subroutine atmos_hydrostatic_buildrho_atmos_1d(dens, temp, pres, pott, qv, qc)
Build up density from lowermost atmosphere (1D)
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
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)
integer, public ks
start point of inner domain: z, local
subroutine atmos_hydrostatic_buildrho_atmos_3d(dens, temp, pres, pott, qv, qc, dz, kref_in)
Build up density from lowermost atmosphere (3D)
module GRID (cartesian)
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)
module profiler
Definition: scale_prof.F90:10
real(rp), public const_eps
small number
Definition: scale_const.F90:36
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)
subroutine atmos_hydrostatic_buildrho_real_3d(dens, temp, pres, pott, qv, qc)
Build up density from surface (3D), not to reverse from TOA.
module PRECISION
subroutine atmos_hydrostatic_buildrho_bytemp_atmos_3d(dens, pott, pres, temp, qv, qc)
Build up density from lowermost atmosphere (3D)
integer, public isb
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
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:66
integer, public jsb
character(len=h_short), public const_thermodyn_type
internal energy type
Definition: scale_const.F90:98
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)
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
subroutine, public atmos_hydrostatic_setup
Setup.
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)
integer, public ja
of whole cells: y, local, with HALO