| 1 | MODULE tendencies |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! tendencies |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Computation and update of PEM ice evolution tendencies. |
|---|
| 8 | ! |
|---|
| 9 | ! AUTHORS & DATE |
|---|
| 10 | ! R. Vandemeulebrouck |
|---|
| 11 | ! L. Lange |
|---|
| 12 | ! JB Clement, 2023-2025 |
|---|
| 13 | ! |
|---|
| 14 | ! NOTES |
|---|
| 15 | ! |
|---|
| 16 | !----------------------------------------------------------------------- |
|---|
| 17 | |
|---|
| 18 | ! DEPENDENCIES |
|---|
| 19 | ! ------------ |
|---|
| 20 | use numerics, only: dp, di, minieps, tol |
|---|
| 21 | |
|---|
| 22 | ! DECLARATION |
|---|
| 23 | ! ----------- |
|---|
| 24 | implicit none |
|---|
| 25 | |
|---|
| 26 | contains |
|---|
| 27 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 28 | |
|---|
| 29 | !======================================================================= |
|---|
| 30 | SUBROUTINE compute_tend(min_ice,d_ice) |
|---|
| 31 | !----------------------------------------------------------------------- |
|---|
| 32 | ! NAME |
|---|
| 33 | ! compute_tend |
|---|
| 34 | ! |
|---|
| 35 | ! DESCRIPTION |
|---|
| 36 | ! Compute initial ice evolution tendencies from PCM data. |
|---|
| 37 | ! |
|---|
| 38 | ! AUTHORS & DATE |
|---|
| 39 | ! R. Vandemeulebrouck |
|---|
| 40 | ! L. Lange |
|---|
| 41 | ! JB Clement, 2023-2025 |
|---|
| 42 | ! |
|---|
| 43 | ! NOTES |
|---|
| 44 | ! Based on minima of ice at each point for the PCM years. |
|---|
| 45 | !----------------------------------------------------------------------- |
|---|
| 46 | |
|---|
| 47 | ! DECLARATION |
|---|
| 48 | ! ----------- |
|---|
| 49 | implicit none |
|---|
| 50 | |
|---|
| 51 | ! ARGUMENTS |
|---|
| 52 | ! --------- |
|---|
| 53 | real(dp), dimension(:,:,:), intent(in) :: min_ice |
|---|
| 54 | real(dp), dimension(:,:), intent(out) :: d_ice |
|---|
| 55 | |
|---|
| 56 | ! CODE |
|---|
| 57 | ! ---- |
|---|
| 58 | ! We compute the difference |
|---|
| 59 | d_ice(:,:) = min_ice(:,:,2) - min_ice(:,:,1) |
|---|
| 60 | |
|---|
| 61 | ! If the difference is too small, then there is no evolution |
|---|
| 62 | where (abs(d_ice) < minieps) d_ice = 0._dp |
|---|
| 63 | |
|---|
| 64 | ! If the minimum over the last year is 0, then we have no perennial ice |
|---|
| 65 | where (abs(min_ice(:,:,2)) < minieps) d_ice = 0._dp |
|---|
| 66 | |
|---|
| 67 | END SUBROUTINE compute_tend |
|---|
| 68 | !======================================================================= |
|---|
| 69 | |
|---|
| 70 | !======================================================================= |
|---|
| 71 | SUBROUTINE evolve_tend_co2(d_co2ice_ini,co2ice,emissivity,q_co2_ts_ini,q_co2_ts,ps_PCM,ps_avg_glob_ini,ps_avg_global,d_co2ice) |
|---|
| 72 | !----------------------------------------------------------------------- |
|---|
| 73 | ! NAME |
|---|
| 74 | ! evolve_tend_co2 |
|---|
| 75 | ! |
|---|
| 76 | ! DESCRIPTION |
|---|
| 77 | ! Recompute CO2 ice tendency based on pressure and atmospheric changes. |
|---|
| 78 | ! |
|---|
| 79 | ! AUTHORS & DATE |
|---|
| 80 | ! L. Lange |
|---|
| 81 | ! JB Clement, 2023-2025 |
|---|
| 82 | ! |
|---|
| 83 | ! NOTES |
|---|
| 84 | ! Adjusts CO2 ice evolution based on Clausius-Clapeyron changes. |
|---|
| 85 | !----------------------------------------------------------------------- |
|---|
| 86 | |
|---|
| 87 | ! DEPENDENCIES |
|---|
| 88 | ! ------------ |
|---|
| 89 | use geometry, only: ngrid, nslope, nday |
|---|
| 90 | use planet, only: alpha_clap_co2, beta_clap_co2, sigmaB, Lco2, m_co2, A, B |
|---|
| 91 | use orbit, only: yr_len_sol, sol_len_s |
|---|
| 92 | use display, only: print_msg |
|---|
| 93 | use utility, only: real2str |
|---|
| 94 | |
|---|
| 95 | ! DECLARATION |
|---|
| 96 | ! ----------- |
|---|
| 97 | implicit none |
|---|
| 98 | |
|---|
| 99 | ! ARGUMENTS |
|---|
| 100 | ! --------- |
|---|
| 101 | real(dp), dimension(:,:), intent(in) :: q_co2_ts, q_co2_ts_ini, ps_PCM, d_co2ice_ini, co2ice, emissivity |
|---|
| 102 | real(dp), intent(in) :: ps_avg_global, ps_avg_glob_ini |
|---|
| 103 | real(dp), dimension(:,:), intent(inout) :: d_co2ice |
|---|
| 104 | |
|---|
| 105 | ! LOCAL VARIABLES |
|---|
| 106 | ! --------------- |
|---|
| 107 | integer(di) :: i, iday, islope |
|---|
| 108 | real(dp) :: coef, avg |
|---|
| 109 | real(dp), dimension(ngrid,nday) :: vmr_co2_ts, vmr_co2_ts_ini |
|---|
| 110 | |
|---|
| 111 | ! CODE |
|---|
| 112 | ! ---- |
|---|
| 113 | call print_msg("> Evolving the CO2 ice tendency according to the new pressure") |
|---|
| 114 | call print_msg('Old CO2 ice tendencies [kg/m2/yr] (min|max): '//real2str(minval(d_co2ice))//' | '//real2str(maxval(d_co2ice))) |
|---|
| 115 | |
|---|
| 116 | ! Compute volume mixing ratio |
|---|
| 117 | vmr_co2_ts(:,:) = q_co2_ts(:,:)/(A*q_co2_ts(:,:) + B)/m_co2 |
|---|
| 118 | vmr_co2_ts_ini(:,:) = q_co2_ts_ini(:,:)/(A*q_co2_ts_ini(:,:) + B)/m_co2 |
|---|
| 119 | |
|---|
| 120 | ! Recompute the tendency |
|---|
| 121 | do i = 1,ngrid |
|---|
| 122 | do islope = 1,nslope |
|---|
| 123 | coef = yr_len_sol*sol_len_s*emissivity(i,islope)*sigmaB/Lco2 |
|---|
| 124 | avg = 0._dp |
|---|
| 125 | if (co2ice(i,islope) > 1.e-4_dp .and. abs(d_co2ice(i,islope)) > 1.e-5_dp) then |
|---|
| 126 | do iday = 1,nday |
|---|
| 127 | avg = avg + (beta_clap_co2/(alpha_clap_co2 - log(q_co2_ts_ini(i,iday)*ps_PCM(i,iday)/100._dp)))**4 & |
|---|
| 128 | - (beta_clap_co2/(alpha_clap_co2 - log(vmr_co2_ts(i,iday)*ps_PCM(i,iday)*(ps_avg_global/ps_avg_glob_ini)/100._dp)))**4 |
|---|
| 129 | end do |
|---|
| 130 | if (avg < 1.e-4_dp) avg = 0._dp |
|---|
| 131 | d_co2ice(i,islope) = d_co2ice_ini(i,islope) - coef*avg/nday |
|---|
| 132 | end if |
|---|
| 133 | end do |
|---|
| 134 | end do |
|---|
| 135 | call print_msg('New CO2 ice tendencies [kg/m2/yr] (min|max): '//real2str( minval(d_co2ice))//' | '//real2str(maxval(d_co2ice))) |
|---|
| 136 | |
|---|
| 137 | END SUBROUTINE evolve_tend_co2 |
|---|
| 138 | !======================================================================= |
|---|
| 139 | |
|---|
| 140 | !======================================================================= |
|---|
| 141 | SUBROUTINE evolve_tend_h2o(h2oice_depth_old,h2oice_depth_new,tsurf,tsoil_ts_old,tsoil_ts_new,d_h2oice) |
|---|
| 142 | !----------------------------------------------------------------------- |
|---|
| 143 | ! NAME |
|---|
| 144 | ! evolve_tend_h2o |
|---|
| 145 | ! |
|---|
| 146 | ! DESCRIPTION |
|---|
| 147 | ! Recompute H2O ice tendency based on soil depth and temperature changes. |
|---|
| 148 | ! |
|---|
| 149 | ! AUTHORS & DATE |
|---|
| 150 | ! JB Clement, 2025 (following E. Vos's work) |
|---|
| 151 | ! |
|---|
| 152 | ! NOTES |
|---|
| 153 | ! |
|---|
| 154 | !----------------------------------------------------------------------- |
|---|
| 155 | |
|---|
| 156 | ! DEPENDENCIES |
|---|
| 157 | ! ------------ |
|---|
| 158 | use soil_temp, only: itp_tsoil |
|---|
| 159 | use subsurface_ice, only: psv |
|---|
| 160 | |
|---|
| 161 | ! DECLARATION |
|---|
| 162 | ! ----------- |
|---|
| 163 | implicit none |
|---|
| 164 | |
|---|
| 165 | ! ARGUMENTS |
|---|
| 166 | ! --------- |
|---|
| 167 | real(dp), intent(in) :: h2oice_depth_old ! Old H2O ice depth |
|---|
| 168 | real(dp), intent(in) :: h2oice_depth_new ! New H2O ice depth |
|---|
| 169 | real(dp), intent(in) :: tsurf ! Surface temperature |
|---|
| 170 | real(dp), dimension(:,:), intent(in) :: tsoil_ts_old ! Old soil temperature time series |
|---|
| 171 | real(dp), dimension(:,:), intent(in) :: tsoil_ts_new ! New soil temperature time series |
|---|
| 172 | real(dp), intent(inout) :: d_h2oice ! Evolution of perennial ice |
|---|
| 173 | |
|---|
| 174 | ! LOCAL VARIABLES |
|---|
| 175 | ! --------------- |
|---|
| 176 | real(dp) :: Rz_old, Rz_new, R_dec, hum_dec, psv_max_old, psv_max_new |
|---|
| 177 | integer(di) :: iday |
|---|
| 178 | real(dp), parameter :: coef_diff = 4.e-4 ! Diffusion coefficient |
|---|
| 179 | real(dp), parameter :: zcdv = 0.0325 ! Drag coefficient |
|---|
| 180 | |
|---|
| 181 | ! CODE |
|---|
| 182 | ! ---- |
|---|
| 183 | ! Higher resistance due to growing lag layer (higher depth) |
|---|
| 184 | Rz_old = h2oice_depth_old*zcdv/coef_diff ! Old resistance from PCM |
|---|
| 185 | Rz_new = h2oice_depth_new*zcdv/coef_diff ! New resistance based on new depth |
|---|
| 186 | R_dec = Rz_old/Rz_new ! Decrease because of resistance |
|---|
| 187 | |
|---|
| 188 | ! The maxmimum of the daily averages over one year for the saturation vapor pressure at the ice table location |
|---|
| 189 | psv_max_old = 0._dp |
|---|
| 190 | psv_max_new = 0._dp |
|---|
| 191 | do iday = 1,size(tsoil_ts_old,2) |
|---|
| 192 | psv_max_old = max(psv_max_old,psv(itp_tsoil(tsoil_ts_old(:,iday),tsurf,h2oice_depth_old))) |
|---|
| 193 | psv_max_new = max(psv_max_new,psv(itp_tsoil(tsoil_ts_new(:,iday),tsurf,h2oice_depth_new))) |
|---|
| 194 | end do |
|---|
| 195 | |
|---|
| 196 | ! Lower humidity due to growing lag layer (higher depth) |
|---|
| 197 | if (abs(psv_max_old) < tol) then |
|---|
| 198 | hum_dec = 1._dp |
|---|
| 199 | else |
|---|
| 200 | hum_dec = psv_max_new/psv_max_old ! Decrease because of lower water vapor pressure at the new depth |
|---|
| 201 | end if |
|---|
| 202 | |
|---|
| 203 | ! Flux correction (decrease) |
|---|
| 204 | d_h2oice = d_h2oice*R_dec*hum_dec |
|---|
| 205 | |
|---|
| 206 | END SUBROUTINE evolve_tend_h2o |
|---|
| 207 | !======================================================================= |
|---|
| 208 | |
|---|
| 209 | END MODULE tendencies |
|---|