source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/watercommon_h.F90 @ 263

Last change on this file since 263 was 227, checked in by milmd, 10 years ago

Last LMDZ version (1315) with OpenMP directives and other stuff

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