30 public :: atmos_hydrometeor_lhv
31 public :: atmos_hydrometeor_lhv_para
32 public :: atmos_hydrometeor_lhs
33 public :: atmos_hydrometeor_lhs_para
34 public :: atmos_hydrometeor_lhf
35 public :: atmos_hydrometeor_lhf_para
36 public :: atmos_hydrometeor_entr
37 public :: atmos_hydrometeor_entr2temp
43 interface atmos_hydrometeor_lhv
46 module procedure atmos_hydrometeor_lhv_2d
47 module procedure atmos_hydrometeor_lhv_3d
48 end interface atmos_hydrometeor_lhv
50 interface atmos_hydrometeor_lhv_para
51 module procedure atmos_hydrometeor_lhv_1d_para
52 end interface atmos_hydrometeor_lhv_para
54 interface atmos_hydrometeor_lhs
55 module procedure atmos_hydrometeor_lhs_0d
57 module procedure atmos_hydrometeor_lhs_2d
58 module procedure atmos_hydrometeor_lhs_3d
59 end interface atmos_hydrometeor_lhs
61 interface atmos_hydrometeor_lhs_para
62 module procedure atmos_hydrometeor_lhs_1d_para
63 end interface atmos_hydrometeor_lhs_para
65 interface atmos_hydrometeor_lhf
66 module procedure atmos_hydrometeor_lhf_0d
68 module procedure atmos_hydrometeor_lhf_2d
69 module procedure atmos_hydrometeor_lhf_3d
70 end interface atmos_hydrometeor_lhf
72 interface atmos_hydrometeor_lhf_para
73 module procedure atmos_hydrometeor_lhf_1d_para
74 end interface atmos_hydrometeor_lhf_para
76 interface atmos_hydrometeor_entr
77 module procedure atmos_hydrometeor_entr_0d
79 module procedure atmos_hydrometeor_entr_3d
80 end interface atmos_hydrometeor_entr
82 interface atmos_hydrometeor_entr2temp
83 module procedure atmos_hydrometeor_entr2temp_0d
84 end interface atmos_hydrometeor_entr2temp
93 integer,
public ::
i_qv = -1
95 integer,
public,
parameter ::
n_hyd = 6
97 integer,
public,
parameter ::
i_hc = 1
98 integer,
public,
parameter ::
i_hr = 2
99 integer,
public,
parameter ::
i_hi = 3
100 integer,
public,
parameter ::
i_hs = 4
101 integer,
public,
parameter ::
i_hg = 5
102 integer,
public,
parameter ::
i_hh = 6
105 (/
"QC",
"QR",
"QI",
"QS",
"QG",
"QH" /)
107 (/
"cloud ",
"rain ",
"ice water",
"snow ",
"graupel ",
"hail " /)
109 (/
"NC",
"NR",
"NI",
"NS",
"NG",
"NH" /)
131 integer,
public ::
qha = 0
132 integer,
public ::
qhs = -1
133 integer,
public ::
qhe = -2
135 integer,
public ::
qla = 0
136 integer,
public ::
qls = -1
137 integer,
public ::
qle = -2
139 integer,
public ::
qia = 0
140 integer,
public ::
qis = -1
141 integer,
public ::
qie = -2
165 real(
rp),
private :: thermodyn_emask = 1.0_rp
168 logical,
private :: initialized = .false.
196 if ( initialized )
return
202 log_info(
"ATMOS_HYDROMETEOR_setup",*)
'Setup'
204 if ( thermodyn_type ==
'EXACT' )
then
216 thermodyn_emask = 1.0_rp
218 elseif( thermodyn_type ==
'SIMPLE' )
then
230 thermodyn_emask = 0.0_rp
233 log_error(
"ATMOS_HYDROMETEOR_setup",*)
'Not appropriate ATMOS_THERMODYN_ENERGY_TYPE. Check!', trim(thermodyn_type)
259 log_info(
"ATMOS_HYDROMETEOR_finalize",*)
'Finalize'
262 initialized = .false.
286 integer,
intent(in) :: nl
287 integer,
intent(in) :: ni
288 character(len=*),
intent(in) :: name(1+nl+ni)
289 character(len=*),
intent(in) :: desc(1+nl+ni)
290 character(len=*),
intent(in) :: unit(1+nl+ni)
292 integer,
intent(out) :: q0
294 logical,
intent(in),
optional :: advc(1+nl+ni)
296 real(
rp) :: cv (1+nl+ni)
297 real(
rp) :: cp (1+nl+ni)
298 real(
rp) :: r (1+nl+ni)
299 real(
rp) :: ei0 (1+nl+ni)
300 logical :: mass (1+nl+ni)
301 logical :: advc_(1+nl+ni)
308 log_error(
"ATMOS_HYDROMETEOR_regist",*)
'tracer for hydrometeor is already registerd'
343 if (
present(advc) )
then
397 real(RP),
intent(in) :: temp
398 real(RP),
intent(out) :: lhv
414 integer,
intent(in) :: KA, KS, KE
416 real(RP),
intent(in) :: temp(KA)
418 real(RP),
intent(out) :: lhv (KA)
431 subroutine atmos_hydrometeor_lhv_1d_para( &
436 integer,
intent(in) :: KA, KS, KE
438 real(RP),
intent(in) :: temp(KA)
440 real(RP),
intent(out) :: lhv (KA)
453 end subroutine atmos_hydrometeor_lhv_1d_para
456 subroutine atmos_hydrometeor_lhv_2d( &
457 IA, IS, IE, JA, JS, JE, &
461 integer,
intent(in) :: IA, IS, IE
462 integer,
intent(in) :: JA, JS, JE
464 real(RP),
intent(in) :: temp(IA,JA)
465 real(RP),
intent(out) :: lhv (IA,JA)
480 end subroutine atmos_hydrometeor_lhv_2d
483 subroutine atmos_hydrometeor_lhv_3d( &
490 integer,
intent(in) :: KA, KS, KE
491 integer,
intent(in) :: IA, IS, IE
492 integer,
intent(in) :: JA, JS, JE
494 real(RP),
intent(in) :: temp(KA,IA,JA)
495 real(RP),
intent(out) :: lhv (KA,IA,JA)
511 end subroutine atmos_hydrometeor_lhv_3d
514 subroutine atmos_hydrometeor_lhs_0d( &
522 real(RP),
intent(in) :: temp
523 real(RP),
intent(out) :: lhs
526 lhs = lhs0 + (
cp_vapor -
cp_ice ) * ( temp - tem00 ) * thermodyn_emask
529 end subroutine atmos_hydrometeor_lhs_0d
539 integer,
intent(in) :: KA, KS, KE
541 real(RP),
intent(in) :: temp(KA)
543 real(RP),
intent(out) :: lhs (KA)
549 call atmos_hydrometeor_lhs( temp(k), lhs(k) )
556 subroutine atmos_hydrometeor_lhs_1d_para( &
561 integer,
intent(in) :: KA, KS, KE
563 real(RP),
intent(in) :: temp(KA)
565 real(RP),
intent(out) :: lhs (KA)
573 call atmos_hydrometeor_lhs( temp(k), lhs(k) )
578 end subroutine atmos_hydrometeor_lhs_1d_para
581 subroutine atmos_hydrometeor_lhs_2d( &
582 IA, IS, IE, JA, JS, JE, &
586 integer,
intent(in) :: IA, IS, IE
587 integer,
intent(in) :: JA, JS, JE
589 real(RP),
intent(in) :: temp(IA,JA)
590 real(RP),
intent(out) :: lhs (IA,JA)
599 call atmos_hydrometeor_lhs( temp(i,j), lhs(i,j) )
605 end subroutine atmos_hydrometeor_lhs_2d
608 subroutine atmos_hydrometeor_lhs_3d( &
609 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
614 integer,
intent(in) :: KA, KS, KE
615 integer,
intent(in) :: IA, IS, IE
616 integer,
intent(in) :: JA, JS, JE
618 real(RP),
intent(in) :: temp(KA,IA,JA)
619 real(RP),
intent(out) :: lhs (KA,IA,JA)
629 call atmos_hydrometeor_lhs( temp(k,i,j), lhs(k,i,j) )
636 end subroutine atmos_hydrometeor_lhs_3d
639 subroutine atmos_hydrometeor_lhf_0d( &
647 real(RP),
intent(in) :: temp
648 real(RP),
intent(out) :: lhf
651 lhf = lhf0 + (
cp_water -
cp_ice ) * ( temp - tem00 ) * thermodyn_emask
654 end subroutine atmos_hydrometeor_lhf_0d
664 integer,
intent(in) :: KA, KS, KE
666 real(RP),
intent(in) :: temp(KA)
667 real(RP),
intent(out) :: lhf (KA)
673 call atmos_hydrometeor_lhf( temp(k), lhf(k) )
680 subroutine atmos_hydrometeor_lhf_1d_para( &
685 integer,
intent(in) :: KA, KS, KE
687 real(RP),
intent(in) :: temp(KA)
688 real(RP),
intent(out) :: lhf (KA)
696 call atmos_hydrometeor_lhf( temp(k), lhf(k) )
701 end subroutine atmos_hydrometeor_lhf_1d_para
704 subroutine atmos_hydrometeor_lhf_2d( &
705 IA, IS, IE, JA, JS, JE, &
709 integer,
intent(in) :: IA, IS, IE
710 integer,
intent(in) :: JA, JS, JE
712 real(RP),
intent(in) :: temp(IA,JA)
713 real(RP),
intent(out) :: lhf (IA,JA)
722 call atmos_hydrometeor_lhf( temp(i,j), lhf(i,j) )
728 end subroutine atmos_hydrometeor_lhf_2d
731 subroutine atmos_hydrometeor_lhf_3d( &
738 integer,
intent(in) :: KA, KS, KE
739 integer,
intent(in) :: IA, IS, IE
740 integer,
intent(in) :: JA, JS, JE
742 real(RP),
intent(in) :: temp(KA,IA,JA)
743 real(RP),
intent(out) :: lhf (KA,IA,JA)
753 call atmos_hydrometeor_lhf( temp(k,i,j), lhf(k,i,j) )
760 end subroutine atmos_hydrometeor_lhf_3d
764 subroutine atmos_hydrometeor_entr_0d( &
780 real(RP),
intent(in) :: TEMP
781 real(RP),
intent(in) :: PRES
782 real(RP),
intent(in) :: QV
783 real(RP),
intent(in) :: QI
784 real(RP),
intent(in) :: Qdry
785 real(RP),
intent(in) :: Rtot
786 real(RP),
intent(in) :: CPtot
788 real(RP),
intent(out) :: entr
790 real(RP) :: pres_dry, pres_vap
793 pres_dry = pres * qdry * rdry / rtot
794 pres_vap = pres * qv * rvap / rtot
796 entr = cptot * log( temp / tem00 ) &
797 - qdry * rdry * log( pres_dry / pre00 ) &
798 - qv * rvap * log( pres_vap / psat0 ) &
799 + ( qv * lhv0 - qi * lhf0 ) / tem00
802 end subroutine atmos_hydrometeor_entr_0d
807 IA, IS, IE, JA, JS, JE, &
813 integer,
intent(in) :: IA, IS, IE
814 integer,
intent(in) :: JA, JS, JE
816 real(RP),
intent(in) :: TEMP (IA,JA)
817 real(RP),
intent(in) :: PRES (IA,JA)
818 real(RP),
intent(in) :: QV (IA,JA)
819 real(RP),
intent(in) :: QI (IA,JA)
820 real(RP),
intent(in) :: Qdry (IA,JA)
821 real(RP),
intent(in) :: Rtot (IA,JA)
822 real(RP),
intent(in) :: CPtot(IA,JA)
824 real(RP),
intent(out) :: entr(IA,JA)
833 call atmos_hydrometeor_entr_0d( &
834 temp(i,j), pres(i,j), &
835 qv(i,j), qi(i,j), qdry(i,j), &
836 rtot(i,j), cptot(i,j), &
846 subroutine atmos_hydrometeor_entr_3d( &
847 KA, KS, KE, IA, IS, IE, JA, JS, JE, &
853 integer,
intent(in) :: KA, KS, KE
854 integer,
intent(in) :: IA, IS, IE
855 integer,
intent(in) :: JA, JS, JE
857 real(RP),
intent(in) :: TEMP (KA,IA,JA)
858 real(RP),
intent(in) :: PRES (KA,IA,JA)
859 real(RP),
intent(in) :: QV (KA,IA,JA)
860 real(RP),
intent(in) :: QI (KA,IA,JA)
861 real(RP),
intent(in) :: Qdry (KA,IA,JA)
862 real(RP),
intent(in) :: Rtot (KA,IA,JA)
863 real(RP),
intent(in) :: CPtot(KA,IA,JA)
865 real(RP),
intent(out) :: entr(KA,IA,JA)
875 call atmos_hydrometeor_entr_0d( &
876 temp(k,i,j), pres(k,i,j), &
877 qv(k,i,j), qi(k,i,j), qdry(k,i,j), &
878 rtot(k,i,j), cptot(k,i,j), &
886 end subroutine atmos_hydrometeor_entr_3d
890 subroutine atmos_hydrometeor_entr2temp_0d( &
906 real(RP),
intent(in) :: entr
907 real(RP),
intent(in) :: pres
908 real(RP),
intent(in) :: qv
909 real(RP),
intent(in) :: qi
910 real(RP),
intent(in) :: qdry
911 real(RP),
intent(in) :: Rtot
912 real(RP),
intent(in) :: CPtot
914 real(RP),
intent(out) :: temp
916 real(RP) :: pres_dry, pres_vap
919 pres_dry = pres * qdry * rdry / rtot
920 pres_vap = pres * qv * rvap / rtot
924 + qdry * rdry * log( pres_dry / pre00 ) &
925 + qv * rvap * log( pres_vap / psat0 ) &
926 - ( qv * lhv0 - qi * lhf0 ) / tem00 ) / cptot )
928 end subroutine atmos_hydrometeor_entr2temp_0d