SCALE-RM
scale_atmos_thermodyn.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
19 
20  use scale_const, only: &
21  rdry => const_rdry, &
22  cpdry => const_cpdry, &
23  cvdry => const_cvdry, &
24  pre00 => const_pre00
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: atmos_thermodyn_setup
34 
35  public :: atmos_thermodyn_qdry
36  public :: atmos_thermodyn_cv
37  public :: atmos_thermodyn_cp
38  public :: atmos_thermodyn_r
39  public :: atmos_thermodyn_specific_heat
40  public :: atmos_thermodyn_rhot2pres
41  public :: atmos_thermodyn_rhot2rhoe
42  public :: atmos_thermodyn_rhoe2rhot
43  public :: atmos_thermodyn_rhot2temp_pres
44  public :: atmos_thermodyn_rhoe2temp_pres
45  public :: atmos_thermodyn_ein2temp_pres
46  public :: atmos_thermodyn_pott2temp_pres
47  public :: atmos_thermodyn_temp_pres2pott
48  public :: atmos_thermodyn_temp_pres2ein
49  public :: atmos_thermodyn_ein_pres2enth
50 
51  interface atmos_thermodyn_qdry
52  module procedure atmos_thermodyn_qdry_0d
53  module procedure atmos_thermodyn_qdry_2d
54  module procedure atmos_thermodyn_qdry_3d
55  end interface atmos_thermodyn_qdry
56 
57  interface atmos_thermodyn_cv
58  module procedure atmos_thermodyn_cv_0d
59  module procedure atmos_thermodyn_cv_3d
60  end interface atmos_thermodyn_cv
61 
62  interface atmos_thermodyn_cp
63  module procedure atmos_thermodyn_cp_0d
64  module procedure atmos_thermodyn_cp_3d
65  end interface atmos_thermodyn_cp
66 
67  interface atmos_thermodyn_r
68  module procedure atmos_thermodyn_r_0d
69  module procedure atmos_thermodyn_r_3d
70  end interface atmos_thermodyn_r
71 
72  interface atmos_thermodyn_specific_heat
73  module procedure atmos_thermodyn_specific_heat_0d
74  module procedure atmos_thermodyn_specific_heat_1d
75  module procedure atmos_thermodyn_specific_heat_2d
76  module procedure atmos_thermodyn_specific_heat_3d
77  end interface atmos_thermodyn_specific_heat
78 
79  interface atmos_thermodyn_rhot2pres
80  module procedure atmos_thermodyn_rhot2pres_0d
81  module procedure atmos_thermodyn_rhot2pres_3d
82  end interface atmos_thermodyn_rhot2pres
83 
84  interface atmos_thermodyn_rhot2rhoe
85  module procedure atmos_thermodyn_rhot2rhoe_0d
86  module procedure atmos_thermodyn_rhot2rhoe_3d
87  end interface atmos_thermodyn_rhot2rhoe
88 
89  interface atmos_thermodyn_rhoe2rhot
90  module procedure atmos_thermodyn_rhoe2rhot_0d
91  module procedure atmos_thermodyn_rhoe2rhot_3d
92  end interface atmos_thermodyn_rhoe2rhot
93 
94  interface atmos_thermodyn_rhoe2temp_pres
95  module procedure atmos_thermodyn_rhoe2temp_pres_0d
96  module procedure atmos_thermodyn_rhoe2temp_pres_2d
97  module procedure atmos_thermodyn_rhoe2temp_pres_3d
98  end interface atmos_thermodyn_rhoe2temp_pres
99 
100  interface atmos_thermodyn_rhot2temp_pres
101  module procedure atmos_thermodyn_rhot2temp_pres_0d
102  module procedure atmos_thermodyn_rhot2temp_pres_1d
103  module procedure atmos_thermodyn_rhot2temp_pres_3d
104  end interface atmos_thermodyn_rhot2temp_pres
105 
106  interface atmos_thermodyn_ein2temp_pres
107  module procedure atmos_thermodyn_ein2temp_pres_0d
108  module procedure atmos_thermodyn_ein2temp_pres_2d
109  module procedure atmos_thermodyn_ein2temp_pres_3d
110  end interface atmos_thermodyn_ein2temp_pres
111 
112  interface atmos_thermodyn_pott2temp_pres
113  module procedure atmos_thermodyn_pott2temp_pres_0d
114  module procedure atmos_thermodyn_pott2temp_pres_3d
115  end interface atmos_thermodyn_pott2temp_pres
116 
117  interface atmos_thermodyn_temp_pres2pott
118  module procedure atmos_thermodyn_temp_pres2pott_0d
119  module procedure atmos_thermodyn_temp_pres2pott_1d
120  module procedure atmos_thermodyn_temp_pres2pott_2d
121  module procedure atmos_thermodyn_temp_pres2pott_3d
122  end interface atmos_thermodyn_temp_pres2pott
123 
124  interface atmos_thermodyn_temp_pres2ein
125  module procedure atmos_thermodyn_temp_pres2ein_0d
126  module procedure atmos_thermodyn_temp_pres2ein_1d
127  module procedure atmos_thermodyn_temp_pres2ein_2d
128  module procedure atmos_thermodyn_temp_pres2ein_3d
129  end interface atmos_thermodyn_temp_pres2ein
130 
131  interface atmos_thermodyn_ein_pres2enth
132  module procedure atmos_thermodyn_ein_pres2enth_0d
133  module procedure atmos_thermodyn_ein_pres2enth_1d
134  module procedure atmos_thermodyn_ein_pres2enth_2d
135  module procedure atmos_thermodyn_ein_pres2enth_3d
136  end interface atmos_thermodyn_ein_pres2enth
137 
138  !-----------------------------------------------------------------------------
139  !
140  !++ Public parameters & variables
141  !
142  !-----------------------------------------------------------------------------
143  !
144  !++ Private procedure
145  !
146  !-----------------------------------------------------------------------------
147  !
148  !++ Private parameters & variables
149  !
150 
151  !-----------------------------------------------------------------------------
152 contains
153  !-----------------------------------------------------------------------------
155  subroutine atmos_thermodyn_setup
156  implicit none
157  !---------------------------------------------------------------------------
158 
159  log_newline
160  log_info("ATMOS_THERMODYN_setup",*) 'Setup'
161  log_info("ATMOS_THERMODYN_setup",*) 'No namelists.'
162 
163  return
164  end subroutine atmos_thermodyn_setup
165 
166  !-----------------------------------------------------------------------------
168 !OCL SERIAL
169 !OCL NOSIMD
170  subroutine atmos_thermodyn_qdry_0d( &
171  QA, &
172  q, q_mass, &
173  qdry )
174  implicit none
175  integer, intent(in) :: QA
176 
177  real(RP), intent(in) :: q(qa)
178  real(RP), intent(in) :: q_mass(qa)
179 
180  real(RP), intent(out) :: qdry
181 
182  integer :: iqw
183  !---------------------------------------------------------------------------
184 
185  qdry = 1.0_rp
186 #ifndef DRY
187  do iqw = 1, qa
188  qdry = qdry - q(iqw)*q_mass(iqw)
189  enddo
190 #endif
191 
192  return
193  end subroutine atmos_thermodyn_qdry_0d
194 
195  !-----------------------------------------------------------------------------
197  subroutine atmos_thermodyn_qdry_2d( &
198  IA, IS, IE, JA, JS, JE, QA, &
199  q, q_mass, &
200  qdry )
201  implicit none
202  integer, intent(in) :: IA, IS, IE
203  integer, intent(in) :: JA, JS, JE
204  integer, intent(in) :: QA
205 
206  real(RP), intent(in) :: q (ia,ja,qa)
207  real(RP), intent(in) :: q_mass(qa)
208  real(RP), intent(out) :: qdry (ia,ja)
209 
210  integer :: i, j
211  !-----------------------------------------------------------------------------
212 
213  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
214  !$omp shared(JS,JE,IS,IE,qdry,q,q_mass,QA)
215  do j = js, je
216  do i = is, ie
217  call atmos_thermodyn_qdry( qa, q(i,j,:), q_mass(:), qdry(i,j) )
218  enddo
219  enddo
220 
221  return
222  end subroutine atmos_thermodyn_qdry_2d
223 
224  !-----------------------------------------------------------------------------
226  subroutine atmos_thermodyn_qdry_3d( &
227  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
228  q, q_mass, &
229  qdry )
230  implicit none
231  integer, intent(in) :: KA, KS, KE
232  integer, intent(in) :: IA, IS, IE
233  integer, intent(in) :: JA, JS, JE
234  integer, intent(in) :: QA
235 
236  real(RP), intent(in) :: q (ka,ia,ja,qa)
237  real(RP), intent(in) :: q_mass(qa)
238  real(RP), intent(out) :: qdry (ka,ia,ja)
239 
240  integer :: k, i, j
241  !-----------------------------------------------------------------------------
242 
243  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
244  !$omp shared(JS,JE,IS,IE,KS,KE,qdry,q,q_mass,QA)
245  do j = js, je
246  do i = is, ie
247  do k = ks, ke
248  call atmos_thermodyn_qdry( qa, q(k,i,j,:), q_mass(:), qdry(k,i,j) )
249  enddo
250  enddo
251  enddo
252  return
253  end subroutine atmos_thermodyn_qdry_3d
254 
255  !-----------------------------------------------------------------------------
257 !OCL SERIAL
258 !OCL NOSIMD
259  subroutine atmos_thermodyn_cv_0d( &
260  QA, &
261  q, CVq, qdry, &
262  CVtot )
263  implicit none
264  integer, intent(in) :: QA
265 
266  real(RP), intent(in) :: q(qa)
267  real(RP), intent(in) :: CVq(qa)
268  real(RP), intent(in) :: qdry
269 
270  real(RP), intent(out) :: CVtot
271 
272  integer :: iqw
273  !---------------------------------------------------------------------------
274 
275  cvtot = qdry * cvdry
276 #ifndef DRY
277  do iqw = 1, qa
278  cvtot = cvtot + q(iqw) * cvq(iqw)
279  enddo
280 #endif
281 
282  return
283  end subroutine atmos_thermodyn_cv_0d
284 
285  !-----------------------------------------------------------------------------
287  subroutine atmos_thermodyn_cv_3d( &
288  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
289  q, CVq, qdry, &
290  CVtot )
291  implicit none
292  integer, intent(in) :: KA, KS, KE
293  integer, intent(in) :: IA, IS, IE
294  integer, intent(in) :: JA, JS, JE
295  integer, intent(in) :: QA
296 
297  real(RP), intent(in) :: q (ka,ia,ja,qa)
298  real(RP), intent(in) :: CVq (qa)
299  real(RP), intent(in) :: qdry (ka,ia,ja)
300 
301  real(RP), intent(out) :: CVtot(ka,ia,ja)
302 
303  integer :: k, i, j
304  !---------------------------------------------------------------------------
305 
306  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
307  !$omp shared(JS,JE,IS,IE,KS,KE,cvtot,qdry,q,CVdry,CVq,QA)
308  do j = js, je
309  do i = is, ie
310  do k = ks, ke
311  call atmos_thermodyn_cv( qa, q(k,i,j,:), cvq(:), qdry(k,i,j), cvtot(k,i,j) )
312  enddo
313  enddo
314  enddo
315 
316  return
317  end subroutine atmos_thermodyn_cv_3d
318 
319  !-----------------------------------------------------------------------------
321 !OCL SERIAL
322 !OCL NOSIMD
323  subroutine atmos_thermodyn_cp_0d( &
324  QA, &
325  q, CPq, qdry, &
326  CPtot )
327  implicit none
328  integer, intent(in) :: QA
329 
330  real(RP), intent(in) :: q(qa)
331  real(RP), intent(in) :: CPq(qa)
332  real(RP), intent(in) :: qdry
333 
334  real(RP), intent(out) :: CPtot
335 
336  integer :: iqw
337  !---------------------------------------------------------------------------
338 
339  cptot = qdry * cpdry
340 #ifndef DRY
341  do iqw = 1, qa
342  cptot = cptot + q(iqw) * cpq(iqw)
343  enddo
344 #endif
345 
346  return
347  end subroutine atmos_thermodyn_cp_0d
348 
349  !-----------------------------------------------------------------------------
351  subroutine atmos_thermodyn_cp_3d( &
352  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
353  q, CPq, qdry, &
354  CPtot )
355  implicit none
356  integer, intent(in) :: KA, KS, KE
357  integer, intent(in) :: IA, IS, IE
358  integer, intent(in) :: JA, JS, JE
359  integer, intent(in) :: QA
360 
361  real(RP), intent(in) :: q (ka,ia,ja,qa)
362  real(RP), intent(in) :: CPq (qa)
363  real(RP), intent(in) :: qdry (ka,ia,ja)
364 
365  real(RP), intent(out) :: CPtot(ka,ia,ja)
366 
367  integer :: k, i, j
368  !---------------------------------------------------------------------------
369 
370  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
371  do j = js, je
372  do i = is, ie
373  do k = ks, ke
374  call atmos_thermodyn_cp( qa, q(k,i,j,:), cpq(:), qdry(k,i,j), cptot(k,i,j) )
375  enddo
376  enddo
377  enddo
378 
379  return
380  end subroutine atmos_thermodyn_cp_3d
381 
382  !-----------------------------------------------------------------------------
384 !OCL SERIAL
385 !OCL NOSIMD
386  subroutine atmos_thermodyn_r_0d( &
387  QA, &
388  q, Rq, qdry, &
389  Rtot )
390  implicit none
391  integer, intent(in) :: QA
392 
393  real(RP), intent(in) :: q(qa)
394  real(RP), intent(in) :: Rq(qa)
395  real(RP), intent(in) :: qdry
396 
397  real(RP), intent(out) :: Rtot
398 
399  integer :: iqw
400  !---------------------------------------------------------------------------
401 
402  rtot = qdry * rdry
403 #ifndef DRY
404  do iqw = 1, qa
405  rtot = rtot + q(iqw) * rq(iqw)
406  enddo
407 #endif
408 
409  return
410  end subroutine atmos_thermodyn_r_0d
411 
412  !-----------------------------------------------------------------------------
414  subroutine atmos_thermodyn_r_3d( &
415  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
416  q, Rq, qdry, &
417  Rtot )
418  implicit none
419  integer, intent(in) :: KA, KS, KE
420  integer, intent(in) :: IA, IS, IE
421  integer, intent(in) :: JA, JS, JE
422  integer, intent(in) :: QA
423 
424  real(RP), intent(in) :: q (ka,ia,ja,qa)
425  real(RP), intent(in) :: Rq (qa)
426  real(RP), intent(in) :: qdry(ka,ia,ja)
427 
428  real(RP), intent(out) :: Rtot(ka,ia,ja)
429 
430  integer :: k, i, j
431  !---------------------------------------------------------------------------
432 
433  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
434  do j = js, je
435  do i = is, ie
436  do k = ks, ke
437  call atmos_thermodyn_r( qa, q(k,i,j,:), rq(:), qdry(k,i,j), rtot(k,i,j) )
438  enddo
439  enddo
440  enddo
441 
442  return
443  end subroutine atmos_thermodyn_r_3d
444 
445  !-----------------------------------------------------------------------------
449 !OCL SERIAL
450 !OCL NOSIMD
451  subroutine atmos_thermodyn_specific_heat_0d( &
452  QA, &
453  q, Mq, Rq, CVq, CPq, &
454  Qdry, Rtot, CVtot, CPtot )
455  integer, intent(in) :: QA
456 
457  real(RP), intent(in) :: q (qa)
458  real(RP), intent(in) :: Mq (qa)
459  real(RP), intent(in) :: Rq (qa)
460  real(RP), intent(in) :: CVq(qa)
461  real(RP), intent(in) :: CPq(qa)
462 
463  real(RP), intent(out) :: Qdry
464  real(RP), intent(out) :: Rtot
465  real(RP), intent(out) :: CVtot
466  real(RP), intent(out) :: CPtot
467 
468  integer :: iqw
469  !---------------------------------------------------------------------------
470 
471  qdry = 1.0_rp
472  rtot = 0.0_rp
473  cvtot = 0.0_rp
474  cptot = 0.0_rp
475  do iqw = 1, qa
476  qdry = qdry - q(iqw) * mq(iqw)
477  rtot = rtot + q(iqw) * rq(iqw)
478  cvtot = cvtot + q(iqw) * cvq(iqw)
479  cptot = cptot + q(iqw) * cpq(iqw)
480  enddo
481  rtot = rtot + qdry * rdry
482  cvtot = cvtot + qdry * cvdry
483  cptot = cptot + qdry * cpdry
484 
485  return
486  end subroutine atmos_thermodyn_specific_heat_0d
487 
488  !-----------------------------------------------------------------------------
492 !OCL SERIAL
493  subroutine atmos_thermodyn_specific_heat_1d( &
494  KA, KS, KE, &
495  QA, &
496  q, Mq, Rq, CVq, CPq, &
497  Qdry, Rtot, CVtot, CPtot )
498  integer, intent(in) :: KA, KS, KE
499  integer, intent(in) :: QA
500 
501  real(RP), intent(in) :: q(ka,qa)
502  real(RP), intent(in) :: Mq (qa)
503  real(RP), intent(in) :: Rq (qa)
504  real(RP), intent(in) :: CVq(qa)
505  real(RP), intent(in) :: CPq(qa)
506 
507  real(RP), intent(out) :: Qdry (ka)
508  real(RP), intent(out) :: Rtot (ka)
509  real(RP), intent(out) :: CVtot(ka)
510  real(RP), intent(out) :: CPtot(ka)
511 
512  integer :: k
513  !---------------------------------------------------------------------------
514 
515  do k = ks, ke
516  call atmos_thermodyn_specific_heat( qa, &
517  q(k,:), mq(:), rq(:), cvq(:), cpq(:), & ! [IN]
518  qdry(k), rtot(k), cvtot(k), cptot(k) ) ! [OUT]
519  end do
520 
521  return
522  end subroutine atmos_thermodyn_specific_heat_1d
523 
524  !-----------------------------------------------------------------------------
528  subroutine atmos_thermodyn_specific_heat_2d( &
529  IA, IS, IE, JA, JS, JE, QA, &
530  q, &
531  Mq, Rq, CVq, CPq, &
532  Qdry, Rtot, CVtot, CPtot )
533  integer, intent(in) :: IA, IS, IE
534  integer, intent(in) :: JA, JS, JE
535  integer, intent(in) :: QA
536 
537  real(RP), intent(in) :: q(ia,ja,qa)
538  real(RP), intent(in) :: Mq (qa)
539  real(RP), intent(in) :: Rq (qa)
540  real(RP), intent(in) :: CVq(qa)
541  real(RP), intent(in) :: CPq(qa)
542 
543  real(RP), intent(out) :: Qdry (ia,ja)
544  real(RP), intent(out) :: Rtot (ia,ja)
545  real(RP), intent(out) :: CVtot(ia,ja)
546  real(RP), intent(out) :: CPtot(ia,ja)
547 
548  integer :: i, j
549  !---------------------------------------------------------------------------
550 
551  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
552  !$omp private(i,j)
553  do j = js, je
554  do i = is, ie
555  call atmos_thermodyn_specific_heat( qa, q(i,j,:), mq(:), rq(:), cvq(:), cpq(:), & ! [IN]
556  qdry(i,j), rtot(i,j), cvtot(i,j), cptot(i,j) ) ! [OUT]
557  end do
558  end do
559 
560  return
561  end subroutine atmos_thermodyn_specific_heat_2d
562 
563  !-----------------------------------------------------------------------------
567  subroutine atmos_thermodyn_specific_heat_3d( &
568  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
569  q, &
570  Mq, Rq, CVq, CPq, &
571  Qdry, Rtot, CVtot, CPtot )
572  integer, intent(in) :: KA, KS, KE
573  integer, intent(in) :: IA, IS, IE
574  integer, intent(in) :: JA, JS, JE
575  integer, intent(in) :: QA
576 
577  real(RP), intent(in) :: q(ka,ia,ja,qa)
578  real(RP), intent(in) :: Mq (qa)
579  real(RP), intent(in) :: Rq (qa)
580  real(RP), intent(in) :: CVq(qa)
581  real(RP), intent(in) :: CPq(qa)
582 
583  real(RP), intent(out) :: Qdry (ka,ia,ja)
584  real(RP), intent(out) :: Rtot (ka,ia,ja)
585  real(RP), intent(out) :: CVtot(ka,ia,ja)
586  real(RP), intent(out) :: CPtot(ka,ia,ja)
587 
588  integer :: k, i, j
589  !---------------------------------------------------------------------------
590 
591  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
592  !$omp private(i,j,k)
593  do j = js, je
594  do i = is, ie
595  do k = ks, ke
596  call atmos_thermodyn_specific_heat( qa, q(k,i,j,:), mq(:), rq(:), cvq(:), cpq(:), & ! [IN]
597  qdry(k,i,j), rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j) ) ! [OUT]
598  end do
599  end do
600  end do
601 
602  return
603  end subroutine atmos_thermodyn_specific_heat_3d
604 
605  !-----------------------------------------------------------------------------
607  subroutine atmos_thermodyn_rhot2pres_0d( &
608  rhot, CVtot, CPtot, Rtot, &
609  pres )
610  implicit none
611  real(RP), intent(in) :: rhot
612  real(RP), intent(in) :: CVtot
613  real(RP), intent(in) :: CPtot
614  real(RP), intent(in) :: Rtot
615 
616  real(RP), intent(out) :: pres
617  !---------------------------------------------------------------------------
618 
619  pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
620 
621  return
622  end subroutine atmos_thermodyn_rhot2pres_0d
623 
624  !-----------------------------------------------------------------------------
626  subroutine atmos_thermodyn_rhot2pres_3d( &
627  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
628  rhot, CVtot, CPtot, Rtot, &
629  pres )
630  implicit none
631  integer, intent(in) :: KA, KS, KE
632  integer, intent(in) :: IA, IS, IE
633  integer, intent(in) :: JA, JS, JE
634 
635  real(RP), intent(in) :: rhot (ka,ia,ja)
636  real(RP), intent(in) :: CVtot(ka,ia,ja)
637  real(RP), intent(in) :: CPtot(ka,ia,ja)
638  real(RP), intent(in) :: Rtot (ka,ia,ja)
639 
640  real(RP), intent(out) :: pres(ka,ia,ja)
641 
642  integer :: k, i, j
643  !---------------------------------------------------------------------------
644 
645  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
646  !$omp shared(JS,JE,IS,IE,KS,KE, &
647  !$omp rhot,CVtot,CPtot,Rtot,pres)
648  do j = js, je
649  do i = is, ie
650  do k = ks, ke
651  call atmos_thermodyn_rhot2pres( rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), pres(k,i,j) )
652  end do
653  end do
654  end do
655 
656  return
657  end subroutine atmos_thermodyn_rhot2pres_3d
658 
659  !-----------------------------------------------------------------------------
661  subroutine atmos_thermodyn_rhot2rhoe_0d( &
662  rhot, CVtot, CPtot, Rtot, &
663  rhoe )
664  implicit none
665  real(RP), intent(in) :: rhot
666  real(RP), intent(in) :: CVtot
667  real(RP), intent(in) :: CPtot
668  real(RP), intent(in) :: Rtot
669 
670  real(RP), intent(out) :: rhoe
671 
672  real(RP) :: pres
673  !---------------------------------------------------------------------------
674 
675  call atmos_thermodyn_rhot2pres( rhot, cvtot, cptot, rtot, pres )
676 
677  rhoe = pres / rtot * cvtot
678 
679  return
680  end subroutine atmos_thermodyn_rhot2rhoe_0d
681 
682  !-----------------------------------------------------------------------------
684  subroutine atmos_thermodyn_rhot2rhoe_3d( &
685  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
686  rhot, CVtot, CPtot, Rtot, &
687  rhoe )
688  implicit none
689  integer, intent(in) :: KA, KS, KE
690  integer, intent(in) :: IA, IS, IE
691  integer, intent(in) :: JA, JS, JE
692 
693  real(RP), intent(in) :: rhot (ka,ia,ja)
694  real(RP), intent(in) :: CVtot(ka,ia,ja)
695  real(RP), intent(in) :: CPtot(ka,ia,ja)
696  real(RP), intent(in) :: Rtot (ka,ia,ja)
697 
698  real(RP), intent(out) :: rhoe(ka,ia,ja)
699 
700  integer :: k, i, j
701  !---------------------------------------------------------------------------
702 
703  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
704  !$omp shared(JS,JE,IS,IE,KS,KE, &
705  !$omp rhot,CVtot,CPtot,Rtot,rhoe)
706  do j = js, je
707  do i = is, ie
708  do k = ks, ke
709  call atmos_thermodyn_rhot2rhoe( rhot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhoe(k,i,j) )
710  enddo
711  enddo
712  enddo
713 
714  return
715  end subroutine atmos_thermodyn_rhot2rhoe_3d
716 
717  !-----------------------------------------------------------------------------
719  subroutine atmos_thermodyn_rhoe2rhot_0d( &
720  rhoe, CVtot, CPtot, Rtot, &
721  rhot )
722  implicit none
723  real(RP), intent(in) :: rhoe
724  real(RP), intent(in) :: CVtot
725  real(RP), intent(in) :: CPtot
726  real(RP), intent(in) :: Rtot
727 
728  real(RP), intent(out) :: rhot
729 
730  real(RP) :: pres
731  !---------------------------------------------------------------------------
732 
733  pres = rhoe * rtot / cvtot
734 
735  rhot = rhoe / cvtot * ( pre00 / pres )**(rtot/cptot)
736 
737  return
738  end subroutine atmos_thermodyn_rhoe2rhot_0d
739 
740  !-----------------------------------------------------------------------------
742  subroutine atmos_thermodyn_rhoe2rhot_3d( &
743  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
744  rhoe, CVtot, CPtot, Rtot, &
745  rhot )
746  implicit none
747  integer, intent(in) :: KA, KS, KE
748  integer, intent(in) :: IA, IS, IE
749  integer, intent(in) :: JA, JS, JE
750 
751  real(RP), intent(in) :: rhoe (ka,ia,ja)
752  real(RP), intent(in) :: CVtot(ka,ia,ja)
753  real(RP), intent(in) :: CPtot(ka,ia,ja)
754  real(RP), intent(in) :: Rtot (ka,ia,ja)
755 
756  real(RP), intent(out) :: rhot(ka,ia,ja)
757 
758  integer :: k, i, j
759  !---------------------------------------------------------------------------
760 
761  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
762  !$omp shared(JS,JE,IS,IE,KS,KE,&
763  !$omp rhoe,CVtot,CPtot,Rtot,rhot)
764  do j = js, je
765  do i = is, ie
766  do k = ks, ke
767  call atmos_thermodyn_rhoe2rhot( rhoe(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhot(k,i,j) )
768  enddo
769  enddo
770  enddo
771 
772  return
773  end subroutine atmos_thermodyn_rhoe2rhot_3d
774 
775  !-----------------------------------------------------------------------------
777  subroutine atmos_thermodyn_rhot2temp_pres_0d( &
778  dens, rhot, Rtot, CVtot, CPtot, &
779  temp, pres )
780  implicit none
781  real(RP), intent(in) :: dens
782  real(RP), intent(in) :: rhot
783  real(RP), intent(in) :: Rtot
784  real(RP), intent(in) :: CVtot
785  real(RP), intent(in) :: CPtot
786 
787  real(RP), intent(out) :: temp
788  real(RP), intent(out) :: pres
789  !---------------------------------------------------------------------------
790 
791  pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
792  temp = pres / ( dens * rtot )
793 
794  return
795  end subroutine atmos_thermodyn_rhot2temp_pres_0d
796 
797  !-----------------------------------------------------------------------------
799 !OCL SERIAL
800  subroutine atmos_thermodyn_rhot2temp_pres_1d( &
801  KA, KS, KE, &
802  dens, rhot, Rtot, CVtot, CPtot, &
803  temp, pres )
804  integer, intent(in) :: KA, KS, KE
805 
806  real(RP), intent(in) :: dens (ka)
807  real(RP), intent(in) :: rhot (ka)
808  real(RP), intent(in) :: Rtot (ka)
809  real(RP), intent(in) :: CVtot(ka)
810  real(RP), intent(in) :: CPtot(ka)
811 
812  real(RP), intent(out) :: temp(ka)
813  real(RP), intent(out) :: pres(ka)
814 
815  integer :: k
816  !---------------------------------------------------------------------------
817 
818  do k = ks, ke
819  call atmos_thermodyn_rhot2temp_pres( dens(k), rhot(k), rtot(k), cvtot(k), cptot(k), &
820  temp(k), pres(k) )
821  enddo
822 
823  return
824  end subroutine atmos_thermodyn_rhot2temp_pres_1d
825 
826  !-----------------------------------------------------------------------------
828  subroutine atmos_thermodyn_rhot2temp_pres_3d( &
829  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
830  dens, rhot, Rtot, CVtot, CPtot, &
831  temp, pres )
832  integer, intent(in) :: KA, KS, KE
833  integer, intent(in) :: IA, IS, IE
834  integer, intent(in) :: JA, JS, JE
835 
836  real(RP), intent(in) :: dens (ka,ia,ja)
837  real(RP), intent(in) :: rhot (ka,ia,ja)
838  real(RP), intent(in) :: Rtot (ka,ia,ja)
839  real(RP), intent(in) :: CVtot(ka,ia,ja)
840  real(RP), intent(in) :: CPtot(ka,ia,ja)
841 
842  real(RP), intent(out) :: temp(ka,ia,ja)
843  real(RP), intent(out) :: pres(ka,ia,ja)
844 
845  integer :: k, i, j
846  !---------------------------------------------------------------------------
847 
848  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
849  !$omp private(i,j,k) &
850  !$omp shared(KS,KE,IS,IE,JS,JE) &
851  !$omp shared(temp,pres,dens,rhot,Rtot,CVtot,CPtot)
852  do j = js, je
853  do i = is, ie
854  do k = ks, ke
855  call atmos_thermodyn_rhot2temp_pres( dens(k,i,j), rhot(k,i,j), rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j), &
856  temp(k,i,j), pres(k,i,j) )
857  enddo
858  enddo
859  enddo
860  return
861  end subroutine atmos_thermodyn_rhot2temp_pres_3d
862 
863  !-----------------------------------------------------------------------------
865 !OCL SERIAL
866 !OCL NOSIMD
867  subroutine atmos_thermodyn_rhoe2temp_pres_0d( &
868  dens, rhoe, CVtot, Rtot, &
869  temp, pres )
870  implicit none
871  real(RP), intent(in) :: dens
872  real(RP), intent(in) :: rhoe
873  real(RP), intent(in) :: CVtot
874  real(RP), intent(in) :: Rtot
875 
876  real(RP), intent(out) :: temp
877  real(RP), intent(out) :: pres
878  !---------------------------------------------------------------------------
879 
880  temp = rhoe / ( dens * cvtot )
881  pres = dens * rtot * temp
882 
883  return
884  end subroutine atmos_thermodyn_rhoe2temp_pres_0d
885 
886  !-----------------------------------------------------------------------------
888  subroutine atmos_thermodyn_rhoe2temp_pres_2d( &
889  IA, IS, IE, JA, JS, JE, &
890  dens, rhoe, CVtot, Rtot, &
891  temp, pres )
892  implicit none
893  integer, intent(in) :: IA, IS, IE
894  integer, intent(in) :: JA, JS, JE
895 
896  real(RP), intent(in) :: dens (ia,ja)
897  real(RP), intent(in) :: rhoe (ia,ja)
898  real(RP), intent(in) :: CVtot(ia,ja)
899  real(RP), intent(in) :: Rtot (ia,ja)
900 
901  real(RP), intent(out) :: temp(ia,ja)
902  real(RP), intent(out) :: pres(ia,ja)
903 
904  integer :: i, j
905  !---------------------------------------------------------------------------
906 
907  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
908  !$omp shared(JS,JE,IS,IE,&
909  !$omp dens,rhoe,CVtot,Rtot,temp,pres)
910  do j = js, je
911  do i = is, ie
912  call atmos_thermodyn_rhoe2temp_pres( dens(i,j), rhoe(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
913  enddo
914  enddo
915 
916  return
917  end subroutine atmos_thermodyn_rhoe2temp_pres_2d
918 
919  !-----------------------------------------------------------------------------
921  subroutine atmos_thermodyn_rhoe2temp_pres_3d( &
922  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
923  dens, rhoe, CVtot, Rtot, &
924  temp, pres )
925  implicit none
926  integer, intent(in) :: KA, KS, KE
927  integer, intent(in) :: IA, IS, IE
928  integer, intent(in) :: JA, JS, JE
929 
930  real(RP), intent(in) :: dens (ka,ia,ja)
931  real(RP), intent(in) :: rhoe (ka,ia,ja)
932  real(RP), intent(in) :: CVtot(ka,ia,ja)
933  real(RP), intent(in) :: Rtot (ka,ia,ja)
934 
935  real(RP), intent(out) :: temp(ka,ia,ja)
936  real(RP), intent(out) :: pres(ka,ia,ja)
937 
938  integer :: k, i, j
939  !---------------------------------------------------------------------------
940 
941  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
942  !$omp shared(JS,JE,IS,IE,KS,KE,&
943  !$omp dens,rhoe,CVtot,Rtot,temp,pres)
944  do j = js, je
945  do i = is, ie
946  do k = ks, ke
947  call atmos_thermodyn_rhoe2temp_pres( dens(k,i,j), rhoe(k,i,j), cvtot(k,i,j), rtot(k,i,j), temp(k,i,j), pres(k,i,j) )
948  enddo
949  enddo
950  enddo
951 
952  return
953  end subroutine atmos_thermodyn_rhoe2temp_pres_3d
954 
955  !-----------------------------------------------------------------------------
957  subroutine atmos_thermodyn_ein2temp_pres_0d( &
958  Ein, dens, CVtot, Rtot, &
959  temp, pres )
960  implicit none
961  real(RP), intent(in) :: Ein ! internal energy
962  real(RP), intent(in) :: dens ! density
963  real(RP), intent(in) :: CVtot ! specific heat
964  real(RP), intent(in) :: Rtot ! gas constant
965 
966  real(RP), intent(out) :: temp ! temperature
967  real(RP), intent(out) :: pres ! pressure
968  !---------------------------------------------------------------------------
969 
970  temp = ein / cvtot
971  pres = dens * rtot * temp
972 
973  return
974  end subroutine atmos_thermodyn_ein2temp_pres_0d
975 
976  !-----------------------------------------------------------------------------
978  subroutine atmos_thermodyn_ein2temp_pres_2d( &
979  IA, IS, IE, JA, JS, JE, &
980  Ein, dens, CVtot, Rtot, &
981  temp, pres )
982  implicit none
983  integer, intent(in) :: IA, IS, IE
984  integer, intent(in) :: JA, JS, JE
985 
986  real(RP), intent(in) :: Ein (ia,ja) ! internal energy
987  real(RP), intent(in) :: dens (ia,ja) ! density
988  real(RP), intent(in) :: CVtot(ia,ja) ! specific heat
989  real(RP), intent(in) :: Rtot (ia,ja) ! gas constant
990 
991  real(RP), intent(out) :: temp(ia,ja) ! temperature
992  real(RP), intent(out) :: pres(ia,ja) ! pressure
993 
994  integer :: i, j
995  !---------------------------------------------------------------------------
996 
997  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
998  do j = js, je
999  do i = is, ie
1000  call atmos_thermodyn_ein2temp_pres( ein(i,j), dens(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
1001  enddo
1002  enddo
1003 
1004  return
1005  end subroutine atmos_thermodyn_ein2temp_pres_2d
1006 
1007  !-----------------------------------------------------------------------------
1009  subroutine atmos_thermodyn_ein2temp_pres_3d( &
1010  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1011  Ein, dens, CVtot, Rtot, &
1012  temp, pres )
1013  implicit none
1014  integer, intent(in) :: KA, KS, KE
1015  integer, intent(in) :: IA, IS, IE
1016  integer, intent(in) :: JA, JS, JE
1017 
1018  real(RP), intent(in) :: Ein (ka,ia,ja) ! internal energy
1019  real(RP), intent(in) :: dens (ka,ia,ja) ! density
1020  real(RP), intent(in) :: CVtot(ka,ia,ja) ! specific heat
1021  real(RP), intent(in) :: Rtot (ka,ia,ja) ! gas constant
1022 
1023  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
1024  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
1025 
1026  integer :: i, j, k
1027  !---------------------------------------------------------------------------
1028 
1029  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1030  do j = js, je
1031  do i = is, ie
1032  do k = ks, ke
1033  call atmos_thermodyn_ein2temp_pres( ein(k,i,j), dens(k,i,j), cvtot(k,i,j), rtot(k,i,j), temp(k,i,j), pres(k,i,j) )
1034  enddo
1035  enddo
1036  enddo
1037 
1038  return
1039  end subroutine atmos_thermodyn_ein2temp_pres_3d
1040 
1041  !-----------------------------------------------------------------------------
1043  subroutine atmos_thermodyn_pott2temp_pres_0d( &
1044  dens, pott, CVtot, CPtot, Rtot, &
1045  temp, pres )
1046  implicit none
1047 
1048  real(RP), intent(in) :: dens ! density
1049  real(RP), intent(in) :: pott ! potential temperature
1050  real(RP), intent(in) :: CVtot ! specific heat
1051  real(RP), intent(in) :: CPtot ! specific heat
1052  real(RP), intent(in) :: Rtot ! gas constant
1053 
1054  real(RP), intent(out) :: temp ! temperature
1055  real(RP), intent(out) :: pres ! pressure
1056  !---------------------------------------------------------------------------
1057 
1058  pres = pre00 * ( dens * pott * rtot / pre00 )**(cptot/cvtot)
1059  temp = pres / ( dens * rtot )
1060 
1061  return
1062  end subroutine atmos_thermodyn_pott2temp_pres_0d
1063 
1064  !-----------------------------------------------------------------------------
1066  subroutine atmos_thermodyn_pott2temp_pres_3d( &
1067  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1068  dens, pott, CVtot, CPtot, Rtot, &
1069  temp, pres )
1070  implicit none
1071  integer, intent(in) :: KA, KS, KE
1072  integer, intent(in) :: IA, IS, IE
1073  integer, intent(in) :: JA, JS, JE
1074 
1075  real(RP), intent(in) :: dens (ka,ia,ja) ! density
1076  real(RP), intent(in) :: pott (ka,ia,ja) ! potential temperature
1077  real(RP), intent(in) :: CVtot(ka,ia,ja) ! specific heat
1078  real(RP), intent(in) :: CPtot(ka,ia,ja) ! specific heat
1079  real(RP), intent(in) :: Rtot (ka,ia,ja) ! gas constant
1080 
1081  real(RP), intent(out) :: temp(ka,ia,ja) ! temperature
1082  real(RP), intent(out) :: pres(ka,ia,ja) ! pressure
1083 
1084  integer :: i, j, k
1085  !---------------------------------------------------------------------------
1086 
1087  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1088  do j = js, je
1089  do i = is, ie
1090  do k = ks, ke
1091  call atmos_thermodyn_pott2temp_pres( dens(k,i,j), pott(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), temp(k,i,j), pres(k,i,j) )
1092  enddo
1093  enddo
1094  enddo
1095 
1096  return
1097  end subroutine atmos_thermodyn_pott2temp_pres_3d
1098 
1099  !-----------------------------------------------------------------------------
1100  subroutine atmos_thermodyn_temp_pres2pott_0d( &
1101  temp, pres, CPtot, Rtot, &
1102  pott )
1103  implicit none
1104  real(RP), intent(in) :: temp ! temperature [K]
1105  real(RP), intent(in) :: pres ! pressure [Pa]
1106  real(RP), intent(in) :: CPtot ! specific heat [J/kg/K]
1107  real(RP), intent(in) :: Rtot ! gas constant [J/kg/K]
1108 
1109  real(RP), intent(out) :: pott ! potential temperature [K]
1110  !---------------------------------------------------------------------------
1111 
1112  pott = temp * ( pre00 / pres )**(rtot/cptot)
1113 
1114  return
1115  end subroutine atmos_thermodyn_temp_pres2pott_0d
1116 
1117  !-----------------------------------------------------------------------------
1118 !OCL SERIAL
1119  subroutine atmos_thermodyn_temp_pres2pott_1d( &
1120  KA, KS, KE, &
1121  temp, pres, CPtot, Rtot, &
1122  pott )
1123  implicit none
1124  integer, intent(in) :: KA, KS, KE
1125 
1126  real(RP), intent(in) :: temp (ka) ! temperature [K]
1127  real(RP), intent(in) :: pres (ka) ! pressure [Pa]
1128  real(RP), intent(in) :: CPtot(ka) ! specific heat [J/kg/K]
1129  real(RP), intent(in) :: Rtot (ka) ! gas constant [J/kg/K]
1130 
1131  real(RP), intent(out) :: pott(ka) ! potential temperature [K]
1132 
1133  integer :: k
1134  !---------------------------------------------------------------------------
1135 
1136  do k = ks, ke
1137  call atmos_thermodyn_temp_pres2pott( temp(k), pres(k), cptot(k), rtot(k), pott(k) )
1138  enddo
1139 
1140  return
1141  end subroutine atmos_thermodyn_temp_pres2pott_1d
1142 
1143  !-----------------------------------------------------------------------------
1144  subroutine atmos_thermodyn_temp_pres2pott_2d( &
1145  IA, IS, IE, JA, JS, JE, &
1146  temp, pres, CPtot, Rtot, &
1147  pott )
1148  implicit none
1149  integer, intent(in) :: IA, IS, IE
1150  integer, intent(in) :: JA, JS, JE
1151 
1152  real(RP), intent(in) :: temp (ia,ja) ! temperature [K]
1153  real(RP), intent(in) :: pres (ia,ja) ! pressure [Pa]
1154  real(RP), intent(in) :: CPtot(ia,ja) ! specific heat [J/kg/K]
1155  real(RP), intent(in) :: Rtot (ia,ja) ! gas constant [J/kg/K]
1156 
1157  real(RP), intent(out) :: pott(ia,ja) ! potential temperature [K]
1158 
1159  integer :: i, j
1160  !---------------------------------------------------------------------------
1161 
1162  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
1163  !$omp private(i,j)
1164  do j = js, je
1165  do i = is, ie
1166  call atmos_thermodyn_temp_pres2pott( temp(i,j), pres(i,j), cptot(i,j), rtot(i,j), pott(i,j) )
1167  enddo
1168  enddo
1169 
1170  return
1171  end subroutine atmos_thermodyn_temp_pres2pott_2d
1172 
1173  !-----------------------------------------------------------------------------
1174  subroutine atmos_thermodyn_temp_pres2pott_3d( &
1175  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1176  temp, pres, CPtot, Rtot, &
1177  pott )
1178  implicit none
1179  integer, intent(in) :: KA, KS, KE
1180  integer, intent(in) :: IA, IS, IE
1181  integer, intent(in) :: JA, JS, JE
1182 
1183  real(RP), intent(in) :: temp (ka,ia,ja) ! temperature [K]
1184  real(RP), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
1185  real(RP), intent(in) :: CPtot(ka,ia,ja) ! specific heat [J/kg/K]
1186  real(RP), intent(in) :: Rtot (ka,ia,ja) ! gas constant [J/kg/K]
1187 
1188  real(RP), intent(out) :: pott(ka,ia,ja) ! potential temperature [K]
1189 
1190  integer :: k, i, j
1191  !---------------------------------------------------------------------------
1192 
1193  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
1194  !$omp private(i,j,k)
1195  do j = js, je
1196  do i = is, ie
1197  do k = ks, ke
1198  call atmos_thermodyn_temp_pres2pott( temp(k,i,j), pres(k,i,j), cptot(k,i,j), rtot(k,i,j), pott(k,i,j) )
1199  enddo
1200  enddo
1201  enddo
1202 
1203  return
1204  end subroutine atmos_thermodyn_temp_pres2pott_3d
1205 
1206  !-----------------------------------------------------------------------------
1207  subroutine atmos_thermodyn_temp_pres2ein_0d( &
1208  temp, pres, CVtot, Rtot, &
1209  ein, dens )
1210  implicit none
1211  real(RP), intent(in) :: temp ! temperature [K]
1212  real(RP), intent(in) :: pres ! pressure [Pa]
1213  real(RP), intent(in) :: CVtot ! specific heat [J/kg/K]
1214  real(RP), intent(in) :: Rtot ! gas constant [J/kg/K]
1215 
1216  real(RP), intent(out) :: ein ! internal energy
1217  real(RP), intent(out) :: dens
1218  !---------------------------------------------------------------------------
1219 
1220  dens = pres / temp / rtot
1221  ein = temp * cvtot
1222 
1223  return
1224  end subroutine atmos_thermodyn_temp_pres2ein_0d
1225 
1226  !-----------------------------------------------------------------------------
1227 !OCL SERIAL
1228  subroutine atmos_thermodyn_temp_pres2ein_1d( &
1229  KA, KS, KE, &
1230  temp, pres, CVtot, Rtot, &
1231  ein, dens )
1232  implicit none
1233  integer, intent(in) :: KA, KS, KE
1234 
1235  real(RP), intent(in) :: temp (ka) ! temperature [K]
1236  real(RP), intent(in) :: pres (ka) ! pressure [Pa]
1237  real(RP), intent(in) :: CVtot(ka) ! specific heat [J/kg/K]
1238  real(RP), intent(in) :: Rtot (ka) ! gas constant [J/kg/K]
1239 
1240  real(RP), intent(out) :: ein (ka) ! internal energy
1241  real(RP), intent(out) :: dens(ka)
1242  integer :: k
1243  !---------------------------------------------------------------------------
1244 
1245  do k = ks, ke
1246  call atmos_thermodyn_temp_pres2ein_0d( temp(k), pres(k), cvtot(k), rtot(k), & ! [IN]
1247  ein(k), dens(k) ) ! [OUT]
1248  end do
1249 
1250  return
1251  end subroutine atmos_thermodyn_temp_pres2ein_1d
1252 
1253  !-----------------------------------------------------------------------------
1254  subroutine atmos_thermodyn_temp_pres2ein_2d( &
1255  IA, IS, IE, JA, JS, JE, &
1256  temp, pres, CVtot, Rtot, &
1257  ein, dens )
1258  implicit none
1259  integer, intent(in) :: IA, IS, IE
1260  integer, intent(in) :: JA, JS, JE
1261 
1262  real(RP), intent(in) :: temp (ia,ja) ! temperature [K]
1263  real(RP), intent(in) :: pres (ia,ja) ! pressure [Pa]
1264  real(RP), intent(in) :: CVtot(ia,ja) ! specific heat [J/kg/K]
1265  real(RP), intent(in) :: Rtot (ia,ja) ! gas constant [J/kg/K]
1266 
1267  real(RP), intent(out) :: ein (ia,ja) ! internal energy
1268  real(RP), intent(out) :: dens(ia,ja)
1269 
1270  integer :: i, j
1271  !---------------------------------------------------------------------------
1272 
1273  !$omp parallel do
1274  do j = js, je
1275  do i = is, ie
1276  call atmos_thermodyn_temp_pres2ein_0d( temp(i,j), pres(i,j), cvtot(i,j), rtot(i,j), & ! [IN]
1277  ein(i,j), dens(i,j) ) ! [OUT]
1278  end do
1279  end do
1280 
1281  return
1282  end subroutine atmos_thermodyn_temp_pres2ein_2d
1283 
1284  !-----------------------------------------------------------------------------
1285  subroutine atmos_thermodyn_temp_pres2ein_3d( &
1286  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1287  temp, pres, CVtot, Rtot, &
1288  ein, dens )
1289  implicit none
1290  integer, intent(in) :: KA, KS, KE
1291  integer, intent(in) :: IA, IS, IE
1292  integer, intent(in) :: JA, JS, JE
1293 
1294  real(RP), intent(in) :: temp (ka,ia,ja) ! temperature [K]
1295  real(RP), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
1296  real(RP), intent(in) :: CVtot(ka,ia,ja) ! specific heat [J/kg/K]
1297  real(RP), intent(in) :: Rtot (ka,ia,ja) ! gas constant [J/kg/K]
1298 
1299  real(RP), intent(out) :: ein (ka,ia,ja) ! internal energy
1300  real(RP), intent(out) :: dens(ka,ia,ja)
1301 
1302  integer :: k, i, j
1303  !---------------------------------------------------------------------------
1304 
1305  !$omp parallel do
1306  do j = js, je
1307  do i = is, ie
1308  do k = ks, ke
1309  call atmos_thermodyn_temp_pres2ein_0d( temp(k,i,j), pres(k,i,j), cvtot(k,i,j), rtot(k,i,j), & ! [IN]
1310  ein(k,i,j), dens(k,i,j) ) ! [OUT]
1311  end do
1312  end do
1313  end do
1314 
1315  return
1316  end subroutine atmos_thermodyn_temp_pres2ein_3d
1317 
1318  !-----------------------------------------------------------------------------
1320  subroutine atmos_thermodyn_ein_pres2enth_0d( &
1321  Ein, PRES, DENS, &
1322  ENTH )
1323  real(RP), intent(in) :: Ein
1324  real(RP), intent(in) :: PRES
1325  real(RP), intent(in) :: DENS
1326 
1327  real(RP), intent(out) :: ENTH
1328 
1329  enth = ein + pres / dens
1330 
1331  return
1332  end subroutine atmos_thermodyn_ein_pres2enth_0d
1333 
1334 !OCL SERIAL
1335  subroutine atmos_thermodyn_ein_pres2enth_1d( &
1336  KA, KS, KE, &
1337  Ein, PRES, DENS, &
1338  ENTH )
1339  integer, intent(in) :: KA, KS, KE
1340 
1341  real(RP), intent(in) :: Ein (ka)
1342  real(RP), intent(in) :: PRES(ka)
1343  real(RP), intent(in) :: DENS(ka)
1344 
1345  real(RP), intent(out) :: ENTH(ka)
1346 
1347  integer :: k
1348 
1349  do k = ks, ke
1350  call atmos_thermodyn_ein_pres2enth_0d( ein(k), pres(k), dens(k), & ! [IN]
1351  enth(k) ) ! [OUT]
1352  end do
1353 
1354  return
1355  end subroutine atmos_thermodyn_ein_pres2enth_1d
1356 
1357  subroutine atmos_thermodyn_ein_pres2enth_2d( &
1358  IA, IS, IE, JA, JS, JE, &
1359  Ein, PRES, DENS, &
1360  ENTH )
1361  integer, intent(in) :: IA, IS, IE
1362  integer, intent(in) :: JA, JS, JE
1363 
1364  real(RP), intent(in) :: Ein (ia,ja)
1365  real(RP), intent(in) :: PRES(ia,ja)
1366  real(RP), intent(in) :: DENS(ia,ja)
1367 
1368  real(RP), intent(out) :: ENTH(ia,ja)
1369 
1370  integer :: i, j
1371 
1372  !$omp parallel do
1373  do j = js, je
1374  do i = is, ie
1375  call atmos_thermodyn_ein_pres2enth_0d( ein(i,j), pres(i,j), dens(i,j), & ! [IN]
1376  enth(i,j) ) ! [OUT]
1377  end do
1378  end do
1379 
1380  return
1381  end subroutine atmos_thermodyn_ein_pres2enth_2d
1382 
1383  subroutine atmos_thermodyn_ein_pres2enth_3d( &
1384  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1385  Ein, PRES, DENS, &
1386  ENTH )
1387  integer, intent(in) :: KA, KS, KE
1388  integer, intent(in) :: IA, IS, IE
1389  integer, intent(in) :: JA, JS, JE
1390 
1391  real(RP), intent(in) :: Ein (ka,ia,ja)
1392  real(RP), intent(in) :: PRES(ka,ia,ja)
1393  real(RP), intent(in) :: DENS(ka,ia,ja)
1394 
1395  real(RP), intent(out) :: ENTH(ka,ia,ja)
1396 
1397  integer :: k, i, j
1398 
1399  !$omp parallel do
1400  do j = js, je
1401  do i = is, ie
1402  do k = ks, ke
1403  call atmos_thermodyn_ein_pres2enth_0d( ein(k,i,j), pres(k,i,j), dens(k,i,j), & ! [IN]
1404  enth(k,i,j) ) ! [OUT]
1405  end do
1406  end do
1407  end do
1408 
1409  return
1410  end subroutine atmos_thermodyn_ein_pres2enth_3d
1411 
1412 end module scale_atmos_thermodyn
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:57
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
subroutine, public atmos_thermodyn_setup
Setup.
module CONSTANT
Definition: scale_const.F90:11
module profiler
Definition: scale_prof.F90:11
module atmosphere / thermodyn
module PRECISION
module STDIO
Definition: scale_io.F90:10