17 #include "macro_thermodyn.h" 18 #include "inc_openmp.h" 45 public :: atmos_thermodyn_qd
46 public :: atmos_thermodyn_cv
47 public :: atmos_thermodyn_cp
48 public :: atmos_thermodyn_r
49 public :: atmos_thermodyn_pott
50 public :: atmos_thermodyn_rhoe
51 public :: atmos_thermodyn_rhot
52 public :: atmos_thermodyn_temp_pres
53 public :: atmos_thermodyn_temp_pres_e
55 interface atmos_thermodyn_qd
56 module procedure atmos_thermodyn_qd_0d
57 module procedure atmos_thermodyn_qd_3d
58 end interface atmos_thermodyn_qd
60 interface atmos_thermodyn_cv
61 module procedure atmos_thermodyn_cv_0d
62 module procedure atmos_thermodyn_cv_3d
63 end interface atmos_thermodyn_cv
65 interface atmos_thermodyn_cp
66 module procedure atmos_thermodyn_cp_0d
67 module procedure atmos_thermodyn_cp_3d
68 end interface atmos_thermodyn_cp
70 interface atmos_thermodyn_r
71 module procedure atmos_thermodyn_r_0d
72 module procedure atmos_thermodyn_r_3d
73 end interface atmos_thermodyn_r
75 interface atmos_thermodyn_rhoe
76 module procedure atmos_thermodyn_rhoe_0d
77 module procedure atmos_thermodyn_rhoe_3d
78 end interface atmos_thermodyn_rhoe
80 interface atmos_thermodyn_rhot
81 module procedure atmos_thermodyn_rhot_0d
82 module procedure atmos_thermodyn_rhot_3d
83 end interface atmos_thermodyn_rhot
85 interface atmos_thermodyn_pott
86 module procedure atmos_thermodyn_pott_0d
87 module procedure atmos_thermodyn_pott_1d
88 module procedure atmos_thermodyn_pott_3d
89 end interface atmos_thermodyn_pott
91 interface atmos_thermodyn_temp_pres
92 module procedure atmos_thermodyn_temp_pres_0d
93 module procedure atmos_thermodyn_temp_pres_3d
94 end interface atmos_thermodyn_temp_pres
96 interface atmos_thermodyn_temp_pres_e
97 module procedure atmos_thermodyn_temp_pres_e_0d
98 module procedure atmos_thermodyn_temp_pres_e_3d
99 end interface atmos_thermodyn_temp_pres_e
126 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[THERMODYN] / Categ[ATMOS SHARE] / Origin[SCALElib]' 134 subroutine atmos_thermodyn_qd_0d( &
140 real(RP),
intent(out) :: qdry
141 real(RP),
intent(in) :: q(
qa)
142 real(RP),
intent(in) :: q_mass(
qa)
150 qdry = qdry - q(iqw)*q_mass(iqw)
155 end subroutine atmos_thermodyn_qd_0d
159 subroutine atmos_thermodyn_qd_3d( &
165 real(RP),
intent(out) :: qdry (
ka,
ia,
ja)
166 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
167 real(RP),
intent(in) :: q_mass(
qa)
169 integer :: k, i, j, iqw
181 calc_qdry( qdry(k,i,j), q, q_mass, k, i, j, iqw )
186 end subroutine atmos_thermodyn_qd_3d
190 subroutine atmos_thermodyn_cv_0d( &
197 real(RP),
intent(out) :: cvtot
198 real(RP),
intent(in) :: q(
qa)
199 real(RP),
intent(in) :: cvq(
qa)
200 real(RP),
intent(in) :: qdry
208 cvtot = cvtot + q(iqw) * cvq(iqw)
213 end subroutine atmos_thermodyn_cv_0d
217 subroutine atmos_thermodyn_cv_3d( &
224 real(RP),
intent(out) :: cvtot(
ka,
ia,
ja)
225 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
226 real(RP),
intent(in) :: cvq (
qa)
227 real(RP),
intent(in) :: qdry (
ka,
ia,
ja)
229 integer :: k, i, j, iqw
241 calc_cv(cvtot(k,i,j), qdry(k,i,j), q, k, i, j, iqw, cvdry, cvq)
247 end subroutine atmos_thermodyn_cv_3d
251 subroutine atmos_thermodyn_cp_0d( &
258 real(RP),
intent(out) :: cptot
259 real(RP),
intent(in) :: q(
qa)
260 real(RP),
intent(in) :: cpq(
qa)
261 real(RP),
intent(in) :: qdry
269 cptot = cptot + q(iqw) * cpq(iqw)
274 end subroutine atmos_thermodyn_cp_0d
278 subroutine atmos_thermodyn_cp_3d( &
285 real(RP),
intent(out) :: cptot(
ka,
ia,
ja)
286 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
287 real(RP),
intent(in) :: cpq (
qa)
288 real(RP),
intent(in) :: qdry (
ka,
ia,
ja)
290 integer :: k, i, j, iqw
301 calc_cp(cptot(k,i,j), qdry(k,i,j), q, k, i, j, iqw, cpdry, cpq)
307 end subroutine atmos_thermodyn_cp_3d
311 subroutine atmos_thermodyn_r_0d( &
318 real(RP),
intent(out) :: rtot
319 real(RP),
intent(in) :: q(
qa)
320 real(RP),
intent(in) :: rq(
qa)
321 real(RP),
intent(in) :: qdry
329 rtot = rtot + q(iqw) * rq(iqw)
334 end subroutine atmos_thermodyn_r_0d
338 subroutine atmos_thermodyn_r_3d( &
345 real(RP),
intent(out) :: rtot(
ka,
ia,
ja)
346 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
347 real(RP),
intent(in) :: rq (
qa)
348 real(RP),
intent(in) :: qdry(
ka,
ia,
ja)
350 integer :: k, i, j, iqw
357 rtot(k,i,j) = qdry(k,i,j) * rdry
359 rtot(k,i,j) = rtot(k,i,j) + q(k,i,j,iqw) * rq(iqw)
366 end subroutine atmos_thermodyn_r_3d
369 subroutine atmos_thermodyn_pott_0d( &
379 real(RP),
intent(out) :: pott
380 real(RP),
intent(in) :: temp
381 real(RP),
intent(in) :: pres
382 real(RP),
intent(in) :: q(
qa)
383 real(RP),
intent(in) :: cvq (
qa)
384 real(RP),
intent(in) :: rq (
qa)
385 real(RP),
intent(in) :: mass(
qa)
388 real(RP) :: rtot, cvtot, rovcp
401 qdry = qdry - q(iqw) * mass(iqw)
402 cvtot = cvtot + q(iqw) * cvq(iqw)
403 rtot = rtot + q(iqw) * rq(iqw)
405 cvtot = cvdry * qdry + cvtot
406 rtot = rdry * qdry + rtot
408 rovcp = rtot / ( cvtot + rtot )
410 pott = temp * ( pre00 / pres )**rovcp
413 end subroutine atmos_thermodyn_pott_0d
416 subroutine atmos_thermodyn_pott_1d( &
426 real(RP),
intent(out) :: pott(
ka)
427 real(RP),
intent(in) :: temp(
ka)
428 real(RP),
intent(in) :: pres(
ka)
429 real(RP),
intent(in) :: q (
ka,
qa)
430 real(RP),
intent(in) :: cvq (
qa)
431 real(RP),
intent(in) :: rq (
qa)
432 real(RP),
intent(in) :: mass(
qa)
435 real(RP) :: rtot, cvtot, rovcp
449 qdry = qdry - q(k,iqw) * mass(iqw)
450 cvtot = cvtot + q(k,iqw) * cvq(iqw)
451 rtot = rtot + q(k,iqw) * rq(iqw)
453 cvtot = cvdry * qdry + cvtot
454 rtot = rdry * qdry + rtot
456 rovcp = rtot / ( cvtot + rtot )
458 pott(k) = temp(k) * ( pre00 / pres(k) )**rovcp
462 end subroutine atmos_thermodyn_pott_1d
465 subroutine atmos_thermodyn_pott_3d( &
475 real(RP),
intent(out) :: pott(
ka,
ia,
ja)
476 real(RP),
intent(in) :: temp(
ka,
ia,
ja)
477 real(RP),
intent(in) :: pres(
ka,
ia,
ja)
478 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
479 real(RP),
intent(in) :: cvq (
qa)
480 real(RP),
intent(in) :: rq (
qa)
481 real(RP),
intent(in) :: mass(
qa)
484 real(RP) :: rtot, cvtot, rovcp
486 integer :: k, i, j, iqw
500 qdry = qdry - q(k,i,j,iqw) * mass(iqw)
501 cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
502 rtot = rtot + q(k,i,j,iqw) * rq(iqw)
504 cvtot = cvdry * qdry + cvtot
505 rtot = rdry * qdry + rtot
507 rovcp = rtot / ( cvtot + rtot )
509 pott(k,i,j) = temp(k,i,j) * ( pre00 / pres(k,i,j) )**rovcp
515 end subroutine atmos_thermodyn_pott_3d
519 subroutine atmos_thermodyn_rhoe_0d( &
528 real(RP),
intent(out) :: rhoe
529 real(RP),
intent(in) :: rhot
530 real(RP),
intent(in) :: q (
qa)
531 real(RP),
intent(in) :: cvq (
qa)
532 real(RP),
intent(in) :: rq (
qa)
533 real(RP),
intent(in) :: mass(
qa)
537 real(RP) :: rtot, cvtot, cpovcv
550 qdry = qdry - q(iqw) * mass(iqw)
551 cvtot = cvtot + q(iqw) * cvq(iqw)
552 rtot = rtot + q(iqw) * rq(iqw)
554 cvtot = cvdry * qdry + cvtot
555 rtot = rdry * qdry + rtot
557 cpovcv = ( cvtot + rtot ) / cvtot
559 pres = pre00 * ( rhot * rtot / pre00 )**cpovcv
561 rhoe = pres / rtot * cvtot
564 end subroutine atmos_thermodyn_rhoe_0d
568 subroutine atmos_thermodyn_rhoe_3d( &
577 real(RP),
intent(out) :: rhoe(
ka,
ia,
ja)
578 real(RP),
intent(in) :: rhot(
ka,
ia,
ja)
579 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
580 real(RP),
intent(in) :: cvq (
qa)
581 real(RP),
intent(in) :: rq (
qa)
582 real(RP),
intent(in) :: mass(
qa)
586 real(RP) :: rtot, cvtot, cpovcv
588 integer :: k, i, j, iqw
604 qdry = qdry - q(k,i,j,iqw) * mass(iqw)
605 cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
606 rtot = rtot + q(k,i,j,iqw) * rq(iqw)
608 cvtot = cvdry * qdry + cvtot
609 rtot = rdry * qdry + rtot
611 cpovcv = ( cvtot + rtot ) / cvtot
613 pres = pre00 * ( rhot(k,i,j) * rtot / pre00 )**cpovcv
615 rhoe(k,i,j) = pres / rtot * cvtot
621 end subroutine atmos_thermodyn_rhoe_3d
625 subroutine atmos_thermodyn_rhot_0d( &
634 real(RP),
intent(out) :: rhot
635 real(RP),
intent(in) :: rhoe
636 real(RP),
intent(in) :: q (
qa)
637 real(RP),
intent(in) :: cvq (
qa)
638 real(RP),
intent(in) :: rq (
qa)
639 real(RP),
intent(in) :: mass(
qa)
643 real(RP) :: rtot, cvtot, rovcp
656 qdry = qdry - q(iqw) * mass(iqw)
657 cvtot = cvtot + q(iqw) * cvq(iqw)
658 rtot = rtot + q(iqw) * rq(iqw)
660 cvtot = cvdry * qdry + cvtot
661 rtot = rdry * qdry + rtot
663 rovcp = rtot / ( cvtot + rtot )
665 pres = rhoe * rtot / cvtot
667 rhot = rhoe / cvtot * ( pre00 / pres )**rovcp
670 end subroutine atmos_thermodyn_rhot_0d
674 subroutine atmos_thermodyn_rhot_3d( &
683 real(RP),
intent(out) :: rhot(
ka,
ia,
ja)
684 real(RP),
intent(in) :: rhoe(
ka,
ia,
ja)
685 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
686 real(RP),
intent(in) :: cvq (
qa)
687 real(RP),
intent(in) :: rq (
qa)
688 real(RP),
intent(in) :: mass(
qa)
692 real(RP) :: rtot, cvtot, rovcp
694 integer :: k, i, j, iqw
711 qdry = qdry - q(k,i,j,iqw) * mass(iqw)
712 cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
713 rtot = rtot + q(k,i,j,iqw) * rq(iqw)
715 cvtot = cvdry * qdry + cvtot
716 rtot = rdry * qdry + rtot
718 rovcp = rtot / ( cvtot + rtot )
720 pres = rhoe(k,i,j) * rtot / cvtot
722 rhot(k,i,j) = rhoe(k,i,j) / cvtot * ( pre00 / pres )**rovcp
727 end subroutine atmos_thermodyn_rhot_3d
731 subroutine atmos_thermodyn_temp_pres_0d( &
742 real(RP),
intent(out) :: temp
743 real(RP),
intent(out) :: pres
744 real(RP),
intent(in) :: dens
745 real(RP),
intent(in) :: rhot
746 real(RP),
intent(in) :: q (
qa)
747 real(RP),
intent(in) :: cvq (
qa)
748 real(RP),
intent(in) :: rq (
qa)
749 real(RP),
intent(in) :: mass(
qa)
752 real(RP) :: rtot, cvtot, cpovcv
765 qdry = qdry - q(iqw) * mass(iqw)
766 cvtot = cvtot + q(iqw) * cvq(iqw)
767 rtot = rtot + q(iqw) * rq(iqw)
769 cvtot = cvdry * qdry + cvtot
770 rtot = rdry * qdry + rtot
772 cpovcv = ( cvtot + rtot ) / cvtot
774 pres = pre00 * ( rhot * rtot / pre00 )**cpovcv
775 temp = pres / ( dens * rtot )
778 end subroutine atmos_thermodyn_temp_pres_0d
782 subroutine atmos_thermodyn_temp_pres_3d( &
793 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
794 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
795 real(RP),
intent(in) :: dens(
ka,
ia,
ja)
796 real(RP),
intent(in) :: rhot(
ka,
ia,
ja)
797 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
798 real(RP),
intent(in) :: cvq (
qa)
799 real(RP),
intent(in) :: rq (
qa)
800 real(RP),
intent(in) :: mass(
qa)
803 real(RP) :: rtot, cvtot, cpovcv
805 integer :: k, i, j, iqw
821 qdry = qdry - q(k,i,j,iqw) * mass(iqw)
822 cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
823 rtot = rtot + q(k,i,j,iqw) * rq(iqw)
825 cvtot = cvdry * qdry + cvtot
826 rtot = rdry * qdry + rtot
828 cpovcv = ( cvtot + rtot ) / cvtot
830 pres(k,i,j) = pre00 * ( rhot(k,i,j) * rtot / pre00 )**cpovcv
831 temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rtot )
836 end subroutine atmos_thermodyn_temp_pres_3d
840 subroutine atmos_thermodyn_temp_pres_e_0d( &
851 real(RP),
intent(out) :: temp
852 real(RP),
intent(out) :: pres
853 real(RP),
intent(in) :: dens
854 real(RP),
intent(in) :: rhoe
855 real(RP),
intent(in) :: q (
qa)
856 real(RP),
intent(in) :: cvq (
qa)
857 real(RP),
intent(in) :: rq (
qa)
858 real(RP),
intent(in) :: mass(
qa)
861 real(RP) :: rtot, cvtot
874 qdry = qdry - q(iqw) * mass(iqw)
875 cvtot = cvtot + q(iqw) * cvq(iqw)
876 rtot = rtot + q(iqw) * rq(iqw)
878 cvtot = cvdry * qdry + cvtot
879 rtot = rdry * qdry + rtot
882 temp = rhoe / ( dens * cvtot )
883 pres = dens * rtot * temp
886 end subroutine atmos_thermodyn_temp_pres_e_0d
890 subroutine atmos_thermodyn_temp_pres_e_3d( &
901 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
902 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
903 real(RP),
intent(in) :: dens(
ka,
ia,
ja)
904 real(RP),
intent(in) :: rhoe(
ka,
ia,
ja)
905 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
906 real(RP),
intent(in) :: cvq (
qa)
907 real(RP),
intent(in) :: rq (
qa)
908 real(RP),
intent(in) :: mass(
qa)
911 real(RP) :: rtot, cvtot
913 integer :: k, i, j, iqw
929 qdry = qdry - q(k,i,j,iqw) * mass(iqw)
930 cvtot = cvtot + q(k,i,j,iqw) * cvq(iqw)
931 rtot = rtot + q(k,i,j,iqw) * rq(iqw)
933 cvtot = cvdry * qdry + cvtot
934 rtot = rdry * qdry + rtot
937 temp(k,i,j) = rhoe(k,i,j) / ( dens(k,i,j) * cvtot )
938 pres(k,i,j) = dens(k,i,j) * rtot * temp(k,i,j)
944 end subroutine atmos_thermodyn_temp_pres_e_3d
949 Ein, dens, qdry, q, &
953 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
954 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
955 real(RP),
intent(in) :: ein (
ka,
ia,
ja)
956 real(RP),
intent(in) :: dens(
ka,
ia,
ja)
957 real(RP),
intent(in) :: qdry(
ka,
ia,
ja)
958 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
959 real(RP),
intent(in) :: cvq (
qa)
960 real(RP),
intent(in) :: rq (
qa)
962 real(RP) :: cv, rmoist
964 integer :: i, j, k, iqw
971 cv = qdry(k,i,j) * cvdry
972 rmoist = qdry(k,i,j) * rdry
974 cv = cv + q(k,i,j,iqw) * cvq(iqw)
975 rmoist = rmoist + q(k,i,j,iqw) * rq(iqw)
978 temp(k,i,j) = ein(k,i,j) / cv
979 pres(k,i,j) = dens(k,i,j) * rmoist * temp(k,i,j)
990 dens, pott, qdry, q, &
994 real(RP),
intent(out) :: temp(
ka,
ia,
ja)
995 real(RP),
intent(out) :: pres(
ka,
ia,
ja)
996 real(RP),
intent(in) :: dens(
ka,
ia,
ja)
997 real(RP),
intent(in) :: pott(
ka,
ia,
ja)
998 real(RP),
intent(in) :: qdry(
ka,
ia,
ja)
999 real(RP),
intent(in) :: q (
ka,
ia,
ja,
qa)
1000 real(RP),
intent(in) :: cpq (
qa)
1001 real(RP),
intent(in) :: rq (
qa)
1003 real(RP) :: rmoist, cp
1005 integer :: i, j, k, iqw
1012 cp = qdry(k,i,j) * cpdry
1013 rmoist = qdry(k,i,j) * rdry
1015 cp = cp + q(k,i,j,iqw) * cpq(iqw)
1016 rmoist = rmoist + q(k,i,j,iqw) * rq(iqw)
1018 calc_pre(pres(k,i,j), dens(k,i,j), pott(k,i,j), rmoist, cp, pre00)
1020 temp(k,i,j) = pres(k,i,j) / ( dens(k,i,j) * rmoist )
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]
logical, public io_l
output log or not? (this process)
subroutine, public atmos_thermodyn_tempre2(temp, pres, dens, pott, qdry, q, CPq, Rq)
integer, public ke
end point of inner domain: z, local
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
integer, public ia
of whole cells: x, local, with HALO
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_pre00
pressure reference [Pa]
subroutine, public atmos_thermodyn_setup
Setup.
integer, public ks
start point of inner domain: z, local
module ATMOSPHERE / Thermodynamics
subroutine, public atmos_thermodyn_tempre(temp, pres, Ein, dens, qdry, q, CVq, Rq)
integer, public io_fid_log
Log file ID.
integer, public ja
of whole cells: y, local, with HALO