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)
596 call atmos_thermodyn_specific_heat( qa, q(k,i,j,:), mq(:), rq(:), cvq(:), cpq(:), &
597 qdry(k,i,j), rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j) )
603 end subroutine atmos_thermodyn_specific_heat_3d
607 subroutine atmos_thermodyn_rhot2pres_0d( &
608 rhot, CVtot, CPtot, Rtot, &
611 real(RP),
intent(in) :: rhot
612 real(RP),
intent(in) :: CVtot
613 real(RP),
intent(in) :: CPtot
614 real(RP),
intent(in) :: Rtot
616 real(RP),
intent(out) :: pres
619 pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
622 end subroutine atmos_thermodyn_rhot2pres_0d
626 subroutine atmos_thermodyn_rhot2pres_3d( &
627 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
628 rhot, CVtot, CPtot, Rtot, &
631 integer,
intent(in) :: KA, KS, KE
632 integer,
intent(in) :: IA, IS, IE
633 integer,
intent(in) :: JA, JS, JE
635 real(RP),
intent(in) :: rhot (ka,ia,ja)
636 real(RP),
intent(in) :: CVtot(ka,ia,ja)
637 real(RP),
intent(in) :: CPtot(ka,ia,ja)
638 real(RP),
intent(in) :: Rtot (ka,ia,ja)
640 real(RP),
intent(out) :: pres(ka,ia,ja)
651 call atmos_thermodyn_rhot2pres( rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), pres(k,i,j) )
657 end subroutine atmos_thermodyn_rhot2pres_3d
661 subroutine atmos_thermodyn_rhot2rhoe_0d( &
662 rhot, CVtot, CPtot, Rtot, &
665 real(RP),
intent(in) :: rhot
666 real(RP),
intent(in) :: CVtot
667 real(RP),
intent(in) :: CPtot
668 real(RP),
intent(in) :: Rtot
670 real(RP),
intent(out) :: rhoe
675 call atmos_thermodyn_rhot2pres( rhot, cvtot, cptot, rtot, pres )
677 rhoe = pres / rtot * cvtot
680 end subroutine atmos_thermodyn_rhot2rhoe_0d
684 subroutine atmos_thermodyn_rhot2rhoe_3d( &
685 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
686 rhot, CVtot, CPtot, Rtot, &
689 integer,
intent(in) :: KA, KS, KE
690 integer,
intent(in) :: IA, IS, IE
691 integer,
intent(in) :: JA, JS, JE
693 real(RP),
intent(in) :: rhot (ka,ia,ja)
694 real(RP),
intent(in) :: CVtot(ka,ia,ja)
695 real(RP),
intent(in) :: CPtot(ka,ia,ja)
696 real(RP),
intent(in) :: Rtot (ka,ia,ja)
698 real(RP),
intent(out) :: rhoe(ka,ia,ja)
709 call atmos_thermodyn_rhot2rhoe( rhot(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhoe(k,i,j) )
715 end subroutine atmos_thermodyn_rhot2rhoe_3d
719 subroutine atmos_thermodyn_rhoe2rhot_0d( &
720 rhoe, CVtot, CPtot, Rtot, &
723 real(RP),
intent(in) :: rhoe
724 real(RP),
intent(in) :: CVtot
725 real(RP),
intent(in) :: CPtot
726 real(RP),
intent(in) :: Rtot
728 real(RP),
intent(out) :: rhot
733 pres = rhoe * rtot / cvtot
735 rhot = rhoe / cvtot * ( pre00 / pres )**(rtot/cptot)
738 end subroutine atmos_thermodyn_rhoe2rhot_0d
742 subroutine atmos_thermodyn_rhoe2rhot_3d( &
743 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
744 rhoe, CVtot, CPtot, Rtot, &
747 integer,
intent(in) :: KA, KS, KE
748 integer,
intent(in) :: IA, IS, IE
749 integer,
intent(in) :: JA, JS, JE
751 real(RP),
intent(in) :: rhoe (ka,ia,ja)
752 real(RP),
intent(in) :: CVtot(ka,ia,ja)
753 real(RP),
intent(in) :: CPtot(ka,ia,ja)
754 real(RP),
intent(in) :: Rtot (ka,ia,ja)
756 real(RP),
intent(out) :: rhot(ka,ia,ja)
767 call atmos_thermodyn_rhoe2rhot( rhoe(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), rhot(k,i,j) )
773 end subroutine atmos_thermodyn_rhoe2rhot_3d
777 subroutine atmos_thermodyn_rhot2temp_pres_0d( &
778 dens, rhot, Rtot, CVtot, CPtot, &
781 real(RP),
intent(in) :: dens
782 real(RP),
intent(in) :: rhot
783 real(RP),
intent(in) :: Rtot
784 real(RP),
intent(in) :: CVtot
785 real(RP),
intent(in) :: CPtot
787 real(RP),
intent(out) :: temp
788 real(RP),
intent(out) :: pres
791 pres = pre00 * ( rhot * rtot / pre00 )**(cptot/cvtot)
792 temp = pres / ( dens * rtot )
795 end subroutine atmos_thermodyn_rhot2temp_pres_0d
800 subroutine atmos_thermodyn_rhot2temp_pres_1d( &
802 dens, rhot, Rtot, CVtot, CPtot, &
804 integer,
intent(in) :: KA, KS, KE
806 real(RP),
intent(in) :: dens (ka)
807 real(RP),
intent(in) :: rhot (ka)
808 real(RP),
intent(in) :: Rtot (ka)
809 real(RP),
intent(in) :: CVtot(ka)
810 real(RP),
intent(in) :: CPtot(ka)
812 real(RP),
intent(out) :: temp(ka)
813 real(RP),
intent(out) :: pres(ka)
819 call atmos_thermodyn_rhot2temp_pres( dens(k), rhot(k), rtot(k), cvtot(k), cptot(k), &
824 end subroutine atmos_thermodyn_rhot2temp_pres_1d
828 subroutine atmos_thermodyn_rhot2temp_pres_3d( &
829 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
830 dens, rhot, Rtot, CVtot, CPtot, &
832 integer,
intent(in) :: KA, KS, KE
833 integer,
intent(in) :: IA, IS, IE
834 integer,
intent(in) :: JA, JS, JE
836 real(RP),
intent(in) :: dens (ka,ia,ja)
837 real(RP),
intent(in) :: rhot (ka,ia,ja)
838 real(RP),
intent(in) :: Rtot (ka,ia,ja)
839 real(RP),
intent(in) :: CVtot(ka,ia,ja)
840 real(RP),
intent(in) :: CPtot(ka,ia,ja)
842 real(RP),
intent(out) :: temp(ka,ia,ja)
843 real(RP),
intent(out) :: pres(ka,ia,ja)
855 call atmos_thermodyn_rhot2temp_pres( dens(k,i,j), rhot(k,i,j), rtot(k,i,j), cvtot(k,i,j), cptot(k,i,j), &
856 temp(k,i,j), pres(k,i,j) )
861 end subroutine atmos_thermodyn_rhot2temp_pres_3d
867 subroutine atmos_thermodyn_rhoe2temp_pres_0d( &
868 dens, rhoe, CVtot, Rtot, &
871 real(RP),
intent(in) :: dens
872 real(RP),
intent(in) :: rhoe
873 real(RP),
intent(in) :: CVtot
874 real(RP),
intent(in) :: Rtot
876 real(RP),
intent(out) :: temp
877 real(RP),
intent(out) :: pres
880 temp = rhoe / ( dens * cvtot )
881 pres = dens * rtot * temp
884 end subroutine atmos_thermodyn_rhoe2temp_pres_0d
888 subroutine atmos_thermodyn_rhoe2temp_pres_2d( &
889 IA, IS, IE, JA, JS, JE, &
890 dens, rhoe, CVtot, Rtot, &
893 integer,
intent(in) :: IA, IS, IE
894 integer,
intent(in) :: JA, JS, JE
896 real(RP),
intent(in) :: dens (ia,ja)
897 real(RP),
intent(in) :: rhoe (ia,ja)
898 real(RP),
intent(in) :: CVtot(ia,ja)
899 real(RP),
intent(in) :: Rtot (ia,ja)
901 real(RP),
intent(out) :: temp(ia,ja)
902 real(RP),
intent(out) :: pres(ia,ja)
912 call atmos_thermodyn_rhoe2temp_pres( dens(i,j), rhoe(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
917 end subroutine atmos_thermodyn_rhoe2temp_pres_2d
921 subroutine atmos_thermodyn_rhoe2temp_pres_3d( &
922 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
923 dens, rhoe, CVtot, Rtot, &
926 integer,
intent(in) :: KA, KS, KE
927 integer,
intent(in) :: IA, IS, IE
928 integer,
intent(in) :: JA, JS, JE
930 real(RP),
intent(in) :: dens (ka,ia,ja)
931 real(RP),
intent(in) :: rhoe (ka,ia,ja)
932 real(RP),
intent(in) :: CVtot(ka,ia,ja)
933 real(RP),
intent(in) :: Rtot (ka,ia,ja)
935 real(RP),
intent(out) :: temp(ka,ia,ja)
936 real(RP),
intent(out) :: pres(ka,ia,ja)
947 call atmos_thermodyn_rhoe2temp_pres( dens(k,i,j), rhoe(k,i,j), cvtot(k,i,j), rtot(k,i,j), temp(k,i,j), pres(k,i,j) )
953 end subroutine atmos_thermodyn_rhoe2temp_pres_3d
957 subroutine atmos_thermodyn_ein2temp_pres_0d( &
958 Ein, dens, CVtot, Rtot, &
961 real(RP),
intent(in) :: Ein
962 real(RP),
intent(in) :: dens
963 real(RP),
intent(in) :: CVtot
964 real(RP),
intent(in) :: Rtot
966 real(RP),
intent(out) :: temp
967 real(RP),
intent(out) :: pres
971 pres = dens * rtot * temp
974 end subroutine atmos_thermodyn_ein2temp_pres_0d
978 subroutine atmos_thermodyn_ein2temp_pres_2d( &
979 IA, IS, IE, JA, JS, JE, &
980 Ein, dens, CVtot, Rtot, &
983 integer,
intent(in) :: IA, IS, IE
984 integer,
intent(in) :: JA, JS, JE
986 real(RP),
intent(in) :: Ein (ia,ja)
987 real(RP),
intent(in) :: dens (ia,ja)
988 real(RP),
intent(in) :: CVtot(ia,ja)
989 real(RP),
intent(in) :: Rtot (ia,ja)
991 real(RP),
intent(out) :: temp(ia,ja)
992 real(RP),
intent(out) :: pres(ia,ja)
1000 call atmos_thermodyn_ein2temp_pres( ein(i,j), dens(i,j), cvtot(i,j), rtot(i,j), temp(i,j), pres(i,j) )
1005 end subroutine atmos_thermodyn_ein2temp_pres_2d
1009 subroutine atmos_thermodyn_ein2temp_pres_3d( &
1010 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1011 Ein, dens, CVtot, Rtot, &
1014 integer,
intent(in) :: KA, KS, KE
1015 integer,
intent(in) :: IA, IS, IE
1016 integer,
intent(in) :: JA, JS, JE
1018 real(RP),
intent(in) :: Ein (ka,ia,ja)
1019 real(RP),
intent(in) :: dens (ka,ia,ja)
1020 real(RP),
intent(in) :: CVtot(ka,ia,ja)
1021 real(RP),
intent(in) :: Rtot (ka,ia,ja)
1023 real(RP),
intent(out) :: temp(ka,ia,ja)
1024 real(RP),
intent(out) :: pres(ka,ia,ja)
1033 call atmos_thermodyn_ein2temp_pres( ein(k,i,j), dens(k,i,j), cvtot(k,i,j), rtot(k,i,j), temp(k,i,j), pres(k,i,j) )
1039 end subroutine atmos_thermodyn_ein2temp_pres_3d
1043 subroutine atmos_thermodyn_pott2temp_pres_0d( &
1044 dens, pott, CVtot, CPtot, Rtot, &
1048 real(RP),
intent(in) :: dens
1049 real(RP),
intent(in) :: pott
1050 real(RP),
intent(in) :: CVtot
1051 real(RP),
intent(in) :: CPtot
1052 real(RP),
intent(in) :: Rtot
1054 real(RP),
intent(out) :: temp
1055 real(RP),
intent(out) :: pres
1058 pres = pre00 * ( dens * pott * rtot / pre00 )**(cptot/cvtot)
1059 temp = pres / ( dens * rtot )
1062 end subroutine atmos_thermodyn_pott2temp_pres_0d
1066 subroutine atmos_thermodyn_pott2temp_pres_3d( &
1067 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1068 dens, pott, CVtot, CPtot, Rtot, &
1071 integer,
intent(in) :: KA, KS, KE
1072 integer,
intent(in) :: IA, IS, IE
1073 integer,
intent(in) :: JA, JS, JE
1075 real(RP),
intent(in) :: dens (ka,ia,ja)
1076 real(RP),
intent(in) :: pott (ka,ia,ja)
1077 real(RP),
intent(in) :: CVtot(ka,ia,ja)
1078 real(RP),
intent(in) :: CPtot(ka,ia,ja)
1079 real(RP),
intent(in) :: Rtot (ka,ia,ja)
1081 real(RP),
intent(out) :: temp(ka,ia,ja)
1082 real(RP),
intent(out) :: pres(ka,ia,ja)
1091 call atmos_thermodyn_pott2temp_pres( dens(k,i,j), pott(k,i,j), cvtot(k,i,j), cptot(k,i,j), rtot(k,i,j), temp(k,i,j), pres(k,i,j) )
1097 end subroutine atmos_thermodyn_pott2temp_pres_3d
1100 subroutine atmos_thermodyn_temp_pres2pott_0d( &
1101 temp, pres, CPtot, Rtot, &
1104 real(RP),
intent(in) :: temp
1105 real(RP),
intent(in) :: pres
1106 real(RP),
intent(in) :: CPtot
1107 real(RP),
intent(in) :: Rtot
1109 real(RP),
intent(out) :: pott
1112 pott = temp * ( pre00 / pres )**(rtot/cptot)
1115 end subroutine atmos_thermodyn_temp_pres2pott_0d
1119 subroutine atmos_thermodyn_temp_pres2pott_1d( &
1121 temp, pres, CPtot, Rtot, &
1124 integer,
intent(in) :: KA, KS, KE
1126 real(RP),
intent(in) :: temp (ka)
1127 real(RP),
intent(in) :: pres (ka)
1128 real(RP),
intent(in) :: CPtot(ka)
1129 real(RP),
intent(in) :: Rtot (ka)
1131 real(RP),
intent(out) :: pott(ka)
1137 call atmos_thermodyn_temp_pres2pott( temp(k), pres(k), cptot(k), rtot(k), pott(k) )
1141 end subroutine atmos_thermodyn_temp_pres2pott_1d
1144 subroutine atmos_thermodyn_temp_pres2pott_2d( &
1145 IA, IS, IE, JA, JS, JE, &
1146 temp, pres, CPtot, Rtot, &
1149 integer,
intent(in) :: IA, IS, IE
1150 integer,
intent(in) :: JA, JS, JE
1152 real(RP),
intent(in) :: temp (ia,ja)
1153 real(RP),
intent(in) :: pres (ia,ja)
1154 real(RP),
intent(in) :: CPtot(ia,ja)
1155 real(RP),
intent(in) :: Rtot (ia,ja)
1157 real(RP),
intent(out) :: pott(ia,ja)
1166 call atmos_thermodyn_temp_pres2pott( temp(i,j), pres(i,j), cptot(i,j), rtot(i,j), pott(i,j) )
1171 end subroutine atmos_thermodyn_temp_pres2pott_2d
1174 subroutine atmos_thermodyn_temp_pres2pott_3d( &
1175 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1176 temp, pres, CPtot, Rtot, &
1179 integer,
intent(in) :: KA, KS, KE
1180 integer,
intent(in) :: IA, IS, IE
1181 integer,
intent(in) :: JA, JS, JE
1183 real(RP),
intent(in) :: temp (ka,ia,ja)
1184 real(RP),
intent(in) :: pres (ka,ia,ja)
1185 real(RP),
intent(in) :: CPtot(ka,ia,ja)
1186 real(RP),
intent(in) :: Rtot (ka,ia,ja)
1188 real(RP),
intent(out) :: pott(ka,ia,ja)
1198 call atmos_thermodyn_temp_pres2pott( temp(k,i,j), pres(k,i,j), cptot(k,i,j), rtot(k,i,j), pott(k,i,j) )
1204 end subroutine atmos_thermodyn_temp_pres2pott_3d
1207 subroutine atmos_thermodyn_temp_pres2ein_0d( &
1208 temp, pres, CVtot, Rtot, &
1211 real(RP),
intent(in) :: temp
1212 real(RP),
intent(in) :: pres
1213 real(RP),
intent(in) :: CVtot
1214 real(RP),
intent(in) :: Rtot
1216 real(RP),
intent(out) :: ein
1217 real(RP),
intent(out) :: dens
1220 dens = pres / temp / rtot
1224 end subroutine atmos_thermodyn_temp_pres2ein_0d
1228 subroutine atmos_thermodyn_temp_pres2ein_1d( &
1230 temp, pres, CVtot, Rtot, &
1233 integer,
intent(in) :: KA, KS, KE
1235 real(RP),
intent(in) :: temp (ka)
1236 real(RP),
intent(in) :: pres (ka)
1237 real(RP),
intent(in) :: CVtot(ka)
1238 real(RP),
intent(in) :: Rtot (ka)
1240 real(RP),
intent(out) :: ein (ka)
1241 real(RP),
intent(out) :: dens(ka)
1246 call atmos_thermodyn_temp_pres2ein_0d( temp(k), pres(k), cvtot(k), rtot(k), &
1251 end subroutine atmos_thermodyn_temp_pres2ein_1d
1254 subroutine atmos_thermodyn_temp_pres2ein_2d( &
1255 IA, IS, IE, JA, JS, JE, &
1256 temp, pres, CVtot, Rtot, &
1259 integer,
intent(in) :: IA, IS, IE
1260 integer,
intent(in) :: JA, JS, JE
1262 real(RP),
intent(in) :: temp (ia,ja)
1263 real(RP),
intent(in) :: pres (ia,ja)
1264 real(RP),
intent(in) :: CVtot(ia,ja)
1265 real(RP),
intent(in) :: Rtot (ia,ja)
1267 real(RP),
intent(out) :: ein (ia,ja)
1268 real(RP),
intent(out) :: dens(ia,ja)
1276 call atmos_thermodyn_temp_pres2ein_0d( temp(i,j), pres(i,j), cvtot(i,j), rtot(i,j), &
1277 ein(i,j), dens(i,j) )
1282 end subroutine atmos_thermodyn_temp_pres2ein_2d
1285 subroutine atmos_thermodyn_temp_pres2ein_3d( &
1286 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1287 temp, pres, CVtot, Rtot, &
1290 integer,
intent(in) :: KA, KS, KE
1291 integer,
intent(in) :: IA, IS, IE
1292 integer,
intent(in) :: JA, JS, JE
1294 real(RP),
intent(in) :: temp (ka,ia,ja)
1295 real(RP),
intent(in) :: pres (ka,ia,ja)
1296 real(RP),
intent(in) :: CVtot(ka,ia,ja)
1297 real(RP),
intent(in) :: Rtot (ka,ia,ja)
1299 real(RP),
intent(out) :: ein (ka,ia,ja)
1300 real(RP),
intent(out) :: dens(ka,ia,ja)
1309 call atmos_thermodyn_temp_pres2ein_0d( temp(k,i,j), pres(k,i,j), cvtot(k,i,j), rtot(k,i,j), &
1310 ein(k,i,j), dens(k,i,j) )
1316 end subroutine atmos_thermodyn_temp_pres2ein_3d
1320 subroutine atmos_thermodyn_ein_pres2enth_0d( &
1323 real(RP),
intent(in) :: Ein
1324 real(RP),
intent(in) :: PRES
1325 real(RP),
intent(in) :: DENS
1327 real(RP),
intent(out) :: ENTH
1329 enth = ein + pres / dens
1332 end subroutine atmos_thermodyn_ein_pres2enth_0d
1335 subroutine atmos_thermodyn_ein_pres2enth_1d( &
1339 integer,
intent(in) :: KA, KS, KE
1341 real(RP),
intent(in) :: Ein (ka)
1342 real(RP),
intent(in) :: PRES(ka)
1343 real(RP),
intent(in) :: DENS(ka)
1345 real(RP),
intent(out) :: ENTH(ka)
1350 call atmos_thermodyn_ein_pres2enth_0d( ein(k), pres(k), dens(k), &
1355 end subroutine atmos_thermodyn_ein_pres2enth_1d
1357 subroutine atmos_thermodyn_ein_pres2enth_2d( &
1358 IA, IS, IE, JA, JS, JE, &
1361 integer,
intent(in) :: IA, IS, IE
1362 integer,
intent(in) :: JA, JS, JE
1364 real(RP),
intent(in) :: Ein (ia,ja)
1365 real(RP),
intent(in) :: PRES(ia,ja)
1366 real(RP),
intent(in) :: DENS(ia,ja)
1368 real(RP),
intent(out) :: ENTH(ia,ja)
1375 call atmos_thermodyn_ein_pres2enth_0d( ein(i,j), pres(i,j), dens(i,j), &
1381 end subroutine atmos_thermodyn_ein_pres2enth_2d
1383 subroutine atmos_thermodyn_ein_pres2enth_3d( &
1384 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1387 integer,
intent(in) :: KA, KS, KE
1388 integer,
intent(in) :: IA, IS, IE
1389 integer,
intent(in) :: JA, JS, JE
1391 real(RP),
intent(in) :: Ein (ka,ia,ja)
1392 real(RP),
intent(in) :: PRES(ka,ia,ja)
1393 real(RP),
intent(in) :: DENS(ka,ia,ja)
1395 real(RP),
intent(out) :: ENTH(ka,ia,ja)
1403 call atmos_thermodyn_ein_pres2enth_0d( ein(k,i,j), pres(k,i,j), dens(k,i,j), &
1410 end subroutine atmos_thermodyn_ein_pres2enth_3d
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
real(rp), public const_pre00
pressure reference [Pa]
subroutine, public atmos_thermodyn_setup
Setup.
module atmosphere / thermodyn