source: trunk/LMDZ.COMMON/libf/evolution/adsorption_mod.F90 @ 2814

Last change on this file since 2814 was 2794, checked in by llange, 3 years ago

MARS PEM:

  • Add a PEMETAT0 that read "startfi_pem.nc"
  • Add the soil in the model: soil temperature, thermal properties, ice table
  • Add a routine that compute CO2 + H2O adsorption
  • Minor corrections in PEM.F90

LL

File size: 11.2 KB
Line 
1  module adsorption_mod 
2  implicit none
3  contains
4
5
6
7  subroutine regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen, ps,q_co2,q_h2o,tsoil_PEM,TI_PEM, theta_h2o_adsorbded,m_h2o_adsorbed)
8
9 use vertical_layers_mod, ONLY: ap,bp
10 use comsoil_h_PEM, only: n_1km
11
12 implicit none
13
14
15 INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen ! size dimension
16 REAL,INTENT(IN) :: ps(ngrid,timelen)                   ! surface pressure (Pa)         
17 REAL,INTENT(IN) :: q_co2(ngrid,timelen)                ! Mass mixing ratio of co2 in the first layer (kg/kg)
18 REAL,INTENT(IN) :: q_h2o(ngrid,timelen)                ! Mass mixing ratio of H2o in the first layer (kg/kg)
19 REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope)      ! Soil Thermal inertia (J/K/^2/s^1/2)
20 REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope)   ! Soil temperature (K)
21
22
23! output
24 REAL,INTENT(OUT) ::  m_h2o_adsorbed(ngrid,nsoil_PEM,nslope)     ! Density of h2o adsorbed (kg/m^3)
25 REAL,INTENT(OUT) :: theta_h2o_adsorbded(ngrid,nsoil_PEM,nslope) ! Fraction of the pores occupied by H2O molecules
26
27
28! constant
29 REAL :: Ko =  1.57e-8            ! Jackosky et al. 1997
30 REAL :: e = 2573.9               ! Jackosky et al. 1997
31 REAL :: mu = 0.48                ! Jackosky et al. 1997
32 REAL ::  inertie_thresold = 800. ! TI > 800 means cementation
33 real :: m_h2o = 18.01528E-3      ! Molecular weight of h2o (kg/mol)
34 real :: m_co2 = 44.01E-3         ! Molecular weight of co2 (kg/mol)
35 real :: m_noco2 = 33.37E-3       ! Molecular weight of non co2 (kg/mol)
36 REAL ::  rho_regolith = 2000.    ! density of the reoglith, Buhler & Piqueux 2021
37 real :: alpha_clapeyron = -6143.7! eq. (2) in Murphy & Koop 2005
38 real :: beta_clapeyron = 29.9074 ! eq. (2) in Murphy & Koop 2005
39 real :: mi = 2.84e-7             ! Mass of h2o per m^2 absorbed Jackosky et al. 1997
40 real :: as = 18.9e3              ! Specific area, Buhler & Piqueux 2021
41
42
43
44
45
46! local variable
47 real :: A,B                                                   ! Used to compute the mean mass above the surface
48 real :: K                                                     ! Used to compute theta
49 real :: p_sat                                                 ! saturated vapor pressure of ice
50 integer ig,iloop, islope,isoil,it                             ! for loops
51 real,allocatable :: mass_mean(:,:)                            ! mean mass above the surface
52 real,allocatable :: zplev_mean(:,:)                           ! pressure above the surface
53 real,allocatable :: pvapor(:,:)                               ! partial pressure above the surface
54 real, allocatable :: pvapor_avg(:)                            ! yearly average vapor pressure above the surface
55
56
57
58! 0. Some initializations
59
60
61   
62     allocate(mass_mean(ngrid,timelen))
63     allocate(zplev_mean(ngrid,timelen))
64     allocate(pvapor(ngrid,timelen))
65     allocate(pvapor_avg(ngrid))
66
67      m_h2o_adsorbed(:,:,:) = 0.
68      theta_h2o_adsorbded(:,:,:) = 0.
69      A =(1/m_co2 - 1/m_noco2)
70      B=1/m_noco2
71! 1. Compute rho surface yearly averaged
72
73!   1.1 Compute the partial pressure of vapor
74!a. the molecular mass into the column
75     do ig = 1,ngrid
76       mass_mean(ig,:) = 1/(A*q_co2(ig,:) +B)
77     enddo
78
79
80! b. pressure level
81     do it = 1,timelen
82       do ig = 1,ngrid
83         zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it)
84       enddo
85     enddo
86
87! c. Vapor pressure
88     pvapor(:,:) = mass_mean(:,:)/m_h2o*q_h2o(:,:)*zplev_mean(:,:)
89     pvapor_avg(:) = sum(pvapor(:,:),2)/timelen
90
91     deallocate(pvapor)
92     deallocate(zplev_mean)
93     deallocate(mass_mean)
94
95
96
97
98
99! 2. we compute the mass of co2 adsorded in each layer of the meshes 
100
101 do ig = 1,ngrid
102  do islope = 1,nslope
103    do iloop = 1,n_1km
104     if(TI_PEM(ig,iloop,islope).lt.inertie_thresold)  then
105       K = Ko*exp(e/tsoil_PEM(ig,iloop,islope))
106       theta_h2o_adsorbded(ig,iloop,islope) = (K*pvapor_avg(ig)/(1+K*pvapor_avg(ig)))**mu
107       m_h2o_adsorbed(ig,iloop,islope) = as*theta_h2o_adsorbded(ig,iloop,islope)*mi*rho_regolith
108     else
109        p_sat =exp(alpha_clapeyron/tsoil_PEM(ig,iloop,islope) +beta_clapeyron) ! we assume fixed temperature in the ice ... not really:q good but ...
110        theta_h2o_adsorbded(ig,iloop,islope) = (K*p_sat/(1+K*p_sat))**mu
111        m_h2o_adsorbed(ig,iloop,islope) =as*theta_h2o_adsorbded(ig,iloop,islope)*mi*rho_regolith
112      endif
113 
114   enddo
115  enddo
116 enddo
117
118 RETURN
119  end subroutine
120
121
122  SUBROUTINE regolith_co2adsorption(ngrid,nslope,nsoil_PEM,timelen,ps,tsoil_PEM,TI_PEM,tend_h2oglaciers,tend_co2glaciers,co2ice,waterice,q_co2,q_h2o,m_co2_completesoil,delta_mreg)
123      use comsoil_h_PEM, only: layer_PEM, mlayer_PEM,n_1km
124      USE comcstfi_h, only: r, cpp, mugaz, g, rcp, pi, rad
125      use comslope_mod, only : subslope_dist,def_slope_mean
126      use vertical_layers_mod, ONLY: ap,bp
127
128      IMPLICIT NONE
129! Input: 
130 INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen             ! size dimension
131 REAL,INTENT(IN) :: ps(ngrid,timelen)                               ! Average surface pressure [Pa]
132 REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope)               ! Mean Soil Temperature [K]
133 REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope)                  ! Mean Thermal Inertia [USI]
134 REAL,INTENT(IN) :: tend_h2oglaciers(ngrid,nslope),tend_co2glaciers(ngrid,nslope) !tendancies on the glaciers ()
135 REAL,INTENT(IN) :: q_co2(ngrid,timelen),q_h2o(ngrid,timelen)       ! Mass mixing ratio of co2 and h2o in the first layer (kg/kg)
136 REAL,INTENT(IN) :: waterice(ngrid,nslope)                          ! water ice at the surface [kg/m^2]
137 REAL,INTENT(IN) :: co2ice(ngrid,nslope)                            ! co2 ice at the surface [kg/m^2]
138
139! Outputs:
140 REAL,INTENT(INOUT) :: m_co2_completesoil(ngrid,nsoil_PEM,nslope)   ! Density of co2 adsorbed (kg/m^3)
141 REAL,INTENT(INOUT) :: delta_mreg(ngrid)                            ! Difference density of co2 adsorbed (kg/m^3)
142
143
144 
145! Constants:
146
147 REAL :: alpha = 7.512e-6 ! Zent & Quinn 1995
148 REAL :: beta =  -1541.5  ! Zent & Quinn 1995
149 REAL ::  inertie_thresold = 800. ! TI > 800 means cementation
150 REAL ::  rho_regolith = 2000. ! density of the reoglith, buhler & piqueux 2021
151 real :: m_co2 = 44.01E-3      ! Molecular weight of co2 (kg/mol)
152 real :: m_noco2 = 33.37E-3    ! Molecular weight of h2o (kg/mol)
153 real :: m_theta = 4.27e-7     ! Mass of co2 per m^2 absorbed
154 real :: as = 18.9e3              ! Specific area, Buhler & Piqueux 2021
155
156! Local         
157 real :: A,B                                             ! Used to compute the mean mass above the surface
158 INTEGER :: ig,islope,iloop,it                           ! for loops
159 REAL :: dm_co2_regolith_slope(ngrid,nsoil_PEM,nslope)   ! elementary mass adsorded per mesh per slope
160 INTEGER :: ispermanent_co2glaciers(ngrid,nslope)        ! Check if the glacier is permanent
161 INTEGER :: ispermanent_h2oglaciers(ngrid,nslope)        ! Check if the glacier is permanent
162 REAL :: deltam_reg_complete(ngrid,n_1km,nslope)         ! Difference in the mass per slope and soil layer (kg/m^3)
163 REAL :: deltam_reg_slope(ngrid,nslope)                  ! Difference in the mass per slope  (kg/m^3)
164 REAL :: m_h2o_adsorbed(ngrid,nsoil_PEM,nslope)          ! Density of CO2 adsorbed (kg/m^3)
165 REAL :: theta_h2o_adsorbed(ngrid,nsoil_PEM,nslope)     ! Fraction of the pores occupied by H2O molecules
166!timelen array are allocated because heavy ...
167 real,allocatable :: mass_mean(:,:)                            ! mean mass above the surface
168 real,allocatable :: zplev_mean(:,:)                           ! pressure above the surface
169 real,allocatable :: pco2(:,:)                                  ! partial pressure above the surface
170 real, allocatable :: pco2_avg(:)                              ! yearly averaged
171
172
173
174! 0. Some initializations
175
176
177
178     allocate(mass_mean(ngrid,timelen))
179     allocate(zplev_mean(ngrid,timelen))
180     allocate(pco2(ngrid,timelen))
181     allocate(pco2_avg(ngrid))
182 
183      m_h2o_adsorbed(:,:,:) = 0.
184      A =(1/m_co2 - 1/m_noco2)
185      B=1/m_noco2
186
187     dm_co2_regolith_slope(:,:,:) = 0
188     delta_mreg(:) = 0.
189
190
191
192!0.1 Look at perenial ice
193  do ig = 1,ngrid
194    do islope = 1,nslope
195     if((abs(tend_h2oglaciers(ig,islope)).gt.1e-5).and.(abs(waterice(ig,islope)).gt.0)) then
196        ispermanent_h2oglaciers(ig,islope) = 1
197     else
198        ispermanent_h2oglaciers(ig,islope) = 0
199     endif
200
201     if((abs(tend_co2glaciers(ig,islope)).gt.1e-5).and.(abs(co2ice(ig,islope)).gt.0)) then
202        ispermanent_co2glaciers(ig,islope) = 1
203     else
204        ispermanent_co2glaciers(ig,islope) = 0
205     endif
206    enddo
207   enddo
208
209
210!   0.2  Compute the partial pressure of CO2
211!a. the molecular mass into the column
212     do ig = 1,ngrid
213       mass_mean(ig,:) = 1/(A*q_co2(ig,:) +B)
214     enddo
215
216
217! b. pressure level
218     do it = 1,timelen
219       do ig = 1,ngrid
220         zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it)
221       enddo
222     enddo
223
224! c. Vapor pressure
225     pco2(:,:) = mass_mean(:,:)/m_co2*q_co2(:,:)*zplev_mean(:,:)
226     pco2_avg(:) = sum(pco2(:,:),2)/timelen
227
228 
229     deallocate(zplev_mean)
230     deallocate(mass_mean)
231     deallocate(pco2)
232
233
234! 1. Compute the fraction of the pores occupied by H2O
235 call regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen, ps,q_co2,q_h2o,tsoil_PEM,TI_PEM,theta_h2o_adsorbed, m_h2o_adsorbed)
236
237
238
239! 2.  we compute the mass of co2 adsorded in each layer of the meshes 
240
241 do ig = 1,ngrid
242  do islope = 1,nslope
243    do iloop = 1,n_1km
244     if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then
245     dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1-theta_h2o_adsorbed(ig,iloop,islope))*alpha*pco2_avg(ig)/ &
246                                             (alpha*pco2_avg(ig)+sqrt(tsoil_PEM(ig,iloop,islope))*exp(beta/tsoil_PEM(ig,iloop,islope)))
247     else
248        if(abs(m_co2_completesoil(ig,iloop,islope)).lt.1-10) then !!! we are at first call
249          dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1-theta_h2o_adsorbed(ig,iloop,islope))*alpha*pco2_avg(ig) &
250                                                 /(alpha*pco2_avg(ig)+sqrt(tsoil_PEM(ig,iloop,islope))*exp(beta/tsoil_PEM(ig,iloop,islope)))
251        else ! no change: permanent ice stick the atoms of CO2
252          dm_co2_regolith_slope(ig,iloop,islope) = m_co2_completesoil(ig,iloop,islope)
253        endif
254     endif
255  enddo
256 enddo
257enddo
258
259
260! 2. Check the exchange between the atmosphere and the regolith
261
262  do ig = 1,ngrid
263   delta_mreg(ig) = 0.
264   do islope = 1,nslope
265    deltam_reg_slope(ig,islope) = 0.
266    do iloop = 1,n_1km
267       if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then
268             deltam_reg_complete(ig,iloop,islope) = (dm_co2_regolith_slope(ig,iloop,islope) - m_co2_completesoil(ig,iloop,islope)) &
269                                                   *(layer_PEM(iloop+1) - layer_PEM(iloop))
270       else ! NO EXCHANGE AS ICE BLOCK THE DYNAMIC!
271             deltam_reg_complete(ig,iloop,islope) = 0.
272       endif
273       deltam_reg_slope(ig,islope) = deltam_reg_slope(ig,islope) + deltam_reg_complete(ig,iloop,islope)
274    enddo
275   delta_mreg(ig) = delta_mreg(ig) + deltam_reg_slope(ig,islope)*subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)
276   enddo
277  enddo
278   m_co2_completesoil(:,:,:) = dm_co2_regolith_slope(:,:,:)
279
280
281!=======================================================================
282      RETURN
283      END
284
285
286   end module
Note: See TracBrowser for help on using the repository browser.