SCALE-RM
scale_atmos_hydrometeor.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  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedure
26  !
27  public :: atmos_hydrometeor_setup
29  public :: atmos_hydrometeor_regist
30  public :: atmos_hydrometeor_lhv
31  public :: atmos_hydrometeor_lhv_para
32  public :: atmos_hydrometeor_lhs
33  public :: atmos_hydrometeor_lhs_para
34  public :: atmos_hydrometeor_lhf
35  public :: atmos_hydrometeor_lhf_para
36  public :: atmos_hydrometeor_entr
37  public :: atmos_hydrometeor_entr2temp
38 
39  interface atmos_hydrometeor_regist
40  module procedure atmos_hydrometeor_regist
41  end interface atmos_hydrometeor_regist
42 
43  interface atmos_hydrometeor_lhv
44  module procedure atmos_hydrometeor_lhv_0d
45  module procedure atmos_hydrometeor_lhv_1d
46  module procedure atmos_hydrometeor_lhv_2d
47  module procedure atmos_hydrometeor_lhv_3d
48  end interface atmos_hydrometeor_lhv
49 
50  interface atmos_hydrometeor_lhv_para
51  module procedure atmos_hydrometeor_lhv_1d_para
52  end interface atmos_hydrometeor_lhv_para
53 
54  interface atmos_hydrometeor_lhs
55  module procedure atmos_hydrometeor_lhs_0d
56  module procedure atmos_hydrometeor_lhs_1d
57  module procedure atmos_hydrometeor_lhs_2d
58  module procedure atmos_hydrometeor_lhs_3d
59  end interface atmos_hydrometeor_lhs
60 
61  interface atmos_hydrometeor_lhs_para
62  module procedure atmos_hydrometeor_lhs_1d_para
63  end interface atmos_hydrometeor_lhs_para
64 
65  interface atmos_hydrometeor_lhf
66  module procedure atmos_hydrometeor_lhf_0d
67  module procedure atmos_hydrometeor_lhf_1d
68  module procedure atmos_hydrometeor_lhf_2d
69  module procedure atmos_hydrometeor_lhf_3d
70  end interface atmos_hydrometeor_lhf
71 
72  interface atmos_hydrometeor_lhf_para
73  module procedure atmos_hydrometeor_lhf_1d_para
74  end interface atmos_hydrometeor_lhf_para
75 
76  interface atmos_hydrometeor_entr
77  module procedure atmos_hydrometeor_entr_0d
78  module procedure atmos_hydrometeor_entr_2d
79  module procedure atmos_hydrometeor_entr_3d
80  end interface atmos_hydrometeor_entr
81 
82  interface atmos_hydrometeor_entr2temp
83  module procedure atmos_hydrometeor_entr2temp_0d
84  end interface atmos_hydrometeor_entr2temp
85 
86  !-----------------------------------------------------------------------------
87  !
88  !++ Public parameters & variables
89  !
90  logical, public :: atmos_hydrometeor_ice_phase
91 
92 
93  integer, public :: i_qv = -1
94 
95  integer, public, parameter :: n_hyd = 6
96 
97  integer, public, parameter :: i_hc = 1
98  integer, public, parameter :: i_hr = 2
99  integer, public, parameter :: i_hi = 3
100  integer, public, parameter :: i_hs = 4
101  integer, public, parameter :: i_hg = 5
102  integer, public, parameter :: i_hh = 6
103 
104  character(len=H_SHORT), public, parameter :: hyd_name(n_hyd) = &
105  (/ "QC", "QR", "QI", "QS", "QG", "QH" /)
106  character(len=H_MID), public, parameter :: hyd_desc(n_hyd) = &
107  (/ "cloud ", "rain ", "ice water", "snow ", "graupel ", "hail " /)
108  character(len=H_SHORT), public, parameter :: num_name(n_hyd) = &
109  (/ "NC", "NR", "NI", "NS", "NG", "NH" /)
110 
111  real(rp), public :: hyd_dens(n_hyd)
112  !$acc declare create(HYD_DENS)
113 
114  logical, public :: atmos_hydrometeor_dry = .true.
115 
116 ! integer, public :: I_NC = -1
117 ! integer, public :: I_NR = -1
118 ! integer, public :: I_NI = -1
119 ! integer, public :: I_NS = -1
120 ! integer, public :: I_NG = -1
121 ! integer, public :: I_NH = -1
122 
123 ! integer, public :: I_QC = -1
124 ! integer, public :: I_QR = -1
125 ! integer, public :: I_QI = -1
126 ! integer, public :: I_QS = -1
127 ! integer, public :: I_QG = -1
128 ! integer, public :: I_QH = -1
129 
130  ! hydrometeor (water + ice)
131  integer, public :: qha = 0
132  integer, public :: qhs = -1
133  integer, public :: qhe = -2
134  ! water
135  integer, public :: qla = 0
136  integer, public :: qls = -1
137  integer, public :: qle = -2
138  ! ice
139  integer, public :: qia = 0
140  integer, public :: qis = -1
141  integer, public :: qie = -2
142  !$acc declare create(QHA, QHS, QHE, QLA, QLS, QLE, QIA, QIS, QIE)
143 
144  real(rp), public :: lhv
145  real(rp), public :: lhs
146  real(rp), public :: lhf
147  !$acc declare create(LHV, LHS, LHF)
148 
149  real(rp), public :: cv_vapor
150  real(rp), public :: cp_vapor
151  real(rp), public :: cv_water
152  real(rp), public :: cp_water
153  real(rp), public :: cv_ice
154  real(rp), public :: cp_ice
155  !$acc declare create(CV_VAPOR,CP_VAPOR,CV_WATER,CP_WATER,CV_ICE,CP_ICE)
156 
157  !-----------------------------------------------------------------------------
158  !
159  !++ Private procedure
160  !
161  !-----------------------------------------------------------------------------
162  !
163  !++ Private parameters & variables
164  !
165  real(rp), private :: thermodyn_emask = 1.0_rp
166  !$acc declare create(THERMODYN_EMASK)
167 
168  logical, private :: initialized = .false.
169 
170  !-----------------------------------------------------------------------------
171 
172 contains
173  !-----------------------------------------------------------------------------
175  subroutine atmos_hydrometeor_setup
176  use scale_const, only: &
177  const_setup, &
178  cpvap => const_cpvap, &
179  cvvap => const_cvvap, &
180  cl => const_cl, &
181  ci => const_ci, &
182  lhv00 => const_lhv00, &
183  lhs00 => const_lhs00, &
184  lhf00 => const_lhf00, &
185  lhv0 => const_lhv0, &
186  lhs0 => const_lhs0, &
187  lhf0 => const_lhf0, &
188  dwatr => const_dwatr, &
189  dice => const_dice, &
190  thermodyn_type => const_thermodyn_type
191  use scale_prc, only: &
192  prc_abort
193  implicit none
194  !---------------------------------------------------------------------------
195 
196  if ( initialized ) return
197  initialized = .true.
198 
199  call const_setup
200 
201  log_newline
202  log_info("ATMOS_HYDROMETEOR_setup",*) 'Setup'
203 
204  if ( thermodyn_type == 'EXACT' ) then
205 
206  cv_vapor = cvvap
207  cp_vapor = cpvap
208  cv_water = cl
210  cv_ice = ci
211  cp_ice = cv_ice
212 
213  lhv = lhv00
214  lhs = lhs00
215  lhf = lhf00
216  thermodyn_emask = 1.0_rp
217 
218  elseif( thermodyn_type == 'SIMPLE' ) then
219 
220  cv_vapor = cvvap
221  cp_vapor = cpvap
222  cv_water = cvvap
224  cv_ice = cvvap
225  cp_ice = cv_ice
226 
227  lhv = lhv0
228  lhs = lhs0
229  lhf = lhf0
230  thermodyn_emask = 0.0_rp
231 
232  else
233  log_error("ATMOS_HYDROMETEOR_setup",*) 'Not appropriate ATMOS_THERMODYN_ENERGY_TYPE. Check!', trim(thermodyn_type)
234  call prc_abort
235  endif
236 
237  hyd_dens(:) = (/ dwatr, & ! HC
238  dwatr, & ! HR
239  dice, & ! HI
240  dice, & ! HS
241  dice, & ! HG
242  dice /) ! HH
243 
244  !$acc update device(HYD_DENS)
245  !$acc update device(LHV, LHS, LHF)
246  !$acc update device(CV_VAPOR,CP_VAPOR,CV_WATER,CP_WATER,CV_ICE,CP_ICE)
247  !$acc update device(THERMODYN_EMASK)
248 
249  return
250  end subroutine atmos_hydrometeor_setup
251 
252  !-----------------------------------------------------------------------------
254  subroutine atmos_hydrometeor_finalize
255  implicit none
256  !---------------------------------------------------------------------------
257 
258  log_newline
259  log_info("ATMOS_HYDROMETEOR_finalize",*) 'Finalize'
260 
261  atmos_hydrometeor_dry = .true.
262  initialized = .false.
263 
264  return
265  end subroutine atmos_hydrometeor_finalize
266 
267  !-----------------------------------------------------------------------------
271  subroutine atmos_hydrometeor_regist( &
272  NL, &
273  NI, &
274  NAME, &
275  DESC, &
276  UNIT, &
277  Q0, &
278  ADVC )
279  use scale_prc, only: &
280  prc_abort
281  use scale_tracer, only: &
283  use scale_const, only: &
284  rvap => const_rvap
285  implicit none
286  integer, intent(in) :: nl
287  integer, intent(in) :: ni
288  character(len=*), intent(in) :: name(1+nl+ni)
289  character(len=*), intent(in) :: desc(1+nl+ni)
290  character(len=*), intent(in) :: unit(1+nl+ni)
291 
292  integer, intent(out) :: q0
293 
294  logical, intent(in), optional :: advc(1+nl+ni)
295 
296  real(rp) :: cv (1+nl+ni)
297  real(rp) :: cp (1+nl+ni)
298  real(rp) :: r (1+nl+ni)
299  real(rp) :: ei0 (1+nl+ni)
300  logical :: mass (1+nl+ni)
301  logical :: advc_(1+nl+ni)
302 
303  integer :: nq
304  integer :: n
305  !---------------------------------------------------------------------------
306 
307  if ( .not. atmos_hydrometeor_dry ) then
308  log_error("ATMOS_HYDROMETEOR_regist",*) 'tracer for hydrometeor is already registerd'
309  call prc_abort
310  endif
311 
312  atmos_hydrometeor_dry = .false.
313 
314  nq = 0
315 
316  ! vapor
317  nq = nq + 1
318  cv(nq) = cv_vapor
319  cp(nq) = cp_vapor
320  r(nq) = rvap
321  ei0(nq) = lhv
322 
323  ! liquid
324  do n = 1, nl
325  nq = nq + 1
326  cv(nq) = cv_water
327  cp(nq) = cp_water
328  r(nq) = 0.0_rp
329  ei0(nq) = 0.0_rp
330  end do
331 
332  ! ice
333  do n = 1, ni
334  nq = nq + 1
335  cv(nq) = cv_ice
336  cp(nq) = cp_ice
337  r(nq) = 0.0_rp
338  ei0(nq) = - lhf
339  end do
340 
341  ! NQ = 1 + NL + NI, vapor + liqid + ice
342 
343  if ( present(advc) ) then
344  advc_(:) = advc(:)
345  else
346  advc_(:) = .true.
347  endif
348 
349  do n = 1, nq
350  mass(n) = .true.
351  end do
352 
353  call tracer_regist( q0, & ! [OUT]
354  nq, & ! [IN]
355  name, & ! [IN]
356  desc, & ! [IN]
357  unit, & ! [IN]
358  cv, cp, r, ei0, & ! [IN], optional
359  advc_, mass ) ! [IN], optional
360 
361  i_qv = q0
362 
363  if ( nq > 1 ) then
364  qhs = q0 + 1
365  qha = nl + ni
366  qhe = qhs + qha - 1
367  endif
368 
369  if ( nl > 0 ) then
370  qls = q0 + 1
371  qla = nl
372  qle = qls + nl - 1
373  endif
374 
375  if ( ni > 0 ) then
376  qis = qle + 1
377  qia = ni
378  qie = qis + ni - 1
379  endif
380 
382 
383  !$acc update device(QHA, QHS, QHE, QLA, QLS, QLE, QIA, QIS, QIE)
384 
385  return
386  end subroutine atmos_hydrometeor_regist
387 
388  !-----------------------------------------------------------------------------
389  subroutine atmos_hydrometeor_lhv_0d( &
390  temp, &
391  lhv )
392  !$acc routine
393  use scale_const, only: &
394  tem00 => const_tem00, &
395  lhv0 => const_lhv0
396  implicit none
397  real(RP), intent(in) :: temp
398  real(RP), intent(out) :: lhv
399  !---------------------------------------------------------------------------
400 
401  lhv = lhv0 + ( cp_vapor - cp_water ) * ( temp - tem00 ) * thermodyn_emask
402 
403  return
404  end subroutine atmos_hydrometeor_lhv_0d
405 
406  !-----------------------------------------------------------------------------
407 !OCL SERIAL
408  subroutine atmos_hydrometeor_lhv_1d( &
409  KA, KS, KE, &
410  temp, &
411  lhv )
412  !$acc routine vector
413  implicit none
414  integer, intent(in) :: KA, KS, KE
415 
416  real(RP), intent(in) :: temp(KA)
417 
418  real(RP), intent(out) :: lhv (KA)
419 
420  integer :: k
421  !---------------------------------------------------------------------------
422 
423  do k = ks, ke
424  call atmos_hydrometeor_lhv_0d( temp(k), lhv(k) )
425  enddo
426 
427  return
428  end subroutine atmos_hydrometeor_lhv_1d
429 
430  !-----------------------------------------------------------------------------
431  subroutine atmos_hydrometeor_lhv_1d_para( &
432  KA, KS, KE, &
433  temp, &
434  lhv )
435  implicit none
436  integer, intent(in) :: KA, KS, KE
437 
438  real(RP), intent(in) :: temp(KA)
439 
440  real(RP), intent(out) :: lhv (KA)
441 
442  integer :: k
443  !---------------------------------------------------------------------------
444 
445  !$omp parallel do
446  !$acc kernels copyin(temp) copyout(lhv)
447  do k = ks, ke
448  call atmos_hydrometeor_lhv_0d( temp(k), lhv(k) )
449  enddo
450  !$acc end kernels
451 
452  return
453  end subroutine atmos_hydrometeor_lhv_1d_para
454 
455  !-----------------------------------------------------------------------------
456  subroutine atmos_hydrometeor_lhv_2d( &
457  IA, IS, IE, JA, JS, JE, &
458  temp, &
459  lhv )
460  implicit none
461  integer, intent(in) :: IA, IS, IE
462  integer, intent(in) :: JA, JS, JE
463 
464  real(RP), intent(in) :: temp(IA,JA)
465  real(RP), intent(out) :: lhv (IA,JA)
466 
467  integer :: i, j
468  !---------------------------------------------------------------------------
469 
470  !$omp parallel do OMP_SCHEDULE_ collapse(2)
471  !$acc kernels copyin(temp) copyout(lhv)
472  do j = js, je
473  do i = is, ie
474  call atmos_hydrometeor_lhv_0d( temp(i,j), lhv(i,j) )
475  enddo
476  enddo
477  !$acc end kernels
478 
479  return
480  end subroutine atmos_hydrometeor_lhv_2d
481 
482  !-----------------------------------------------------------------------------
483  subroutine atmos_hydrometeor_lhv_3d( &
484  KA, KS, KE, &
485  IA, IS, IE, &
486  JA, JS, JE, &
487  temp, &
488  lhv )
489  implicit none
490  integer, intent(in) :: KA, KS, KE
491  integer, intent(in) :: IA, IS, IE
492  integer, intent(in) :: JA, JS, JE
493 
494  real(RP), intent(in) :: temp(KA,IA,JA)
495  real(RP), intent(out) :: lhv (KA,IA,JA)
496 
497  integer :: k, i, j
498  !---------------------------------------------------------------------------
499 
500  !$acc kernels copyin(temp) copyout(lhv)
501  do j = js, je
502  do i = is, ie
503  do k = ks, ke
504  call atmos_hydrometeor_lhv_0d( temp(k,i,j), lhv(k,i,j) )
505  enddo
506  enddo
507  enddo
508  !$acc end kernels
509 
510  return
511  end subroutine atmos_hydrometeor_lhv_3d
512 
513  !-----------------------------------------------------------------------------
514  subroutine atmos_hydrometeor_lhs_0d( &
515  temp, &
516  lhs )
517  !$acc routine
518  use scale_const, only: &
519  tem00 => const_tem00, &
520  lhs0 => const_lhs0
521  implicit none
522  real(RP), intent(in) :: temp
523  real(RP), intent(out) :: lhs
524  !---------------------------------------------------------------------------
525 
526  lhs = lhs0 + ( cp_vapor - cp_ice ) * ( temp - tem00 ) * thermodyn_emask
527 
528  return
529  end subroutine atmos_hydrometeor_lhs_0d
530 
531  !-----------------------------------------------------------------------------
532 !OCL SERIAL
533  subroutine atmos_hydrometeor_lhs_1d( &
534  KA, KS, KE, &
535  temp, &
536  lhs )
537  !$acc routine vector
538  implicit none
539  integer, intent(in) :: KA, KS, KE
540 
541  real(RP), intent(in) :: temp(KA)
542 
543  real(RP), intent(out) :: lhs (KA)
544 
545  integer :: k
546  !---------------------------------------------------------------------------
547 
548  do k = ks, ke
549  call atmos_hydrometeor_lhs( temp(k), lhs(k) )
550  enddo
551 
552  return
553  end subroutine atmos_hydrometeor_lhs_1d
554 
555  !-----------------------------------------------------------------------------
556  subroutine atmos_hydrometeor_lhs_1d_para( &
557  KA, KS, KE, &
558  temp, &
559  lhs )
560  implicit none
561  integer, intent(in) :: KA, KS, KE
562 
563  real(RP), intent(in) :: temp(KA)
564 
565  real(RP), intent(out) :: lhs (KA)
566 
567  integer :: k
568  !---------------------------------------------------------------------------
569 
570  !$omp parallel do
571  !$acc kernels copyin(temp) copyout(lhs)
572  do k = ks, ke
573  call atmos_hydrometeor_lhs( temp(k), lhs(k) )
574  enddo
575  !$acc end kernels
576 
577  return
578  end subroutine atmos_hydrometeor_lhs_1d_para
579 
580  !-----------------------------------------------------------------------------
581  subroutine atmos_hydrometeor_lhs_2d( &
582  IA, IS, IE, JA, JS, JE, &
583  temp, &
584  lhs )
585  implicit none
586  integer, intent(in) :: IA, IS, IE
587  integer, intent(in) :: JA, JS, JE
588 
589  real(RP), intent(in) :: temp(IA,JA)
590  real(RP), intent(out) :: lhs (IA,JA)
591 
592  integer :: i, j
593  !---------------------------------------------------------------------------
594 
595  !$omp parallel do OMP_SCHEDULE_ collapse(2)
596  !$acc kernels copyin(temp) copyout(lhs)
597  do j = js, je
598  do i = is, ie
599  call atmos_hydrometeor_lhs( temp(i,j), lhs(i,j) )
600  enddo
601  enddo
602  !$acc end kernels
603 
604  return
605  end subroutine atmos_hydrometeor_lhs_2d
606 
607  !-----------------------------------------------------------------------------
608  subroutine atmos_hydrometeor_lhs_3d( &
609  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
610  temp, &
611  lhs )
612  implicit none
613 
614  integer, intent(in) :: KA, KS, KE
615  integer, intent(in) :: IA, IS, IE
616  integer, intent(in) :: JA, JS, JE
617 
618  real(RP), intent(in) :: temp(KA,IA,JA)
619  real(RP), intent(out) :: lhs (KA,IA,JA)
620 
621  integer :: k, i, j
622  !---------------------------------------------------------------------------
623 
624  !$omp parallel do OMP_SCHEDULE_ collapse(2)
625  !$acc kernels copyin(temp) copyout(lhs)
626  do j = js, je
627  do i = is, ie
628  do k = ks, ke
629  call atmos_hydrometeor_lhs( temp(k,i,j), lhs(k,i,j) )
630  enddo
631  enddo
632  enddo
633  !$acc end kernels
634 
635  return
636  end subroutine atmos_hydrometeor_lhs_3d
637 
638  !-----------------------------------------------------------------------------
639  subroutine atmos_hydrometeor_lhf_0d( &
640  temp, &
641  lhf )
642  !$acc routine
643  use scale_const, only: &
644  tem00 => const_tem00, &
645  lhf0 => const_lhf0
646  implicit none
647  real(RP), intent(in) :: temp
648  real(RP), intent(out) :: lhf
649  !---------------------------------------------------------------------------
650 
651  lhf = lhf0 + ( cp_water - cp_ice ) * ( temp - tem00 ) * thermodyn_emask
652 
653  return
654  end subroutine atmos_hydrometeor_lhf_0d
655 
656  !-----------------------------------------------------------------------------
657 !OCL SERIAL
658  subroutine atmos_hydrometeor_lhf_1d( &
659  KA, KS, KE, &
660  temp, &
661  lhf )
662  !$acc routine vector
663  implicit none
664  integer, intent(in) :: KA, KS, KE
665 
666  real(RP), intent(in) :: temp(KA)
667  real(RP), intent(out) :: lhf (KA)
668 
669  integer :: k
670  !---------------------------------------------------------------------------
671 
672  do k = ks, ke
673  call atmos_hydrometeor_lhf( temp(k), lhf(k) )
674  enddo
675 
676  return
677  end subroutine atmos_hydrometeor_lhf_1d
678 
679  !-----------------------------------------------------------------------------
680  subroutine atmos_hydrometeor_lhf_1d_para( &
681  KA, KS, KE, &
682  temp, &
683  lhf )
684  implicit none
685  integer, intent(in) :: KA, KS, KE
686 
687  real(RP), intent(in) :: temp(KA)
688  real(RP), intent(out) :: lhf (KA)
689 
690  integer :: k
691  !---------------------------------------------------------------------------
692 
693  !$omp parallel do
694  !$acc kernels copyin(temp) copyout(lhf)
695  do k = ks, ke
696  call atmos_hydrometeor_lhf( temp(k), lhf(k) )
697  enddo
698  !$acc end kernels
699 
700  return
701  end subroutine atmos_hydrometeor_lhf_1d_para
702 
703  !-----------------------------------------------------------------------------
704  subroutine atmos_hydrometeor_lhf_2d( &
705  IA, IS, IE, JA, JS, JE, &
706  temp, &
707  lhf )
708  implicit none
709  integer, intent(in) :: IA, IS, IE
710  integer, intent(in) :: JA, JS, JE
711 
712  real(RP), intent(in) :: temp(IA,JA)
713  real(RP), intent(out) :: lhf (IA,JA)
714 
715  integer :: i, j
716  !---------------------------------------------------------------------------
717 
718  !$omp parallel do OMP_SCHEDULE_ collapse(2)
719  !$acc kernels copyin(temp) copyout(lhf)
720  do j = js, je
721  do i = is, ie
722  call atmos_hydrometeor_lhf( temp(i,j), lhf(i,j) )
723  enddo
724  enddo
725  !$acc end kernels
726 
727  return
728  end subroutine atmos_hydrometeor_lhf_2d
729 
730  !-----------------------------------------------------------------------------
731  subroutine atmos_hydrometeor_lhf_3d( &
732  KA, KS, KE, &
733  IA, IS, IE, &
734  JA, JS, JE, &
735  temp, &
736  lhf )
737  implicit none
738  integer, intent(in) :: KA, KS, KE
739  integer, intent(in) :: IA, IS, IE
740  integer, intent(in) :: JA, JS, JE
741 
742  real(RP), intent(in) :: temp(KA,IA,JA)
743  real(RP), intent(out) :: lhf (KA,IA,JA)
744 
745  integer :: k, i, j
746  !---------------------------------------------------------------------------
747 
748  !$omp parallel do OMP_SCHEDULE_ collapse(2)
749  !$acc kernels copyin(temp) copyout(lhf)
750  do j = js, je
751  do i = is, ie
752  do k = ks, ke
753  call atmos_hydrometeor_lhf( temp(k,i,j), lhf(k,i,j) )
754  enddo
755  enddo
756  enddo
757  !$acc end kernels
758 
759  return
760  end subroutine atmos_hydrometeor_lhf_3d
761 
762  !-----------------------------------------------------------------------------
764  subroutine atmos_hydrometeor_entr_0d( &
765  TEMP, PRES, &
766  QV, QI, Qdry, &
767  Rtot, CPtot, &
768  entr )
769  !$acc routine
770  use scale_const, only: &
771  pre00 => const_pre00, &
772  tem00 => const_tem00, &
773  rdry => const_rdry, &
774  rvap => const_rvap, &
775  lhv0 => const_lhv0, &
776  lhf0 => const_lhf0, &
777  psat0 => const_psat0
778  implicit none
779 
780  real(RP), intent(in) :: TEMP
781  real(RP), intent(in) :: PRES
782  real(RP), intent(in) :: QV
783  real(RP), intent(in) :: QI
784  real(RP), intent(in) :: Qdry
785  real(RP), intent(in) :: Rtot
786  real(RP), intent(in) :: CPtot
787 
788  real(RP), intent(out) :: entr
789 
790  real(RP) :: pres_dry, pres_vap
791  !---------------------------------------------------------------------------
792 
793  pres_dry = pres * qdry * rdry / rtot
794  pres_vap = pres * qv * rvap / rtot
795 
796  entr = cptot * log( temp / tem00 ) &
797  - qdry * rdry * log( pres_dry / pre00 ) &
798  - qv * rvap * log( pres_vap / psat0 ) &
799  + ( qv * lhv0 - qi * lhf0 ) / tem00
800 
801  return
802  end subroutine atmos_hydrometeor_entr_0d
803 
804  !-----------------------------------------------------------------------------
806  subroutine atmos_hydrometeor_entr_2d( &
807  IA, IS, IE, JA, JS, JE, &
808  TEMP, PRES, &
809  QV, QI, Qdry, &
810  Rtot, CPtot, &
811  entr )
812  implicit none
813  integer, intent(in) :: IA, IS, IE
814  integer, intent(in) :: JA, JS, JE
815 
816  real(RP), intent(in) :: TEMP (IA,JA)
817  real(RP), intent(in) :: PRES (IA,JA)
818  real(RP), intent(in) :: QV (IA,JA)
819  real(RP), intent(in) :: QI (IA,JA)
820  real(RP), intent(in) :: Qdry (IA,JA)
821  real(RP), intent(in) :: Rtot (IA,JA)
822  real(RP), intent(in) :: CPtot(IA,JA)
823 
824  real(RP), intent(out) :: entr(IA,JA)
825 
826  integer :: i, j
827  !---------------------------------------------------------------------------
828 
829  !$omp parallel do OMP_SCHEDULE_ collapse(2)
830  !$acc kernels copyin(temp,pres,qv,qi,qdry,rtot,cptot) copyout(entr)
831  do j = js, je
832  do i = is, ie
833  call atmos_hydrometeor_entr_0d( &
834  temp(i,j), pres(i,j), & ! [IN]
835  qv(i,j), qi(i,j), qdry(i,j), & ! [IN]
836  rtot(i,j), cptot(i,j), & ! [IN]
837  entr(i,j) ) ! [OUT]
838  enddo
839  enddo
840  !$acc end kernels
841 
842  return
843  end subroutine atmos_hydrometeor_entr_2d
844  !-----------------------------------------------------------------------------
846  subroutine atmos_hydrometeor_entr_3d( &
847  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
848  TEMP, PRES, &
849  QV, QI, Qdry, &
850  Rtot, CPtot, &
851  entr )
852  implicit none
853  integer, intent(in) :: KA, KS, KE
854  integer, intent(in) :: IA, IS, IE
855  integer, intent(in) :: JA, JS, JE
856 
857  real(RP), intent(in) :: TEMP (KA,IA,JA)
858  real(RP), intent(in) :: PRES (KA,IA,JA)
859  real(RP), intent(in) :: QV (KA,IA,JA)
860  real(RP), intent(in) :: QI (KA,IA,JA)
861  real(RP), intent(in) :: Qdry (KA,IA,JA)
862  real(RP), intent(in) :: Rtot (KA,IA,JA)
863  real(RP), intent(in) :: CPtot(KA,IA,JA)
864 
865  real(RP), intent(out) :: entr(KA,IA,JA)
866 
867  integer :: k, i, j
868  !---------------------------------------------------------------------------
869 
870  !$omp parallel do OMP_SCHEDULE_ collapse(2)
871  !$acc kernels copyin(temp,pres,qv,qi,qdry,rtot,cptot) copyout(entr)
872  do j = js, je
873  do i = is, ie
874  do k = ks, ke
875  call atmos_hydrometeor_entr_0d( &
876  temp(k,i,j), pres(k,i,j), & ! [IN]
877  qv(k,i,j), qi(k,i,j), qdry(k,i,j), & ! [IN]
878  rtot(k,i,j), cptot(k,i,j), & ! [IN]
879  entr(k,i,j) ) ! [OUT]
880  enddo
881  enddo
882  enddo
883  !$acc end kernels
884 
885  return
886  end subroutine atmos_hydrometeor_entr_3d
887 
888  !-----------------------------------------------------------------------------
890  subroutine atmos_hydrometeor_entr2temp_0d( &
891  entr, pres, &
892  qv, qi, qdry, &
893  Rtot, CPtot, &
894  temp )
895  !$acc routine
896  use scale_const, only: &
897  pre00 => const_pre00, &
898  tem00 => const_tem00, &
899  rdry => const_rdry, &
900  rvap => const_rvap, &
901  lhv0 => const_lhv0, &
902  lhf0 => const_lhf0, &
903  psat0 => const_psat0
904  implicit none
905 
906  real(RP), intent(in) :: entr
907  real(RP), intent(in) :: pres
908  real(RP), intent(in) :: qv
909  real(RP), intent(in) :: qi
910  real(RP), intent(in) :: qdry
911  real(RP), intent(in) :: Rtot
912  real(RP), intent(in) :: CPtot
913 
914  real(RP), intent(out) :: temp
915 
916  real(RP) :: pres_dry, pres_vap
917  !---------------------------------------------------------------------------
918 
919  pres_dry = pres * qdry * rdry / rtot
920  pres_vap = pres * qv * rvap / rtot
921 
922  temp = tem00 &
923  * exp( ( entr &
924  + qdry * rdry * log( pres_dry / pre00 ) &
925  + qv * rvap * log( pres_vap / psat0 ) &
926  - ( qv * lhv0 - qi * lhf0 ) / tem00 ) / cptot )
927  return
928  end subroutine atmos_hydrometeor_entr2temp_0d
929 
930 end module scale_atmos_hydrometeor
scale_atmos_hydrometeor::num_name
character(len=h_short), dimension(n_hyd), parameter, public num_name
Definition: scale_atmos_hydrometeor.F90:108
scale_const::const_lhv0
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:82
scale_atmos_hydrometeor::hyd_desc
character(len=h_mid), dimension(n_hyd), parameter, public hyd_desc
Definition: scale_atmos_hydrometeor.F90:106
scale_atmos_hydrometeor::atmos_hydrometeor_lhf_1d
subroutine atmos_hydrometeor_lhf_1d(KA, KS, KE, temp, lhf)
Definition: scale_atmos_hydrometeor.F90:662
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:98
scale_atmos_hydrometeor::qia
integer, public qia
Definition: scale_atmos_hydrometeor.F90:139
scale_atmos_hydrometeor::cp_water
real(rp), public cp_water
CP for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:152
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
scale_atmos_hydrometeor::hyd_name
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
Definition: scale_atmos_hydrometeor.F90:104
scale_const::const_lhs00
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
Definition: scale_const.F90:85
scale_const::const_lhv00
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:83
scale_atmos_hydrometeor::qhs
integer, public qhs
Definition: scale_atmos_hydrometeor.F90:132
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_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_hydrometeor::atmos_hydrometeor_entr_2d
subroutine atmos_hydrometeor_entr_2d(IA, IS, IE, JA, JS, JE, TEMP, PRES, QV, QI, Qdry, Rtot, CPtot, entr)
calc temp, pres, q -> entropy (2D)
Definition: scale_atmos_hydrometeor.F90:812
scale_atmos_hydrometeor::i_hh
integer, parameter, public i_hh
hail
Definition: scale_atmos_hydrometeor.F90:102
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
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::qle
integer, public qle
Definition: scale_atmos_hydrometeor.F90:137
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_hydrometeor::cv_vapor
real(rp), public cv_vapor
CV for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:149
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_hydrometeor::qhe
integer, public qhe
Definition: scale_atmos_hydrometeor.F90:133
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_hydrometeor::atmos_hydrometeor_lhs_1d
subroutine atmos_hydrometeor_lhs_1d(KA, KS, KE, temp, lhs)
Definition: scale_atmos_hydrometeor.F90:537
scale_atmos_hydrometeor::qis
integer, public qis
Definition: scale_atmos_hydrometeor.F90:140
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_hydrometeor::qha
integer, public qha
Definition: scale_atmos_hydrometeor.F90:131
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_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_tracer::tracer_regist
subroutine, public tracer_regist(QS, NQ, NAME, DESC, UNIT, CV, CP, R, ENGI0, ADVC, MASS)
Regist tracer.
Definition: scale_tracer.F90:68
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
scale_const::const_ci
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:72
scale_const::const_dwatr
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:89
scale_atmos_hydrometeor::atmos_hydrometeor_ice_phase
logical, public atmos_hydrometeor_ice_phase
Definition: scale_atmos_hydrometeor.F90:90
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_atmos_hydrometeor::atmos_hydrometeor_lhv_0d
subroutine atmos_hydrometeor_lhv_0d(temp, lhv)
Definition: scale_atmos_hydrometeor.F90:392
scale_const::const_cvvap
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:70
scale_const::const_lhf0
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:86
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:146
scale_atmos_hydrometeor::atmos_hydrometeor_regist
subroutine, public atmos_hydrometeor_regist(NL, NI, NAME, DESC, UNIT, Q0, ADVC)
ATMOS_HYDROMETEOR_regist Regist tracer.
Definition: scale_atmos_hydrometeor.F90:279
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
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_hydrometeor::lhs
real(rp), public lhs
latent heat of sublimation for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:145
scale_atmos_hydrometeor::atmos_hydrometeor_finalize
subroutine, public atmos_hydrometeor_finalize
Finalize.
Definition: scale_atmos_hydrometeor.F90:255
scale_const::const_dice
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:90
scale_atmos_hydrometeor::qla
integer, public qla
Definition: scale_atmos_hydrometeor.F90:135
scale_atmos_hydrometeor::qie
integer, public qie
Definition: scale_atmos_hydrometeor.F90:141
scale_atmos_hydrometeor::qls
integer, public qls
Definition: scale_atmos_hydrometeor.F90:136
scale_atmos_hydrometeor::atmos_hydrometeor_setup
subroutine, public atmos_hydrometeor_setup
Setup.
Definition: scale_atmos_hydrometeor.F90:176
scale_atmos_hydrometeor::hyd_dens
real(rp), dimension(n_hyd), public hyd_dens
Definition: scale_atmos_hydrometeor.F90:111
scale_atmos_hydrometeor::atmos_hydrometeor_lhv_1d
subroutine atmos_hydrometeor_lhv_1d(KA, KS, KE, temp, lhv)
Definition: scale_atmos_hydrometeor.F90:412
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
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_const::const_setup
subroutine, public const_setup
Setup.
Definition: scale_const.F90:128
scale_atmos_hydrometeor::cp_ice
real(rp), public cp_ice
CP for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:154
scale_const::const_lhf00
real(rp), public const_lhf00
latent heat of fusion at 0K [J/kg]
Definition: scale_const.F90:87
scale_atmos_hydrometeor::i_hg
integer, parameter, public i_hg
graupel
Definition: scale_atmos_hydrometeor.F90:101