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

Last change on this file since 2841 was 2835, checked in by romain.vande, 3 years ago

Mars PEM:
Introduction of the possibility to follow an orbital forcing.
Introduction of new control parameters.
Cleaning of the PEM (removing unused files, add comments and new files)

A file named run_PEM.def can be added to the run.def. It contains the following variables:

_ evol_orbit_pem: Boolean. Do you want to follow an orbital forcing predefined (read in ob_ex_lsp.asc for example)? (default=false)
_ year_bp_ini: Integer. Number of year before present to start the pem run if evol_orbit_pem=.true. , default=0
_ Max_iter_pem: Integer. Maximal number of iteration if none of the stopping criterion is reached and if evol_orbit_pem=.false., default=99999999
_ dt_pem: Integer. Time step of the PEM in year, default=1
_ alpha_criterion: Real. Acceptance rate of sublimating ice surface change, default=0.2
_ soil_pem: Boolean. Do you want to run with subsurface physical processes in the PEM? default=.true.

RV

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