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