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

Last change on this file since 1243 was 1123, checked in by jleconte, 11 years ago

correct a copy/paste bug from last commit and gfortran bug in watercommon

File size: 9.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         ! Ehouarn: gfortran says: Error: Result of EXP underflows its kind,
65         !          so set psat to the smallest possible value instead
66            psat=tiny(psat)
67         else                 
68            psat = Pref_solid_liquid*Exp(r3sublimation*(T-T_h2O_ice_liq)/(T-Trefsublimation)) ! solid / vapour
69         endif
70         if(psat.gt.p) then
71            qsat=1.
72         else
73            qsat=epsi*psat/(p-(1.-epsi)*psat)
74         endif
75         return
76      end subroutine Psat_water
77
78
79
80
81!==================================================================
82      subroutine Lcpdqsat_water(T,p,psat,qsat,dqsat,dlnpsat)
83
84         implicit none
85
86!==================================================================
87!     Purpose
88!     -------
89!     Compute dqsat=L/cp*d (q_sat)/d T and dlnpsat=L/cp d(ln Psat)/d T
90!     for a given temperature (K)!
91!     Based on the Tetens formula from L.Li physical parametrization manual
92!
93!     Authors
94!     -------
95!     Jeremy Leconte (2012)
96!
97!==================================================================
98
99!        input
100         real T, p, psat, qsat
101 
102!        output
103         real dqsat,dlnpsat
104
105! JL12 variables for tetens formula
106         real,parameter :: Pref_solid_liquid=611.14
107         real,parameter :: Trefvaporization=35.86
108         real,parameter :: Tmin=8.
109         real,parameter :: Trefsublimation=7.66
110         real,parameter :: r3vaporization=17.269
111         real,parameter :: r3sublimation=21.875
112
113         real :: dummy
114
115         if (psat.gt.p) then
116            dqsat=0.
117            return
118         endif
119
120         if (T.gt.T_h2O_ice_liq) then
121            dummy = r3vaporization*(T_h2O_ice_liq-Trefvaporization)/(T-Trefvaporization)**2  ! liquid / vapour
122         else if (T.lt.Tmin) then
123            print*, "careful, T<Tmin in Lcp psat water"
124            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(Tmin-Trefsublimation)**2  ! solid / vapour
125         else               
126            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(T-Trefsublimation)**2  ! solid / vapour
127         endif
128
129         dqsat=RLVTT/RCPD*qsat*(p/(p-(1.-epsi)*psat))*dummy
130         dlnpsat=RLVTT/RCPD*dummy
131         return
132      end subroutine Lcpdqsat_water
133
134
135
136
137!==================================================================
138      subroutine Tsat_water(p,Tsat)
139
140         implicit none
141
142!==================================================================
143!     Purpose
144!     -------
145!     Compute the saturation temperature
146!     for a given pressure (Pa)
147!     Based on the Tetens formula from L.Li physical parametrization manual
148!
149!     Authors
150!     -------
151!     Jeremy Leconte (2012)
152!
153!==================================================================
154
155!        input
156         real p
157 
158!        output
159         real Tsat
160
161! JL12 variables for tetens formula
162         real,parameter :: Pref_solid_liquid=611.14
163         real,parameter :: Trefvaporization=35.86
164         real,parameter :: Trefsublimation=7.66
165         real,parameter :: r3vaporization=17.269
166         real,parameter :: r3sublimation=21.875
167
168         if (p.lt.Pref_solid_liquid) then ! solid / vapour
169            Tsat =(T_h2O_ice_liq*r3sublimation- Trefsublimation*Log(p/Pref_solid_liquid))/(r3sublimation-Log(p/Pref_solid_liquid))
170         else                 ! liquid / vapour
171            Tsat =(T_h2O_ice_liq*r3vaporization- Trefvaporization*Log(p/Pref_solid_liquid))/(r3vaporization-Log(p/Pref_solid_liquid))
172         endif
173
174         return
175      end subroutine Tsat_water
176
177!==================================================================
178      subroutine Psat_waterDP(T,p,psat,qsat)
179
180         implicit none
181
182!==================================================================
183!     Purpose
184!     -------
185!     Compute the saturation vapor pressure and mass mixing ratio at saturation (kg/kg)
186!     for a given pressure (Pa) and temperature (K)
187!     DOUBLE PRECISION
188!     Based on the Tetens formula from L.Li physical parametrization manual
189!
190!     Authors
191!     -------
192!     Jeremy Leconte (2012)
193!
194!==================================================================
195
196!        input
197         double precision, intent(in) :: T, p
198 
199!        output
200         double precision psat,qsat
201
202! JL12 variables for tetens formula
203         double precision,parameter :: Pref_solid_liquid=611.14d0
204         double precision,parameter :: Trefvaporization=35.86D0
205         double precision,parameter :: Trefsublimation=7.66d0
206         double precision,parameter :: Tmin=8.d0
207         double precision,parameter :: r3vaporization=17.269d0
208         double precision,parameter :: r3sublimation=21.875d0
209
210! checked vs. old watersat data 14/05/2012 by JL.
211
212         if (T.gt.T_h2O_ice_liq) then
213            psat = Pref_solid_liquid*Exp(r3vaporization*(T-T_h2O_ice_liq)/(T-Trefvaporization)) ! liquid / vapour
214         else if (T.lt.Tmin) then
215            print*, "careful, T<Tmin in psat water"
216         !   psat = Pref_solid_liquid*Exp(r3sublimation*(Tmin-T_h2O_ice_liq)/(Tmin-Trefsublimation)) ! min psat 
217         ! Ehouarn: gfortran says: Error: Result of EXP underflows its kind,
218         !          so set psat to the smallest possible value instead
219            psat=tiny(psat)
220         else                 
221            psat = Pref_solid_liquid*Exp(r3sublimation*(T-T_h2O_ice_liq)/(T-Trefsublimation)) ! solid / vapour
222         endif
223         if(psat.gt.p) then
224            qsat=1.d0
225         else
226            qsat=epsi*psat/(p-(1.d0-epsi)*psat)
227         endif
228         return
229      end subroutine Psat_waterDP
230
231
232
233
234!==================================================================
235      subroutine Lcpdqsat_waterDP(T,p,psat,qsat,dqsat,dlnpsat)
236
237         implicit none
238
239!==================================================================
240!     Purpose
241!     -------
242!     Compute dqsat=L/cp*d (q_sat)/d T and dlnpsat=L/cp d(ln Psat)/d T
243!     for a given temperature (K)!
244!     Based on the Tetens formula from L.Li physical parametrization manual
245!
246!     Authors
247!     -------
248!     Jeremy Leconte (2012)
249!
250!==================================================================
251
252!        input
253         double precision T, p, psat, qsat
254 
255!        output
256         double precision dqsat,dlnpsat
257
258! JL12 variables for tetens formula
259         double precision,parameter :: Pref_solid_liquid=611.14d0
260         double precision,parameter :: Trefvaporization=35.86d0
261         double precision,parameter :: Tmin=8.d0
262         double precision,parameter :: Trefsublimation=7.66d0
263         double precision,parameter :: r3vaporization=17.269d0
264         double precision,parameter :: r3sublimation=21.875d0
265
266         double precision :: dummy
267
268         if (psat.gt.p) then
269            dqsat=0.d0
270            return
271         endif
272
273         if (T.gt.T_h2O_ice_liq) then
274            dummy = r3vaporization*(T_h2O_ice_liq-Trefvaporization)/(T-Trefvaporization)**2  ! liquid / vapour
275         else if (T.lt.Tmin) then
276            print*, "careful, T<Tmin in Lcp psat water"
277            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(Tmin-Trefsublimation)**2  ! solid / vapour
278         else               
279            dummy = r3sublimation*(T_h2O_ice_liq-Trefsublimation)/(T-Trefsublimation)**2  ! solid / vapour
280         endif
281
282         dqsat=RLVTT/RCPD*qsat*(p/(p-(1.d0-epsi)*psat))*dummy
283         dlnpsat=RLVTT/RCPD*dummy
284         return
285      end subroutine Lcpdqsat_waterDP
286
287
288end module watercommon_h
Note: See TracBrowser for help on using the repository browser.