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