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  !$acc routine seq
175  implicit none
176  integer, intent(in) :: qa
177 
178  real(rp), intent(in) :: q(qa)
179  real(rp), intent(in) :: q_mass(qa)
180 
181  real(rp), intent(out) :: qdry
182 
183  integer :: iqw
184  !---------------------------------------------------------------------------
185 
186  qdry = 1.0_rp
187  do iqw = 1, qa
188  qdry = qdry - q(iqw)*q_mass(iqw)
189  enddo
190 
191  return
192  end subroutine atmos_thermodyn_qdry_0d
193 
194  !-----------------------------------------------------------------------------
196  subroutine atmos_thermodyn_qdry_2d( &
197  IA, IS, IE, JA, JS, JE, QA, &
198  q, q_mass, &
199  qdry )
200  implicit none
201  integer, intent(in) :: ia, is, ie
202  integer, intent(in) :: ja, js, je
203  integer, intent(in) :: qa
204 
205  real(rp), intent(in) :: q (ia,ja,qa)
206  real(rp), intent(in) :: q_mass(qa)
207  real(rp), intent(out) :: qdry (ia,ja)
208 
209  integer :: i, j, iqw
210  !-----------------------------------------------------------------------------
211 
212  !$omp parallel do default(none) OMP_SCHEDULE_ &
213  !$omp private(i,j,iqw) &
214  !$omp shared(IS,IE,JS,JE,QA,qdry,q,q_mass)
215  !$acc kernels copyin(q, q_mass) copyout(qdry)
216  do j = js, je
217  do i = is, ie
218  qdry(i,j) = 1.0_rp
219  end do
220  !$acc loop seq
221  do iqw = 1, qa
222  do i = is, ie
223  qdry(i,j) = qdry(i,j) - q(i,j,iqw) * q_mass(iqw)
224  enddo
225  enddo
226  enddo
227  !$acc end kernels
228 
229  return
230  end subroutine atmos_thermodyn_qdry_2d
231 
232  !-----------------------------------------------------------------------------
234  subroutine atmos_thermodyn_qdry_3d( &
235  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
236  q, q_mass, &
237  qdry )
238  implicit none
239  integer, intent(in) :: ka, ks, ke
240  integer, intent(in) :: ia, is, ie
241  integer, intent(in) :: ja, js, je
242  integer, intent(in) :: qa
243 
244  real(rp), intent(in) :: q (ka,ia,ja,qa)
245  real(rp), intent(in) :: q_mass(qa)
246  real(rp), intent(out) :: qdry (ka,ia,ja)
247 
248  integer :: k, i, j, iqw
249  !-----------------------------------------------------------------------------
250 
251  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
252  !$omp private(i,j,k,iqw) &
253  !$omp shared(KS,KE,IS,IE,JS,JE,QA,qdry,q,q_mass)
254  !$acc kernels copyin(q, q_mass) copyout(qdry)
255  do j = js, je
256  do i = is, ie
257  do k = ks, ke
258  qdry(k,i,j) = 1.0_rp
259  end do
260  !$acc loop seq
261  do iqw = 1, qa
262  do k = ks, ke
263  qdry(k,i,j) = qdry(k,i,j) - q(k,i,j,iqw) * q_mass(iqw)
264  enddo
265  enddo
266  enddo
267  enddo
268  !$acc end kernels
269 
270  return
271  end subroutine atmos_thermodyn_qdry_3d
272 
273  !-----------------------------------------------------------------------------
275 !OCL SERIAL
276 !OCL NOSIMD
277  subroutine atmos_thermodyn_cv_0d( &
278  QA, &
279  q, CVq, qdry, &
280  CVtot )
281  !$acc routine seq
282  implicit none
283  integer, intent(in) :: qa
284 
285  real(rp), intent(in) :: q(qa)
286  real(rp), intent(in) :: cvq(qa)
287  real(rp), intent(in) :: qdry
288 
289  real(rp), intent(out) :: cvtot
290 
291  integer :: iqw
292  !---------------------------------------------------------------------------
293 
294  cvtot = qdry * cvdry
295  do iqw = 1, qa
296  cvtot = cvtot + q(iqw) * cvq(iqw)
297  enddo
298 
299  return
300  end subroutine atmos_thermodyn_cv_0d
301 
302  !-----------------------------------------------------------------------------
304  subroutine atmos_thermodyn_cv_3d( &
305  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
306  q, CVq, qdry, &
307  CVtot )
308  implicit none
309  integer, intent(in) :: ka, ks, ke
310  integer, intent(in) :: ia, is, ie
311  integer, intent(in) :: ja, js, je
312  integer, intent(in) :: qa
313 
314  real(rp), intent(in) :: q (ka,ia,ja,qa)
315  real(rp), intent(in) :: cvq (qa)
316  real(rp), intent(in) :: qdry (ka,ia,ja)
317 
318  real(rp), intent(out) :: cvtot(ka,ia,ja)
319 
320  integer :: k, i, j, iqw
321  !---------------------------------------------------------------------------
322 
323  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
324  !$omp private(i,j,k) &
325  !$omp shared(KS,KE,IS,IE,JS,JE,QA,cvtot,qdry,q,CVdry,CVq)
326  !$acc kernels copyin(q, CVq, qdry) copyout(CVtot)
327  do j = js, je
328  do i = is, ie
329  do k = ks, ke
330  cvtot(k,i,j) = qdry(k,i,j) * cvdry
331  end do
332  !$acc loop seq
333  do iqw = 1, qa
334  do k = ks, ke
335  cvtot(k,i,j) = cvtot(k,i,j) + q(k,i,j,iqw) * cvq(iqw)
336  enddo
337  enddo
338  enddo
339  enddo
340  !$acc end kernels
341 
342  return
343  end subroutine atmos_thermodyn_cv_3d
344 
345  !-----------------------------------------------------------------------------
347 !OCL SERIAL
348 !OCL NOSIMD
349  subroutine atmos_thermodyn_cp_0d( &
350  QA, &
351  q, CPq, qdry, &
352  CPtot )
353  !$acc routine seq
354  implicit none
355  integer, intent(in) :: qa
356 
357  real(rp), intent(in) :: q(qa)
358  real(rp), intent(in) :: cpq(qa)
359  real(rp), intent(in) :: qdry
360 
361  real(rp), intent(out) :: cptot
362 
363  integer :: iqw
364  !---------------------------------------------------------------------------
365 
366  cptot = qdry * cpdry
367  do iqw = 1, qa
368  cptot = cptot + q(iqw) * cpq(iqw)
369  enddo
370 
371  return
372  end subroutine atmos_thermodyn_cp_0d
373 
374  !-----------------------------------------------------------------------------
376  subroutine atmos_thermodyn_cp_3d( &
377  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
378  q, CPq, qdry, &
379  CPtot )
380  implicit none
381  integer, intent(in) :: ka, ks, ke
382  integer, intent(in) :: ia, is, ie
383  integer, intent(in) :: ja, js, je
384  integer, intent(in) :: qa
385 
386  real(rp), intent(in) :: q (ka,ia,ja,qa)
387  real(rp), intent(in) :: cpq (qa)
388  real(rp), intent(in) :: qdry (ka,ia,ja)
389 
390  real(rp), intent(out) :: cptot(ka,ia,ja)
391 
392  integer :: k, i, j, iqw
393  !---------------------------------------------------------------------------
394 
395  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
396  !$omp private(i,j,k,iqw)
397  !$acc kernels copyin(q, CPq, qdry) copyout(CPtot)
398  do j = js, je
399  do i = is, ie
400  do k = ks, ke
401  cptot(k,i,j) = qdry(k,i,j) * cpdry
402  end do
403  !$acc loop seq
404  do iqw = 1, qa
405  do k = ks, ke
406  cptot(k,i,j) = cptot(k,i,j) + q(k,i,j,iqw) * cpq(iqw)
407  enddo
408  enddo
409  enddo
410  enddo
411  !$acc end kernels
412 
413  return
414  end subroutine atmos_thermodyn_cp_3d
415 
416  !-----------------------------------------------------------------------------
418 !OCL SERIAL
419 !OCL NOSIMD
420  subroutine atmos_thermodyn_r_0d( &
421  QA, &
422  q, Rq, qdry, &
423  Rtot )
424  !$acc routine seq
425  implicit none
426  integer, intent(in) :: qa
427 
428  real(rp), intent(in) :: q(qa)
429  real(rp), intent(in) :: rq(qa)
430  real(rp), intent(in) :: qdry
431 
432  real(rp), intent(out) :: rtot
433 
434  integer :: iqw
435  !---------------------------------------------------------------------------
436 
437  rtot = qdry * rdry
438  do iqw = 1, qa
439  rtot = rtot + q(iqw) * rq(iqw)
440  enddo
441 
442  return
443  end subroutine atmos_thermodyn_r_0d
444 
445  !-----------------------------------------------------------------------------
447  subroutine atmos_thermodyn_r_3d( &
448  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
449  q, Rq, qdry, &
450  Rtot )
451  implicit none
452  integer, intent(in) :: ka, ks, ke
453  integer, intent(in) :: ia, is, ie
454  integer, intent(in) :: ja, js, je
455  integer, intent(in) :: qa
456 
457  real(rp), intent(in) :: q (ka,ia,ja,qa)
458  real(rp), intent(in) :: rq (qa)
459  real(rp), intent(in) :: qdry(ka,ia,ja)
460 
461  real(rp), intent(out) :: rtot(ka,ia,ja)
462 
463  integer :: k, i, j, iqw
464  !---------------------------------------------------------------------------
465 
466  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
467  !$omp private(i,j,k,iqw)
468  !$acc kernels copyin(q, Rq, qdry) copyout(Rtot)
469  do j = js, je
470  do i = is, ie
471  do k = ks, ke
472  rtot(k,i,j) = qdry(k,i,j) * rdry
473  end do
474  !$acc loop seq
475  do iqw = 1, qa
476  do k = ks, ke
477  rtot(k,i,j) = rtot(k,i,j) + q(k,i,j,iqw) * rq(iqw)
478  enddo
479  enddo
480  enddo
481  enddo
482  !$acc end kernels
483 
484  return
485  end subroutine atmos_thermodyn_r_3d
486 
487  !-----------------------------------------------------------------------------
491 !OCL SERIAL
492 !OCL NOSIMD
493  subroutine atmos_thermodyn_specific_heat_0d( &
494  QA, &
495  q, Mq, Rq, CVq, CPq, &
496  Qdry, Rtot, CVtot, CPtot )
497  !$acc routine seq
498  integer, intent(in) :: qa
499 
500  real(rp), intent(in) :: q (qa)
501  real(rp), intent(in) :: mq (qa)
502  real(rp), intent(in) :: rq (qa)
503  real(rp), intent(in) :: cvq(qa)
504  real(rp), intent(in) :: cpq(qa)
505 
506  real(rp), intent(out) :: qdry
507  real(rp), intent(out) :: rtot
508  real(rp), intent(out) :: cvtot
509  real(rp), intent(out) :: cptot
510 
511  integer :: iqw
512  !---------------------------------------------------------------------------
513 
514  qdry = 1.0_rp
515  rtot = 0.0_rp
516  cvtot = 0.0_rp
517  cptot = 0.0_rp
518  do iqw = 1, qa
519  qdry = qdry - q(iqw) * mq(iqw)
520  rtot = rtot + q(iqw) * rq(iqw)
521  cvtot = cvtot + q(iqw) * cvq(iqw)
522  cptot = cptot + q(iqw) * cpq(iqw)
523  enddo
524  rtot = rtot + qdry * rdry
525  cvtot = cvtot + qdry * cvdry
526  cptot = cptot + qdry * cpdry
527 
528  return
529  end subroutine atmos_thermodyn_specific_heat_0d
530 
531  !-----------------------------------------------------------------------------
535 !OCL SERIAL
536  subroutine atmos_thermodyn_specific_heat_1d( &
537  KA, KS, KE, &
538  QA, &
539  q, Mq, Rq, CVq, CPq, &
540  Qdry, Rtot, CVtot, CPtot )
541  !$acc routine vector
542  integer, intent(in) :: ka, ks, ke
543  integer, intent(in) :: qa
544 
545  real(rp), intent(in) :: q(ka,qa)
546  real(rp), intent(in) :: mq (qa)
547  real(rp), intent(in) :: rq (qa)
548  real(rp), intent(in) :: cvq(qa)
549  real(rp), intent(in) :: cpq(qa)
550 
551  real(rp), intent(out) :: qdry (ka)
552  real(rp), intent(out) :: rtot (ka)
553  real(rp), intent(out) :: cvtot(ka)
554  real(rp), intent(out) :: cptot(ka)
555 
556  integer :: k, iqw
557  !---------------------------------------------------------------------------
558 
559  do k = ks, ke
560  qdry(k) = 1.0_rp
561  rtot(k) = 0.0_rp
562  cvtot(k) = 0.0_rp
563  cptot(k) = 0.0_rp
564  end do
565  do iqw = 1, qa
566  do k = ks, ke
567  qdry(k) = qdry(k) - q(k,iqw) * mq(iqw)
568  rtot(k) = rtot(k) + q(k,iqw) * rq(iqw)
569  cvtot(k) = cvtot(k) + q(k,iqw) * cvq(iqw)
570  cptot(k) = cptot(k) + q(k,iqw) * cpq(iqw)
571  enddo
572  enddo
573  do k = ks, ke
574  rtot(k) = rtot(k) + qdry(k) * rdry
575  cvtot(k) = cvtot(k) + qdry(k) * cvdry
576  cptot(k) = cptot(k) + qdry(k) * cpdry
577  end do
578 
579  return
580  end subroutine atmos_thermodyn_specific_heat_1d
581 
582  !-----------------------------------------------------------------------------
586  subroutine atmos_thermodyn_specific_heat_2d( &
587  IA, IS, IE, JA, JS, JE, QA, &
588  q, &
589  Mq, Rq, CVq, CPq, &
590  Qdry, Rtot, CVtot, CPtot )
591  integer, intent(in) :: ia, is, ie
592  integer, intent(in) :: ja, js, je
593  integer, intent(in) :: qa
594 
595  real(rp), intent(in) :: q(ia,ja,qa)
596  real(rp), intent(in) :: mq (qa)
597  real(rp), intent(in) :: rq (qa)
598  real(rp), intent(in) :: cvq(qa)
599  real(rp), intent(in) :: cpq(qa)
600 
601  real(rp), intent(out) :: qdry (ia,ja)
602  real(rp), intent(out) :: rtot (ia,ja)
603  real(rp), intent(out) :: cvtot(ia,ja)
604  real(rp), intent(out) :: cptot(ia,ja)
605 
606  integer :: i, j, iqw
607  !---------------------------------------------------------------------------
608 
609  !$omp parallel do OMP_SCHEDULE_ &
610  !$omp private(i,j,iqw)
611  !$acc kernels copyin(q, Mq, Rq, CVq, CPq) copyout(Qdry, Rtot, CVtot, CPtot)
612  do j = js, je
613  do i = is, ie
614  qdry(i,j) = 1.0_rp
615  rtot(i,j) = 0.0_rp
616  cvtot(i,j) = 0.0_rp
617  cptot(i,j) = 0.0_rp
618  end do
619  !$acc loop seq
620  do iqw = 1, qa
621  do i = is, ie
622  qdry(i,j) = qdry(i,j) - q(i,j,iqw) * mq(iqw)
623  rtot(i,j) = rtot(i,j) + q(i,j,iqw) * rq(iqw)
624  cvtot(i,j) = cvtot(i,j) + q(i,j,iqw) * cvq(iqw)
625  cptot(i,j) = cptot(i,j) + q(i,j,iqw) * cpq(iqw)
626  enddo
627  enddo
628  do i = is, ie
629  rtot(i,j) = rtot(i,j) + qdry(i,j) * rdry
630  cvtot(i,j) = cvtot(i,j) + qdry(i,j) * cvdry
631  cptot(i,j) = cptot(i,j) + qdry(i,j) * cpdry
632  end do
633  end do
634  !$acc end kernels
635 
636  return
637  end subroutine atmos_thermodyn_specific_heat_2d
638 
639  !-----------------------------------------------------------------------------
643  subroutine atmos_thermodyn_specific_heat_3d( &
644  KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
645  q, &
646  Mq, Rq, CVq, CPq, &
647  Qdry, Rtot, CVtot, CPtot )
648  integer, intent(in) :: ka, ks, ke
649  integer, intent(in) :: ia, is, ie
650  integer, intent(in) :: ja, js, je
651  integer, intent(in) :: qa
652 
653  real(rp), intent(in) :: q(ka,ia,ja,qa)
654  real(rp), intent(in) :: mq (qa)
655  real(rp), intent(in) :: rq (qa)
656  real(rp), intent(in) :: cvq(qa)
657  real(rp), intent(in) :: cpq(qa)
658 
659  real(rp), intent(out) :: qdry (ka,ia,ja)
660  real(rp), intent(out) :: rtot (ka,ia,ja)
661  real(rp), intent(out) :: cvtot(ka,ia,ja)
662  real(rp), intent(out) :: cptot(ka,ia,ja)
663 
664  integer :: k, i, j, iqw
665  !---------------------------------------------------------------------------
666 
667  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
668  !$omp private(i,j,k,iqw)
669  !$acc kernels copyin(q, Mq, Rq, CVq, CPq) copyout(Qdry, Rtot, CVtot, CPtot)
670  do j = js, je
671  do i = is, ie
672  do k = ks, ke
673  qdry(k,i,j) = 1.0_rp
674  rtot(k,i,j) = 0.0_rp
675  cvtot(k,i,j) = 0.0_rp
676  cptot(k,i,j) = 0.0_rp
677  end do
678  !$acc loop seq
679  do iqw = 1, qa
680  do k = ks, ke
681  qdry(k,i,j) = qdry(k,i,j) - q(k,i,j,iqw) * mq(iqw)
682  rtot(k,i,j) = rtot(k,i,j) + q(k,i,j,iqw) * rq(iqw)
683  cvtot(k,i,j) = cvtot(k,i,j) + q(k,i,j,iqw) * cvq(iqw)
684  cptot(k,i,j) = cptot(k,i,j) + q(k,i,j,iqw) * cpq(iqw)
685  enddo
686  enddo
687  do k = ks, ke
688  rtot(k,i,j) = rtot(k,i,j) + qdry(k,i,j) * rdry
689  cvtot(k,i,j) = cvtot(k,i,j) + qdry(k,i,j) * cvdry
690  cptot(k,i,j) = cptot(k,i,j) + qdry(k,i,j) * cpdry
691  end do
692  end do
693  end do
694  !$acc end kernels
695 
696  return
697  end subroutine atmos_thermodyn_specific_heat_3d
698 
699  !-----------------------------------------------------------------------------
701  subroutine atmos_thermodyn_rhot2pres_0d( &
702  rhot, CVtot, CPtot, Rtot, &
703  pres )
704  !$acc routine seq
705  implicit none
706  real(rp), intent(in) :: rhot
707  real(rp), intent(in) :: cvtot
708  real(rp), intent(in) :: cptot
709  real(rp), intent(in) :: rtot
710 
711  real(rp), intent(out) :: pres
712  !---------------------------------------------------------------------------
713 
714  pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
715 
716  return
717  end subroutine atmos_thermodyn_rhot2pres_0d
718 
719  !-----------------------------------------------------------------------------
721  subroutine atmos_thermodyn_rhot2pres_3d( &
722  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
723  rhot, CVtot, CPtot, Rtot, &
724  pres )
725  implicit none
726  integer, intent(in) :: ka, ks, ke
727  integer, intent(in) :: ia, is, ie
728  integer, intent(in) :: ja, js, je
729 
730  real(rp), intent(in) :: rhot (ka,ia,ja)
731  real(rp), intent(in) :: cvtot(ka,ia,ja)
732  real(rp), intent(in) :: cptot(ka,ia,ja)
733  real(rp), intent(in) :: rtot (ka,ia,ja)
734 
735  real(rp), intent(out) :: pres(ka,ia,ja)
736 
737  integer :: k, i, j
738  !---------------------------------------------------------------------------
739 
740  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
741  !$omp shared(JS,JE,IS,IE,KS,KE, &
742  !$omp rhot,CVtot,CPtot,Rtot,pres)
743  !$acc kernels copyin(rhot, CVtot, CPtot, Rtot) copyout(pres)
744  do j = js, je
745  do i = is, ie
746  do k = ks, ke
747  call atmos_thermodyn_rhot2pres( rhot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), pres(k,i,j) )
748  end do
749  end do
750  end do
751  !$acc end kernels
752 
753  return
754  end subroutine atmos_thermodyn_rhot2pres_3d
755 
756  !-----------------------------------------------------------------------------
758  subroutine atmos_thermodyn_rhot2rhoe_0d( &
759  rhot, CVtot, CPtot, Rtot, &
760  rhoe )
761  !$acc routine seq
762  implicit none
763  real(rp), intent(in) :: rhot
764  real(rp), intent(in) :: cvtot
765  real(rp), intent(in) :: cptot
766  real(rp), intent(in) :: rtot
767 
768  real(rp), intent(out) :: rhoe
769 
770  real(rp) :: pres
771  !---------------------------------------------------------------------------
772 
773  call atmos_thermodyn_rhot2pres( rhot, cvtot, cptot, rtot, pres )
774 
775  rhoe = pres / rtot * cvtot
776 
777  return
778  end subroutine atmos_thermodyn_rhot2rhoe_0d
779 
780  !-----------------------------------------------------------------------------
782  subroutine atmos_thermodyn_rhot2rhoe_3d( &
783  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
784  rhot, CVtot, CPtot, Rtot, &
785  rhoe )
786  implicit none
787  integer, intent(in) :: ka, ks, ke
788  integer, intent(in) :: ia, is, ie
789  integer, intent(in) :: ja, js, je
790 
791  real(rp), intent(in) :: rhot (ka,ia,ja)
792  real(rp), intent(in) :: cvtot(ka,ia,ja)
793  real(rp), intent(in) :: cptot(ka,ia,ja)
794  real(rp), intent(in) :: rtot (ka,ia,ja)
795 
796  real(rp), intent(out) :: rhoe(ka,ia,ja)
797 
798  integer :: k, i, j
799  !---------------------------------------------------------------------------
800 
801  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
802  !$omp shared(JS,JE,IS,IE,KS,KE, &
803  !$omp rhot,CVtot,CPtot,Rtot,rhoe)
804  !$acc kernels copyin(rhot, CVtot, CPtot, Rtot) copyout(rhoe)
805  do j = js, je
806  do i = is, ie
807  do k = ks, ke
808  call atmos_thermodyn_rhot2rhoe( rhot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhoe(k,i,j) )
809  enddo
810  enddo
811  enddo
812  !$acc end kernels
813 
814  return
815  end subroutine atmos_thermodyn_rhot2rhoe_3d
816 
817  !-----------------------------------------------------------------------------
819  subroutine atmos_thermodyn_rhoe2rhot_0d( &
820  rhoe, CVtot, CPtot, Rtot, &
821  rhot )
822  !$acc routine seq
823  implicit none
824  real(rp), intent(in) :: rhoe
825  real(rp), intent(in) :: cvtot
826  real(rp), intent(in) :: cptot
827  real(rp), intent(in) :: rtot
828 
829  real(rp), intent(out) :: rhot
830 
831  real(rp) :: pres
832  !---------------------------------------------------------------------------
833 
834  pres = rhoe * rtot / cvtot
835 
836  rhot = rhoe / cvtot * ( pre00 / pres )**(rtot/cptot)
837 
838  return
839  end subroutine atmos_thermodyn_rhoe2rhot_0d
840 
841  !-----------------------------------------------------------------------------
843  subroutine atmos_thermodyn_rhoe2rhot_3d( &
844  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
845  rhoe, CVtot, CPtot, Rtot, &
846  rhot )
847  implicit none
848  integer, intent(in) :: ka, ks, ke
849  integer, intent(in) :: ia, is, ie
850  integer, intent(in) :: ja, js, je
851 
852  real(rp), intent(in) :: rhoe (ka,ia,ja)
853  real(rp), intent(in) :: cvtot(ka,ia,ja)
854  real(rp), intent(in) :: cptot(ka,ia,ja)
855  real(rp), intent(in) :: rtot (ka,ia,ja)
856 
857  real(rp), intent(out) :: rhot(ka,ia,ja)
858 
859  integer :: k, i, j
860  !---------------------------------------------------------------------------
861 
862  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
863  !$omp shared(JS,JE,IS,IE,KS,KE,&
864  !$omp rhoe,CVtot,CPtot,Rtot,rhot)
865  !$acc kernels copyin(rhoe, CVtot, CPtot, Rtot) copyout(rhot)
866  do j = js, je
867  do i = is, ie
868  do k = ks, ke
869  call atmos_thermodyn_rhoe2rhot( rhoe(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhot(k,i,j) )
870  enddo
871  enddo
872  enddo
873  !$acc end kernels
874 
875  return
876  end subroutine atmos_thermodyn_rhoe2rhot_3d
877 
878  !-----------------------------------------------------------------------------
880  subroutine atmos_thermodyn_rhot2temp_pres_0d( &
881  dens, rhot, Rtot, CVtot, CPtot, &
882  temp, pres )
883  !$acc routine seq
884  implicit none
885  real(rp), intent(in) :: dens
886  real(rp), intent(in) :: rhot
887  real(rp), intent(in) :: rtot
888  real(rp), intent(in) :: cvtot
889  real(rp), intent(in) :: cptot
890 
891  real(rp), intent(out) :: temp
892  real(rp), intent(out) :: pres
893  !---------------------------------------------------------------------------
894 
895  pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
896  temp = pres / ( dens * rtot )
897 
898  return
899  end subroutine atmos_thermodyn_rhot2temp_pres_0d
900 
901  !-----------------------------------------------------------------------------
903 !OCL SERIAL
904  subroutine atmos_thermodyn_rhot2temp_pres_1d( &
905  KA, KS, KE, &
906  dens, rhot, Rtot, CVtot, CPtot, &
907  temp, pres )
908  !$acc routine vector
909  integer, intent(in) :: ka, ks, ke
910 
911  real(rp), intent(in) :: dens (ka)
912  real(rp), intent(in) :: rhot (ka)
913  real(rp), intent(in) :: rtot (ka)
914  real(rp), intent(in) :: cvtot(ka)
915  real(rp), intent(in) :: cptot(ka)
916 
917  real(rp), intent(out) :: temp(ka)
918  real(rp), intent(out) :: pres(ka)
919 
920  integer :: k
921  !---------------------------------------------------------------------------
922 
923  do k = ks, ke
924  call atmos_thermodyn_rhot2temp_pres( dens(k), rhot(k), rtot(k), cvtot(k), cptot(k), &
925  temp(k), pres(k) )
926  enddo
927 
928  return
929  end subroutine atmos_thermodyn_rhot2temp_pres_1d
930 
931  !-----------------------------------------------------------------------------
933  subroutine atmos_thermodyn_rhot2temp_pres_3d( &
934  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
935  dens, rhot, Rtot, CVtot, CPtot, &
936  temp, pres )
937  integer, intent(in) :: ka, ks, ke
938  integer, intent(in) :: ia, is, ie
939  integer, intent(in) :: ja, js, je
940 
941  real(rp), intent(in) :: dens (ka,ia,ja)
942  real(rp), intent(in) :: rhot (ka,ia,ja)
943  real(rp), intent(in) :: rtot (ka,ia,ja)
944  real(rp), intent(in) :: cvtot(ka,ia,ja)
945  real(rp), intent(in) :: cptot(ka,ia,ja)
946 
947  real(rp), intent(out) :: temp(ka,ia,ja)
948  real(rp), intent(out) :: pres(ka,ia,ja)
949 
950  integer :: k, i, j
951  !---------------------------------------------------------------------------
952 
953  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
954  !$omp private(i,j,k) &
955  !$omp shared(KS,KE,IS,IE,JS,JE) &
956  !$omp shared(temp,pres,dens,rhot,Rtot,CVtot,CPtot)
957  !$acc kernels copyin(dens, rhot, Rtot, CVtot, CPtot) copyout(temp, pres)
958  do j = js, je
959  do i = is, ie
960  do k = ks, ke
961  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), &
962  temp(k,i,j), pres(k,i,j) )
963  enddo
964  enddo
965  enddo
966  !$acc end kernels
967 
968  return
969  end subroutine atmos_thermodyn_rhot2temp_pres_3d
970 
971  !-----------------------------------------------------------------------------
973 !OCL SERIAL
974  subroutine atmos_thermodyn_rhoe2temp_pres_0d( &
975  dens, rhoe, CVtot, Rtot, &
976  temp, pres )
977  !$acc routine seq
978  implicit none
979  real(rp), intent(in) :: dens
980  real(rp), intent(in) :: rhoe
981  real(rp), intent(in) :: cvtot
982  real(rp), intent(in) :: rtot
983 
984  real(rp), intent(out) :: temp
985  real(rp), intent(out) :: pres
986  !---------------------------------------------------------------------------
987 
988  temp = rhoe / ( dens * cvtot )
989  pres = dens * rtot * temp
990 
991  return
992  end subroutine atmos_thermodyn_rhoe2temp_pres_0d
993 
994  !-----------------------------------------------------------------------------
996  subroutine atmos_thermodyn_rhoe2temp_pres_2d( &
997  IA, IS, IE, JA, JS, JE, &
998  dens, rhoe, CVtot, Rtot, &
999  temp, pres )
1000  implicit none
1001  integer, intent(in) :: ia, is, ie
1002  integer, intent(in) :: ja, js, je
1003 
1004  real(rp), intent(in) :: dens (ia,ja)
1005  real(rp), intent(in) :: rhoe (ia,ja)
1006  real(rp), intent(in) :: cvtot(ia,ja)
1007  real(rp), intent(in) :: rtot (ia,ja)
1008 
1009  real(rp), intent(out) :: temp(ia,ja)
1010  real(rp), intent(out) :: pres(ia,ja)
1011 
1012  integer :: i, j
1013  !---------------------------------------------------------------------------
1014 
1015  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ &
1016  !$omp shared(JS,JE,IS,IE,&
1017  !$omp dens,rhoe,CVtot,Rtot,temp,pres)
1018  !$acc kernels copyin(dens, rhoe, CVtot, Rtot) copyout(temp, pres)
1019  do j = js, je
1020  do i = is, ie
1021 ! call ATMOS_THERMODYN_rhoe2temp_pres_0D( dens(i,j), rhoe(i,j), CVtot(i,j), Rtot(i,j), temp(i,j), pres(i,j) )
1022  temp(i,j) = rhoe(i,j) / ( dens(i,j) * cvtot(i,j) )
1023  pres(i,j) = dens(i,j) * rtot(i,j) * temp(i,j)
1024  enddo
1025  enddo
1026  !$acc end kernels
1027 
1028  return
1029  end subroutine atmos_thermodyn_rhoe2temp_pres_2d
1030 
1031  !-----------------------------------------------------------------------------
1033  subroutine atmos_thermodyn_rhoe2temp_pres_3d( &
1034  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1035  dens, rhoe, CVtot, Rtot, &
1036  temp, pres )
1037  implicit none
1038  integer, intent(in) :: ka, ks, ke
1039  integer, intent(in) :: ia, is, ie
1040  integer, intent(in) :: ja, js, je
1041 
1042  real(rp), intent(in) :: dens (ka,ia,ja)
1043  real(rp), intent(in) :: rhoe (ka,ia,ja)
1044  real(rp), intent(in) :: cvtot(ka,ia,ja)
1045  real(rp), intent(in) :: rtot (ka,ia,ja)
1046 
1047  real(rp), intent(out) :: temp(ka,ia,ja)
1048  real(rp), intent(out) :: pres(ka,ia,ja)
1049 
1050  integer :: k, i, j
1051  !---------------------------------------------------------------------------
1052 
1053  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1054  !$omp shared(JS,JE,IS,IE,KS,KE,&
1055  !$omp dens,rhoe,CVtot,Rtot,temp,pres)
1056  !$acc kernels copyin(dens, rhoe, CVtot, Rtot) copyout(temp, pres)
1057  do j = js, je
1058  do i = is, ie
1059  do k = ks, ke
1060 ! 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) )
1061  temp(k,i,j) = rhoe(k,i,j) / ( dens(k,i,j) * cvtot(k,i,j) )
1062  pres(k,i,j) = dens(k,i,j) * rtot(k,i,j) * temp(k,i,j)
1063  enddo
1064  enddo
1065  enddo
1066  !$acc end kernels
1067 
1068  return
1069  end subroutine atmos_thermodyn_rhoe2temp_pres_3d
1070 
1071  !-----------------------------------------------------------------------------
1073  subroutine atmos_thermodyn_ein2temp_pres_0d( &
1074  Ein, dens, CVtot, Rtot, &
1075  temp, pres )
1076  !$acc routine seq
1077  implicit none
1078  real(rp), intent(in) :: ein ! internal energy
1079  real(rp), intent(in) :: dens ! density
1080  real(rp), intent(in) :: cvtot ! specific heat
1081  real(rp), intent(in) :: rtot ! gas constant
1082 
1083  real(rp), intent(out) :: temp ! temperature
1084  real(rp), intent(out) :: pres ! pressure
1085  !---------------------------------------------------------------------------
1086 
1087  temp = ein / cvtot
1088  pres = dens * rtot * temp
1089 
1090  return
1091  end subroutine atmos_thermodyn_ein2temp_pres_0d
1092 
1093  !-----------------------------------------------------------------------------
1095  subroutine atmos_thermodyn_ein2temp_pres_2d( &
1096  IA, IS, IE, JA, JS, JE, &
1097  Ein, dens, CVtot, Rtot, &
1098  temp, pres )
1099  implicit none
1100  integer, intent(in) :: ia, is, ie
1101  integer, intent(in) :: ja, js, je
1102 
1103  real(rp), intent(in) :: ein (ia,ja) ! internal energy
1104  real(rp), intent(in) :: dens (ia,ja) ! density
1105  real(rp), intent(in) :: cvtot(ia,ja) ! specific heat
1106  real(rp), intent(in) :: rtot (ia,ja) ! gas constant
1107 
1108  real(rp), intent(out) :: temp(ia,ja) ! temperature
1109  real(rp), intent(out) :: pres(ia,ja) ! pressure
1110 
1111  integer :: i, j
1112  !---------------------------------------------------------------------------
1113 
1114  !$omp parallel do private(i,j) OMP_SCHEDULE_
1115  !$acc kernels copyin(Ein, dens, CVtot, Rtot) copyout(temp, pres)
1116  do j = js, je
1117  do i = is, ie
1118  call atmos_thermodyn_ein2temp_pres( ein(i,j), dens(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
1119  enddo
1120  enddo
1121  !$acc end kernels
1122 
1123  return
1124  end subroutine atmos_thermodyn_ein2temp_pres_2d
1125 
1126  !-----------------------------------------------------------------------------
1128  subroutine atmos_thermodyn_ein2temp_pres_3d( &
1129  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1130  Ein, dens, CVtot, Rtot, &
1131  temp, pres )
1132  implicit none
1133  integer, intent(in) :: ka, ks, ke
1134  integer, intent(in) :: ia, is, ie
1135  integer, intent(in) :: ja, js, je
1136 
1137  real(rp), intent(in) :: ein (ka,ia,ja) ! internal energy
1138  real(rp), intent(in) :: dens (ka,ia,ja) ! density
1139  real(rp), intent(in) :: cvtot(ka,ia,ja) ! specific heat
1140  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant
1141 
1142  real(rp), intent(out) :: temp(ka,ia,ja) ! temperature
1143  real(rp), intent(out) :: pres(ka,ia,ja) ! pressure
1144 
1145  integer :: i, j, k
1146  !---------------------------------------------------------------------------
1147 
1148  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1149  !$acc kernels copyin(Ein, dens, CVtot, Rtot) copyout(temp, pres)
1150  do j = js, je
1151  do i = is, ie
1152  do k = ks, ke
1153  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) )
1154  enddo
1155  enddo
1156  enddo
1157  !$acc end kernels
1158 
1159  return
1160  end subroutine atmos_thermodyn_ein2temp_pres_3d
1161 
1162  !-----------------------------------------------------------------------------
1164  subroutine atmos_thermodyn_pott2temp_pres_0d( &
1165  dens, pott, CVtot, CPtot, Rtot, &
1166  temp, pres )
1167  !$acc routine seq
1168  implicit none
1169 
1170  real(rp), intent(in) :: dens ! density
1171  real(rp), intent(in) :: pott ! potential temperature
1172  real(rp), intent(in) :: cvtot ! specific heat
1173  real(rp), intent(in) :: cptot ! specific heat
1174  real(rp), intent(in) :: rtot ! gas constant
1175 
1176  real(rp), intent(out) :: temp ! temperature
1177  real(rp), intent(out) :: pres ! pressure
1178  !---------------------------------------------------------------------------
1179 
1180  pres = pre00 * ( dens * pott * rtot / pre00 )**(cptot/cvtot)
1181  temp = pres / ( dens * rtot )
1182 
1183  return
1184  end subroutine atmos_thermodyn_pott2temp_pres_0d
1185 
1186  !-----------------------------------------------------------------------------
1188  subroutine atmos_thermodyn_pott2temp_pres_3d( &
1189  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1190  dens, pott, CVtot, CPtot, Rtot, &
1191  temp, pres )
1192  implicit none
1193  integer, intent(in) :: ka, ks, ke
1194  integer, intent(in) :: ia, is, ie
1195  integer, intent(in) :: ja, js, je
1196 
1197  real(rp), intent(in) :: dens (ka,ia,ja) ! density
1198  real(rp), intent(in) :: pott (ka,ia,ja) ! potential temperature
1199  real(rp), intent(in) :: cvtot(ka,ia,ja) ! specific heat
1200  real(rp), intent(in) :: cptot(ka,ia,ja) ! specific heat
1201  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant
1202 
1203  real(rp), intent(out) :: temp(ka,ia,ja) ! temperature
1204  real(rp), intent(out) :: pres(ka,ia,ja) ! pressure
1205 
1206  integer :: i, j, k
1207  !---------------------------------------------------------------------------
1208 
1209  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1210  !$acc kernels copyin(dens, pott, CVtot, CPtot, Rtot) copyout(temp, pres)
1211  do j = js, je
1212  do i = is, ie
1213  do k = ks, ke
1214  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) )
1215  enddo
1216  enddo
1217  enddo
1218  !$acc end kernels
1219 
1220  return
1221  end subroutine atmos_thermodyn_pott2temp_pres_3d
1222 
1223  !-----------------------------------------------------------------------------
1224  subroutine atmos_thermodyn_temp_pres2pott_0d( &
1225  temp, pres, CPtot, Rtot, &
1226  pott )
1227  !$acc routine seq
1228  implicit none
1229  real(rp), intent(in) :: temp ! temperature [K]
1230  real(rp), intent(in) :: pres ! pressure [Pa]
1231  real(rp), intent(in) :: cptot ! specific heat [J/kg/K]
1232  real(rp), intent(in) :: rtot ! gas constant [J/kg/K]
1233 
1234  real(rp), intent(out) :: pott ! potential temperature [K]
1235  !---------------------------------------------------------------------------
1236 
1237  pott = temp * ( pre00 / pres )**(rtot/cptot)
1238 
1239  return
1240  end subroutine atmos_thermodyn_temp_pres2pott_0d
1241 
1242  !-----------------------------------------------------------------------------
1243 !OCL SERIAL
1244  subroutine atmos_thermodyn_temp_pres2pott_1d( &
1245  KA, KS, KE, &
1246  temp, pres, CPtot, Rtot, &
1247  pott )
1248  !$acc routine vector
1249  implicit none
1250  integer, intent(in) :: ka, ks, ke
1251 
1252  real(rp), intent(in) :: temp (ka) ! temperature [K]
1253  real(rp), intent(in) :: pres (ka) ! pressure [Pa]
1254  real(rp), intent(in) :: cptot(ka) ! specific heat [J/kg/K]
1255  real(rp), intent(in) :: rtot (ka) ! gas constant [J/kg/K]
1256 
1257  real(rp), intent(out) :: pott(ka) ! potential temperature [K]
1258 
1259  integer :: k
1260  !---------------------------------------------------------------------------
1261 
1262  do k = ks, ke
1263  call atmos_thermodyn_temp_pres2pott( temp(k), pres(k), cptot(k), rtot(k), pott(k) )
1264  enddo
1265 
1266  return
1267  end subroutine atmos_thermodyn_temp_pres2pott_1d
1268 
1269  !-----------------------------------------------------------------------------
1270  subroutine atmos_thermodyn_temp_pres2pott_2d( &
1271  IA, IS, IE, JA, JS, JE, &
1272  temp, pres, CPtot, Rtot, &
1273  pott )
1274  implicit none
1275  integer, intent(in) :: ia, is, ie
1276  integer, intent(in) :: ja, js, je
1277 
1278  real(rp), intent(in) :: temp (ia,ja) ! temperature [K]
1279  real(rp), intent(in) :: pres (ia,ja) ! pressure [Pa]
1280  real(rp), intent(in) :: cptot(ia,ja) ! specific heat [J/kg/K]
1281  real(rp), intent(in) :: rtot (ia,ja) ! gas constant [J/kg/K]
1282 
1283  real(rp), intent(out) :: pott(ia,ja) ! potential temperature [K]
1284 
1285  integer :: i, j
1286  !---------------------------------------------------------------------------
1287 
1288  !$omp parallel do OMP_SCHEDULE_ &
1289  !$omp private(i,j)
1290  !$acc kernels copyin(temp, pres, CPtot, Rtot) copyout(pott)
1291  do j = js, je
1292  do i = is, ie
1293  call atmos_thermodyn_temp_pres2pott( temp(i,j), pres(i,j), cptot(i,j), rtot(i,j), pott(i,j) )
1294  enddo
1295  enddo
1296  !$acc end kernels
1297 
1298  return
1299  end subroutine atmos_thermodyn_temp_pres2pott_2d
1300 
1301  !-----------------------------------------------------------------------------
1302  subroutine atmos_thermodyn_temp_pres2pott_3d( &
1303  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1304  temp, pres, CPtot, Rtot, &
1305  pott )
1306  implicit none
1307  integer, intent(in) :: ka, ks, ke
1308  integer, intent(in) :: ia, is, ie
1309  integer, intent(in) :: ja, js, je
1310 
1311  real(rp), intent(in) :: temp (ka,ia,ja) ! temperature [K]
1312  real(rp), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
1313  real(rp), intent(in) :: cptot(ka,ia,ja) ! specific heat [J/kg/K]
1314  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant [J/kg/K]
1315 
1316  real(rp), intent(out) :: pott(ka,ia,ja) ! potential temperature [K]
1317 
1318  integer :: k, i, j
1319  !---------------------------------------------------------------------------
1320 
1321  !$omp parallel do OMP_SCHEDULE_ collapse(2) &
1322  !$omp private(i,j,k)
1323  !$acc kernels copyin(temp, pres, CPtot, Rtot) copyout(pott)
1324  do j = js, je
1325  do i = is, ie
1326  do k = ks, ke
1327  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) )
1328  enddo
1329  enddo
1330  enddo
1331  !$acc end kernels
1332 
1333  return
1334  end subroutine atmos_thermodyn_temp_pres2pott_3d
1335 
1336  !-----------------------------------------------------------------------------
1337  subroutine atmos_thermodyn_temp_pres2ein_0d( &
1338  temp, pres, CVtot, Rtot, &
1339  ein, dens )
1340  !$acc routine seq
1341  implicit none
1342  real(rp), intent(in) :: temp ! temperature [K]
1343  real(rp), intent(in) :: pres ! pressure [Pa]
1344  real(rp), intent(in) :: cvtot ! specific heat [J/kg/K]
1345  real(rp), intent(in) :: rtot ! gas constant [J/kg/K]
1346 
1347  real(rp), intent(out) :: ein ! internal energy
1348  real(rp), intent(out) :: dens
1349  !---------------------------------------------------------------------------
1350 
1351  dens = pres / temp / rtot
1352  ein = temp * cvtot
1353 
1354  return
1355  end subroutine atmos_thermodyn_temp_pres2ein_0d
1356 
1357  !-----------------------------------------------------------------------------
1358 !OCL SERIAL
1359  subroutine atmos_thermodyn_temp_pres2ein_1d( &
1360  KA, KS, KE, &
1361  temp, pres, CVtot, Rtot, &
1362  ein, dens )
1363  !$acc routine vector
1364  implicit none
1365  integer, intent(in) :: ka, ks, ke
1366 
1367  real(rp), intent(in) :: temp (ka) ! temperature [K]
1368  real(rp), intent(in) :: pres (ka) ! pressure [Pa]
1369  real(rp), intent(in) :: cvtot(ka) ! specific heat [J/kg/K]
1370  real(rp), intent(in) :: rtot (ka) ! gas constant [J/kg/K]
1371 
1372  real(rp), intent(out) :: ein (ka) ! internal energy
1373  real(rp), intent(out) :: dens(ka)
1374  integer :: k
1375  !---------------------------------------------------------------------------
1376 
1377  do k = ks, ke
1378  call atmos_thermodyn_temp_pres2ein_0d( temp(k), pres(k), cvtot(k), rtot(k), & ! [IN]
1379  ein(k), dens(k) ) ! [OUT]
1380  end do
1381 
1382  return
1383  end subroutine atmos_thermodyn_temp_pres2ein_1d
1384 
1385  !-----------------------------------------------------------------------------
1386  subroutine atmos_thermodyn_temp_pres2ein_2d( &
1387  IA, IS, IE, JA, JS, JE, &
1388  temp, pres, CVtot, Rtot, &
1389  ein, dens )
1390  implicit none
1391  integer, intent(in) :: ia, is, ie
1392  integer, intent(in) :: ja, js, je
1393 
1394  real(rp), intent(in) :: temp (ia,ja) ! temperature [K]
1395  real(rp), intent(in) :: pres (ia,ja) ! pressure [Pa]
1396  real(rp), intent(in) :: cvtot(ia,ja) ! specific heat [J/kg/K]
1397  real(rp), intent(in) :: rtot (ia,ja) ! gas constant [J/kg/K]
1398 
1399  real(rp), intent(out) :: ein (ia,ja) ! internal energy
1400  real(rp), intent(out) :: dens(ia,ja)
1401 
1402  integer :: i, j
1403  !---------------------------------------------------------------------------
1404 
1405  !$omp parallel do
1406  !$acc kernels copyin(temp, pres, CVtot, Rtot) copyout(ein, dens)
1407  do j = js, je
1408  do i = is, ie
1409  call atmos_thermodyn_temp_pres2ein_0d( temp(i,j), pres(i,j), cvtot(i,j), rtot(i,j), & ! [IN]
1410  ein(i,j), dens(i,j) ) ! [OUT]
1411  end do
1412  end do
1413  !$acc end kernels
1414 
1415  return
1416  end subroutine atmos_thermodyn_temp_pres2ein_2d
1417 
1418  !-----------------------------------------------------------------------------
1419  subroutine atmos_thermodyn_temp_pres2ein_3d( &
1420  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1421  temp, pres, CVtot, Rtot, &
1422  ein, dens )
1423  implicit none
1424  integer, intent(in) :: ka, ks, ke
1425  integer, intent(in) :: ia, is, ie
1426  integer, intent(in) :: ja, js, je
1427 
1428  real(rp), intent(in) :: temp (ka,ia,ja) ! temperature [K]
1429  real(rp), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
1430  real(rp), intent(in) :: cvtot(ka,ia,ja) ! specific heat [J/kg/K]
1431  real(rp), intent(in) :: rtot (ka,ia,ja) ! gas constant [J/kg/K]
1432 
1433  real(rp), intent(out) :: ein (ka,ia,ja) ! internal energy
1434  real(rp), intent(out) :: dens(ka,ia,ja)
1435 
1436  integer :: k, i, j
1437  !---------------------------------------------------------------------------
1438 
1439  !$omp parallel do
1440  !$acc kernels copyin(temp, pres, CVtot, Rtot) copyout(ein, dens)
1441  do j = js, je
1442  do i = is, ie
1443  do k = ks, ke
1444  call atmos_thermodyn_temp_pres2ein_0d( temp(k,i,j), pres(k,i,j), cvtot(k,i,j), rtot(k,i,j), & ! [IN]
1445  ein(k,i,j), dens(k,i,j) ) ! [OUT]
1446  end do
1447  end do
1448  end do
1449  !$acc end kernels
1450 
1451  return
1452  end subroutine atmos_thermodyn_temp_pres2ein_3d
1453 
1454  !-----------------------------------------------------------------------------
1456  subroutine atmos_thermodyn_ein_pres2enth_0d( &
1457  Ein, PRES, DENS, &
1458  ENTH )
1459  !$acc routine seq
1460  implicit none
1461  real(rp), intent(in) :: ein
1462  real(rp), intent(in) :: pres
1463  real(rp), intent(in) :: dens
1464 
1465  real(rp), intent(out) :: enth
1466 
1467  enth = ein + pres / dens
1468 
1469  return
1470  end subroutine atmos_thermodyn_ein_pres2enth_0d
1471 
1472 !OCL SERIAL
1473  subroutine atmos_thermodyn_ein_pres2enth_1d( &
1474  KA, KS, KE, &
1475  Ein, PRES, DENS, &
1476  ENTH )
1477  !$acc routine vector
1478  implicit none
1479  integer, intent(in) :: ka, ks, ke
1480 
1481  real(rp), intent(in) :: ein (ka)
1482  real(rp), intent(in) :: pres(ka)
1483  real(rp), intent(in) :: dens(ka)
1484 
1485  real(rp), intent(out) :: enth(ka)
1486 
1487  integer :: k
1488 
1489  do k = ks, ke
1490  call atmos_thermodyn_ein_pres2enth_0d( ein(k), pres(k), dens(k), & ! [IN]
1491  enth(k) ) ! [OUT]
1492  end do
1493 
1494  return
1495  end subroutine atmos_thermodyn_ein_pres2enth_1d
1496 
1497  subroutine atmos_thermodyn_ein_pres2enth_2d( &
1498  IA, IS, IE, JA, JS, JE, &
1499  Ein, PRES, DENS, &
1500  ENTH )
1501  implicit none
1502  integer, intent(in) :: ia, is, ie
1503  integer, intent(in) :: ja, js, je
1504 
1505  real(rp), intent(in) :: ein (ia,ja)
1506  real(rp), intent(in) :: pres(ia,ja)
1507  real(rp), intent(in) :: dens(ia,ja)
1508 
1509  real(rp), intent(out) :: enth(ia,ja)
1510 
1511  integer :: i, j
1512 
1513  !$omp parallel do
1514  !$acc kernels copyin(Ein, PRES, DENS) copyout(ENTH)
1515  do j = js, je
1516  do i = is, ie
1517  call atmos_thermodyn_ein_pres2enth_0d( ein(i,j), pres(i,j), dens(i,j), & ! [IN]
1518  enth(i,j) ) ! [OUT]
1519  end do
1520  end do
1521  !$acc end kernels
1522 
1523  return
1524  end subroutine atmos_thermodyn_ein_pres2enth_2d
1525 
1526  subroutine atmos_thermodyn_ein_pres2enth_3d( &
1527  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1528  Ein, PRES, DENS, &
1529  ENTH )
1530  implicit none
1531  integer, intent(in) :: ka, ks, ke
1532  integer, intent(in) :: ia, is, ie
1533  integer, intent(in) :: ja, js, je
1534 
1535  real(rp), intent(in) :: ein (ka,ia,ja)
1536  real(rp), intent(in) :: pres(ka,ia,ja)
1537  real(rp), intent(in) :: dens(ka,ia,ja)
1538 
1539  real(rp), intent(out) :: enth(ka,ia,ja)
1540 
1541  integer :: k, i, j
1542 
1543  !$omp parallel do
1544  !$acc kernels copyin(Ein, PRES, DENS) copyout(ENTH)
1545  do j = js, je
1546  do i = is, ie
1547  do k = ks, ke
1548  call atmos_thermodyn_ein_pres2enth_0d( ein(k,i,j), pres(k,i,j), dens(k,i,j), & ! [IN]
1549  enth(k,i,j) ) ! [OUT]
1550  end do
1551  end do
1552  end do
1553  !$acc end kernels
1554 
1555  return
1556  end subroutine atmos_thermodyn_ein_pres2enth_3d
1557 
1558 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:61
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
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:59
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:97