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_all_para
46  public :: atmos_saturation_pres2qsat_liq
47  public :: atmos_saturation_pres2qsat_liq_para
48  public :: atmos_saturation_pres2qsat_ice
49  public :: atmos_saturation_pres2qsat_ice_para
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 
57  public :: atmos_saturation_dqs_dtem_dens_liq
58  public :: atmos_saturation_dqs_dtem_dens_ice
59  public :: atmos_saturation_dqs_dtem_dpre_liq
60  public :: atmos_saturation_dqs_dtem_dpre_ice
61 
62  public :: atmos_saturation_tdew_liq
63 
64  public :: atmos_saturation_pote
65 
66  public :: atmos_saturation_moist_conversion_dens_liq
67  public :: atmos_saturation_moist_conversion_dens_all
68  public :: atmos_saturation_moist_conversion_pres_liq
69 
70  interface atmos_saturation_alpha
71  module procedure atmos_saturation_alpha_0d
72  module procedure atmos_saturation_alpha_1d
73  module procedure atmos_saturation_alpha_3d
74  end interface atmos_saturation_alpha
75 
76  interface atmos_saturation_psat_all
77  module procedure atmos_saturation_psat_all_0d
78  module procedure atmos_saturation_psat_all_1d
79  module procedure atmos_saturation_psat_all_2d
80  module procedure atmos_saturation_psat_all_3d
81  end interface atmos_saturation_psat_all
82  interface atmos_saturation_psat_liq
83  module procedure atmos_saturation_psat_liq_0d
84  module procedure atmos_saturation_psat_liq_1d
85  module procedure atmos_saturation_psat_liq_2d
86  module procedure atmos_saturation_psat_liq_3d
87  end interface atmos_saturation_psat_liq
88  interface atmos_saturation_psat_ice
89  module procedure atmos_saturation_psat_ice_0d
90  module procedure atmos_saturation_psat_ice_1d
91  module procedure atmos_saturation_psat_ice_2d
92  module procedure atmos_saturation_psat_ice_3d
93  end interface atmos_saturation_psat_ice
94 
95  interface atmos_saturation_psat2qsat_pres
96  module procedure atmos_saturation_psat2qsat_pres_0d
97  module procedure atmos_saturation_psat2qsat_pres_qdry_0d
98  end interface atmos_saturation_psat2qsat_pres
99  interface atmos_saturation_psat2qsat_dens
100  module procedure atmos_saturation_psat2qsat_dens_0d
101  end interface atmos_saturation_psat2qsat_dens
102 
103  interface atmos_saturation_pres2qsat_all
104  module procedure atmos_saturation_pres2qsat_all_0d
105  module procedure atmos_saturation_pres2qsat_all_1d
106  module procedure atmos_saturation_pres2qsat_all_2d
107  module procedure atmos_saturation_pres2qsat_all_3d
108  module procedure atmos_saturation_pres2qsat_all_qdry_0d
109  module procedure atmos_saturation_pres2qsat_all_qdry_1d
110  module procedure atmos_saturation_pres2qsat_all_qdry_2d
111  module procedure atmos_saturation_pres2qsat_all_qdry_3d
112  end interface atmos_saturation_pres2qsat_all
113  interface atmos_saturation_pres2qsat_all_para
114  module procedure atmos_saturation_pres2qsat_all_1d_para
115  module procedure atmos_saturation_pres2qsat_all_qdry_1d_para
116  end interface atmos_saturation_pres2qsat_all_para
117  interface atmos_saturation_pres2qsat_liq
118  module procedure atmos_saturation_pres2qsat_liq_0d
119  module procedure atmos_saturation_pres2qsat_liq_1d
120  module procedure atmos_saturation_pres2qsat_liq_3d
121  module procedure atmos_saturation_pres2qsat_liq_qdry_0d
122  module procedure atmos_saturation_pres2qsat_liq_qdry_1d
123  module procedure atmos_saturation_pres2qsat_liq_qdry_3d
124  end interface atmos_saturation_pres2qsat_liq
125  interface atmos_saturation_pres2qsat_liq_para
126  module procedure atmos_saturation_pres2qsat_liq_1d_para
127  module procedure atmos_saturation_pres2qsat_liq_qdry_1d_para
128  end interface atmos_saturation_pres2qsat_liq_para
129  interface atmos_saturation_pres2qsat_ice
130  module procedure atmos_saturation_pres2qsat_ice_0d
131  module procedure atmos_saturation_pres2qsat_ice_1d
132  module procedure atmos_saturation_pres2qsat_ice_3d
133  module procedure atmos_saturation_pres2qsat_ice_qdry_0d
134  module procedure atmos_saturation_pres2qsat_ice_qdry_1d
135  module procedure atmos_saturation_pres2qsat_ice_qdry_3d
136  end interface atmos_saturation_pres2qsat_ice
137  interface atmos_saturation_pres2qsat_ice_para
138  module procedure atmos_saturation_pres2qsat_ice_1d_para
139  module procedure atmos_saturation_pres2qsat_ice_qdry_1d_para
140  end interface atmos_saturation_pres2qsat_ice_para
141 
142  interface atmos_saturation_dens2qsat_all
143  module procedure atmos_saturation_dens2qsat_all_0d
144  module procedure atmos_saturation_dens2qsat_all_1d
145  module procedure atmos_saturation_dens2qsat_all_3d
146  end interface atmos_saturation_dens2qsat_all
147  interface atmos_saturation_dens2qsat_liq
148  module procedure atmos_saturation_dens2qsat_liq_0d
149  module procedure atmos_saturation_dens2qsat_liq_1d
150  module procedure atmos_saturation_dens2qsat_liq_3d
151  end interface atmos_saturation_dens2qsat_liq
152  interface atmos_saturation_dens2qsat_ice
153  module procedure atmos_saturation_dens2qsat_ice_0d
154  module procedure atmos_saturation_dens2qsat_ice_1d
155  module procedure atmos_saturation_dens2qsat_ice_3d
156  end interface atmos_saturation_dens2qsat_ice
157 
158  interface atmos_saturation_dalphadt
159  module procedure atmos_saturation_dalphadt_0d
160  module procedure atmos_saturation_dalphadt_1d
161  module procedure atmos_saturation_dalphadt_3d
162  end interface atmos_saturation_dalphadt
163 
164  interface atmos_saturation_dqs_dtem_dens_liq
165  module procedure atmos_saturation_dqs_dtem_dens_liq_0d
167  module procedure atmos_saturation_dqs_dtem_dens_liq_3d
168  end interface atmos_saturation_dqs_dtem_dens_liq
169  interface atmos_saturation_dqs_dtem_dens_ice
170  module procedure atmos_saturation_dqs_dtem_dens_ice_0d
172  module procedure atmos_saturation_dqs_dtem_dens_ice_3d
173  end interface atmos_saturation_dqs_dtem_dens_ice
174  interface atmos_saturation_dqs_dtem_dpre_liq
175  module procedure atmos_saturation_dqs_dtem_dpre_liq_0d
177  module procedure atmos_saturation_dqs_dtem_dpre_liq_3d
178  end interface atmos_saturation_dqs_dtem_dpre_liq
179  interface atmos_saturation_dqs_dtem_dpre_ice
180  module procedure atmos_saturation_dqs_dtem_dpre_ice_0d
182  module procedure atmos_saturation_dqs_dtem_dpre_ice_3d
183  end interface atmos_saturation_dqs_dtem_dpre_ice
184 
185  interface atmos_saturation_tdew_liq
186  module procedure atmos_saturation_tdew_liq_0d
187  module procedure atmos_saturation_tdew_liq_1d
188  module procedure atmos_saturation_tdew_liq_3d
189  end interface atmos_saturation_tdew_liq
190 
191  interface atmos_saturation_pote
192  module procedure atmos_saturation_pote_0d
193  module procedure atmos_saturation_pote_1d
194  module procedure atmos_saturation_pote_3d
195  end interface atmos_saturation_pote
196 
197  interface atmos_saturation_moist_conversion_dens_liq
198  module procedure atmos_saturation_moist_conversion_dens_liq_0d
199  end interface atmos_saturation_moist_conversion_dens_liq
200  interface atmos_saturation_moist_conversion_dens_all
202  end interface atmos_saturation_moist_conversion_dens_all
203  interface atmos_saturation_moist_conversion_pres_liq
205  end interface atmos_saturation_moist_conversion_pres_liq
206 
207  !-----------------------------------------------------------------------------
208  !
209  !++ Public parameters & variables
210  !
211  !-----------------------------------------------------------------------------
212  !
213  !++ Private procedure
214  !
215  !-----------------------------------------------------------------------------
216  !
217  !++ Private parameters & variables
218  !
219  real(RP), private, parameter :: TEM_MIN = 10.0_rp
220 
221  real(RP), private :: ATMOS_SATURATION_ULIMIT_TEMP = 273.15_rp
222  real(RP), private :: ATMOS_SATURATION_LLIMIT_TEMP = 233.15_rp
223  !$acc declare create(ATMOS_SATURATION_ULIMIT_TEMP, ATMOS_SATURATION_LLIMIT_TEMP)
224 
225  real(RP), private :: RTEM00
226  real(RP), private :: dalphadT_const
227  real(RP), private :: psat_min_liq
228  real(RP), private :: psat_min_ice
229  !$acc declare create(RTEM00, dalphadT_const, psat_min_liq, psat_min_ice)
230 
231  real(RP), private :: CPovR_liq
232  real(RP), private :: CPovR_ice
233  real(RP), private :: CVovR_liq
234  real(RP), private :: CVovR_ice
235  real(RP), private :: LovR_liq
236  real(RP), private :: LovR_ice
237  !$acc declare create(CPovR_liq, CPovR_ice, CVovR_liq, CVovR_ice, LovR_liq, LovR_ice)
238 
239  !-----------------------------------------------------------------------------
240 contains
241  !-----------------------------------------------------------------------------
243  subroutine atmos_saturation_setup
244  use scale_prc, only: &
245  prc_abort
246  use scale_const, only: &
247  cpvap => const_cpvap, &
248  cvvap => const_cvvap, &
249  cl => const_cl, &
250  ci => const_ci, &
251  lhv00 => const_lhv00, &
252  lhs00 => const_lhs00, &
253  lhv0 => const_lhv0, &
254  lhs0 => const_lhs0, &
256  use scale_atmos_hydrometeor, only: &
258  implicit none
259 
260  namelist / param_atmos_saturation / &
261  atmos_saturation_ulimit_temp, &
262  atmos_saturation_llimit_temp
263 
264  integer :: ierr
265  !---------------------------------------------------------------------------
266 
267  log_newline
268  log_info("ATMOS_SATURATION_setup",*) 'Setup'
269 
270  !--- read namelist
271  rewind(io_fid_conf)
272  read(io_fid_conf,nml=param_atmos_saturation,iostat=ierr)
273  if( ierr < 0 ) then !--- missing
274  log_info("ATMOS_SATURATION_setup",*) 'Not found namelist. Default used.'
275  elseif( ierr > 0 ) then !--- fatal error
276  log_error("ATMOS_SATURATION_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SATURATION. Check!'
277  call prc_abort
278  endif
279  log_nml(param_atmos_saturation)
280 
281  rtem00 = 1.0_rp / tem00
282 
283  if ( const_thermodyn_type == 'EXACT' ) then
284 
285  cpovr_liq = ( cpvap - cl ) / rvap
286  cpovr_ice = ( cpvap - ci ) / rvap
287  cvovr_liq = ( cvvap - cl ) / rvap
288  cvovr_ice = ( cvvap - ci ) / rvap
289 
290  lovr_liq = lhv00 / rvap
291  lovr_ice = lhs00 / rvap
292 
293  elseif( const_thermodyn_type == 'SIMPLE' ) then
294 
295  cpovr_liq = 0.0_rp
296  cpovr_ice = 0.0_rp
297  cvovr_liq = 0.0_rp
298  cvovr_ice = 0.0_rp
299 
300  lovr_liq = lhv0 / rvap
301  lovr_ice = lhs0 / rvap
302 
303  endif
304 
305  dalphadt_const = 1.0_rp / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
306 
307  log_newline
308  log_info("ATMOS_SATURATION_setup",'(1x,A,F7.2,A,F7.2)') 'Temperature range for liquid/ice mixture : ', &
309  atmos_saturation_llimit_temp, ' - ', &
310  atmos_saturation_ulimit_temp
311 
312  !$acc update device(ATMOS_SATURATION_ULIMIT_TEMP, ATMOS_SATURATION_LLIMIT_TEMP)
313  !$acc update device(RTEM00, dalphadT_const)
314  !$acc update device(CPovR_liq, CPovR_ice, CVovR_liq, CVovR_ice, LovR_liq, LovR_ice)
315 
317 
318  call atmos_saturation_psat_liq( tem_min, psat_min_liq ) ! [IN], [OUT]
319  call atmos_saturation_psat_ice( tem_min, psat_min_ice ) ! [IN], [OUT]
320 
321  !$acc update device(psat_min_liq, psat_min_ice)
322 
323  return
324  end subroutine atmos_saturation_setup
325 
326  !-----------------------------------------------------------------------------
328  subroutine atmos_saturation_alpha_0d( &
329  temp, &
330  alpha )
331  !$acc routine
332  implicit none
333 
334  real(RP), intent(in) :: temp
335  real(RP), intent(out) :: alpha
336  !---------------------------------------------------------------------------
337 
338  alpha = ( temp - atmos_saturation_llimit_temp ) &
339  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
340 
341  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
342 
343  return
344  end subroutine atmos_saturation_alpha_0d
345 
346  !-----------------------------------------------------------------------------
348 !OCL SERIAL
349  subroutine atmos_saturation_alpha_1d( &
350  KA, KS, KE, &
351  temp, &
352  alpha )
353  !$acc routine vector
354  implicit none
355  integer, intent(in) :: KA, KS, KE
356 
357  real(RP), intent(in) :: temp (KA)
358 
359  real(RP), intent(out) :: alpha(KA)
360 
361  integer :: k
362  !---------------------------------------------------------------------------
363 
364  do k = ks, ke
365  call atmos_saturation_alpha_0d( temp(k), & ! [IN]
366  alpha(k) ) ! [OUT]
367  enddo
368 
369  return
370  end subroutine atmos_saturation_alpha_1d
371 
372  !-----------------------------------------------------------------------------
374  subroutine atmos_saturation_alpha_3d( &
375  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
376  temp, &
377  alpha )
378  implicit none
379  integer, intent(in) :: KA, KS, KE
380  integer, intent(in) :: IA, IS, IE
381  integer, intent(in) :: JA, JS, JE
382 
383  real(RP), intent(in) :: temp (KA,IA,JA)
384  real(RP), intent(out) :: alpha(KA,IA,JA)
385 
386  integer :: k, i, j
387  !---------------------------------------------------------------------------
388 
389  !$omp parallel do OMP_SCHEDULE_ collapse(2)
390  !$acc kernels copyin(temp) copyout(alpha)
391  do j = js, je
392  do i = is, ie
393  do k = ks, ke
394  call atmos_saturation_alpha_0d( temp(k,i,j), & ! [IN]
395  alpha(k,i,j) ) ! [OUT]
396 
397  enddo
398  enddo
399  enddo
400  !$acc end kernels
401 
402  return
403  end subroutine atmos_saturation_alpha_3d
404 
405  !-----------------------------------------------------------------------------
407  subroutine atmos_saturation_psat_all_0d( &
408  temp, &
409  psat )
410  !$acc routine
411  implicit none
412 
413  real(RP), intent(in) :: temp
414  real(RP), intent(out) :: psat
415 
416  real(RP) :: alpha, psatl, psati
417  !---------------------------------------------------------------------------
418 
419  call atmos_saturation_alpha ( temp, alpha )
420  call atmos_saturation_psat_liq( temp, psatl )
421  call atmos_saturation_psat_ice( temp, psati )
422 
423  psat = psatl * ( alpha ) &
424  + psati * ( 1.0_rp - alpha )
425 
426  return
427  end subroutine atmos_saturation_psat_all_0d
428 
429  !-----------------------------------------------------------------------------
431 !OCL SERIAL
432  subroutine atmos_saturation_psat_all_1d( &
433  KA, KS, KE, &
434  temp, &
435  psat )
436  !$acc routine vector
437  implicit none
438  integer, intent(in) :: KA, KS, KE
439 
440  real(RP), intent(in) :: temp(KA)
441  real(RP), intent(out) :: psat(KA)
442 
443  real(RP) :: alpha, psatl, psati
444  integer :: k
445  !---------------------------------------------------------------------------
446 
447  !$acc loop private(alpha,psatl,psati)
448  do k = ks, ke
449 ! call ATMOS_SATURATION_psat_all_0D( temp(k), & ! [IN]
450 ! psat(k) ) ! [OUT]
451  call atmos_saturation_alpha ( temp(k), alpha )
452  call atmos_saturation_psat_liq( temp(k), psatl )
453  call atmos_saturation_psat_ice( temp(k), psati )
454 
455  psat(k) = psatl * ( alpha ) &
456  + psati * ( 1.0_rp - alpha )
457  enddo
458 
459  return
460  end subroutine atmos_saturation_psat_all_1d
461 
462  !-----------------------------------------------------------------------------
464  subroutine atmos_saturation_psat_all_2d( &
465  IA, IS, IE, &
466  JA, JS, JE, &
467  temp, &
468  psat )
469  implicit none
470  integer, intent(in) :: IA, IS, IE
471  integer, intent(in) :: JA, JS, JE
472 
473  real(RP), intent(in) :: temp(IA,JA)
474 
475  real(RP), intent(out) :: psat(IA,JA)
476 
477  real(RP) :: alpha, psatl, psati
478  integer :: i, j
479  !---------------------------------------------------------------------------
480 
481  !$omp parallel do OMP_SCHEDULE_ &
482  !$omp private(alpha,psatl,psati)
483  !$acc kernels copyin(temp) copyout(psat)
484  do j = js, je
485  !$acc loop private(alpha,psatl,psati)
486  do i = is, ie
487 ! call ATMOS_SATURATION_psat_all_0D( temp(i,j), & ! [IN]
488 ! psat(i,j) ) ! [OUT]
489  call atmos_saturation_alpha ( temp(i,j), alpha )
490  call atmos_saturation_psat_liq( temp(i,j), psatl )
491  call atmos_saturation_psat_ice( temp(i,j), psati )
492 
493  psat(i,j) = psatl * ( alpha ) &
494  + psati * ( 1.0_rp - alpha )
495  enddo
496  enddo
497  !$acc end kernels
498 
499  return
500  end subroutine atmos_saturation_psat_all_2d
501 
502  !-----------------------------------------------------------------------------
504  subroutine atmos_saturation_psat_all_3d( &
505  KA, KS, KE, &
506  IA, IS, IE, &
507  JA, JS, JE, &
508  temp, &
509  psat )
510  implicit none
511 
512  integer, intent(in) :: KA, KS, KE
513  integer, intent(in) :: IA, IS, IE
514  integer, intent(in) :: JA, JS, JE
515 
516  real(RP), intent(in) :: temp(KA,IA,JA)
517  real(RP), intent(out) :: psat(KA,IA,JA)
518 
519  real(RP) :: alpha, psatl, psati
520  integer :: k, i, j
521  !---------------------------------------------------------------------------
522 
523  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
524  !$omp private(alpha,psatl,psati)
525  !$acc kernels copyin(temp) copyout(psat)
526  do j = js, je
527  do i = is, ie
528  !$acc loop private(alpha,psatl,psati)
529  do k = ks, ke
530 ! call ATMOS_SATURATION_psat_all_0D( temp(k,i,j), & ! [IN]
531 ! psat(k,i,j) ) ! [OUT]
532  call atmos_saturation_alpha ( temp(k,i,j), alpha )
533  call atmos_saturation_psat_liq( temp(k,i,j), psatl )
534  call atmos_saturation_psat_ice( temp(k,i,j), psati )
535 
536  psat(k,i,j) = psatl * ( alpha ) &
537  + psati * ( 1.0_rp - alpha )
538  enddo
539  enddo
540  enddo
541  !$acc end kernels
542 
543  return
544  end subroutine atmos_saturation_psat_all_3d
545 
546  !-----------------------------------------------------------------------------
548  subroutine atmos_saturation_psat_liq_0d( &
549  temp, &
550  psat )
551  !$acc routine
552  implicit none
553 
554  real(RP), intent(in) :: temp
555  real(RP), intent(out) :: psat
556  !---------------------------------------------------------------------------
557 
558 ! psat = PSAT0 * ( temp * RTEM00 )**CPovR_liq &
559 ! * exp( LovR_liq * ( RTEM00 - 1.0_RP/temp ) )
560  psat = psat0 * exp( log( temp * rtem00 ) * cpovr_liq &
561  + lovr_liq * ( rtem00 - 1.0_rp/temp ) )
562 
563  return
564  end subroutine atmos_saturation_psat_liq_0d
565 
566  !-----------------------------------------------------------------------------
568 !OCL SERIAL
569  subroutine atmos_saturation_psat_liq_1d( &
570  KA, KS, KE, &
571  temp, &
572  psat )
573  !$acc routine vector
574  implicit none
575  integer, intent(in) :: KA, KS, KE
576 
577  real(RP), intent(in) :: temp(KA)
578  real(RP), intent(out) :: psat(KA)
579 
580  integer :: k
581  !---------------------------------------------------------------------------
582 
583  do k = ks, ke
584  call atmos_saturation_psat_liq_0d( temp(k), &
585  psat(k) )
586  enddo
587 
588  return
589  end subroutine atmos_saturation_psat_liq_1d
590 
591  !-----------------------------------------------------------------------------
593  subroutine atmos_saturation_psat_liq_2d( &
594  IA, IS, IE, &
595  JA, JS, JE, &
596  temp, &
597  psat )
598  implicit none
599 
600  integer, intent(in) :: IA, IS, IE
601  integer, intent(in) :: JA, JS, JE
602 
603  real(RP), intent(in) :: temp(IA,JA)
604  real(RP), intent(out) :: psat(IA,JA)
605 
606  integer :: i, j
607  !---------------------------------------------------------------------------
608 
609  !$omp parallel do OMP_SCHEDULE_
610  !$acc kernels copyin(temp) copyout(psat)
611  do j = js, je
612  do i = is, ie
613  call atmos_saturation_psat_liq_0d( temp(i,j), & ! [IN]
614  psat(i,j) ) ! [OUT]
615  enddo
616  enddo
617  !$acc end kernels
618 
619  return
620  end subroutine atmos_saturation_psat_liq_2d
621 
622  !-----------------------------------------------------------------------------
624  subroutine atmos_saturation_psat_liq_3d( &
625  KA, KS, KE, &
626  IA, IS, IE, &
627  JA, JS, JE, &
628  temp, &
629  psat )
630  implicit none
631 
632  integer, intent(in) :: KA, KS, KE
633  integer, intent(in) :: IA, IS, IE
634  integer, intent(in) :: JA, JS, JE
635 
636  real(RP), intent(in) :: temp(KA,IA,JA)
637  real(RP), intent(out) :: psat(KA,IA,JA)
638 
639  integer :: k, i, j
640  !---------------------------------------------------------------------------
641 
642  !$omp parallel do OMP_SCHEDULE_ collapse(2)
643  !$acc kernels copyin(temp) copyout(psat)
644  do j = js, je
645  do i = is, ie
646  do k = ks, ke
647  call atmos_saturation_psat_liq_0d( temp(k,i,j), & ! [IN]
648  psat(k,i,j) ) ! [OUT]
649  enddo
650  enddo
651  enddo
652  !$acc end kernels
653 
654  return
655  end subroutine atmos_saturation_psat_liq_3d
656 
657  !-----------------------------------------------------------------------------
659  subroutine atmos_saturation_psat_ice_0d( &
660  temp, &
661  psat )
662  !$acc routine
663  implicit none
664 
665  real(RP), intent(in) :: temp
666 
667  real(RP), intent(out) :: psat
668  !---------------------------------------------------------------------------
669 
670 ! psat = PSAT0 * ( temp * RTEM00 )**CPovR_ice &
671 ! * exp( LovR_ice * ( RTEM00 - 1.0_RP/temp ) )
672  psat = psat0 * exp( log( temp * rtem00 ) * cpovr_ice &
673  + lovr_ice * ( rtem00 - 1.0_rp/temp ) )
674 
675  return
676  end subroutine atmos_saturation_psat_ice_0d
677 
678  !-----------------------------------------------------------------------------
680 !OCL SERIAL
681  subroutine atmos_saturation_psat_ice_1d( &
682  KA, KS, KE, &
683  temp, &
684  psat )
685  !$acc routine vector
686  implicit none
687  integer, intent(in) :: KA, KS, KE
688 
689  real(RP), intent(in) :: temp(KA)
690  real(RP), intent(out) :: psat(KA)
691 
692  integer :: k
693  !---------------------------------------------------------------------------
694 
695  do k = ks, ke
696  call atmos_saturation_psat_ice_0d( temp(k), & ! [IN]
697  psat(k) ) ! [OUT]
698  enddo
699 
700  return
701  end subroutine atmos_saturation_psat_ice_1d
702 
703  !-----------------------------------------------------------------------------
705  subroutine atmos_saturation_psat_ice_2d( &
706  IA, IS, IE, &
707  JA, JS, JE, &
708  temp, &
709  psat )
710  implicit none
711 
712  integer, intent(in) :: IA, IS, IE
713  integer, intent(in) :: JA, JS, JE
714 
715  real(RP), intent(in) :: temp(IA,JA)
716  real(RP), intent(out) :: psat(IA,JA)
717 
718  integer :: i, j
719  !---------------------------------------------------------------------------
720 
721  !$omp parallel do OMP_SCHEDULE_
722  !$acc kernels copyin(temp) copyout(psat)
723  do j = js, je
724  do i = is, ie
725  call atmos_saturation_psat_ice_0d( temp(i,j), & ! [IN]
726  psat(i,j) ) ! [OUT]
727  enddo
728  enddo
729  !$acc end kernels
730 
731  return
732  end subroutine atmos_saturation_psat_ice_2d
733 
734  !-----------------------------------------------------------------------------
736  subroutine atmos_saturation_psat_ice_3d( &
737  KA, KS, KE, &
738  IA, IS, IE, &
739  JA, JS, JE, &
740  temp, &
741  psat )
742  implicit none
743 
744  integer, intent(in) :: KA, KS, KE
745  integer, intent(in) :: IA, IS, IE
746  integer, intent(in) :: JA, JS, JE
747 
748  real(RP), intent(in) :: temp(KA,IA,JA)
749  real(RP), intent(out) :: psat(KA,IA,JA)
750 
751  integer :: k, i, j
752  !---------------------------------------------------------------------------
753 
754  !$omp parallel do OMP_SCHEDULE_ collapse(2)
755  !$acc kernels copyin(temp) copyout(psat)
756  do j = js, je
757  do i = is, ie
758  do k = ks, ke
759  call atmos_saturation_psat_ice_0d( temp(k,i,j), & ! [IN]
760  psat(k,i,j) ) ! [OUT]
761  enddo
762  enddo
763  enddo
764  !$acc end kernels
765 
766  return
767  end subroutine atmos_saturation_psat_ice_3d
768 
769 
770  !-----------------------------------------------------------------------------
772  subroutine atmos_saturation_psat2qsat_pres_0d( &
773  psat, pres, &
774  qsat )
775  !$acc routine
776  implicit none
777 
778  real(RP), intent(in) :: psat
779  real(RP), intent(in) :: pres
780 
781  real(RP), intent(out) :: qsat
782 
783  ! qdry is assumed to be 1 - qsat
784  qsat = epsvap * psat / ( pres - ( 1.0_rp-epsvap ) * psat )
785 
786  return
787  end subroutine atmos_saturation_psat2qsat_pres_0d
788 
789  !-----------------------------------------------------------------------------
791  subroutine atmos_saturation_psat2qsat_pres_qdry_0d( &
792  psat, pres, qdry, &
793  qsat )
794  !$acc routine
795  implicit none
796 
797  real(RP), intent(in) :: psat
798  real(RP), intent(in) :: pres
799  real(RP), intent(in) :: qdry
800 
801  real(RP), intent(out) :: qsat
802 
803  qsat = epsvap * qdry * psat / ( pres - psat )
804 
805  return
806  end subroutine atmos_saturation_psat2qsat_pres_qdry_0d
807 
808  !-----------------------------------------------------------------------------
810  subroutine atmos_saturation_pres2qsat_all_0d( &
811  temp, pres, &
812  qsat )
813  !$acc routine
814  implicit none
815 
816  real(RP), intent(in) :: temp
817  real(RP), intent(in) :: pres
818 
819  real(RP), intent(out) :: qsat
820 
821  real(RP) :: psat
822  !---------------------------------------------------------------------------
823 
824  call atmos_saturation_psat_all_0d( temp, psat ) ! [IN], [OUT]
825 
826  call atmos_saturation_psat2qsat_pres_0d( psat, pres, & ! [IN]
827  qsat ) ! [OUT]
828 
829  return
830  end subroutine atmos_saturation_pres2qsat_all_0d
831 
832  !-----------------------------------------------------------------------------
834  subroutine atmos_saturation_pres2qsat_all_qdry_0d( &
835  temp, pres, qdry, &
836  qsat )
837  !$acc routine
838  implicit none
839 
840  real(RP), intent(in) :: temp
841  real(RP), intent(in) :: pres
842  real(RP), intent(in) :: qdry
843 
844  real(RP), intent(out) :: qsat
845 
846  real(RP) :: psat
847  !---------------------------------------------------------------------------
848 
849  call atmos_saturation_psat_all_0d( temp, psat ) ! [IN], [OUT]
850 
851  call atmos_saturation_psat2qsat_pres_qdry_0d( psat, pres, qdry, & ! [IN]
852  qsat ) ! [OUT]
853 
854  return
855  end subroutine atmos_saturation_pres2qsat_all_qdry_0d
856 
857  !-----------------------------------------------------------------------------
859 !OCL SERIAL
860  subroutine atmos_saturation_pres2qsat_all_1d( &
861  KA, KS, KE, &
862  temp, pres, &
863  qsat )
864  !$acc routine vector
865  implicit none
866  integer, intent(in) :: KA, KS, KE
867 
868  real(RP), intent(in) :: temp(KA)
869  real(RP), intent(in) :: pres(KA)
870 
871  real(RP), intent(out) :: qsat(KA)
872 
873  integer :: k
874  !---------------------------------------------------------------------------
875 
876  do k = ks, ke
877  call atmos_saturation_pres2qsat_all_0d( temp(k), pres(k), & ! [IN]
878  qsat(k) ) ! [OUT]
879  enddo
880 
881  return
882  end subroutine atmos_saturation_pres2qsat_all_1d
883 
884  !-----------------------------------------------------------------------------
886  subroutine atmos_saturation_pres2qsat_all_1d_para( &
887  KA, KS, KE, &
888  temp, pres, &
889  qsat )
890  implicit none
891  integer, intent(in) :: KA, KS, KE
892 
893  real(RP), intent(in) :: temp(KA)
894  real(RP), intent(in) :: pres(KA)
895 
896  real(RP), intent(out) :: qsat(KA)
897 
898  integer :: k
899  !---------------------------------------------------------------------------
900 
901  !$omp parallel do OMP_SCHEDULE_
902  !$acc kernels copyin(temp,pres) copyout(qsat)
903  do k = ks, ke
904  call atmos_saturation_pres2qsat_all_0d( temp(k), pres(k), & ! [IN]
905  qsat(k) ) ! [OUT]
906  enddo
907  !$acc end kernels
908 
909  return
910  end subroutine atmos_saturation_pres2qsat_all_1d_para
911 
912  !-----------------------------------------------------------------------------
914 !OCL SERIAL
915  subroutine atmos_saturation_pres2qsat_all_qdry_1d( &
916  KA, KS, KE, &
917  temp, pres, qdry, &
918  qsat )
919  !$acc routine vector
920  implicit none
921  integer, intent(in) :: KA, KS, KE
922 
923  real(RP), intent(in) :: temp(KA)
924  real(RP), intent(in) :: pres(KA)
925  real(RP), intent(in) :: qdry(KA)
926 
927  real(RP), intent(out) :: qsat(KA)
928 
929  integer :: k
930  !---------------------------------------------------------------------------
931 
932  do k = ks, ke
933  call atmos_saturation_pres2qsat_all_qdry_0d( temp(k), pres(k), qdry(k), & ! [IN]
934  qsat(k) ) ! [OUT]
935  enddo
936 
937  return
938  end subroutine atmos_saturation_pres2qsat_all_qdry_1d
939 
940  !-----------------------------------------------------------------------------
942  subroutine atmos_saturation_pres2qsat_all_qdry_1d_para( &
943  KA, KS, KE, &
944  temp, pres, qdry, &
945  qsat )
946  implicit none
947  integer, intent(in) :: KA, KS, KE
948 
949  real(RP), intent(in) :: temp(KA)
950  real(RP), intent(in) :: pres(KA)
951  real(RP), intent(in) :: qdry(KA)
952 
953  real(RP), intent(out) :: qsat(KA)
954 
955  integer :: k
956  !---------------------------------------------------------------------------
957 
958  !$omp parallel do OMP_SCHEDULE_
959  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
960  do k = ks, ke
961  call atmos_saturation_pres2qsat_all_qdry_0d( temp(k), pres(k), qdry(k), & ! [IN]
962  qsat(k) ) ! [OUT]
963  enddo
964  !$acc end kernels
965 
966  return
967  end subroutine atmos_saturation_pres2qsat_all_qdry_1d_para
968 
969  !-----------------------------------------------------------------------------
971  subroutine atmos_saturation_pres2qsat_all_2d( &
972  IA, IS, IE, JA, JS, JE, &
973  temp, pres, &
974  qsat )
975  implicit none
976  integer, intent(in) :: IA, IS, IE
977  integer, intent(in) :: JA, JS, JE
978 
979  real(RP), intent(in) :: temp(IA,JA)
980  real(RP), intent(in) :: pres(IA,JA)
981 
982  real(RP), intent(out) :: qsat(IA,JA)
983 
984  integer :: i, j
985  !---------------------------------------------------------------------------
986 
987  !$omp parallel do OMP_SCHEDULE_
988  !$acc kernels copyin(temp,pres) copyout(qsat)
989  do j = js, je
990  do i = is, ie
991  call atmos_saturation_pres2qsat_all_0d( temp(i,j), pres(i,j), &
992  qsat(i,j) )
993  enddo
994  enddo
995  !$acc end kernels
996 
997  return
998  end subroutine atmos_saturation_pres2qsat_all_2d
999 
1000  !-----------------------------------------------------------------------------
1002  subroutine atmos_saturation_pres2qsat_all_qdry_2d( &
1003  IA, IS, IE, JA, JS, JE, &
1004  temp, pres, qdry, &
1005  qsat )
1006  implicit none
1007  integer, intent(in) :: IA, IS, IE
1008  integer, intent(in) :: JA, JS, JE
1009 
1010  real(RP), intent(in) :: temp(IA,JA)
1011  real(RP), intent(in) :: pres(IA,JA)
1012  real(RP), intent(in) :: qdry(IA,JA)
1013 
1014  real(RP), intent(out) :: qsat(IA,JA)
1015 
1016  integer :: i, j
1017  !---------------------------------------------------------------------------
1018 
1019  !$omp parallel do OMP_SCHEDULE_
1020  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
1021  do j = js, je
1022  do i = is, ie
1023  call atmos_saturation_pres2qsat_all_qdry_0d( temp(i,j), pres(i,j), qdry(i,j), &
1024  qsat(i,j) )
1025  enddo
1026  enddo
1027  !$acc end kernels
1028 
1029  return
1030  end subroutine atmos_saturation_pres2qsat_all_qdry_2d
1031 
1032  !-----------------------------------------------------------------------------
1034  subroutine atmos_saturation_pres2qsat_all_3d( &
1035  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1036  temp, pres, &
1037  qsat )
1038  implicit none
1039  integer, intent(in) :: KA, KS, KE
1040  integer, intent(in) :: IA, IS, IE
1041  integer, intent(in) :: JA, JS, JE
1042 
1043  real(RP), intent(in) :: temp(KA,IA,JA)
1044  real(RP), intent(in) :: pres(KA,IA,JA)
1045 
1046  real(RP), intent(out) :: qsat(KA,IA,JA)
1047 
1048  integer :: k, i, j
1049  !---------------------------------------------------------------------------
1050 
1051  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1052  !$acc kernels copyin(temp,pres) copyout(qsat)
1053  do j = js, je
1054  do i = is, ie
1055  do k = ks, ke
1056  call atmos_saturation_pres2qsat_all_0d( temp(k,i,j), pres(k,i,j), & ! [IN]
1057  qsat(k,i,j) ) ! [OUT]
1058  enddo
1059  enddo
1060  enddo
1061  !$acc end kernels
1062 
1063  return
1064  end subroutine atmos_saturation_pres2qsat_all_3d
1065 
1066  !-----------------------------------------------------------------------------
1068  subroutine atmos_saturation_pres2qsat_all_qdry_3d( &
1069  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1070  temp, pres, qdry, &
1071  qsat )
1072  implicit none
1073  integer, intent(in) :: KA, KS, KE
1074  integer, intent(in) :: IA, IS, IE
1075  integer, intent(in) :: JA, JS, JE
1076 
1077  real(RP), intent(in) :: temp(KA,IA,JA)
1078  real(RP), intent(in) :: pres(KA,IA,JA)
1079  real(RP), intent(in) :: qdry(KA,IA,JA)
1080 
1081  real(RP), intent(out) :: qsat(KA,IA,JA)
1082 
1083  integer :: k, i, j
1084  !---------------------------------------------------------------------------
1085 
1086  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1087  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
1088  do j = js, je
1089  do i = is, ie
1090  do k = ks, ke
1091  call atmos_saturation_pres2qsat_all_qdry_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1092  qsat(k,i,j) ) ! [OUT]
1093  enddo
1094  enddo
1095  enddo
1096  !$acc end kernels
1097 
1098  return
1099  end subroutine atmos_saturation_pres2qsat_all_qdry_3d
1100 
1101  !-----------------------------------------------------------------------------
1103  subroutine atmos_saturation_pres2qsat_liq_0d( &
1104  temp, pres, &
1105  qsat )
1106  !$acc routine
1107  implicit none
1108 
1109  real(RP), intent(in) :: temp
1110  real(RP), intent(in) :: pres
1111 
1112  real(RP), intent(out) :: qsat
1113 
1114  real(RP) :: psat
1115  !---------------------------------------------------------------------------
1116 
1117  call atmos_saturation_psat_liq( temp, psat ) ! [IN], [OUT]
1118  call atmos_saturation_psat2qsat_pres_0d( psat, pres, & ! [IN]
1119  qsat ) ! [OUT]
1120 
1121  return
1122  end subroutine atmos_saturation_pres2qsat_liq_0d
1123 
1124  !-----------------------------------------------------------------------------
1126  subroutine atmos_saturation_pres2qsat_liq_qdry_0d( &
1127  temp, pres, qdry, &
1128  qsat )
1129  !$acc routine
1130  implicit none
1131 
1132  real(RP), intent(in) :: temp
1133  real(RP), intent(in) :: pres
1134  real(RP), intent(in) :: qdry
1135 
1136  real(RP), intent(out) :: qsat
1137 
1138  real(RP) :: psat
1139  !---------------------------------------------------------------------------
1140 
1141  call atmos_saturation_psat_liq( temp, psat ) ! [IN], [OUT]
1142  call atmos_saturation_psat2qsat_pres_qdry_0d( psat, pres, qdry, & ! [IN]
1143  qsat ) ! [OUT]
1144 
1145  return
1146  end subroutine atmos_saturation_pres2qsat_liq_qdry_0d
1147 
1148  !-----------------------------------------------------------------------------
1150 !OCL SERIAL
1151  subroutine atmos_saturation_pres2qsat_liq_1d( &
1152  KA, KS, KE, &
1153  temp, pres, &
1154  qsat )
1155  !$acc routine vector
1156  implicit none
1157  integer, intent(in) :: KA, KS, KE
1158 
1159  real(RP), intent(in) :: temp(KA)
1160  real(RP), intent(in) :: pres(KA)
1161 
1162  real(RP), intent(out) :: qsat(KA)
1163 
1164  integer :: k
1165  !---------------------------------------------------------------------------
1166 
1167  do k = ks, ke
1168  call atmos_saturation_pres2qsat_liq_0d( temp(k), pres(k), & ! [IN]
1169  qsat(k) ) ! [OUT]
1170  enddo
1171 
1172  return
1173  end subroutine atmos_saturation_pres2qsat_liq_1d
1174 
1175  !-----------------------------------------------------------------------------
1177  subroutine atmos_saturation_pres2qsat_liq_1d_para( &
1178  KA, KS, KE, &
1179  temp, pres, &
1180  qsat )
1181  implicit none
1182  integer, intent(in) :: KA, KS, KE
1183 
1184  real(RP), intent(in) :: temp(KA)
1185  real(RP), intent(in) :: pres(KA)
1186 
1187  real(RP), intent(out) :: qsat(KA)
1188 
1189  integer :: k
1190  !---------------------------------------------------------------------------
1191 
1192  !$omp parallel do OMP_SCHEDULE_
1193  !$acc kernels copyin(temp,pres) copyout(qsat)
1194  do k = ks, ke
1195  call atmos_saturation_pres2qsat_liq_0d( temp(k), pres(k), & ! [IN]
1196  qsat(k) ) ! [OUT]
1197  enddo
1198  !$acc end kernels
1199 
1200  return
1201  end subroutine atmos_saturation_pres2qsat_liq_1d_para
1202 
1203  !-----------------------------------------------------------------------------
1205 !OCL SERIAL
1206  subroutine atmos_saturation_pres2qsat_liq_qdry_1d( &
1207  KA, KS, KE, &
1208  temp, pres, qdry, &
1209  qsat )
1210  !$acc routine vector
1211  implicit none
1212  integer, intent(in) :: KA, KS, KE
1213 
1214  real(RP), intent(in) :: temp(KA)
1215  real(RP), intent(in) :: pres(KA)
1216  real(RP), intent(in) :: qdry(KA)
1217 
1218  real(RP), intent(out) :: qsat(KA)
1219 
1220  integer :: k
1221  !---------------------------------------------------------------------------
1222 
1223  do k = ks, ke
1224  call atmos_saturation_pres2qsat_liq_qdry_0d( temp(k), pres(k), qdry(k), & ! [IN]
1225  qsat(k) ) ! [OUT]
1226  enddo
1227 
1228  return
1229  end subroutine atmos_saturation_pres2qsat_liq_qdry_1d
1230 
1231  !-----------------------------------------------------------------------------
1233  subroutine atmos_saturation_pres2qsat_liq_qdry_1d_para( &
1234  KA, KS, KE, &
1235  temp, pres, qdry, &
1236  qsat )
1237  implicit none
1238  integer, intent(in) :: KA, KS, KE
1239 
1240  real(RP), intent(in) :: temp(KA)
1241  real(RP), intent(in) :: pres(KA)
1242  real(RP), intent(in) :: qdry(KA)
1243 
1244  real(RP), intent(out) :: qsat(KA)
1245 
1246  integer :: k
1247  !---------------------------------------------------------------------------
1248 
1249  !$omp parallel do OMP_SCHEDULE_
1250  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
1251  do k = ks, ke
1252  call atmos_saturation_pres2qsat_liq_qdry_0d( temp(k), pres(k), qdry(k), & ! [IN]
1253  qsat(k) ) ! [OUT]
1254  enddo
1255  !$acc end kernels
1256 
1257  return
1258  end subroutine atmos_saturation_pres2qsat_liq_qdry_1d_para
1259 
1260  !-----------------------------------------------------------------------------
1262  subroutine atmos_saturation_pres2qsat_liq_3d( &
1263  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1264  temp, pres, &
1265  qsat )
1266  implicit none
1267  integer, intent(in) :: KA, KS, KE
1268  integer, intent(in) :: IA, IS, IE
1269  integer, intent(in) :: JA, JS, JE
1270 
1271  real(RP), intent(in) :: temp(KA,IA,JA)
1272  real(RP), intent(in) :: pres(KA,IA,JA)
1273 
1274  real(RP), intent(out) :: qsat(KA,IA,JA)
1275 
1276  integer :: k, i, j
1277  !---------------------------------------------------------------------------
1278 
1279  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1280  !$acc kernels copyin(temp,pres) copyout(qsat)
1281  do j = js, je
1282  do i = is, ie
1283  do k = ks, ke
1284  call atmos_saturation_pres2qsat_liq_0d( temp(k,i,j), pres(k,i,j), & ! [IN]
1285  qsat(k,i,j) ) ! [OUT]
1286  enddo
1287  enddo
1288  enddo
1289  !$acc end kernels
1290 
1291  return
1292  end subroutine atmos_saturation_pres2qsat_liq_3d
1293 
1294  !-----------------------------------------------------------------------------
1296  subroutine atmos_saturation_pres2qsat_liq_qdry_3d( &
1297  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1298  temp, pres, qdry, &
1299  qsat )
1300  implicit none
1301  integer, intent(in) :: KA, KS, KE
1302  integer, intent(in) :: IA, IS, IE
1303  integer, intent(in) :: JA, JS, JE
1304 
1305  real(RP), intent(in) :: temp(KA,IA,JA)
1306  real(RP), intent(in) :: pres(KA,IA,JA)
1307  real(RP), intent(in) :: qdry(KA,IA,JA)
1308 
1309  real(RP), intent(out) :: qsat(KA,IA,JA)
1310 
1311  integer :: k, i, j
1312  !---------------------------------------------------------------------------
1313 
1314  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1315  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
1316  do j = js, je
1317  do i = is, ie
1318  do k = ks, ke
1319  call atmos_saturation_pres2qsat_liq_qdry_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1320  qsat(k,i,j) ) ! [OUT]
1321  enddo
1322  enddo
1323  enddo
1324  !$acc end kernels
1325 
1326  return
1327  end subroutine atmos_saturation_pres2qsat_liq_qdry_3d
1328 
1329  !-----------------------------------------------------------------------------
1331  subroutine atmos_saturation_pres2qsat_ice_0d( &
1332  temp, pres, &
1333  qsat )
1334  !$acc routine
1335  implicit none
1336 
1337  real(RP), intent(in) :: temp
1338  real(RP), intent(in) :: pres
1339 
1340  real(RP), intent(out) :: qsat
1341 
1342  real(RP) :: psat
1343  !---------------------------------------------------------------------------
1344 
1345  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
1346  call atmos_saturation_psat2qsat_pres_0d( psat, pres, & ! [IN]
1347  qsat ) ! [OUT]
1348 
1349  return
1350  end subroutine atmos_saturation_pres2qsat_ice_0d
1351 
1352  !-----------------------------------------------------------------------------
1354  subroutine atmos_saturation_pres2qsat_ice_qdry_0d( &
1355  temp, pres, qdry, &
1356  qsat )
1357  !$acc routine
1358  implicit none
1359 
1360  real(RP), intent(in) :: temp
1361  real(RP), intent(in) :: pres
1362  real(RP), intent(in) :: qdry
1363 
1364  real(RP), intent(out) :: qsat
1365 
1366  real(RP) :: psat
1367  !---------------------------------------------------------------------------
1368 
1369  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
1370  call atmos_saturation_psat2qsat_pres_qdry_0d( psat, pres, qdry, & ! [IN]
1371  qsat ) ! [OUT]
1372 
1373  return
1374  end subroutine atmos_saturation_pres2qsat_ice_qdry_0d
1375 
1376  !-----------------------------------------------------------------------------
1378 !OCL SERIAL
1379  subroutine atmos_saturation_pres2qsat_ice_1d( &
1380  KA, KS, KE, &
1381  temp, pres, &
1382  qsat )
1383  !$acc routine vector
1384  implicit none
1385  integer, intent(in) :: KA, KS, KE
1386 
1387  real(RP), intent(in) :: temp(KA)
1388  real(RP), intent(in) :: pres(KA)
1389 
1390  real(RP), intent(out) :: qsat(KA)
1391 
1392  integer :: k
1393  !---------------------------------------------------------------------------
1394 
1395  do k = ks, ke
1396  call atmos_saturation_pres2qsat_ice_0d( temp(k), pres(k), & ! [IN]
1397  qsat(k) ) ! [OUT]
1398  enddo
1399 
1400  return
1401  end subroutine atmos_saturation_pres2qsat_ice_1d
1402 
1403  !-----------------------------------------------------------------------------
1405  subroutine atmos_saturation_pres2qsat_ice_1d_para( &
1406  KA, KS, KE, &
1407  temp, pres, &
1408  qsat )
1409  implicit none
1410  integer, intent(in) :: KA, KS, KE
1411 
1412  real(RP), intent(in) :: temp(KA)
1413  real(RP), intent(in) :: pres(KA)
1414 
1415  real(RP), intent(out) :: qsat(KA)
1416 
1417  integer :: k
1418  !---------------------------------------------------------------------------
1419 
1420  !$omp parallel do OMP_SCHEDULE_
1421  !$acc kernels copyin(temp,pres) copyout(qsat)
1422  do k = ks, ke
1423  call atmos_saturation_pres2qsat_ice_0d( temp(k), pres(k), & ! [IN]
1424  qsat(k) ) ! [OUT]
1425  enddo
1426  !$acc end kernels
1427 
1428  return
1429  end subroutine atmos_saturation_pres2qsat_ice_1d_para
1430 
1431  !-----------------------------------------------------------------------------
1433 !OCL SERIAL
1434  subroutine atmos_saturation_pres2qsat_ice_qdry_1d( &
1435  KA, KS, KE, &
1436  temp, pres, qdry, &
1437  qsat )
1438  !$acc routine vector
1439  implicit none
1440  integer, intent(in) :: KA, KS, KE
1441 
1442  real(RP), intent(in) :: temp(KA)
1443  real(RP), intent(in) :: pres(KA)
1444  real(RP), intent(in) :: qdry(KA)
1445 
1446  real(RP), intent(out) :: qsat(KA)
1447 
1448  integer :: k
1449  !---------------------------------------------------------------------------
1450 
1451  do k = ks, ke
1452  call atmos_saturation_pres2qsat_ice_qdry_0d( temp(k), pres(k), qdry(k), & ! [IN]
1453  qsat(k) ) ! [OUT]
1454  enddo
1455 
1456  return
1457  end subroutine atmos_saturation_pres2qsat_ice_qdry_1d
1458 
1459  !-----------------------------------------------------------------------------
1461  subroutine atmos_saturation_pres2qsat_ice_qdry_1d_para( &
1462  KA, KS, KE, &
1463  temp, pres, qdry, &
1464  qsat )
1465  implicit none
1466  integer, intent(in) :: KA, KS, KE
1467 
1468  real(RP), intent(in) :: temp(KA)
1469  real(RP), intent(in) :: pres(KA)
1470  real(RP), intent(in) :: qdry(KA)
1471 
1472  real(RP), intent(out) :: qsat(KA)
1473 
1474  integer :: k
1475  !---------------------------------------------------------------------------
1476 
1477  !$omp parallel do OMP_SCHEDULE_
1478  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
1479  do k = ks, ke
1480  call atmos_saturation_pres2qsat_ice_qdry_0d( temp(k), pres(k), qdry(k), & ! [IN]
1481  qsat(k) ) ! [OUT]
1482  enddo
1483  !$acc end kernels
1484 
1485  return
1486  end subroutine atmos_saturation_pres2qsat_ice_qdry_1d_para
1487 
1488  !-----------------------------------------------------------------------------
1490  subroutine atmos_saturation_pres2qsat_ice_3d( &
1491  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1492  temp, pres, &
1493  qsat )
1494  implicit none
1495  integer, intent(in) :: KA, KS, KE
1496  integer, intent(in) :: IA, IS, IE
1497  integer, intent(in) :: JA, JS, JE
1498 
1499  real(RP), intent(in) :: temp(KA,IA,JA)
1500  real(RP), intent(in) :: pres(KA,IA,JA)
1501 
1502  real(RP), intent(out) :: qsat(KA,IA,JA)
1503 
1504  integer :: k, i, j
1505  !---------------------------------------------------------------------------
1506 
1507  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1508  !$acc kernels copyin(temp,pres) copyout(qsat)
1509  do j = js, je
1510  do i = is, ie
1511  do k = ks, ke
1512  call atmos_saturation_pres2qsat_ice_0d( temp(k,i,j), pres(k,i,j), & ! [IN]
1513  qsat(k,i,j) ) ! [OUT]
1514  enddo
1515  enddo
1516  enddo
1517  !$acc end kernels
1518 
1519  return
1520  end subroutine atmos_saturation_pres2qsat_ice_3d
1521 
1522  !-----------------------------------------------------------------------------
1524  subroutine atmos_saturation_pres2qsat_ice_qdry_3d( &
1525  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1526  temp, pres, qdry, &
1527  qsat )
1528  implicit none
1529  integer, intent(in) :: KA, KS, KE
1530  integer, intent(in) :: IA, IS, IE
1531  integer, intent(in) :: JA, JS, JE
1532 
1533  real(RP), intent(in) :: temp(KA,IA,JA)
1534  real(RP), intent(in) :: pres(KA,IA,JA)
1535  real(RP), intent(in) :: qdry(KA,IA,JA)
1536 
1537  real(RP), intent(out) :: qsat(KA,IA,JA)
1538 
1539  integer :: k, i, j
1540  !---------------------------------------------------------------------------
1541 
1542  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1543  !$acc kernels copyin(temp,pres,qdry) copyout(qsat)
1544  do j = js, je
1545  do i = is, ie
1546  do k = ks, ke
1547  call atmos_saturation_pres2qsat_ice_qdry_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1548  qsat(k,i,j) ) ! [OUT]
1549  enddo
1550  enddo
1551  enddo
1552  !$acc end kernels
1553 
1554  return
1555  end subroutine atmos_saturation_pres2qsat_ice_qdry_3d
1556 
1557  !-----------------------------------------------------------------------------
1559  subroutine atmos_saturation_psat2qsat_dens_0d( &
1560  psat, temp, dens, &
1561  qsat )
1562  !$acc routine
1563  implicit none
1564  real(RP), intent(in) :: psat
1565  real(RP), intent(in) :: temp
1566  real(RP), intent(in) :: dens
1567 
1568  real(RP), intent(out) :: qsat
1569 
1570  qsat = psat / ( dens * rvap * temp )
1571 
1572  return
1573  end subroutine atmos_saturation_psat2qsat_dens_0d
1574 
1575  !-----------------------------------------------------------------------------
1577  subroutine atmos_saturation_dens2qsat_all_0d( &
1578  temp, &
1579  dens, &
1580  qsat )
1581  !$acc routine
1582  implicit none
1583  real(RP), intent(in) :: temp
1584  real(RP), intent(in) :: dens
1585 
1586  real(RP), intent(out) :: qsat
1587 
1588  real(RP) :: psat
1589  !---------------------------------------------------------------------------
1590 
1591  call atmos_saturation_psat_all( temp, psat ) ! [IN], [OUT]
1592  call atmos_saturation_psat2qsat_dens_0d( psat, temp, dens, & ! [IN]
1593  qsat ) ! [OUT]
1594 
1595  return
1596  end subroutine atmos_saturation_dens2qsat_all_0d
1597 
1598  !-----------------------------------------------------------------------------
1600 !OCL SERIAL
1601  subroutine atmos_saturation_dens2qsat_all_1d( &
1602  KA, KS, KE, &
1603  temp, dens, &
1604  qsat )
1605  !$acc routine vector
1606  implicit none
1607  integer, intent(in) :: KA, KS, KE
1608 
1609  real(RP), intent(in) :: temp(KA)
1610  real(RP), intent(in) :: dens(KA)
1611 
1612  real(RP), intent(out) :: qsat(KA)
1613 
1614  integer :: k
1615  !---------------------------------------------------------------------------
1616 
1617  do k = ks, ke
1618  call atmos_saturation_dens2qsat_all_0d( temp(k), dens(k), & ! [IN]
1619  qsat(k) ) ! [OUT]
1620  enddo
1621 
1622  return
1623  end subroutine atmos_saturation_dens2qsat_all_1d
1624 
1625  !-----------------------------------------------------------------------------
1627  subroutine atmos_saturation_dens2qsat_all_3d( &
1628  KA, KS, KE, &
1629  IA, IS, IE, &
1630  JA, JS, JE, &
1631  temp, &
1632  dens, &
1633  qsat )
1634  implicit none
1635  integer, intent(in) :: KA, KS, KE
1636  integer, intent(in) :: IA, IS, IE
1637  integer, intent(in) :: JA, JS, JE
1638 
1639  real(RP), intent(in) :: temp(KA,IA,JA)
1640  real(RP), intent(in) :: dens(KA,IA,JA)
1641 
1642  real(RP), intent(out) :: qsat(KA,IA,JA)
1643 
1644  integer :: k, i, j
1645  !---------------------------------------------------------------------------
1646 
1647  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1648  !$acc kernels copyin(temp,dens) copyout(qsat)
1649  do j = js, je
1650  do i = is, ie
1651  do k = ks, ke
1652  call atmos_saturation_dens2qsat_all_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1653  qsat(k,i,j) ) ! [OUT]
1654  enddo
1655  enddo
1656  enddo
1657  !$acc end kernels
1658 
1659  return
1660  end subroutine atmos_saturation_dens2qsat_all_3d
1661 
1662  !-----------------------------------------------------------------------------
1664  subroutine atmos_saturation_dens2qsat_liq_0d( &
1665  temp, &
1666  dens, &
1667  qsat )
1668  !$acc routine
1669  implicit none
1670  real(RP), intent(in) :: temp
1671  real(RP), intent(in) :: dens
1672 
1673  real(RP), intent(out) :: qsat
1674 
1675  real(RP) :: psat
1676  !---------------------------------------------------------------------------
1677 
1678  call atmos_saturation_psat_liq_0d( temp, psat ) ! [IN], [OUT]
1679  call atmos_saturation_psat2qsat_dens_0d( psat, temp, dens, & ! [IN]
1680  qsat ) ! [OUT]
1681 
1682  return
1683  end subroutine atmos_saturation_dens2qsat_liq_0d
1684 
1685  !-----------------------------------------------------------------------------
1687 !OCL SERIAL
1688  subroutine atmos_saturation_dens2qsat_liq_1d( &
1689  KA, KS, KE, &
1690  temp, dens, &
1691  qsat )
1692  !$acc routine vector
1693  implicit none
1694  integer, intent(in) :: KA, KS, KE
1695 
1696  real(RP), intent(in) :: temp(KA)
1697  real(RP), intent(in) :: dens(KA)
1698 
1699  real(RP), intent(out) :: qsat(KA)
1700 
1701  integer :: k
1702  !---------------------------------------------------------------------------
1703 
1704  do k = ks, ke
1705  call atmos_saturation_dens2qsat_liq_0d( temp(k), dens(k), & ! [IN]
1706  qsat(k) ) ! [OUT]
1707  enddo
1708 
1709  return
1710  end subroutine atmos_saturation_dens2qsat_liq_1d
1711 
1712  !-----------------------------------------------------------------------------
1714  subroutine atmos_saturation_dens2qsat_liq_3d( &
1715  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1716  temp, dens, &
1717  qsat )
1718  implicit none
1719  integer, intent(in) :: KA, KS, KE
1720  integer, intent(in) :: IA, IS, IE
1721  integer, intent(in) :: JA, JS, JE
1722 
1723  real(RP), intent(in) :: temp(KA,IA,JA)
1724  real(RP), intent(in) :: dens(KA,IA,JA)
1725 
1726  real(RP), intent(out) :: qsat(KA,IA,JA)
1727 
1728  integer :: k, i, j
1729  !---------------------------------------------------------------------------
1730 
1731  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1732  !$acc kernels copyin(temp,dens) copyout(qsat)
1733  do j = js, je
1734  do i = is, ie
1735  do k = ks, ke
1736  call atmos_saturation_dens2qsat_liq_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1737  qsat(k,i,j) ) ! [OUT]
1738  enddo
1739  enddo
1740  enddo
1741  !$acc end kernels
1742 
1743  return
1744  end subroutine atmos_saturation_dens2qsat_liq_3d
1745 
1746  !-----------------------------------------------------------------------------
1748  subroutine atmos_saturation_dens2qsat_ice_0d( &
1749  temp, &
1750  dens, &
1751  qsat )
1752  !$acc routine
1753  implicit none
1754  real(RP), intent(in) :: temp
1755  real(RP), intent(in) :: dens
1756 
1757  real(RP), intent(out) :: qsat
1758 
1759  real(RP) :: psat
1760  !---------------------------------------------------------------------------
1761 
1762  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
1763  call atmos_saturation_psat2qsat_dens( psat, temp, dens, & ! [IN]
1764  qsat ) ! [OUT]
1765  return
1766  end subroutine atmos_saturation_dens2qsat_ice_0d
1767 
1768  !-----------------------------------------------------------------------------
1770 !OCL SERIAL
1771  subroutine atmos_saturation_dens2qsat_ice_1d( &
1772  KA, KS, KE, &
1773  temp, dens, &
1774  qsat )
1775  !$acc routine vector
1776  implicit none
1777  integer, intent(in) :: KA, KS, KE
1778 
1779  real(RP), intent(in) :: temp(KA)
1780  real(RP), intent(in) :: dens(KA)
1781 
1782  real(RP), intent(out) :: qsat(KA)
1783 
1784  integer :: k
1785  !---------------------------------------------------------------------------
1786 
1787  do k = ks, ke
1788  call atmos_saturation_dens2qsat_ice_0d( temp(k), dens(k), & ! [IN]
1789  qsat(k) ) ! [OUT]
1790  enddo
1791 
1792  return
1793  end subroutine atmos_saturation_dens2qsat_ice_1d
1794 
1795  !-----------------------------------------------------------------------------
1797  subroutine atmos_saturation_dens2qsat_ice_3d( &
1798  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1799  temp, dens, &
1800  qsat )
1801  implicit none
1802  integer, intent(in) :: KA, KS, KE
1803  integer, intent(in) :: IA, IS, IE
1804  integer, intent(in) :: JA, JS, JE
1805 
1806  real(RP), intent(in) :: temp(KA,IA,JA)
1807  real(RP), intent(in) :: dens(KA,IA,JA)
1808 
1809  real(RP), intent(out) :: qsat(KA,IA,JA)
1810 
1811  integer :: k, i, j
1812  !---------------------------------------------------------------------------
1813 
1814  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1815  !$acc kernels copyin(temp,dens) copyout(qsat)
1816  do j = js, je
1817  do i = is, ie
1818  do k = ks, ke
1819  call atmos_saturation_dens2qsat_ice_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1820  qsat(k,i,j) ) ! [OUT]
1821  enddo
1822  enddo
1823  enddo
1824  !$acc end kernels
1825 
1826  return
1827  end subroutine atmos_saturation_dens2qsat_ice_3d
1828 
1829  !-----------------------------------------------------------------------------
1831  subroutine atmos_saturation_dalphadt_0d( &
1832  temp, &
1833  dalpha_dT )
1834  !$acc routine
1835  implicit none
1836  real(RP), intent(in) :: temp
1837 
1838  real(RP), intent(out) :: dalpha_dT
1839 
1840  real(RP) :: lim1, lim2
1841  !---------------------------------------------------------------------------
1842 
1843  ! if Tup < temp, dalpha/dT = 0 (no slope)
1844  lim1 = 0.5_rp + sign( 0.5_rp, atmos_saturation_ulimit_temp - temp )
1845  ! if Tdn > temp, dalpha/dT = 0 (no slope)
1846  lim2 = 0.5_rp + sign( 0.5_rp, temp - atmos_saturation_llimit_temp )
1847 
1848  dalpha_dt = dalphadt_const * lim1 * lim2
1849 
1850  return
1851  end subroutine atmos_saturation_dalphadt_0d
1852 
1853  !-----------------------------------------------------------------------------
1855 !OCL SERIAL
1856  subroutine atmos_saturation_dalphadt_1d( &
1857  KA, KS, KE, &
1858  temp, &
1859  dalpha_dT )
1860  !$acc routine vector
1861  implicit none
1862  integer, intent(in) :: KA, KS, KE
1863 
1864  real(RP), intent(in) :: temp (KA)
1865 
1866  real(RP), intent(out) :: dalpha_dT(KA)
1867 
1868  integer :: k
1869  !---------------------------------------------------------------------------
1870 
1871  do k = ks, ke
1872  call atmos_saturation_dalphadt_0d( temp(k), dalpha_dt(k) ) ! [IN], [OUT]
1873  enddo
1874 
1875  return
1876  end subroutine atmos_saturation_dalphadt_1d
1877 
1878  !-----------------------------------------------------------------------------
1880  subroutine atmos_saturation_dalphadt_3d( &
1881  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1882  temp, &
1883  dalpha_dT )
1884  implicit none
1885  integer, intent(in) :: KA, KS, KE
1886  integer, intent(in) :: IA, IS, IE
1887  integer, intent(in) :: JA, JS, JE
1888 
1889  real(RP), intent(in) :: temp (KA,IA,JA)
1890 
1891  real(RP), intent(out) :: dalpha_dT(KA,IA,JA)
1892 
1893  integer :: k, i, j
1894  !---------------------------------------------------------------------------
1895 
1896  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1897  !$acc kernels copyin(temp) copyout(dalpha_dT)
1898  do j = js, je
1899  do i = is, ie
1900  do k = ks, ke
1901  call atmos_saturation_dalphadt_0d( temp(k,i,j), dalpha_dt(k,i,j) ) ! [IN], [OUT]
1902  enddo
1903  enddo
1904  enddo
1905  !$acc end kernels
1906 
1907  return
1908  end subroutine atmos_saturation_dalphadt_3d
1909 
1910  !-----------------------------------------------------------------------------
1911  ! (d qsw/d T)_{rho}: partial difference of qsat_water at constant density
1912  subroutine atmos_saturation_dqs_dtem_dens_liq_0d( &
1913  temp, dens, &
1914  dqsdtem, qsat )
1915  !$acc routine
1916  use scale_atmos_hydrometeor, only: &
1917  hydrometeor_lhv => atmos_hydrometeor_lhv
1918  implicit none
1919 
1920  real(RP), intent(in) :: temp
1921  real(RP), intent(in) :: dens
1922 
1923  real(RP), intent(out) :: dqsdtem
1924  real(RP), intent(out), optional :: qsat
1925 
1926  real(RP) :: LHV ! latent heat of vaporization [J/kg]
1927  real(RP) :: psat
1928  !---------------------------------------------------------------------------
1929 
1930  call hydrometeor_lhv( temp, lhv )
1931  call atmos_saturation_psat_liq( temp, psat ) ! [IN], [OUT]
1932 
1933  dqsdtem = psat / ( dens* rvap * temp**2 ) &
1934  * ( lhv / ( rvap * temp ) - 1.0_rp )
1935 
1936  if ( present(qsat) ) &
1937  call atmos_saturation_psat2qsat_dens( psat, temp, dens, & ! [IN]
1938  qsat ) ! [OUT]
1939 
1940  return
1941  end subroutine atmos_saturation_dqs_dtem_dens_liq_0d
1942 
1943  !-----------------------------------------------------------------------------
1944 !OCL SERIAL
1946  KA, KS, KE, &
1947  temp, dens, &
1948  dqsdtem )
1949  !$acc routine vector
1950  implicit none
1951  integer, intent(in) :: KA, KS, KE
1952 
1953  real(RP), intent(in) :: temp (KA)
1954  real(RP), intent(in) :: dens (KA)
1955 
1956  real(RP), intent(out) :: dqsdtem(KA)
1957 
1958  integer :: k
1959  !---------------------------------------------------------------------------
1960 
1961  do k = ks, ke
1962  call atmos_saturation_dqs_dtem_dens_liq_0d( temp(k), dens(k), & ! [IN]
1963  dqsdtem(k) ) ! [OUT]
1964  enddo
1965 
1966  return
1968 
1969  !-----------------------------------------------------------------------------
1970  ! (d qsw/d T)_{rho}: partial difference of qsat_water at constant density
1971  subroutine atmos_saturation_dqs_dtem_dens_liq_3d( &
1972  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1973  temp, dens, &
1974  dqsdtem )
1975  implicit none
1976  integer, intent(in) :: KA, KS, KE
1977  integer, intent(in) :: IA, IS, IE
1978  integer, intent(in) :: JA, JS, JE
1979 
1980  real(RP), intent(in) :: temp (KA,IA,JA)
1981  real(RP), intent(in) :: dens (KA,IA,JA)
1982 
1983  real(RP), intent(out) :: dqsdtem(KA,IA,JA)
1984 
1985  integer :: k, i, j
1986  !---------------------------------------------------------------------------
1987 
1988  !$omp parallel do OMP_SCHEDULE_ collapse(2)
1989  !$acc kernels copyin(temp,dens) copyout(dqsdtem)
1990  do j = js, je
1991  do i = is, ie
1992  do k = ks, ke
1993  call atmos_saturation_dqs_dtem_dens_liq_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
1994  dqsdtem(k,i,j) ) ! [OUT]
1995  enddo
1996  enddo
1997  enddo
1998  !$acc end kernels
1999 
2000  return
2001  end subroutine atmos_saturation_dqs_dtem_dens_liq_3d
2002 
2003  !-----------------------------------------------------------------------------
2004  ! (d qsi/d T)_{rho}: partial difference of qsat_ice at constant density
2005  subroutine atmos_saturation_dqs_dtem_dens_ice_0d( &
2006  temp, dens, &
2007  dqsdtem, qsat )
2008  !$acc routine
2009  use scale_atmos_hydrometeor, only: &
2010  hydrometeor_lhs => atmos_hydrometeor_lhs
2011  implicit none
2012 
2013  real(RP), intent(in) :: temp
2014  real(RP), intent(in) :: dens
2015 
2016  real(RP), intent(out) :: dqsdtem
2017  real(RP), intent(out), optional :: qsat
2018 
2019  real(RP) :: LHS ! latent heat of sublimation [J/kg]
2020  real(RP) :: psat
2021  !---------------------------------------------------------------------------
2022 
2023  call hydrometeor_lhs( temp, lhs ) ! [IN], [OUT]
2024  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
2025 
2026  dqsdtem = psat / ( dens * rvap * temp**2 ) &
2027  * ( lhs / ( rvap * temp ) - 1.0_rp )
2028 
2029  if ( present(qsat) ) &
2030  call atmos_saturation_psat2qsat_dens( psat, temp, dens, & ! [IN]
2031  qsat ) ! [OUT]
2032 
2033  return
2034  end subroutine atmos_saturation_dqs_dtem_dens_ice_0d
2035 
2036  !-----------------------------------------------------------------------------
2037 !OCL SERIAL
2039  KA, KS, KE, &
2040  temp, dens, &
2041  dqsdtem )
2042  !$acc routine vector
2043  implicit none
2044  integer, intent(in) :: KA, KS, KE
2045 
2046  real(RP), intent(in) :: temp (KA)
2047  real(RP), intent(in) :: dens (KA)
2048 
2049  real(RP), intent(out) :: dqsdtem(KA)
2050 
2051  integer :: k
2052  !---------------------------------------------------------------------------
2053 
2054  do k = ks, ke
2055  call atmos_saturation_dqs_dtem_dens_ice_0d( temp(k), dens(k), & ! [IN]
2056  dqsdtem(k) ) ! [OUT]
2057  enddo
2058 
2059  return
2061 
2062  !-----------------------------------------------------------------------------
2063  ! (d qsi/d T)_{rho}: partial difference of qsat_ice at constant density
2064  subroutine atmos_saturation_dqs_dtem_dens_ice_3d( &
2065  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2066  temp, dens, &
2067  dqsdtem )
2068  implicit none
2069  integer, intent(in) :: KA, KS, KE
2070  integer, intent(in) :: IA, IS, IE
2071  integer, intent(in) :: JA, JS, JE
2072 
2073  real(RP), intent(in) :: temp (KA,IA,JA)
2074  real(RP), intent(in) :: dens (KA,IA,JA)
2075 
2076  real(RP), intent(out) :: dqsdtem(KA,IA,JA)
2077 
2078  integer :: k, i, j
2079  !---------------------------------------------------------------------------
2080 
2081  !$omp parallel do OMP_SCHEDULE_ collapse(2)
2082  !$acc kernels copyin(temp,dens) copyout(dqsdtem)
2083  do j = js, je
2084  do i = is, ie
2085  do k = ks, ke
2086  call atmos_saturation_dqs_dtem_dens_ice_0d( temp(k,i,j), dens(k,i,j), & ! [IN]
2087  dqsdtem(k,i,j) ) ! [OUT]
2088  enddo
2089  enddo
2090  enddo
2091  !$acc end kernels
2092 
2093  return
2094  end subroutine atmos_saturation_dqs_dtem_dens_ice_3d
2095 
2096  !-----------------------------------------------------------------------------
2097  ! (d qsi/d T)_{rho}: partial difference of qsat_all at constant density
2098  subroutine atmos_saturation_dqs_dtem_dens_all_0d( &
2099  temp, dens, &
2100  dqsat_dT, &
2101  qsat, qsat_liq, qsat_ice, &
2102  alpha )
2103  !$acc routine
2104  implicit none
2105 
2106  real(RP), intent(in) :: temp
2107  real(RP), intent(in) :: dens
2108 
2109  real(RP), intent(out) :: dqsat_dT
2110 
2111  real(RP), intent(out), optional :: qsat
2112  real(RP), intent(out), optional :: qsat_liq
2113  real(RP), intent(out), optional :: qsat_ice
2114  real(RP), intent(out), optional :: alpha
2115 
2116  real(RP) :: qsat_liq_, qsat_ice_, alpha_
2117  real(RP) :: dqsat_dT_liq, dqsat_dT_ice, dalpha_dT
2118 
2119  !---------------------------------------------------------------------------
2120 
2121  call atmos_saturation_dqs_dtem_dens_liq_0d( temp, dens, & ! [IN]
2122  dqsat_dt_liq, qsat_liq_ ) ! [OUT]
2123  call atmos_saturation_dqs_dtem_dens_ice_0d( temp, dens, & ! [IN]
2124  dqsat_dt_ice, qsat_ice_ ) ! [OUT]
2125  call atmos_saturation_alpha ( temp, alpha_ ) ! [IN], [OUT]
2126  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
2127 
2128  dqsat_dt = qsat_liq_ * dalpha_dt + dqsat_dt_liq * ( alpha_ ) &
2129  - qsat_ice_ * dalpha_dt + dqsat_dt_ice * ( 1.0_rp-alpha_ )
2130 
2131  if ( present(qsat) ) qsat = qsat_liq_ * alpha_ + qsat_ice_ * ( 1.0_rp-alpha_ )
2132  if ( present(qsat_liq) ) qsat_liq = qsat_liq_
2133  if ( present(qsat_ice) ) qsat_ice = qsat_ice_
2134  if ( present(alpha ) ) alpha = alpha_
2135 
2136  return
2137  end subroutine atmos_saturation_dqs_dtem_dens_all_0d
2138 
2139  !-----------------------------------------------------------------------------
2140  ! (d qs/d T)_{p} and (d qs/d p)_{T}
2141  subroutine atmos_saturation_dqs_dtem_dpre_liq_0d( &
2142  temp, pres, qdry, &
2143  dqsat_dT, dqsat_dP, &
2144  qsat, psat )
2145  !$acc routine
2146  use scale_atmos_hydrometeor, only: &
2147  hydrometeor_lhv => atmos_hydrometeor_lhv
2148  implicit none
2149 
2150  real(RP), intent(in) :: temp
2151  real(RP), intent(in) :: pres
2152  real(RP), intent(in) :: qdry
2153 
2154  real(RP), intent(out) :: dqsat_dT
2155  real(RP), intent(out) :: dqsat_dP
2156 
2157  real(RP), intent(out), optional :: qsat
2158  real(RP), intent(out), optional :: psat
2159 
2160  real(RP) :: LHV ! latent heat of vaporization [J/kg]
2161  real(RP) :: psat_
2162  real(RP) :: den1, den2 ! denominator
2163 
2164  !---------------------------------------------------------------------------
2165 
2166  call hydrometeor_lhv( temp, lhv )
2167  call atmos_saturation_psat_liq( temp, psat_ ) ! [IN], [OUT]
2168 
2169  den1 = ( pres - (1.0_rp-epsvap) * psat_ )**2
2170  den2 = den1 * rvap * temp**2
2171 
2172  dqsat_dp = - epsvap * psat_ / den1
2173  dqsat_dt = epsvap * psat_ / den2 * lhv * pres
2174 
2175  if ( present(qsat) ) &
2176  call atmos_saturation_psat2qsat_pres_qdry_0d( psat_, pres, qdry, & ! [IN]
2177  qsat ) ! [OUT]
2178  if ( present(psat) ) psat = psat_
2179 
2180  return
2181  end subroutine atmos_saturation_dqs_dtem_dpre_liq_0d
2182 
2183  !-----------------------------------------------------------------------------
2184 !OCL SERIAL
2186  KA, KS, KE, &
2187  temp, pres, qdry, &
2188  dqsat_dT, dqsat_dP )
2189  !$acc routine vector
2190  implicit none
2191  integer, intent(in) :: KA, KS, KE
2192 
2193  real(RP), intent(in) :: temp (KA)
2194  real(RP), intent(in) :: pres (KA)
2195  real(RP), intent(in) :: qdry (KA)
2196 
2197  real(RP), intent(out) :: dqsat_dT(KA)
2198  real(RP), intent(out) :: dqsat_dP(KA)
2199 
2200  integer :: k
2201  !---------------------------------------------------------------------------
2202 
2203  do k = ks, ke
2204  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp(k), pres(k), qdry(k), & ! [IN]
2205  dqsat_dt(k), dqsat_dp(k) ) ! [OUT]
2206  enddo
2207 
2208  return
2210 
2211  !-----------------------------------------------------------------------------
2212  ! (d qs/d T)_{p} and (d qs/d p)_{T}
2213  subroutine atmos_saturation_dqs_dtem_dpre_liq_3d( &
2214  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2215  temp, pres, qdry, &
2216  dqsat_dT, dqsat_dP )
2217  implicit none
2218  integer, intent(in) :: KA, KS, KE
2219  integer, intent(in) :: IA, IS, IE
2220  integer, intent(in) :: JA, JS, JE
2221 
2222  real(RP), intent(in) :: temp (KA,IA,JA)
2223  real(RP), intent(in) :: pres (KA,IA,JA)
2224  real(RP), intent(in) :: qdry (KA,IA,JA)
2225 
2226  real(RP), intent(out) :: dqsat_dT(KA,IA,JA)
2227  real(RP), intent(out) :: dqsat_dP(KA,IA,JA)
2228 
2229  integer :: k, i, j
2230  !---------------------------------------------------------------------------
2231 
2232  !$omp parallel do OMP_SCHEDULE_ collapse(2)
2233  !$acc kernels copyin(temp,pres,qdry) copyout(dqsat_dT,dqsat_dP)
2234  do j = js, je
2235  do i = is, ie
2236  do k = ks, ke
2237  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
2238  dqsat_dt(k,i,j), dqsat_dp(k,i,j) ) ! [OUT]
2239  enddo
2240  enddo
2241  enddo
2242  !$acc end kernels
2243 
2244  return
2245  end subroutine atmos_saturation_dqs_dtem_dpre_liq_3d
2246 
2247  !-----------------------------------------------------------------------------
2248  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
2249  !-----------------------------------------------------------------------------
2250  subroutine atmos_saturation_dqs_dtem_dpre_ice_0d( &
2251  temp, pres, qdry, &
2252  dqsat_dT, dqsat_dP, &
2253  qsat )
2254  !$acc routine
2255  use scale_atmos_hydrometeor, only: &
2256  hydrometeor_lhs => atmos_hydrometeor_lhs
2257  implicit none
2258 
2259  real(RP), intent(in) :: temp
2260  real(RP), intent(in) :: pres
2261  real(RP), intent(in) :: qdry
2262 
2263  real(RP), intent(out) :: dqsat_dT
2264  real(RP), intent(out) :: dqsat_dP
2265 
2266  real(RP), intent(out), optional :: qsat
2267 
2268  real(RP) :: LHS ! latent heat of sublimation [J/kg]
2269  real(RP) :: psat
2270  real(RP) :: den1, den2 ! denominator
2271 
2272  !---------------------------------------------------------------------------
2273 
2274  call hydrometeor_lhs( temp, lhs ) ! [IN], [OUT]
2275  call atmos_saturation_psat_ice( temp, psat ) ! [IN], [OUT]
2276 
2277  den1 = ( pres - (1.0_rp-epsvap) * psat )**2
2278  den2 = den1 * rvap * temp**2
2279 
2280  dqsat_dp = - epsvap * psat / den1
2281  dqsat_dt = epsvap * psat / den2 * lhs * pres
2282 
2283  if ( present(qsat) ) &
2284  call atmos_saturation_psat2qsat_pres_qdry_0d( psat, pres, qdry, & ! [IN]
2285  qsat ) ! [OUT]
2286  return
2287  end subroutine atmos_saturation_dqs_dtem_dpre_ice_0d
2288 
2289  !-----------------------------------------------------------------------------
2290 !OCL SERIAL
2292  KA, KS, KE, &
2293  temp, pres, qdry, &
2294  dqsat_dT, dqsat_dP )
2295  !$acc routine vector
2296  implicit none
2297  integer, intent(in) :: KA, KS, KE
2298 
2299  real(RP), intent(in) :: temp (KA)
2300  real(RP), intent(in) :: pres (KA)
2301  real(RP), intent(in) :: qdry (KA)
2302 
2303  real(RP), intent(out) :: dqsat_dT(KA)
2304  real(RP), intent(out) :: dqsat_dP(KA)
2305 
2306  integer :: k
2307  !---------------------------------------------------------------------------
2308 
2309  do k = ks, ke
2310  call atmos_saturation_dqs_dtem_dpre_ice( temp(k), pres(k), qdry(k), & ! [IN]
2311  dqsat_dt(k), dqsat_dp(k) ) ! [OUT]
2312  enddo
2313 
2314  return
2316 
2317  !-----------------------------------------------------------------------------
2318  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
2319  !-----------------------------------------------------------------------------
2320  subroutine atmos_saturation_dqs_dtem_dpre_ice_3d( &
2321  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2322  temp, pres, qdry, &
2323  dqsat_dT, dqsat_dP )
2324  implicit none
2325  integer, intent(in) :: KA, KS, KE
2326  integer, intent(in) :: IA, IS, IE
2327  integer, intent(in) :: JA, JS, JE
2328 
2329  real(RP), intent(in) :: temp (KA,IA,JA)
2330  real(RP), intent(in) :: pres (KA,IA,JA)
2331  real(RP), intent(in) :: qdry (KA,IA,JA)
2332 
2333  real(RP), intent(out) :: dqsat_dT(KA,IA,JA)
2334  real(RP), intent(out) :: dqsat_dP(KA,IA,JA)
2335 
2336  integer :: k, i, j
2337  !---------------------------------------------------------------------------
2338 
2339  !$omp parallel do OMP_SCHEDULE_ collapse(2)
2340  !$acc kernels copyin(temp,pres,qdry) copyout(dqsat_dT,dqsat_dP)
2341  do j = js, je
2342  do i = is, ie
2343  do k = ks, ke
2344  call atmos_saturation_dqs_dtem_dpre_ice( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
2345  dqsat_dt(k,i,j), dqsat_dp(k,i,j) ) ! [OUT]
2346  enddo
2347  enddo
2348  enddo
2349  !$acc end kernels
2350 
2351  return
2352  end subroutine atmos_saturation_dqs_dtem_dpre_ice_3d
2353 
2354  !-----------------------------------------------------------------------------
2355  ! (d qsi/d T)_{p} and (d qs/d p)_{T}
2356  !-----------------------------------------------------------------------------
2357  subroutine atmos_saturation_dqs_dtem_dpre_all_0d( &
2358  temp, pres, qdry, &
2359  dqsat_dT, dqsat_dP, &
2360  qsat, qsat_liq, qsat_ice, &
2361  alpha )
2362  !$acc routine
2363  implicit none
2364 
2365  real(RP), intent(in) :: temp
2366  real(RP), intent(in) :: pres
2367  real(RP), intent(in) :: qdry
2368 
2369  real(RP), intent(out) :: dqsat_dT
2370  real(RP), intent(out) :: dqsat_dP
2371 
2372  real(RP), intent(out), optional :: qsat
2373  real(RP), intent(out), optional :: qsat_liq
2374  real(RP), intent(out), optional :: qsat_ice
2375  real(RP), intent(out), optional :: alpha
2376 
2377  real(RP) :: qsat_liq_, qsat_ice_, alpha_
2378  real(RP) :: dqsat_dT_liq, dqsat_dT_ice
2379  real(RP) :: dqsat_dP_liq, dqsat_dP_ice
2380  real(RP) :: dalpha_dT
2381 
2382  !---------------------------------------------------------------------------
2383 
2384  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp, pres, qdry, & ! [IN]
2385  dqsat_dt_liq, dqsat_dp_liq, & ! [OUT]
2386  qsat_liq_ ) ! [OUT]
2387  call atmos_saturation_dqs_dtem_dpre_ice_0d( temp, pres, qdry, & ! [IN]
2388  dqsat_dt_ice, dqsat_dp_ice, & ! [OUT]
2389  qsat_ice_ ) ! [OUT]
2390  call atmos_saturation_alpha ( temp, alpha_ ) ! [IN], [OUT]
2391  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
2392 
2393  dqsat_dt = qsat_liq_ * dalpha_dt + dqsat_dt_liq * ( alpha_ ) &
2394  - qsat_ice_ * dalpha_dt + dqsat_dt_ice * ( 1.0_rp-alpha_ )
2395  dqsat_dp = dqsat_dp_liq * ( alpha_ ) &
2396  + dqsat_dp_ice * ( 1.0_rp-alpha_ )
2397 
2398  if ( present(qsat) ) qsat = qsat_liq * alpha_ + qsat_ice * ( 1.0_rp-alpha_ )
2399  if ( present(qsat_liq) ) qsat_liq = qsat_liq_
2400  if ( present(qsat_ice) ) qsat_ice = qsat_ice_
2401  if ( present(alpha ) ) alpha = alpha_
2402 
2403  return
2404  end subroutine atmos_saturation_dqs_dtem_dpre_all_0d
2405 
2406  !-----------------------------------------------------------------------------
2408  !-----------------------------------------------------------------------------
2409  subroutine atmos_saturation_tdew_liq_0d( &
2410  DENS, TEMP, QV, &
2411  Tdew, &
2412  converged )
2413  !$acc routine
2414  use scale_const, only: &
2415  undef => const_undef
2416  use scale_atmos_hydrometeor, only: &
2417  atmos_hydrometeor_lhv
2418  real(RP), intent(in) :: DENS
2419  real(RP), intent(in) :: TEMP
2420  real(RP), intent(in) :: QV
2421 
2422  real(RP), intent(out) :: Tdew
2423  logical, intent(out) :: converged
2424 
2425  real(RP), parameter :: A = 17.625_rp
2426  real(RP), parameter :: B = 243.04_rp
2427  real(RP), parameter :: C = 610.94_rp
2428  integer, parameter :: itelim = 100
2429  real(RP), parameter :: criteria = 0.1_rp**(2+rp/2)
2430 
2431  real(RP) :: lhv
2432  real(RP) :: pvap, psat
2433  real(RP) :: dpsat_dT
2434  real(RP) :: dTdew
2435 
2436  integer :: ite
2437  !---------------------------------------------------------------------------
2438 
2439  pvap = dens * qv * rvap * temp
2440 
2441  if ( pvap < psat_min_liq ) then
2442  converged = .true.
2443  tdew = undef
2444  return
2445  end if
2446 
2447  ! first guess is calculated by Alduchov and Eskridge (1996)
2448  ! See Lawrence (2005) BAMS
2449  tdew = b * log( pvap / c ) / ( a - log( pvap / c ) ) + tem00
2450  converged = .false.
2451  !$acc loop seq
2452  do ite = 1, itelim
2453 
2454  call atmos_saturation_psat_liq( tdew, psat ) ! [IN], [OUT]
2455  call atmos_hydrometeor_lhv( tdew, lhv )
2456 
2457  dpsat_dt = psat * lhv / ( rvap * tdew**2 )
2458  dtdew = ( psat - pvap ) / dpsat_dt
2459  if ( abs(dtdew) < criteria .or. abs(psat-pvap) < criteria ) then
2460  converged = .true.
2461  exit
2462  end if
2463 
2464  tdew = tdew - dtdew
2465  end do
2466 #ifndef _OPENACC
2467  if( .not. converged ) then
2468  log_warn("ATMOS_SATURATION_tdew_liq_0D",*) dens, temp, qv, pvap, tdew, dtdew, dpsat_dt
2469  endif
2470 #endif
2471 
2472  return
2473  end subroutine atmos_saturation_tdew_liq_0d
2474 
2475 !OCL SERIAL
2476  subroutine atmos_saturation_tdew_liq_1d( &
2477  KA, KS, KE, &
2478  DENS, TEMP, QV, &
2479  Tdew, &
2480  converged )
2481  !$acc routine vector
2482  integer, intent(in) :: KA, KS, KE
2483 
2484  real(RP), intent(in) :: DENS(KA)
2485  real(RP), intent(in) :: TEMP(KA)
2486  real(RP), intent(in) :: QV (KA)
2487 
2488  real(RP), intent(out) :: Tdew(KA)
2489  logical, intent(out) :: converged
2490 
2491  integer :: k
2492  !---------------------------------------------------------------------------
2493 
2494  do k = ks, ke
2495  call atmos_saturation_tdew_liq_0d( dens(k), temp(k), qv(k), & ! [IN]
2496  tdew(k), converged ) ! [OUT]
2497  if ( .not. converged ) exit
2498  end do
2499 
2500  end subroutine atmos_saturation_tdew_liq_1d
2501 
2502  subroutine atmos_saturation_tdew_liq_3d( &
2503  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2504  DENS, TEMP, QV, &
2505  Tdew )
2506  use scale_prc, only: &
2507  prc_abort
2508  integer, intent(in) :: KA, KS, KE
2509  integer, intent(in) :: IA, IS, IE
2510  integer, intent(in) :: JA, JS, JE
2511 
2512  real(RP), intent(in) :: DENS(KA,IA,JA)
2513  real(RP), intent(in) :: TEMP(KA,IA,JA)
2514  real(RP), intent(in) :: QV (KA,IA,JA)
2515 
2516  real(RP), intent(out) :: Tdew(KA,IA,JA)
2517 
2518  logical :: converged, error
2519  integer :: k, i, j
2520  !---------------------------------------------------------------------------
2521 
2522  error = .false.
2523  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
2524  !$omp private(converged)
2525  !$acc kernels copyin(DENS,TEMP,QV) copyout(Tdew)
2526  !$acc loop reduction(.or.:error)
2527  do j = js, je
2528  !$acc loop reduction(.or.:error)
2529  do i = is, ie
2530  !$acc loop private(converged) reduction(.or.:error)
2531  do k = ks, ke
2532  call atmos_saturation_tdew_liq_0d( dens(k,i,j), temp(k,i,j), qv(k,i,j), & ! [IN]
2533  tdew(k,i,j), converged ) ! [OUT]
2534  if ( .not. converged ) then
2535  error = .true.
2536 #ifndef _OPENACC
2537  log_error("ATMOS_SATURATION_tdew_liq_3D",*) 'not converged! ', k,i,j
2538  exit
2539 #endif
2540  end if
2541  end do
2542  end do
2543  end do
2544  !$acc end kernels
2545 
2546  if ( error ) call prc_abort
2547 
2548  end subroutine atmos_saturation_tdew_liq_3d
2549 
2550  !-----------------------------------------------------------------------------
2557  !-----------------------------------------------------------------------------
2558  subroutine atmos_saturation_pote_0d( &
2559  DENS, POTT, TEMP, QV, &
2560  POTE )
2561  !$acc routine
2562  use scale_const, only: &
2563  rdry => const_rdry, &
2564  cpdry => const_cpdry
2565  use scale_atmos_hydrometeor, only: &
2566  atmos_hydrometeor_lhv
2567  real(RP), intent(in) :: DENS
2568  real(RP), intent(in) :: POTT
2569  real(RP), intent(in) :: TEMP
2570  real(RP), intent(in) :: QV
2571 
2572  real(RP), intent(out) :: POTE
2573 
2574  real(RP) :: TL
2575  real(RP) :: Pv
2576  real(RP) :: LHV
2577 
2578  pv = dens * qv * rvap * temp
2579  tl = 55.0_rp + 2840.0_rp / ( cpdry / rdry * log(temp) - log(pv) - 4.805_rp )
2580  call atmos_hydrometeor_lhv( temp, lhv ) ! [IN], [OUT]
2581 
2582  pote = pott * exp( lhv * qv / ( cpdry * tl ) &
2583  * 1.0784_rp * ( 1.0_rp + 0.810_rp * qv ) )
2584 
2585  return
2586  end subroutine atmos_saturation_pote_0d
2587 
2588 !OCL SERIAL
2589  subroutine atmos_saturation_pote_1d( &
2590  KA, KS, KE, &
2591  DENS, POTT, TEMP, QV, &
2592  POTE )
2593  !$acc routine vector
2594  integer, intent(in) :: KA, KS, KE
2595 
2596  real(RP), intent(in) :: DENS(KA)
2597  real(RP), intent(in) :: POTT(KA)
2598  real(RP), intent(in) :: TEMP(KA)
2599  real(RP), intent(in) :: QV (KA)
2600 
2601  real(RP), intent(out) :: POTE(KA)
2602 
2603  integer :: k
2604 
2605  do k = ks, ke
2606  call atmos_saturation_pote_0d( dens(k), pott(k), temp(k), qv(k), & ! [IN]
2607  pote(k) ) ! [OUT]
2608  end do
2609 
2610  return
2611  end subroutine atmos_saturation_pote_1d
2612 
2613  subroutine atmos_saturation_pote_3d( &
2614  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2615  DENS, POTT, TEMP, QV, &
2616  POTE )
2617  integer, intent(in) :: KA, KS, KE
2618  integer, intent(in) :: IA, IS, IE
2619  integer, intent(in) :: JA, JS, JE
2620 
2621  real(RP), intent(in) :: DENS(KA,IA,JA)
2622  real(RP), intent(in) :: POTT(KA,IA,JA)
2623  real(RP), intent(in) :: TEMP(KA,IA,JA)
2624  real(RP), intent(in) :: QV (KA,IA,JA)
2625 
2626  real(RP), intent(out) :: POTE(KA,IA,JA)
2627 
2628  integer :: k, i, j
2629 
2630  !$omp parallel do
2631  !$acc kernels copyin(DENS,POTT,TEMP,QV) copyout(POTE)
2632  do j = js, je
2633  do i = is, ie
2634  do k = ks, ke
2635  call atmos_saturation_pote_0d( &
2636  dens(k,i,j), pott(k,i,j), temp(k,i,j), qv(k,i,j), & ! [IN]
2637  pote(k,i,j) ) ! [OUT]
2638  end do
2639  end do
2640  end do
2641  !$acc end kernels
2642 
2643  return
2644  end subroutine atmos_saturation_pote_3d
2645 
2646  !-----------------------------------------------------------------------------
2648  !-----------------------------------------------------------------------------
2649 !OCL SERIAL
2650  subroutine atmos_saturation_moist_conversion_dens_liq_0d( &
2651  DENS, Emoist0, &
2652  TEMP, QV, QC, CPtot, CVtot, &
2653  converged )
2654  !$acc routine
2655  use scale_atmos_hydrometeor, only: &
2656  cp_vapor, &
2657  cp_water, &
2658  cv_vapor, &
2659  cv_water, &
2660  lhv
2661  implicit none
2662 
2663  real(RP), intent(in) :: DENS
2664  real(RP), intent(in) :: Emoist0
2665 
2666  real(RP), intent(inout) :: TEMP
2667  real(RP), intent(inout) :: QV
2668  real(RP), intent(inout) :: QC
2669  real(RP), intent(inout) :: CPtot
2670  real(RP), intent(inout) :: CVtot
2671 
2672  logical, intent(out) :: converged
2673 
2674  ! working
2675  real(RP) :: QSUM
2676  real(RP) :: TEMP0, QV0, QC0
2677  real(RP) :: CPtot0, CVtot0
2678  real(RP) :: qsat
2679  real(RP) :: Emoist ! moist internal energy
2680 
2681  ! d(X)/dT
2682  real(RP) :: dqsatl_dT
2683  real(RP) :: dqc_dT
2684  real(RP) :: dCVtot_dT
2685  real(RP) :: dEmoist_dT
2686  real(RP) :: dtemp
2687 
2688  integer, parameter :: itelim = 100
2689  real(RP), parameter :: dtemp_criteria = 0.1_rp**(2+rp/2)
2690  integer :: ite
2691  !---------------------------------------------------------------------------
2692 
2693  temp0 = temp
2694  cptot0 = cptot
2695  cvtot0 = cvtot
2696  qv0 = qv
2697  qc0 = qc
2698  qsum = qv + qc
2699 
2700  ! Check if saturation condition would not be met after evaporating all cloud water and ice virtually.
2701  qv = qsum
2702  qc = 0.0_rp
2703  cptot = cptot0 + qc0 * ( cp_vapor - cp_water )
2704  cvtot = cvtot0 + qc0 * ( cv_vapor - cv_water )
2705  temp = ( emoist0 - lhv * qv ) / cvtot
2706 
2707  call atmos_saturation_dens2qsat_liq( temp, dens, & ! [IN]
2708  qsat ) ! [OUT]
2709 
2710  if ( qsum <= qsat ) then
2711  converged = .true.
2712  return
2713  end if
2714 
2715  ! Iteration
2716  qv = qv0
2717  qc = qc0
2718  temp = temp0
2719  cptot = cptot0
2720  cvtot = cvtot0
2721 
2722  converged = .false.
2723  !$acc loop seq
2724  do ite = 1, itelim
2725 
2726  ! get qsat and dX/dT under the temp
2727  call atmos_saturation_dqs_dtem_dens_liq_0d( temp, dens, & ! [IN]
2728  dqsatl_dt, qsat ) ! [OUT]
2729 
2730  dqc_dt = - dqsatl_dt
2731 
2732  dcvtot_dt = dqsatl_dt * cv_vapor &
2733  + dqc_dt * cv_water
2734 
2735  demoist_dt = temp * dcvtot_dt + cvtot + dqsatl_dt * lhv
2736 
2737  ! diagnose quantities with the qsat
2738  qc = qsum - qsat
2739  cvtot = cvtot0 + ( cv_vapor - cv_water ) * ( qsat - qv )
2740  emoist = temp * cvtot + qsat * lhv
2741 
2742 
2743  ! update temp by the newtonian method
2744  dtemp = ( emoist - emoist0 ) / demoist_dt
2745  temp = temp - dtemp
2746 
2747  if ( abs(dtemp) < dtemp_criteria ) then
2748  converged = .true.
2749  exit
2750  endif
2751 
2752  if( temp*0.0_rp /= 0.0_rp ) exit
2753  enddo
2754 
2755  cptot = cptot + ( cp_vapor - cp_water ) * ( qsat - qv )
2756 
2757  qv = qsat
2758 
2759  return
2760  end subroutine atmos_saturation_moist_conversion_dens_liq_0d
2761 
2762  !-----------------------------------------------------------------------------
2764  !-----------------------------------------------------------------------------
2765 !OCL SERIAL
2767  DENS, Emoist0, &
2768  TEMP, QV, QC, QI, CPtot, CVtot, &
2769  converged )
2770  !$acc routine
2771  use scale_atmos_hydrometeor, only: &
2772  cp_vapor, &
2773  cp_water, &
2774  cp_ice, &
2775  cv_vapor, &
2776  cv_water, &
2777  cv_ice, &
2778  lhv, &
2779  lhf
2780  implicit none
2781  real(RP), intent(in) :: DENS
2782  real(RP), intent(in) :: Emoist0
2783 
2784  real(RP), intent(inout) :: TEMP
2785  real(RP), intent(inout) :: QV
2786  real(RP), intent(inout) :: QC
2787  real(RP), intent(inout) :: QI
2788  real(RP), intent(inout) :: CPtot
2789  real(RP), intent(inout) :: CVtot
2790 
2791  logical, intent(out) :: converged
2792 
2793  ! working
2794  real(RP) :: QSUM
2795  real(RP) :: TEMP0
2796  real(RP) :: QV0, QC0, QI0
2797  real(RP) :: CPtot0, CVtot0
2798  real(RP) :: alpha
2799  real(RP) :: qsat
2800  real(RP) :: Emoist ! moist internal energy
2801 
2802  ! d(X)/dT
2803  real(RP) :: dalpha_dT
2804  real(RP) :: dqsat_dT
2805  real(RP) :: dqc_dT, dqi_dT
2806  real(RP) :: dCVtot_dT
2807  real(RP) :: dEmoist_dT
2808  real(RP) :: dtemp
2809 
2810  integer, parameter :: itelim = 100
2811  real(RP), parameter :: dtemp_criteria = 0.1_rp**(2+rp/2)
2812  integer :: ite
2813  !---------------------------------------------------------------------------
2814 
2815  temp0 = temp
2816  cptot0 = cptot
2817  cvtot0 = cvtot
2818  qv0 = qv
2819  qc0 = qc
2820  qi0 = qi
2821  qsum = qv + qc + qi
2822 
2823  ! Check if saturation condition would not be met after evaporating all cloud water and ice virtually.
2824  qv = qsum
2825  qc = 0.0_rp
2826  qi = 0.0_rp
2827  cptot = cptot0 + qc0 * ( cp_vapor - cp_water ) + qi0 * ( cp_vapor - cp_ice )
2828  cvtot = cvtot0 + qc0 * ( cv_vapor - cv_water ) + qi0 * ( cv_vapor - cv_ice )
2829  temp = ( emoist0 - lhv * qv ) / cvtot
2830 
2831  call atmos_saturation_dens2qsat_all( temp, dens, & ! [IN]
2832  qsat ) ! [OUT]
2833 
2834  if ( qsum <= qsat ) then
2835  converged = .true.
2836  return
2837  end if
2838 
2839  ! Iteration
2840  qv = qv0
2841  qc = qc0
2842  qi = qi0
2843  temp = temp0
2844  cptot = cptot0
2845  cvtot = cvtot0
2846 
2847  converged = .false.
2848  !$acc loop seq
2849  do ite = 1, itelim
2850 
2851  ! dX/dT
2852  call atmos_saturation_dqs_dtem_dens_all_0d( temp, dens, & ! [IN]
2853  dqsat_dt, & ! [OUT]
2854  qsat=qsat, alpha=alpha ) ! [OUT]
2855  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
2856 
2857  dqc_dt = ( qsum - qv ) * dalpha_dt - dqsat_dt * ( alpha )
2858  dqi_dt = -( qsum - qv ) * dalpha_dt - dqsat_dt * ( 1.0_rp-alpha )
2859 
2860  dcvtot_dt = dqsat_dt * cv_vapor &
2861  + dqc_dt * cv_water &
2862  + dqi_dt * cv_ice
2863 
2864  demoist_dt = temp * dcvtot_dt + cvtot + dqsat_dt * lhv - dqi_dt * lhf
2865 
2866  ! Saturation
2867  qv = qsat
2868  qc = ( qsum - qsat ) * ( alpha )
2869  qi = ( qsum - qsat ) * ( 1.0_rp - alpha )
2870 
2871  cvtot = cvtot0 &
2872  + cv_vapor * ( qv - qv0 ) &
2873  + cv_water * ( qc - qc0 ) &
2874  + cv_ice * ( qi - qi0 )
2875 
2876 
2877 ! Emoist = temp * CVtot + qv * LHV - qi * LHF
2878  emoist = temp * cvtot + qsat * lhv - qi * lhf
2879 
2880  ! update temp by the newtonian method
2881  dtemp = ( emoist - emoist0 ) / demoist_dt
2882  temp = temp - dtemp
2883 
2884  if ( abs(dtemp) < dtemp_criteria ) then
2885  converged = .true.
2886  exit
2887  endif
2888 
2889  if( temp*0.0_rp /= 0.0_rp ) exit
2890  enddo
2891 
2892  cptot = cptot &
2893  + cp_vapor * ( qv - qv0 ) &
2894  + cp_water * ( qc - qc0 ) &
2895  + cp_ice * ( qi - qi0 )
2896 
2897  return
2899 
2900  !-----------------------------------------------------------------------------
2902  !-----------------------------------------------------------------------------
2904  PRES, Entr, Qdry, &
2905  QV, QC, &
2906  Rtot, CPtot, &
2907  TEMP, &
2908  converged )
2909  !$acc routine
2910  use scale_const, only: &
2911  eps => const_eps, &
2912  lhv0 => const_lhv0
2913  use scale_atmos_hydrometeor, only: &
2914  atmos_hydrometeor_entr, &
2915  atmos_hydrometeor_entr2temp, &
2916  cp_vapor, &
2917  cp_water
2918 
2919  real(RP), intent(in) :: PRES
2920  real(RP), intent(in) :: ENTR
2921  real(RP), intent(in) :: Qdry
2922 
2923  real(RP), intent(inout) :: QV
2924  real(RP), intent(inout) :: QC
2925  real(RP), intent(inout) :: CPtot
2926  real(RP), intent(inout) :: Rtot
2927 
2928  real(RP), intent(out) :: TEMP
2929  logical, intent(out) :: converged
2930 
2931  real(RP), parameter :: TEMMIN = 0.1_rp
2932  real(RP), parameter :: criteria = 1.e-8_rp
2933  integer, parameter :: itelim = 100
2934  integer :: ite
2935 
2936  real(RP) :: Qsum
2937  real(RP) :: qsat, psat
2938 
2939  real(RP) :: TEMP_ite
2940  real(RP) :: QV_ite
2941  real(RP) :: ENTR_ite
2942  real(RP) :: Rtot_ite
2943  real(RP) :: CPtot_ite
2944 
2945  real(RP) :: dqsat_dT, dqsat_dP
2946 
2947  real(RP) :: TEMP_prev
2948  real(RP) :: dENTR_dT
2949  !---------------------------------------------------------------------------
2950 
2951  qsum = qv + qc
2952 
2953  ! all liquid evaporates
2954  rtot = rtot + qc * rvap
2955  cptot = cptot + qc * ( cp_vapor - cp_water )
2956  qv = qsum
2957  qc = 0.0_rp
2958 
2959  call atmos_hydrometeor_entr2temp( entr, pres, qsum, 0.0_rp, qdry, & ! [IN]
2960  rtot, cptot, & ! [IN]
2961  temp ) ! [OUT]
2962 
2963  call atmos_saturation_pres2qsat_liq( temp, pres, qdry, & ! [IN]
2964  qsat ) ! [OUT]
2965  if ( qsum <= qsat ) then
2966  ! unsaturated
2967 
2968  converged = .true.
2969 
2970  return
2971 
2972  else
2973  ! saturated
2974 
2975  temp_ite = temp
2976 
2977  !$acc loop seq
2978  do ite = 1, itelim
2979 
2980  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp_ite, pres, qdry, & ! [IN]
2981  dqsat_dt, dqsat_dp, & ! [OUT]
2982  qsat=qsat, psat=psat ) ! [OUT]
2983 
2984  qv_ite = min( qsum, qsat )
2985 
2986  rtot_ite = rtot - ( qsum - qv_ite ) * rvap
2987  cptot_ite = cptot - ( qsum - qv_ite ) * ( cp_vapor - cp_water )
2988 
2989  dentr_dt = cptot_ite / temp_ite &
2990  + ( ( cp_vapor - cp_water ) * log( temp_ite / tem00 ) &
2991  - rvap * log( psat/ psat0 ) &
2992  + lhv0 / tem00 &
2993  ) * dqsat_dt
2994 
2995  call atmos_hydrometeor_entr( temp_ite, pres, & ! [IN]
2996  qv_ite, 0.0_rp, qdry, & ! [IN] ! QI = 0
2997  rtot_ite, cptot_ite, & ! [IN]
2998  entr_ite ) ! [OUT]
2999 
3000  temp_prev = temp_ite
3001  temp_ite = temp_ite - ( entr_ite - entr ) / max( dentr_dt, eps )
3002  temp_ite = max( temp_ite, temmin )
3003 
3004  if( abs(temp_ite-temp_prev) < criteria ) then
3005  converged = .true.
3006  exit
3007  end if
3008 
3009  enddo
3010 
3011  end if
3012 
3013  qv = qv_ite
3014  qc = qsum - qv
3015  cptot = cptot_ite
3016  rtot = rtot_ite
3017  temp = temp_ite
3018 
3019  return
3021 
3022 end module scale_atmos_saturation
scale_const::const_lhv0
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:82
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_hydrometeor::cp_water
real(rp), public cp_water
CP for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:152
scale_const::const_lhs00
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
Definition: scale_const.F90:85
scale_atmos_saturation::atmos_saturation_moist_conversion_pres_liq_0d
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.
Definition: scale_atmos_saturation.F90:2909
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:75
scale_const::const_lhv00
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:83
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:68
scale_atmos_saturation::atmos_saturation_dqs_dtem_dpre_ice_1d
subroutine atmos_saturation_dqs_dtem_dpre_ice_1d(KA, KS, KE, temp, pres, qdry, dqsat_dT, dqsat_dP)
Definition: scale_atmos_saturation.F90:2295
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_saturation::atmos_saturation_moist_conversion_dens_all_0d
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)
Definition: scale_atmos_saturation.F90:2770
scale_const::const_thermodyn_type
character(len=h_short), public const_thermodyn_type
internal energy type
Definition: scale_const.F90:111
scale_const::const_cpvap
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:69
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_hydrometeor::cv_vapor
real(rp), public cv_vapor
CV for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:149
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_saturation::atmos_saturation_tdew_liq_1d
subroutine atmos_saturation_tdew_liq_1d(KA, KS, KE, DENS, TEMP, QV, Tdew, converged)
Definition: scale_atmos_saturation.F90:2481
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_hydrometeor::lhv
real(rp), public lhv
latent heat of vaporization for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:144
scale_atmos_saturation::atmos_saturation_dqs_dtem_dens_ice_1d
subroutine atmos_saturation_dqs_dtem_dens_ice_1d(KA, KS, KE, temp, dens, dqsdtem)
Definition: scale_atmos_saturation.F90:2042
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_const::const_psat0
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:88
scale_atmos_saturation::atmos_saturation_setup
subroutine, public atmos_saturation_setup
Setup.
Definition: scale_atmos_saturation.F90:244
scale_const::const_ci
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:72
scale_atmos_saturation::atmos_saturation_pote_0d
subroutine atmos_saturation_pote_0d(DENS, POTT, TEMP, QV, POTE)
calculate equivalent potential temperature Bolton, D., 1980: The computation of equivalent potential ...
Definition: scale_atmos_saturation.F90:2561
scale_const::const_tem00
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:99
scale_const::const_cl
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:71
scale_const::const_cvvap
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:70
scale_atmos_saturation::atmos_saturation_pote_1d
subroutine atmos_saturation_pote_1d(KA, KS, KE, DENS, POTT, TEMP, QV, POTE)
Definition: scale_atmos_saturation.F90:2593
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:146
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
scale_atmos_saturation::atmos_saturation_alpha_0d
subroutine atmos_saturation_alpha_0d(temp, alpha)
calc liquid/ice separation factor (0D)
Definition: scale_atmos_saturation.F90:331
scale_atmos_saturation::atmos_saturation_dqs_dtem_dpre_liq_1d
subroutine atmos_saturation_dqs_dtem_dpre_liq_1d(KA, KS, KE, temp, pres, qdry, dqsat_dT, dqsat_dP)
Definition: scale_atmos_saturation.F90:2189
scale_const::const_lhs0
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:84
scale_atmos_saturation::atmos_saturation_dqs_dtem_dens_liq_1d
subroutine atmos_saturation_dqs_dtem_dens_liq_1d(KA, KS, KE, temp, dens, dqsdtem)
Definition: scale_atmos_saturation.F90:1949
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_hydrometeor::atmos_hydrometeor_setup
subroutine, public atmos_hydrometeor_setup
Setup.
Definition: scale_atmos_hydrometeor.F90:176
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_hydrometeor::cp_vapor
real(rp), public cp_vapor
CP for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:150
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:151
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:153
scale_atmos_hydrometeor::cp_ice
real(rp), public cp_ice
CP for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:154