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
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
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
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
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
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
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
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
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
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
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
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
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
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
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
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
160 log_info(
"ATMOS_THERMODYN_setup",*)
'Setup'
161 log_info(
"ATMOS_THERMODYN_setup",*)
'No namelists.'
170 subroutine atmos_thermodyn_qdry_0d( &
175 integer,
intent(in) :: qa
177 real(
rp),
intent(in) :: q(qa)
178 real(
rp),
intent(in) :: q_mass(qa)
180 real(
rp),
intent(out) :: qdry
188 qdry = qdry - q(iqw)*q_mass(iqw)
193 end subroutine atmos_thermodyn_qdry_0d
197 subroutine atmos_thermodyn_qdry_2d( &
198 IA, IS, IE, JA, JS, JE, QA, &
202 integer,
intent(in) :: ia, is, ie
203 integer,
intent(in) :: ja, js, je
204 integer,
intent(in) :: qa
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)
217 call atmos_thermodyn_qdry( qa, q(i,j,:), q_mass(:), qdry(i,j) )
222 end subroutine atmos_thermodyn_qdry_2d
226 subroutine atmos_thermodyn_qdry_3d( &
227 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
248 call atmos_thermodyn_qdry( qa, q(k,i,j,:), q_mass(:), qdry(k,i,j) )
253 end subroutine atmos_thermodyn_qdry_3d
259 subroutine atmos_thermodyn_cv_0d( &
264 integer,
intent(in) :: qa
266 real(
rp),
intent(in) :: q(qa)
267 real(
rp),
intent(in) :: cvq(qa)
268 real(
rp),
intent(in) :: qdry
270 real(
rp),
intent(out) :: cvtot
278 cvtot = cvtot + q(iqw) * cvq(iqw)
283 end subroutine atmos_thermodyn_cv_0d
287 subroutine atmos_thermodyn_cv_3d( &
288 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
301 real(
rp),
intent(out) :: cvtot(ka,ia,ja)
311 call atmos_thermodyn_cv( qa, q(k,i,j,:), cvq(:), qdry(k,i,j), cvtot(k,i,j) )
317 end subroutine atmos_thermodyn_cv_3d
323 subroutine atmos_thermodyn_cp_0d( &
328 integer,
intent(in) :: qa
330 real(
rp),
intent(in) :: q(qa)
331 real(
rp),
intent(in) :: cpq(qa)
332 real(
rp),
intent(in) :: qdry
334 real(
rp),
intent(out) :: cptot
342 cptot = cptot + q(iqw) * cpq(iqw)
347 end subroutine atmos_thermodyn_cp_0d
351 subroutine atmos_thermodyn_cp_3d( &
352 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
365 real(
rp),
intent(out) :: cptot(ka,ia,ja)
374 call atmos_thermodyn_cp( qa, q(k,i,j,:), cpq(:), qdry(k,i,j), cptot(k,i,j) )
380 end subroutine atmos_thermodyn_cp_3d
386 subroutine atmos_thermodyn_r_0d( &
391 integer,
intent(in) :: qa
393 real(
rp),
intent(in) :: q(qa)
394 real(
rp),
intent(in) :: rq(qa)
395 real(
rp),
intent(in) :: qdry
397 real(
rp),
intent(out) :: rtot
405 rtot = rtot + q(iqw) * rq(iqw)
410 end subroutine atmos_thermodyn_r_0d
414 subroutine atmos_thermodyn_r_3d( &
415 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
428 real(
rp),
intent(out) :: rtot(ka,ia,ja)
437 call atmos_thermodyn_r( qa, q(k,i,j,:), rq(:), qdry(k,i,j), rtot(k,i,j) )
443 end subroutine atmos_thermodyn_r_3d
451 subroutine atmos_thermodyn_specific_heat_0d( &
453 q, Mq, Rq, CVq, CPq, &
454 Qdry, Rtot, CVtot, CPtot )
455 integer,
intent(in) :: qa
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)
463 real(
rp),
intent(out) :: qdry
464 real(
rp),
intent(out) :: rtot
465 real(
rp),
intent(out) :: cvtot
466 real(
rp),
intent(out) :: cptot
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)
481 rtot = rtot + qdry * rdry
482 cvtot = cvtot + qdry * cvdry
483 cptot = cptot + qdry * cpdry
486 end subroutine atmos_thermodyn_specific_heat_0d
493 subroutine atmos_thermodyn_specific_heat_1d( &
496 q, Mq, Rq, CVq, CPq, &
497 Qdry, Rtot, CVtot, CPtot )
498 integer,
intent(in) :: ka, ks, ke
499 integer,
intent(in) :: qa
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)
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)
516 call atmos_thermodyn_specific_heat( qa, &
517 q(k,:), mq(:), rq(:), cvq(:), cpq(:), &
518 qdry(k), rtot(k), cvtot(k), cptot(k) )
522 end subroutine atmos_thermodyn_specific_heat_1d
528 subroutine atmos_thermodyn_specific_heat_2d( &
529 IA, IS, IE, JA, JS, JE, QA, &
532 Qdry, Rtot, CVtot, CPtot )
533 integer,
intent(in) :: ia, is, ie
534 integer,
intent(in) :: ja, js, je
535 integer,
intent(in) :: qa
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)
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)
555 call atmos_thermodyn_specific_heat( qa, q(i,j,:), mq(:), rq(:), cvq(:), cpq(:), &
556 qdry(i,j), rtot(i,j), cvtot(i,j), cptot(i,j) )
561 end subroutine atmos_thermodyn_specific_heat_2d
567 subroutine atmos_thermodyn_specific_heat_3d( &
568 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
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)
610 call atmos_thermodyn_specific_heat( qa, q(k,i,j,:), mq(:), rq(:), cvq(:), cpq(:), &
611 qdry(k,i,j), rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j) )
618 end subroutine atmos_thermodyn_specific_heat_3d
622 subroutine atmos_thermodyn_rhot2pres_0d( &
623 rhot, CVtot, CPtot, Rtot, &
626 real(
rp),
intent(in) :: rhot
627 real(
rp),
intent(in) :: cvtot
628 real(
rp),
intent(in) :: cptot
629 real(
rp),
intent(in) :: rtot
631 real(
rp),
intent(out) :: pres
634 pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
637 end subroutine atmos_thermodyn_rhot2pres_0d
641 subroutine atmos_thermodyn_rhot2pres_3d( &
642 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
643 rhot, CVtot, CPtot, Rtot, &
646 integer,
intent(in) :: ka, ks, ke
647 integer,
intent(in) :: ia, is, ie
648 integer,
intent(in) :: ja, js, je
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)
655 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
672 end subroutine atmos_thermodyn_rhot2pres_3d
676 subroutine atmos_thermodyn_rhot2rhoe_0d( &
677 rhot, CVtot, CPtot, Rtot, &
680 real(
rp),
intent(in) :: rhot
681 real(
rp),
intent(in) :: cvtot
682 real(
rp),
intent(in) :: cptot
683 real(
rp),
intent(in) :: rtot
685 real(
rp),
intent(out) :: rhoe
690 call atmos_thermodyn_rhot2pres( rhot, cvtot, cptot, rtot, pres )
692 rhoe = pres / rtot * cvtot
695 end subroutine atmos_thermodyn_rhot2rhoe_0d
699 subroutine atmos_thermodyn_rhot2rhoe_3d( &
700 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
701 rhot, CVtot, CPtot, Rtot, &
704 integer,
intent(in) :: ka, ks, ke
705 integer,
intent(in) :: ia, is, ie
706 integer,
intent(in) :: ja, js, je
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)
713 real(
rp),
intent(out) :: rhoe(ka,ia,ja)
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) )
730 end subroutine atmos_thermodyn_rhot2rhoe_3d
734 subroutine atmos_thermodyn_rhoe2rhot_0d( &
735 rhoe, CVtot, CPtot, Rtot, &
738 real(
rp),
intent(in) :: rhoe
739 real(
rp),
intent(in) :: cvtot
740 real(
rp),
intent(in) :: cptot
741 real(
rp),
intent(in) :: rtot
743 real(
rp),
intent(out) :: rhot
748 pres = rhoe * rtot / cvtot
750 rhot = rhoe / cvtot * ( pre00 / pres )**(rtot/cptot)
753 end subroutine atmos_thermodyn_rhoe2rhot_0d
757 subroutine atmos_thermodyn_rhoe2rhot_3d( &
758 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
759 rhoe, CVtot, CPtot, Rtot, &
762 integer,
intent(in) :: ka, ks, ke
763 integer,
intent(in) :: ia, is, ie
764 integer,
intent(in) :: ja, js, je
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)
771 real(
rp),
intent(out) :: rhot(ka,ia,ja)
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) )
788 end subroutine atmos_thermodyn_rhoe2rhot_3d
792 subroutine atmos_thermodyn_rhot2temp_pres_0d( &
793 dens, rhot, Rtot, CVtot, CPtot, &
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
802 real(
rp),
intent(out) :: temp
803 real(
rp),
intent(out) :: pres
806 pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
807 temp = pres / ( dens * rtot )
810 end subroutine atmos_thermodyn_rhot2temp_pres_0d
815 subroutine atmos_thermodyn_rhot2temp_pres_1d( &
817 dens, rhot, Rtot, CVtot, CPtot, &
819 integer,
intent(in) :: ka, ks, ke
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)
827 real(
rp),
intent(out) :: temp(ka)
828 real(
rp),
intent(out) :: pres(ka)
834 call atmos_thermodyn_rhot2temp_pres( dens(k), rhot(k), rtot(k), cvtot(k), cptot(k), &
839 end subroutine atmos_thermodyn_rhot2temp_pres_1d
843 subroutine atmos_thermodyn_rhot2temp_pres_3d( &
844 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
845 dens, rhot, Rtot, CVtot, CPtot, &
847 integer,
intent(in) :: ka, ks, ke
848 integer,
intent(in) :: ia, is, ie
849 integer,
intent(in) :: ja, js, je
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)
857 real(
rp),
intent(out) :: temp(ka,ia,ja)
858 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
876 end subroutine atmos_thermodyn_rhot2temp_pres_3d
882 subroutine atmos_thermodyn_rhoe2temp_pres_0d( &
883 dens, rhoe, CVtot, Rtot, &
886 real(
rp),
intent(in) :: dens
887 real(
rp),
intent(in) :: rhoe
888 real(
rp),
intent(in) :: cvtot
889 real(
rp),
intent(in) :: rtot
891 real(
rp),
intent(out) :: temp
892 real(
rp),
intent(out) :: pres
895 temp = rhoe / ( dens * cvtot )
896 pres = dens * rtot * temp
899 end subroutine atmos_thermodyn_rhoe2temp_pres_0d
903 subroutine atmos_thermodyn_rhoe2temp_pres_2d( &
904 IA, IS, IE, JA, JS, JE, &
905 dens, rhoe, CVtot, Rtot, &
908 integer,
intent(in) :: ia, is, ie
909 integer,
intent(in) :: ja, js, je
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)
916 real(
rp),
intent(out) :: temp(ia,ja)
917 real(
rp),
intent(out) :: pres(ia,ja)
927 call atmos_thermodyn_rhoe2temp_pres( dens(i,j), rhoe(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
932 end subroutine atmos_thermodyn_rhoe2temp_pres_2d
936 subroutine atmos_thermodyn_rhoe2temp_pres_3d( &
937 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
938 dens, rhoe, CVtot, Rtot, &
941 integer,
intent(in) :: ka, ks, ke
942 integer,
intent(in) :: ia, is, ie
943 integer,
intent(in) :: ja, js, je
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)
950 real(
rp),
intent(out) :: temp(ka,ia,ja)
951 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
968 end subroutine atmos_thermodyn_rhoe2temp_pres_3d
972 subroutine atmos_thermodyn_ein2temp_pres_0d( &
973 Ein, dens, CVtot, Rtot, &
976 real(
rp),
intent(in) :: ein
977 real(
rp),
intent(in) :: dens
978 real(
rp),
intent(in) :: cvtot
979 real(
rp),
intent(in) :: rtot
981 real(
rp),
intent(out) :: temp
982 real(
rp),
intent(out) :: pres
986 pres = dens * rtot * temp
989 end subroutine atmos_thermodyn_ein2temp_pres_0d
993 subroutine atmos_thermodyn_ein2temp_pres_2d( &
994 IA, IS, IE, JA, JS, JE, &
995 Ein, dens, CVtot, Rtot, &
998 integer,
intent(in) :: ia, is, ie
999 integer,
intent(in) :: ja, js, je
1001 real(
rp),
intent(in) :: ein (ia,ja)
1002 real(
rp),
intent(in) :: dens (ia,ja)
1003 real(
rp),
intent(in) :: cvtot(ia,ja)
1004 real(
rp),
intent(in) :: rtot (ia,ja)
1006 real(
rp),
intent(out) :: temp(ia,ja)
1007 real(
rp),
intent(out) :: pres(ia,ja)
1015 call atmos_thermodyn_ein2temp_pres( ein(i,j), dens(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
1020 end subroutine atmos_thermodyn_ein2temp_pres_2d
1024 subroutine atmos_thermodyn_ein2temp_pres_3d( &
1025 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1026 Ein, dens, CVtot, Rtot, &
1029 integer,
intent(in) :: ka, ks, ke
1030 integer,
intent(in) :: ia, is, ie
1031 integer,
intent(in) :: ja, js, je
1033 real(
rp),
intent(in) :: ein (ka,ia,ja)
1034 real(
rp),
intent(in) :: dens (ka,ia,ja)
1035 real(
rp),
intent(in) :: cvtot(ka,ia,ja)
1036 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1038 real(
rp),
intent(out) :: temp(ka,ia,ja)
1039 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
1054 end subroutine atmos_thermodyn_ein2temp_pres_3d
1058 subroutine atmos_thermodyn_pott2temp_pres_0d( &
1059 dens, pott, CVtot, CPtot, Rtot, &
1063 real(
rp),
intent(in) :: dens
1064 real(
rp),
intent(in) :: pott
1065 real(
rp),
intent(in) :: cvtot
1066 real(
rp),
intent(in) :: cptot
1067 real(
rp),
intent(in) :: rtot
1069 real(
rp),
intent(out) :: temp
1070 real(
rp),
intent(out) :: pres
1073 pres = pre00 * ( dens * pott * rtot / pre00 )**(cptot/cvtot)
1074 temp = pres / ( dens * rtot )
1077 end subroutine atmos_thermodyn_pott2temp_pres_0d
1081 subroutine atmos_thermodyn_pott2temp_pres_3d( &
1082 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1083 dens, pott, CVtot, CPtot, Rtot, &
1086 integer,
intent(in) :: ka, ks, ke
1087 integer,
intent(in) :: ia, is, ie
1088 integer,
intent(in) :: ja, js, je
1090 real(
rp),
intent(in) :: dens (ka,ia,ja)
1091 real(
rp),
intent(in) :: pott (ka,ia,ja)
1092 real(
rp),
intent(in) :: cvtot(ka,ia,ja)
1093 real(
rp),
intent(in) :: cptot(ka,ia,ja)
1094 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1096 real(
rp),
intent(out) :: temp(ka,ia,ja)
1097 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
1112 end subroutine atmos_thermodyn_pott2temp_pres_3d
1115 subroutine atmos_thermodyn_temp_pres2pott_0d( &
1116 temp, pres, CPtot, Rtot, &
1119 real(
rp),
intent(in) :: temp
1120 real(
rp),
intent(in) :: pres
1121 real(
rp),
intent(in) :: cptot
1122 real(
rp),
intent(in) :: rtot
1124 real(
rp),
intent(out) :: pott
1127 pott = temp * ( pre00 / pres )**(rtot/cptot)
1130 end subroutine atmos_thermodyn_temp_pres2pott_0d
1134 subroutine atmos_thermodyn_temp_pres2pott_1d( &
1136 temp, pres, CPtot, Rtot, &
1139 integer,
intent(in) :: ka, ks, ke
1141 real(
rp),
intent(in) :: temp (ka)
1142 real(
rp),
intent(in) :: pres (ka)
1143 real(
rp),
intent(in) :: cptot(ka)
1144 real(
rp),
intent(in) :: rtot (ka)
1146 real(
rp),
intent(out) :: pott(ka)
1152 call atmos_thermodyn_temp_pres2pott( temp(k), pres(k), cptot(k), rtot(k), pott(k) )
1156 end subroutine atmos_thermodyn_temp_pres2pott_1d
1159 subroutine atmos_thermodyn_temp_pres2pott_2d( &
1160 IA, IS, IE, JA, JS, JE, &
1161 temp, pres, CPtot, Rtot, &
1164 integer,
intent(in) :: ia, is, ie
1165 integer,
intent(in) :: ja, js, je
1167 real(
rp),
intent(in) :: temp (ia,ja)
1168 real(
rp),
intent(in) :: pres (ia,ja)
1169 real(
rp),
intent(in) :: cptot(ia,ja)
1170 real(
rp),
intent(in) :: rtot (ia,ja)
1172 real(
rp),
intent(out) :: pott(ia,ja)
1181 call atmos_thermodyn_temp_pres2pott( temp(i,j), pres(i,j), cptot(i,j), rtot(i,j), pott(i,j) )
1186 end subroutine atmos_thermodyn_temp_pres2pott_2d
1189 subroutine atmos_thermodyn_temp_pres2pott_3d( &
1190 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1191 temp, pres, CPtot, Rtot, &
1194 integer,
intent(in) :: ka, ks, ke
1195 integer,
intent(in) :: ia, is, ie
1196 integer,
intent(in) :: ja, js, je
1198 real(
rp),
intent(in) :: temp (ka,ia,ja)
1199 real(
rp),
intent(in) :: pres (ka,ia,ja)
1200 real(
rp),
intent(in) :: cptot(ka,ia,ja)
1201 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1203 real(
rp),
intent(out) :: pott(ka,ia,ja)
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) )
1219 end subroutine atmos_thermodyn_temp_pres2pott_3d
1222 subroutine atmos_thermodyn_temp_pres2ein_0d( &
1223 temp, pres, CVtot, Rtot, &
1226 real(
rp),
intent(in) :: temp
1227 real(
rp),
intent(in) :: pres
1228 real(
rp),
intent(in) :: cvtot
1229 real(
rp),
intent(in) :: rtot
1231 real(
rp),
intent(out) :: ein
1232 real(
rp),
intent(out) :: dens
1235 dens = pres / temp / rtot
1239 end subroutine atmos_thermodyn_temp_pres2ein_0d
1243 subroutine atmos_thermodyn_temp_pres2ein_1d( &
1245 temp, pres, CVtot, Rtot, &
1248 integer,
intent(in) :: ka, ks, ke
1250 real(
rp),
intent(in) :: temp (ka)
1251 real(
rp),
intent(in) :: pres (ka)
1252 real(
rp),
intent(in) :: cvtot(ka)
1253 real(
rp),
intent(in) :: rtot (ka)
1255 real(
rp),
intent(out) :: ein (ka)
1256 real(
rp),
intent(out) :: dens(ka)
1261 call atmos_thermodyn_temp_pres2ein_0d( temp(k), pres(k), cvtot(k), rtot(k), &
1266 end subroutine atmos_thermodyn_temp_pres2ein_1d
1269 subroutine atmos_thermodyn_temp_pres2ein_2d( &
1270 IA, IS, IE, JA, JS, JE, &
1271 temp, pres, CVtot, Rtot, &
1274 integer,
intent(in) :: ia, is, ie
1275 integer,
intent(in) :: ja, js, je
1277 real(
rp),
intent(in) :: temp (ia,ja)
1278 real(
rp),
intent(in) :: pres (ia,ja)
1279 real(
rp),
intent(in) :: cvtot(ia,ja)
1280 real(
rp),
intent(in) :: rtot (ia,ja)
1282 real(
rp),
intent(out) :: ein (ia,ja)
1283 real(
rp),
intent(out) :: dens(ia,ja)
1291 call atmos_thermodyn_temp_pres2ein_0d( temp(i,j), pres(i,j), cvtot(i,j), rtot(i,j), &
1292 ein(i,j), dens(i,j) )
1297 end subroutine atmos_thermodyn_temp_pres2ein_2d
1300 subroutine atmos_thermodyn_temp_pres2ein_3d( &
1301 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1302 temp, pres, CVtot, Rtot, &
1305 integer,
intent(in) :: ka, ks, ke
1306 integer,
intent(in) :: ia, is, ie
1307 integer,
intent(in) :: ja, js, je
1309 real(
rp),
intent(in) :: temp (ka,ia,ja)
1310 real(
rp),
intent(in) :: pres (ka,ia,ja)
1311 real(
rp),
intent(in) :: cvtot(ka,ia,ja)
1312 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1314 real(
rp),
intent(out) :: ein (ka,ia,ja)
1315 real(
rp),
intent(out) :: dens(ka,ia,ja)
1324 call atmos_thermodyn_temp_pres2ein_0d( temp(k,i,j), pres(k,i,j), cvtot(k,i,j), rtot(k,i,j), &
1325 ein(k,i,j), dens(k,i,j) )
1331 end subroutine atmos_thermodyn_temp_pres2ein_3d
1335 subroutine atmos_thermodyn_ein_pres2enth_0d( &
1338 real(
rp),
intent(in) :: ein
1339 real(
rp),
intent(in) :: pres
1340 real(
rp),
intent(in) :: dens
1342 real(
rp),
intent(out) :: enth
1344 enth = ein + pres / dens
1347 end subroutine atmos_thermodyn_ein_pres2enth_0d
1350 subroutine atmos_thermodyn_ein_pres2enth_1d( &
1354 integer,
intent(in) :: ka, ks, ke
1356 real(
rp),
intent(in) :: ein (ka)
1357 real(
rp),
intent(in) :: pres(ka)
1358 real(
rp),
intent(in) :: dens(ka)
1360 real(
rp),
intent(out) :: enth(ka)
1365 call atmos_thermodyn_ein_pres2enth_0d( ein(k), pres(k), dens(k), &
1370 end subroutine atmos_thermodyn_ein_pres2enth_1d
1372 subroutine atmos_thermodyn_ein_pres2enth_2d( &
1373 IA, IS, IE, JA, JS, JE, &
1376 integer,
intent(in) :: ia, is, ie
1377 integer,
intent(in) :: ja, js, je
1379 real(
rp),
intent(in) :: ein (ia,ja)
1380 real(
rp),
intent(in) :: pres(ia,ja)
1381 real(
rp),
intent(in) :: dens(ia,ja)
1383 real(
rp),
intent(out) :: enth(ia,ja)
1390 call atmos_thermodyn_ein_pres2enth_0d( ein(i,j), pres(i,j), dens(i,j), &
1396 end subroutine atmos_thermodyn_ein_pres2enth_2d
1398 subroutine atmos_thermodyn_ein_pres2enth_3d( &
1399 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1402 integer,
intent(in) :: ka, ks, ke
1403 integer,
intent(in) :: ia, is, ie
1404 integer,
intent(in) :: ja, js, je
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)
1410 real(
rp),
intent(out) :: enth(ka,ia,ja)
1418 call atmos_thermodyn_ein_pres2enth_0d( ein(k,i,j), pres(k,i,j), dens(k,i,j), &
1425 end subroutine atmos_thermodyn_ein_pres2enth_3d