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

Last change on this file since 2950 was 2944, checked in by llange, 3 years ago

PEM
Introduce module 'constants_marspem': Constants that are used are now put in this module, and load in each subroutine.
LL

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