module adsorption_mod implicit none contains subroutine regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen, ps,q_co2,q_h2o,tsoil_PEM,TI_PEM, theta_h2o_adsorbded,m_h2o_adsorbed) use vertical_layers_mod, ONLY: ap,bp use comsoil_h_PEM, only: n_1km implicit none INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen ! size dimension REAL,INTENT(IN) :: ps(ngrid,timelen) ! surface pressure (Pa) REAL,INTENT(IN) :: q_co2(ngrid,timelen) ! Mass mixing ratio of co2 in the first layer (kg/kg) REAL,INTENT(IN) :: q_h2o(ngrid,timelen) ! Mass mixing ratio of H2o in the first layer (kg/kg) REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope) ! Soil Thermal inertia (J/K/^2/s^1/2) REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope) ! Soil temperature (K) ! output REAL,INTENT(OUT) :: m_h2o_adsorbed(ngrid,nsoil_PEM,nslope) ! Density of h2o adsorbed (kg/m^3) REAL,INTENT(OUT) :: theta_h2o_adsorbded(ngrid,nsoil_PEM,nslope) ! Fraction of the pores occupied by H2O molecules ! constant REAL :: Ko = 1.57e-8 ! Jackosky et al. 1997 REAL :: e = 2573.9 ! Jackosky et al. 1997 REAL :: mu = 0.48 ! Jackosky et al. 1997 REAL :: inertie_thresold = 800. ! TI > 800 means cementation real :: m_h2o = 18.01528E-3 ! Molecular weight of h2o (kg/mol) real :: m_co2 = 44.01E-3 ! Molecular weight of co2 (kg/mol) real :: m_noco2 = 33.37E-3 ! Molecular weight of non co2 (kg/mol) REAL :: rho_regolith = 2000. ! density of the reoglith, Buhler & Piqueux 2021 real :: alpha_clapeyron = -6143.7! eq. (2) in Murphy & Koop 2005 real :: beta_clapeyron = 29.9074 ! eq. (2) in Murphy & Koop 2005 real :: mi = 2.84e-7 ! Mass of h2o per m^2 absorbed Jackosky et al. 1997 real :: as = 18.9e3 ! Specific area, Buhler & Piqueux 2021 ! local variable real :: A,B ! Used to compute the mean mass above the surface real :: K ! Used to compute theta real :: p_sat ! saturated vapor pressure of ice integer ig,iloop, islope,isoil,it ! for loops real,allocatable :: mass_mean(:,:) ! mean mass above the surface real,allocatable :: zplev_mean(:,:) ! pressure above the surface real,allocatable :: pvapor(:,:) ! partial pressure above the surface real, allocatable :: pvapor_avg(:) ! yearly average vapor pressure above the surface ! 0. Some initializations allocate(mass_mean(ngrid,timelen)) allocate(zplev_mean(ngrid,timelen)) allocate(pvapor(ngrid,timelen)) allocate(pvapor_avg(ngrid)) m_h2o_adsorbed(:,:,:) = 0. theta_h2o_adsorbded(:,:,:) = 0. A =(1/m_co2 - 1/m_noco2) B=1/m_noco2 ! 1. Compute rho surface yearly averaged ! 1.1 Compute the partial pressure of vapor !a. the molecular mass into the column do ig = 1,ngrid mass_mean(ig,:) = 1/(A*q_co2(ig,:) +B) enddo ! b. pressure level do it = 1,timelen do ig = 1,ngrid zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it) enddo enddo ! c. Vapor pressure pvapor(:,:) = mass_mean(:,:)/m_h2o*q_h2o(:,:)*zplev_mean(:,:) pvapor_avg(:) = sum(pvapor(:,:),2)/timelen deallocate(pvapor) deallocate(zplev_mean) deallocate(mass_mean) ! 2. we compute the mass of co2 adsorded in each layer of the meshes do ig = 1,ngrid do islope = 1,nslope do iloop = 1,n_1km if(TI_PEM(ig,iloop,islope).lt.inertie_thresold) then K = Ko*exp(e/tsoil_PEM(ig,iloop,islope)) theta_h2o_adsorbded(ig,iloop,islope) = (K*pvapor_avg(ig)/(1+K*pvapor_avg(ig)))**mu m_h2o_adsorbed(ig,iloop,islope) = as*theta_h2o_adsorbded(ig,iloop,islope)*mi*rho_regolith else p_sat =exp(alpha_clapeyron/tsoil_PEM(ig,iloop,islope) +beta_clapeyron) ! we assume fixed temperature in the ice ... not really:q good but ... theta_h2o_adsorbded(ig,iloop,islope) = (K*p_sat/(1+K*p_sat))**mu m_h2o_adsorbed(ig,iloop,islope) =as*theta_h2o_adsorbded(ig,iloop,islope)*mi*rho_regolith endif enddo enddo enddo RETURN end subroutine 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) use comsoil_h_PEM, only: layer_PEM, mlayer_PEM,n_1km USE comcstfi_h, only: r, cpp, mugaz, g, rcp, pi, rad use comslope_mod, only : subslope_dist,def_slope_mean use vertical_layers_mod, ONLY: ap,bp IMPLICIT NONE ! Input: INTEGER,INTENT(IN) :: ngrid, nslope, nsoil_PEM,timelen ! size dimension REAL,INTENT(IN) :: ps(ngrid,timelen) ! Average surface pressure [Pa] REAL,INTENT(IN) :: tsoil_PEM(ngrid,nsoil_PEM,nslope) ! Mean Soil Temperature [K] REAL,INTENT(IN) :: TI_PEM(ngrid,nsoil_PEM,nslope) ! Mean Thermal Inertia [USI] REAL,INTENT(IN) :: tend_h2oglaciers(ngrid,nslope),tend_co2glaciers(ngrid,nslope) !tendancies on the glaciers () REAL,INTENT(IN) :: q_co2(ngrid,timelen),q_h2o(ngrid,timelen) ! Mass mixing ratio of co2 and h2o in the first layer (kg/kg) REAL,INTENT(IN) :: waterice(ngrid,nslope) ! water ice at the surface [kg/m^2] REAL,INTENT(IN) :: co2ice(ngrid,nslope) ! co2 ice at the surface [kg/m^2] ! Outputs: REAL,INTENT(INOUT) :: m_co2_completesoil(ngrid,nsoil_PEM,nslope) ! Density of co2 adsorbed (kg/m^3) REAL,INTENT(INOUT) :: delta_mreg(ngrid) ! Difference density of co2 adsorbed (kg/m^3) ! Constants: REAL :: alpha = 7.512e-6 ! Zent & Quinn 1995 REAL :: beta = -1541.5 ! Zent & Quinn 1995 REAL :: inertie_thresold = 800. ! TI > 800 means cementation REAL :: rho_regolith = 2000. ! density of the reoglith, buhler & piqueux 2021 real :: m_co2 = 44.01E-3 ! Molecular weight of co2 (kg/mol) real :: m_noco2 = 33.37E-3 ! Molecular weight of h2o (kg/mol) real :: m_theta = 4.27e-7 ! Mass of co2 per m^2 absorbed real :: as = 18.9e3 ! Specific area, Buhler & Piqueux 2021 ! Local real :: A,B ! Used to compute the mean mass above the surface INTEGER :: ig,islope,iloop,it ! for loops REAL :: dm_co2_regolith_slope(ngrid,nsoil_PEM,nslope) ! elementary mass adsorded per mesh per slope INTEGER :: ispermanent_co2glaciers(ngrid,nslope) ! Check if the glacier is permanent INTEGER :: ispermanent_h2oglaciers(ngrid,nslope) ! Check if the glacier is permanent REAL :: deltam_reg_complete(ngrid,n_1km,nslope) ! Difference in the mass per slope and soil layer (kg/m^3) REAL :: deltam_reg_slope(ngrid,nslope) ! Difference in the mass per slope (kg/m^3) REAL :: m_h2o_adsorbed(ngrid,nsoil_PEM,nslope) ! Density of CO2 adsorbed (kg/m^3) REAL :: theta_h2o_adsorbed(ngrid,nsoil_PEM,nslope) ! Fraction of the pores occupied by H2O molecules !timelen array are allocated because heavy ... real,allocatable :: mass_mean(:,:) ! mean mass above the surface real,allocatable :: zplev_mean(:,:) ! pressure above the surface real,allocatable :: pco2(:,:) ! partial pressure above the surface real, allocatable :: pco2_avg(:) ! yearly averaged ! 0. Some initializations allocate(mass_mean(ngrid,timelen)) allocate(zplev_mean(ngrid,timelen)) allocate(pco2(ngrid,timelen)) allocate(pco2_avg(ngrid)) m_h2o_adsorbed(:,:,:) = 0. A =(1/m_co2 - 1/m_noco2) B=1/m_noco2 dm_co2_regolith_slope(:,:,:) = 0 delta_mreg(:) = 0. !0.1 Look at perenial ice do ig = 1,ngrid do islope = 1,nslope if((abs(tend_h2oglaciers(ig,islope)).gt.1e-5).and.(abs(waterice(ig,islope)).gt.0)) then ispermanent_h2oglaciers(ig,islope) = 1 else ispermanent_h2oglaciers(ig,islope) = 0 endif if((abs(tend_co2glaciers(ig,islope)).gt.1e-5).and.(abs(co2ice(ig,islope)).gt.0)) then ispermanent_co2glaciers(ig,islope) = 1 else ispermanent_co2glaciers(ig,islope) = 0 endif enddo enddo ! 0.2 Compute the partial pressure of CO2 !a. the molecular mass into the column do ig = 1,ngrid mass_mean(ig,:) = 1/(A*q_co2(ig,:) +B) enddo ! b. pressure level do it = 1,timelen do ig = 1,ngrid zplev_mean(ig,it) = ap(1) + bp(1)*ps(ig,it) enddo enddo ! c. Vapor pressure pco2(:,:) = mass_mean(:,:)/m_co2*q_co2(:,:)*zplev_mean(:,:) pco2_avg(:) = sum(pco2(:,:),2)/timelen deallocate(zplev_mean) deallocate(mass_mean) deallocate(pco2) ! 1. Compute the fraction of the pores occupied by H2O call regolith_h2oadsorption(ngrid,nslope,nsoil_PEM,timelen, ps,q_co2,q_h2o,tsoil_PEM,TI_PEM,theta_h2o_adsorbed, m_h2o_adsorbed) ! 2. we compute the mass of co2 adsorded in each layer of the meshes do ig = 1,ngrid do islope = 1,nslope do iloop = 1,n_1km if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1-theta_h2o_adsorbed(ig,iloop,islope))*alpha*pco2_avg(ig)/ & (alpha*pco2_avg(ig)+sqrt(tsoil_PEM(ig,iloop,islope))*exp(beta/tsoil_PEM(ig,iloop,islope))) else if(abs(m_co2_completesoil(ig,iloop,islope)).lt.1-10) then !!! we are at first call dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1-theta_h2o_adsorbed(ig,iloop,islope))*alpha*pco2_avg(ig) & /(alpha*pco2_avg(ig)+sqrt(tsoil_PEM(ig,iloop,islope))*exp(beta/tsoil_PEM(ig,iloop,islope))) else ! no change: permanent ice stick the atoms of CO2 dm_co2_regolith_slope(ig,iloop,islope) = m_co2_completesoil(ig,iloop,islope) endif endif enddo enddo enddo ! 2. Check the exchange between the atmosphere and the regolith do ig = 1,ngrid delta_mreg(ig) = 0. do islope = 1,nslope deltam_reg_slope(ig,islope) = 0. do iloop = 1,n_1km if((TI_PEM(ig,iloop,islope).lt.inertie_thresold).and.(ispermanent_h2oglaciers(ig,islope).eq.0).and.(ispermanent_co2glaciers(ig,islope).eq.0)) then deltam_reg_complete(ig,iloop,islope) = (dm_co2_regolith_slope(ig,iloop,islope) - m_co2_completesoil(ig,iloop,islope)) & *(layer_PEM(iloop+1) - layer_PEM(iloop)) else ! NO EXCHANGE AS ICE BLOCK THE DYNAMIC! deltam_reg_complete(ig,iloop,islope) = 0. endif deltam_reg_slope(ig,islope) = deltam_reg_slope(ig,islope) + deltam_reg_complete(ig,iloop,islope) enddo delta_mreg(ig) = delta_mreg(ig) + deltam_reg_slope(ig,islope)*subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.) enddo enddo m_co2_completesoil(:,:,:) = dm_co2_regolith_slope(:,:,:) !======================================================================= RETURN END end module