SCALE-RM
scale_atmos_sub_thermodyn.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
16 !-------------------------------------------------------------------------------
17 #include "macro_thermodyn.h"
18 #include "inc_openmp.h"
20  !-----------------------------------------------------------------------------
21  !
22  !++ used modules
23  !
24  use scale_precision
25  use scale_stdio
26  use scale_prof
28  use scale_tracer
29 
30  use scale_const, only: &
31  eps => const_eps, &
32  rdry => const_rdry, &
33  cpdry => const_cpdry, &
34  cpvap => const_cpvap, &
35  cvdry => const_cvdry, &
36  cl => const_cl, &
37  ci => const_ci, &
38  rvap => const_rvap, &
39  lhv0 => const_lhv0, &
40  lhs0 => const_lhs0, &
41  lhf0 => const_lhf0, &
42  psat0 => const_psat0, &
43  pre00 => const_pre00, &
44  tem00 => const_tem00
45  !-----------------------------------------------------------------------------
46  implicit none
47  private
48 
49  !-----------------------------------------------------------------------------
50  !
51  !++ Public procedure
52  !
53  public :: atmos_thermodyn_setup
54 
55  public :: atmos_thermodyn_qd
56  public :: atmos_thermodyn_cv
57  public :: atmos_thermodyn_cp
58  public :: atmos_thermodyn_pott
59  public :: atmos_thermodyn_rhoe
60  public :: atmos_thermodyn_rhot
61  public :: atmos_thermodyn_temp_pres
62  public :: atmos_thermodyn_temp_pres_e
63  public :: atmos_thermodyn_templhv
64  public :: atmos_thermodyn_templhs
65  public :: atmos_thermodyn_templhf
66  public :: atmos_thermodyn_entr
67 
68  interface atmos_thermodyn_qd
69  module procedure atmos_thermodyn_qd_0d
70  module procedure atmos_thermodyn_qd_3d
71  end interface atmos_thermodyn_qd
72 
73  interface atmos_thermodyn_cp
74  module procedure atmos_thermodyn_cp_0d
75  module procedure atmos_thermodyn_cp_3d
76  end interface atmos_thermodyn_cp
77 
78  interface atmos_thermodyn_cv
79  module procedure atmos_thermodyn_cv_0d
80  module procedure atmos_thermodyn_cv_3d
81  end interface atmos_thermodyn_cv
82 
83  interface atmos_thermodyn_rhoe
84  module procedure atmos_thermodyn_rhoe_0d
85  module procedure atmos_thermodyn_rhoe_3d
86  end interface atmos_thermodyn_rhoe
87 
88  interface atmos_thermodyn_rhot
89  module procedure atmos_thermodyn_rhot_0d
90  module procedure atmos_thermodyn_rhot_3d
91  end interface atmos_thermodyn_rhot
92 
93  interface atmos_thermodyn_pott
94  module procedure atmos_thermodyn_pott_0d
95  module procedure atmos_thermodyn_pott_1d
96  module procedure atmos_thermodyn_pott_3d
97  end interface atmos_thermodyn_pott
98 
99  interface atmos_thermodyn_temp_pres
100  module procedure atmos_thermodyn_temp_pres_0d
101  module procedure atmos_thermodyn_temp_pres_3d
102  end interface atmos_thermodyn_temp_pres
103 
104  interface atmos_thermodyn_temp_pres_e
105  module procedure atmos_thermodyn_temp_pres_e_0d
106  module procedure atmos_thermodyn_temp_pres_e_3d
107  end interface atmos_thermodyn_temp_pres_e
108 
109  interface atmos_thermodyn_templhv
110  module procedure atmos_thermodyn_templhv_0d
111  module procedure atmos_thermodyn_templhv_2d
112  module procedure atmos_thermodyn_templhv_3d
113  end interface atmos_thermodyn_templhv
114 
115  interface atmos_thermodyn_templhs
116  module procedure atmos_thermodyn_templhs_0d
117  module procedure atmos_thermodyn_templhs_2d
118  module procedure atmos_thermodyn_templhs_3d
119  end interface atmos_thermodyn_templhs
120 
121  interface atmos_thermodyn_templhf
122  module procedure atmos_thermodyn_templhf_0d
123  module procedure atmos_thermodyn_templhf_2d
124  module procedure atmos_thermodyn_templhf_3d
125  end interface atmos_thermodyn_templhf
126 
127  interface atmos_thermodyn_entr
128  module procedure atmos_thermodyn_entr_0d
129  module procedure atmos_thermodyn_entr_2d
130  module procedure atmos_thermodyn_entr_3d
131  end interface atmos_thermodyn_entr
132 
133  public :: atmos_thermodyn_tempre
134  public :: atmos_thermodyn_tempre2
135 
136  !-----------------------------------------------------------------------------
137  !
138  !++ Public parameters & variables
139  !
140  real(RP), public :: thermodyn_emask = 0.0_rp
141 
142  real(RP), public, allocatable :: aq_cp(:)
143  real(RP), public, allocatable :: aq_cv(:)
144 
145  !-----------------------------------------------------------------------------
146  !
147  !++ Private procedure
148  !
149  !-----------------------------------------------------------------------------
150  !
151  !++ Private parameters & variables
152  !
153 
154  !-----------------------------------------------------------------------------
155 contains
156  !-----------------------------------------------------------------------------
158  subroutine atmos_thermodyn_setup
159  use scale_process, only: &
161  use scale_const, only: &
162  cpvap => const_cpvap, &
163  cvvap => const_cvvap, &
164  cl => const_cl, &
165  ci => const_ci, &
166  thermodyn_type => const_thermodyn_type
167  implicit none
168 
169  integer :: n
170  !---------------------------------------------------------------------------
171 
172  if( io_l ) write(io_fid_log,*)
173  if( io_l ) write(io_fid_log,*) '++++++ Module[THERMODYN] / Categ[ATMOS SHARE] / Origin[SCALElib]'
174 
175  allocate( aq_cp(qqa) )
176  allocate( aq_cv(qqa) )
177 
178  if ( thermodyn_type == 'EXACT' ) then
179  thermodyn_emask = 1.0_rp
180 
181  aq_cp(i_qv) = cpvap
182  aq_cv(i_qv) = cvvap
183 
184  if ( qws /= 0 ) then
185  do n = qws, qwe
186  aq_cp(n) = cl
187  aq_cv(n) = cl
188  enddo
189  endif
190 
191  if ( qis /= 0 ) then
192  do n = qis, qie
193  aq_cp(n) = ci
194  aq_cv(n) = ci
195  enddo
196  endif
197  elseif( thermodyn_type == 'SIMPLE' ) then
198  thermodyn_emask = 0.0_rp
199 
200  aq_cp(i_qv) = cpdry
201  aq_cv(i_qv) = cvdry
202 
203  if ( qws /= 0 ) then
204  do n = qws, qwe
205  aq_cp(n) = cvdry
206  aq_cv(n) = cvdry
207  enddo
208  endif
209 
210  if ( qis /= 0 ) then
211  do n = qis, qie
212  aq_cp(n) = cvdry
213  aq_cv(n) = cvdry
214  enddo
215  endif
216  endif
217 
218  return
219  end subroutine atmos_thermodyn_setup
220 
221  !-----------------------------------------------------------------------------
223  subroutine atmos_thermodyn_qd_0d( &
224  qdry, &
225  q )
226  implicit none
227 
228  real(RP), intent(out) :: qdry
229  real(RP), intent(in) :: q(qa)
230 
231  integer :: iqw
232  !-----------------------------------------------------------------------------
233 
234  qdry = 1.0_rp
235 #ifndef DRY
236  do iqw = qqs, qqe
237  qdry = qdry - q(iqw)
238  enddo
239 #endif
240 
241  return
242  end subroutine atmos_thermodyn_qd_0d
243 
244  !-----------------------------------------------------------------------------
246  subroutine atmos_thermodyn_qd_3d( &
247  qdry, &
248  q )
249  implicit none
250 
251  real(RP), intent(out) :: qdry(ka,ia,ja)
252  real(RP), intent(in) :: q (ka,ia,ja,qa)
253 
254  integer :: k, i, j, iqw
255  !-----------------------------------------------------------------------------
256 
257  !$omp parallel do private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2)
258  do j = jsb, jeb
259  do i = isb, ieb
260  do k = ks, ke
261 
262 ! qdry(k,i,j) = 1.0_RP
263 ! do iqw = QQS, QQE
264 ! qdry(k,i,j) = qdry(k,i,j) - q(k,i,j,iqw)
265 ! enddo
266  calc_qdry( qdry(k,i,j), q, k, i, j, iqw )
267 
268  enddo
269  enddo
270  enddo
271 
272  return
273  end subroutine atmos_thermodyn_qd_3d
274 
275  !-----------------------------------------------------------------------------
277  subroutine atmos_thermodyn_cp_0d( &
278  CPtot, &
279  q, &
280  qdry )
281  implicit none
282 
283  real(RP), intent(out) :: CPtot
284  real(RP), intent(in) :: q(qa)
285  real(RP), intent(in) :: qdry
286 
287  integer :: iqw
288  !---------------------------------------------------------------------------
289 
290  cptot = qdry * cpdry
291 #ifndef DRY
292  do iqw = qqs, qqe
293  cptot = cptot + q(iqw) * aq_cp(iqw)
294  enddo
295 #endif
296 
297  return
298  end subroutine atmos_thermodyn_cp_0d
299 
300  !-----------------------------------------------------------------------------
302  subroutine atmos_thermodyn_cp_3d( &
303  CPtot, &
304  q, &
305  qdry )
306  implicit none
307 
308  real(RP), intent(out) :: CPtot(ka,ia,ja)
309  real(RP), intent(in) :: q (ka,ia,ja,qa)
310  real(RP), intent(in) :: qdry (ka,ia,ja)
311 
312  integer :: k, i, j, iqw
313  !---------------------------------------------------------------------------
314 
315  !$omp parallel do private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2)
316  do j = jsb, jeb
317  do i = isb, ieb
318  do k = ks, ke
319 
320 ! CPtot(k,i,j) = qdry(k,i,j) * CPdry
321 ! do iqw = QQS, QQE
322 ! CPtot(k,i,j) = CPtot(k,i,j) + q(k,i,j,iqw) * AQ_CP(iqw)
323 ! enddo
324 
325  calc_cp(cptot(k,i,j), qdry(k,i,j), q, k, i, j, iqw, cpdry, aq_cp)
326 
327  enddo
328  enddo
329  enddo
330 
331  return
332  end subroutine atmos_thermodyn_cp_3d
333 
334  !-----------------------------------------------------------------------------
336  subroutine atmos_thermodyn_cv_0d( &
337  CVtot, &
338  q, &
339  qdry )
340  implicit none
341 
342  real(RP), intent(out) :: CVtot
343  real(RP), intent(in) :: q(qa)
344  real(RP), intent(in) :: qdry
345 
346  integer :: iqw
347  !---------------------------------------------------------------------------
348 
349  cvtot = qdry * cvdry
350 #ifndef DRY
351  do iqw = qqs, qqe
352  cvtot = cvtot + q(iqw) * aq_cv(iqw)
353  enddo
354 #endif
355 
356  return
357  end subroutine atmos_thermodyn_cv_0d
358 
359  !-----------------------------------------------------------------------------
361  subroutine atmos_thermodyn_cv_3d( &
362  CVtot, &
363  q, &
364  qdry )
365  implicit none
366 
367  real(RP), intent(out) :: CVtot(ka,ia,ja)
368  real(RP), intent(in) :: q (ka,ia,ja,qa)
369  real(RP), intent(in) :: qdry (ka,ia,ja)
370 
371  integer :: k, i, j, iqw
372  !---------------------------------------------------------------------------
373 
374  !$omp parallel do private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2)
375  do j = jsb, jeb
376  do i = isb, ieb
377  do k = ks, ke
378 
379 ! CVtot(k,i,j) = qdry(k,i,j) * CVdry
380 ! do iqw = QQS, QQE
381 ! CVtot(k,i,j) = CVtot(k,i,j) + q(k,i,j,iqw) * AQ_CV(iqw)
382 ! enddo
383 
384  calc_cv(cvtot(k,i,j), qdry(k,i,j), q, k, i, j, iqw, cvdry, aq_cv)
385 
386  enddo
387  enddo
388  enddo
389 
390  return
391  end subroutine atmos_thermodyn_cv_3d
392 
393  !-----------------------------------------------------------------------------
394  subroutine atmos_thermodyn_pott_0d( &
395  pott, &
396  temp, &
397  pres, &
398  q )
399  implicit none
400 
401  real(RP), intent(out) :: pott ! potential temperature [K]
402  real(RP), intent(in) :: temp ! temperature [K]
403  real(RP), intent(in) :: pres ! pressure [Pa]
404  real(RP), intent(in) :: q(qa) ! mass concentration [kg/kg]
405 
406  real(RP) :: qdry
407  real(RP) :: Rtot, CVtot, RovCP
408 
409  integer :: iqw
410  !---------------------------------------------------------------------------
411 
412 #ifdef DRY
413  cvtot = cvdry
414  rtot = rdry
415 #else
416  qdry = 1.0_rp
417  cvtot = 0.0_rp
418  do iqw = qqs, qqe
419  qdry = qdry - q(iqw)
420  cvtot = cvtot + q(iqw) * aq_cv(iqw)
421  enddo
422  cvtot = cvdry * qdry + cvtot
423  rtot = rdry * qdry + rvap * q(i_qv)
424 #endif
425 
426  rovcp = rtot / ( cvtot + rtot )
427 
428  pott = temp * ( pre00 / pres )**rovcp
429 
430  return
431  end subroutine atmos_thermodyn_pott_0d
432 
433  !-----------------------------------------------------------------------------
434  subroutine atmos_thermodyn_pott_1d( &
435  pott, &
436  temp, &
437  pres, &
438  q )
439  implicit none
440 
441  real(RP), intent(out) :: pott(ka) ! potential temperature [K]
442  real(RP), intent(in) :: temp(ka) ! temperature [K]
443  real(RP), intent(in) :: pres(ka) ! pressure [Pa]
444  real(RP), intent(in) :: q (ka,qa) ! mass concentration [kg/kg]
445 
446  real(RP) :: qdry
447  real(RP) :: Rtot, CVtot, RovCP
448 
449  integer :: k, iqw
450  !---------------------------------------------------------------------------
451 
452  do k = ks, ke
453 #ifdef DRY
454  cvtot = cvdry
455  rtot = rdry
456 #else
457  qdry = 1.0_rp
458  cvtot = 0.0_rp
459  do iqw = qqs, qqe
460  qdry = qdry - q(k,iqw)
461  cvtot = cvtot + q(k,iqw) * aq_cv(iqw)
462  enddo
463  cvtot = cvdry * qdry + cvtot
464  rtot = rdry * qdry + rvap * q(k,i_qv)
465 #endif
466 
467  rovcp = rtot / ( cvtot + rtot )
468 
469  pott(k) = temp(k) * ( pre00 / pres(k) )**rovcp
470  enddo
471 
472  return
473  end subroutine atmos_thermodyn_pott_1d
474 
475  !-----------------------------------------------------------------------------
476  subroutine atmos_thermodyn_pott_3d( &
477  pott, &
478  temp, &
479  pres, &
480  q )
481  implicit none
482 
483  real(RP), intent(out) :: pott(ka,ia,ja) ! potential temperature [K]
484  real(RP), intent(in) :: temp(ka,ia,ja) ! temperature [K]
485  real(RP), intent(in) :: pres(ka,ia,ja) ! pressure [Pa]
486  real(RP), intent(in) :: q (ka,ia,ja,qa) ! mass concentration [kg/kg]
487 
488  real(RP) :: qdry
489  real(RP) :: Rtot, CVtot, RovCP
490 
491  integer :: k, i, j, iqw
492  !---------------------------------------------------------------------------
493 
494  do j = jsb, jeb
495  do i = isb, ieb
496  do k = ks, ke
497 #ifdef DRY
498  cvtot = cvdry
499  rtot = rdry
500 #else
501  qdry = 1.0_rp
502  cvtot = 0.0_rp
503  do iqw = qqs, qqe
504  qdry = qdry - q(k,i,j,iqw)
505  cvtot = cvtot + q(k,i,j,iqw) * aq_cv(iqw)
506  enddo
507  cvtot = cvdry * qdry + cvtot
508  rtot = rdry * qdry + rvap * q(k,i,j,i_qv)
509 #endif
510 
511  rovcp = rtot / ( cvtot + rtot )
512 
513  pott(k,i,j) = temp(k,i,j) * ( pre00 / pres(k,i,j) )**rovcp
514  enddo
515  enddo
516  enddo
517 
518  return
519  end subroutine atmos_thermodyn_pott_3d
520 
521  !-----------------------------------------------------------------------------
523  subroutine atmos_thermodyn_rhoe_0d( &
524  rhoe, &
525  rhot, &
526  q )
527  implicit none
528 
529  real(RP), intent(out) :: rhoe
530  real(RP), intent(in) :: rhot
531  real(RP), intent(in) :: q(qa)
532 
533  real(RP) :: qdry
534  real(RP) :: pres
535  real(RP) :: Rtot, CVtot, CPovCV
536 
537  integer :: iqw
538  !---------------------------------------------------------------------------
539 
540 #ifdef DRY
541  cvtot = cvdry
542  rtot = rdry
543 #else
544  qdry = 1.0_rp
545  cvtot = 0.0_rp
546  do iqw = qqs, qqe
547  qdry = qdry - q(iqw)
548  cvtot = cvtot + q(iqw) * aq_cv(iqw)
549  enddo
550  cvtot = cvdry * qdry + cvtot
551  rtot = rdry * qdry + rvap * q(i_qv)
552 #endif
553 
554  cpovcv = ( cvtot + rtot ) / cvtot
555 
556  pres = pre00 * ( rhot * rtot / pre00 )**cpovcv
557 
558  rhoe = pres / rtot * cvtot
559 
560  return
561  end subroutine atmos_thermodyn_rhoe_0d
562 
563  !-----------------------------------------------------------------------------
565  subroutine atmos_thermodyn_rhoe_3d( &
566  rhoe, &
567  rhot, &
568  q )
569  implicit none
570 
571  real(RP), intent(out) :: rhoe(ka,ia,ja)
572  real(RP), intent(in) :: rhot(ka,ia,ja)
573  real(RP), intent(in) :: q (ka,ia,ja,qa)
574 
575  real(RP) :: qdry
576  real(RP) :: pres
577  real(RP) :: Rtot, CVtot, CPovCV
578 
579  integer :: k, i, j, iqw
580  !---------------------------------------------------------------------------
581 
582  !$omp parallel do private(i,j,k,iqw,qdry,pres,Rtot,CVtot,CPovCV) OMP_SCHEDULE_ collapse(2)
583  do j = jsb, jeb
584  do i = isb, ieb
585  do k = ks, ke
586 #ifdef DRY
587  cvtot = cvdry
588  rtot = rdry
589 #else
590  qdry = 1.0_rp
591  cvtot = 0.0_rp
592  do iqw = qqs, qqe
593  qdry = qdry - q(k,i,j,iqw)
594  cvtot = cvtot + q(k,i,j,iqw) * aq_cv(iqw)
595  enddo
596  cvtot = cvdry * qdry + cvtot
597  rtot = rdry * qdry + rvap * q(k,i,j,i_qv)
598 #endif
599 
600 
601  cpovcv = ( cvtot + rtot ) / cvtot
602 
603  pres = pre00 * ( rhot(k,i,j) * rtot / pre00 )**cpovcv
604 
605  rhoe(k,i,j) = pres / rtot * cvtot
606  enddo
607  enddo
608  enddo
609 
610  return
611  end subroutine atmos_thermodyn_rhoe_3d
612 
613  !-----------------------------------------------------------------------------
615  subroutine atmos_thermodyn_rhot_0d( &
616  rhot, &
617  rhoe, &
618  q )
619  implicit none
620 
621  real(RP), intent(out) :: rhot
622  real(RP), intent(in) :: rhoe
623  real(RP), intent(in) :: q(qa)
624 
625  real(RP) :: qdry
626  real(RP) :: pres
627  real(RP) :: Rtot, CVtot, RovCP
628 
629  integer :: iqw
630  !---------------------------------------------------------------------------
631 
632 #ifdef DRY
633  cvtot = cvdry
634  rtot = rdry
635 #else
636  qdry = 1.0_rp
637  cvtot = 0.0_rp
638  do iqw = qqs, qqe
639  qdry = qdry - q(iqw)
640  cvtot = cvtot + q(iqw) * aq_cv(iqw)
641  enddo
642  cvtot = cvdry * qdry + cvtot
643  rtot = rdry * qdry + rvap * q(i_qv)
644 #endif
645 
646  rovcp = rtot / ( cvtot + rtot )
647 
648  pres = rhoe * rtot / cvtot
649 
650  rhot = rhoe / cvtot * ( pre00 / pres )**rovcp
651 
652  return
653  end subroutine atmos_thermodyn_rhot_0d
654 
655  !-----------------------------------------------------------------------------
657  subroutine atmos_thermodyn_rhot_3d( &
658  rhot, &
659  rhoe, &
660  q )
661  implicit none
662 
663  real(RP), intent(out) :: rhot(ka,ia,ja)
664  real(RP), intent(in) :: rhoe(ka,ia,ja)
665  real(RP), intent(in) :: q (ka,ia,ja,qa)
666 
667  real(RP) :: qdry
668  real(RP) :: pres
669  real(RP) :: Rtot, CVtot, RovCP
670 
671  integer :: k, i, j, iqw
672  !---------------------------------------------------------------------------
673 
674  !$omp parallel do private(i,j,k,iqw,qdry,pres,Rtot,CVtot,RovCP) OMP_SCHEDULE_ collapse(2)
675  do j = jsb, jeb
676  do i = isb, ieb
677  do k = ks, ke
678 
679 #ifdef DRY
680  cvtot = cvdry
681  rtot = rdry
682 #else
683  qdry = 1.0_rp
684  cvtot = 0.0_rp
685  do iqw = qqs, qqe
686  qdry = qdry - q(k,i,j,iqw)
687  cvtot = cvtot + q(k,i,j,iqw) * aq_cv(iqw)
688  enddo
689  cvtot = cvdry * qdry + cvtot
690  rtot = rdry * qdry + rvap * q(k,i,j,i_qv)
691 #endif
692 
693  rovcp = rtot / ( cvtot + rtot )
694 
695  pres = rhoe(k,i,j) * rtot / cvtot
696 
697  rhot(k,i,j) = rhoe(k,i,j) / cvtot * ( pre00 / pres )**rovcp
698  enddo
699  enddo
700  enddo
701 
702  return
703  end subroutine atmos_thermodyn_rhot_3d
704 
705  !-----------------------------------------------------------------------------
707  subroutine atmos_thermodyn_temp_pres_0d( &
708  temp, &
709  pres, &
710  dens, &
711  rhot, &
712  q )
713  implicit none
714 
715  real(RP), intent(out) :: temp
716  real(RP), intent(out) :: pres
717  real(RP), intent(in) :: dens
718  real(RP), intent(in) :: rhot
719  real(RP), intent(in) :: q(qa)
720 
721  real(RP) :: qdry
722  real(RP) :: Rtot, CVtot, CPovCV
723 
724  integer :: iqw
725  !---------------------------------------------------------------------------
726 
727 #ifdef DRY
728  cvtot = cvdry
729  rtot = rdry
730 #else
731  qdry = 1.0_rp
732  cvtot = 0.0_rp
733  do iqw = qqs, qqe
734  qdry = qdry - q(iqw)
735  cvtot = cvtot + q(iqw) * aq_cv(iqw)
736  enddo
737  cvtot = cvdry * qdry + cvtot
738  rtot = rdry * qdry + rvap * q(i_qv)
739 #endif
740 
741  cpovcv = ( cvtot + rtot ) / cvtot
742 
743  pres = pre00 * ( rhot * rtot / pre00 )**cpovcv
744  temp = pres / ( dens * rtot )
745 
746  return
747  end subroutine atmos_thermodyn_temp_pres_0d
748 
749  !-----------------------------------------------------------------------------
751  subroutine atmos_thermodyn_temp_pres_3d( &
752  temp, &
753  pres, &
754  dens, &
755  rhot, &
756  q )
757  implicit none
758 
759  real(RP), intent(out) :: temp(ka,ia,ja)
760  real(RP), intent(out) :: pres(ka,ia,ja)
761  real(RP), intent(in) :: dens(ka,ia,ja)
762  real(RP), intent(in) :: rhot(ka,ia,ja)
763  real(RP), intent(in) :: q (ka,ia,ja,qa)
764 
765  real(RP) :: qdry
766  real(RP) :: Rtot, CVtot, CPovCV
767 
768  integer :: k, i, j, iqw
769  !---------------------------------------------------------------------------
770 
771  !$omp parallel do private(i,j,k,iqw,qdry,Rtot,CVtot,CPovCV) OMP_SCHEDULE_ collapse(2)
772  do j = jsb, jeb
773  do i = isb, ieb
774  do k = ks, ke
775 #ifdef DRY
776  cvtot = cvdry
777  rtot = rdry
778 #else
779  qdry = 1.0_rp
780  cvtot = 0.0_rp
781  do iqw = qqs, qqe
782  qdry = qdry - q(k,i,j,iqw)
783  cvtot = cvtot + q(k,i,j,iqw) * aq_cv(iqw)
784  enddo
785  cvtot = cvdry * qdry + cvtot
786  rtot = rdry * qdry + rvap * q(k,i,j,i_qv)
787 #endif
788 
789  cpovcv = ( cvtot + rtot ) / cvtot
790 
791  pres(k,i,j) = pre00 * ( rhot(k,i,j) * rtot / pre00 )**cpovcv
792  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot )
793  enddo
794  enddo
795  enddo
796 
797  return
798  end subroutine atmos_thermodyn_temp_pres_3d
799 
800  !-----------------------------------------------------------------------------
802  subroutine atmos_thermodyn_temp_pres_e_0d( &
803  temp, &
804  pres, &
805  dens, &
806  rhoe, &
807  q )
808  implicit none
809 
810  real(RP), intent(out) :: temp
811  real(RP), intent(out) :: pres
812  real(RP), intent(in) :: dens
813  real(RP), intent(in) :: rhoe
814  real(RP), intent(in) :: q(qa)
815 
816  real(RP) :: qdry
817  real(RP) :: Rtot, CVtot
818 
819  integer :: iqw
820  !---------------------------------------------------------------------------
821 
822 #ifdef DRY
823  cvtot = cvdry
824  rtot = rdry
825 #else
826  qdry = 1.0_rp
827  cvtot = 0.0_rp
828  do iqw = qqs, qqe
829  qdry = qdry - q(iqw)
830  cvtot = cvtot + q(iqw) * aq_cv(iqw)
831  enddo
832  cvtot = cvdry * qdry + cvtot
833  rtot = rdry * qdry + rvap * q(i_qv)
834 #endif
835 
836  temp = rhoe / ( dens * cvtot )
837  pres = dens * rtot * temp
838 
839  return
840  end subroutine atmos_thermodyn_temp_pres_e_0d
841 
842  !-----------------------------------------------------------------------------
844  subroutine atmos_thermodyn_temp_pres_e_3d( &
845  temp, &
846  pres, &
847  dens, &
848  rhoe, &
849  q )
850  implicit none
851 
852  real(RP), intent(out) :: temp(ka,ia,ja)
853  real(RP), intent(out) :: pres(ka,ia,ja)
854  real(RP), intent(in) :: dens(ka,ia,ja)
855  real(RP), intent(in) :: rhoe(ka,ia,ja)
856  real(RP), intent(in) :: q (ka,ia,ja,qa)
857 
858  real(RP) :: qdry
859  real(RP) :: Rtot, CVtot
860 
861  integer :: k, i, j, iqw
862  !---------------------------------------------------------------------------
863 
864  !$omp parallel do private(i,j,k,iqw,qdry,Rtot,CVtot) OMP_SCHEDULE_ collapse(2)
865  do j = jsb, jeb
866  do i = isb, ieb
867  do k = ks, ke
868 #ifdef DRY
869  cvtot = cvdry
870  rtot = rdry
871 #else
872  qdry = 1.0_rp
873  cvtot = 0.0_rp
874  do iqw = qqs, qqe
875  qdry = qdry - q(k,i,j,iqw)
876  cvtot = cvtot + q(k,i,j,iqw) * aq_cv(iqw)
877  enddo
878  cvtot = cvdry * qdry + cvtot
879  rtot = rdry * qdry + rvap * q(k,i,j,i_qv)
880 #endif
881 
882  temp(k,i,j) = rhoe(k,i,j) / ( dens(k,i,j) * cvtot )
883  pres(k,i,j) = dens(k,i,j) * rtot * temp(k,i,j)
884  enddo
885  enddo
886  enddo
887 
888  return
889  end subroutine atmos_thermodyn_temp_pres_e_3d
890 
891  !-----------------------------------------------------------------------------
892  subroutine atmos_thermodyn_templhv_0d( &
893  lhv, &
894  temp )
895  implicit none
896 
897  real(RP), intent(out) :: lhv ! potential temperature [K]
898  real(RP), intent(in) :: temp ! temperature [K]
899  !---------------------------------------------------------------------------
900 
901  lhv = lhv0 + ( cpvap - cl ) * ( temp - tem00 ) * thermodyn_emask
902 
903  return
904  end subroutine atmos_thermodyn_templhv_0d
905 
906  !-----------------------------------------------------------------------------
907  subroutine atmos_thermodyn_templhv_2d( &
908  lhv, &
909  temp )
910  implicit none
911 
912  real(RP), intent(out) :: lhv (ia,ja) ! potential temperature [K]
913  real(RP), intent(in) :: temp(ia,ja) ! temperature [K]
914 
915  integer :: i, j
916  !---------------------------------------------------------------------------
917 
918  do j = jsb, jeb
919  do i = isb, ieb
920  lhv(i,j) = lhv0 + ( cpvap - cl ) * ( temp(i,j) - tem00 ) * thermodyn_emask
921  enddo
922  enddo
923 
924  return
925  end subroutine atmos_thermodyn_templhv_2d
926 
927  !-----------------------------------------------------------------------------
928  subroutine atmos_thermodyn_templhv_3d( &
929  lhv, &
930  temp )
931  implicit none
932 
933  real(RP), intent(out) :: lhv (ka,ia,ja) ! potential temperature [K]
934  real(RP), intent(in) :: temp(ka,ia,ja) ! temperature [K]
935 
936  integer :: k, i, j
937  !---------------------------------------------------------------------------
938 
939  do j = jsb, jeb
940  do i = isb, ieb
941  do k = ks, ke
942  lhv(k,i,j) = lhv0 + ( cpvap - cl ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
943  enddo
944  enddo
945  enddo
946 
947  return
948  end subroutine atmos_thermodyn_templhv_3d
949 
950  !-----------------------------------------------------------------------------
951  subroutine atmos_thermodyn_templhs_0d( &
952  lhs, &
953  temp )
954  implicit none
955 
956  real(RP), intent(out) :: lhs ! potential temperature [K]
957  real(RP), intent(in) :: temp ! temperature [K]
958  !---------------------------------------------------------------------------
959 
960  lhs = lhs0 + ( cpvap - ci ) * ( temp - tem00 ) * thermodyn_emask
961 
962  return
963  end subroutine atmos_thermodyn_templhs_0d
964 
965  !-----------------------------------------------------------------------------
966  subroutine atmos_thermodyn_templhs_2d( &
967  lhs, &
968  temp )
969  implicit none
970 
971  real(RP), intent(out) :: lhs (ia,ja) ! potential temperature [K]
972  real(RP), intent(in) :: temp(ia,ja) ! temperature [K]
973 
974  integer :: i, j
975  !---------------------------------------------------------------------------
976 
977  do j = jsb, jeb
978  do i = isb, ieb
979  lhs(i,j) = lhs0 + ( cpvap - ci ) * ( temp(i,j) - tem00 ) * thermodyn_emask
980  enddo
981  enddo
982 
983  return
984  end subroutine atmos_thermodyn_templhs_2d
985 
986  !-----------------------------------------------------------------------------
987  subroutine atmos_thermodyn_templhs_3d( &
988  lhs, &
989  temp )
990  implicit none
991 
992  real(RP), intent(out) :: lhs (ka,ia,ja) ! potential temperature [K]
993  real(RP), intent(in) :: temp(ka,ia,ja) ! temperature [K]
994 
995  integer :: k, i, j
996  !---------------------------------------------------------------------------
997 
998  do j = jsb, jeb
999  do i = isb, ieb
1000  do k = ks, ke
1001  lhs(k,i,j) = lhs0 + ( cpvap - ci ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
1002  enddo
1003  enddo
1004  enddo
1005 
1006  return
1007  end subroutine atmos_thermodyn_templhs_3d
1008 
1009  !-----------------------------------------------------------------------------
1010  subroutine atmos_thermodyn_templhf_0d( &
1011  lhf, &
1012  temp )
1013  implicit none
1014 
1015  real(RP), intent(out) :: lhf ! potential temperature [K]
1016  real(RP), intent(in) :: temp ! temperature [K]
1017  !---------------------------------------------------------------------------
1018 
1019  lhf = lhf0 + ( cl - ci ) * ( temp - tem00 ) * thermodyn_emask
1020 
1021  return
1022  end subroutine atmos_thermodyn_templhf_0d
1023 
1024  !-----------------------------------------------------------------------------
1025  subroutine atmos_thermodyn_templhf_2d( &
1026  lhf, &
1027  temp )
1028  implicit none
1029 
1030  real(RP), intent(out) :: lhf (ia,ja) ! potential temperature [K]
1031  real(RP), intent(in) :: temp(ia,ja) ! temperature [K]
1032 
1033  integer :: i, j
1034  !---------------------------------------------------------------------------
1035 
1036  do j = jsb, jeb
1037  do i = isb, ieb
1038  lhf(i,j) = lhf0 + ( cl - ci ) * ( temp(i,j) - tem00 ) * thermodyn_emask
1039  enddo
1040  enddo
1041 
1042  return
1043  end subroutine atmos_thermodyn_templhf_2d
1044 
1045  !-----------------------------------------------------------------------------
1046  subroutine atmos_thermodyn_templhf_3d( &
1047  lhf, &
1048  temp )
1049  implicit none
1050 
1051  real(RP), intent(out) :: lhf (ka,ia,ja) ! potential temperature [K]
1052  real(RP), intent(in) :: temp(ka,ia,ja) ! temperature [K]
1053 
1054  integer :: k, i, j
1055  !---------------------------------------------------------------------------
1056 
1057  do j = jsb, jeb
1058  do i = isb, ieb
1059  do k = ks, ke
1060  lhf(k,i,j) = lhf0 + ( cl - ci ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
1061  enddo
1062  enddo
1063  enddo
1064 
1065  return
1066  end subroutine atmos_thermodyn_templhf_3d
1067 
1068  !-----------------------------------------------------------------------------
1070  subroutine atmos_thermodyn_entr_0d( &
1071  entr, &
1072  temp, &
1073  pres, &
1074  q )
1075  implicit none
1076 
1077  real(RP), intent(out) :: entr
1078  real(RP), intent(in) :: temp
1079  real(RP), intent(in) :: pres
1080  real(RP), intent(in) :: q(qa)
1081 
1082  real(RP) :: lhv, lhf
1083  real(RP) :: qdry, Rtot
1084  real(RP) :: logT_T0, pres_dry, pres_vap
1085  real(RP) :: sq
1086 
1087  integer :: iqw
1088  !---------------------------------------------------------------------------
1089 
1090  lhv = lhv0 + ( cpvap - cl ) * ( temp - tem00 ) * thermodyn_emask
1091  lhf = lhf0 + ( cl - ci ) * ( temp - tem00 ) * thermodyn_emask
1092 
1093  qdry = 1.0_rp
1094  do iqw = qqs, qqe
1095  qdry = qdry - q(iqw)
1096  enddo
1097  rtot = rdry * qdry + rvap * q(i_qv)
1098 
1099  logt_t0 = log( temp / tem00 )
1100 
1101  ! dry air + vapor
1102  pres_dry = max( pres * qdry * rdry / rtot, eps )
1103  pres_vap = max( pres * q(i_qv) * rvap / rtot, eps )
1104 
1105  entr = qdry * cpdry * logt_t0 &
1106  - qdry * rdry * log( pres_dry / pre00 ) &
1107  + q(i_qv) * cpvap * logt_t0 &
1108  - q(i_qv) * rvap * log( pres_vap / psat0 ) &
1109  + q(i_qv) * lhv / tem00
1110 
1111  ! hydrometeors
1112  do iqw = qqs, qqe
1113  sq = 0.0_rp
1114 
1115  if( iqw == i_qc ) sq = q(i_qc) * cl * logt_t0
1116 
1117  if( iqw == i_qr ) sq = q(i_qr) * cl * logt_t0
1118 
1119  if( iqw == i_qi ) sq = q(i_qi) * ci * logt_t0 &
1120  - q(i_qi) * lhf / tem00
1121 
1122  if( iqw == i_qs ) sq = q(i_qs) * ci * logt_t0 &
1123  - q(i_qs) * lhf / tem00
1124 
1125  if( iqw == i_qg ) sq = q(i_qg) * ci * logt_t0 &
1126  - q(i_qg) * lhf / tem00
1127 
1128  entr = entr + sq
1129  enddo
1130 
1131  return
1132  end subroutine atmos_thermodyn_entr_0d
1133 
1134  !-----------------------------------------------------------------------------
1136  subroutine atmos_thermodyn_entr_2d( &
1137  entr, &
1138  temp, &
1139  pres, &
1140  q )
1141  implicit none
1142 
1143  real(RP), intent(out) :: entr(ia,ja)
1144  real(RP), intent(in) :: temp(ia,ja)
1145  real(RP), intent(in) :: pres(ia,ja)
1146  real(RP), intent(in) :: q (ia,ja,qa)
1147 
1148  real(RP) :: lhv, lhf
1149  real(RP) :: qdry, Rtot
1150  real(RP) :: logT_T0, pres_dry, pres_vap
1151  real(RP) :: sq
1152 
1153  integer :: i, j, iqw
1154  !---------------------------------------------------------------------------
1155 
1156  ! dry air + vapor
1157  do j = jsb, jeb
1158  do i = isb, ieb
1159  lhv = lhv0 + ( cpvap - cl ) * ( temp(i,j) - tem00 ) * thermodyn_emask
1160 
1161  qdry = 1.0_rp
1162  do iqw = qqs, qqe
1163  qdry = qdry - q(i,j,iqw)
1164  enddo
1165  rtot = rdry * qdry + rvap * q(i,j,i_qv)
1166 
1167  logt_t0 = log( temp(i,j) / tem00 )
1168 
1169  pres_dry = max( pres(i,j) * qdry * rdry / rtot, eps )
1170  pres_vap = max( pres(i,j) * q(i,j,i_qv) * rvap / rtot, eps )
1171 
1172  entr(i,j) = qdry * cpdry * logt_t0 &
1173  - qdry * rdry * log( pres_dry / pre00 ) &
1174  + q(i,j,i_qv) * cpvap * logt_t0 &
1175  - q(i,j,i_qv) * rvap * log( pres_vap / psat0 ) &
1176  + q(i,j,i_qv) * lhv / tem00
1177  enddo
1178  enddo
1179 
1180  ! hydrometeors
1181  do j = jsb, jeb
1182  do i = isb, ieb
1183  lhf = lhf0 + ( cl - ci ) * ( temp(i,j) - tem00 ) * thermodyn_emask
1184 
1185  logt_t0 = log( temp(i,j) / tem00 )
1186 
1187  do iqw = qqs, qqe
1188  sq = 0.0_rp
1189 
1190  if( iqw == i_qc ) sq = q(i,j,i_qc) * cl * logt_t0
1191  if( iqw == i_qr ) sq = q(i,j,i_qr) * cl * logt_t0
1192  if( iqw == i_qi ) sq = q(i,j,i_qi) * ( ci * logt_t0 - lhf / tem00 )
1193  if( iqw == i_qs ) sq = q(i,j,i_qs) * ( ci * logt_t0 - lhf / tem00 )
1194  if( iqw == i_qg ) sq = q(i,j,i_qg) * ( ci * logt_t0 - lhf / tem00 )
1195 
1196  entr(i,j) = entr(i,j) + sq
1197  enddo
1198  enddo
1199  enddo
1200 
1201  return
1202  end subroutine atmos_thermodyn_entr_2d
1203 
1204  !-----------------------------------------------------------------------------
1206  subroutine atmos_thermodyn_entr_3d( &
1207  entr, &
1208  temp, &
1209  pres, &
1210  q )
1211  implicit none
1212 
1213  real(RP), intent(out) :: entr(ka,ia,ja)
1214  real(RP), intent(in) :: temp(ka,ia,ja)
1215  real(RP), intent(in) :: pres(ka,ia,ja)
1216  real(RP), intent(in) :: q (ka,ia,ja,qa)
1217 
1218  real(RP) :: lhv, lhf
1219  real(RP) :: qdry, Rtot
1220  real(RP) :: logT_T0, pres_dry, pres_vap
1221  real(RP) :: sq
1222 
1223  integer :: k, i, j, iqw
1224  !---------------------------------------------------------------------------
1225 
1226  ! dry air + vapor
1227  do j = jsb, jeb
1228  do i = isb, ieb
1229  do k = ks, ke
1230  lhv = lhv0 + ( cpvap - cl ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
1231 
1232  qdry = 1.0_rp
1233  do iqw = qqs, qqe
1234  qdry = qdry - q(k,i,j,iqw)
1235  enddo
1236  rtot = rdry * qdry + rvap * q(k,i,j,i_qv)
1237 
1238  logt_t0 = log( temp(k,i,j) / tem00 )
1239 
1240  pres_dry = max( pres(k,i,j) * qdry * rdry / rtot, eps )
1241  pres_vap = max( pres(k,i,j) * q(k,i,j,i_qv) * rvap / rtot, eps )
1242 
1243  entr(k,i,j) = qdry * cpdry * logt_t0 &
1244  - qdry * rdry * log( pres_dry / pre00 ) &
1245  + q(k,i,j,i_qv) * cpvap * logt_t0 &
1246  - q(k,i,j,i_qv) * rvap * log( pres_vap / psat0 ) &
1247  + q(k,i,j,i_qv) * lhv / tem00
1248  enddo
1249  enddo
1250  enddo
1251 
1252  ! hydrometeors
1253  do j = jsb, jeb
1254  do i = isb, ieb
1255  do k = ks, ke
1256  lhf = lhf0 + ( cl - ci ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
1257 
1258  logt_t0 = log( temp(k,i,j) / tem00 )
1259 
1260  do iqw = qqs, qqe
1261  sq = 0.0_rp
1262 
1263  if( iqw == i_qc ) sq = q(k,i,j,i_qc) * cl * logt_t0
1264  if( iqw == i_qr ) sq = q(k,i,j,i_qr) * cl * logt_t0
1265  if( iqw == i_qi ) sq = q(k,i,j,i_qi) * ( ci * logt_t0 - lhf / tem00 )
1266  if( iqw == i_qs ) sq = q(k,i,j,i_qs) * ( ci * logt_t0 - lhf / tem00 )
1267  if( iqw == i_qg ) sq = q(k,i,j,i_qg) * ( ci * logt_t0 - lhf / tem00 )
1268 
1269  entr(k,i,j) = entr(k,i,j) + sq
1270  enddo
1271  enddo
1272  enddo
1273  enddo
1274 
1275  return
1276  end subroutine atmos_thermodyn_entr_3d
1277 
1278  !-----------------------------------------------------------------------------
1279  subroutine atmos_thermodyn_tempre( &
1280  temp, pres, &
1281  Ein, dens, qdry, q )
1282  implicit none
1283 
1284  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
1285  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
1286  real(RP), intent(in) :: Ein (ka,ia,ja) ! internal energy
1287  real(RP), intent(in) :: dens(ka,ia,ja) ! density
1288  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry concentration
1289  real(RP), intent(in) :: q (ka,ia,ja,qa) ! water concentration
1290 
1291  real(RP) :: cv, Rmoist
1292 
1293  integer :: i, j, k, iqw
1294  !---------------------------------------------------------------------------
1295 
1296  !$omp parallel do private(i,j,k,iqw,cv,Rmoist) OMP_SCHEDULE_ collapse(2)
1297  do j = jsb, jeb
1298  do i = isb, ieb
1299  do k = ks, ke
1300 
1301  calc_cv(cv, qdry(k,i,j), q, k, i, j, iqw, cvdry, aq_cv)
1302  calc_r(rmoist, q(k,i,j,i_qv), qdry(k,i,j), rdry, rvap)
1303 
1304  temp(k,i,j) = ein(k,i,j) / cv
1305 
1306  pres(k,i,j) = dens(k,i,j) * rmoist * temp(k,i,j)
1307 
1308  enddo
1309  enddo
1310  enddo
1311 
1312  return
1313  end subroutine atmos_thermodyn_tempre
1314 
1315  !-----------------------------------------------------------------------------
1316  subroutine atmos_thermodyn_tempre2( &
1317  temp, pres, &
1318  dens, pott, qdry, q )
1319  implicit none
1320 
1321  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
1322  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
1323  real(RP), intent(in) :: dens(ka,ia,ja) ! density
1324  real(RP), intent(in) :: pott(ka,ia,ja) ! potential temperature
1325  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry concentration
1326  real(RP), intent(in) :: q (ka,ia,ja,qa) ! water concentration
1327 
1328  real(RP) :: Rmoist, cp
1329 
1330  integer :: i, j, k, iqw
1331  !---------------------------------------------------------------------------
1332 
1333  !$omp parallel do private(i,j,k,iqw,cp,Rmoist) OMP_SCHEDULE_ collapse(2)
1334  do j = jsb, jeb
1335  do i = isb, ieb
1336  do k = ks, ke
1337 
1338  calc_cp(cp, qdry(k,i,j), q, k, i, j, iqw, cpdry, aq_cp)
1339  calc_r(rmoist, q(k,i,j,i_qv), qdry(k,i,j), rdry, rvap)
1340  calc_pre(pres(k,i,j), dens(k,i,j), pott(k,i,j), rmoist, cp, pre00)
1341 
1342  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rmoist )
1343 
1344  enddo
1345  enddo
1346  enddo
1347 
1348  return
1349  end subroutine atmos_thermodyn_tempre2
1350 
1351 end module scale_atmos_thermodyn
integer, public qie
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
subroutine, public prc_mpistop
Abort MPI.
integer, public jeb
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:69
real(rp), dimension(:), allocatable, public aq_cp
CP for each hydrometeors [J/kg/K].
integer, public qqe
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public qwe
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), public thermodyn_emask
=0: SIMPLE, 1: EXACT
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:95
integer, public qa
subroutine atmos_thermodyn_qd_0d(qdry, q)
calc dry air mass (0D)
integer, public qws
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:67
integer, public qis
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:84
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
subroutine, public atmos_thermodyn_tempre(temp, pres, Ein, dens, qdry, q)
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
integer, public i_qv
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:80
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
integer, public qqa
module PROCESS
subroutine, public atmos_thermodyn_setup
Setup.
subroutine, public log(type, message)
Definition: dc_log.f90:133
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:), allocatable, public aq_cv
CV for each hydrometeors [J/kg/K].
integer, public ks
start point of inner domain: z, local
integer, public i_qs
module profiler
Definition: scale_prof.F90:10
integer, public i_qi
real(rp), public const_eps
small number
Definition: scale_const.F90:36
module ATMOSPHERE / Thermodynamics
subroutine, public atmos_thermodyn_tempre2(temp, pres, dens, pott, qdry, q)
integer, public qqs
integer, public i_qg
module PRECISION
integer, public isb
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 i_qr
integer, public i_qc
integer, public ja
of y whole cells (local, with HALO)