[3083] | 1 | subroutine evapCH4(ngrid,nlayer,ptimestep,pplev,zzlay,zzlev, & |
---|
| 2 | u,v,tsurf,pqCH4, & |
---|
| 3 | tankCH4,pdqCH4,dtsurfevap) |
---|
| 4 | |
---|
| 5 | use radcommon_h, only: gzlat ! Gravity of the planet [m.s-2] |
---|
| 6 | use planete_mod, only: z0 ! Surface roughness of the planet [m] |
---|
| 7 | use geometry_mod, only: latitude_deg ! Latitude grid of the planet [°] |
---|
| 8 | |
---|
| 9 | implicit 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 | ! |
---|
[3090] | 46 | ! Author |
---|
| 47 | ! ------ |
---|
[3083] | 48 | ! B. de Batz de Trenquelléon (10/2022) |
---|
| 49 | ! |
---|
| 50 | !==================================================================== |
---|
| 51 | |
---|
| 52 | |
---|
| 53 | !------------------------------------ |
---|
| 54 | ! 0. DECLARATIONS |
---|
| 55 | !------------------------------------ |
---|
| 56 | |
---|
| 57 | ! Inputs : |
---|
| 58 | !--------- |
---|
| 59 | integer, intent(in) :: ngrid |
---|
| 60 | integer, intent(in) :: nlayer |
---|
| 61 | real, intent(in) :: ptimestep |
---|
| 62 | real, intent(in) :: pplev(ngrid,nlayer+1) |
---|
| 63 | real, intent(in) :: zzlay(ngrid,nlayer) |
---|
| 64 | real, intent(in) :: zzlev(ngrid,nlayer+1) |
---|
| 65 | real, intent(in) :: u(ngrid,nlayer) |
---|
| 66 | real, intent(in) :: v(ngrid,nlayer) |
---|
| 67 | real, intent(in) :: tsurf(ngrid) |
---|
| 68 | real, intent(in) :: pqCH4(ngrid,nlayer) |
---|
| 69 | |
---|
| 70 | ! Outputs : |
---|
| 71 | !---------- |
---|
| 72 | real, intent(out) :: tankCH4(ngrid) |
---|
| 73 | real, intent(out) :: pdqCH4(ngrid) |
---|
| 74 | real, intent(out) :: dtsurfevap(ngrid) |
---|
| 75 | |
---|
| 76 | ! Parameters : |
---|
| 77 | !------------- |
---|
[3318] | 78 | real, parameter :: karman = 0.4 ! Karman constant [-] |
---|
| 79 | real, parameter :: humCH4 = 0.4 ! Imposed surface humidity for CH4 [-] |
---|
[3083] | 80 | |
---|
[3318] | 81 | real, parameter :: Flnp = 0.08 ! Fraction occupied by lakes (North Pole) |
---|
| 82 | real, parameter :: Flsp = 0.02 ! Fraction occupied by lakes (South Pole) |
---|
| 83 | real, parameter :: Flml = 0.02 ! Fraction not infiltrated into the ground (Mid latitudes) |
---|
[3083] | 84 | |
---|
| 85 | real, parameter :: mmolair = 28.e-3 ! Molar mass of air [kg.mol-1] |
---|
| 86 | real, parameter :: mmolCH4 = 16.e-3 ! Molar mass of CH4 [kg.mol-1] |
---|
| 87 | real, parameter :: rhoiCH4 = 425. ! Density of ice of CH4 [kg.m-3] |
---|
| 88 | |
---|
| 89 | real, parameter :: TcCH4 = 190.56 ! Critical point of CH4 [K] |
---|
| 90 | real, parameter :: Cplake = 2689992 ! Surface heat capacity for hydrocarbon lakes [J.m-2.K-1] (Tokano 2005) |
---|
| 91 | |
---|
| 92 | ! Local variables : |
---|
| 93 | !------------------ |
---|
| 94 | integer :: ig |
---|
| 95 | |
---|
| 96 | ! Initialisation : |
---|
[3318] | 97 | real*8 :: Tlake ! Lake temperature [K] |
---|
[3083] | 98 | real*8 :: rhoair ! Density of air [kg.m-3] |
---|
| 99 | real*8 :: ws ! Horizontal wind at the surface [m.s-1] |
---|
| 100 | real*8 :: Cd ! Turbulent term [-] |
---|
| 101 | real*8 :: qsatCH4 ! Saturation profile of CH4 [mol.mol-1] |
---|
| 102 | |
---|
| 103 | ! Flux : |
---|
| 104 | real*8 :: flux ! Surface flux [kg.m-2.s-1] |
---|
| 105 | real*8 :: fluxCH4 ! Surface flux fo CH4 [kg.m-2.s-1.mol.mol-1] |
---|
| 106 | |
---|
| 107 | ! Variations of CH4 : |
---|
| 108 | real*8 :: newpqCH4 ! New molar fraction of CH4 in the first layer [mol.mol-1] |
---|
| 109 | real*8 :: dtankCH4 ! Trend of CH4's tank [m] |
---|
| 110 | |
---|
| 111 | ! Latente heat : |
---|
| 112 | real*8 :: ftm, LvCH4 ! Variables for latente heat [-, J.kg-1] |
---|
| 113 | |
---|
| 114 | |
---|
| 115 | !------------------------------------ |
---|
| 116 | ! 1. INITIALISATION |
---|
| 117 | !------------------------------------ |
---|
| 118 | |
---|
| 119 | do 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) : |
---|
[3318] | 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) |
---|
[3090] | 133 | ! CH4 : 0.80 * qsat because of dissolution in N2 |
---|
| 134 | qsatCH4 = 0.80 * qsatCH4 |
---|
[3083] | 135 | |
---|
| 136 | ! Flux at the surface [kg.m-2.s-1] : |
---|
| 137 | flux = rhoair * Cd * ws |
---|
[3318] | 138 | |
---|
| 139 | ! Surface humidity : |
---|
| 140 | qsatCH4 = humCH4 * qsatCH4 |
---|
[3083] | 141 | |
---|
| 142 | ! <North Pole> |
---|
[3318] | 143 | if (REAL(latitude_deg(ig)) .ge. 70.0) then |
---|
[3083] | 144 | flux = Flnp * flux |
---|
| 145 | ! <South Pole> |
---|
[3318] | 146 | else if (REAL(latitude_deg(ig)) .le. -70.0) then |
---|
[3083] | 147 | flux = Flsp * flux |
---|
| 148 | ! <Mid Latitudes> |
---|
| 149 | else |
---|
[3318] | 150 | flux = Flml * flux |
---|
[3083] | 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 | |
---|
[3318] | 163 | ! Flux of CH4 at the surface [kg.m-2.s-1.mol.mol-1] : |
---|
| 164 | fluxCH4 = flux * (qsatCH4 - pqCH4(ig,1)) |
---|
| 165 | |
---|
[3083] | 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 / rhoiCH4 |
---|
| 174 | |
---|
| 175 | ! New tank depth of CH4 [m] |
---|
[3318] | 176 | if ((tankCH4(ig) + dtankCH4) .ge. 1.e-30) then |
---|
[3083] | 177 | tankCH4(ig) = tankCH4(ig) + dtankCH4 |
---|
| 178 | else |
---|
[3318] | 179 | !write(*,*) 'Evaporation CH4 : Empty lakes...' |
---|
[3083] | 180 | newpqCH4 = tankCH4(ig) / (rhoair * (zzlev(ig,2) - zzlev(ig,1)) * mmolCH4 / mmolair / rhoiCH4) + pqCH4(ig,1) |
---|
| 181 | tankCH4(ig) = 1.e-30 |
---|
| 182 | endif |
---|
| 183 | |
---|
| 184 | ! Trend of CH4's molar fraction in the first layer [mol.mol-1.s-1] |
---|
| 185 | pdqCH4(ig) = (newpqCH4 - pqCH4(ig,1)) / ptimestep |
---|
| 186 | |
---|
| 187 | |
---|
| 188 | !------------------------------------ |
---|
| 189 | ! 3. LATENTE HEAT |
---|
| 190 | !------------------------------------ |
---|
| 191 | |
---|
[3318] | 192 | ftm = (1. - Tlake / TcCH4) |
---|
[3083] | 193 | if(ftm .le. 1.e-3) then |
---|
| 194 | ftm = 1.e-3 |
---|
| 195 | endif |
---|
| 196 | |
---|
| 197 | ! Latente heat of CH4 vaporisation [J.kg-1.mol.mol-1] |
---|
[3318] | 198 | LvCH4 = 8.314 * TcCH4 * (7.08 * ftm**0.354 + 10.95 * 1.1e-2 * ftm**0.456) / mmolCH4 |
---|
[3083] | 199 | |
---|
| 200 | ! Evaporation heating rate [K.s-1] |
---|
| 201 | dtsurfevap(ig) = - (fluxCH4 * LvCH4) / Cplake |
---|
| 202 | |
---|
| 203 | enddo ! End of main loop on the horizontal grid |
---|
| 204 | |
---|
| 205 | return |
---|
| 206 | |
---|
| 207 | end subroutine evapCH4 |
---|