SCALE-RM
scale_atmos_saturation.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
20 
21  use scale_const, only: &
22  rvap => const_rvap, &
23  psat0 => const_psat0, &
24  epsvap => const_epsvap, &
25  tem00 => const_tem00
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: atmos_saturation_setup
34 
35  public :: atmos_saturation_alpha
36 
37  public :: atmos_saturation_psat_all
38  public :: atmos_saturation_psat_liq
39  public :: atmos_saturation_psat_ice
40 
41  public :: atmos_saturation_psat2qsat_pres
42  public :: atmos_saturation_psat2qsat_dens
43 
44  public :: atmos_saturation_pres2qsat_all
45  public :: atmos_saturation_pres2qsat_liq
46  public :: atmos_saturation_pres2qsat_ice
47 
48  public :: atmos_saturation_dens2qsat_all
49  public :: atmos_saturation_dens2qsat_liq
50  public :: atmos_saturation_dens2qsat_ice
51 
52  public :: atmos_saturation_dalphadt
53 
54  public :: atmos_saturation_dqs_dtem_dens_liq
55  public :: atmos_saturation_dqs_dtem_dens_ice
56  public :: atmos_saturation_dqs_dtem_dpre_liq
57  public :: atmos_saturation_dqs_dtem_dpre_ice
58 
59  public :: atmos_saturation_tdew_liq
60 
61  public :: atmos_saturation_pote
62 
63  public :: atmos_saturation_moist_conversion_dens_liq
64  public :: atmos_saturation_moist_conversion_dens_all
65  public :: atmos_saturation_moist_conversion_pres_liq
66 
67  interface atmos_saturation_alpha
68  module procedure atmos_saturation_alpha_0d
69  module procedure atmos_saturation_alpha_1d
70  module procedure atmos_saturation_alpha_3d
71  end interface atmos_saturation_alpha
72 
73  interface atmos_saturation_psat_all
74  module procedure atmos_saturation_psat_all_0d
75  module procedure atmos_saturation_psat_all_1d
76  module procedure atmos_saturation_psat_all_2d
77  module procedure atmos_saturation_psat_all_3d
78  end interface atmos_saturation_psat_all
79  interface atmos_saturation_psat_liq
80  module procedure atmos_saturation_psat_liq_0d
81  module procedure atmos_saturation_psat_liq_1d
82  module procedure atmos_saturation_psat_liq_2d
83  module procedure atmos_saturation_psat_liq_3d
84  end interface atmos_saturation_psat_liq
85  interface atmos_saturation_psat_ice
86  module procedure atmos_saturation_psat_ice_0d
87  module procedure atmos_saturation_psat_ice_1d
88  module procedure atmos_saturation_psat_ice_2d
89  module procedure atmos_saturation_psat_ice_3d
90  end interface atmos_saturation_psat_ice
91 
92  interface atmos_saturation_psat2qsat_pres
93  module procedure atmos_saturation_psat2qsat_pres_0d
94  end interface atmos_saturation_psat2qsat_pres
95  interface atmos_saturation_psat2qsat_dens
96  module procedure atmos_saturation_psat2qsat_dens_0d
97  end interface atmos_saturation_psat2qsat_dens
98 
99  interface atmos_saturation_pres2qsat_all
100  module procedure atmos_saturation_pres2qsat_all_0d
101  module procedure atmos_saturation_pres2qsat_all_1d
102  module procedure atmos_saturation_pres2qsat_all_2d
103  module procedure atmos_saturation_pres2qsat_all_3d
104  end interface atmos_saturation_pres2qsat_all
105  interface atmos_saturation_pres2qsat_liq
106  module procedure atmos_saturation_pres2qsat_liq_0d
107  module procedure atmos_saturation_pres2qsat_liq_1d
108  module procedure atmos_saturation_pres2qsat_liq_3d
109  end interface atmos_saturation_pres2qsat_liq
110  interface atmos_saturation_pres2qsat_ice
111  module procedure atmos_saturation_pres2qsat_ice_0d
112  module procedure atmos_saturation_pres2qsat_ice_1d
113  module procedure atmos_saturation_pres2qsat_ice_3d
114  end interface atmos_saturation_pres2qsat_ice
115 
116  interface atmos_saturation_dens2qsat_all
117  module procedure atmos_saturation_dens2qsat_all_0d
118  module procedure atmos_saturation_dens2qsat_all_1d
119  module procedure atmos_saturation_dens2qsat_all_3d
120  end interface atmos_saturation_dens2qsat_all
121  interface atmos_saturation_dens2qsat_liq
122  module procedure atmos_saturation_dens2qsat_liq_0d
123  module procedure atmos_saturation_dens2qsat_liq_1d
124  module procedure atmos_saturation_dens2qsat_liq_3d
125  end interface atmos_saturation_dens2qsat_liq
126  interface atmos_saturation_dens2qsat_ice
127  module procedure atmos_saturation_dens2qsat_ice_0d
128  module procedure atmos_saturation_dens2qsat_ice_1d
129  module procedure atmos_saturation_dens2qsat_ice_3d
130  end interface atmos_saturation_dens2qsat_ice
131 
132  interface atmos_saturation_dalphadt
133  module procedure atmos_saturation_dalphadt_0d
134  module procedure atmos_saturation_dalphadt_1d
135  module procedure atmos_saturation_dalphadt_3d
136  end interface atmos_saturation_dalphadt
137 
138  interface atmos_saturation_dqs_dtem_dens_liq
139  module procedure atmos_saturation_dqs_dtem_dens_liq_0d
141  end interface atmos_saturation_dqs_dtem_dens_liq
142  interface atmos_saturation_dqs_dtem_dens_ice
143  module procedure atmos_saturation_dqs_dtem_dens_ice_0d
145  end interface atmos_saturation_dqs_dtem_dens_ice
146  interface atmos_saturation_dqs_dtem_dpre_liq
147  module procedure atmos_saturation_dqs_dtem_dpre_liq_0d
149  end interface atmos_saturation_dqs_dtem_dpre_liq
150  interface atmos_saturation_dqs_dtem_dpre_ice
151  module procedure atmos_saturation_dqs_dtem_dpre_ice_0d
153  end interface atmos_saturation_dqs_dtem_dpre_ice
154 
155  interface atmos_saturation_tdew_liq
156  module procedure atmos_saturation_tdew_liq_0d
157  module procedure atmos_saturation_tdew_liq_1d
158  module procedure atmos_saturation_tdew_liq_3d
159  end interface atmos_saturation_tdew_liq
160 
161  interface atmos_saturation_pote
162  module procedure atmos_saturation_pote_0d
163  module procedure atmos_saturation_pote_1d
164  module procedure atmos_saturation_pote_3d
165  end interface atmos_saturation_pote
166 
167  interface atmos_saturation_moist_conversion_dens_liq
168  module procedure atmos_saturation_moist_conversion_dens_liq_0d
169  end interface atmos_saturation_moist_conversion_dens_liq
170  interface atmos_saturation_moist_conversion_dens_all
172  end interface atmos_saturation_moist_conversion_dens_all
173  interface atmos_saturation_moist_conversion_pres_liq
175  end interface atmos_saturation_moist_conversion_pres_liq
176 
177  !-----------------------------------------------------------------------------
178  !
179  !++ Public parameters & variables
180  !
181  !-----------------------------------------------------------------------------
182  !
183  !++ Private procedure
184  !
185  !-----------------------------------------------------------------------------
186  !
187  !++ Private parameters & variables
188  !
189  real(RP), private, parameter :: tem_min = 10.0_rp
190 
191  real(RP), private :: atmos_saturation_ulimit_temp = 273.15_rp
192  real(RP), private :: atmos_saturation_llimit_temp = 233.15_rp
193 
194  real(RP), private :: rtem00
195  real(RP), private :: dalphadt_const
196  real(RP), private :: psat_min_liq
197  real(RP), private :: psat_min_ice
198 
199  real(RP), private :: cpovr_liq
200  real(RP), private :: cpovr_ice
201  real(RP), private :: cvovr_liq
202  real(RP), private :: cvovr_ice
203  real(RP), private :: lovr_liq
204  real(RP), private :: lovr_ice
205 
206  !-----------------------------------------------------------------------------
207 contains
208  !-----------------------------------------------------------------------------
210  subroutine atmos_saturation_setup
211  use scale_prc, only: &
212  prc_abort
213  use scale_const, only: &
214  cpvap => const_cpvap, &
215  cvvap => const_cvvap, &
216  cl => const_cl, &
217  ci => const_ci, &
218  lhv00 => const_lhv00, &
219  lhs00 => const_lhs00, &
220  lhv0 => const_lhv0, &
221  lhs0 => const_lhs0, &
223  use scale_atmos_hydrometeor, only: &
225  implicit none
226 
227  namelist / param_atmos_saturation / &
228  atmos_saturation_ulimit_temp, &
229  atmos_saturation_llimit_temp
230 
231  integer :: ierr
232  !---------------------------------------------------------------------------
233 
234  log_newline
235  log_info("ATMOS_SATURATION_setup",*) 'Setup'
236 
237  !--- read namelist
238  rewind(io_fid_conf)
239  read(io_fid_conf,nml=param_atmos_saturation,iostat=ierr)
240  if( ierr < 0 ) then !--- missing
241  log_info("ATMOS_SATURATION_setup",*) 'Not found namelist. Default used.'
242  elseif( ierr > 0 ) then !--- fatal error
243  log_error("ATMOS_SATURATION_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SATURATION. Check!'
244  call prc_abort
245  endif
246  log_nml(param_atmos_saturation)
247 
248  rtem00 = 1.0_rp / tem00
249 
250  if ( const_thermodyn_type == 'EXACT' ) then
251 
252  cpovr_liq = ( cpvap - cl ) / rvap
253  cpovr_ice = ( cpvap - ci ) / rvap
254  cvovr_liq = ( cvvap - cl ) / rvap
255  cvovr_ice = ( cvvap - ci ) / rvap
256 
257  lovr_liq = lhv00 / rvap
258  lovr_ice = lhs00 / rvap
259 
260  elseif( const_thermodyn_type == 'SIMPLE' ) then
261 
262  cpovr_liq = 0.0_rp
263  cpovr_ice = 0.0_rp
264  cvovr_liq = 0.0_rp
265  cvovr_ice = 0.0_rp
266 
267  lovr_liq = lhv0 / rvap
268  lovr_ice = lhs0 / rvap
269 
270  endif
271 
272  dalphadt_const = 1.0_rp / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
273 
274  log_newline
275  log_info("ATMOS_SATURATION_setup",'(1x,A,F7.2,A,F7.2)') 'Temperature range for liquid/ice mixture : ', &
276  atmos_saturation_llimit_temp, ' - ', &
277  atmos_saturation_ulimit_temp
278 
280 
281  call atmos_saturation_psat_liq( tem_min, psat_min_liq ) ! [IN], [OUT]
282  call atmos_saturation_psat_ice( tem_min, psat_min_ice ) ! [IN], [OUT]
283 
284  return
285  end subroutine atmos_saturation_setup
286 
287  !-----------------------------------------------------------------------------
289  subroutine atmos_saturation_alpha_0d( &
290  temp, &
291  alpha )
292  implicit none
293 
294  real(RP), intent(in) :: temp
295  real(RP), intent(out) :: alpha
296  !---------------------------------------------------------------------------
297 
298  alpha = ( temp - atmos_saturation_llimit_temp ) &
299  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
300 
301  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
302 
303  return
304  end subroutine atmos_saturation_alpha_0d
305 
306  !-----------------------------------------------------------------------------
308 !OCL SERIAL
309  subroutine atmos_saturation_alpha_1d( &
310  KA, KS, KE, &
311  temp, &
312  alpha )
313  implicit none
314  integer, intent(in) :: KA, KS, KE
315 
316  real(RP), intent(in) :: temp (ka)
317 
318  real(RP), intent(out) :: alpha(ka)
319 
320  integer :: k
321  !---------------------------------------------------------------------------
322 
323  do k = ks, ke
324  call atmos_saturation_alpha_0d( temp(k), & ! [IN]
325  alpha(k) ) ! [OUT]
326  enddo
327 
328  return
329  end subroutine atmos_saturation_alpha_1d
330 
331  !-----------------------------------------------------------------------------
333  subroutine atmos_saturation_alpha_3d( &
334  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
335  temp, &
336  alpha )
337  implicit none
338  integer, intent(in) :: KA, KS, KE
339  integer, intent(in) :: IA, IS, IE
340  integer, intent(in) :: JA, JS, JE
341 
342  real(RP), intent(in) :: temp (ka,ia,ja)
343  real(RP), intent(out) :: alpha(ka,ia,ja)
344 
345  integer :: k, i, j
346  !---------------------------------------------------------------------------
347 
348  !$omp parallel do OMP_SCHEDULE_ collapse(2)
349  do j = js, je
350  do i = is, ie
351  do k = ks, ke
352  call atmos_saturation_alpha_0d( temp(k,i,j), & ! [IN]
353  alpha(k,i,j) ) ! [OUT]
354 
355  enddo
356  enddo
357  enddo
358 
359  return
360  end subroutine atmos_saturation_alpha_3d
361 
362  !-----------------------------------------------------------------------------
364  subroutine atmos_saturation_psat_all_0d( &
365  temp, &
366  psat )
367  implicit none
368 
369  real(RP), intent(in) :: temp
370  real(RP), intent(out) :: psat
371 
372  real(RP) :: alpha, psatl, psati
373  !---------------------------------------------------------------------------
374 
375  call atmos_saturation_alpha ( temp, alpha )
376  call atmos_saturation_psat_liq( temp, psatl )
377  call atmos_saturation_psat_ice( temp, psati )
378 
379  psat = psatl * ( alpha ) &
380  + psati * ( 1.0_rp - alpha )
381 
382  return
383  end subroutine atmos_saturation_psat_all_0d
384 
385  !-----------------------------------------------------------------------------
387 !OCL SERIAL
388  subroutine atmos_saturation_psat_all_1d( &
389  KA, KS, KE, &
390  temp, &
391  psat )
392  implicit none
393  integer, intent(in) :: KA, KS, KE
394 
395  real(RP), intent(in) :: temp(ka)
396  real(RP), intent(out) :: psat(ka)
397 
398  integer :: k
399  !---------------------------------------------------------------------------
400 
401  do k = ks, ke
402  call atmos_saturation_psat_all_0d( temp(k), & ! [IN]
403  psat(k) ) ! [OUT]
404  enddo
405 
406  return
407  end subroutine atmos_saturation_psat_all_1d
408 
409  !-----------------------------------------------------------------------------
411  subroutine atmos_saturation_psat_all_2d( &
412  IA, IS, IE, &
413  JA, JS, JE, &
414  temp, &
415  psat )
416  implicit none
417  integer, intent(in) :: IA, IS, IE
418  integer, intent(in) :: JA, JS, JE
419 
420  real(RP), intent(in) :: temp(ia,ja)
421 
422  real(RP), intent(out) :: psat(ia,ja)
423 
424  integer :: i, j
425  !---------------------------------------------------------------------------
426 
427  !$omp parallel do OMP_SCHEDULE_ collapse(2)
428  do j = js, je
429  do i = is, ie
430  call atmos_saturation_psat_all_0d( temp(i,j), & ! [IN]
431  psat(i,j) ) ! [OUT]
432  enddo
433  enddo
434 
435  return
436  end subroutine atmos_saturation_psat_all_2d
437 
438  !-----------------------------------------------------------------------------
440  subroutine atmos_saturation_psat_all_3d( &
441  KA, KS, KE, &
442  IA, IS, IE, &
443  JA, JS, JE, &
444  temp, &
445  psat )
446  implicit none
447 
448  integer, intent(in) :: KA, KS, KE
449  integer, intent(in) :: IA, IS, IE
450  integer, intent(in) :: JA, JS, JE
451 
452  real(RP), intent(in) :: temp(ka,ia,ja)
453  real(RP), intent(out) :: psat(ka,ia,ja)
454 
455  integer :: k, i, j
456  !---------------------------------------------------------------------------
457 
458  !$omp parallel do OMP_SCHEDULE_ collapse(2)
459  do j = js, je
460  do i = is, ie
461  do k = ks, ke
462  call atmos_saturation_psat_all_0d( temp(k,i,j), & ! [IN]
463  psat(k,i,j) ) ! [OUT]
464  enddo
465  enddo
466  enddo
467 
468  return
469  end subroutine atmos_saturation_psat_all_3d
470 
471  !-----------------------------------------------------------------------------
473  subroutine atmos_saturation_psat_liq_0d( &
474  temp, &
475  psat )
476  implicit none
477 
478  real(RP), intent(in) :: temp
479  real(RP), intent(out) :: psat
480  !---------------------------------------------------------------------------
481 
482  psat = psat0 * ( temp * rtem00 )**cpovr_liq &
483  * exp( lovr_liq * ( rtem00 - 1.0_rp/temp ) )
484 
485  return
486  end subroutine atmos_saturation_psat_liq_0d
487 
488  !-----------------------------------------------------------------------------
490 !OCL SERIAL
491  subroutine atmos_saturation_psat_liq_1d( &
492  KA, KS, KE, &
493  temp, &
494  psat )
495  implicit none
496  integer, intent(in) :: KA, KS, KE
497 
498  real(RP), intent(in) :: temp(ka)
499  real(RP), intent(out) :: psat(ka)
500 
501  integer :: k
502  !---------------------------------------------------------------------------
503 
504  do k = ks, ke
505  call atmos_saturation_psat_liq_0d( temp(k), &
506  psat(k) )
507  enddo
508 
509  return
510  end subroutine atmos_saturation_psat_liq_1d
511 
512  !-----------------------------------------------------------------------------
514  subroutine atmos_saturation_psat_liq_2d( &
515  IA, IS, IE, &
516  JA, JS, JE, &
517  temp, &
518  psat )
519  implicit none
520 
521  integer, intent(in) :: IA, IS, IE
522  integer, intent(in) :: JA, JS, JE
523 
524  real(RP), intent(in) :: temp(ia,ja)
525  real(RP), intent(out) :: psat(ia,ja)
526 
527  integer :: i, j
528  !---------------------------------------------------------------------------
529 
530  !$omp parallel do OMP_SCHEDULE_ collapse(2)
531  do j = js, je
532  do i = is, ie
533  call atmos_saturation_psat_liq_0d( temp(i,j), & ! [IN]
534  psat(i,j) ) ! [OUT]
535  enddo
536  enddo
537 
538  return
539  end subroutine atmos_saturation_psat_liq_2d
540 
541  !-----------------------------------------------------------------------------
543  subroutine atmos_saturation_psat_liq_3d( &
544  KA, KS, KE, &
545  IA, IS, IE, &
546  JA, JS, JE, &
547  temp, &
548  psat )
549  implicit none
550 
551  integer, intent(in) :: KA, KS, KE
552  integer, intent(in) :: IA, IS, IE
553  integer, intent(in) :: JA, JS, JE
554 
555  real(RP), intent(in) :: temp(ka,ia,ja)
556  real(RP), intent(out) :: psat(ka,ia,ja)
557 
558  integer :: k, i, j
559  !---------------------------------------------------------------------------
560 
561  !$omp parallel do OMP_SCHEDULE_ collapse(2)
562  do j = js, je
563  do i = is, ie
564  do k = ks, ke
565  call atmos_saturation_psat_liq_0d( temp(k,i,j), & ! [IN]
566  psat(k,i,j) ) ! [OUT]
567  enddo
568  enddo
569  enddo
570 
571  return
572  end subroutine atmos_saturation_psat_liq_3d
573 
574  !-----------------------------------------------------------------------------
576  subroutine atmos_saturation_psat_ice_0d( &
577  temp, &
578  psat )
579  implicit none
580 
581  real(RP), intent(in) :: temp
582 
583  real(RP), intent(out) :: psat
584  !---------------------------------------------------------------------------
585 
586  psat = psat0 * ( temp * rtem00 )**cpovr_ice &
587  * exp( lovr_ice * ( rtem00 - 1.0_rp/temp ) )
588 
589  return
590  end subroutine atmos_saturation_psat_ice_0d
591 
592  !-----------------------------------------------------------------------------
594 !OCL SERIAL
595  subroutine atmos_saturation_psat_ice_1d( &
596  KA, KS, KE, &
597  temp, &
598  psat )
599  implicit none
600  integer, intent(in) :: KA, KS, KE
601 
602  real(RP), intent(in) :: temp(ka)
603  real(RP), intent(out) :: psat(ka)
604 
605  integer :: k
606  !---------------------------------------------------------------------------
607 
608  do k = ks, ke
609  call atmos_saturation_psat_ice_0d( temp(k), & ! [IN]
610  psat(k) ) ! [OUT]
611  enddo
612 
613  return
614  end subroutine atmos_saturation_psat_ice_1d
615 
616  !-----------------------------------------------------------------------------
618  subroutine atmos_saturation_psat_ice_2d( &
619  IA, IS, IE, &
620  JA, JS, JE, &
621  temp, &
622  psat )
623  implicit none
624 
625  integer, intent(in) :: IA, IS, IE
626  integer, intent(in) :: JA, JS, JE
627 
628  real(RP), intent(in) :: temp(ia,ja)
629  real(RP), intent(out) :: psat(ia,ja)
630 
631  integer :: i, j
632  !---------------------------------------------------------------------------
633 
634  !$omp parallel do OMP_SCHEDULE_ collapse(2)
635  do j = js, je
636  do i = is, ie
637  call atmos_saturation_psat_ice_0d( temp(i,j), & ! [IN]
638  psat(i,j) ) ! [OUT]
639  enddo
640  enddo
641 
642  return
643  end subroutine atmos_saturation_psat_ice_2d
644 
645  !-----------------------------------------------------------------------------
647  subroutine atmos_saturation_psat_ice_3d( &
648  KA, KS, KE, &
649  IA, IS, IE, &
650  JA, JS, JE, &
651  temp, &
652  psat )
653  implicit none
654 
655  integer, intent(in) :: KA, KS, KE
656  integer, intent(in) :: IA, IS, IE
657  integer, intent(in) :: JA, JS, JE
658 
659  real(RP), intent(in) :: temp(ka,ia,ja)
660  real(RP), intent(out) :: psat(ka,ia,ja)
661 
662  integer :: k, i, j
663  !---------------------------------------------------------------------------
664 
665  !$omp parallel do OMP_SCHEDULE_ collapse(2)
666  do j = js, je
667  do i = is, ie
668  do k = ks, ke
669  call atmos_saturation_psat_ice_0d( temp(k,i,j), & ! [IN]
670  psat(k,i,j) ) ! [OUT]
671  enddo
672  enddo
673  enddo
674 
675  return
676  end subroutine atmos_saturation_psat_ice_3d
677 
678 
679  !-----------------------------------------------------------------------------
681  subroutine atmos_saturation_psat2qsat_pres_0d( &
682  psat, pres, qdry, &
683  qsat )
684  implicit none
685 
686  real(RP), intent(in) :: psat
687  real(RP), intent(in) :: pres
688  real(RP), intent(in) :: qdry
689 
690  real(RP), intent(out) :: qsat
691 
692  ! ! qdry is assumed to be 1 - qsat
693  ! qsat = EPSvap * psat / ( pres - ( 1.0_RP-EPSvap ) * psat )
694 
695  qsat = epsvap * qdry * psat / ( pres - psat )
696 
697  return
698  end subroutine atmos_saturation_psat2qsat_pres_0d
699  !-----------------------------------------------------------------------------
701  subroutine atmos_saturation_pres2qsat_all_0d( &
702  temp, pres, qdry, &
703  qsat )
704  implicit none
705 
706  real(RP), intent(in) :: temp
707  real(RP), intent(in) :: pres
708  real(RP), intent(in) :: qdry
709 
710  real(RP), intent(out) :: qsat
711 
712  real(RP) :: psat
713  !---------------------------------------------------------------------------
714 
715  call atmos_saturation_psat_all_0d( temp, psat ) ! [IN], [OUT]
716  call atmos_saturation_psat2qsat_pres_0d( psat, pres, qdry, & ! [IN]
717  qsat ) ! [OUT]
718 
719  return
720  end subroutine atmos_saturation_pres2qsat_all_0d
721 
722  !-----------------------------------------------------------------------------
724 !OCL SERIAL
725  subroutine atmos_saturation_pres2qsat_all_1d( &
726  KA, KS, KE, &
727  temp, pres, qdry, &
728  qsat )
729  implicit none
730  integer, intent(in) :: KA, KS, KE
731 
732  real(RP), intent(in) :: temp(ka)
733  real(RP), intent(in) :: pres(ka)
734  real(RP), intent(in) :: qdry(ka)
735 
736  real(RP), intent(out) :: qsat(ka)
737 
738  integer :: k
739  !---------------------------------------------------------------------------
740 
741  do k = ks, ke
742  call atmos_saturation_pres2qsat_all_0d( temp(k), pres(k), qdry(k), & ! [IN]
743  qsat(k) ) ! [OUT]
744  enddo
745 
746  return
747  end subroutine atmos_saturation_pres2qsat_all_1d
748 
749  !-----------------------------------------------------------------------------
751  subroutine atmos_saturation_pres2qsat_all_2d( &
752  IA, IS, IE, JA, JS, JE, &
753  temp, pres, qdry, &
754  qsat )
755  implicit none
756  integer, intent(in) :: IA, IS, IE
757  integer, intent(in) :: JA, JS, JE
758 
759  real(RP), intent(in) :: temp(ia,ja)
760  real(RP), intent(in) :: pres(ia,ja)
761  real(RP), intent(in) :: qdry(ia,ja)
762 
763  real(RP), intent(out) :: qsat(ia,ja)
764 
765  integer :: i, j
766  !---------------------------------------------------------------------------
767 
768  !$omp parallel do OMP_SCHEDULE_ collapse(2)
769  do j = js, je
770  do i = is, ie
771  call atmos_saturation_pres2qsat_all_0d( temp(i,j), pres(i,j), qdry(i,j), &
772  qsat(i,j) )
773  enddo
774  enddo
775 
776  return
777  end subroutine atmos_saturation_pres2qsat_all_2d
778 
779  !-----------------------------------------------------------------------------
781  subroutine atmos_saturation_pres2qsat_all_3d( &
782  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
783  temp, pres, qdry, &
784  qsat )
785  implicit none
786  integer, intent(in) :: KA, KS, KE
787  integer, intent(in) :: IA, IS, IE
788  integer, intent(in) :: JA, JS, JE
789 
790  real(RP), intent(in) :: temp(ka,ia,ja)
791  real(RP), intent(in) :: pres(ka,ia,ja)
792  real(RP), intent(in) :: qdry(ka,ia,ja)
793 
794  real(RP), intent(out) :: qsat(ka,ia,ja)
795 
796  integer :: k, i, j
797  !---------------------------------------------------------------------------
798 
799  !$omp parallel do OMP_SCHEDULE_ collapse(2)
800  do j = js, je
801  do i = is, ie
802  do k = ks, ke
803  call atmos_saturation_pres2qsat_all_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
804  qsat(k,i,j) ) ! [OUT]
805  enddo
806  enddo
807  enddo
808 
809  return
810  end subroutine atmos_saturation_pres2qsat_all_3d
811 
812  !-----------------------------------------------------------------------------
814  subroutine atmos_saturation_pres2qsat_liq_0d( &
815  temp, pres, qdry, &
816  qsat )
817  implicit none
818 
819  real(RP), intent(in) :: temp
820  real(RP), intent(in) :: pres
821  real(RP), intent(in) :: qdry
822 
823  real(RP), intent(out) :: qsat
824 
825  real(RP) :: psat
826  !---------------------------------------------------------------------------
827 
828  call atmos_saturation_psat_liq( temp, psat ) ! [IN], [OUT]
829  call atmos_saturation_psat2qsat_pres( psat, pres, qdry, & ! [IN]
830  qsat ) ! [OUT]
831 
832  return
833  end subroutine atmos_saturation_pres2qsat_liq_0d
834 
835  !-----------------------------------------------------------------------------
837 !OCL SERIAL
838  subroutine atmos_saturation_pres2qsat_liq_1d( &
839  KA, KS, KE, &
840  temp, pres, qdry, &
841  qsat )
842  implicit none
843  integer, intent(in) :: KA, KS, KE
844 
845  real(RP), intent(in) :: temp(ka)
846  real(RP), intent(in) :: pres(ka)
847  real(RP), intent(in) :: qdry(ka)
848 
849  real(RP), intent(out) :: qsat(ka)
850 
851  integer :: k
852  !---------------------------------------------------------------------------
853 
854  do k = ks, ke
855  call atmos_saturation_pres2qsat_liq_0d( temp(k), pres(k), qdry(k), & ! [IN]
856  qsat(k) ) ! [OUT]
857  enddo
858 
859  return
860  end subroutine atmos_saturation_pres2qsat_liq_1d
861 
862  !-----------------------------------------------------------------------------
864  subroutine atmos_saturation_pres2qsat_liq_3d( &
865  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
866  temp, pres, qdry, &
867  qsat )
868  implicit none
869  integer, intent(in) :: KA, KS, KE
870  integer, intent(in) :: IA, IS, IE
871  integer, intent(in) :: JA, JS, JE
872 
873  real(RP), intent(in) :: temp(ka,ia,ja)
874  real(RP), intent(in) :: pres(ka,ia,ja)
875  real(RP), intent(in) :: qdry(ka,ia,ja)
876 
877  real(RP), intent(out) :: qsat(ka,ia,ja)
878 
879  integer :: k, i, j
880  !---------------------------------------------------------------------------
881 
882  !$omp parallel do OMP_SCHEDULE_ collapse(2)
883  do j = js, je
884  do i = is, ie
885  do k = ks, ke
886  call atmos_saturation_pres2qsat_liq_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
887  qsat(k,i,j) ) ! [OUT]
888  enddo
889  enddo
890  enddo
891 
892  return
893  end subroutine atmos_saturation_pres2qsat_liq_3d
894 
895  !-----------------------------------------------------------------------------
897  subroutine atmos_saturation_pres2qsat_ice_0d( &
898  temp, pres, qdry, &
899  qsat )
900  implicit none
901 
902  real(RP), intent(in) :: temp
903  real(RP), intent(in) :: pres
904  real(RP), intent(in) :: qdry
905 
906  real(RP), intent(out) :: qsat
907 
908  real(RP) :: psat
909  !---------------------------------------------------------------------------
910 
911  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
912  call atmos_saturation_psat2qsat_pres( psat, pres, qdry, & ! [IN]
913  qsat ) ! [OUT]
914 
915  return
916  end subroutine atmos_saturation_pres2qsat_ice_0d
917 
918  !-----------------------------------------------------------------------------
920 !OCL SERIAL
921  subroutine atmos_saturation_pres2qsat_ice_1d( &
922  KA, KS, KE, &
923  temp, pres, qdry, &
924  qsat )
925  implicit none
926  integer, intent(in) :: KA, KS, KE
927 
928  real(RP), intent(in) :: temp(ka)
929  real(RP), intent(in) :: pres(ka)
930  real(RP), intent(in) :: qdry(ka)
931 
932  real(RP), intent(out) :: qsat(ka)
933 
934  integer :: k
935  !---------------------------------------------------------------------------
936 
937  do k = ks, ke
938  call atmos_saturation_pres2qsat_ice_0d( temp(k), pres(k), qdry(k), & ! [IN]
939  qsat(k) ) ! [OUT]
940  enddo
941 
942  return
943  end subroutine atmos_saturation_pres2qsat_ice_1d
944 
945  !-----------------------------------------------------------------------------
947  subroutine atmos_saturation_pres2qsat_ice_3d( &
948  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
949  temp, pres, qdry, &
950  qsat )
951  implicit none
952  integer, intent(in) :: KA, KS, KE
953  integer, intent(in) :: IA, IS, IE
954  integer, intent(in) :: JA, JS, JE
955 
956  real(RP), intent(in) :: temp(ka,ia,ja)
957  real(RP), intent(in) :: pres(ka,ia,ja)
958  real(RP), intent(in) :: qdry(ka,ia,ja)
959 
960  real(RP), intent(out) :: qsat(ka,ia,ja)
961 
962  integer :: k, i, j
963  !---------------------------------------------------------------------------
964 
965  !$omp parallel do OMP_SCHEDULE_ collapse(2)
966  do j = js, je
967  do i = is, ie
968  do k = ks, ke
969  call atmos_saturation_pres2qsat_ice_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
970  qsat(k,i,j) ) ! [OUT]
971  enddo
972  enddo
973  enddo
974 
975  return
976  end subroutine atmos_saturation_pres2qsat_ice_3d
977 
978  !-----------------------------------------------------------------------------
980  subroutine atmos_saturation_psat2qsat_dens_0d( &
981  psat, temp, dens, &
982  qsat )
983  implicit none
984  real(RP), intent(in) :: psat
985  real(RP), intent(in) :: temp
986  real(RP), intent(in) :: dens
987 
988  real(RP), intent(out) :: qsat
989 
990  qsat = psat / ( dens * rvap * temp )
991 
992  return
993  end subroutine atmos_saturation_psat2qsat_dens_0d
994 
995  !-----------------------------------------------------------------------------
997  subroutine atmos_saturation_dens2qsat_all_0d( &
998  temp, &
999  dens, &
1000  qsat )
1001  implicit none
1002  real(RP), intent(in) :: temp
1003  real(RP), intent(in) :: dens
1004 
1005  real(RP), intent(out) :: qsat
1006 
1007  real(RP) :: psat
1008  !---------------------------------------------------------------------------
1009 
1010  call atmos_saturation_psat_all( temp, psat ) ! [IN], [OUT]
1011  call atmos_saturation_psat2qsat_dens_0d( psat, temp, dens, & ! [IN]
1012  qsat ) ! [OUT]
1013 
1014  return
1015  end subroutine atmos_saturation_dens2qsat_all_0d
1016 
1017  !-----------------------------------------------------------------------------
1019 !OCL SERIAL
1020  subroutine atmos_saturation_dens2qsat_all_1d( &
1021  KA, KS, KE, &
1022  temp, dens, &
1023  qsat )
1024  implicit none
1025  integer, intent(in) :: KA, KS, KE
1026 
1027  real(RP), intent(in) :: temp(ka)
1028  real(RP), intent(in) :: dens(ka)
1029 
1030  real(RP), intent(out) :: qsat(ka)
1031 
1032  integer :: k
1033  !---------------------------------------------------------------------------
1034 
1035  do k = ks, ke
1036  call atmos_saturation_dens2qsat_all_0d( temp(k), dens(k), & ! [IN]
1037  qsat(k) ) ! [OUT]
1038  enddo
1039 
1040  return
1041  end subroutine atmos_saturation_dens2qsat_all_1d
1042 
1043  !-----------------------------------------------------------------------------
1045  subroutine atmos_saturation_dens2qsat_all_3d( &
1046  KA, KS, KE, &
1047  IA, IS, IE, &
1048  JA, JS, JE, &
1049  temp, &
1050  dens, &
1051  qsat )
1052  implicit none
1053  integer, intent(in) :: KA, KS, KE
1054  integer, intent(in) :: IA, IS, IE
1055  integer, intent(in) :: JA, JS, JE
1056 
1057  real(RP), intent(in) :: temp(ka,ia,ja)
1058  real(RP), intent(in) :: dens(ka,ia,ja)
1059 
1060  real(RP), intent(out) :: qsat(ka,ia,ja)
1061 
1062  integer :: k, i, j
1063  !---------------------------------------------------------------------------
1064 
1065  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1066  do j = js, je
1067  do i = is, ie
1068  do k = ks, ke
1069  call atmos_saturation_dens2qsat_all_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1070  qsat(k,i,j) ) ! [OUT]
1071  enddo
1072  enddo
1073  enddo
1074 
1075  return
1076  end subroutine atmos_saturation_dens2qsat_all_3d
1077 
1078  !-----------------------------------------------------------------------------
1080  subroutine atmos_saturation_dens2qsat_liq_0d( &
1081  temp, &
1082  dens, &
1083  qsat )
1084  implicit none
1085  real(RP), intent(in) :: temp
1086  real(RP), intent(in) :: dens
1087 
1088  real(RP), intent(out) :: qsat
1089 
1090  real(RP) :: psat
1091  !---------------------------------------------------------------------------
1092 
1093  call atmos_saturation_psat_liq_0d( temp, psat ) ! [IN], [OUT]
1094  call atmos_saturation_psat2qsat_dens_0d( psat, temp, dens, & ! [IN]
1095  qsat ) ! [OUT]
1096 
1097  return
1098  end subroutine atmos_saturation_dens2qsat_liq_0d
1099 
1100  !-----------------------------------------------------------------------------
1102 !OCL SERIAL
1103  subroutine atmos_saturation_dens2qsat_liq_1d( &
1104  KA, KS, KE, &
1105  temp, dens, &
1106  qsat )
1107  implicit none
1108  integer, intent(in) :: KA, KS, KE
1109 
1110  real(RP), intent(in) :: temp(ka)
1111  real(RP), intent(in) :: dens(ka)
1112 
1113  real(RP), intent(out) :: qsat(ka)
1114 
1115  integer :: k
1116  !---------------------------------------------------------------------------
1117 
1118  do k = ks, ke
1119  call atmos_saturation_dens2qsat_liq_0d( temp(k), dens(k), & ! [IN]
1120  qsat(k) ) ! [OUT]
1121  enddo
1122 
1123  return
1124  end subroutine atmos_saturation_dens2qsat_liq_1d
1125 
1126  !-----------------------------------------------------------------------------
1128  subroutine atmos_saturation_dens2qsat_liq_3d( &
1129  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1130  temp, dens, &
1131  qsat )
1132  implicit none
1133  integer, intent(in) :: KA, KS, KE
1134  integer, intent(in) :: IA, IS, IE
1135  integer, intent(in) :: JA, JS, JE
1136 
1137  real(RP), intent(in) :: temp(ka,ia,ja)
1138  real(RP), intent(in) :: dens(ka,ia,ja)
1139 
1140  real(RP), intent(out) :: qsat(ka,ia,ja)
1141 
1142  integer :: k, i, j
1143  !---------------------------------------------------------------------------
1144 
1145  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1146  do j = js, je
1147  do i = is, ie
1148  do k = ks, ke
1149  call atmos_saturation_dens2qsat_liq_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1150  qsat(k,i,j) ) ! [OUT]
1151  enddo
1152  enddo
1153  enddo
1154 
1155  return
1156  end subroutine atmos_saturation_dens2qsat_liq_3d
1157 
1158  !-----------------------------------------------------------------------------
1160  subroutine atmos_saturation_dens2qsat_ice_0d( &
1161  temp, &
1162  dens, &
1163  qsat )
1164  implicit none
1165  real(RP), intent(in) :: temp
1166  real(RP), intent(in) :: dens
1167 
1168  real(RP), intent(out) :: qsat
1169 
1170  real(RP) :: psat
1171  !---------------------------------------------------------------------------
1172 
1173  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
1174  call atmos_saturation_psat2qsat_dens( psat, temp, dens, & ! [IN]
1175  qsat ) ! [OUT]
1176  return
1177  end subroutine atmos_saturation_dens2qsat_ice_0d
1178 
1179  !-----------------------------------------------------------------------------
1181 !OCL SERIAL
1182  subroutine atmos_saturation_dens2qsat_ice_1d( &
1183  KA, KS, KE, &
1184  temp, dens, &
1185  qsat )
1186  implicit none
1187  integer, intent(in) :: KA, KS, KE
1188 
1189  real(RP), intent(in) :: temp(ka)
1190  real(RP), intent(in) :: dens(ka)
1191 
1192  real(RP), intent(out) :: qsat(ka)
1193 
1194  integer :: k
1195  !---------------------------------------------------------------------------
1196 
1197  do k = ks, ke
1198  call atmos_saturation_dens2qsat_ice_0d( temp(k), dens(k), & ! [IN]
1199  qsat(k) ) ! [OUT]
1200  enddo
1201 
1202  return
1203  end subroutine atmos_saturation_dens2qsat_ice_1d
1204 
1205  !-----------------------------------------------------------------------------
1207  subroutine atmos_saturation_dens2qsat_ice_3d( &
1208  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1209  temp, dens, &
1210  qsat )
1211  implicit none
1212  integer, intent(in) :: KA, KS, KE
1213  integer, intent(in) :: IA, IS, IE
1214  integer, intent(in) :: JA, JS, JE
1215 
1216  real(RP), intent(in) :: temp(ka,ia,ja)
1217  real(RP), intent(in) :: dens(ka,ia,ja)
1218 
1219  real(RP), intent(out) :: qsat(ka,ia,ja)
1220 
1221  integer :: k, i, j
1222  !---------------------------------------------------------------------------
1223 
1224  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1225  do j = js, je
1226  do i = is, ie
1227  do k = ks, ke
1228  call atmos_saturation_dens2qsat_ice_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1229  qsat(k,i,j) ) ! [OUT]
1230  enddo
1231  enddo
1232  enddo
1233 
1234  return
1235  end subroutine atmos_saturation_dens2qsat_ice_3d
1236 
1237  !-----------------------------------------------------------------------------
1239  subroutine atmos_saturation_dalphadt_0d( &
1240  temp, &
1241  dalpha_dT )
1242  implicit none
1243  real(RP), intent(in) :: temp
1244 
1245  real(RP), intent(out) :: dalpha_dT
1246 
1247  real(RP) :: lim1, lim2
1248  !---------------------------------------------------------------------------
1249 
1250  ! if Tup < temp, dalpha/dT = 0 (no slope)
1251  lim1 = 0.5_rp + sign( 0.5_rp, atmos_saturation_ulimit_temp - temp )
1252  ! if Tdn > temp, dalpha/dT = 0 (no slope)
1253  lim2 = 0.5_rp + sign( 0.5_rp, temp - atmos_saturation_llimit_temp )
1254 
1255  dalpha_dt = dalphadt_const * lim1 * lim2
1256 
1257  return
1258  end subroutine atmos_saturation_dalphadt_0d
1259 
1260  !-----------------------------------------------------------------------------
1262 !OCL SERIAL
1263  subroutine atmos_saturation_dalphadt_1d( &
1264  KA, KS, KE, &
1265  temp, &
1266  dalpha_dT )
1267  implicit none
1268  integer, intent(in) :: KA, KS, KE
1269 
1270  real(RP), intent(in) :: temp (ka)
1271 
1272  real(RP), intent(out) :: dalpha_dT(ka)
1273 
1274  integer :: k
1275  !---------------------------------------------------------------------------
1276 
1277  do k = ks, ke
1278  call atmos_saturation_dalphadt_0d( temp(k), dalpha_dt(k) ) ! [IN], [OUT]
1279  enddo
1280 
1281  return
1282  end subroutine atmos_saturation_dalphadt_1d
1283 
1284  !-----------------------------------------------------------------------------
1286  subroutine atmos_saturation_dalphadt_3d( &
1287  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1288  temp, &
1289  dalpha_dT )
1290  implicit none
1291  integer, intent(in) :: KA, KS, KE
1292  integer, intent(in) :: IA, IS, IE
1293  integer, intent(in) :: JA, JS, JE
1294 
1295  real(RP), intent(in) :: temp (ka,ia,ja)
1296 
1297  real(RP), intent(out) :: dalpha_dT(ka,ia,ja)
1298 
1299  integer :: k, i, j
1300  !---------------------------------------------------------------------------
1301 
1302  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1303  do j = js, je
1304  do i = is, ie
1305  do k = ks, ke
1306  call atmos_saturation_dalphadt_0d( temp(k,i,j), dalpha_dt(k,i,j) ) ! [IN], [OUT]
1307  enddo
1308  enddo
1309  enddo
1310 
1311  return
1312  end subroutine atmos_saturation_dalphadt_3d
1313 
1314  !-----------------------------------------------------------------------------
1315  ! (d qsw/d T)_{rho}: partial difference of qsat_water at constant density
1316  subroutine atmos_saturation_dqs_dtem_dens_liq_0d( &
1317  temp, dens, &
1318  dqsdtem, qsat )
1319  use scale_atmos_hydrometeor, only: &
1320  hydrometeor_lhv => atmos_hydrometeor_lhv
1321  implicit none
1322 
1323  real(RP), intent(in) :: temp
1324  real(RP), intent(in) :: dens
1325 
1326  real(RP), intent(out) :: dqsdtem
1327  real(RP), intent(out), optional :: qsat
1328 
1329  real(RP) :: LHV ! latent heat of vaporization [J/kg]
1330  real(RP) :: psat
1331  !---------------------------------------------------------------------------
1332 
1333  call hydrometeor_lhv( temp, lhv )
1334  call atmos_saturation_psat_liq( temp, psat ) ! [IN], [OUT]
1335 
1336  dqsdtem = psat / ( dens* rvap * temp**2 ) &
1337  * ( lhv / ( rvap * temp ) - 1.0_rp )
1338 
1339  if ( present(qsat) ) &
1340  call atmos_saturation_psat2qsat_dens( psat, temp, dens, & ! [IN]
1341  qsat ) ! [OUT]
1342 
1343  return
1344  end subroutine atmos_saturation_dqs_dtem_dens_liq_0d
1345 
1346  !-----------------------------------------------------------------------------
1347  ! (d qsw/d T)_{rho}: partial difference of qsat_water at constant density
1349  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1350  temp, dens, &
1351  dqsdtem )
1352  implicit none
1353  integer, intent(in) :: KA, KS, KE
1354  integer, intent(in) :: IA, IS, IE
1355  integer, intent(in) :: JA, JS, JE
1356 
1357  real(RP), intent(in) :: temp (ka,ia,ja)
1358  real(RP), intent(in) :: dens (ka,ia,ja)
1359 
1360  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1361 
1362  integer :: k, i, j
1363  !---------------------------------------------------------------------------
1364 
1365  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1366  do j = js, je
1367  do i = is, ie
1368  do k = ks, ke
1369  call atmos_saturation_dqs_dtem_dens_liq_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1370  dqsdtem(k,i,j) ) ! [OUT]
1371  enddo
1372  enddo
1373  enddo
1374 
1375  return
1377 
1378  !-----------------------------------------------------------------------------
1379  ! (d qsi/d T)_{rho}: partial difference of qsat_ice at constant density
1380  subroutine atmos_saturation_dqs_dtem_dens_ice_0d( &
1381  temp, dens, &
1382  dqsdtem, qsat )
1383  use scale_atmos_hydrometeor, only: &
1384  hydrometeor_lhs => atmos_hydrometeor_lhs
1385  implicit none
1386 
1387  real(RP), intent(in) :: temp
1388  real(RP), intent(in) :: dens
1389 
1390  real(RP), intent(out) :: dqsdtem
1391  real(RP), intent(out), optional :: qsat
1392 
1393  real(RP) :: LHS ! latent heat of sublimation [J/kg]
1394  real(RP) :: psat
1395 
1396  !---------------------------------------------------------------------------
1397 
1398  call hydrometeor_lhs( temp, lhs ) ! [IN], [OUT]
1399  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
1400 
1401  dqsdtem = psat / ( dens * rvap * temp**2 ) &
1402  * ( lhs / ( rvap * temp ) - 1.0_rp )
1403 
1404  if ( present(qsat) ) &
1405  call atmos_saturation_psat2qsat_dens( psat, temp, dens, & ! [IN]
1406  qsat ) ! [OUT]
1407 
1408  return
1409  end subroutine atmos_saturation_dqs_dtem_dens_ice_0d
1410 
1411  !-----------------------------------------------------------------------------
1412  ! (d qsi/d T)_{rho}: partial difference of qsat_ice at constant density
1414  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1415  temp, dens, &
1416  dqsdtem )
1417  implicit none
1418  integer, intent(in) :: KA, KS, KE
1419  integer, intent(in) :: IA, IS, IE
1420  integer, intent(in) :: JA, JS, JE
1421 
1422  real(RP), intent(in) :: temp (ka,ia,ja)
1423  real(RP), intent(in) :: dens (ka,ia,ja)
1424 
1425  real(RP), intent(out) :: dqsdtem(ka,ia,ja)
1426 
1427  integer :: k, i, j
1428  !---------------------------------------------------------------------------
1429 
1430  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1431  do j = js, je
1432  do i = is, ie
1433  do k = ks, ke
1434  call atmos_saturation_dqs_dtem_dens_ice_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1435  dqsdtem(k,i,j) ) ! [OUT]
1436  enddo
1437  enddo
1438  enddo
1439 
1440  return
1442 
1443  !-----------------------------------------------------------------------------
1444  ! (d qsi/d T)_{rho}: partial difference of qsat_all at constant density
1445  subroutine atmos_saturation_dqs_dtem_dens_all_0d( &
1446  temp, dens, &
1447  dqsat_dT, &
1448  qsat, qsat_liq, qsat_ice, &
1449  alpha )
1450  implicit none
1451 
1452  real(RP), intent(in) :: temp
1453  real(RP), intent(in) :: dens
1454 
1455  real(RP), intent(out) :: dqsat_dT
1456 
1457  real(RP), intent(out), optional :: qsat
1458  real(RP), intent(out), optional :: qsat_liq
1459  real(RP), intent(out), optional :: qsat_ice
1460  real(RP), intent(out), optional :: alpha
1461 
1462  real(RP) :: qsat_liq_, qsat_ice_, alpha_
1463  real(RP) :: dqsat_dT_liq, dqsat_dT_ice, dalpha_dT
1464 
1465  !---------------------------------------------------------------------------
1466 
1467  call atmos_saturation_dqs_dtem_dens_liq_0d( temp, dens, & ! [IN]
1468  dqsat_dt_liq, qsat_liq_ ) ! [OUT]
1469  call atmos_saturation_dqs_dtem_dens_ice_0d( temp, dens, & ! [IN]
1470  dqsat_dt_ice, qsat_ice_ ) ! [OUT]
1471  call atmos_saturation_alpha ( temp, alpha_ ) ! [IN], [OUT]
1472  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
1473 
1474  dqsat_dt = qsat_liq_ * dalpha_dt + dqsat_dt_liq * ( alpha_ ) &
1475  - qsat_ice_ * dalpha_dt + dqsat_dt_ice * ( 1.0_rp-alpha_ )
1476 
1477  if ( present(qsat) ) qsat = qsat_liq_ * alpha_ + qsat_ice_ * ( 1.0_rp-alpha_ )
1478  if ( present(qsat_liq) ) qsat_liq = qsat_liq_
1479  if ( present(qsat_ice) ) qsat_ice = qsat_ice_
1480  if ( present(alpha ) ) alpha = alpha_
1481 
1482  return
1483  end subroutine atmos_saturation_dqs_dtem_dens_all_0d
1484 
1485  !-----------------------------------------------------------------------------
1486  ! (d qs/d T)_{p} and (d qs/d p)_{T}
1487  subroutine atmos_saturation_dqs_dtem_dpre_liq_0d( &
1488  temp, pres, qdry, &
1489  dqsat_dT, dqsat_dP, &
1490  qsat, psat )
1491  use scale_atmos_hydrometeor, only: &
1492  hydrometeor_lhv => atmos_hydrometeor_lhv
1493  implicit none
1494 
1495  real(RP), intent(in) :: temp
1496  real(RP), intent(in) :: pres
1497  real(RP), intent(in) :: qdry
1498 
1499  real(RP), intent(out) :: dqsat_dT
1500  real(RP), intent(out) :: dqsat_dP
1501 
1502  real(RP), intent(out), optional :: qsat
1503  real(RP), intent(out), optional :: psat
1504 
1505  real(RP) :: LHV ! latent heat of vaporization [J/kg]
1506  real(RP) :: psat_
1507  real(RP) :: den1, den2 ! denominator
1508 
1509  !---------------------------------------------------------------------------
1510 
1511  call hydrometeor_lhv( temp, lhv )
1512  call atmos_saturation_psat_liq( temp, psat_ ) ! [IN], [OUT]
1513 
1514  den1 = ( pres - (1.0_rp-epsvap) * psat_ )**2
1515  den2 = den1 * rvap * temp**2
1516 
1517  dqsat_dp = - epsvap * psat_ / den1
1518  dqsat_dt = epsvap * psat_ / den2 * lhv * pres
1519 
1520  if ( present(qsat) ) &
1521  call atmos_saturation_psat2qsat_pres( psat_, pres, qdry, & ! [IN]
1522  qsat ) ! [OUT]
1523  if ( present(psat) ) psat = psat_
1524 
1525  return
1526  end subroutine atmos_saturation_dqs_dtem_dpre_liq_0d
1527 
1528  !-----------------------------------------------------------------------------
1529  ! (d qs/d T)_{p} and (d qs/d p)_{T}
1531  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1532  temp, pres, qdry, &
1533  dqsat_dT, dqsat_dP )
1534  implicit none
1535  integer, intent(in) :: KA, KS, KE
1536  integer, intent(in) :: IA, IS, IE
1537  integer, intent(in) :: JA, JS, JE
1538 
1539  real(RP), intent(in) :: temp (ka,ia,ja)
1540  real(RP), intent(in) :: pres (ka,ia,ja)
1541  real(RP), intent(in) :: qdry (ka,ia,ja)
1542 
1543  real(RP), intent(out) :: dqsat_dT(ka,ia,ja)
1544  real(RP), intent(out) :: dqsat_dP(ka,ia,ja)
1545 
1546  integer :: k, i, j
1547  !---------------------------------------------------------------------------
1548 
1549  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1550  do j = js, je
1551  do i = is, ie
1552  do k = ks, ke
1553  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1554  dqsat_dt(k,i,j), dqsat_dp(k,i,j) ) ! [OUT]
1555  enddo
1556  enddo
1557  enddo
1558 
1559  return
1561 
1562  !-----------------------------------------------------------------------------
1563  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
1564  !-----------------------------------------------------------------------------
1565  subroutine atmos_saturation_dqs_dtem_dpre_ice_0d( &
1566  temp, pres, qdry, &
1567  dqsat_dT, dqsat_dP, &
1568  qsat )
1569  use scale_atmos_hydrometeor, only: &
1570  hydrometeor_lhs => atmos_hydrometeor_lhs
1571  implicit none
1572 
1573  real(RP), intent(in) :: temp
1574  real(RP), intent(in) :: pres
1575  real(RP), intent(in) :: qdry
1576 
1577  real(RP), intent(out) :: dqsat_dT
1578  real(RP), intent(out) :: dqsat_dP
1579 
1580  real(RP), intent(out), optional :: qsat
1581 
1582  real(RP) :: LHS ! latent heat of sublimation [J/kg]
1583  real(RP) :: psat
1584  real(RP) :: den1, den2 ! denominator
1585 
1586  !---------------------------------------------------------------------------
1587 
1588  call hydrometeor_lhs( temp, lhs ) ! [IN], [OUT]
1589  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
1590 
1591  den1 = ( pres - (1.0_rp-epsvap) * psat )**2
1592  den2 = den1 * rvap * temp**2
1593 
1594  dqsat_dp = - epsvap * psat / den1
1595  dqsat_dt = epsvap * psat / den2 * lhs * pres
1596 
1597  if ( present(qsat) ) &
1598  call atmos_saturation_psat2qsat_pres( psat, pres, qdry, & ! [IN]
1599  qsat ) ! [OUT]
1600  return
1601  end subroutine atmos_saturation_dqs_dtem_dpre_ice_0d
1602 
1603  !-----------------------------------------------------------------------------
1604  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
1605  !-----------------------------------------------------------------------------
1607  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1608  temp, pres, qdry, &
1609  dqsat_dT, dqsat_dP )
1610  implicit none
1611  integer, intent(in) :: KA, KS, KE
1612  integer, intent(in) :: IA, IS, IE
1613  integer, intent(in) :: JA, JS, JE
1614 
1615  real(RP), intent(in) :: temp (ka,ia,ja)
1616  real(RP), intent(in) :: pres (ka,ia,ja)
1617  real(RP), intent(in) :: qdry (ka,ia,ja)
1618 
1619  real(RP), intent(out) :: dqsat_dT(ka,ia,ja)
1620  real(RP), intent(out) :: dqsat_dP(ka,ia,ja)
1621 
1622  integer :: k, i, j
1623  !---------------------------------------------------------------------------
1624 
1625  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1626  do j = js, je
1627  do i = is, ie
1628  do k = ks, ke
1629  call atmos_saturation_dqs_dtem_dpre_ice( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1630  dqsat_dt(k,i,j), dqsat_dp(k,i,j) ) ! [OUT]
1631  enddo
1632  enddo
1633  enddo
1634 
1635  return
1637 
1638  !-----------------------------------------------------------------------------
1639  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
1640  !-----------------------------------------------------------------------------
1641  subroutine atmos_saturation_dqs_dtem_dpre_all_0d( &
1642  temp, pres, qdry, &
1643  dqsat_dT, dqsat_dP, &
1644  qsat, qsat_liq, qsat_ice, &
1645  alpha )
1646  implicit none
1647 
1648  real(RP), intent(in) :: temp
1649  real(RP), intent(in) :: pres
1650  real(RP), intent(in) :: qdry
1651 
1652  real(RP), intent(out) :: dqsat_dT
1653  real(RP), intent(out) :: dqsat_dP
1654 
1655  real(RP), intent(out), optional :: qsat
1656  real(RP), intent(out), optional :: qsat_liq
1657  real(RP), intent(out), optional :: qsat_ice
1658  real(RP), intent(out), optional :: alpha
1659 
1660  real(RP) :: qsat_liq_, qsat_ice_, alpha_
1661  real(RP) :: dqsat_dT_liq, dqsat_dT_ice
1662  real(RP) :: dqsat_dP_liq, dqsat_dP_ice
1663  real(RP) :: dalpha_dT
1664 
1665  !---------------------------------------------------------------------------
1666 
1667  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp, pres, qdry, & ! [IN]
1668  dqsat_dt_liq, dqsat_dp_liq, & ! [OUT]
1669  qsat_liq_ ) ! [OUT]
1670  call atmos_saturation_dqs_dtem_dpre_ice_0d( temp, pres, qdry, & ! [IN]
1671  dqsat_dt_ice, dqsat_dp_ice, & ! [OUT]
1672  qsat_ice_ ) ! [OUT]
1673  call atmos_saturation_alpha ( temp, alpha_ ) ! [IN], [OUT]
1674  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
1675 
1676  dqsat_dt = qsat_liq_ * dalpha_dt + dqsat_dt_liq * ( alpha_ ) &
1677  - qsat_ice_ * dalpha_dt + dqsat_dt_ice * ( 1.0_rp-alpha_ )
1678  dqsat_dp = dqsat_dp_liq * ( alpha_ ) &
1679  + dqsat_dp_ice * ( 1.0_rp-alpha_ )
1680 
1681  if ( present(qsat) ) qsat = qsat_liq * alpha_ + qsat_ice * ( 1.0_rp-alpha_ )
1682  if ( present(qsat_liq) ) qsat_liq = qsat_liq_
1683  if ( present(qsat_ice) ) qsat_ice = qsat_ice_
1684  if ( present(alpha ) ) alpha = alpha_
1685 
1686  return
1687  end subroutine atmos_saturation_dqs_dtem_dpre_all_0d
1688 
1689  !-----------------------------------------------------------------------------
1691  !-----------------------------------------------------------------------------
1692  subroutine atmos_saturation_tdew_liq_0d( &
1693  DENS, TEMP, QV, &
1694  Tdew, &
1695  converged )
1696  use scale_const, only: &
1697  undef => const_undef
1698  use scale_atmos_hydrometeor, only: &
1699  atmos_hydrometeor_lhv
1700  real(RP), intent(in) :: DENS
1701  real(RP), intent(in) :: TEMP
1702  real(RP), intent(in) :: QV
1703 
1704  real(RP), intent(out) :: Tdew
1705  logical, intent(out) :: converged
1706 
1707  real(RP), parameter :: A = 17.625_rp
1708  real(RP), parameter :: B = 243.04_rp
1709  real(RP), parameter :: C = 610.94_rp
1710  integer, parameter :: itelim = 100
1711  real(RP), parameter :: criteria = 0.1_rp**(2+rp/2)
1712 
1713  real(RP) :: lhv
1714  real(RP) :: pvap, psat
1715  real(RP) :: dpsat_dT
1716  real(RP) :: dTdew
1717 
1718  integer :: ite
1719  !---------------------------------------------------------------------------
1720 
1721  pvap = dens * qv * rvap * temp
1722 
1723  if ( pvap < psat_min_liq ) then
1724  converged = .true.
1725  tdew = undef
1726  return
1727  end if
1728 
1729  ! first guess is calculated by Alduchov and Eskridge (1996)
1730  ! See Lawrence (2005) BAMS
1731  tdew = b * log( pvap / c ) / ( a - log( pvap / c ) ) + tem00
1732  converged = .false.
1733  do ite = 1, itelim
1734 
1735  call atmos_saturation_psat_liq( tdew, psat ) ! [IN], [OUT]
1736  call atmos_hydrometeor_lhv( tdew, lhv )
1737 
1738  dpsat_dt = psat * lhv / ( rvap * tdew**2 )
1739  dtdew = ( psat - pvap ) / dpsat_dt
1740  if ( dtdew < criteria ) then
1741  converged = .true.
1742  exit
1743  end if
1744 
1745  tdew = tdew - dtdew
1746  end do
1747  if( .not. converged ) then
1748  log_warn("ATMOS_SATURATION_tdew_liq_0D",*) dens, temp, qv, pvap, tdew, dtdew, dpsat_dt
1749  endif
1750 
1751  return
1752  end subroutine atmos_saturation_tdew_liq_0d
1753 
1754 !OCL SERIAL
1755  subroutine atmos_saturation_tdew_liq_1d( &
1756  KA, KS, KE, &
1757  DENS, TEMP, QV, &
1758  Tdew, &
1759  converged )
1760  integer, intent(in) :: KA, KS, KE
1761 
1762  real(RP), intent(in) :: DENS(ka)
1763  real(RP), intent(in) :: TEMP(ka)
1764  real(RP), intent(in) :: QV (ka)
1765 
1766  real(RP), intent(out) :: Tdew(ka)
1767  logical, intent(out) :: converged
1768 
1769  integer :: k
1770  !---------------------------------------------------------------------------
1771 
1772  do k = ks, ke
1773  call atmos_saturation_tdew_liq_0d( dens(k), temp(k), qv(k), & ! [IN]
1774  tdew(k), converged ) ! [OUT]
1775  if ( .not. converged ) exit
1776  end do
1777 
1778  end subroutine atmos_saturation_tdew_liq_1d
1779 
1780  subroutine atmos_saturation_tdew_liq_3d( &
1781  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1782  DENS, TEMP, QV, &
1783  Tdew )
1784  use scale_prc, only: &
1785  prc_abort
1786  integer, intent(in) :: KA, KS, KE
1787  integer, intent(in) :: IA, IS, IE
1788  integer, intent(in) :: JA, JS, JE
1789 
1790  real(RP), intent(in) :: DENS(ka,ia,ja)
1791  real(RP), intent(in) :: TEMP(ka,ia,ja)
1792  real(RP), intent(in) :: QV (ka,ia,ja)
1793 
1794  real(RP), intent(out) :: Tdew(ka,ia,ja)
1795 
1796  logical :: converged, error
1797  integer :: k, i, j
1798  !---------------------------------------------------------------------------
1799 
1800  error = .false.
1801  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
1802  !$omp private(converged)
1803  do j = js, je
1804  do i = is, ie
1805  do k = ks, ke
1806  call atmos_saturation_tdew_liq_0d( dens(k,i,j), temp(k,i,j), qv(k,i,j), & ! [IN]
1807  tdew(k,i,j), converged ) ! [OUT]
1808  if ( .not. converged ) then
1809  log_error("ATMOS_SATURATION_tdew_liq_3D",*) 'not converged! ', k,i,j
1810  error = .true.
1811  exit
1812  end if
1813  end do
1814  end do
1815  end do
1816 
1817  if ( error ) call prc_abort
1818 
1819  end subroutine atmos_saturation_tdew_liq_3d
1820 
1821  !-----------------------------------------------------------------------------
1828  !-----------------------------------------------------------------------------
1829  subroutine atmos_saturation_pote_0d( &
1830  DENS, POTT, TEMP, QV, &
1831  POTE )
1832  use scale_const, only: &
1833  rdry => const_rdry, &
1834  cpdry => const_cpdry
1835  use scale_atmos_hydrometeor, only: &
1836  atmos_hydrometeor_lhv
1837  real(RP), intent(in) :: DENS
1838  real(RP), intent(in) :: POTT
1839  real(RP), intent(in) :: TEMP
1840  real(RP), intent(in) :: QV
1841 
1842  real(RP), intent(out) :: POTE
1843 
1844  real(RP) :: TL
1845  real(RP) :: Pv
1846  real(RP) :: LHV
1847 
1848  pv = dens * qv * rvap * temp
1849  tl = 55.0_rp + 2840.0_rp / ( cpdry / rdry * log(temp) - log(pv) - 4.805_rp )
1850  call atmos_hydrometeor_lhv( temp, lhv ) ! [IN], [OUT]
1851 
1852  pote = pott * exp( lhv * qv / ( cpdry * temp ) &
1853  * 1.0784_rp * ( 1.0_rp + 0.810_rp * qv ) )
1854 
1855  return
1856  end subroutine atmos_saturation_pote_0d
1857 
1858 !OCL SERIAL
1859  subroutine atmos_saturation_pote_1d( &
1860  KA, KS, KE, &
1861  DENS, POTT, TEMP, QV, &
1862  POTE )
1863  integer, intent(in) :: KA, KS, KE
1864 
1865  real(RP), intent(in) :: DENS(ka)
1866  real(RP), intent(in) :: POTT(ka)
1867  real(RP), intent(in) :: TEMP(ka)
1868  real(RP), intent(in) :: QV (ka)
1869 
1870  real(RP), intent(out) :: POTE(ka)
1871 
1872  integer :: k
1873 
1874  do k = ks, ke
1875  call atmos_saturation_pote_0d( dens(k), pott(k), temp(k), qv(k), & ! [IN]
1876  pote(k) ) ! [OUT]
1877  end do
1878 
1879  return
1880  end subroutine atmos_saturation_pote_1d
1881 
1882  subroutine atmos_saturation_pote_3d( &
1883  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1884  DENS, POTT, TEMP, QV, &
1885  POTE )
1886  integer, intent(in) :: KA, KS, KE
1887  integer, intent(in) :: IA, IS, IE
1888  integer, intent(in) :: JA, JS, JE
1889 
1890  real(RP), intent(in) :: DENS(ka,ia,ja)
1891  real(RP), intent(in) :: POTT(ka,ia,ja)
1892  real(RP), intent(in) :: TEMP(ka,ia,ja)
1893  real(RP), intent(in) :: QV (ka,ia,ja)
1894 
1895  real(RP), intent(out) :: POTE(ka,ia,ja)
1896 
1897  integer :: k, i, j
1898 
1899  !$omp parallel do
1900  do j = js, je
1901  do i = is, ie
1902  do k = ks, ke
1903  call atmos_saturation_pote_0d( &
1904  dens(k,i,j), pott(k,i,j), temp(k,i,j), qv(k,i,j), & ! [IN]
1905  pote(k,i,j) ) ! [OUT]
1906  end do
1907  end do
1908  end do
1909 
1910  return
1911  end subroutine atmos_saturation_pote_3d
1912 
1913  !-----------------------------------------------------------------------------
1915  !-----------------------------------------------------------------------------
1916 !OCL SERIAL
1917  subroutine atmos_saturation_moist_conversion_dens_liq_0d( &
1918  DENS, Emoist0, &
1919  TEMP, QV, QC, CPtot, CVtot, &
1920  converged )
1921  use scale_prc, only: &
1922  prc_abort
1923  use scale_atmos_hydrometeor, only: &
1924  cp_vapor, &
1925  cp_water, &
1926  cv_vapor, &
1927  cv_water, &
1928  lhv
1929  implicit none
1930 
1931  real(RP), intent(in) :: DENS
1932  real(RP), intent(in) :: Emoist0
1933 
1934  real(RP), intent(inout) :: TEMP
1935  real(RP), intent(inout) :: QV
1936  real(RP), intent(inout) :: QC
1937  real(RP), intent(inout) :: CPtot
1938  real(RP), intent(inout) :: CVtot
1939 
1940  logical, intent(out) :: converged
1941 
1942  ! working
1943  real(RP) :: QSUM
1944  real(RP) :: CVtot0
1945  real(RP) :: qsat
1946  real(RP) :: Emoist ! moist internal energy
1947 
1948  ! d(X)/dT
1949  real(RP) :: dqsatl_dT
1950  real(RP) :: dqc_dT
1951  real(RP) :: dCVtot_dT
1952  real(RP) :: dEmoist_dT
1953  real(RP) :: dtemp
1954 
1955  integer, parameter :: itelim = 100
1956  real(RP), parameter :: dtemp_criteria = 0.1_rp**(2+rp/2)
1957  integer :: ite
1958  !---------------------------------------------------------------------------
1959 
1960  qsum = qv + qc
1961 
1962  call atmos_saturation_dens2qsat_liq( temp, dens, & ! [IN]
1963  qsat ) ! [OUT]
1964 
1965  if ( qsum <= qsat ) then
1966  qv = qsum
1967  cptot = cptot + qc * ( cp_vapor - cp_water )
1968  cvtot = cvtot + qc * ( cv_vapor - cv_water )
1969  qc = 0.0_rp
1970  temp = ( emoist0 - lhv * qv ) / cvtot
1971  converged = .true.
1972  return
1973  end if
1974 
1975  cvtot0 = cvtot
1976 
1977  converged = .false.
1978  do ite = 1, itelim
1979 
1980  ! get qsat and dX/dT under the temp
1981  call atmos_saturation_dqs_dtem_dens_liq_0d( temp, dens, & ! [IN]
1982  dqsatl_dt, qsat ) ! [OUT]
1983 
1984  dqc_dt = - dqsatl_dt
1985 
1986  dcvtot_dt = dqsatl_dt * cv_vapor &
1987  + dqc_dt * cv_water
1988 
1989  demoist_dt = temp * dcvtot_dt + cvtot + dqsatl_dt * lhv
1990 
1991  ! diagnose quantities with the qsat
1992  qc = qsum - qsat
1993  cvtot = cvtot0 + ( cv_vapor - cv_water ) * ( qsat - qv )
1994  emoist = temp * cvtot + qsat * lhv
1995 
1996 
1997  ! update temp by the newtonian method
1998  dtemp = ( emoist - emoist0 ) / demoist_dt
1999  temp = temp - dtemp
2000 
2001  if ( abs(dtemp) < dtemp_criteria ) then
2002  converged = .true.
2003  exit
2004  endif
2005 
2006  if( temp*0.0_rp /= 0.0_rp ) exit
2007  enddo
2008 
2009  cptot = cptot + ( cp_vapor - cp_water ) * ( qsat - qv )
2010 
2011  qv = qsat
2012 
2013  return
2014  end subroutine atmos_saturation_moist_conversion_dens_liq_0d
2015 
2016  !-----------------------------------------------------------------------------
2018  !-----------------------------------------------------------------------------
2019 !OCL SERIAL
2021  DENS, Emoist0, &
2022  TEMP, QV, QC, QI, CPtot, CVtot, &
2023  converged )
2024  use scale_prc, only: &
2025  prc_abort
2026  use scale_atmos_hydrometeor, only: &
2027  cp_vapor, &
2028  cp_water, &
2029  cp_ice, &
2030  cv_vapor, &
2031  cv_water, &
2032  cv_ice, &
2033  lhv, &
2034  lhf
2035  implicit none
2036 
2037  real(RP), intent(in) :: DENS
2038  real(RP), intent(in) :: Emoist0
2039 
2040  real(RP), intent(inout) :: TEMP
2041  real(RP), intent(inout) :: QV
2042  real(RP), intent(inout) :: QC
2043  real(RP), intent(inout) :: QI
2044  real(RP), intent(inout) :: CPtot
2045  real(RP), intent(inout) :: CVtot
2046 
2047  logical, intent(out) :: converged
2048 
2049  ! working
2050  real(RP) :: QSUM
2051  real(RP) :: TEMP0
2052  real(RP) :: QV0, QC0, QI0
2053  real(RP) :: CVtot0
2054  real(RP) :: alpha
2055  real(RP) :: qsat
2056  real(RP) :: Emoist ! moist internal energy
2057 
2058  ! d(X)/dT
2059  real(RP) :: dalpha_dT
2060  real(RP) :: dqsat_dT
2061  real(RP) :: dqc_dT, dqi_dT
2062  real(RP) :: dCVtot_dT
2063  real(RP) :: dEmoist_dT
2064  real(RP) :: dtemp
2065 
2066  integer, parameter :: itelim = 100
2067  real(RP), parameter :: dtemp_criteria = 0.1_rp**(2+rp/2)
2068  integer :: ite
2069  !---------------------------------------------------------------------------
2070 
2071  temp0 = temp
2072  cvtot0 = cvtot
2073  qv0 = qv
2074  qc0 = qc
2075  qi0 = qi
2076  qsum = qv + qc + qi
2077 
2078  call atmos_saturation_dens2qsat_all( temp, dens, & ! [IN]
2079  qsat ) ! [OUT]
2080 
2081  if ( qsum <= qsat ) then
2082  qv = qsum
2083  cptot = cptot + qc * ( cp_vapor - cp_water ) + qi * ( cp_vapor - cp_ice )
2084  cvtot = cvtot + qc * ( cv_vapor - cv_water ) + qi * ( cv_vapor - cv_ice )
2085  qc = 0.0_rp
2086  qi = 0.0_rp
2087  temp = ( emoist0 - lhv * qv ) / cvtot
2088  converged = .true.
2089  return
2090  end if
2091 
2092  converged = .false.
2093  do ite = 1, itelim
2094 
2095  ! dX/dT
2096  call atmos_saturation_dqs_dtem_dens_all_0d( temp, dens, & ! [IN]
2097  dqsat_dt, & ! [OUT]
2098  qsat=qsat, alpha=alpha ) ! [OUT]
2099  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
2100 
2101  dqc_dt = ( qsum - qv ) * dalpha_dt - dqsat_dt * ( alpha )
2102  dqi_dt = -( qsum - qv ) * dalpha_dt - dqsat_dt * ( 1.0_rp-alpha )
2103 
2104  dcvtot_dt = dqsat_dt * cv_vapor &
2105  + dqc_dt * cv_water &
2106  + dqi_dt * cv_ice
2107 
2108  demoist_dt = temp * dcvtot_dt + cvtot + dqsat_dt * lhv - dqi_dt * lhf
2109 
2110  ! Saturation
2111  qv = qsat
2112  qc = ( qsum - qsat ) * ( alpha )
2113  qi = ( qsum - qsat ) * ( 1.0_rp - alpha )
2114 
2115  cvtot = cvtot0 &
2116  + cv_vapor * ( qv - qv0 ) &
2117  + cv_water * ( qc - qc0 ) &
2118  + cv_ice * ( qi - qi0 )
2119 
2120  emoist = temp * cvtot + qv * lhv - qi * lhf
2121 
2122  ! update temp by the newtonian method
2123  dtemp = ( emoist - emoist0 ) / demoist_dt
2124  temp = temp - dtemp
2125 
2126  if ( abs(dtemp) < dtemp_criteria ) then
2127  converged = .true.
2128  exit
2129  endif
2130 
2131  if( temp*0.0_rp /= 0.0_rp ) exit
2132  enddo
2133 
2134  cptot = cptot &
2135  + cp_vapor * ( qv - qv0 ) &
2136  + cp_water * ( qc - qc0 ) &
2137  + cp_ice * ( qi - qi0 )
2138 
2139  return
2141 
2142  !-----------------------------------------------------------------------------
2144  !-----------------------------------------------------------------------------
2146  PRES, Entr, Qdry, &
2147  QV, QC, &
2148  Rtot, CPtot, &
2149  TEMP, &
2150  converged )
2151  use scale_const, only: &
2152  eps => const_eps, &
2153  lhv0 => const_lhv0
2154  use scale_atmos_hydrometeor, only: &
2155  atmos_hydrometeor_entr, &
2156  atmos_hydrometeor_entr2temp, &
2157  cp_vapor, &
2158  cp_water
2159 
2160  real(RP), intent(in) :: PRES
2161  real(RP), intent(in) :: ENTR
2162  real(RP), intent(in) :: Qdry
2163 
2164  real(RP), intent(inout) :: QV
2165  real(RP), intent(inout) :: QC
2166  real(RP), intent(inout) :: CPtot
2167  real(RP), intent(inout) :: Rtot
2168 
2169  real(RP), intent(out) :: TEMP
2170  logical, intent(out) :: converged
2171 
2172  real(RP), parameter :: TEMMIN = 0.1_rp
2173  real(RP), parameter :: criteria = 1.e-8_rp
2174  integer, parameter :: itelim = 100
2175  integer :: ite
2176 
2177  real(RP) :: Qsum
2178  real(RP) :: qsat, psat
2179 
2180  real(RP) :: TEMP_ite
2181  real(RP) :: QV_ite
2182  real(RP) :: ENTR_ite
2183  real(RP) :: Rtot_ite
2184  real(RP) :: CPtot_ite
2185 
2186  real(RP) :: dqsat_dT, dqsat_dP
2187 
2188  real(RP) :: TEMP_prev
2189  real(RP) :: dENTR_dT
2190  !---------------------------------------------------------------------------
2191 
2192  qsum = qv + qc
2193 
2194  ! all liquid evaporates
2195  rtot = rtot + qc * rvap
2196  cptot = cptot + qc * ( cp_vapor - cp_water )
2197  qv = qsum
2198  qc = 0.0_rp
2199 
2200  call atmos_hydrometeor_entr2temp( entr, pres, qsum, 0.0_rp, qdry, & ! [IN]
2201  rtot, cptot, & ! [IN]
2202  temp ) ! [OUT]
2203 
2204  call atmos_saturation_pres2qsat_liq( temp, pres, qdry, & ! [IN]
2205  qsat ) ! [OUT]
2206  if ( qsum <= qsat ) then
2207  ! unsaturated
2208 
2209  converged = .true.
2210 
2211  return
2212 
2213  else
2214  ! saturated
2215 
2216  temp_ite = temp
2217 
2218  do ite = 1, itelim
2219 
2220  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp_ite, pres, qdry, & ! [IN]
2221  dqsat_dt, dqsat_dp, & ! [OUT]
2222  qsat=qsat, psat=psat ) ! [OUT]
2223 
2224  qv_ite = min( qsum, qsat )
2225 
2226  rtot_ite = rtot - ( qsum - qv_ite ) * rvap
2227  cptot_ite = cptot - ( qsum - qv_ite ) * ( cp_vapor - cp_water )
2228 
2229  dentr_dt = cptot_ite / temp_ite &
2230  + ( ( cp_vapor - cp_water ) * log( temp_ite / tem00 ) &
2231  - rvap * log( psat/ psat0 ) &
2232  + lhv0 / tem00 &
2233  ) * dqsat_dt
2234 
2235  call atmos_hydrometeor_entr( temp_ite, pres, & ! [IN]
2236  qv_ite, 0.0_rp, qdry, & ! [IN] ! QI = 0
2237  rtot_ite, cptot_ite, & ! [IN]
2238  entr_ite ) ! [OUT]
2239 
2240  temp_prev = temp_ite
2241  temp_ite = temp_ite - ( entr_ite - entr ) / max( dentr_dt, eps )
2242  temp_ite = max( temp_ite, temmin )
2243 
2244  if( abs(temp_ite-temp_prev) < criteria ) then
2245  converged = .true.
2246  exit
2247  end if
2248 
2249  enddo
2250 
2251  end if
2252 
2253  qv = qv_ite
2254  qc = qsum - qv
2255  cptot = cptot_ite
2256  rtot = rtot_ite
2257  temp = temp_ite
2258 
2259  return
2261 
2262 end module scale_atmos_saturation
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:81
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
module atmosphere / saturation
subroutine, public atmos_saturation_setup
Setup.
real(rp), public cv_ice
CV for ice [J/kg/K].
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:67
real(rp), public cp_ice
CP for ice [J/kg/K].
subroutine atmos_saturation_dqs_dtem_dens_ice_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, temp, dens, dqsdtem)
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:66
subroutine atmos_saturation_pote_1d(KA, KS, KE, DENS, POTT, TEMP, QV, POTE)
subroutine atmos_saturation_pote_0d(DENS, POTT, TEMP, QV, POTE)
calculate equivalent potential temperature Bolton, D., 1980: The computation of equivalent potential ...
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
subroutine atmos_saturation_alpha_0d(temp, alpha)
calc liquid/ice separation factor (0D)
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:65
subroutine atmos_saturation_dqs_dtem_dpre_ice_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, temp, pres, qdry, dqsat_dT, dqsat_dP)
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
real(rp), public const_undef
Definition: scale_const.F90:41
subroutine atmos_saturation_dqs_dtem_dens_liq_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, temp, dens, dqsdtem)
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:77
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:75
module PROCESS
Definition: scale_prc.F90:11
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
Definition: scale_const.F90:78
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:76
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:69
real(rp), public lhf
latent heat of fusion for use [J/kg]
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:63
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
real(rp), public lhv
latent heat of vaporization for use [J/kg]
module CONSTANT
Definition: scale_const.F90:11
subroutine atmos_saturation_moist_conversion_pres_liq_0d(PRES, Entr, Qdry, QV, QC, Rtot, CPtot, TEMP, converged)
Iterative moist conversion for liquid water at constant pressure.
subroutine atmos_saturation_dqs_dtem_dpre_liq_3d(KA, KS, KE, IA, IS, IE, JA, JS, JE, temp, pres, qdry, dqsat_dT, dqsat_dP)
module profiler
Definition: scale_prof.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
subroutine, public atmos_hydrometeor_setup
Setup.
module PRECISION
subroutine atmos_saturation_moist_conversion_dens_all_0d(DENS, Emoist0, TEMP, QV, QC, QI, CPtot, CVtot, converged)
Iterative moist conversion (liquid/ice mixture) at constant density (volume)
module STDIO
Definition: scale_io.F90:10
subroutine atmos_saturation_tdew_liq_1d(KA, KS, KE, DENS, TEMP, QV, Tdew, converged)
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:64
character(len=h_short), public const_thermodyn_type
internal energy type
Definition: scale_const.F90:96
real(rp), public cp_vapor
CP for vapor [J/kg/K].
integer, parameter, public rp
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].