source: trunk/LMDZ.COMMON/libf/evolution/tendencies.F90 @ 4092

Last change on this file since 4092 was 4074, checked in by jbclement, 2 weeks ago

PEM:

  • Correct management of H2O ice tendency in 1D when there is not enough ice anymore.
  • Clean initialization of allocatable module arrays (especially needed when no slope)
  • One more renaming for consistency + few small updates thoughout the code.

JBC

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