Changeset 4117 for trunk/LMDZ.COMMON


Ignore:
Timestamp:
Mar 10, 2026, 4:52:43 PM (8 days ago)
Author:
jbclement
Message:

PEM:

  • Remove 'd_h2oice' and 'd_co2ice' arguments from adsorption and soil thermal inertia routines (not needed and possible bugs).
  • Move 'read_callphys' call from 'read_start1D' to 'read_rundef' so that 'CO2cond_ps' is initialized no matter what.
  • Explicit initialization of 'delta_h2o_ads'/'delta_co2_ads' before 'read_startpem'.
  • Add time index argument to 'put_var_nc' calls to write "restart.nc".
  • Increase temporary string buffer size for 'real2str_qp'.
  • Fixed some typos.

JBC

Location:
trunk/LMDZ.COMMON/libf/evolution
Files:
13 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.COMMON/libf/evolution/atmosphere.F90

    r4110 r4117  
    181181! ----
    182182CO2cond_ps_PCM = CO2cond_ps_PCM_in
    183 call print_msg('set_CO2cond_ps_PCM = '//real2str(CO2cond_ps_PCM),LVL_NFO)
     183call print_msg('CO2cond_ps_PCM = '//real2str(CO2cond_ps_PCM),LVL_NFO)
    184184if (CO2cond_ps_PCM < 0._dp .or. CO2cond_ps_PCM > 1._dp) call stop_clean(__FILE__,__LINE__,'Value for ''CO2cond_ps'' is not between 0 and 1!',1)
    185185
  • trunk/LMDZ.COMMON/libf/evolution/changelog.txt

    r4110 r4117  
    899899- Addition of workflow safety checks for required executables, dependencies (e.g. 'ncdump'), input files and callphys keys.
    900900- Renaming of PEM starting and diagnostic files ("startevol.nc" into "startevo.nc", "diagevol.nc" into "diagevo.nc").
     901
     902== 10/03/2026 == JBC
     903- Remove 'd_h2oice' and 'd_co2ice' arguments from adsorption and soil thermal inertia routines (not needed and possible bugs).
     904- Move 'read_callphys' call from 'read_start1D' to 'read_rundef' so that 'CO2cond_ps' is initialized no matter what.
     905- Explicit initialization of 'delta_h2o_ads'/'delta_co2_ads' before 'read_startpem'.
     906- Add time index argument to 'put_var_nc' calls to write "restart.nc".
     907- Increase temporary string buffer size for 'real2str_qp'.
     908- Fixed some typos.
  • trunk/LMDZ.COMMON/libf/evolution/clim_state_init.F90

    r4110 r4117  
    125125use tracers,    only: nq, set_q_PCM
    126126use stoppage,   only: stop_clean
    127 use config,     only: read_callphys
    128127use display,    only: print_msg, LVL_NFO
    129128
     
    174173call set_ap(ap_tmp)
    175174call set_bp(bp_tmp)
    176 
    177 ! Read the "callphys.def"
    178 call read_callphys() ! To get 'CO2cond_ps'
    179175
    180176END SUBROUTINE read_start1D
     
    300296
    301297!=======================================================================
    302 SUBROUTINE read_startpem(tsurf_avg_yr1,tsurf_avg_yr2,ps_avg_glob,ps_ts,q_co2_ts,q_h2o_ts,h2o_surfdensity_avg,d_h2oice,d_co2ice,h2o_ice,co2_ice, &
    303                          tsoil_avg,h2o_soildensity_avg,icetable_depth,icetable_thickness,ice_porefilling,layerings_map,                         &
     298SUBROUTINE read_startpem(tsurf_avg_yr1,tsurf_avg_yr2,ps_avg_glob,ps_ts,q_co2_ts,q_h2o_ts,h2o_surfdensity_avg,h2o_ice,co2_ice, &
     299                         tsoil_avg,h2o_soildensity_avg,icetable_depth,icetable_thickness,ice_porefilling,layerings_map,       &
    304300                         h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    305301!-----------------------------------------------------------------------
     
    351347real(dp),       dimension(:,:),   intent(in)    :: q_co2_ts, q_h2o_ts           ! MMR of CO2 and H2O
    352348real(dp),       dimension(:,:),   intent(in)    :: h2o_surfdensity_avg          ! Average of surface water density
    353 real(dp),       dimension(:,:),   intent(in)    :: d_h2oice, d_co2ice           ! Tendencies on ice
    354349real(dp),       dimension(:,:),   intent(inout) :: h2o_ice, co2_ice             ! Surface ice
    355350real(dp),       dimension(:,:,:), intent(inout) :: tsoil_avg                    ! Soil temperature
     
    471466        if (icetable_equilibrium) then
    472467            call computeice_table_equilibrium(is_h2o_perice_PCM,h2o_surfdensity_avg,h2o_soildensity_avg,TI(:,1,:),icetable_depth,icetable_thickness)
    473             call update_soil_TI(d_h2oice,h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
     468            call update_soil_TI(h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
    474469            do islope = 1,nslope
    475470                call ini_tsoil_pem(ngrid,nsoil,TI(:,:,islope),tsurf_avg_yr2(:,islope),tsoil_avg(:,:,islope))
    476471            end do
    477472        else if (icetable_dynamic) then
    478             call update_soil_TI(d_h2oice,h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
     473            call update_soil_TI(h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
    479474            do islope = 1,nslope
    480475                call ini_tsoil_pem(ngrid,nsoil,TI(:,:,islope),tsurf_avg_yr2(:,islope),tsoil_avg(:,:,islope))
     
    486481            h2o_ads_reg = 0._dp
    487482            co2_ads_reg = 0._dp
    488             call evolve_regolith_adsorption(d_h2oice,d_co2ice,h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    489             delta_co2_ads = 0._dp
    490             delta_h2o_ads = 0._dp
     483            call evolve_regolith_adsorption(h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    491484        end if ! do_sorption
    492485    end if ! do_soil
     
    565558        call get_var_nc('co2_ads_reg',co2_ads_reg(:,:,:))
    566559        call get_var_nc('h2o_ads_reg',h2o_ads_reg(:,:,:))
    567         call evolve_regolith_adsorption(d_h2oice,d_co2ice,h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
     560        call evolve_regolith_adsorption(h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    568561    end if ! do_sorption
    569562end if ! do_soil
  • trunk/LMDZ.COMMON/libf/evolution/clim_state_rec.F90

    r4110 r4117  
    9595! Surface pressure
    9696call vect2dyngrd(nlon + 1,nlat,ngrid,ps4PCM,ps_dyn)
    97 call put_var_nc('ps',ps_dyn)
     97call put_var_nc('ps',ps_dyn,1)
    9898
    9999! Potential temperature
    100100call vect2dyngrd(nlon + 1,nlat,ngrid,teta4PCM,var_dyn)
    101 call put_var_nc('teta',var_dyn)
     101call put_var_nc('teta',var_dyn,1)
    102102
    103103! Air mass
    104104call vect2dyngrd(nlon + 1,nlat,ngrid,air_mass4PCM,var_dyn)
    105 call put_var_nc('masse',var_dyn)
     105call put_var_nc('masse',var_dyn,1)
    106106
    107107! Tracers
     
    110110        call vect2dyngrd(nlon + 1,nlat,ngrid,q4PCM(:,l,i),var_dyn(:,:,l))
    111111    end do
    112     call put_var_nc(qnames(i),var_dyn)
     112    call put_var_nc(qnames(i),var_dyn,1)
    113113end do
    114114
  • trunk/LMDZ.COMMON/libf/evolution/config.F90

    r4110 r4117  
    208208call set_atmosphere_config(hybrid)
    209209
     210! Read the "callphys.def"
     211call read_callphys() ! To get 'CO2cond_ps'
     212
    210213! Check incompatibilities
    211214call check_config_incompatibility()
  • trunk/LMDZ.COMMON/libf/evolution/deftank/pem_workflow_lib.sh

    r4110 r4117  
    2828    if ! command -v "$1" >/dev/null 2>&1; then
    2929        echo "Error: required command '$1' not found."
    30         echo "Please install or load the corresonding library."
     30        echo "Please install or load the corresponding library."
    3131        abort_workflow
    3232    fi
  • trunk/LMDZ.COMMON/libf/evolution/display.F90

    r4110 r4117  
    195195call print_msg('',LVL_NFO)
    196196call print_msg('****** PEM final information *******',LVL_NFO)
    197 write(msg,'(a,i0,a,f16.4,a)') ' + The run PEM #',i_pem_run,' achieved ', n_yr_run, ' Planetary years, completed in '//format_duration(dur_secs)//'.'
     197write(msg,'(a,i0,a,f16.4,a)') '+ The run PEM #',i_pem_run,' achieved ', n_yr_run, ' Planetary years, completed in '//format_duration(dur_secs)//'.'
    198198call print_msg(msg,LVL_NFO)
    199 write(msg,'(a,f16.4,a,f16.4,a)') ' + The workflow has achieved ', n_yr_sim, ' Planetary years =', n_yr_sim*r_plnt2earth_yr, ' Earth years.'
     199write(msg,'(a,f16.4,a,f16.4,a)') '+ The workflow has achieved ', n_yr_sim, ' Planetary years =', n_yr_sim*r_plnt2earth_yr, ' Earth years.'
    200200call print_msg(msg,LVL_NFO)
    201 write(msg,'(a,f16.4,a)') ' + The reached date is now ', (pem_ini_date + n_yr_sim)*r_plnt2earth_yr, ' Earth years.'
     201write(msg,'(a,f16.4,a)') '+ The reached date is now ', (pem_ini_date + n_yr_sim)*r_plnt2earth_yr, ' Earth years.'
    202202call print_msg(msg,LVL_NFO)
    203203call print_msg('+ PEM: so far, so good!',LVL_NFO)
  • trunk/LMDZ.COMMON/libf/evolution/orbit.F90

    r4110 r4117  
    9292max_change_ecc = max_change_ecc_in
    9393max_change_lsp = max_change_lsp_in
    94 call print_msg('evo_orbit     = '//bool2str(evo_orbit),LVL_NFO)
    95 call print_msg('evo_obl       = '//bool2str(evo_obl),LVL_NFO)
    96 call print_msg('evo_ecc       = '//bool2str(evo_ecc),LVL_NFO)
    97 call print_msg('evo_lsp       = '//bool2str(evo_lsp),LVL_NFO)
     94call print_msg('evo_orbit      = '//bool2str(evo_orbit),LVL_NFO)
     95call print_msg('evo_obl        = '//bool2str(evo_obl),LVL_NFO)
     96call print_msg('evo_ecc        = '//bool2str(evo_ecc),LVL_NFO)
     97call print_msg('evo_lsp        = '//bool2str(evo_lsp),LVL_NFO)
    9898call print_msg('max_change_obl = '//real2str(max_change_obl),LVL_NFO)
    9999call print_msg('max_change_ecc = '//real2str(max_change_ecc),LVL_NFO)
  • trunk/LMDZ.COMMON/libf/evolution/pem.F90

    r4110 r4117  
    208208allocate(h2o_ads_reg(ngrid,nsoil,nslope),co2_ads_reg(ngrid,nsoil,nslope),delta_h2o_ads(ngrid),delta_co2_ads(ngrid))
    209209allocate(layerings_map(ngrid,nslope))
    210 call read_startpem(tsurf_avg_yr1,tsurf_avg,ps_avg_glob_ini,ps_ts,q_co2_ts,q_h2o_ts,h2o_surfdensity_avg,d_h2oice,d_co2ice,h2o_ice,co2_ice, &
    211                    tsoil_avg,h2o_soildensity_avg,icetable_depth,icetable_thickness,ice_porefilling,layerings_map,                         &
     210delta_h2o_ads(:) = 0._dp
     211delta_co2_ads(:) = 0._dp
     212call read_startpem(tsurf_avg_yr1,tsurf_avg,ps_avg_glob_ini,ps_ts,q_co2_ts,q_h2o_ts,h2o_surfdensity_avg,h2o_ice,co2_ice, &
     213                   tsoil_avg,h2o_soildensity_avg,icetable_depth,icetable_thickness,ice_porefilling,layerings_map,       &
    212214                   h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    213215deallocate(tsurf_avg_yr1)
     
    342344
    343345        ! Update soil thermal properties according to ice table and soil temperature
    344         call update_soil_TI(d_h2oice,h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
     346        call update_soil_TI(h2o_ice,ps_avg_glob,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
    345347
    346348        ! Evolve adsorbed species
    347349        if (do_sorption) then
    348             call evolve_regolith_adsorption(d_h2oice,d_co2ice,h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
     350            call evolve_regolith_adsorption(h2o_ice,co2_ice,tsoil_avg,TI,ps_ts,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    349351            call compute_totmass_adsorbed(h2o_ads_reg,co2_ads_reg,totmass_adsco2,totmass_adsh2o)
    350352        else
  • trunk/LMDZ.COMMON/libf/evolution/soil_therm_inertia.F90

    r4110 r4117  
    8585
    8686!=======================================================================
    87 SUBROUTINE update_soil_TI(tendencies_waterice,waterice,p_avg_new,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
     87SUBROUTINE update_soil_TI(h2o_ice,p_avg_new,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI)
    8888!-----------------------------------------------------------------------
    8989! NAME
     
    115115! ARGUMENTS
    116116! ---------
    117 real(dp),                   intent(in)    :: p_avg_new            ! Global average surface pressure [Pa]
    118 real(dp), dimension(:,:),   intent(in)    :: tendencies_waterice  ! Tendencies on the water ice [kg/m^2/year]
    119 real(dp), dimension(:,:),   intent(in)    :: waterice             ! Surface Water ice [kg/m^2]
    120 real(dp), dimension(:,:),   intent(in)    :: icetable_depth       ! Ice table depth [m]
    121 real(dp), dimension(:,:),   intent(in)    :: icetable_thickness   ! Ice table thickness [m]
    122 real(dp), dimension(:,:,:), intent(in)    :: ice_porefilling      ! Ice porefilling [m^3/m^3]
     117real(dp),                   intent(in)    :: p_avg_new          ! Global average surface pressure [Pa]
     118real(dp), dimension(:,:),   intent(in)    :: h2o_ice            ! Surface Water ice [kg/m^2]
     119real(dp), dimension(:,:),   intent(in)    :: icetable_depth     ! Ice table depth [m]
     120real(dp), dimension(:,:),   intent(in)    :: icetable_thickness ! Ice table thickness [m]
     121real(dp), dimension(:,:,:), intent(in)    :: ice_porefilling    ! Ice porefilling [m^3/m^3]
    123122logical(k4),                intent(in)    :: icetable_equilibrium, icetable_dynamic ! Computing method for ice table
    124 real(dp), dimension(:,:,:), intent(inout) :: TI                   ! Soil Thermal Inertia [J/m^2/K/s^1/2]
     123real(dp), dimension(:,:,:), intent(inout) :: TI                 ! Soil Thermal Inertia [J/m^2/K/s^1/2]
    125124
    126125! LOCAL VARIABLES
     
    141140    regolith_inertia(:,islope) = inertiedat(:,1)
    142141    do ig = 1,ngrid
    143         if (tendencies_waterice(ig,islope) < -1.e-5_dp .and. abs(waterice(ig,islope)) < minieps) regolith_inertia(ig,islope) = TI_regolith_avg
     142        if (abs(h2o_ice(ig,islope)) < minieps) regolith_inertia(ig,islope) = TI_regolith_avg
    144143        if (reg_thprop_dependp) then
    145144            if (TI(ig,1,islope) < reg_inertie_thresold) then
     
    153152end do
    154153
    155 ! 2. Build new Thermal Inertia
     154! 2. Build new thermal inertia
    156155do islope = 1,nslope
    157156    do ig = 1,ngrid
  • trunk/LMDZ.COMMON/libf/evolution/sorption.F90

    r4110 r4117  
    6969
    7070!=======================================================================
    71 SUBROUTINE evolve_regolith_adsorption(d_h2oice,d_co2ice,h2oice,co2ice,tsoil,TI,ps,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
     71SUBROUTINE evolve_regolith_adsorption(h2o_ice,co2_ice,tsoil,TI,ps,q_h2o_ts,q_co2_ts,h2o_ads_reg,co2_ads_reg,delta_h2o_ads,delta_co2_ads)
    7272!-----------------------------------------------------------------------
    7373! NAME
     
    9696! ARGUMENTS
    9797! ---------
    98 real(dp), dimension(:,:),   intent(in)    :: d_h2oice, d_co2ice ! tendencies on the glaciers [1]
    99 real(dp), dimension(:,:),   intent(in)    :: h2oice             ! H2O ice at the surface [kg/m^2]
    100 real(dp), dimension(:,:),   intent(in)    :: co2ice             ! CO2 ice at the surface [kg/m^2]
     98real(dp), dimension(:,:),   intent(in)    :: h2o_ice            ! H2O ice at the surface [kg/m^2]
     99real(dp), dimension(:,:),   intent(in)    :: co2_ice            ! CO2 ice at the surface [kg/m^2]
    101100real(dp), dimension(:,:,:), intent(in)    :: TI                 ! Soil Thermal inertia [J/K/^2/s^1/2]
    102101real(dp), dimension(:,:,:), intent(in)    :: tsoil              ! Soil temperature [K]
     
    114113! CODE
    115114! ----
    116 call evolve_h2o_ads(d_h2oice,d_co2ice,h2oice,co2ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,theta_h2o_ads,h2o_ads_reg,delta_h2o_ads)
    117 call evolve_co2_ads(d_h2oice,d_co2ice,h2oice,co2ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,co2_ads_reg,delta_co2_ads)
     115call evolve_h2o_ads(h2o_ice,co2_ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,theta_h2o_ads,h2o_ads_reg,delta_h2o_ads)
     116call evolve_co2_ads(h2o_ice,co2_ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,co2_ads_reg,delta_co2_ads)
    118117
    119118END SUBROUTINE evolve_regolith_adsorption
     
    121120
    122121!=======================================================================
    123 SUBROUTINE evolve_h2o_ads(d_h2oice,d_co2ice,h2oice,co2ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,theta_h2o_ads,h2o_ads_reg,delta_h2o_ads)
     122SUBROUTINE evolve_h2o_ads(h2o_ice,co2_ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,theta_h2o_ads,h2o_ads_reg,delta_h2o_ads)
    124123!-----------------------------------------------------------------------
    125124! NAME
     
    152151! ARGUMENTS
    153152! ---------
    154 real(dp), dimension(:,:),   intent(in)    :: d_h2oice, d_co2ice ! Tendencies on the glaciers ()
    155 real(dp), dimension(:,:),   intent(in)    :: h2oice             ! H2O ice at the surface [kg/m^2]
    156 real(dp), dimension(:,:),   intent(in)    :: co2ice             ! CO2 ice at the surface [kg/m^2]
    157 real(dp), dimension(:,:),   intent(in)    :: ps                 ! Surface pressure [Pa]
    158 real(dp), dimension(:,:),   intent(in)    :: q_co2_ts           ! Mass mixing ratio of CO2 in the first layer [kg/kg]
    159 real(dp), dimension(:,:),   intent(in)    :: q_h2o_ts           ! Mass mixing ratio of H2O in the first layer [kg/kg]
    160 real(dp), dimension(:,:,:), intent(in)    :: TI                 ! Soil Thermal inertia [J/K/^2/s^1/2]
    161 real(dp), dimension(:,:,:), intent(in)    :: tsoil              ! Soil temperature [K]
    162 real(dp), dimension(:,:,:), intent(inout) :: h2o_ads_reg        ! Density of H2O adsorbed [kg/m^3]
    163 real(dp), dimension(:,:,:), intent(out)   :: theta_h2o_ads      ! Fraction of the pores occupied by H2O molecules
    164 real(dp), dimension(:),     intent(out)   :: delta_h2o_ads      ! Difference density of H2O adsorbed [kg/m^3]
     153real(dp), dimension(:,:),   intent(in)    :: h2o_ice       ! H2O ice at the surface [kg/m^2]
     154real(dp), dimension(:,:),   intent(in)    :: co2_ice       ! CO2 ice at the surface [kg/m^2]
     155real(dp), dimension(:,:),   intent(in)    :: ps            ! Surface pressure [Pa]
     156real(dp), dimension(:,:),   intent(in)    :: q_co2_ts      ! Mass mixing ratio of CO2 in the first layer [kg/kg]
     157real(dp), dimension(:,:),   intent(in)    :: q_h2o_ts      ! Mass mixing ratio of H2O in the first layer [kg/kg]
     158real(dp), dimension(:,:,:), intent(in)    :: TI            ! Soil Thermal inertia [J/K/^2/s^1/2]
     159real(dp), dimension(:,:,:), intent(in)    :: tsoil         ! Soil temperature [K]
     160real(dp), dimension(:,:,:), intent(inout) :: h2o_ads_reg   ! Density of H2O adsorbed [kg/m^3]
     161real(dp), dimension(:,:,:), intent(out)   :: theta_h2o_ads ! Fraction of the pores occupied by H2O molecules
     162real(dp), dimension(:),     intent(out)   :: delta_h2o_ads ! Difference density of H2O adsorbed [kg/m^3]
    165163
    166164! LOCAL VARIABLES
     
    177175real(dp)                                        :: K                       ! Used to compute theta
    178176integer(di)                                     :: ig, iloop, islope, it   ! For loops
    179 logical(k4),  dimension(ngrid,nslope)               :: ispermanent_co2glaciers ! Check if the CO2 glacier is permanent
    180 logical(k4),  dimension(ngrid,nslope)               :: ispermanent_h2oglaciers ! Check if the H2O glacier is permanent
    181 real(dp), dimension(ngrid,nslope)               :: deltam_reg_slope        ! Difference density of H2O adsorbed per slope [kg/m^3]
    182 real(dp), dimension(ngrid,nsoil,nslope)         :: dm_h2o_regolith_slope   ! Elementary H2O mass adsorded per mesh per slope
     177real(dp),    dimension(ngrid,nslope)            :: deltam_reg_slope        ! Difference density of H2O adsorbed per slope [kg/m^3]
     178real(dp),    dimension(ngrid,nsoil,nslope)      :: dm_h2o_regolith_slope   ! Elementary H2O mass adsorded per mesh per slope
    183179real(dp)                                        :: A, B                    ! Used to compute the mean mass above the surface
    184180real(dp)                                        :: p_sat                   ! Saturated vapor pressure of ice
    185 real(dp), dimension(:,:), allocatable           :: mass_mean               ! Mean mass above the surface
    186 real(dp), dimension(:,:), allocatable           :: zplev_mean              ! Pressure above the surface
    187 real(dp), dimension(:,:), allocatable           :: pvapor                  ! Partial pressure above the surface
    188 real(dp), dimension(:)  , allocatable           :: pvapor_avg              ! Yearly averaged
     181real(dp),    dimension(:,:), allocatable        :: mass_mean               ! Mean mass above the surface
     182real(dp),    dimension(:,:), allocatable        :: zplev_mean              ! Pressure above the surface
     183real(dp),    dimension(:,:), allocatable        :: pvapor                  ! Partial pressure above the surface
     184real(dp),    dimension(:)  , allocatable        :: pvapor_avg              ! Yearly averaged
    189185
    190186! CODE
     
    196192theta_h2o_ads = 0._dp
    197193dm_h2o_regolith_slope = 0._dp
    198 ispermanent_h2oglaciers = .false.
    199 ispermanent_co2glaciers = .false.
    200 
    201 ! 0.1 Look at perennial ice
    202 where (abs(d_h2oice(:,:)) > 1.e-5 .and. h2oice(:,:) > 0._dp) ispermanent_h2oglaciers(:,:) = .true.
    203 where (abs(d_co2ice(:,:)) > 1.e-5 .and. co2ice(:,:) > 0._dp) ispermanent_co2glaciers(:,:) = .true.
    204 
    205 ! 0.2 Compute the partial pressure of vapor
     194
     195! 0. Compute the partial pressure of vapor
    206196! a. the molecular mass into the column
    207197do ig = 1,ngrid
     
    243233        deltam_reg_slope(ig,islope) = 0._dp
    244234        do iloop = 1,index_breccia
    245             if (TI(ig,iloop,islope) < inertie_thresold .and. .not. ispermanent_h2oglaciers(ig,islope) .and. .not. ispermanent_co2glaciers(ig,islope)) then
     235            if (TI(ig,iloop,islope) < inertie_thresold .and. abs(h2o_ice(ig,islope)) < minieps .and. abs(co2_ice(ig,islope)) < minieps) then
    246236                if (iloop == 1) then
    247237                    deltam_reg_complete(ig,iloop,islope) = (dm_h2o_regolith_slope(ig,iloop,islope) - h2o_ads_reg(ig,iloop,islope))*(layer(iloop))
     
    263253
    264254!=======================================================================
    265 SUBROUTINE evolve_co2_ads(d_h2oice,d_co2ice,h2oice,co2ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,co2_ads_reg,delta_co2_ads)
     255SUBROUTINE evolve_co2_ads(h2o_ice,co2_ice,ps,q_h2o_ts,q_co2_ts,tsoil,TI,co2_ads_reg,delta_co2_ads)
    266256
    267257!-----------------------------------------------------------------------
     
    298288real(dp), dimension(:,:,:), intent(in)    :: tsoil              ! Mean Soil Temperature [K]
    299289real(dp), dimension(:,:,:), intent(in)    :: TI                 ! Mean Thermal Inertia [USI]
    300 real(dp), dimension(:,:),   intent(in)    :: d_h2oice, d_co2ice ! Tendencies on the glaciers ()
    301290real(dp), dimension(:,:),   intent(in)    :: q_co2_ts, q_h2o_ts ! Mass mixing ratio of CO2 and H2O in the first layer [kg/kg]
    302 real(dp), dimension(:,:),   intent(in)    :: h2oice             ! Water ice at the surface [kg/m^2]
    303 real(dp), dimension(:,:),   intent(in)    :: co2ice             ! CO2 ice at the surface [kg/m^2]
     291real(dp), dimension(:,:),   intent(in)    :: h2o_ice            ! Water ice at the surface [kg/m^2]
     292real(dp), dimension(:,:),   intent(in)    :: co2_ice            ! CO2 ice at the surface [kg/m^2]
    304293real(dp), dimension(:,:,:), intent(inout) :: co2_ads_reg        ! Density of CO2 adsorbed [kg/m^3]
    305294real(dp), dimension(:),     intent(out)   :: delta_co2_ads      ! Difference density of CO2 adsorbed [kg/m^3]
     
    317306integer(di)                                        :: ig, islope, iloop, it   ! For loops
    318307real(dp),    dimension(ngrid,nsoil,nslope)         :: dm_co2_regolith_slope   ! Elementary mass adsorded per mesh per slope
    319 logical(k4), dimension(ngrid,nslope)               :: ispermanent_co2glaciers ! Check if the CO2 glacier is permanent
    320 logical(k4), dimension(ngrid,nslope)               :: ispermanent_h2oglaciers ! Check if the H2O glacier is permanent
    321308real(dp),    dimension(ngrid,index_breccia,nslope) :: deltam_reg_complete     ! Difference in the mass per slope and soil layer [kg/m^3]
    322309real(dp),    dimension(ngrid,nslope)               :: deltam_reg_slope        ! Difference in the mass per slope  [kg/m^3]
     
    339326dm_co2_regolith_slope = 0._dp
    340327delta_co2_ads = 0._dp
    341 ispermanent_h2oglaciers = .false.
    342 ispermanent_co2glaciers = .false.
    343 
    344 ! 0.1 Look at perennial ice
    345 where (abs(d_h2oice(:,:)) > 1.e-5 .and. h2oice(:,:) > 0._dp) ispermanent_h2oglaciers(:,:) = .true.
    346 where (abs(d_co2ice(:,:)) > 1.e-5 .and. co2ice(:,:) > 0._dp) ispermanent_co2glaciers(:,:) = .true.
    347 
    348 ! 0.2  Compute the partial pressure of CO2
     328
     329! 0.  Compute the partial pressure of CO2
    349330! a. the molecular mass into the column
    350331do ig = 1,ngrid
     
    366347
    367348! 1. Compute the fraction of the pores occupied by H2O
    368 call evolve_h2o_ads(d_h2oice,d_co2ice,h2oice,co2ice,ps,q_co2_ts,q_h2o_ts,tsoil,TI,theta_h2o_ads,m_h2o_adsorbed,delta_mh2o)
     349call evolve_h2o_ads(h2o_ice,co2_ice,ps,q_co2_ts,q_h2o_ts,tsoil,TI,theta_h2o_ads,m_h2o_adsorbed,delta_mh2o)
    369350
    370351! 2. we compute the mass of CO2 adsorded in each layer of the meshes
     
    372353    do islope = 1,nslope
    373354        do iloop = 1,index_breccia
    374             if (TI(ig,iloop,islope) < inertie_thresold .and. .not. ispermanent_h2oglaciers(ig,islope) .and. .not. ispermanent_co2glaciers(ig,islope)) then
     355            if (TI(ig,iloop,islope) < inertie_thresold .and. abs(h2o_ice(ig,islope)) < minieps .and. abs(co2_ice(ig,islope)) < minieps) then
    375356                dm_co2_regolith_slope(ig,iloop,islope) = as*rho_regolith*m_theta*(1._dp - theta_h2o_ads(ig,iloop,islope))*alpha*pco2_avg(ig)/ &
    376357                                                         (alpha*pco2_avg(ig) + sqrt(tsoil(ig,iloop,islope))*exp(beta/tsoil(ig,iloop,islope)))
     
    393374        deltam_reg_slope(ig,islope) = 0._dp
    394375        do iloop = 1,index_breccia
    395             if (TI(ig,iloop,islope) < inertie_thresold .and. .not. ispermanent_h2oglaciers(ig,islope) .and. .not. ispermanent_co2glaciers(ig,islope)) then
     376            if (TI(ig,iloop,islope) < inertie_thresold .and. abs(h2o_ice(ig,islope)) < minieps .and. abs(co2_ice(ig,islope)) < minieps) then
    396377                if (iloop == 1) then
    397378                    deltam_reg_complete(ig,iloop,islope) = (dm_co2_regolith_slope(ig,iloop,islope) - co2_ads_reg(ig,iloop,islope))*(layer(iloop))
  • trunk/LMDZ.COMMON/libf/evolution/utility.F90

    r4110 r4117  
    209209! ---------------
    210210!~ integer(di)               :: len_trimmed
    211 character(32)             :: str_tmp
     211character(128)            :: str_tmp
    212212character(:), allocatable :: str
    213213
  • trunk/LMDZ.COMMON/libf/evolution/workflow_status.F90

    r4110 r4117  
    7676close(funit)
    7777
    78 call print_msg('Current PEM run: '//int2str(i_pem_run),LVL_NFO)
    79 call print_msg('Input PCM runs : '//int2str(i_pcm_run - 1)//' and '//int2str(i_pcm_run),LVL_NFO)
     78call print_msg('Current PEM run: #'//int2str(i_pem_run),LVL_NFO)
     79call print_msg('Input PCM runs : #'//int2str(i_pcm_run - 1)//' and #'//int2str(i_pcm_run),LVL_NFO)
    8080
    8181END SUBROUTINE read_workflow_status
Note: See TracChangeset for help on using the changeset viewer.