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

Last change on this file since 2937 was 2895, checked in by llange, 3 years ago

PEM
Soil temperature initialisation has been updated
Conf_PEM improved by adding some options to the users (thermal regolith depend on the pressure, depth of the subsurface layers, etc.)
Minor edits then (+ svn update with RV had some issues, so there are some "artefact changes" ...)
LL

File size: 18.6 KB
Line 
1  module adsorption_mod 
2  implicit none
3  LOGICAL adsorption_pem ! True by default, to compute adsorption/desorption. Read in  pem.def
4  real, save, allocatable :: co2_adsorbded_phys(:,:,:)  ! co2 that is in the regolith [kg/m^2]
5  real, save, allocatable :: h2o_adsorbded_phys(:,:,:)  ! h2o that is in the regolith [kg/m^2]
6  contains
7
8!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
9!!!
10!!! Purpose: Compute CO2 and H2O adsorption, following the methods from Zent & Quinn 1995, Jackosky et al., 1997
11!!!
12!!! Author: LL, 01/2023
13!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15
16
17  subroutine ini_adsorption_h_PEM(ngrid,nslope,nsoilmx_PEM)
18
19  implicit none
20  integer,intent(in) :: ngrid ! number of atmospheric columns
21  integer,intent(in) :: nslope ! number of slope within a mesh
22  integer,intent(in) :: nsoilmx_PEM ! number of soil layer in the PEM
23    allocate(co2_adsorbded_phys(ngrid,nsoilmx_PEM,nslope))
24    allocate(h2o_adsorbded_phys(ngrid,nsoilmx_PEM,nslope))
25  end subroutine ini_adsorption_h_PEM
26
27  subroutine end_adsorption_h_PEM
28
29  implicit none
30    if (allocated(co2_adsorbded_phys)) deallocate(co2_adsorbded_phys)
31    if (allocated(h2o_adsorbded_phys)) deallocate(h2o_adsorbded_phys)
32  end subroutine end_adsorption_h_PEM
33
34
35
36
37
38
39  subroutine regolith_adsorption(ngrid,nslope,nsoil_PEM,timelen,tend_h2oglaciers,tend_co2glaciers,waterice,co2ice,tsoil_PEM,TI_PEM,ps,q_co2,q_h2o, &
40                                m_h2o_completesoil,delta_mh2oreg, m_co2_completesoil,delta_mco2reg)
41#ifndef CPP_STD
42! inputs
43 INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen ! size dimension: physics x subslope x soil x timeseries
44 REAL,INTENT(IN) :: tend_h2oglaciers(ngrid,nslope),tend_co2glaciers(ngrid,nslope) !tendancies on the glaciers [1]
45 REAL,INTENT(IN) :: waterice(ngrid,nslope)              ! water ice at the surface [kg/m^2]
46 REAL,INTENT(IN) :: co2ice(ngrid,nslope)                ! co2 ice at the surface [kg/m^2]
47 REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope)      ! Soil Thermal inertia (J/K/^2/s^1/2)
48 REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope)   ! Soil temperature (K)
49 REAL,INTENT(IN) :: ps(ngrid,timelen)                   ! Average surface pressure [Pa]
50 REAL,INTENT(IN) :: q_co2(ngrid,timelen)                ! Mass mixing ratio of co2 in the first layer (kg/kg)
51 REAL,INTENT(IN) :: q_h2o(ngrid,timelen)                ! Mass mixing ratio of H2o in the first layer (kg/kg)
52
53! outputs
54 REAL,INTENT(INOUT) :: m_h2o_completesoil(ngrid,nsoil_PEM,nslope) ! Density of h2o adsorbed (kg/m^3)(ngrid,nsoil_PEM,nslope)     
55 REAL,INTENT(OUT) :: delta_mh2oreg(ngrid)                         ! Difference density of h2o adsorbed (kg/m^3)
56
57 REAL,INTENT(INOUT) :: m_co2_completesoil(ngrid,nsoil_PEM,nslope)   ! Density of co2 adsorbed (kg/m^3)
58 REAL,INTENT(OUT) :: delta_mco2reg(ngrid)                            ! Difference density of co2 adsorbed (kg/m^3)
59 
60! local variables
61 REAL :: theta_h2o_adsorbded(ngrid,nsoil_PEM,nslope) ! Fraction of the pores occupied by H2O molecules
62! -------------
63
64! Compute H2O adsorption, then CO2 adsorption
65
66  call regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen,tend_h2oglaciers,tend_co2glaciers,waterice,co2ice,ps,q_co2,q_h2o,tsoil_PEM,TI_PEM, &
67                                   theta_h2o_adsorbded,m_h2o_completesoil,delta_mh2oreg)
68
69
70  call regolith_co2adsorption(ngrid,nslope,nsoil_PEM,timelen,tend_h2oglaciers,tend_co2glaciers,waterice,co2ice,ps,q_co2,q_h2o, &
71                                                          tsoil_PEM,TI_PEM,m_co2_completesoil,delta_mco2reg)
72
73#endif
74   RETURN
75  end subroutine
76
77!------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
78
79  subroutine regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen,tend_h2oglaciers,tend_co2glaciers,waterice,co2ice,ps,q_co2,q_h2o,tsoil_PEM,TI_PEM, &
80                                   theta_h2o_adsorbded,m_h2o_completesoil,delta_mreg)
81
82!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
83!!!
84!!! Purpose: Compute H2O adsorption, following the methods from Jackosky et al., 1997
85!!!
86!!! Author: LL, 01/2023
87!!!
88!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
89
90#ifndef CPP_STD
91      use comsoil_h_PEM, only: layer_PEM, mlayer_PEM,index_breccia,index_bedrock
92      USE comcstfi_h, only:  pi
93      use comslope_mod, only : subslope_dist,def_slope_mean
94      use vertical_layers_mod, ONLY: ap,bp
95#endif
96
97 implicit none
98! inputs
99 INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen ! size dimension
100 REAL,INTENT(IN) :: tend_h2oglaciers(ngrid,nslope),tend_co2glaciers(ngrid,nslope) !tendancies on the glaciers ()
101 REAL,INTENT(IN) :: waterice(ngrid,nslope)              ! water ice at the surface [kg/m^2]
102 REAL,INTENT(IN) :: co2ice(ngrid,nslope)                ! co2 ice at the surface [kg/m^2]
103 REAL,INTENT(IN) :: ps(ngrid,timelen)                   ! surface pressure (Pa)         
104 REAL,INTENT(IN) :: q_co2(ngrid,timelen)                ! Mass mixing ratio of co2 in the first layer (kg/kg)
105 REAL,INTENT(IN) :: q_h2o(ngrid,timelen)                ! Mass mixing ratio of H2o in the first layer (kg/kg)
106 REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope)      ! Soil Thermal inertia (J/K/^2/s^1/2)
107 REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope)   ! Soil temperature (K)
108
109! outputs
110 REAL,INTENT(INOUT) :: m_h2o_completesoil(ngrid,nsoil_PEM,nslope) ! Density of h2o adsorbed (kg/m^3)(ngrid,nsoil_PEM,nslope)     
111 REAL,INTENT(OUT) :: theta_h2o_adsorbded(ngrid,nsoil_PEM,nslope) ! Fraction of the pores occupied by H2O molecules
112 REAL,INTENT(OUT) :: delta_mreg(ngrid)                         ! Difference density of h2o adsorbed (kg/m^3)
113
114! constants
115 REAL :: Ko =  1.57e-8            ! Jackosky et al. 1997
116 REAL :: e = 2573.9               ! Jackosky et al. 1997
117 REAL :: mu = 0.48                ! Jackosky et al. 1997
118 real :: m_theta = 2.84e-7        ! Mass of h2o per m^2 absorbed Jackosky et al. 1997
119! real :: as = 18.9e3             ! Specific area, Buhler & Piqueux 2021
120 real :: as = 9.48e4              ! Specific area, Zent
121 real :: as_breccia = 1.7e4       ! Specific area of basalt, Zent
122 real ::  inertie_thresold = 800. ! TI > 800 means cementation
123 real :: m_h2o = 18.01528E-3      ! Molecular weight of h2o (kg/mol)
124 real :: m_co2 = 44.01E-3         ! Molecular weight of co2 (kg/mol)
125 real :: m_noco2 = 33.37E-3       ! Molecular weight of non co2 (kg/mol)
126 real ::  rho_regolith = 2000.    ! density of the regolith (Zent 1995, Buhler & piqueux 2021)
127 real :: alpha_clapeyron = -6143.7! eq. (2) in Murphy & Koop 2005
128 real :: beta_clapeyron = 28.9074 ! eq. (2) in Murphy & Koop 2005
129! local variables
130#ifndef CPP_STD
131 REAL :: deltam_reg_complete(ngrid,index_breccia,nslope)         ! Difference in the mass per slope and soil layer (kg/m^3)
132#endif
133 real :: K                        ! Used to compute theta
134 integer ig,iloop, islope,isoil,it   ! for loops
135 INTEGER :: ispermanent_co2glaciers(ngrid,nslope)        ! Check if the co2 glacier is permanent
136 INTEGER :: ispermanent_h2oglaciers(ngrid,nslope)        ! Check if the h2o glacier is permanent
137 REAL :: deltam_reg_slope(ngrid,nslope)  ! Difference density of h2o adsorbed per slope (kg/m^3)
138 REAL :: dm_h2o_regolith_slope(ngrid,nsoil_PEM,nslope)   ! elementary h2o mass adsorded per mesh per slope
139 real :: A,B                                                   ! Used to compute the mean mass above the surface
140 real :: p_sat                                                 ! saturated vapor pressure of ice
141 real,allocatable :: mass_mean(:,:)                            ! mean mass above the surface
142 real,allocatable :: zplev_mean(:,:)                           ! pressure above the surface
143 real,allocatable :: pvapor(:,:)                               ! partial pressure above the surface
144 real, allocatable :: pvapor_avg(:)                            ! yearly
145
146! 0. Some initializations
147#ifndef CPP_STD
148
149     allocate(mass_mean(ngrid,timelen))
150     allocate(zplev_mean(ngrid,timelen))
151     allocate(pvapor(ngrid,timelen))
152     allocate(pvapor_avg(ngrid))
153     A =(1/m_co2 - 1/m_noco2)
154     B=1/m_noco2
155     theta_h2o_adsorbded(:,:,:) = 0.
156     dm_h2o_regolith_slope(:,:,:) = 0.
157
158!0.1 Look at perenial ice
159  do ig = 1,ngrid
160    do islope = 1,nslope
161     if((abs(tend_h2oglaciers(ig,islope)).gt.1e-5).and.(abs(waterice(ig,islope)).gt.0)) then
162        ispermanent_h2oglaciers(ig,islope) = 1
163     else
164        ispermanent_h2oglaciers(ig,islope) = 0
165     endif
166
167     if((abs(tend_co2glaciers(ig,islope)).gt.1e-5).and.(abs(co2ice(ig,islope)).gt.0)) then
168        ispermanent_co2glaciers(ig,islope) = 1
169     else
170        ispermanent_co2glaciers(ig,islope) = 0
171     endif
172    enddo
173   enddo
174
175!   0.2 Compute the partial pressure of vapor
176!a. the molecular mass into the column
177     do ig = 1,ngrid
178       mass_mean(ig,:) = 1/(A*q_co2(ig,:) +B)
179     enddo
180
181
182! b. pressure level
183     do it = 1,timelen
184       do ig = 1,ngrid
185         zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it)
186       enddo
187     enddo
188! c. Vapor pressure
189     pvapor(:,:) = mass_mean(:,:)/m_h2o*q_h2o(:,:)*zplev_mean(:,:)
190     pvapor_avg(:) = sum(pvapor(:,:),2)/timelen
191     deallocate(pvapor)
192     deallocate(zplev_mean)
193     deallocate(mass_mean)
194
195! 1. we compute the mass of H2O adsorded in each layer of the meshes 
196
197 do ig = 1,ngrid
198  do islope = 1,nslope
199    do iloop = 1,index_breccia
200      K = Ko*exp(e/tsoil_PEM(ig,iloop,islope))
201      if(TI_PEM(ig,iloop,islope).lt.inertie_thresold)  then
202        theta_h2o_adsorbded(ig,iloop,islope) = (K*pvapor_avg(ig)/(1+K*pvapor_avg(ig)))**mu
203      else
204        p_sat =exp(alpha_clapeyron/tsoil_PEM(ig,iloop,islope) +beta_clapeyron) ! we assume fixed temperature in the ice ... not really:q good but ...
205        theta_h2o_adsorbded(ig,iloop,islope) = (K*p_sat/(1+K*p_sat))**mu
206      endif
207        dm_h2o_regolith_slope(ig,iloop,islope) = as*theta_h2o_adsorbded(ig,iloop,islope)*m_theta*rho_regolith
208   enddo
209  enddo
210 enddo
211
212! 2. Check the exchange between the atmosphere and the regolith
213
214  do ig = 1,ngrid
215   delta_mreg(ig) = 0.
216   do islope = 1,nslope
217    deltam_reg_slope(ig,islope) = 0.
218    do iloop = 1,index_breccia
219       if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then
220             deltam_reg_complete(ig,iloop,islope) = (dm_h2o_regolith_slope(ig,iloop,islope) - m_h2o_completesoil(ig,iloop,islope)) &
221                                                   *(layer_PEM(iloop+1) - layer_PEM(iloop))
222       else ! NO EXCHANGE AS ICE BLOCK THE DYNAMIC!
223             deltam_reg_complete(ig,iloop,islope) = 0.
224       endif
225       deltam_reg_slope(ig,islope) = deltam_reg_slope(ig,islope) + deltam_reg_complete(ig,iloop,islope)
226    enddo
227   delta_mreg(ig) = delta_mreg(ig) + deltam_reg_slope(ig,islope)*subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)
228   enddo
229  enddo
230   m_h2o_completesoil(:,:,:) = dm_h2o_regolith_slope(:,:,:)
231
232
233 RETURN
234#endif
235  end subroutine
236
237!------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
238!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
239!!!
240!!! Purpose: Compute CO2  following the methods from Zent & Quinn 1995
241!!!
242!!! Author: LL, 01/2023
243!!!
244!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
245
246  SUBROUTINE regolith_co2adsorption(ngrid,nslope,nsoil_PEM,timelen,tend_h2oglaciers,tend_co2glaciers,waterice,co2ice,ps,q_co2,q_h2o,&
247                                   tsoil_PEM,TI_PEM,m_co2_completesoil,delta_mreg)
248#ifndef CPP_STD
249      use comsoil_h_PEM, only: layer_PEM, mlayer_PEM,index_breccia,index_bedrock,index_breccia
250      USE comcstfi_h, only: pi
251      use comslope_mod, only : subslope_dist,def_slope_mean
252      use vertical_layers_mod, ONLY: ap,bp
253#endif
254
255      IMPLICIT NONE
256! Inputs: 
257 INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen             ! size dimension
258 REAL,INTENT(IN) :: ps(ngrid,timelen)                               ! Average surface pressure [Pa]
259 REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope)               ! Mean Soil Temperature [K]
260 REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope)                  ! Mean Thermal Inertia [USI]
261 REAL,INTENT(IN) :: tend_h2oglaciers(ngrid,nslope),tend_co2glaciers(ngrid,nslope) !tendancies on the glaciers ()
262 REAL,INTENT(IN) :: q_co2(ngrid,timelen),q_h2o(ngrid,timelen)       ! Mass mixing ratio of co2 and h2o in the first layer (kg/kg)
263 REAL,INTENT(IN) :: waterice(ngrid,nslope)                          ! water ice at the surface [kg/m^2]
264 REAL,INTENT(IN) :: co2ice(ngrid,nslope)                            ! co2 ice at the surface [kg/m^2]
265
266! Outputs:
267 REAL,INTENT(INOUT) :: m_co2_completesoil(ngrid,nsoil_PEM,nslope)   ! Density of co2 adsorbed (kg/m^3)
268 REAL,INTENT(OUT) :: delta_mreg(ngrid)                              ! Difference density of co2 adsorbed (kg/m^3)
269 
270! Constants:
271
272 REAL :: alpha = 7.512e-6 ! Zent & Quinn 1995
273 REAL :: beta =  -1541.5  ! Zent & Quinn 1995
274 REAL ::  inertie_thresold = 800. ! TI > 800 means cementation
275 REAL ::  rho_regolith = 2000. ! density of the regolith, buhler & piqueux 2021
276 real :: m_co2 = 44.01E-3      ! Molecular weight of co2 (kg/mol)
277 real :: m_noco2 = 33.37E-3    ! Molecular weight of h2o (kg/mol)
278 real :: m_theta = 4.27e-7     ! Mass of co2 per m^2 absorbed
279! real :: as = 18.9e3              ! Specific area, Buhler & Piqueux 2021
280 real :: as = 9.48e4           ! same as previous but from zent
281 real :: as_breccia = 1.7e4       ! Specific area of basalt, Zent 1997
282! Local         
283 real :: A,B                                             ! Used to compute the mean mass above the surface
284 INTEGER :: ig,islope,iloop,it                           ! for loops
285 REAL :: dm_co2_regolith_slope(ngrid,nsoil_PEM,nslope)   ! elementary mass adsorded per mesh per slope
286 INTEGER :: ispermanent_co2glaciers(ngrid,nslope)        ! Check if the co2 glacier is permanent
287 INTEGER :: ispermanent_h2oglaciers(ngrid,nslope)        ! Check if the h2o glacier is permanent
288#ifndef CPP_STD
289 REAL :: deltam_reg_complete(ngrid,index_breccia,nslope)         ! Difference in the mass per slope and soil layer (kg/m^3)
290#endif
291 REAL :: deltam_reg_slope(ngrid,nslope)                  ! Difference in the mass per slope  (kg/m^3)
292 REAL :: m_h2o_adsorbed(ngrid,nsoil_PEM,nslope)          ! Density of CO2 adsorbed (kg/m^3)
293 REAL :: theta_h2o_adsorbed(ngrid,nsoil_PEM,nslope)     ! Fraction of the pores occupied by H2O molecules
294 REAL :: delta_mh2o(ngrid)                              ! Difference density of h2o adsorbed (kg/m^3)
295!timelen array are allocated because heavy ...
296 real,allocatable :: mass_mean(:,:)                            ! mean mass above the surface
297 real,allocatable :: zplev_mean(:,:)                           ! pressure above the surface
298 real,allocatable :: pco2(:,:)                                  ! partial pressure above the surface
299 real, allocatable :: pco2_avg(:)                              ! yearly averaged
300
301! 0. Some initializations
302
303     allocate(mass_mean(ngrid,timelen))
304     allocate(zplev_mean(ngrid,timelen))
305     allocate(pco2(ngrid,timelen))
306     allocate(pco2_avg(ngrid))
307
308#ifndef CPP_STD
309 
310      m_h2o_adsorbed(:,:,:) = 0.
311      A =(1/m_co2 - 1/m_noco2)
312      B=1/m_noco2
313
314     dm_co2_regolith_slope(:,:,:) = 0
315     delta_mreg(:) = 0.
316
317!0.1 Look at perenial ice
318  do ig = 1,ngrid
319    do islope = 1,nslope
320     if((abs(tend_h2oglaciers(ig,islope)).gt.1e-5).and.(abs(waterice(ig,islope)).gt.0)) then
321        ispermanent_h2oglaciers(ig,islope) = 1
322     else
323        ispermanent_h2oglaciers(ig,islope) = 0
324     endif
325
326     if((abs(tend_co2glaciers(ig,islope)).gt.1e-5).and.(abs(co2ice(ig,islope)).gt.0)) then
327        ispermanent_co2glaciers(ig,islope) = 1
328     else
329        ispermanent_co2glaciers(ig,islope) = 0
330     endif
331    enddo
332   enddo
333
334!   0.2  Compute the partial pressure of CO2
335!a. the molecular mass into the column
336     do ig = 1,ngrid
337       mass_mean(ig,:) = 1./(A*q_co2(ig,:) +B)
338     enddo
339
340! b. pressure level
341     do it = 1,timelen
342       do ig = 1,ngrid
343         zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it)
344       enddo
345     enddo
346
347! c. Vapor pressure
348     pco2(:,:) = mass_mean(:,:)/m_co2*q_co2(:,:)*zplev_mean(:,:)
349     pco2_avg(:) = sum(pco2(:,:),2)/timelen
350
351     deallocate(zplev_mean)
352     deallocate(mass_mean)
353     deallocate(pco2)
354
355
356! 1. Compute the fraction of the pores occupied by H2O
357
358 call regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen,tend_h2oglaciers,tend_co2glaciers,waterice,co2ice,ps,q_co2,q_h2o,tsoil_PEM,TI_PEM, &
359                                   theta_h2o_adsorbed, m_h2o_adsorbed,delta_mh2o)
360
361
362
363! 2.  we compute the mass of co2 adsorded in each layer of the meshes 
364
365 do ig = 1,ngrid
366  do islope = 1,nslope
367    do iloop = 1,index_breccia
368     if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then
369     dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1-theta_h2o_adsorbed(ig,iloop,islope))*alpha*pco2_avg(ig)/ &
370                                             (alpha*pco2_avg(ig)+sqrt(tsoil_PEM(ig,iloop,islope))*exp(beta/tsoil_PEM(ig,iloop,islope)))
371     else
372        if(abs(m_co2_completesoil(ig,iloop,islope)).lt.(1e-10)) then !!! we are at first call
373          dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1-theta_h2o_adsorbed(ig,iloop,islope))*alpha*pco2_avg(ig) &
374                                                 /(alpha*pco2_avg(ig)+sqrt(tsoil_PEM(ig,iloop,islope))*exp(beta/tsoil_PEM(ig,iloop,islope)))
375        else ! no change: permanent ice stick the atoms of CO2
376          dm_co2_regolith_slope(ig,iloop,islope) = m_co2_completesoil(ig,iloop,islope)
377        endif
378     endif
379    enddo
380 enddo
381enddo
382
383! 3. Check the exchange between the atmosphere and the regolith
384
385  do ig = 1,ngrid
386   delta_mreg(ig) = 0.
387   do islope = 1,nslope
388    deltam_reg_slope(ig,islope) = 0.
389    do iloop = 1,index_breccia
390       if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then
391             deltam_reg_complete(ig,iloop,islope) = (dm_co2_regolith_slope(ig,iloop,islope) - m_co2_completesoil(ig,iloop,islope)) &
392                                                   *(layer_PEM(iloop+1) - layer_PEM(iloop))
393       else ! NO EXCHANGE AS ICE BLOCK THE DYNAMIC!
394             deltam_reg_complete(ig,iloop,islope) = 0.
395       endif
396       deltam_reg_slope(ig,islope) = deltam_reg_slope(ig,islope) + deltam_reg_complete(ig,iloop,islope)
397    enddo
398   delta_mreg(ig) = delta_mreg(ig) + deltam_reg_slope(ig,islope)*subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)
399   enddo
400  enddo
401  m_co2_completesoil(:,:,:) = dm_co2_regolith_slope(:,:,:)
402
403!=======================================================================
404      RETURN
405#endif
406      END
407
408
409   end module
Note: See TracBrowser for help on using the repository browser.