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

Last change on this file since 3026 was 2985, checked in by romain.vande, 18 months ago

Generic PEM :

Adapt PEM to run with the generic model.
(CPP_STD keyword to exclude some part of the code at the compilation)

RV

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