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  rdry => const_rdry, &
32  cpdry => const_cpdry, &
33  cvdry => const_cvdry, &
34  pre00 => const_pre00
35  !-----------------------------------------------------------------------------
36  implicit none
37  private
38 
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public procedure
42  !
43  public :: atmos_thermodyn_setup
44 
45  public :: atmos_thermodyn_qd
46  public :: atmos_thermodyn_cv
47  public :: atmos_thermodyn_cp
48  public :: atmos_thermodyn_r
49  public :: atmos_thermodyn_pott
50  public :: atmos_thermodyn_rhoe
51  public :: atmos_thermodyn_rhot
52  public :: atmos_thermodyn_temp_pres
53  public :: atmos_thermodyn_temp_pres_e
54 
55  interface atmos_thermodyn_qd
56  module procedure atmos_thermodyn_qd_0d
57  module procedure atmos_thermodyn_qd_3d
58  end interface atmos_thermodyn_qd
59 
60  interface atmos_thermodyn_cv
61  module procedure atmos_thermodyn_cv_0d
62  module procedure atmos_thermodyn_cv_3d
63  end interface atmos_thermodyn_cv
64 
65  interface atmos_thermodyn_cp
66  module procedure atmos_thermodyn_cp_0d
67  module procedure atmos_thermodyn_cp_3d
68  end interface atmos_thermodyn_cp
69 
70  interface atmos_thermodyn_r
71  module procedure atmos_thermodyn_r_0d
72  module procedure atmos_thermodyn_r_3d
73  end interface atmos_thermodyn_r
74 
75  interface atmos_thermodyn_rhoe
76  module procedure atmos_thermodyn_rhoe_0d
77  module procedure atmos_thermodyn_rhoe_3d
78  end interface atmos_thermodyn_rhoe
79 
80  interface atmos_thermodyn_rhot
81  module procedure atmos_thermodyn_rhot_0d
82  module procedure atmos_thermodyn_rhot_3d
83  end interface atmos_thermodyn_rhot
84 
85  interface atmos_thermodyn_pott
86  module procedure atmos_thermodyn_pott_0d
87  module procedure atmos_thermodyn_pott_1d
88  module procedure atmos_thermodyn_pott_3d
89  end interface atmos_thermodyn_pott
90 
91  interface atmos_thermodyn_temp_pres
92  module procedure atmos_thermodyn_temp_pres_0d
93  module procedure atmos_thermodyn_temp_pres_3d
94  end interface atmos_thermodyn_temp_pres
95 
96  interface atmos_thermodyn_temp_pres_e
97  module procedure atmos_thermodyn_temp_pres_e_0d
98  module procedure atmos_thermodyn_temp_pres_e_3d
99  end interface atmos_thermodyn_temp_pres_e
100 
101  public :: atmos_thermodyn_tempre
102  public :: atmos_thermodyn_tempre2
103 
104  !-----------------------------------------------------------------------------
105  !
106  !++ Public parameters & variables
107  !
108  !-----------------------------------------------------------------------------
109  !
110  !++ Private procedure
111  !
112  !-----------------------------------------------------------------------------
113  !
114  !++ Private parameters & variables
115  !
116 
117  !-----------------------------------------------------------------------------
118 contains
119  !-----------------------------------------------------------------------------
121  subroutine atmos_thermodyn_setup
122  implicit none
123  !---------------------------------------------------------------------------
124 
125  if( io_l ) write(io_fid_log,*)
126  if( io_l ) write(io_fid_log,*) '++++++ Module[THERMODYN] / Categ[ATMOS SHARE] / Origin[SCALElib]'
127  if( io_l ) write(io_fid_log,*) '*** No namelists.'
128 
129  return
130  end subroutine atmos_thermodyn_setup
131 
132  !-----------------------------------------------------------------------------
134  subroutine atmos_thermodyn_qd_0d( &
135  qdry, &
136  q, &
137  q_mass )
138  implicit none
139 
140  real(RP), intent(out) :: qdry
141  real(RP), intent(in) :: q(qa)
142  real(RP), intent(in) :: q_mass(qa)
143 
144  integer :: iqw
145  !-----------------------------------------------------------------------------
146 
147  qdry = 1.0_rp
148 #ifndef DRY
149  do iqw = 1, qa
150  qdry = qdry - q(iqw)*q_mass(iqw)
151  enddo
152 #endif
153 
154  return
155  end subroutine atmos_thermodyn_qd_0d
156 
157  !-----------------------------------------------------------------------------
159  subroutine atmos_thermodyn_qd_3d( &
160  qdry, &
161  q, &
162  q_mass )
163  implicit none
164 
165  real(RP), intent(out) :: qdry (ka,ia,ja)
166  real(RP), intent(in) :: q (ka,ia,ja,qa)
167  real(RP), intent(in) :: q_mass(qa)
168 
169  integer :: k, i, j, iqw
170  !-----------------------------------------------------------------------------
171 
172  !$omp parallel do default(none) private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2) &
173  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,qdry,q,q_mass,QA)
174  do j = jsb, jeb
175  do i = isb, ieb
176  do k = ks, ke
177 ! qdry(k,i,j) = 1.0_RP
178 ! do iqw = 1, QA
179 ! qdry(k,i,j) = qdry(k,i,j) - q(k,i,j,iqw) * q_mass(iqw)
180 ! enddo
181  calc_qdry( qdry(k,i,j), q, q_mass, k, i, j, iqw )
182  enddo
183  enddo
184  enddo
185  return
186  end subroutine atmos_thermodyn_qd_3d
187 
188  !-----------------------------------------------------------------------------
190  subroutine atmos_thermodyn_cv_0d( &
191  CVtot, &
192  q, &
193  CVq, &
194  qdry )
195  implicit none
196 
197  real(RP), intent(out) :: cvtot
198  real(RP), intent(in) :: q(qa)
199  real(RP), intent(in) :: cvq(qa)
200  real(RP), intent(in) :: qdry
201 
202  integer :: iqw
203  !---------------------------------------------------------------------------
204 
205  cvtot = qdry * cvdry
206 #ifndef DRY
207  do iqw = 1, qa
208  cvtot = cvtot + q(iqw) * cvq(iqw)
209  enddo
210 #endif
211 
212  return
213  end subroutine atmos_thermodyn_cv_0d
214 
215  !-----------------------------------------------------------------------------
217  subroutine atmos_thermodyn_cv_3d( &
218  CVtot, &
219  q, &
220  CVq, &
221  qdry )
222  implicit none
223 
224  real(RP), intent(out) :: cvtot(ka,ia,ja)
225  real(RP), intent(in) :: q (ka,ia,ja,qa)
226  real(RP), intent(in) :: cvq (qa)
227  real(RP), intent(in) :: qdry (ka,ia,ja)
228 
229  integer :: k, i, j, iqw
230  !---------------------------------------------------------------------------
231 
232  !$omp parallel do default(none) private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2) &
233  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,cvtot,qdry,q,CVdry,CVq,QA)
234  do j = jsb, jeb
235  do i = isb, ieb
236  do k = ks, ke
237 ! CVtot(k,i,j) = qdry(k,i,j) * CVdry
238 ! do iqw = 1, QA
239 ! CVtot(k,i,j) = CVtot(k,i,j) + q(k,i,j,iqw) * CVq(iqw)
240 ! enddo
241  calc_cv(cvtot(k,i,j), qdry(k,i,j), q, k, i, j, iqw, cvdry, cvq)
242  enddo
243  enddo
244  enddo
245 
246  return
247  end subroutine atmos_thermodyn_cv_3d
248 
249  !-----------------------------------------------------------------------------
251  subroutine atmos_thermodyn_cp_0d( &
252  CPtot, &
253  q, &
254  CPq, &
255  qdry )
256  implicit none
257 
258  real(RP), intent(out) :: cptot
259  real(RP), intent(in) :: q(qa)
260  real(RP), intent(in) :: cpq(qa)
261  real(RP), intent(in) :: qdry
262 
263  integer :: iqw
264  !---------------------------------------------------------------------------
265 
266  cptot = qdry * cpdry
267 #ifndef DRY
268  do iqw = 1, qa
269  cptot = cptot + q(iqw) * cpq(iqw)
270  enddo
271 #endif
272 
273  return
274  end subroutine atmos_thermodyn_cp_0d
275 
276  !-----------------------------------------------------------------------------
278  subroutine atmos_thermodyn_cp_3d( &
279  CPtot, &
280  q, &
281  CPq, &
282  qdry )
283  implicit none
284 
285  real(RP), intent(out) :: cptot(ka,ia,ja)
286  real(RP), intent(in) :: q (ka,ia,ja,qa)
287  real(RP), intent(in) :: cpq (qa)
288  real(RP), intent(in) :: qdry (ka,ia,ja)
289 
290  integer :: k, i, j, iqw
291  !---------------------------------------------------------------------------
292 
293  !$omp parallel do private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2)
294  do j = jsb, jeb
295  do i = isb, ieb
296  do k = ks, ke
297 ! CPtot(k,i,j) = qdry(k,i,j) * CPdry
298 ! do iqw = 1, QA
299 ! CPtot(k,i,j) = CPtot(k,i,j) + q(k,i,j,iqw) * CPq(k,i,j,iqw)
300 ! enddo
301  calc_cp(cptot(k,i,j), qdry(k,i,j), q, k, i, j, iqw, cpdry, cpq)
302  enddo
303  enddo
304  enddo
305 
306  return
307  end subroutine atmos_thermodyn_cp_3d
308 
309  !-----------------------------------------------------------------------------
311  subroutine atmos_thermodyn_r_0d( &
312  Rtot, &
313  q, &
314  Rq, &
315  qdry )
316  implicit none
317 
318  real(RP), intent(out) :: rtot
319  real(RP), intent(in) :: q(qa)
320  real(RP), intent(in) :: rq(qa)
321  real(RP), intent(in) :: qdry
322 
323  integer :: iqw
324  !---------------------------------------------------------------------------
325 
326  rtot = qdry * rdry
327 #ifndef DRY
328  do iqw = 1, qa
329  rtot = rtot + q(iqw) * rq(iqw)
330  enddo
331 #endif
332 
333  return
334  end subroutine atmos_thermodyn_r_0d
335 
336  !-----------------------------------------------------------------------------
338  subroutine atmos_thermodyn_r_3d( &
339  Rtot, &
340  q, &
341  Rq, &
342  qdry )
343  implicit none
344 
345  real(RP), intent(out) :: rtot(ka,ia,ja)
346  real(RP), intent(in) :: q (ka,ia,ja,qa)
347  real(RP), intent(in) :: rq (qa)
348  real(RP), intent(in) :: qdry(ka,ia,ja)
349 
350  integer :: k, i, j, iqw
351  !---------------------------------------------------------------------------
352 
353  !$omp parallel do private(i,j,k,iqw) OMP_SCHEDULE_ collapse(2)
354  do j = jsb, jeb
355  do i = isb, ieb
356  do k = ks, ke
357  rtot(k,i,j) = qdry(k,i,j) * rdry
358  do iqw = 1, qa
359  rtot(k,i,j) = rtot(k,i,j) + q(k,i,j,iqw) * rq(iqw)
360  enddo
361  enddo
362  enddo
363  enddo
364 
365  return
366  end subroutine atmos_thermodyn_r_3d
367 
368  !-----------------------------------------------------------------------------
369  subroutine atmos_thermodyn_pott_0d( &
370  pott, &
371  temp, &
372  pres, &
373  q, &
374  CVq, &
375  Rq, &
376  mass )
377  implicit none
378 
379  real(RP), intent(out) :: pott ! potential temperature [K]
380  real(RP), intent(in) :: temp ! temperature [K]
381  real(RP), intent(in) :: pres ! pressure [Pa]
382  real(RP), intent(in) :: q(qa) ! mass concentration [kg/kg]
383  real(RP), intent(in) :: cvq (qa) ! specific heat [J/kg/K]
384  real(RP), intent(in) :: rq (qa) ! gas constant [J/kg/K]
385  real(RP), intent(in) :: mass(qa) !
386 
387  real(RP) :: qdry
388  real(RP) :: rtot, cvtot, rovcp
389 
390  integer :: iqw
391  !---------------------------------------------------------------------------
392 
393 #ifdef DRY
394  cvtot = cvdry
395  rtot = rdry
396 #else
397  qdry = 1.0_rp
398  cvtot = 0.0_rp
399  rtot = 0.0_rp
400  do iqw = 1, qa
401  qdry = qdry - q(iqw) * mass(iqw)
402  cvtot = cvtot + q(iqw) * cvq(iqw)
403  rtot = rtot + q(iqw) * rq(iqw)
404  enddo
405  cvtot = cvdry * qdry + cvtot
406  rtot = rdry * qdry + rtot
407 #endif
408  rovcp = rtot / ( cvtot + rtot )
409 
410  pott = temp * ( pre00 / pres )**rovcp
411 
412  return
413  end subroutine atmos_thermodyn_pott_0d
414 
415  !-----------------------------------------------------------------------------
416  subroutine atmos_thermodyn_pott_1d( &
417  pott, &
418  temp, &
419  pres, &
420  q, &
421  CVq, &
422  Rq, &
423  mass )
424  implicit none
425 
426  real(RP), intent(out) :: pott(ka) ! potential temperature [K]
427  real(RP), intent(in) :: temp(ka) ! temperature [K]
428  real(RP), intent(in) :: pres(ka) ! pressure [Pa]
429  real(RP), intent(in) :: q (ka,qa) ! mass concentration [kg/kg]
430  real(RP), intent(in) :: cvq (qa) ! specific heat [J/kg/K]
431  real(RP), intent(in) :: rq (qa) ! gas constant [J/kg/K]
432  real(RP), intent(in) :: mass(qa) !
433 
434  real(RP) :: qdry
435  real(RP) :: rtot, cvtot, rovcp
436 
437  integer :: k, iqw
438  !---------------------------------------------------------------------------
439 
440  do k = ks, ke
441 #ifdef DRY
442  cvtot = cvdry
443  rtot = rdry
444 #else
445  qdry = 1.0_rp
446  cvtot = 0.0_rp
447  rtot = 0.0_rp
448  do iqw = 1, qa
449  qdry = qdry - q(k,iqw) * mass(iqw)
450  cvtot = cvtot + q(k,iqw) * cvq(iqw)
451  rtot = rtot + q(k,iqw) * rq(iqw)
452  enddo
453  cvtot = cvdry * qdry + cvtot
454  rtot = rdry * qdry + rtot
455 #endif
456  rovcp = rtot / ( cvtot + rtot )
457 
458  pott(k) = temp(k) * ( pre00 / pres(k) )**rovcp
459  enddo
460 
461  return
462  end subroutine atmos_thermodyn_pott_1d
463 
464  !-----------------------------------------------------------------------------
465  subroutine atmos_thermodyn_pott_3d( &
466  pott, &
467  temp, &
468  pres, &
469  q, &
470  CVq, &
471  Rq, &
472  mass )
473  implicit none
474 
475  real(RP), intent(out) :: pott(ka,ia,ja) ! potential temperature [K]
476  real(RP), intent(in) :: temp(ka,ia,ja) ! temperature [K]
477  real(RP), intent(in) :: pres(ka,ia,ja) ! pressure [Pa]
478  real(RP), intent(in) :: q (ka,ia,ja,qa) ! mass concentration [kg/kg]
479  real(RP), intent(in) :: cvq (qa) ! specific heat [J/kg/K]
480  real(RP), intent(in) :: rq (qa) ! gas constant [J/kg/K]
481  real(RP), intent(in) :: mass(qa) !
482 
483  real(RP) :: qdry
484  real(RP) :: rtot, cvtot, rovcp
485 
486  integer :: k, i, j, iqw
487  !---------------------------------------------------------------------------
488 
489  do j = jsb, jeb
490  do i = isb, ieb
491  do k = ks, ke
492 #ifdef DRY
493  cvtot = cvdry
494  rtot = rdry
495 #else
496  qdry = 1.0_rp
497  cvtot = 0.0_rp
498  rtot = 0.0_rp
499  do iqw = 1, qa
500  qdry = qdry - q(k,i,j,iqw) * mass(iqw)
501  cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
502  rtot = rtot + q(k,i,j,iqw) * rq(iqw)
503  enddo
504  cvtot = cvdry * qdry + cvtot
505  rtot = rdry * qdry + rtot
506 #endif
507  rovcp = rtot / ( cvtot + rtot )
508 
509  pott(k,i,j) = temp(k,i,j) * ( pre00 / pres(k,i,j) )**rovcp
510  enddo
511  enddo
512  enddo
513 
514  return
515  end subroutine atmos_thermodyn_pott_3d
516 
517  !-----------------------------------------------------------------------------
519  subroutine atmos_thermodyn_rhoe_0d( &
520  rhoe, &
521  rhot, &
522  q, &
523  CVq, &
524  Rq, &
525  mass )
526  implicit none
527 
528  real(RP), intent(out) :: rhoe
529  real(RP), intent(in) :: rhot
530  real(RP), intent(in) :: q (qa)
531  real(RP), intent(in) :: cvq (qa)
532  real(RP), intent(in) :: rq (qa)
533  real(RP), intent(in) :: mass(qa)
534 
535  real(RP) :: qdry
536  real(RP) :: pres
537  real(RP) :: rtot, cvtot, cpovcv
538 
539  integer :: iqw
540  !---------------------------------------------------------------------------
541 
542 #ifdef DRY
543  cvtot = cvdry
544  rtot = rdry
545 #else
546  qdry = 1.0_rp
547  cvtot = 0.0_rp
548  rtot = 0.0_rp
549  do iqw = 1, qa
550  qdry = qdry - q(iqw) * mass(iqw)
551  cvtot = cvtot + q(iqw) * cvq(iqw)
552  rtot = rtot + q(iqw) * rq(iqw)
553  enddo
554  cvtot = cvdry * qdry + cvtot
555  rtot = rdry * qdry + rtot
556 #endif
557  cpovcv = ( cvtot + rtot ) / cvtot
558 
559  pres = pre00 * ( rhot * rtot / pre00 )**cpovcv
560 
561  rhoe = pres / rtot * cvtot
562 
563  return
564  end subroutine atmos_thermodyn_rhoe_0d
565 
566  !-----------------------------------------------------------------------------
568  subroutine atmos_thermodyn_rhoe_3d( &
569  rhoe, &
570  rhot, &
571  q, &
572  CVq, &
573  Rq, &
574  mass )
575  implicit none
576 
577  real(RP), intent(out) :: rhoe(ka,ia,ja)
578  real(RP), intent(in) :: rhot(ka,ia,ja)
579  real(RP), intent(in) :: q (ka,ia,ja,qa)
580  real(RP), intent(in) :: cvq (qa)
581  real(RP), intent(in) :: rq (qa)
582  real(RP), intent(in) :: mass(qa)
583 
584  real(RP) :: qdry
585  real(RP) :: pres
586  real(RP) :: rtot, cvtot, cpovcv
587 
588  integer :: k, i, j, iqw
589  !---------------------------------------------------------------------------
590 
591  !$omp parallel do default(none) private(i,j,k,iqw,qdry,pres,Rtot,CVtot,CPovCV) OMP_SCHEDULE_ collapse(2) &
592  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,CVdry,Rdry,QA,q,mass,CVq,Rq,rhoe,rhot,PRE00)
593  do j = jsb, jeb
594  do i = isb, ieb
595  do k = ks, ke
596 #ifdef DRY
597  cvtot = cvdry
598  rtot = rdry
599 #else
600  qdry = 1.0_rp
601  cvtot = 0.0_rp
602  rtot = 0.0_rp
603  do iqw = 1, qa
604  qdry = qdry - q(k,i,j,iqw) * mass(iqw)
605  cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
606  rtot = rtot + q(k,i,j,iqw) * rq(iqw)
607  enddo
608  cvtot = cvdry * qdry + cvtot
609  rtot = rdry * qdry + rtot
610 #endif
611  cpovcv = ( cvtot + rtot ) / cvtot
612 
613  pres = pre00 * ( rhot(k,i,j) * rtot / pre00 )**cpovcv
614 
615  rhoe(k,i,j) = pres / rtot * cvtot
616  enddo
617  enddo
618  enddo
619 
620  return
621  end subroutine atmos_thermodyn_rhoe_3d
622 
623  !-----------------------------------------------------------------------------
625  subroutine atmos_thermodyn_rhot_0d( &
626  rhot, &
627  rhoe, &
628  q, &
629  CVq, &
630  Rq, &
631  mass )
632  implicit none
633 
634  real(RP), intent(out) :: rhot
635  real(RP), intent(in) :: rhoe
636  real(RP), intent(in) :: q (qa)
637  real(RP), intent(in) :: cvq (qa)
638  real(RP), intent(in) :: rq (qa)
639  real(RP), intent(in) :: mass(qa)
640 
641  real(RP) :: qdry
642  real(RP) :: pres
643  real(RP) :: rtot, cvtot, rovcp
644 
645  integer :: iqw
646  !---------------------------------------------------------------------------
647 
648 #ifdef DRY
649  cvtot = cvdry
650  rtot = rdry
651 #else
652  qdry = 1.0_rp
653  cvtot = 0.0_rp
654  rtot = 0.0_rp
655  do iqw = 1, qa
656  qdry = qdry - q(iqw) * mass(iqw)
657  cvtot = cvtot + q(iqw) * cvq(iqw)
658  rtot = rtot + q(iqw) * rq(iqw)
659  enddo
660  cvtot = cvdry * qdry + cvtot
661  rtot = rdry * qdry + rtot
662 #endif
663  rovcp = rtot / ( cvtot + rtot )
664 
665  pres = rhoe * rtot / cvtot
666 
667  rhot = rhoe / cvtot * ( pre00 / pres )**rovcp
668 
669  return
670  end subroutine atmos_thermodyn_rhot_0d
671 
672  !-----------------------------------------------------------------------------
674  subroutine atmos_thermodyn_rhot_3d( &
675  rhot, &
676  rhoe, &
677  q, &
678  CVq, &
679  Rq, &
680  mass )
681  implicit none
682 
683  real(RP), intent(out) :: rhot(ka,ia,ja)
684  real(RP), intent(in) :: rhoe(ka,ia,ja)
685  real(RP), intent(in) :: q (ka,ia,ja,qa)
686  real(RP), intent(in) :: cvq (qa)
687  real(RP), intent(in) :: rq (qa)
688  real(RP), intent(in) :: mass(qa)
689 
690  real(RP) :: qdry
691  real(RP) :: pres
692  real(RP) :: rtot, cvtot, rovcp
693 
694  integer :: k, i, j, iqw
695  !---------------------------------------------------------------------------
696 
697  !$omp parallel do default(none) private(i,j,k,iqw,qdry,pres,Rtot,CVtot,RovCP) OMP_SCHEDULE_ collapse(2) &
698  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,CVdry,Rdry,QA,q,mass,CVq,Rq,rhoe,rhot,PRE00)
699  do j = jsb, jeb
700  do i = isb, ieb
701  do k = ks, ke
702 
703 #ifdef DRY
704  cvtot = cvdry
705  rtot = rdry
706 #else
707  qdry = 1.0_rp
708  cvtot = 0.0_rp
709  rtot = 0.0_rp
710  do iqw = 1, qa
711  qdry = qdry - q(k,i,j,iqw) * mass(iqw)
712  cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
713  rtot = rtot + q(k,i,j,iqw) * rq(iqw)
714  enddo
715  cvtot = cvdry * qdry + cvtot
716  rtot = rdry * qdry + rtot
717 #endif
718  rovcp = rtot / ( cvtot + rtot )
719 
720  pres = rhoe(k,i,j) * rtot / cvtot
721 
722  rhot(k,i,j) = rhoe(k,i,j) / cvtot * ( pre00 / pres )**rovcp
723  enddo
724  enddo
725  enddo
726  return
727  end subroutine atmos_thermodyn_rhot_3d
728 
729  !-----------------------------------------------------------------------------
731  subroutine atmos_thermodyn_temp_pres_0d( &
732  temp, &
733  pres, &
734  dens, &
735  rhot, &
736  q, &
737  CVq, &
738  Rq, &
739  mass )
740  implicit none
741 
742  real(RP), intent(out) :: temp
743  real(RP), intent(out) :: pres
744  real(RP), intent(in) :: dens
745  real(RP), intent(in) :: rhot
746  real(RP), intent(in) :: q (qa)
747  real(RP), intent(in) :: cvq (qa)
748  real(RP), intent(in) :: rq (qa)
749  real(RP), intent(in) :: mass(qa)
750 
751  real(RP) :: qdry
752  real(RP) :: rtot, cvtot, cpovcv
753 
754  integer :: iqw
755  !---------------------------------------------------------------------------
756 
757 #ifdef DRY
758  cvtot = cvdry
759  rtot = rdry
760 #else
761  qdry = 1.0_rp
762  cvtot = 0.0_rp
763  rtot = 0.0_rp
764  do iqw = 1, qa
765  qdry = qdry - q(iqw) * mass(iqw)
766  cvtot = cvtot + q(iqw) * cvq(iqw)
767  rtot = rtot + q(iqw) * rq(iqw)
768  enddo
769  cvtot = cvdry * qdry + cvtot
770  rtot = rdry * qdry + rtot
771 #endif
772  cpovcv = ( cvtot + rtot ) / cvtot
773 
774  pres = pre00 * ( rhot * rtot / pre00 )**cpovcv
775  temp = pres / ( dens * rtot )
776 
777  return
778  end subroutine atmos_thermodyn_temp_pres_0d
779 
780  !-----------------------------------------------------------------------------
782  subroutine atmos_thermodyn_temp_pres_3d( &
783  temp, &
784  pres, &
785  dens, &
786  rhot, &
787  q, &
788  CVq, &
789  Rq, &
790  mass )
791  implicit none
792 
793  real(RP), intent(out) :: temp(ka,ia,ja)
794  real(RP), intent(out) :: pres(ka,ia,ja)
795  real(RP), intent(in) :: dens(ka,ia,ja)
796  real(RP), intent(in) :: rhot(ka,ia,ja)
797  real(RP), intent(in) :: q (ka,ia,ja,qa)
798  real(RP), intent(in) :: cvq (qa)
799  real(RP), intent(in) :: rq (qa)
800  real(RP), intent(in) :: mass(qa)
801 
802  real(RP) :: qdry
803  real(RP) :: rtot, cvtot, cpovcv
804 
805  integer :: k, i, j, iqw
806  !---------------------------------------------------------------------------
807 
808  !$omp parallel do default(none) private(i,j,k,iqw,qdry,Rtot,CVtot,CPovCV) OMP_SCHEDULE_ collapse(2) &
809  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,CVdry,Rdry,QA,q,mass,CVq,Rq,temp,rhot,dens,pres,PRE00)
810  do j = jsb, jeb
811  do i = isb, ieb
812  do k = ks, ke
813 #ifdef DRY
814  cvtot = cvdry
815  rtot = rdry
816 #else
817  qdry = 1.0_rp
818  cvtot = 0.0_rp
819  rtot = 0.0_rp
820  do iqw = 1, qa
821  qdry = qdry - q(k,i,j,iqw) * mass(iqw)
822  cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
823  rtot = rtot + q(k,i,j,iqw) * rq(iqw)
824  enddo
825  cvtot = cvdry * qdry + cvtot
826  rtot = rdry * qdry + rtot
827 #endif
828  cpovcv = ( cvtot + rtot ) / cvtot
829 
830  pres(k,i,j) = pre00 * ( rhot(k,i,j) * rtot / pre00 )**cpovcv
831  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot )
832  enddo
833  enddo
834  enddo
835  return
836  end subroutine atmos_thermodyn_temp_pres_3d
837 
838  !-----------------------------------------------------------------------------
840  subroutine atmos_thermodyn_temp_pres_e_0d( &
841  temp, &
842  pres, &
843  dens, &
844  rhoe, &
845  q, &
846  CVq, &
847  Rq, &
848  mass )
849  implicit none
850 
851  real(RP), intent(out) :: temp
852  real(RP), intent(out) :: pres
853  real(RP), intent(in) :: dens
854  real(RP), intent(in) :: rhoe
855  real(RP), intent(in) :: q (qa)
856  real(RP), intent(in) :: cvq (qa)
857  real(RP), intent(in) :: rq (qa)
858  real(RP), intent(in) :: mass(qa)
859 
860  real(RP) :: qdry
861  real(RP) :: rtot, cvtot
862 
863  integer :: iqw
864  !---------------------------------------------------------------------------
865 
866 #ifdef DRY
867  cvtot = cvdry
868  rtot = rdry
869 #else
870  qdry = 1.0_rp
871  cvtot = 0.0_rp
872  rtot = 0.0_rp
873  do iqw = 1, qa
874  qdry = qdry - q(iqw) * mass(iqw)
875  cvtot = cvtot + q(iqw) * cvq(iqw)
876  rtot = rtot + q(iqw) * rq(iqw)
877  enddo
878  cvtot = cvdry * qdry + cvtot
879  rtot = rdry * qdry + rtot
880 #endif
881 
882  temp = rhoe / ( dens * cvtot )
883  pres = dens * rtot * temp
884 
885  return
886  end subroutine atmos_thermodyn_temp_pres_e_0d
887 
888  !-----------------------------------------------------------------------------
890  subroutine atmos_thermodyn_temp_pres_e_3d( &
891  temp, &
892  pres, &
893  dens, &
894  rhoe, &
895  q, &
896  CVq, &
897  Rq, &
898  mass )
899  implicit none
900 
901  real(RP), intent(out) :: temp(ka,ia,ja)
902  real(RP), intent(out) :: pres(ka,ia,ja)
903  real(RP), intent(in) :: dens(ka,ia,ja)
904  real(RP), intent(in) :: rhoe(ka,ia,ja)
905  real(RP), intent(in) :: q (ka,ia,ja,qa)
906  real(RP), intent(in) :: cvq (qa)
907  real(RP), intent(in) :: rq (qa)
908  real(RP), intent(in) :: mass(qa)
909 
910  real(RP) :: qdry
911  real(RP) :: rtot, cvtot
912 
913  integer :: k, i, j, iqw
914  !---------------------------------------------------------------------------
915 
916  !$omp parallel do default(none) private(i,j,k,iqw,qdry,Rtot,CVtot) OMP_SCHEDULE_ collapse(2) &
917  !$omp shared(JSB,JEB,ISB,IEB,KS,KE,CVdry,Rdry,QA,q,mass,CVq,Rq,temp,rhoe,dens,pres)
918  do j = jsb, jeb
919  do i = isb, ieb
920  do k = ks, ke
921 #ifdef DRY
922  cvtot = cvdry
923  rtot = rdry
924 #else
925  qdry = 1.0_rp
926  cvtot = 0.0_rp
927  rtot = 0.0_rp
928  do iqw = 1, qa
929  qdry = qdry - q(k,i,j,iqw) * mass(iqw)
930  cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
931  rtot = rtot + q(k,i,j,iqw) * rq(iqw)
932  enddo
933  cvtot = cvdry * qdry + cvtot
934  rtot = rdry * qdry + rtot
935 #endif
936 
937  temp(k,i,j) = rhoe(k,i,j) / ( dens(k,i,j) * cvtot )
938  pres(k,i,j) = dens(k,i,j) * rtot * temp(k,i,j)
939  enddo
940  enddo
941  enddo
942 
943  return
944  end subroutine atmos_thermodyn_temp_pres_e_3d
945 
946  !-----------------------------------------------------------------------------
947  subroutine atmos_thermodyn_tempre( &
948  temp, pres, &
949  Ein, dens, qdry, q, &
950  CVq, Rq )
951  implicit none
952 
953  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
954  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
955  real(RP), intent(in) :: ein (ka,ia,ja) ! internal energy
956  real(RP), intent(in) :: dens(ka,ia,ja) ! density
957  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry concentration
958  real(RP), intent(in) :: q (ka,ia,ja,qa) ! water concentration
959  real(RP), intent(in) :: cvq (qa) ! specific heat
960  real(RP), intent(in) :: rq (qa) ! gas constant
961 
962  real(RP) :: cv, rmoist
963 
964  integer :: i, j, k, iqw
965  !---------------------------------------------------------------------------
966 
967  !$omp parallel do private(i,j,k,iqw,cv,Rmoist) OMP_SCHEDULE_ collapse(2)
968  do j = jsb, jeb
969  do i = isb, ieb
970  do k = ks, ke
971  cv = qdry(k,i,j) * cvdry
972  rmoist = qdry(k,i,j) * rdry
973  do iqw =1, qa
974  cv = cv + q(k,i,j,iqw) * cvq(iqw)
975  rmoist = rmoist + q(k,i,j,iqw) * rq(iqw)
976  enddo
977 
978  temp(k,i,j) = ein(k,i,j) / cv
979  pres(k,i,j) = dens(k,i,j) * rmoist * temp(k,i,j)
980  enddo
981  enddo
982  enddo
983 
984  return
985  end subroutine atmos_thermodyn_tempre
986 
987  !-----------------------------------------------------------------------------
988  subroutine atmos_thermodyn_tempre2( &
989  temp, pres, &
990  dens, pott, qdry, q, &
991  CPq, Rq )
992  implicit none
993 
994  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
995  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
996  real(RP), intent(in) :: dens(ka,ia,ja) ! density
997  real(RP), intent(in) :: pott(ka,ia,ja) ! potential temperature
998  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry concentration
999  real(RP), intent(in) :: q (ka,ia,ja,qa) ! water concentration
1000  real(RP), intent(in) :: cpq (qa) ! specific heat
1001  real(RP), intent(in) :: rq (qa) ! gas constant
1002 
1003  real(RP) :: rmoist, cp
1004 
1005  integer :: i, j, k, iqw
1006  !---------------------------------------------------------------------------
1007 
1008  !$omp parallel do private(i,j,k,iqw,cp,Rmoist) OMP_SCHEDULE_ collapse(2)
1009  do j = jsb, jeb
1010  do i = isb, ieb
1011  do k = ks, ke
1012  cp = qdry(k,i,j) * cpdry
1013  rmoist = qdry(k,i,j) * rdry
1014  do iqw = 1, qa
1015  cp = cp + q(k,i,j,iqw) * cpq(iqw)
1016  rmoist = rmoist + q(k,i,j,iqw) * rq(iqw)
1017  enddo
1018  calc_pre(pres(k,i,j), dens(k,i,j), pott(k,i,j), rmoist, cp, pre00)
1019 
1020  temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rmoist )
1021  enddo
1022  enddo
1023  enddo
1024 
1025  return
1026  end subroutine atmos_thermodyn_tempre2
1027 
1028 end module scale_atmos_thermodyn
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
integer, public jeb
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
subroutine, public atmos_thermodyn_tempre2(temp, pres, dens, pott, qdry, q, CPq, Rq)
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
integer, public qa
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
integer, public ieb
module grid index
module TRACER
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
subroutine, public atmos_thermodyn_setup
Setup.
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
module ATMOSPHERE / Thermodynamics
module PRECISION
integer, public isb
subroutine, public atmos_thermodyn_tempre(temp, pres, Ein, dens, qdry, q, CVq, Rq)
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public jsb
integer, public ja
of whole cells: y, local, with HALO