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  if ( qa == 0 ) then
592  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
593  !$omp private(i,j,k)
594  do j = js, je
595  do i = is, ie
596  do k = ks, ke
597  qdry(k,i,j) = 1.0_rp
598  rtot(k,i,j) = rdry
599  cvtot(k,i,j) = cvdry
600  cptot(k,i,j) = cpdry
601  end do
602  end do
603  end do
604  else
605  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
606  !$omp private(i,j,k)
607  do j = js, je
608  do i = is, ie
609  do k = ks, ke
610  call atmos_thermodyn_specific_heat( qa, q(k,i,j,:), mq(:), rq(:), cvq(:), cpq(:), & ! [IN]
611  qdry(k,i,j), rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j) ) ! [OUT]
612  end do
613  end do
614  end do
615  end if
616 
617  return
618  end subroutine atmos_thermodyn_specific_heat_3d
619 
620  !-----------------------------------------------------------------------------
622  subroutine atmos_thermodyn_rhot2pres_0d( &
623  rhot, CVtot, CPtot, Rtot, &
624  pres )
625  implicit none
626  real(rp), intent(in) :: rhot
627  real(rp), intent(in) :: cvtot
628  real(rp), intent(in) :: cptot
629  real(rp), intent(in) :: rtot
630 
631  real(rp), intent(out) :: pres
632  !---------------------------------------------------------------------------
633 
634  pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
635 
636  return
637  end subroutine atmos_thermodyn_rhot2pres_0d
638 
639  !-----------------------------------------------------------------------------
641  subroutine atmos_thermodyn_rhot2pres_3d( &
642  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
643  rhot, CVtot, CPtot, Rtot, &
644  pres )
645  implicit none
646  integer, intent(in) :: ka, ks, ke
647  integer, intent(in) :: ia, is, ie
648  integer, intent(in) :: ja, js, je
649 
650  real(rp), intent(in) :: rhot (ka,ia,ja)
651  real(rp), intent(in) :: cvtot(ka,ia,ja)
652  real(rp), intent(in) :: cptot(ka,ia,ja)
653  real(rp), intent(in) :: rtot (ka,ia,ja)
654 
655  real(rp), intent(out) :: pres(ka,ia,ja)
656 
657  integer :: k, i, j
658  !---------------------------------------------------------------------------
659 
660  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
661  !$omp shared(JS,JE,IS,IE,KS,KE, &
662  !$omp rhot,CVtot,CPtot,Rtot,pres)
663  do j = js, je
664  do i = is, ie
665  do k = ks, ke
666  call atmos_thermodyn_rhot2pres( rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), pres(k,i,j) )
667  end do
668  end do
669  end do
670 
671  return
672  end subroutine atmos_thermodyn_rhot2pres_3d
673 
674  !-----------------------------------------------------------------------------
676  subroutine atmos_thermodyn_rhot2rhoe_0d( &
677  rhot, CVtot, CPtot, Rtot, &
678  rhoe )
679  implicit none
680  real(rp), intent(in) :: rhot
681  real(rp), intent(in) :: cvtot
682  real(rp), intent(in) :: cptot
683  real(rp), intent(in) :: rtot
684 
685  real(rp), intent(out) :: rhoe
686 
687  real(rp) :: pres
688  !---------------------------------------------------------------------------
689 
690  call atmos_thermodyn_rhot2pres( rhot, cvtot, cptot, rtot, pres )
691 
692  rhoe = pres / rtot * cvtot
693 
694  return
695  end subroutine atmos_thermodyn_rhot2rhoe_0d
696 
697  !-----------------------------------------------------------------------------
699  subroutine atmos_thermodyn_rhot2rhoe_3d( &
700  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
701  rhot, CVtot, CPtot, Rtot, &
702  rhoe )
703  implicit none
704  integer, intent(in) :: ka, ks, ke
705  integer, intent(in) :: ia, is, ie
706  integer, intent(in) :: ja, js, je
707 
708  real(rp), intent(in) :: rhot (ka,ia,ja)
709  real(rp), intent(in) :: cvtot(ka,ia,ja)
710  real(rp), intent(in) :: cptot(ka,ia,ja)
711  real(rp), intent(in) :: rtot (ka,ia,ja)
712 
713  real(rp), intent(out) :: rhoe(ka,ia,ja)
714 
715  integer :: k, i, j
716  !---------------------------------------------------------------------------
717 
718  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
719  !$omp shared(JS,JE,IS,IE,KS,KE, &
720  !$omp rhot,CVtot,CPtot,Rtot,rhoe)
721  do j = js, je
722  do i = is, ie
723  do k = ks, ke
724  call atmos_thermodyn_rhot2rhoe( rhot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhoe(k,i,j) )
725  enddo
726  enddo
727  enddo
728 
729  return
730  end subroutine atmos_thermodyn_rhot2rhoe_3d
731 
732  !-----------------------------------------------------------------------------
734  subroutine atmos_thermodyn_rhoe2rhot_0d( &
735  rhoe, CVtot, CPtot, Rtot, &
736  rhot )
737  implicit none
738  real(rp), intent(in) :: rhoe
739  real(rp), intent(in) :: cvtot
740  real(rp), intent(in) :: cptot
741  real(rp), intent(in) :: rtot
742 
743  real(rp), intent(out) :: rhot
744 
745  real(rp) :: pres
746  !---------------------------------------------------------------------------
747 
748  pres = rhoe * rtot / cvtot
749 
750  rhot = rhoe / cvtot * ( pre00 / pres )**(rtot/cptot)
751 
752  return
753  end subroutine atmos_thermodyn_rhoe2rhot_0d
754 
755  !-----------------------------------------------------------------------------
757  subroutine atmos_thermodyn_rhoe2rhot_3d( &
758  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
759  rhoe, CVtot, CPtot, Rtot, &
760  rhot )
761  implicit none
762  integer, intent(in) :: ka, ks, ke
763  integer, intent(in) :: ia, is, ie
764  integer, intent(in) :: ja, js, je
765 
766  real(rp), intent(in) :: rhoe (ka,ia,ja)
767  real(rp), intent(in) :: cvtot(ka,ia,ja)
768  real(rp), intent(in) :: cptot(ka,ia,ja)
769  real(rp), intent(in) :: rtot (ka,ia,ja)
770 
771  real(rp), intent(out) :: rhot(ka,ia,ja)
772 
773  integer :: k, i, j
774  !---------------------------------------------------------------------------
775 
776  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
777  !$omp shared(JS,JE,IS,IE,KS,KE,&
778  !$omp rhoe,CVtot,CPtot,Rtot,rhot)
779  do j = js, je
780  do i = is, ie
781  do k = ks, ke
782  call atmos_thermodyn_rhoe2rhot( rhoe(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhot(k,i,j) )
783  enddo
784  enddo
785  enddo
786 
787  return
788  end subroutine atmos_thermodyn_rhoe2rhot_3d
789 
790  !-----------------------------------------------------------------------------
792  subroutine atmos_thermodyn_rhot2temp_pres_0d( &
793  dens, rhot, Rtot, CVtot, CPtot, &
794  temp, pres )
795  implicit none
796  real(rp), intent(in) :: dens
797  real(rp), intent(in) :: rhot
798  real(rp), intent(in) :: rtot
799  real(rp), intent(in) :: cvtot
800  real(rp), intent(in) :: cptot
801 
802  real(rp), intent(out) :: temp
803  real(rp), intent(out) :: pres
804  !---------------------------------------------------------------------------
805 
806  pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
807  temp = pres / ( dens * rtot )
808 
809  return
810  end subroutine atmos_thermodyn_rhot2temp_pres_0d
811 
812  !-----------------------------------------------------------------------------
814 !OCL SERIAL
815  subroutine atmos_thermodyn_rhot2temp_pres_1d( &
816  KA, KS, KE, &
817  dens, rhot, Rtot, CVtot, CPtot, &
818  temp, pres )
819  integer, intent(in) :: ka, ks, ke
820 
821  real(rp), intent(in) :: dens (ka)
822  real(rp), intent(in) :: rhot (ka)
823  real(rp), intent(in) :: rtot (ka)
824  real(rp), intent(in) :: cvtot(ka)
825  real(rp), intent(in) :: cptot(ka)
826 
827  real(rp), intent(out) :: temp(ka)
828  real(rp), intent(out) :: pres(ka)
829 
830  integer :: k
831  !---------------------------------------------------------------------------
832 
833  do k = ks, ke
834  call atmos_thermodyn_rhot2temp_pres( dens(k), rhot(k), rtot(k), cvtot(k), cptot(k), &
835  temp(k), pres(k) )
836  enddo
837 
838  return
839  end subroutine atmos_thermodyn_rhot2temp_pres_1d
840 
841  !-----------------------------------------------------------------------------
843  subroutine atmos_thermodyn_rhot2temp_pres_3d( &
844  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
845  dens, rhot, Rtot, CVtot, CPtot, &
846  temp, pres )
847  integer, intent(in) :: ka, ks, ke
848  integer, intent(in) :: ia, is, ie
849  integer, intent(in) :: ja, js, je
850 
851  real(rp), intent(in) :: dens (ka,ia,ja)
852  real(rp), intent(in) :: rhot (ka,ia,ja)
853  real(rp), intent(in) :: rtot (ka,ia,ja)
854  real(rp), intent(in) :: cvtot(ka,ia,ja)
855  real(rp), intent(in) :: cptot(ka,ia,ja)
856 
857  real(rp), intent(out) :: temp(ka,ia,ja)
858  real(rp), intent(out) :: pres(ka,ia,ja)
859 
860  integer :: k, i, j
861  !---------------------------------------------------------------------------
862 
863  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
864  !$omp private(i,j,k) &
865  !$omp shared(KS,KE,IS,IE,JS,JE) &
866  !$omp shared(temp,pres,dens,rhot,Rtot,CVtot,CPtot)
867  do j = js, je
868  do i = is, ie
869  do k = ks, ke
870  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), &
871  temp(k,i,j), pres(k,i,j) )
872  enddo
873  enddo
874  enddo
875  return
876  end subroutine atmos_thermodyn_rhot2temp_pres_3d
877 
878  !-----------------------------------------------------------------------------
880 !OCL SERIAL
881 !OCL NOSIMD
882  subroutine atmos_thermodyn_rhoe2temp_pres_0d( &
883  dens, rhoe, CVtot, Rtot, &
884  temp, pres )
885  implicit none
886  real(rp), intent(in) :: dens
887  real(rp), intent(in) :: rhoe
888  real(rp), intent(in) :: cvtot
889  real(rp), intent(in) :: rtot
890 
891  real(rp), intent(out) :: temp
892  real(rp), intent(out) :: pres
893  !---------------------------------------------------------------------------
894 
895  temp = rhoe / ( dens * cvtot )
896  pres = dens * rtot * temp
897 
898  return
899  end subroutine atmos_thermodyn_rhoe2temp_pres_0d
900 
901  !-----------------------------------------------------------------------------
903  subroutine atmos_thermodyn_rhoe2temp_pres_2d( &
904  IA, IS, IE, JA, JS, JE, &
905  dens, rhoe, CVtot, Rtot, &
906  temp, pres )
907  implicit none
908  integer, intent(in) :: ia, is, ie
909  integer, intent(in) :: ja, js, je
910 
911  real(rp), intent(in) :: dens (ia,ja)
912  real(rp), intent(in) :: rhoe (ia,ja)
913  real(rp), intent(in) :: cvtot(ia,ja)
914  real(rp), intent(in) :: rtot (ia,ja)
915 
916  real(rp), intent(out) :: temp(ia,ja)
917  real(rp), intent(out) :: pres(ia,ja)
918 
919  integer :: i, j
920  !---------------------------------------------------------------------------
921 
922  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
923  !$omp shared(JS,JE,IS,IE,&
924  !$omp dens,rhoe,CVtot,Rtot,temp,pres)
925  do j = js, je
926  do i = is, ie
927  call atmos_thermodyn_rhoe2temp_pres( dens(i,j), rhoe(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
928  enddo
929  enddo
930 
931  return
932  end subroutine atmos_thermodyn_rhoe2temp_pres_2d
933 
934  !-----------------------------------------------------------------------------
936  subroutine atmos_thermodyn_rhoe2temp_pres_3d( &
937  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
938  dens, rhoe, CVtot, Rtot, &
939  temp, pres )
940  implicit none
941  integer, intent(in) :: ka, ks, ke
942  integer, intent(in) :: ia, is, ie
943  integer, intent(in) :: ja, js, je
944 
945  real(rp), intent(in) :: dens (ka,ia,ja)
946  real(rp), intent(in) :: rhoe (ka,ia,ja)
947  real(rp), intent(in) :: cvtot(ka,ia,ja)
948  real(rp), intent(in) :: rtot (ka,ia,ja)
949 
950  real(rp), intent(out) :: temp(ka,ia,ja)
951  real(rp), intent(out) :: pres(ka,ia,ja)
952 
953  integer :: k, i, j
954  !---------------------------------------------------------------------------
955 
956  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
957  !$omp shared(JS,JE,IS,IE,KS,KE,&
958  !$omp dens,rhoe,CVtot,Rtot,temp,pres)
959  do j = js, je
960  do i = is, ie
961  do k = ks, ke
962  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) )
963  enddo
964  enddo
965  enddo
966 
967  return
968  end subroutine atmos_thermodyn_rhoe2temp_pres_3d
969 
970  !-----------------------------------------------------------------------------
972  subroutine atmos_thermodyn_ein2temp_pres_0d( &
973  Ein, dens, CVtot, Rtot, &
974  temp, pres )
975  implicit none
976  real(rp), intent(in) :: ein ! internal energy
977  real(rp), intent(in) :: dens ! density
978  real(rp), intent(in) :: cvtot ! specific heat
979  real(rp), intent(in) :: rtot ! gas constant
980 
981  real(rp), intent(out) :: temp ! temperature
982  real(rp), intent(out) :: pres ! pressure
983  !---------------------------------------------------------------------------
984 
985  temp = ein / cvtot
986  pres = dens * rtot * temp
987 
988  return
989  end subroutine atmos_thermodyn_ein2temp_pres_0d
990 
991  !-----------------------------------------------------------------------------
993  subroutine atmos_thermodyn_ein2temp_pres_2d( &
994  IA, IS, IE, JA, JS, JE, &
995  Ein, dens, CVtot, Rtot, &
996  temp, pres )
997  implicit none
998  integer, intent(in) :: ia, is, ie
999  integer, intent(in) :: ja, js, je
1000 
1001  real(rp), intent(in) :: ein (ia,ja) ! internal energy
1002  real(rp), intent(in) :: dens (ia,ja) ! density
1003  real(rp), intent(in) :: cvtot(ia,ja) ! specific heat
1004  real(rp), intent(in) :: rtot (ia,ja) ! gas constant
1005 
1006  real(rp), intent(out) :: temp(ia,ja) ! temperature
1007  real(rp), intent(out) :: pres(ia,ja) ! pressure
1008 
1009  integer :: i, j
1010  !---------------------------------------------------------------------------
1011 
1012  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
1013  do j = js, je
1014  do i = is, ie
1015  call atmos_thermodyn_ein2temp_pres( ein(i,j), dens(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
1016  enddo
1017  enddo
1018 
1019  return
1020  end subroutine atmos_thermodyn_ein2temp_pres_2d
1021 
1022  !-----------------------------------------------------------------------------
1024  subroutine atmos_thermodyn_ein2temp_pres_3d( &
1025  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1026  Ein, dens, CVtot, Rtot, &
1027  temp, pres )
1028  implicit none
1029  integer, intent(in) :: ka, ks, ke
1030  integer, intent(in) :: ia, is, ie
1031  integer, intent(in) :: ja, js, je
1032 
1033  real(rp), intent(in) :: ein (ka,ia,ja) ! internal energy
1034  real(rp), intent(in) :: dens (ka,ia,ja) ! density
1035  real(rp), intent(in) :: cvtot(ka,ia,ja) ! specific heat
1036  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant
1037 
1038  real(rp), intent(out) :: temp(ka,ia,ja) ! temperature
1039  real(rp), intent(out) :: pres(ka,ia,ja) ! pressure
1040 
1041  integer :: i, j, k
1042  !---------------------------------------------------------------------------
1043 
1044  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1045  do j = js, je
1046  do i = is, ie
1047  do k = ks, ke
1048  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) )
1049  enddo
1050  enddo
1051  enddo
1052 
1053  return
1054  end subroutine atmos_thermodyn_ein2temp_pres_3d
1055 
1056  !-----------------------------------------------------------------------------
1058  subroutine atmos_thermodyn_pott2temp_pres_0d( &
1059  dens, pott, CVtot, CPtot, Rtot, &
1060  temp, pres )
1061  implicit none
1062 
1063  real(rp), intent(in) :: dens ! density
1064  real(rp), intent(in) :: pott ! potential temperature
1065  real(rp), intent(in) :: cvtot ! specific heat
1066  real(rp), intent(in) :: cptot ! specific heat
1067  real(rp), intent(in) :: rtot ! gas constant
1068 
1069  real(rp), intent(out) :: temp ! temperature
1070  real(rp), intent(out) :: pres ! pressure
1071  !---------------------------------------------------------------------------
1072 
1073  pres = pre00 * ( dens * pott * rtot / pre00 )**(cptot/cvtot)
1074  temp = pres / ( dens * rtot )
1075 
1076  return
1077  end subroutine atmos_thermodyn_pott2temp_pres_0d
1078 
1079  !-----------------------------------------------------------------------------
1081  subroutine atmos_thermodyn_pott2temp_pres_3d( &
1082  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1083  dens, pott, CVtot, CPtot, Rtot, &
1084  temp, pres )
1085  implicit none
1086  integer, intent(in) :: ka, ks, ke
1087  integer, intent(in) :: ia, is, ie
1088  integer, intent(in) :: ja, js, je
1089 
1090  real(rp), intent(in) :: dens (ka,ia,ja) ! density
1091  real(rp), intent(in) :: pott (ka,ia,ja) ! potential temperature
1092  real(rp), intent(in) :: cvtot(ka,ia,ja) ! specific heat
1093  real(rp), intent(in) :: cptot(ka,ia,ja) ! specific heat
1094  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant
1095 
1096  real(rp), intent(out) :: temp(ka,ia,ja) ! temperature
1097  real(rp), intent(out) :: pres(ka,ia,ja) ! pressure
1098 
1099  integer :: i, j, k
1100  !---------------------------------------------------------------------------
1101 
1102  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1103  do j = js, je
1104  do i = is, ie
1105  do k = ks, ke
1106  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) )
1107  enddo
1108  enddo
1109  enddo
1110 
1111  return
1112  end subroutine atmos_thermodyn_pott2temp_pres_3d
1113 
1114  !-----------------------------------------------------------------------------
1115  subroutine atmos_thermodyn_temp_pres2pott_0d( &
1116  temp, pres, CPtot, Rtot, &
1117  pott )
1118  implicit none
1119  real(rp), intent(in) :: temp ! temperature [K]
1120  real(rp), intent(in) :: pres ! pressure [Pa]
1121  real(rp), intent(in) :: cptot ! specific heat [J/kg/K]
1122  real(rp), intent(in) :: rtot ! gas constant [J/kg/K]
1123 
1124  real(rp), intent(out) :: pott ! potential temperature [K]
1125  !---------------------------------------------------------------------------
1126 
1127  pott = temp * ( pre00 / pres )**(rtot/cptot)
1128 
1129  return
1130  end subroutine atmos_thermodyn_temp_pres2pott_0d
1131 
1132  !-----------------------------------------------------------------------------
1133 !OCL SERIAL
1134  subroutine atmos_thermodyn_temp_pres2pott_1d( &
1135  KA, KS, KE, &
1136  temp, pres, CPtot, Rtot, &
1137  pott )
1138  implicit none
1139  integer, intent(in) :: ka, ks, ke
1140 
1141  real(rp), intent(in) :: temp (ka) ! temperature [K]
1142  real(rp), intent(in) :: pres (ka) ! pressure [Pa]
1143  real(rp), intent(in) :: cptot(ka) ! specific heat [J/kg/K]
1144  real(rp), intent(in) :: rtot (ka) ! gas constant [J/kg/K]
1145 
1146  real(rp), intent(out) :: pott(ka) ! potential temperature [K]
1147 
1148  integer :: k
1149  !---------------------------------------------------------------------------
1150 
1151  do k = ks, ke
1152  call atmos_thermodyn_temp_pres2pott( temp(k), pres(k), cptot(k), rtot(k), pott(k) )
1153  enddo
1154 
1155  return
1156  end subroutine atmos_thermodyn_temp_pres2pott_1d
1157 
1158  !-----------------------------------------------------------------------------
1159  subroutine atmos_thermodyn_temp_pres2pott_2d( &
1160  IA, IS, IE, JA, JS, JE, &
1161  temp, pres, CPtot, Rtot, &
1162  pott )
1163  implicit none
1164  integer, intent(in) :: ia, is, ie
1165  integer, intent(in) :: ja, js, je
1166 
1167  real(rp), intent(in) :: temp (ia,ja) ! temperature [K]
1168  real(rp), intent(in) :: pres (ia,ja) ! pressure [Pa]
1169  real(rp), intent(in) :: cptot(ia,ja) ! specific heat [J/kg/K]
1170  real(rp), intent(in) :: rtot (ia,ja) ! gas constant [J/kg/K]
1171 
1172  real(rp), intent(out) :: pott(ia,ja) ! potential temperature [K]
1173 
1174  integer :: i, j
1175  !---------------------------------------------------------------------------
1176 
1177  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
1178  !$omp private(i,j)
1179  do j = js, je
1180  do i = is, ie
1181  call atmos_thermodyn_temp_pres2pott( temp(i,j), pres(i,j), cptot(i,j), rtot(i,j), pott(i,j) )
1182  enddo
1183  enddo
1184 
1185  return
1186  end subroutine atmos_thermodyn_temp_pres2pott_2d
1187 
1188  !-----------------------------------------------------------------------------
1189  subroutine atmos_thermodyn_temp_pres2pott_3d( &
1190  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1191  temp, pres, CPtot, Rtot, &
1192  pott )
1193  implicit none
1194  integer, intent(in) :: ka, ks, ke
1195  integer, intent(in) :: ia, is, ie
1196  integer, intent(in) :: ja, js, je
1197 
1198  real(rp), intent(in) :: temp (ka,ia,ja) ! temperature [K]
1199  real(rp), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
1200  real(rp), intent(in) :: cptot(ka,ia,ja) ! specific heat [J/kg/K]
1201  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant [J/kg/K]
1202 
1203  real(rp), intent(out) :: pott(ka,ia,ja) ! potential temperature [K]
1204 
1205  integer :: k, i, j
1206  !---------------------------------------------------------------------------
1207 
1208  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
1209  !$omp private(i,j,k)
1210  do j = js, je
1211  do i = is, ie
1212  do k = ks, ke
1213  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) )
1214  enddo
1215  enddo
1216  enddo
1217 
1218  return
1219  end subroutine atmos_thermodyn_temp_pres2pott_3d
1220 
1221  !-----------------------------------------------------------------------------
1222  subroutine atmos_thermodyn_temp_pres2ein_0d( &
1223  temp, pres, CVtot, Rtot, &
1224  ein, dens )
1225  implicit none
1226  real(rp), intent(in) :: temp ! temperature [K]
1227  real(rp), intent(in) :: pres ! pressure [Pa]
1228  real(rp), intent(in) :: cvtot ! specific heat [J/kg/K]
1229  real(rp), intent(in) :: rtot ! gas constant [J/kg/K]
1230 
1231  real(rp), intent(out) :: ein ! internal energy
1232  real(rp), intent(out) :: dens
1233  !---------------------------------------------------------------------------
1234 
1235  dens = pres / temp / rtot
1236  ein = temp * cvtot
1237 
1238  return
1239  end subroutine atmos_thermodyn_temp_pres2ein_0d
1240 
1241  !-----------------------------------------------------------------------------
1242 !OCL SERIAL
1243  subroutine atmos_thermodyn_temp_pres2ein_1d( &
1244  KA, KS, KE, &
1245  temp, pres, CVtot, Rtot, &
1246  ein, dens )
1247  implicit none
1248  integer, intent(in) :: ka, ks, ke
1249 
1250  real(rp), intent(in) :: temp (ka) ! temperature [K]
1251  real(rp), intent(in) :: pres (ka) ! pressure [Pa]
1252  real(rp), intent(in) :: cvtot(ka) ! specific heat [J/kg/K]
1253  real(rp), intent(in) :: rtot (ka) ! gas constant [J/kg/K]
1254 
1255  real(rp), intent(out) :: ein (ka) ! internal energy
1256  real(rp), intent(out) :: dens(ka)
1257  integer :: k
1258  !---------------------------------------------------------------------------
1259 
1260  do k = ks, ke
1261  call atmos_thermodyn_temp_pres2ein_0d( temp(k), pres(k), cvtot(k), rtot(k), & ! [IN]
1262  ein(k), dens(k) ) ! [OUT]
1263  end do
1264 
1265  return
1266  end subroutine atmos_thermodyn_temp_pres2ein_1d
1267 
1268  !-----------------------------------------------------------------------------
1269  subroutine atmos_thermodyn_temp_pres2ein_2d( &
1270  IA, IS, IE, JA, JS, JE, &
1271  temp, pres, CVtot, Rtot, &
1272  ein, dens )
1273  implicit none
1274  integer, intent(in) :: ia, is, ie
1275  integer, intent(in) :: ja, js, je
1276 
1277  real(rp), intent(in) :: temp (ia,ja) ! temperature [K]
1278  real(rp), intent(in) :: pres (ia,ja) ! pressure [Pa]
1279  real(rp), intent(in) :: cvtot(ia,ja) ! specific heat [J/kg/K]
1280  real(rp), intent(in) :: rtot (ia,ja) ! gas constant [J/kg/K]
1281 
1282  real(rp), intent(out) :: ein (ia,ja) ! internal energy
1283  real(rp), intent(out) :: dens(ia,ja)
1284 
1285  integer :: i, j
1286  !---------------------------------------------------------------------------
1287 
1288  !$omp parallel do
1289  do j = js, je
1290  do i = is, ie
1291  call atmos_thermodyn_temp_pres2ein_0d( temp(i,j), pres(i,j), cvtot(i,j), rtot(i,j), & ! [IN]
1292  ein(i,j), dens(i,j) ) ! [OUT]
1293  end do
1294  end do
1295 
1296  return
1297  end subroutine atmos_thermodyn_temp_pres2ein_2d
1298 
1299  !-----------------------------------------------------------------------------
1300  subroutine atmos_thermodyn_temp_pres2ein_3d( &
1301  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1302  temp, pres, CVtot, Rtot, &
1303  ein, dens )
1304  implicit none
1305  integer, intent(in) :: ka, ks, ke
1306  integer, intent(in) :: ia, is, ie
1307  integer, intent(in) :: ja, js, je
1308 
1309  real(rp), intent(in) :: temp (ka,ia,ja) ! temperature [K]
1310  real(rp), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
1311  real(rp), intent(in) :: cvtot(ka,ia,ja) ! specific heat [J/kg/K]
1312  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant [J/kg/K]
1313 
1314  real(rp), intent(out) :: ein (ka,ia,ja) ! internal energy
1315  real(rp), intent(out) :: dens(ka,ia,ja)
1316 
1317  integer :: k, i, j
1318  !---------------------------------------------------------------------------
1319 
1320  !$omp parallel do
1321  do j = js, je
1322  do i = is, ie
1323  do k = ks, ke
1324  call atmos_thermodyn_temp_pres2ein_0d( temp(k,i,j), pres(k,i,j), cvtot(k,i,j), rtot(k,i,j), & ! [IN]
1325  ein(k,i,j), dens(k,i,j) ) ! [OUT]
1326  end do
1327  end do
1328  end do
1329 
1330  return
1331  end subroutine atmos_thermodyn_temp_pres2ein_3d
1332 
1333  !-----------------------------------------------------------------------------
1335  subroutine atmos_thermodyn_ein_pres2enth_0d( &
1336  Ein, PRES, DENS, &
1337  ENTH )
1338  real(rp), intent(in) :: ein
1339  real(rp), intent(in) :: pres
1340  real(rp), intent(in) :: dens
1341 
1342  real(rp), intent(out) :: enth
1343 
1344  enth = ein + pres / dens
1345 
1346  return
1347  end subroutine atmos_thermodyn_ein_pres2enth_0d
1348 
1349 !OCL SERIAL
1350  subroutine atmos_thermodyn_ein_pres2enth_1d( &
1351  KA, KS, KE, &
1352  Ein, PRES, DENS, &
1353  ENTH )
1354  integer, intent(in) :: ka, ks, ke
1355 
1356  real(rp), intent(in) :: ein (ka)
1357  real(rp), intent(in) :: pres(ka)
1358  real(rp), intent(in) :: dens(ka)
1359 
1360  real(rp), intent(out) :: enth(ka)
1361 
1362  integer :: k
1363 
1364  do k = ks, ke
1365  call atmos_thermodyn_ein_pres2enth_0d( ein(k), pres(k), dens(k), & ! [IN]
1366  enth(k) ) ! [OUT]
1367  end do
1368 
1369  return
1370  end subroutine atmos_thermodyn_ein_pres2enth_1d
1371 
1372  subroutine atmos_thermodyn_ein_pres2enth_2d( &
1373  IA, IS, IE, JA, JS, JE, &
1374  Ein, PRES, DENS, &
1375  ENTH )
1376  integer, intent(in) :: ia, is, ie
1377  integer, intent(in) :: ja, js, je
1378 
1379  real(rp), intent(in) :: ein (ia,ja)
1380  real(rp), intent(in) :: pres(ia,ja)
1381  real(rp), intent(in) :: dens(ia,ja)
1382 
1383  real(rp), intent(out) :: enth(ia,ja)
1384 
1385  integer :: i, j
1386 
1387  !$omp parallel do
1388  do j = js, je
1389  do i = is, ie
1390  call atmos_thermodyn_ein_pres2enth_0d( ein(i,j), pres(i,j), dens(i,j), & ! [IN]
1391  enth(i,j) ) ! [OUT]
1392  end do
1393  end do
1394 
1395  return
1396  end subroutine atmos_thermodyn_ein_pres2enth_2d
1397 
1398  subroutine atmos_thermodyn_ein_pres2enth_3d( &
1399  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1400  Ein, PRES, DENS, &
1401  ENTH )
1402  integer, intent(in) :: ka, ks, ke
1403  integer, intent(in) :: ia, is, ie
1404  integer, intent(in) :: ja, js, je
1405 
1406  real(rp), intent(in) :: ein (ka,ia,ja)
1407  real(rp), intent(in) :: pres(ka,ia,ja)
1408  real(rp), intent(in) :: dens(ka,ia,ja)
1409 
1410  real(rp), intent(out) :: enth(ka,ia,ja)
1411 
1412  integer :: k, i, j
1413 
1414  !$omp parallel do
1415  do j = js, je
1416  do i = is, ie
1417  do k = ks, ke
1418  call atmos_thermodyn_ein_pres2enth_0d( ein(k,i,j), pres(k,i,j), dens(k,i,j), & ! [IN]
1419  enth(k,i,j) ) ! [OUT]
1420  end do
1421  end do
1422  end do
1423 
1424  return
1425  end subroutine atmos_thermodyn_ein_pres2enth_3d
1426 
1427 end module scale_atmos_thermodyn
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cvdry
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:57
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_thermodyn
module atmosphere / thermodyn
Definition: scale_atmos_thermodyn.F90:11
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
scale_atmos_thermodyn::atmos_thermodyn_setup
subroutine, public atmos_thermodyn_setup
Setup.
Definition: scale_atmos_thermodyn.F90:156
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88