SCALE-RM
scale_atmos_sub_saturation.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
15 #include "inc_openmp.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_stdio
23  use scale_prof
25  use scale_tracer
26 
27  use scale_const, only: &
28  rvap => const_rvap, &
29  psat0 => const_psat0, &
30  epsvap => const_epsvap, &
31  tem00 => const_tem00
32  !-----------------------------------------------------------------------------
33  implicit none
34  private
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public procedure
38  !
39  public :: atmos_saturation_setup
40 
41  public :: atmos_saturation_alpha
42 
43  public :: atmos_saturation_psat_all
44  public :: atmos_saturation_psat_liq
45  public :: atmos_saturation_psat_ice
46 
47  public :: atmos_saturation_pres2qsat_all
48  public :: atmos_saturation_pres2qsat_liq
49  public :: atmos_saturation_pres2qsat_ice
50 
51  public :: atmos_saturation_dens2qsat_all
52  public :: atmos_saturation_dens2qsat_liq
53  public :: atmos_saturation_dens2qsat_ice
54 
55  public :: atmos_saturation_dalphadt
56 
61 
62  interface atmos_saturation_alpha
63  module procedure atmos_saturation_alpha_0d
64  module procedure atmos_saturation_alpha_1d
65  module procedure atmos_saturation_alpha_3d
66  end interface atmos_saturation_alpha
67 
68  interface atmos_saturation_psat_all
69  module procedure atmos_saturation_psat_all_0d
70  module procedure atmos_saturation_psat_all_1d
71  module procedure atmos_saturation_psat_all_3d
72  end interface atmos_saturation_psat_all
73  interface atmos_saturation_psat_liq
74  module procedure atmos_saturation_psat_liq_0d
75  module procedure atmos_saturation_psat_liq_1d
76  module procedure atmos_saturation_psat_liq_3d
77  end interface atmos_saturation_psat_liq
78  interface atmos_saturation_psat_ice
79  module procedure atmos_saturation_psat_ice_0d
80  module procedure atmos_saturation_psat_ice_1d
81  module procedure atmos_saturation_psat_ice_3d
82  end interface atmos_saturation_psat_ice
83 
84  interface atmos_saturation_pres2qsat_all
85  module procedure atmos_saturation_pres2qsat_all_0d
86  module procedure atmos_saturation_pres2qsat_all_1d
87  module procedure atmos_saturation_pres2qsat_all_2d
88  module procedure atmos_saturation_pres2qsat_all_3d
89  module procedure atmos_saturation_pres2qsat_all_3d_k
90  end interface atmos_saturation_pres2qsat_all
91  interface atmos_saturation_pres2qsat_liq
92  module procedure atmos_saturation_pres2qsat_liq_0d
93  module procedure atmos_saturation_pres2qsat_liq_1d
94  module procedure atmos_saturation_pres2qsat_liq_3d
95  end interface atmos_saturation_pres2qsat_liq
96  interface atmos_saturation_pres2qsat_ice
97  module procedure atmos_saturation_pres2qsat_ice_0d
98  module procedure atmos_saturation_pres2qsat_ice_1d
99  module procedure atmos_saturation_pres2qsat_ice_3d
100  end interface atmos_saturation_pres2qsat_ice
101 
102  interface atmos_saturation_dens2qsat_all
103  module procedure atmos_saturation_dens2qsat_all_0d
104  module procedure atmos_saturation_dens2qsat_all_1d
105  module procedure atmos_saturation_dens2qsat_all_3d
106  end interface atmos_saturation_dens2qsat_all
107  interface atmos_saturation_dens2qsat_liq
108  module procedure atmos_saturation_dens2qsat_liq_0d
109  module procedure atmos_saturation_dens2qsat_liq_1d
110  module procedure atmos_saturation_dens2qsat_liq_3d
111  end interface atmos_saturation_dens2qsat_liq
112  interface atmos_saturation_dens2qsat_ice
113  module procedure atmos_saturation_dens2qsat_ice_0d
114  module procedure atmos_saturation_dens2qsat_ice_1d
115  module procedure atmos_saturation_dens2qsat_ice_3d
116  end interface atmos_saturation_dens2qsat_ice
117 
118  interface atmos_saturation_dalphadt
119  module procedure atmos_saturation_dalphadt_0d
120  module procedure atmos_saturation_dalphadt_1d
121  module procedure atmos_saturation_dalphadt_3d
122  end interface atmos_saturation_dalphadt
123 
124  !-----------------------------------------------------------------------------
125  !
126  !++ Public parameters & variables
127  !
128  real(RP), public :: cpovr_liq
129  real(RP), public :: cpovr_ice
130  real(RP), public :: cvovr_liq
131  real(RP), public :: cvovr_ice
132  real(RP), public :: lovr_liq
133  real(RP), public :: lovr_ice
134 
135  !-----------------------------------------------------------------------------
136  !
137  !++ Private procedure
138  !
139  !-----------------------------------------------------------------------------
140  !
141  !++ Private parameters & variables
142  !
143  real(RP), private, parameter :: tem_min = 10.0_rp
144 
145  real(RP), private, save :: atmos_saturation_ulimit_temp = 273.15_rp
146  real(RP), private, save :: atmos_saturation_llimit_temp = 233.15_rp
147 
148  real(RP), private, save :: rtem00
149  real(RP), private, save :: dalphadt_const
150 
151  !-----------------------------------------------------------------------------
152 contains
153  !-----------------------------------------------------------------------------
155  subroutine atmos_saturation_setup
156  use scale_process, only: &
158  use scale_const, only: &
159  cpvap => const_cpvap, &
160  cvvap => const_cvvap, &
161  cl => const_cl, &
162  ci => const_ci, &
163  lhv00 => const_lhv00, &
164  lhs00 => const_lhs00, &
165  lhv0 => const_lhv0, &
166  lhs0 => const_lhs0, &
168  implicit none
169 
170  namelist / param_atmos_saturation / &
171  atmos_saturation_ulimit_temp, &
172  atmos_saturation_llimit_temp
173 
174  integer :: ierr
175  !---------------------------------------------------------------------------
176 
177  if( io_l ) write(io_fid_log,*)
178  if( io_l ) write(io_fid_log,*) '++++++ Module[SATURATION] / Categ[ATMOS SHARE] / Origin[SCALElib]'
179 
180  !--- read namelist
181  rewind(io_fid_conf)
182  read(io_fid_conf,nml=param_atmos_saturation,iostat=ierr)
183  if( ierr < 0 ) then !--- missing
184  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
185  elseif( ierr > 0 ) then !--- fatal error
186  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_SATURATION. Check!'
187  call prc_mpistop
188  endif
189  if( io_nml ) write(io_fid_nml,nml=param_atmos_saturation)
190 
191  rtem00 = 1.0_rp / tem00
192 
193  if ( const_thermodyn_type == 'EXACT' ) then
194 
195  cpovr_liq = ( cpvap - cl ) / rvap
196  cpovr_ice = ( cpvap - ci ) / rvap
197  cvovr_liq = ( cvvap - cl ) / rvap
198  cvovr_ice = ( cvvap - ci ) / rvap
199 
200  lovr_liq = lhv00 / rvap
201  lovr_ice = lhs00 / rvap
202 
203  elseif( const_thermodyn_type == 'SIMPLE' ) then
204 
205  cpovr_liq = 0.0_rp
206  cpovr_ice = 0.0_rp
207  cvovr_liq = 0.0_rp
208  cvovr_ice = 0.0_rp
209 
210  lovr_liq = lhv0 / rvap
211  lovr_ice = lhs0 / rvap
212 
213  endif
214 
215  dalphadt_const = 1.0_rp / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
216 
217  if( io_l ) write(io_fid_log,*)
218  if( io_l ) write(io_fid_log,'(1x,A,F7.2,A,F7.2)') '*** Temperature range for liquid/ice mixture : ', &
219  atmos_saturation_llimit_temp, ' - ', &
220  atmos_saturation_ulimit_temp
221 
222  return
223  end subroutine atmos_saturation_setup
224 
225  !-----------------------------------------------------------------------------
227  subroutine atmos_saturation_alpha_0d( &
228  alpha, &
229  temp )
230  implicit none
231 
232  real(RP), intent(out) :: alpha
233  real(RP), intent(in) :: temp
234  !---------------------------------------------------------------------------
235 
236  alpha = ( temp - atmos_saturation_llimit_temp ) &
237  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
238 
239  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
240 
241  return
242  end subroutine atmos_saturation_alpha_0d
243 
244  !-----------------------------------------------------------------------------
246  subroutine atmos_saturation_alpha_1d( &
247  alpha, &
248  temp )
249  implicit none
250 
251  real(RP), intent(out) :: alpha(KA)
252  real(RP), intent(in) :: temp (KA)
253 
254  integer :: k
255  !---------------------------------------------------------------------------
256 
257  do k = ks, ke
258 
259  alpha(k) = ( temp(k) - atmos_saturation_llimit_temp ) &
260  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
261  alpha(k) = min( max( alpha(k), 0.0_rp ), 1.0_rp )
262 
263  enddo
264 
265  return
266  end subroutine atmos_saturation_alpha_1d
267 
268  !-----------------------------------------------------------------------------
270  subroutine atmos_saturation_alpha_3d( &
271  alpha, &
272  temp )
273  implicit none
274 
275  real(RP), intent(out) :: alpha(KA,IA,JA)
276  real(RP), intent(in) :: temp (KA,IA,JA)
277 
278  integer :: k, i, j
279  !---------------------------------------------------------------------------
280 
281  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
282  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,alpha,temp,ATMOS_SATURATION_LLIMIT_TEMP,ATMOS_SATURATION_ULIMIT_TEMP)
283  do j = jsb, jeb
284  do i = isb, ieb
285  do k = ks, ke
286 
287  alpha(k,i,j) = ( temp(k,i,j) - atmos_saturation_llimit_temp ) &
288  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
289  alpha(k,i,j) = min( max( alpha(k,i,j), 0.0_rp ), 1.0_rp )
290 
291  enddo
292  enddo
293  enddo
294 
295  return
296  end subroutine atmos_saturation_alpha_3d
297 
298  !-----------------------------------------------------------------------------
300  subroutine atmos_saturation_psat_all_0d( &
301  psat, &
302  temp )
303  implicit none
304 
305  real(RP), intent(out) :: psat
306  real(RP), intent(in) :: temp
307 
308  real(RP) :: alpha, psatl, psati
309  !---------------------------------------------------------------------------
310 
311  call atmos_saturation_alpha ( alpha, temp )
312  call atmos_saturation_psat_liq( psatl, temp )
313  call atmos_saturation_psat_ice( psati, temp )
314 
315  psat = psatl * ( alpha ) &
316  + psati * ( 1.0_rp - alpha )
317 
318  return
319  end subroutine atmos_saturation_psat_all_0d
320 
321  !-----------------------------------------------------------------------------
323  subroutine atmos_saturation_psat_all_1d( &
324  psat, &
325  temp )
326  implicit none
327 
328  real(RP), intent(out) :: psat(KA)
329  real(RP), intent(in) :: temp(KA)
330 
331  real(RP) :: alpha(KA), psatl(KA), psati(KA)
332 
333  integer :: k
334  !---------------------------------------------------------------------------
335 
336  call atmos_saturation_alpha ( alpha(:), temp(:) )
337  call atmos_saturation_psat_liq( psatl(:), temp(:) )
338  call atmos_saturation_psat_ice( psati(:), temp(:) )
339 
340  do k = ks, ke
341  psat(k) = psatl(k) * ( alpha(k) ) &
342  + psati(k) * ( 1.0_rp - alpha(k) )
343  enddo
344 
345  return
346  end subroutine atmos_saturation_psat_all_1d
347 
348  !-----------------------------------------------------------------------------
350  subroutine atmos_saturation_psat_all_3d( &
351  psat, &
352  temp )
353  implicit none
354 
355  real(RP), intent(out) :: psat(KA,IA,JA)
356  real(RP), intent(in) :: temp(KA,IA,JA)
357 
358  real(RP) :: alpha, psatl, psati
359 
360  integer :: k, i, j
361  !---------------------------------------------------------------------------
362 
363  !$omp parallel do private(i,j,k,alpha,psatl,psati) OMP_SCHEDULE_ collapse(2)
364  do j = jsb, jeb
365  do i = isb, ieb
366  do k = ks, ke
367 
368  alpha = ( temp(k,i,j) - atmos_saturation_llimit_temp ) &
369  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
370  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
371 
372  psatl = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_liq &
373  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
374 
375  psati = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_ice &
376  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
377 
378  psat(k,i,j) = psatl * ( alpha ) &
379  + psati * ( 1.0_rp - alpha )
380 
381  enddo
382  enddo
383  enddo
384 
385  return
386  end subroutine atmos_saturation_psat_all_3d
387 
388  !-----------------------------------------------------------------------------
390  subroutine atmos_saturation_psat_liq_0d( &
391  psat, &
392  temp )
393  implicit none
394 
395  real(RP), intent(out) :: psat
396  real(RP), intent(in) :: temp
397  !---------------------------------------------------------------------------
398 
399  psat = psat0 * ( temp * rtem00 )**cpovr_liq &
400  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp ) )
401 
402  return
403  end subroutine atmos_saturation_psat_liq_0d
404 
405  !-----------------------------------------------------------------------------
407  subroutine atmos_saturation_psat_liq_1d( &
408  psat, &
409  temp )
410  implicit none
411 
412  real(RP), intent(out) :: psat(KA)
413  real(RP), intent(in) :: temp(KA)
414 
415  integer :: k
416  !---------------------------------------------------------------------------
417 
418  do k = ks, ke
419  psat(k) = psat0 * ( temp(k) * rtem00 )**cpovr_liq &
420  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k) ) )
421  enddo
422 
423  return
424  end subroutine atmos_saturation_psat_liq_1d
425 
426  !-----------------------------------------------------------------------------
428  subroutine atmos_saturation_psat_liq_3d( &
429  psat, &
430  temp )
431  implicit none
432 
433  real(RP), intent(out) :: psat(KA,IA,JA)
434  real(RP), intent(in) :: temp(KA,IA,JA)
435 
436  integer :: k, i, j
437  !---------------------------------------------------------------------------
438 
439  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
440  do j = jsb, jeb
441  do i = isb, ieb
442  do k = ks, ke
443  psat(k,i,j) = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_liq &
444  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
445  enddo
446  enddo
447  enddo
448 
449  return
450  end subroutine atmos_saturation_psat_liq_3d
451 
452  !-----------------------------------------------------------------------------
454  subroutine atmos_saturation_psat_ice_0d( &
455  psat, &
456  temp )
457  implicit none
458 
459  real(RP), intent(out) :: psat
460  real(RP), intent(in) :: temp
461  !---------------------------------------------------------------------------
462 
463  psat = psat0 * ( temp * rtem00 )**cpovr_ice &
464  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp ) )
465 
466  return
467  end subroutine atmos_saturation_psat_ice_0d
468 
469  !-----------------------------------------------------------------------------
471  subroutine atmos_saturation_psat_ice_1d( &
472  psat, &
473  temp )
474  implicit none
475 
476  real(RP), intent(out) :: psat(KA)
477  real(RP), intent(in) :: temp(KA)
478 
479  integer :: k
480  !---------------------------------------------------------------------------
481 
482  do k = ks, ke
483  psat(k) = psat0 * ( temp(k) * rtem00 )**cpovr_ice &
484  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k) ) )
485  enddo
486 
487  return
488  end subroutine atmos_saturation_psat_ice_1d
489 
490  !-----------------------------------------------------------------------------
492  subroutine atmos_saturation_psat_ice_3d( &
493  psat, &
494  temp )
495  implicit none
496 
497  real(RP), intent(out) :: psat(KA,IA,JA)
498  real(RP), intent(in) :: temp(KA,IA,JA)
499 
500  integer :: k, i, j
501  !---------------------------------------------------------------------------
502 
503  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
504  do j = jsb, jeb
505  do i = isb, ieb
506  do k = ks, ke
507  psat(k,i,j) = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_ice &
508  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
509  enddo
510  enddo
511  enddo
512 
513  return
514  end subroutine atmos_saturation_psat_ice_3d
515 
516  !-----------------------------------------------------------------------------
518  subroutine atmos_saturation_pres2qsat_all_0d( &
519  qsat, &
520  temp, &
521  pres )
522  implicit none
523 
524  real(RP), intent(out) :: qsat
525  real(RP), intent(in) :: temp
526  real(RP), intent(in) :: pres
527 
528  real(RP) :: psat
529  !---------------------------------------------------------------------------
530 
531  call atmos_saturation_psat_all( psat, temp )
532 
533  qsat = epsvap * psat / ( pres - ( 1.0_rp-epsvap ) * psat )
534 
535  return
536  end subroutine atmos_saturation_pres2qsat_all_0d
537 
538  !-----------------------------------------------------------------------------
540  subroutine atmos_saturation_pres2qsat_all_1d( &
541  qsat, &
542  temp, &
543  pres )
544  implicit none
545 
546  real(RP), intent(out) :: qsat(KA)
547  real(RP), intent(in) :: temp(KA)
548  real(RP), intent(in) :: pres(KA)
549 
550  real(RP) :: psat(KA)
551 
552  integer :: k
553  !---------------------------------------------------------------------------
554 
555  call atmos_saturation_psat_all( psat(:), temp(:) )
556 
557  do k = ks, ke
558  qsat(k) = epsvap * psat(k) / ( pres(k) - ( 1.0_rp-epsvap ) * psat(k) )
559  enddo
560 
561  return
562  end subroutine atmos_saturation_pres2qsat_all_1d
563 
564  !-----------------------------------------------------------------------------
566  subroutine atmos_saturation_pres2qsat_all_2d( &
567  qsat, &
568  temp, &
569  pres )
570  implicit none
571 
572  real(RP), intent(out) :: qsat(IA,JA)
573  real(RP), intent(in) :: temp(IA,JA)
574  real(RP), intent(in) :: pres(IA,JA)
575 
576  real(RP) :: psat
577 
578  integer :: i, j
579  !---------------------------------------------------------------------------
580 
581  do j = jsb, jeb
582  do i = isb, ieb
583  call atmos_saturation_psat_all( psat, temp(i,j) )
584 
585  qsat(i,j) = epsvap * psat / ( pres(i,j) - ( 1.0_rp-epsvap ) * psat )
586  enddo
587  enddo
588 
589  return
590  end subroutine atmos_saturation_pres2qsat_all_2d
591 
592  !-----------------------------------------------------------------------------
594  subroutine atmos_saturation_pres2qsat_all_3d( &
595  qsat, &
596  temp, &
597  pres )
598  implicit none
599 
600  real(RP), intent(out) :: qsat(KA,IA,JA)
601  real(RP), intent(in) :: temp(KA,IA,JA)
602  real(RP), intent(in) :: pres(KA,IA,JA)
603 
604  real(RP) :: alpha, psatl, psati
605  real(RP) :: psat
606 
607  integer :: k, i, j
608  !---------------------------------------------------------------------------
609 
610  !$omp parallel do private(i,j,k,alpha,psatl,psati,psat) OMP_SCHEDULE_ collapse(2)
611  do j = jsb, jeb
612  do i = isb, ieb
613  do k = ks, ke
614 
615  alpha = ( temp(k,i,j) - atmos_saturation_llimit_temp ) &
616  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
617  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
618 
619  psatl = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_liq &
620  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
621 
622  psati = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_ice &
623  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
624 
625  psat = psatl * ( alpha ) &
626  + psati * ( 1.0_rp - alpha )
627 
628  qsat(k,i,j) = epsvap * psat / ( pres(k,i,j) - ( 1.0_rp-epsvap ) * psat )
629 
630  enddo
631  enddo
632  enddo
633 
634  return
635  end subroutine atmos_saturation_pres2qsat_all_3d
636 
637  !-----------------------------------------------------------------------------
639  subroutine atmos_saturation_pres2qsat_all_3d_k( &
640  qsat, &
641  temp, &
642  pres, &
643  knum )
644  implicit none
645 
646  integer, intent(in) :: knum
647  real(RP), intent(out) :: qsat(knum,IA,JA)
648  real(RP), intent(in) :: temp(knum,IA,JA)
649  real(RP), intent(in) :: pres(knum,IA,JA)
650 
651  real(RP) :: psat
652 
653  integer :: k, i, j
654  !---------------------------------------------------------------------------
655 
656  !$omp parallel do private(i,j,k,psat) OMP_SCHEDULE_ collapse(2)
657  do j = jsb, jeb
658  do i = isb, ieb
659  do k = 1, knum
660  call atmos_saturation_psat_all( psat, temp(k,i,j) )
661 
662  qsat(k,i,j) = epsvap * psat / ( pres(k,i,j) - ( 1.0_rp-epsvap ) * psat )
663  enddo
664  enddo
665  enddo
666 
667  return
668  end subroutine atmos_saturation_pres2qsat_all_3d_k
669 
670  !-----------------------------------------------------------------------------
672  subroutine atmos_saturation_pres2qsat_liq_0d( &
673  qsat, &
674  temp, &
675  pres )
676  implicit none
677 
678  real(RP), intent(out) :: qsat
679  real(RP), intent(in) :: temp
680  real(RP), intent(in) :: pres
681 
682  real(RP) :: psat
683  !---------------------------------------------------------------------------
684 
685  call atmos_saturation_psat_liq( psat, temp )
686 
687  qsat = epsvap * psat / ( pres - ( 1.0_rp-epsvap ) * psat )
688 
689  return
690  end subroutine atmos_saturation_pres2qsat_liq_0d
691 
692  !-----------------------------------------------------------------------------
694  subroutine atmos_saturation_pres2qsat_liq_1d( &
695  qsat, &
696  temp, &
697  pres )
698  implicit none
699 
700  real(RP), intent(out) :: qsat(KA)
701  real(RP), intent(in) :: temp(KA)
702  real(RP), intent(in) :: pres(KA)
703 
704  real(RP) :: psat(KA)
705 
706  integer :: k
707  !---------------------------------------------------------------------------
708 
709  call atmos_saturation_psat_liq( psat(:), temp(:) )
710 
711  do k = ks, ke
712  qsat(k) = epsvap * psat(k) / ( pres(k) - ( 1.0_rp-epsvap ) * psat(k) )
713  enddo
714 
715  return
716  end subroutine atmos_saturation_pres2qsat_liq_1d
717 
718  !-----------------------------------------------------------------------------
720  subroutine atmos_saturation_pres2qsat_liq_3d( &
721  qsat, &
722  temp, &
723  pres )
724  implicit none
725 
726  real(RP), intent(out) :: qsat(KA,IA,JA)
727  real(RP), intent(in) :: temp(KA,IA,JA)
728  real(RP), intent(in) :: pres(KA,IA,JA)
729 
730  real(RP) :: psat
731 
732  integer :: k, i, j
733  !---------------------------------------------------------------------------
734 
735  !$omp parallel do private(i,j,k,psat) OMP_SCHEDULE_ collapse(2)
736  do j = jsb, jeb
737  do i = isb, ieb
738  do k = ks, ke
739  psat = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_liq &
740  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
741 
742  qsat(k,i,j) = epsvap * psat / ( pres(k,i,j) - ( 1.0_rp-epsvap ) * psat )
743  enddo
744  enddo
745  enddo
746 
747  return
748  end subroutine atmos_saturation_pres2qsat_liq_3d
749 
750  !-----------------------------------------------------------------------------
752  subroutine atmos_saturation_pres2qsat_ice_0d( &
753  qsat, &
754  temp, &
755  pres )
756  implicit none
757 
758  real(RP), intent(out) :: qsat
759  real(RP), intent(in) :: temp
760  real(RP), intent(in) :: pres
761 
762  real(RP) :: psat
763  !---------------------------------------------------------------------------
764 
765  call atmos_saturation_psat_ice( psat, temp )
766 
767  qsat = epsvap * psat / ( pres - ( 1.0_rp-epsvap ) * psat )
768 
769  return
770  end subroutine atmos_saturation_pres2qsat_ice_0d
771 
772  !-----------------------------------------------------------------------------
774  subroutine atmos_saturation_pres2qsat_ice_1d( &
775  qsat, &
776  temp, &
777  pres )
778  implicit none
779 
780  real(RP), intent(out) :: qsat(KA)
781  real(RP), intent(in) :: temp(KA)
782  real(RP), intent(in) :: pres(KA)
783 
784  real(RP) :: psat(KA)
785 
786  integer :: k
787  !---------------------------------------------------------------------------
788 
789  call atmos_saturation_psat_ice( psat(:), temp(:) )
790 
791  do k = ks, ke
792  qsat(k) = epsvap * psat(k) / ( pres(k) - ( 1.0_rp-epsvap ) * psat(k) )
793  enddo
794 
795  return
796  end subroutine atmos_saturation_pres2qsat_ice_1d
797 
798  !-----------------------------------------------------------------------------
800  subroutine atmos_saturation_pres2qsat_ice_3d( &
801  qsat, &
802  temp, &
803  pres )
804  implicit none
805 
806  real(RP), intent(out) :: qsat(KA,IA,JA)
807  real(RP), intent(in) :: temp(KA,IA,JA)
808  real(RP), intent(in) :: pres(KA,IA,JA)
809 
810  real(RP) :: psat
811 
812  integer :: k, i, j
813  !---------------------------------------------------------------------------
814 
815  !$omp parallel do private(i,j,k,psat) OMP_SCHEDULE_ collapse(2)
816  do j = jsb, jeb
817  do i = isb, ieb
818  do k = ks, ke
819  psat = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_ice &
820  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
821 
822  qsat(k,i,j) = epsvap * psat / ( pres(k,i,j) - ( 1.0_rp-epsvap ) * psat )
823  enddo
824  enddo
825  enddo
826 
827  return
828  end subroutine atmos_saturation_pres2qsat_ice_3d
829 
830  !-----------------------------------------------------------------------------
832  subroutine atmos_saturation_dens2qsat_all_0d( &
833  qsat, &
834  temp, &
835  dens )
836  implicit none
837 
838  real(RP), intent(out) :: qsat
839  real(RP), intent(in) :: temp
840  real(RP), intent(in) :: dens
841 
842  real(RP) :: psat
843  !---------------------------------------------------------------------------
844 
845  call atmos_saturation_psat_all( psat, temp )
846 
847  qsat = psat / ( dens * rvap * temp )
848 
849  return
850  end subroutine atmos_saturation_dens2qsat_all_0d
851 
852  !-----------------------------------------------------------------------------
854  subroutine atmos_saturation_dens2qsat_all_1d( &
855  qsat, &
856  temp, &
857  dens )
858  implicit none
859 
860  real(RP), intent(out) :: qsat(KA)
861  real(RP), intent(in) :: temp(KA)
862  real(RP), intent(in) :: dens(KA)
863 
864  real(RP) :: psat(KA)
865 
866  integer :: k
867  !---------------------------------------------------------------------------
868 
869  call atmos_saturation_psat_all( psat(:), temp(:) )
870 
871  do k = ks, ke
872  qsat(k) = psat(k) / ( dens(k) * rvap * temp(k) )
873  enddo
874 
875  return
876  end subroutine atmos_saturation_dens2qsat_all_1d
877 
878  !-----------------------------------------------------------------------------
880  subroutine atmos_saturation_dens2qsat_all_3d( &
881  qsat, &
882  temp, &
883  dens )
884  implicit none
885 
886  real(RP), intent(out) :: qsat(KA,IA,JA)
887  real(RP), intent(in) :: temp(KA,IA,JA)
888  real(RP), intent(in) :: dens(KA,IA,JA)
889 
890  real(RP) :: alpha, psatl, psati
891  real(RP) :: psat
892 
893  integer :: k, i, j
894  !---------------------------------------------------------------------------
895 
896  !$omp parallel do default(none) private(i,j,k,alpha,psatl,psati,psat) OMP_SCHEDULE_ collapse(2) &
897  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,temp,ATMOS_SATURATION_LLIMIT_TEMP,RTEM00,CPovR_liq,LovR_liq, &
898  !$omp qsat,dens,LovR_ice,CPovR_ice,ATMOS_SATURATION_ULIMIT_TEMP)
899  do j = jsb, jeb
900  do i = isb, ieb
901  do k = ks, ke
902 
903  alpha = ( temp(k,i,j) - atmos_saturation_llimit_temp ) &
904  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
905  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
906 
907  psatl = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_liq &
908  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
909 
910  psati = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_ice &
911  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
912 
913  psat = psatl * ( alpha ) &
914  + psati * ( 1.0_rp - alpha )
915 
916  qsat(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) )
917 
918  enddo
919  enddo
920  enddo
921  return
922  end subroutine atmos_saturation_dens2qsat_all_3d
923 
924  !-----------------------------------------------------------------------------
926  subroutine atmos_saturation_dens2qsat_liq_0d( &
927  qsat, &
928  temp, &
929  dens )
930  implicit none
931 
932  real(RP), intent(out) :: qsat
933  real(RP), intent(in) :: temp
934  real(RP), intent(in) :: dens
935 
936  real(RP) :: psat
937  !---------------------------------------------------------------------------
938 
939  call atmos_saturation_psat_liq( psat, temp )
940 
941  qsat = psat / ( dens * rvap * temp )
942 
943  return
944  end subroutine atmos_saturation_dens2qsat_liq_0d
945 
946  !-----------------------------------------------------------------------------
948  subroutine atmos_saturation_dens2qsat_liq_1d( &
949  qsat, &
950  temp, &
951  dens )
952  implicit none
953 
954  real(RP), intent(out) :: qsat(KA)
955  real(RP), intent(in) :: temp(KA)
956  real(RP), intent(in) :: dens(KA)
957 
958  real(RP) :: psat(KA)
959 
960  integer :: k
961  !---------------------------------------------------------------------------
962 
963  call atmos_saturation_psat_liq( psat(:), temp(:) )
964 
965  do k = ks, ke
966  qsat(k) = psat(k) / ( dens(k) * rvap * temp(k) )
967  enddo
968 
969  return
970  end subroutine atmos_saturation_dens2qsat_liq_1d
971 
972  !-----------------------------------------------------------------------------
974  subroutine atmos_saturation_dens2qsat_liq_3d( &
975  qsat, &
976  temp, &
977  dens )
978  implicit none
979 
980  real(RP), intent(out) :: qsat(KA,IA,JA)
981  real(RP), intent(in) :: temp(KA,IA,JA)
982  real(RP), intent(in) :: dens(KA,IA,JA)
983 
984  real(RP) :: psat
985 
986  integer :: k, i, j
987  !---------------------------------------------------------------------------
988 
989  !$omp parallel do default(none) private(i,j,k,psat) OMP_SCHEDULE_ collapse(2) &
990  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,temp,RTEM00,CPovR_liq,LovR_liq,qsat,dens)
991  do j = jsb, jeb
992  do i = isb, ieb
993  do k = ks, ke
994  psat = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_liq &
995  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
996 
997  qsat(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) )
998  enddo
999  enddo
1000  enddo
1001 
1002  return
1003  end subroutine atmos_saturation_dens2qsat_liq_3d
1004 
1005  !-----------------------------------------------------------------------------
1007  subroutine atmos_saturation_dens2qsat_ice_0d( &
1008  qsat, &
1009  temp, &
1010  dens )
1011  implicit none
1012 
1013  real(RP), intent(out) :: qsat
1014  real(RP), intent(in) :: temp
1015  real(RP), intent(in) :: dens
1016 
1017  real(RP) :: psat
1018  !---------------------------------------------------------------------------
1019 
1020  call atmos_saturation_psat_ice( psat, temp )
1021 
1022  qsat = psat / ( dens * rvap * temp )
1023 
1024  return
1025  end subroutine atmos_saturation_dens2qsat_ice_0d
1026 
1027  !-----------------------------------------------------------------------------
1029  subroutine atmos_saturation_dens2qsat_ice_1d( &
1030  qsat, &
1031  temp, &
1032  dens )
1033  implicit none
1034 
1035  real(RP), intent(out) :: qsat(KA)
1036  real(RP), intent(in) :: temp(KA)
1037  real(RP), intent(in) :: dens(KA)
1038 
1039  real(RP) :: psat(KA)
1040 
1041  integer :: k
1042  !---------------------------------------------------------------------------
1043 
1044  call atmos_saturation_psat_ice( psat(:), temp(:) )
1045 
1046  do k = ks, ke
1047  qsat(k) = psat(k) / ( dens(k) * rvap * temp(k) )
1048  enddo
1049 
1050  return
1051  end subroutine atmos_saturation_dens2qsat_ice_1d
1052 
1053  !-----------------------------------------------------------------------------
1055  subroutine atmos_saturation_dens2qsat_ice_3d( &
1056  qsat, &
1057  temp, &
1058  dens )
1059  implicit none
1060 
1061  real(RP), intent(out) :: qsat(KA,IA,JA)
1062  real(RP), intent(in) :: temp(KA,IA,JA)
1063  real(RP), intent(in) :: dens(KA,IA,JA)
1064 
1065  real(RP) :: psat
1066 
1067  integer :: k, i, j
1068  !---------------------------------------------------------------------------
1069 
1070  !$omp parallel do default(none) private(i,j,k,psat) OMP_SCHEDULE_ collapse(2) &
1071  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,temp,RTEM00,CPovR_ice,LovR_ice,qsat,dens)
1072  do j = jsb, jeb
1073  do i = isb, ieb
1074  do k = ks, ke
1075  psat = psat0 * ( temp(k,i,j) * rtem00 )**cpovr_ice &
1076  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp(k,i,j) ) )
1077 
1078  qsat(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) )
1079  enddo
1080  enddo
1081  enddo
1082 
1083  return
1084  end subroutine atmos_saturation_dens2qsat_ice_3d
1085 
1086  !-----------------------------------------------------------------------------
1088  subroutine atmos_saturation_dalphadt_0d( &
1089  dalpha_dT, &
1090  temp )
1091  implicit none
1092 
1093  real(RP), intent(out) :: dalpha_dT
1094  real(RP), intent(in) :: temp
1095 
1096  real(RP) :: lim1, lim2
1097  !---------------------------------------------------------------------------
1098 
1099  ! if Tup < temp, dalpha/dT = 0 (no slope)
1100  lim1 = 0.5_rp + sign( 0.5_rp, atmos_saturation_ulimit_temp - temp )
1101  ! if Tdn > temp, dalpha/dT = 0 (no slope)
1102  lim2 = 0.5_rp + sign( 0.5_rp, temp - atmos_saturation_llimit_temp )
1103 
1104  dalpha_dt = dalphadt_const * lim1 * lim2
1105 
1106  return
1107  end subroutine atmos_saturation_dalphadt_0d
1108 
1109  !-----------------------------------------------------------------------------
1111  subroutine atmos_saturation_dalphadt_1d( &
1112  dalpha_dT, &
1113  temp )
1114  implicit none
1115 
1116  real(RP), intent(out) :: dalpha_dT(KA)
1117  real(RP), intent(in) :: temp (KA)
1118 
1119  real(RP) :: lim1, lim2
1120 
1121  integer :: k
1122  !---------------------------------------------------------------------------
1123 
1124  do k = ks, ke
1125 
1126  ! if Tup < temp(k), dalpha/dT = 0 (no slope)
1127  lim1 = 0.5_rp + sign( 0.5_rp, atmos_saturation_ulimit_temp - temp(k) )
1128  ! if Tdn > temp(k), dalpha/dT = 0 (no slope)
1129  lim2 = 0.5_rp + sign( 0.5_rp, temp(k) - atmos_saturation_llimit_temp )
1130 
1131  dalpha_dt(k) = dalphadt_const * lim1 * lim2
1132 
1133  enddo
1134 
1135  return
1136  end subroutine atmos_saturation_dalphadt_1d
1137 
1138  !-----------------------------------------------------------------------------
1140  subroutine atmos_saturation_dalphadt_3d( &
1141  dalpha_dT, &
1142  temp )
1143  implicit none
1144 
1145  real(RP), intent(out) :: dalpha_dT(KA,IA,JA)
1146  real(RP), intent(in) :: temp (KA,IA,JA)
1147 
1148  real(RP) :: lim1, lim2
1149 
1150  integer :: k, i, j
1151  !---------------------------------------------------------------------------
1152 
1153  !$omp parallel do private(i,j,k,lim1,lim2) OMP_SCHEDULE_ collapse(2)
1154  do j = jsb, jeb
1155  do i = isb, ieb
1156  do k = ks, ke
1157 
1158  ! if Tup < temp(k,i,j), dalpha/dT = 0 (no slope)
1159  lim1 = 0.5_rp + sign( 0.5_rp, atmos_saturation_ulimit_temp - temp(k,i,j) )
1160  ! if Tdn > temp(k,i,j), dalpha/dT = 0 (no slope)
1161  lim2 = 0.5_rp + sign( 0.5_rp, temp(k,i,j) - atmos_saturation_llimit_temp )
1162 
1163  dalpha_dt(k,i,j) = dalphadt_const * lim1 * lim2
1164 
1165  enddo
1166  enddo
1167  enddo
1168 
1169  return
1170  end subroutine atmos_saturation_dalphadt_3d
1171 
1172  !-----------------------------------------------------------------------------
1173  ! (d qsw/d T)_{rho}: partial difference of qsat_water
1174  subroutine atmos_saturation_dqsw_dtem_rho( &
1175  dqsdtem, &
1176  temp, &
1177  dens )
1179  hydrometeor_lhv => atmos_hydrometeor_lhv
1180  implicit none
1181 
1182  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1183  real(RP), intent(in) :: temp (ka,ia,ja)
1184  real(RP), intent(in) :: dens (ka,ia,ja)
1185 
1186  real(RP) :: lhv(ka,ia,ja) ! latent heat of vaporization [J/kg]
1187  real(RP) :: psat
1188  real(RP) :: tem
1189 
1190  integer :: k, i, j
1191  !---------------------------------------------------------------------------
1192 
1193  call hydrometeor_lhv( lhv(:,:,:), temp(:,:,:) )
1194 
1195  !$omp parallel do private(i,j,k,TEM,psat) OMP_SCHEDULE_ collapse(2)
1196  do j = jsb, jeb
1197  do i = isb, ieb
1198  do k = ks, ke
1199  tem = max( temp(k,i,j), tem_min )
1200 
1201  psat = psat0 * ( tem * rtem00 )**cpovr_liq &
1202  * exp( lovr_liq * ( rtem00 - 1.0_rp/tem ) )
1203 
1204  dqsdtem(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) * temp(k,i,j) ) &
1205  * ( lhv(k,i,j) / ( rvap * temp(k,i,j) ) - 1.0_rp )
1206  enddo
1207  enddo
1208  enddo
1209 
1210  return
1211  end subroutine atmos_saturation_dqsw_dtem_rho
1212 
1213  !-----------------------------------------------------------------------------
1214  ! (d qsi/d T)_{rho}: partial difference of qsat_ice
1215  subroutine atmos_saturation_dqsi_dtem_rho( &
1216  dqsdtem, &
1217  temp, &
1218  dens )
1220  hydrometeor_lhs => atmos_hydrometeor_lhs
1221  implicit none
1222 
1223  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1224  real(RP), intent(in) :: temp (ka,ia,ja)
1225  real(RP), intent(in) :: dens (ka,ia,ja)
1226 
1227  real(RP) :: lhs(ka,ia,ja) ! latent heat of sublimation [J/kg]
1228  real(RP) :: psat
1229  real(RP) :: tem
1230 
1231  integer :: k, i, j
1232  !---------------------------------------------------------------------------
1233 
1234  call hydrometeor_lhs( lhs(:,:,:), temp(:,:,:) )
1235 
1236  !$omp parallel do private(i,j,k,TEM,psat) OMP_SCHEDULE_ collapse(2)
1237  do j = jsb, jeb
1238  do i = isb, ieb
1239  do k = ks, ke
1240  tem = max( temp(k,i,j), tem_min )
1241 
1242  psat = psat0 * ( tem * rtem00 )**cpovr_ice &
1243  * exp( lovr_ice * ( rtem00 - 1.0_rp/tem ) )
1244 
1245  dqsdtem(k,i,j) = psat / ( dens(k,i,j) * rvap * temp(k,i,j) * temp(k,i,j) ) &
1246  * ( lhs(k,i,j) / ( rvap * temp(k,i,j) ) - 1.0_rp )
1247  enddo
1248  enddo
1249  enddo
1250 
1251  return
1252  end subroutine atmos_saturation_dqsi_dtem_rho
1253 
1254  !-----------------------------------------------------------------------------
1255  ! (d qs/d T)_{p} and (d qs/d p)_{T}
1256  subroutine atmos_saturation_dqsw_dtem_dpre( &
1257  dqsdtem, &
1258  dqsdpre, &
1259  temp, &
1260  pres )
1262  hydrometeor_lhv => atmos_hydrometeor_lhv
1263  implicit none
1264 
1265  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1266  real(RP), intent(out) :: dqsdpre(ka,ia,ja)
1267  real(RP), intent(in) :: temp (ka,ia,ja)
1268  real(RP), intent(in) :: pres (ka,ia,ja)
1269 
1270  real(RP) :: lhv(ka,ia,ja) ! latent heat of vaporization [J/kg]
1271  real(RP) :: psat
1272  real(RP) :: tem
1273  real(RP) :: den1, den2 ! denominator
1274 
1275  integer :: k, i, j
1276  !---------------------------------------------------------------------------
1277 
1278  call hydrometeor_lhv( lhv(:,:,:), temp(:,:,:) )
1279 
1280  !$omp parallel do private(i,j,k,TEM,psat,den1,den2) OMP_SCHEDULE_ collapse(2)
1281  do j = jsb, jeb
1282  do i = isb, ieb
1283  do k = ks, ke
1284  tem = max( temp(k,i,j), tem_min )
1285 
1286  psat = psat0 * ( tem * rtem00 )**cpovr_liq &
1287  * exp( lovr_liq * ( rtem00 - 1.0_rp/tem ) )
1288 
1289  den1 = ( pres(k,i,j) - (1.0_rp-epsvap) * psat ) &
1290  * ( pres(k,i,j) - (1.0_rp-epsvap) * psat )
1291  den2 = den1 * rvap * temp(k,i,j) * temp(k,i,j)
1292 
1293  dqsdpre(k,i,j) = - epsvap * psat / den1
1294  dqsdtem(k,i,j) = epsvap * psat / den2 * lhv(k,i,j) * pres(k,i,j)
1295  enddo
1296  enddo
1297  enddo
1298 
1299  return
1300  end subroutine atmos_saturation_dqsw_dtem_dpre
1301 
1302  !-----------------------------------------------------------------------------
1303  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
1304  !-----------------------------------------------------------------------------
1305  subroutine atmos_saturation_dqsi_dtem_dpre( &
1306  dqsdtem, &
1307  dqsdpre, &
1308  temp, &
1309  pres )
1311  hydrometeor_lhs => atmos_hydrometeor_lhs
1312  implicit none
1313 
1314  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1315  real(RP), intent(out) :: dqsdpre(ka,ia,ja)
1316  real(RP), intent(in) :: temp (ka,ia,ja)
1317  real(RP), intent(in) :: pres (ka,ia,ja)
1318 
1319  real(RP) :: lhs(ka,ia,ja) ! latent heat of sublimation [J/kg]
1320  real(RP) :: psat
1321  real(RP) :: tem
1322  real(RP) :: den1, den2 ! denominator
1323 
1324  integer :: k, i, j
1325  !---------------------------------------------------------------------------
1326 
1327  call hydrometeor_lhs( lhs(:,:,:), temp(:,:,:) )
1328 
1329  !$omp parallel do private(i,j,k,TEM,psat,den1,den2) OMP_SCHEDULE_ collapse(2)
1330  do j = jsb, jeb
1331  do i = isb, ieb
1332  do k = ks, ke
1333  tem = max( temp(k,i,j), tem_min )
1334 
1335  psat = psat0 * ( tem * rtem00 )**cpovr_ice &
1336  * exp( lovr_ice * ( rtem00 - 1.0_rp/tem ) )
1337 
1338  den1 = ( pres(k,i,j) - (1.0_rp-epsvap) * psat ) &
1339  * ( pres(k,i,j) - (1.0_rp-epsvap) * psat )
1340  den2 = den1 * rvap * temp(k,i,j) * temp(k,i,j)
1341 
1342  dqsdpre(k,i,j) = - epsvap * psat / den1
1343  dqsdtem(k,i,j) = epsvap * psat / den2 * lhs(k,i,j) * pres(k,i,j)
1344  enddo
1345  enddo
1346  enddo
1347 
1348  return
1349  end subroutine atmos_saturation_dqsi_dtem_dpre
1350 
1351 end module scale_atmos_saturation
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:83
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
subroutine, public atmos_saturation_setup
Setup.
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:69
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
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), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:67
subroutine, public atmos_saturation_dqsw_dtem_rho(dqsdtem, temp, dens)
subroutine, public atmos_saturation_dqsw_dtem_dpre(dqsdtem, dqsdpre, temp, pres)
integer, public ieb
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:79
module grid index
subroutine, public atmos_saturation_dqsi_dtem_rho(dqsdtem, temp, dens)
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
integer, public ka
of whole cells: z, local, with HALO
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
subroutine atmos_saturation_alpha_0d(alpha, temp)
calc liquid/ice separation factor (0D)
module PROCESS
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
Definition: scale_const.F90:80
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:78
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:71
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
subroutine, public atmos_saturation_dqsi_dtem_dpre(dqsdtem, dqsdpre, temp, pres)
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
module PRECISION
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
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
integer, public ja
of whole cells: y, local, with HALO