source: trunk/LMDZ.TITAN/libf/phytitan/evapCH4.F90 @ 3529

Last change on this file since 3529 was 3497, checked in by debatzbr, 3 weeks ago

Add AC6H6 condensation in the microphysics

File size: 6.5 KB
Line 
1subroutine evapCH4(ngrid,nlayer,ptimestep,pplev,zzlay,zzlev, &
2                   u,v,tsurf,pqCH4, &
3                   tankCH4,pdqCH4,dtsurfevap)
4
5use radcommon_h, only: gzlat         ! Gravity of the planet [m.s-2]
6use planete_mod, only: z0            ! Surface roughness of the planet [m]
7use geometry_mod, only: latitude_deg ! Latitude grid of the planet [°]
8
9implicit none
10
11!====================================================================
12!
13!   Purpose
14!   -------
15!   Surface flux for methane evaporation.
16!   The routine calculates the surface flux of methane
17!   And the latente heat of methane evaporation.
18!
19!   Implicit scheme.
20!   We start by adding to variables x the physical tendencies already computed.
21!   We resolve the equation :
22!   x(t+1) =  x(t) + A * (dx/dt) * dt
23!
24!
25!   INPUT
26!   -----
27!   ngrid     = Number of grid points in the chunk [-]
28!   nlayer    = Number of vertical layers [-]
29!   ptimestep = Time step [s]
30!   pplev     = Intermediate pressure levels [Pa]
31!   zzlay     = Altitude at the middle of the atmospheric layers (ref : local surf) [m]
32!   zzlev     = Altitude at the atmospheric layer boundaries (ref : local surf) [m]
33!   u         = Zonal component of the wind [m.s-1]
34!   v         = Meridional component of the wind [m.s-1]
35!   tsurf     = Surface temperature [K]
36!   pqCH4     = Molar fraction of CH4 [mol.mol-1]
37!
38!
39!   OUTPUT
40!   ------
41!   tankCH4    = Depth of surface methane tank [m]
42!   pdqCH4     = Molar fraction tendency of CH4 at the surface [mol.mol-1.s-1]
43!   dtsurfevap = Evaporation heating rate at the surface [K.s-1]
44!
45!
46!   Author
47!   ------
48!   B. de Batz de Trenquelléon (10/2022)
49!
50!====================================================================   
51
52
53!------------------------------------
54! 0. DECLARATIONS
55!------------------------------------
56
57! Inputs :
58!---------
59integer, intent(in) :: ngrid
60integer, intent(in) :: nlayer
61real, intent(in)    :: ptimestep
62real, intent(in)    :: pplev(ngrid,nlayer+1)
63real, intent(in)    :: zzlay(ngrid,nlayer)
64real, intent(in)    :: zzlev(ngrid,nlayer+1)
65real, intent(in)    :: u(ngrid,nlayer)
66real, intent(in)    :: v(ngrid,nlayer)
67real, intent(in)    :: tsurf(ngrid)
68real, intent(in)    :: pqCH4(ngrid,nlayer)
69
70! Outputs :
71!----------
72real, intent(out) :: tankCH4(ngrid)
73real, intent(out) :: pdqCH4(ngrid)
74real, intent(out) :: dtsurfevap(ngrid)
75
76! Parameters :
77!-------------
78real, parameter :: karman = 0.4 ! Karman constant [-]
79real, parameter :: humCH4 = 0.4 ! Imposed surface humidity for CH4 [-]
80
81real, parameter :: Flnp = 0.10  ! Fraction occupied by lakes (North Pole)
82real, parameter :: Flsp = 0.01  ! Fraction occupied by lakes (South Pole)
83real, parameter :: Flml = 1.    ! Fraction not infiltrated into the ground (Mid latitudes)
84
85real, parameter :: mmolair = 28.e-3 ! Molar mass of air [kg.mol-1]
86real, parameter :: mmolCH4 = 16.e-3 ! Molar mass of CH4 [kg.mol-1]
87real, parameter :: rholCH4 = 425.   ! Density of liquid CH4 [kg.m-3]
88
89real, parameter :: TcCH4  = 190.56  ! Critical point of CH4 [K]
90real, parameter :: Cplake = 2689992 ! Surface heat capacity for hydrocarbon lakes [J.m-2.K-1] (Tokano 2005)
91
92! Local variables :
93!------------------
94integer :: ig
95
96! Initialisation :
97real*8  :: Tlake      ! Lake temperature [K]
98real*8  :: rhoair     ! Density of air [kg.m-3]
99real*8  :: ws         ! Horizontal wind at the surface [m.s-1]
100real*8  :: Cd         ! Turbulent term [-]
101real*8  :: qsatCH4    ! Saturation profile of CH4 [mol.mol-1]
102
103! Flux :
104real*8  :: flux       ! Surface flux [kg.m-2.s-1]
105real*8  :: fluxCH4    ! Surface flux fo CH4 [kg.m-2.s-1.mol.mol-1]
106
107! Variations of CH4 :
108real*8  :: newpqCH4   ! New molar fraction of CH4 in the first layer [mol.mol-1]
109real*8  :: dtankCH4   ! Trend of CH4's tank [m]
110
111! Latente heat :
112real*8  :: ftm, LvCH4 ! Variables for latente heat [-, J.kg-1]
113
114
115!------------------------------------
116! 1. INITIALISATION
117!------------------------------------
118
119do ig = 1, ngrid ! Main loop on the horizontal grid
120
121  ! Density of the first layer of the atmosphere [kg.m-3]
122  rhoair = (pplev(ig,1) - pplev(ig,2)) / gzlat(ig,1) / (zzlev(ig,2) - zzlev(ig,1))
123
124  ! Horizontal winds at the surface [m.s-1]
125  ws = sqrt(u(ig,1)*u(ig,1) + v(ig,1)*v(ig,1)) * (10. / zzlay(ig,1))**0.2
126
127  ! Source of turbulent kinetic energy at the surface [-] (Forget et al. 1999)
128  Cd = (karman / log(1. + zzlay(ig,1)/z0))**2
129
130  ! Saturation profile of CH4 [mol.mol-1] (Fray and Schmidt 2009) :
131  Tlake = tsurf(ig) - 7 ! Lakes are 2-7K less than surface.
132  qsatCH4 = (1.0e5 / pplev(ig,1)) * exp(1.051e1 - 1.110e3/Tlake - 4.341e3/Tlake**2 + 1.035e5/Tlake**3 - 7.910e5/Tlake**4)
133  ! CH4 : 0.70 * qsat because of dissolution in N2
134  qsatCH4 = 0.70 * qsatCH4
135
136  ! Flux at the surface [kg.m-2.s-1] :
137  flux = rhoair * Cd * ws
138 
139  ! Surface humidity :
140  qsatCH4 = humCH4 * qsatCH4
141
142  ! <North Pole>
143  if (REAL(latitude_deg(ig)) .ge. 70.0) then
144    flux = Flnp * flux
145  ! <South Pole>
146  else if (REAL(latitude_deg(ig)) .le. -70.0) then
147    flux = Flsp * flux
148  ! <Mid Latitudes>
149  else
150    flux = Flml * flux
151  endif
152
153  ! Empty tank ?
154  if (tankCH4(ig) .le. 1.e-30) then
155    flux = 0.0
156    tankCH4(ig) = 1.e-30
157  endif
158
159!------------------------------------
160! 2. IMPLICIT SCHEME
161!------------------------------------
162
163  ! Flux of CH4 at the surface [kg.m-2.s-1.mol.mol-1] :
164  fluxCH4 = flux * (qsatCH4 - pqCH4(ig,1))
165
166  ! Flux at the surface [kg.m-2.s-1] --> [s-1] :
167  flux = flux / rhoair / (zzlev(ig,2) - zzlev(ig,1))
168
169  ! New molar fraction of CH4 in the first layer [mol.mol-1] :
170  newpqCH4 = (pqCH4(ig,1) + flux * ptimestep * qsatCH4) / (1. + flux * ptimestep)
171
172  ! Trend of CH4's tank [m]
173  dtankCH4 = - (newpqCH4 - pqCH4(ig,1)) * rhoair * (zzlev(ig,2) - zzlev(ig,1)) * mmolCH4 / mmolair / rholCH4
174
175  ! New tank depth of CH4 [m]
176  if ((tankCH4(ig) + dtankCH4) .ge. 1.e-30) then
177    tankCH4(ig) = tankCH4(ig) + dtankCH4
178  else
179    newpqCH4 = tankCH4(ig) / (rhoair * (zzlev(ig,2) - zzlev(ig,1)) * mmolCH4 / mmolair / rholCH4) + pqCH4(ig,1)
180    tankCH4(ig) = 1.e-30
181  endif
182
183  ! Trend of CH4's molar fraction in the first layer [mol.mol-1.s-1]
184  pdqCH4(ig) = (newpqCH4 - pqCH4(ig,1)) / ptimestep
185
186
187!------------------------------------
188! 3. LATENTE HEAT
189!------------------------------------
190
191  ftm = (1. - Tlake / TcCH4)
192  if(ftm .le. 1.e-3) then
193    ftm = 1.e-3
194  endif
195 
196  ! Latente heat of CH4 vaporisation [J.kg-1.mol.mol-1]
197  LvCH4 = 8.314 * TcCH4 * (7.08 * ftm**0.354 + 10.95 * 1.1e-2 * ftm**0.456) / mmolCH4
198
199  ! Evaporation heating rate [K.s-1]
200  dtsurfevap(ig) = - (fluxCH4 * LvCH4) / Cplake
201
202enddo ! End of main loop on the horizontal grid
203
204return
205
206end subroutine evapCH4
Note: See TracBrowser for help on using the repository browser.