source: trunk/LMDZ.GENERIC/libf/phystd/watercommon_h.F90 @ 848

Last change on this file since 848 was 786, checked in by jleconte, 12 years ago

19/09/2012 == JL

  • Correction in largescale to improve robustness when large water vapor amount
  • Correction in soil_setting to allow change of the number of subsurface layers


File size: 5.5 KB
Line 
1module watercommon_h
2
3      implicit none
4
5      real, parameter :: T_coup = 234.0
6      real, parameter :: T_h2O_ice_liq = 273.16
7      real, parameter :: T_h2O_ice_clouds = T_h2O_ice_liq-15.
8      real, parameter :: mH2O = 18.01528   
9
10      ! benjamin additions
11      real, parameter :: RLVTT = 2.257E+6 ! Latent heat of vaporization (J kg-1)
12      real, parameter :: RLSTT = 2.257E+6 ! 2.591E+6 in reality ! Latent heat of sublimation (J kg-1)
13
14      real, parameter :: RLFTT = 3.334E+5 ! Latent heat of fusion (J kg-1) ! entails an energy sink but better description of albedo
15      real, parameter :: rhowater = 1.0E+3 ! mass of water (kg/m^3)
16      real, parameter :: rhowaterice = 9.2E+2 ! mass of water (kg/m^3)
17      real, parameter :: capcal_h2o_liq = 4181.3 ! specific heat capacity of liquid water J/kg/K
18      real, parameter :: mx_eau_sol = 150 ! mass of water (kg/m^2)
19
20      real, save :: epsi, RCPD, RCPV, RV, RVTMP2
21     
22      contains
23
24     
25!==================================================================
26      subroutine Psat_water(T,p,psat,qsat)
27
28         implicit none
29
30!==================================================================
31!     Purpose
32!     -------
33!     Compute the saturation vapor pressure and mass mixing ratio at saturation (kg/kg)
34!     for a given pressure (Pa) and temperature (K)
35!     Based on the Tetens formula from L.Li physical parametrization manual
36!
37!     Authors
38!     -------
39!     Jeremy Leconte (2012)
40!
41!==================================================================
42
43!        input
44         real, intent(in) :: T, p
45 
46!        output
47         real psat,qsat
48
49! JL12 variables for tetens formula
50         real,parameter :: Pref_solid_liquid=611.14
51         real,parameter :: Trefvaporization=35.86
52         real,parameter :: Trefsublimation=7.66
53         real,parameter :: Tmin=8.
54         real,parameter :: r3vaporization=17.269
55         real,parameter :: r3sublimation=21.875
56
57! checked vs. old watersat data 14/05/2012 by JL.
58
59         if (T.gt.T_h2O_ice_liq) then
60            psat = Pref_solid_liquid*Exp(r3vaporization*(T-T_h2O_ice_liq)/(T-Trefvaporization)) ! liquid / vapour
61         else if (T.lt.Tmin) then
62            print*, "careful, T<Tmin in psat water"
63            psat = Pref_solid_liquid*Exp(r3sublimation*(Tmin-T_h2O_ice_liq)/(Tmin-Trefsublimation)) ! min psat 
64         else                 
65            psat = Pref_solid_liquid*Exp(r3sublimation*(T-T_h2O_ice_liq)/(T-Trefsublimation)) ! solid / vapour
66         endif
67         if(psat.gt.p) then
68            qsat=1.
69         else
70            qsat=epsi*psat/(p-(1.-epsi)*psat)
71         endif
72         return
73      end subroutine Psat_water
74
75
76
77
78!==================================================================
79      subroutine Lcpdqsat_water(T,p,psat,qsat,dqsat)
80
81         implicit none
82
83!==================================================================
84!     Purpose
85!     -------
86!     Compute L/cp*d (q_sat)/d T
87!     for a given temperature (K)!
88!     Based on the Tetens formula from L.Li physical parametrization manual
89!
90!     Authors
91!     -------
92!     Jeremy Leconte (2012)
93!
94!==================================================================
95
96!        input
97         real T, p, psat, qsat
98 
99!        output
100         real dqsat
101
102! JL12 variables for tetens formula
103         real,parameter :: Pref_solid_liquid=611.14
104         real,parameter :: Trefvaporization=35.86
105         real,parameter :: Tmin=8.
106         real,parameter :: Trefsublimation=7.66
107         real,parameter :: r3vaporization=17.269
108         real,parameter :: r3sublimation=21.875
109
110         real :: dummy
111
112         if (psat.gt.p) then
113            dqsat=0.
114            return
115         endif
116
117         if (T.gt.T_h2O_ice_liq) then
118            dummy = r3vaporization*(T_h2O_ice_liq-Trefvaporization)/(T-Trefvaporization)**2  ! liquid / vapour
119         else if (T.lt.Tmin) then
120            print*, "careful, T<Tmin in Lcp psat water"
121            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(Tmin-Trefsublimation)**2  ! solid / vapour
122         else               
123            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(T-Trefsublimation)**2  ! solid / vapour
124         endif
125
126         dqsat=RLVTT/RCPD*qsat*(p/(p-(1.-epsi)*psat))*dummy
127         return
128      end subroutine Lcpdqsat_water
129
130
131
132
133!==================================================================
134      subroutine Tsat_water(p,Tsat)
135
136         implicit none
137
138!==================================================================
139!     Purpose
140!     -------
141!     Compute the saturation temperature
142!     for a given pressure (Pa)
143!     Based on the Tetens formula from L.Li physical parametrization manual
144!
145!     Authors
146!     -------
147!     Jeremy Leconte (2012)
148!
149!==================================================================
150
151!        input
152         real p
153 
154!        output
155         real Tsat
156
157! JL12 variables for tetens formula
158         real,parameter :: Pref_solid_liquid=611.14
159         real,parameter :: Trefvaporization=35.86
160         real,parameter :: Trefsublimation=7.66
161         real,parameter :: r3vaporization=17.269
162         real,parameter :: r3sublimation=21.875
163
164         if (p.lt.Pref_solid_liquid) then ! solid / vapour
165            Tsat =(T_h2O_ice_liq*r3sublimation- Trefsublimation*Log(p/Pref_solid_liquid))/(r3sublimation-Log(p/Pref_solid_liquid))
166         else                 ! liquid / vapour
167            Tsat =(T_h2O_ice_liq*r3vaporization- Trefvaporization*Log(p/Pref_solid_liquid))/(r3vaporization-Log(p/Pref_solid_liquid))
168         endif
169
170         return
171      end subroutine Tsat_water
172
173
174end module watercommon_h
Note: See TracBrowser for help on using the repository browser.