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( &
176 integer,
intent(in) :: qa
178 real(
rp),
intent(in) :: q(qa)
179 real(
rp),
intent(in) :: q_mass(qa)
181 real(
rp),
intent(out) :: qdry
188 qdry = qdry - q(iqw)*q_mass(iqw)
192 end subroutine atmos_thermodyn_qdry_0d
196 subroutine atmos_thermodyn_qdry_2d( &
197 IA, IS, IE, JA, JS, JE, QA, &
201 integer,
intent(in) :: ia, is, ie
202 integer,
intent(in) :: ja, js, je
203 integer,
intent(in) :: qa
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)
223 qdry(i,j) = qdry(i,j) - q(i,j,iqw) * q_mass(iqw)
230 end subroutine atmos_thermodyn_qdry_2d
234 subroutine atmos_thermodyn_qdry_3d( &
235 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
248 integer :: k, i, j, iqw
263 qdry(k,i,j) = qdry(k,i,j) - q(k,i,j,iqw) * q_mass(iqw)
271 end subroutine atmos_thermodyn_qdry_3d
277 subroutine atmos_thermodyn_cv_0d( &
283 integer,
intent(in) :: qa
285 real(
rp),
intent(in) :: q(qa)
286 real(
rp),
intent(in) :: cvq(qa)
287 real(
rp),
intent(in) :: qdry
289 real(
rp),
intent(out) :: cvtot
296 cvtot = cvtot + q(iqw) * cvq(iqw)
300 end subroutine atmos_thermodyn_cv_0d
304 subroutine atmos_thermodyn_cv_3d( &
305 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
318 real(
rp),
intent(out) :: cvtot(ka,ia,ja)
320 integer :: k, i, j, iqw
330 cvtot(k,i,j) = qdry(k,i,j) * cvdry
335 cvtot(k,i,j) = cvtot(k,i,j) + q(k,i,j,iqw) * cvq(iqw)
343 end subroutine atmos_thermodyn_cv_3d
349 subroutine atmos_thermodyn_cp_0d( &
355 integer,
intent(in) :: qa
357 real(
rp),
intent(in) :: q(qa)
358 real(
rp),
intent(in) :: cpq(qa)
359 real(
rp),
intent(in) :: qdry
361 real(
rp),
intent(out) :: cptot
368 cptot = cptot + q(iqw) * cpq(iqw)
372 end subroutine atmos_thermodyn_cp_0d
376 subroutine atmos_thermodyn_cp_3d( &
377 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
390 real(
rp),
intent(out) :: cptot(ka,ia,ja)
392 integer :: k, i, j, iqw
401 cptot(k,i,j) = qdry(k,i,j) * cpdry
406 cptot(k,i,j) = cptot(k,i,j) + q(k,i,j,iqw) * cpq(iqw)
414 end subroutine atmos_thermodyn_cp_3d
420 subroutine atmos_thermodyn_r_0d( &
426 integer,
intent(in) :: qa
428 real(
rp),
intent(in) :: q(qa)
429 real(
rp),
intent(in) :: rq(qa)
430 real(
rp),
intent(in) :: qdry
432 real(
rp),
intent(out) :: rtot
439 rtot = rtot + q(iqw) * rq(iqw)
443 end subroutine atmos_thermodyn_r_0d
447 subroutine atmos_thermodyn_r_3d( &
448 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
461 real(
rp),
intent(out) :: rtot(ka,ia,ja)
463 integer :: k, i, j, iqw
472 rtot(k,i,j) = qdry(k,i,j) * rdry
477 rtot(k,i,j) = rtot(k,i,j) + q(k,i,j,iqw) * rq(iqw)
485 end subroutine atmos_thermodyn_r_3d
493 subroutine atmos_thermodyn_specific_heat_0d( &
495 q, Mq, Rq, CVq, CPq, &
496 Qdry, Rtot, CVtot, CPtot )
498 integer,
intent(in) :: qa
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)
506 real(
rp),
intent(out) :: qdry
507 real(
rp),
intent(out) :: rtot
508 real(
rp),
intent(out) :: cvtot
509 real(
rp),
intent(out) :: cptot
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)
524 rtot = rtot + qdry * rdry
525 cvtot = cvtot + qdry * cvdry
526 cptot = cptot + qdry * cpdry
529 end subroutine atmos_thermodyn_specific_heat_0d
536 subroutine atmos_thermodyn_specific_heat_1d( &
539 q, Mq, Rq, CVq, CPq, &
540 Qdry, Rtot, CVtot, CPtot )
542 integer,
intent(in) :: ka, ks, ke
543 integer,
intent(in) :: qa
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)
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)
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)
574 rtot(k) = rtot(k) + qdry(k) * rdry
575 cvtot(k) = cvtot(k) + qdry(k) * cvdry
576 cptot(k) = cptot(k) + qdry(k) * cpdry
580 end subroutine atmos_thermodyn_specific_heat_1d
586 subroutine atmos_thermodyn_specific_heat_2d( &
587 IA, IS, IE, JA, JS, JE, QA, &
590 Qdry, Rtot, CVtot, CPtot )
591 integer,
intent(in) :: ia, is, ie
592 integer,
intent(in) :: ja, js, je
593 integer,
intent(in) :: qa
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)
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)
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)
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
637 end subroutine atmos_thermodyn_specific_heat_2d
643 subroutine atmos_thermodyn_specific_heat_3d( &
644 KA, KS, KE, IA, IS, IE, JA, JS, JE, QA, &
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
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)
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)
664 integer :: k, i, j, iqw
675 cvtot(k,i,j) = 0.0_rp
676 cptot(k,i,j) = 0.0_rp
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)
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
697 end subroutine atmos_thermodyn_specific_heat_3d
701 subroutine atmos_thermodyn_rhot2pres_0d( &
702 rhot, CVtot, CPtot, Rtot, &
706 real(
rp),
intent(in) :: rhot
707 real(
rp),
intent(in) :: cvtot
708 real(
rp),
intent(in) :: cptot
709 real(
rp),
intent(in) :: rtot
711 real(
rp),
intent(out) :: pres
714 pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
717 end subroutine atmos_thermodyn_rhot2pres_0d
721 subroutine atmos_thermodyn_rhot2pres_3d( &
722 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
723 rhot, CVtot, CPtot, Rtot, &
726 integer,
intent(in) :: ka, ks, ke
727 integer,
intent(in) :: ia, is, ie
728 integer,
intent(in) :: ja, js, je
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)
735 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
754 end subroutine atmos_thermodyn_rhot2pres_3d
758 subroutine atmos_thermodyn_rhot2rhoe_0d( &
759 rhot, CVtot, CPtot, Rtot, &
763 real(
rp),
intent(in) :: rhot
764 real(
rp),
intent(in) :: cvtot
765 real(
rp),
intent(in) :: cptot
766 real(
rp),
intent(in) :: rtot
768 real(
rp),
intent(out) :: rhoe
773 call atmos_thermodyn_rhot2pres( rhot, cvtot, cptot, rtot, pres )
775 rhoe = pres / rtot * cvtot
778 end subroutine atmos_thermodyn_rhot2rhoe_0d
782 subroutine atmos_thermodyn_rhot2rhoe_3d( &
783 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
784 rhot, CVtot, CPtot, Rtot, &
787 integer,
intent(in) :: ka, ks, ke
788 integer,
intent(in) :: ia, is, ie
789 integer,
intent(in) :: ja, js, je
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)
796 real(
rp),
intent(out) :: rhoe(ka,ia,ja)
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) )
815 end subroutine atmos_thermodyn_rhot2rhoe_3d
819 subroutine atmos_thermodyn_rhoe2rhot_0d( &
820 rhoe, CVtot, CPtot, Rtot, &
824 real(
rp),
intent(in) :: rhoe
825 real(
rp),
intent(in) :: cvtot
826 real(
rp),
intent(in) :: cptot
827 real(
rp),
intent(in) :: rtot
829 real(
rp),
intent(out) :: rhot
834 pres = rhoe * rtot / cvtot
836 rhot = rhoe / cvtot * ( pre00 / pres )**(rtot/cptot)
839 end subroutine atmos_thermodyn_rhoe2rhot_0d
843 subroutine atmos_thermodyn_rhoe2rhot_3d( &
844 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
845 rhoe, CVtot, CPtot, Rtot, &
848 integer,
intent(in) :: ka, ks, ke
849 integer,
intent(in) :: ia, is, ie
850 integer,
intent(in) :: ja, js, je
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)
857 real(
rp),
intent(out) :: rhot(ka,ia,ja)
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) )
876 end subroutine atmos_thermodyn_rhoe2rhot_3d
880 subroutine atmos_thermodyn_rhot2temp_pres_0d( &
881 dens, rhot, Rtot, CVtot, CPtot, &
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
891 real(
rp),
intent(out) :: temp
892 real(
rp),
intent(out) :: pres
895 pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
896 temp = pres / ( dens * rtot )
899 end subroutine atmos_thermodyn_rhot2temp_pres_0d
904 subroutine atmos_thermodyn_rhot2temp_pres_1d( &
906 dens, rhot, Rtot, CVtot, CPtot, &
909 integer,
intent(in) :: ka, ks, ke
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)
917 real(
rp),
intent(out) :: temp(ka)
918 real(
rp),
intent(out) :: pres(ka)
924 call atmos_thermodyn_rhot2temp_pres( dens(k), rhot(k), rtot(k), cvtot(k), cptot(k), &
929 end subroutine atmos_thermodyn_rhot2temp_pres_1d
933 subroutine atmos_thermodyn_rhot2temp_pres_3d( &
934 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
935 dens, rhot, Rtot, CVtot, CPtot, &
937 integer,
intent(in) :: ka, ks, ke
938 integer,
intent(in) :: ia, is, ie
939 integer,
intent(in) :: ja, js, je
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)
947 real(
rp),
intent(out) :: temp(ka,ia,ja)
948 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
969 end subroutine atmos_thermodyn_rhot2temp_pres_3d
974 subroutine atmos_thermodyn_rhoe2temp_pres_0d( &
975 dens, rhoe, CVtot, Rtot, &
979 real(
rp),
intent(in) :: dens
980 real(
rp),
intent(in) :: rhoe
981 real(
rp),
intent(in) :: cvtot
982 real(
rp),
intent(in) :: rtot
984 real(
rp),
intent(out) :: temp
985 real(
rp),
intent(out) :: pres
988 temp = rhoe / ( dens * cvtot )
989 pres = dens * rtot * temp
992 end subroutine atmos_thermodyn_rhoe2temp_pres_0d
996 subroutine atmos_thermodyn_rhoe2temp_pres_2d( &
997 IA, IS, IE, JA, JS, JE, &
998 dens, rhoe, CVtot, Rtot, &
1001 integer,
intent(in) :: ia, is, ie
1002 integer,
intent(in) :: ja, js, je
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)
1009 real(
rp),
intent(out) :: temp(ia,ja)
1010 real(
rp),
intent(out) :: pres(ia,ja)
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)
1029 end subroutine atmos_thermodyn_rhoe2temp_pres_2d
1033 subroutine atmos_thermodyn_rhoe2temp_pres_3d( &
1034 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1035 dens, rhoe, CVtot, Rtot, &
1038 integer,
intent(in) :: ka, ks, ke
1039 integer,
intent(in) :: ia, is, ie
1040 integer,
intent(in) :: ja, js, je
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)
1047 real(
rp),
intent(out) :: temp(ka,ia,ja)
1048 real(
rp),
intent(out) :: pres(ka,ia,ja)
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)
1069 end subroutine atmos_thermodyn_rhoe2temp_pres_3d
1073 subroutine atmos_thermodyn_ein2temp_pres_0d( &
1074 Ein, dens, CVtot, Rtot, &
1078 real(
rp),
intent(in) :: ein
1079 real(
rp),
intent(in) :: dens
1080 real(
rp),
intent(in) :: cvtot
1081 real(
rp),
intent(in) :: rtot
1083 real(
rp),
intent(out) :: temp
1084 real(
rp),
intent(out) :: pres
1088 pres = dens * rtot * temp
1091 end subroutine atmos_thermodyn_ein2temp_pres_0d
1095 subroutine atmos_thermodyn_ein2temp_pres_2d( &
1096 IA, IS, IE, JA, JS, JE, &
1097 Ein, dens, CVtot, Rtot, &
1100 integer,
intent(in) :: ia, is, ie
1101 integer,
intent(in) :: ja, js, je
1103 real(
rp),
intent(in) :: ein (ia,ja)
1104 real(
rp),
intent(in) :: dens (ia,ja)
1105 real(
rp),
intent(in) :: cvtot(ia,ja)
1106 real(
rp),
intent(in) :: rtot (ia,ja)
1108 real(
rp),
intent(out) :: temp(ia,ja)
1109 real(
rp),
intent(out) :: pres(ia,ja)
1118 call atmos_thermodyn_ein2temp_pres( ein(i,j), dens(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
1124 end subroutine atmos_thermodyn_ein2temp_pres_2d
1128 subroutine atmos_thermodyn_ein2temp_pres_3d( &
1129 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1130 Ein, dens, CVtot, Rtot, &
1133 integer,
intent(in) :: ka, ks, ke
1134 integer,
intent(in) :: ia, is, ie
1135 integer,
intent(in) :: ja, js, je
1137 real(
rp),
intent(in) :: ein (ka,ia,ja)
1138 real(
rp),
intent(in) :: dens (ka,ia,ja)
1139 real(
rp),
intent(in) :: cvtot(ka,ia,ja)
1140 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1142 real(
rp),
intent(out) :: temp(ka,ia,ja)
1143 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
1160 end subroutine atmos_thermodyn_ein2temp_pres_3d
1164 subroutine atmos_thermodyn_pott2temp_pres_0d( &
1165 dens, pott, CVtot, CPtot, Rtot, &
1170 real(
rp),
intent(in) :: dens
1171 real(
rp),
intent(in) :: pott
1172 real(
rp),
intent(in) :: cvtot
1173 real(
rp),
intent(in) :: cptot
1174 real(
rp),
intent(in) :: rtot
1176 real(
rp),
intent(out) :: temp
1177 real(
rp),
intent(out) :: pres
1180 pres = pre00 * ( dens * pott * rtot / pre00 )**(cptot/cvtot)
1181 temp = pres / ( dens * rtot )
1184 end subroutine atmos_thermodyn_pott2temp_pres_0d
1188 subroutine atmos_thermodyn_pott2temp_pres_3d( &
1189 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1190 dens, pott, CVtot, CPtot, Rtot, &
1193 integer,
intent(in) :: ka, ks, ke
1194 integer,
intent(in) :: ia, is, ie
1195 integer,
intent(in) :: ja, js, je
1197 real(
rp),
intent(in) :: dens (ka,ia,ja)
1198 real(
rp),
intent(in) :: pott (ka,ia,ja)
1199 real(
rp),
intent(in) :: cvtot(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) :: temp(ka,ia,ja)
1204 real(
rp),
intent(out) :: pres(ka,ia,ja)
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) )
1221 end subroutine atmos_thermodyn_pott2temp_pres_3d
1224 subroutine atmos_thermodyn_temp_pres2pott_0d( &
1225 temp, pres, CPtot, Rtot, &
1229 real(
rp),
intent(in) :: temp
1230 real(
rp),
intent(in) :: pres
1231 real(
rp),
intent(in) :: cptot
1232 real(
rp),
intent(in) :: rtot
1234 real(
rp),
intent(out) :: pott
1237 pott = temp * ( pre00 / pres )**(rtot/cptot)
1240 end subroutine atmos_thermodyn_temp_pres2pott_0d
1244 subroutine atmos_thermodyn_temp_pres2pott_1d( &
1246 temp, pres, CPtot, Rtot, &
1250 integer,
intent(in) :: ka, ks, ke
1252 real(
rp),
intent(in) :: temp (ka)
1253 real(
rp),
intent(in) :: pres (ka)
1254 real(
rp),
intent(in) :: cptot(ka)
1255 real(
rp),
intent(in) :: rtot (ka)
1257 real(
rp),
intent(out) :: pott(ka)
1263 call atmos_thermodyn_temp_pres2pott( temp(k), pres(k), cptot(k), rtot(k), pott(k) )
1267 end subroutine atmos_thermodyn_temp_pres2pott_1d
1270 subroutine atmos_thermodyn_temp_pres2pott_2d( &
1271 IA, IS, IE, JA, JS, JE, &
1272 temp, pres, CPtot, Rtot, &
1275 integer,
intent(in) :: ia, is, ie
1276 integer,
intent(in) :: ja, js, je
1278 real(
rp),
intent(in) :: temp (ia,ja)
1279 real(
rp),
intent(in) :: pres (ia,ja)
1280 real(
rp),
intent(in) :: cptot(ia,ja)
1281 real(
rp),
intent(in) :: rtot (ia,ja)
1283 real(
rp),
intent(out) :: pott(ia,ja)
1293 call atmos_thermodyn_temp_pres2pott( temp(i,j), pres(i,j), cptot(i,j), rtot(i,j), pott(i,j) )
1299 end subroutine atmos_thermodyn_temp_pres2pott_2d
1302 subroutine atmos_thermodyn_temp_pres2pott_3d( &
1303 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1304 temp, pres, CPtot, Rtot, &
1307 integer,
intent(in) :: ka, ks, ke
1308 integer,
intent(in) :: ia, is, ie
1309 integer,
intent(in) :: ja, js, je
1311 real(
rp),
intent(in) :: temp (ka,ia,ja)
1312 real(
rp),
intent(in) :: pres (ka,ia,ja)
1313 real(
rp),
intent(in) :: cptot(ka,ia,ja)
1314 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1316 real(
rp),
intent(out) :: pott(ka,ia,ja)
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) )
1334 end subroutine atmos_thermodyn_temp_pres2pott_3d
1337 subroutine atmos_thermodyn_temp_pres2ein_0d( &
1338 temp, pres, CVtot, Rtot, &
1342 real(
rp),
intent(in) :: temp
1343 real(
rp),
intent(in) :: pres
1344 real(
rp),
intent(in) :: cvtot
1345 real(
rp),
intent(in) :: rtot
1347 real(
rp),
intent(out) :: ein
1348 real(
rp),
intent(out) :: dens
1351 dens = pres / temp / rtot
1355 end subroutine atmos_thermodyn_temp_pres2ein_0d
1359 subroutine atmos_thermodyn_temp_pres2ein_1d( &
1361 temp, pres, CVtot, Rtot, &
1365 integer,
intent(in) :: ka, ks, ke
1367 real(
rp),
intent(in) :: temp (ka)
1368 real(
rp),
intent(in) :: pres (ka)
1369 real(
rp),
intent(in) :: cvtot(ka)
1370 real(
rp),
intent(in) :: rtot (ka)
1372 real(
rp),
intent(out) :: ein (ka)
1373 real(
rp),
intent(out) :: dens(ka)
1378 call atmos_thermodyn_temp_pres2ein_0d( temp(k), pres(k), cvtot(k), rtot(k), &
1383 end subroutine atmos_thermodyn_temp_pres2ein_1d
1386 subroutine atmos_thermodyn_temp_pres2ein_2d( &
1387 IA, IS, IE, JA, JS, JE, &
1388 temp, pres, CVtot, Rtot, &
1391 integer,
intent(in) :: ia, is, ie
1392 integer,
intent(in) :: ja, js, je
1394 real(
rp),
intent(in) :: temp (ia,ja)
1395 real(
rp),
intent(in) :: pres (ia,ja)
1396 real(
rp),
intent(in) :: cvtot(ia,ja)
1397 real(
rp),
intent(in) :: rtot (ia,ja)
1399 real(
rp),
intent(out) :: ein (ia,ja)
1400 real(
rp),
intent(out) :: dens(ia,ja)
1409 call atmos_thermodyn_temp_pres2ein_0d( temp(i,j), pres(i,j), cvtot(i,j), rtot(i,j), &
1410 ein(i,j), dens(i,j) )
1416 end subroutine atmos_thermodyn_temp_pres2ein_2d
1419 subroutine atmos_thermodyn_temp_pres2ein_3d( &
1420 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1421 temp, pres, CVtot, Rtot, &
1424 integer,
intent(in) :: ka, ks, ke
1425 integer,
intent(in) :: ia, is, ie
1426 integer,
intent(in) :: ja, js, je
1428 real(
rp),
intent(in) :: temp (ka,ia,ja)
1429 real(
rp),
intent(in) :: pres (ka,ia,ja)
1430 real(
rp),
intent(in) :: cvtot(ka,ia,ja)
1431 real(
rp),
intent(in) :: rtot (ka,ia,ja)
1433 real(
rp),
intent(out) :: ein (ka,ia,ja)
1434 real(
rp),
intent(out) :: dens(ka,ia,ja)
1444 call atmos_thermodyn_temp_pres2ein_0d( temp(k,i,j), pres(k,i,j), cvtot(k,i,j), rtot(k,i,j), &
1445 ein(k,i,j), dens(k,i,j) )
1452 end subroutine atmos_thermodyn_temp_pres2ein_3d
1456 subroutine atmos_thermodyn_ein_pres2enth_0d( &
1461 real(
rp),
intent(in) :: ein
1462 real(
rp),
intent(in) :: pres
1463 real(
rp),
intent(in) :: dens
1465 real(
rp),
intent(out) :: enth
1467 enth = ein + pres / dens
1470 end subroutine atmos_thermodyn_ein_pres2enth_0d
1473 subroutine atmos_thermodyn_ein_pres2enth_1d( &
1479 integer,
intent(in) :: ka, ks, ke
1481 real(
rp),
intent(in) :: ein (ka)
1482 real(
rp),
intent(in) :: pres(ka)
1483 real(
rp),
intent(in) :: dens(ka)
1485 real(
rp),
intent(out) :: enth(ka)
1490 call atmos_thermodyn_ein_pres2enth_0d( ein(k), pres(k), dens(k), &
1495 end subroutine atmos_thermodyn_ein_pres2enth_1d
1497 subroutine atmos_thermodyn_ein_pres2enth_2d( &
1498 IA, IS, IE, JA, JS, JE, &
1502 integer,
intent(in) :: ia, is, ie
1503 integer,
intent(in) :: ja, js, je
1505 real(
rp),
intent(in) :: ein (ia,ja)
1506 real(
rp),
intent(in) :: pres(ia,ja)
1507 real(
rp),
intent(in) :: dens(ia,ja)
1509 real(
rp),
intent(out) :: enth(ia,ja)
1517 call atmos_thermodyn_ein_pres2enth_0d( ein(i,j), pres(i,j), dens(i,j), &
1524 end subroutine atmos_thermodyn_ein_pres2enth_2d
1526 subroutine atmos_thermodyn_ein_pres2enth_3d( &
1527 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1531 integer,
intent(in) :: ka, ks, ke
1532 integer,
intent(in) :: ia, is, ie
1533 integer,
intent(in) :: ja, js, je
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)
1539 real(
rp),
intent(out) :: enth(ka,ia,ja)
1548 call atmos_thermodyn_ein_pres2enth_0d( ein(k,i,j), pres(k,i,j), dens(k,i,j), &
1556 end subroutine atmos_thermodyn_ein_pres2enth_3d