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

Last change on this file since 2019 was 1993, checked in by jleconte, 7 years ago

29/08/2018 == JL

-watersat was used only in vdifc and thus it was not consistent with other routines (turbdiff, rain, largescale...) which used Psat_water from watercommon.
This is now harmonized. ALl routines use Psat_water. Watersat.F has been removed, but the routine is now in watercommon for archival purpose. It is not used anymore.
-also changed the number of chars for tname in the dyn3D/infotrac.F90 to be able to run rcm1d.

File size: 9.1 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      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
213!==================================================================
214subroutine watersat(T,p,qsat)
215  implicit none
216
217!==================================================================
218!     Purpose
219!     -------
220!     Compute the water mass mixing ratio at saturation (kg/kg)
221!     for a given pressure (Pa) and temperature (K)
222!     A replacement for the old watersat.F in the Martian GCM.
223!     Based on FCTTRE.h in the LMDTERRE model.
224!
225!     JL18 watersat was used only in vdifc and thus it was not consistent with other routines (turbdiff, rain, largescale...)
226!        which used Psat_water. This is now harmonized
227!        we put it here for archival purpose, but it is not used anymore.
228!
229!     Authors
230!     -------
231!     Robin Wordsworth (2010)
232!
233!==================================================================
234
235!   input
236  real T, p
237 
238!   output
239  real qsat
240
241! checked vs. NIST data 22/06/2010 by RW.
242! / by p gives partial pressure
243! x by epsi converts to mass mixing ratio
244
245  if (T.lt.T_h2O_ice_liq) then ! solid / vapour
246     qsat = 100.0 * 10**(2.07023 - 0.00320991             &
247          * T - 2484.896 / T + 3.56654 * alog10(T))
248  else                 ! liquid / vapour
249     qsat = 100.0 * 10**(23.8319 - 2948.964 / T - 5.028  &
250          * alog10(T) - 29810.16 * exp( -0.0699382 * T)  &
251          + 25.21935 * exp(-2999.924/T))
252  endif
253!  qsat=epsi*qsat/p
254  if(qsat.gt.p) then
255     qsat=1.
256  else
257     qsat=epsi*qsat/(p-(1.-epsi)*qsat)
258  endif
259  return
260end subroutine watersat
261
262subroutine watersat_grad(T,qsat,dqsat)
263
264  implicit none
265
266!==================================================================
267!     Purpose
268!     -------
269!     Compute the L/cp*d (q_sat)/d T
270!     for a given temperature (K)
271!
272!     JL18: see watersat
273!
274!     Authors
275!     -------
276!     Robin Wordsworth (2010)
277!
278!==================================================================
279
280!   input
281  real T,qsat
282 
283!   output
284  real dqsat
285
286!  if (T.lt.T_coup) then ! solid / vapour !why use T_coup?????????? JL12
287  if (T.lt.T_h2O_ice_liq) then ! solid / vapour
288     dqsat = RLVTT/RCPD*qsat*(3.56654/T             &
289          +2484.896*LOG(10.)/T**2                   &
290          -0.00320991*LOG(10.))
291  else                 ! liquid / vapour
292     dqsat = RLVTT/RCPD*qsat*LOG(10.)*              &
293          (2948.964/T**2-5.028/LOG(10.)/T           &
294          +25.21935*2999.924/T**2*EXP(-2999.924/T)  &
295          +29810.16*0.0699382*EXP(-0.0699382*T))
296  end if
297
298  return
299end subroutine watersat_grad
300
301
302
303end module watercommon_h
Note: See TracBrowser for help on using the repository browser.