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

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

PEM:
Major refactor following the previous ones (r3989 and r3991) completing the large structural reorganization and cleanup of the PEM codebase. This revision introduces newly designed modules, standardizes interfaces with explicit ini/end APIs and adds native NetCDF I/O together with explicit PCM/PEM adapters. In detail:

  • Some PEM models were corrected or improved:
    • Frost/perennial ice semantics are clarified via renaming;
    • Soil temperature remapping clarified, notably by removing the rescaling of temperature deviation;
    • Geothermal flux for the PCM is computed based on the PEM state;
  • New explicit PEM/PCM adapters ("set_*"/"build4PCM_*") to decouple PEM internal representation from PCM file layouts and reconstruct consistent fields returned to the PCM;
  • New explicit build/teardown routines that centralize allocation and initialization ordering, reducing accidental use of uninitialized data and making the model lifecycle explicit;
  • Add native read/write helpers for NetCDF that centralize all low-level NetCDF interactions with major improvements (and more simplicity) compared to legacy PEM/PCM I/O (see the modules "io_netcdf" and "output"). They support reading, creation and writing of "diagevol.nc" (renamed from "diagpem.nc") and start/restart files;
  • Provide well-focused modules ("numerics"/"maths"/"utility"/"display") to host commonly-used primitives:
    • "numerics" defines numerical types and constants for reproducibility, portability across compilers and future transitions (e.g. quadruple precision experiments);
    • "display" provides a single controlled interface for runtime messages, status output and diagnostics, avoiding direct 'print'/'write' to enable silent mode, log redirection, and MPI-safe output in the future.
    • "utility" (new module) hosts generic helpers used throughout the code (e.g. "int2str" or "real2str");
  • Add modules "clim_state_init"/"clim_state_rec" which provide robust read/write logic for "start/startfi/startpem", including 1D fallbacks, mesh conversions and dimension checks. NetCDF file creation is centralized and explicit. Restart files are now self-consistent and future-proof, requiring changes only to affected variables;
  • Add module "atmosphere" which computes pressure fields, reconstructs potential temperature and air mass. It also holds the whole logic to define sigma or hybrid coordinates for altitudes;
  • Add module "geometry" to centrilize dimensions logic and grid conversions routines (including 2 new ones "dyngrd2vect"/"vect2dyngrd");
  • Add module "slopes" to isolate slopes handling;
  • Add module "surface" to isolate surface management. Notably, albedo and emissivity are now fully reconstructed following the PCM settings;
  • Add module "allocation" to check the dimension initialization and centrilize allocation/deallocation;
  • Finalize module decomposition and renaming to consolidate domain-based modules, purpose-based routines and physics/process-based variables;
  • The main program now drives a clearer sequence of conceptual steps (initialization / reading / evolution / update / build / writing) and fails explicitly instead of silently defaulting;
  • Ice table logic is made restart-safe;
  • 'Open'/'read' intrinsic logic is made safe and automatic;
  • Improve discoverability and standardize the data handling (private vs protected vs public);
  • Apply consistent documentation/header style already introduced;
  • Update deftank scripts to reflect new names and launch wrappers;

This revision is a structural milestone aiming to be behavior-preserving where possible. It has been tested via compilation and short integration runs. However, due to extensive renames, moves, and API changes, full validation is still ongoing.
Note: the revision includes one (possibly two) easter egg hidden in the code for future archaeologists of the PEM. No physics were harmed.
JBC

File size: 6.7 KB
Line 
1MODULE 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! ------------
20use numerics, only: dp, di, minieps, tol
21
22! DECLARATION
23! -----------
24implicit none
25
26contains
27!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28
29!=======================================================================
30SUBROUTINE 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! -----------
49implicit none
50
51! ARGUMENTS
52! ---------
53real(dp), dimension(:,:,:), intent(in)  :: min_ice
54real(dp), dimension(:,:),   intent(out) :: d_ice
55
56! CODE
57! ----
58! We compute the difference
59d_ice(:,:) = min_ice(:,:,2) - min_ice(:,:,1)
60
61! If the difference is too small, then there is no evolution
62where (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
65where (abs(min_ice(:,:,2)) < minieps) d_ice = 0._dp
66
67END SUBROUTINE compute_tend
68!=======================================================================
69
70!=======================================================================
71SUBROUTINE 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! ------------
89use geometry, only: ngrid, nslope, nday
90use planet,   only: alpha_clap_co2, beta_clap_co2, sigmaB, Lco2, m_co2, A, B
91use orbit,    only: yr_len_sol, sol_len_s
92use display,  only: print_msg
93use utility,  only: real2str
94
95! DECLARATION
96! -----------
97implicit none
98
99! ARGUMENTS
100! ---------
101real(dp), dimension(:,:), intent(in)    :: q_co2_ts, q_co2_ts_ini, ps_PCM, d_co2ice_ini, co2ice, emissivity
102real(dp),                 intent(in)    :: ps_avg_global, ps_avg_glob_ini
103real(dp), dimension(:,:), intent(inout) :: d_co2ice
104
105! LOCAL VARIABLES
106! ---------------
107integer(di)                     :: i, iday, islope
108real(dp)                        :: coef, avg
109real(dp), dimension(ngrid,nday) :: vmr_co2_ts, vmr_co2_ts_ini
110
111! CODE
112! ----
113call print_msg("> Evolving the CO2 ice tendency according to the new pressure")
114call 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
117vmr_co2_ts(:,:) = q_co2_ts(:,:)/(A*q_co2_ts(:,:) + B)/m_co2
118vmr_co2_ts_ini(:,:) = q_co2_ts_ini(:,:)/(A*q_co2_ts_ini(:,:) + B)/m_co2
119
120! Recompute the tendency
121do 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
134end do
135call print_msg('New CO2 ice tendencies [kg/m2/yr] (min|max): '//real2str( minval(d_co2ice))//' | '//real2str(maxval(d_co2ice)))
136
137END SUBROUTINE evolve_tend_co2
138!=======================================================================
139
140!=======================================================================
141SUBROUTINE 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! ------------
158use soil_temp,      only: itp_tsoil
159use subsurface_ice, only: psv
160
161! DECLARATION
162! -----------
163implicit none
164
165! ARGUMENTS
166! ---------
167real(dp),                 intent(in)    :: h2oice_depth_old ! Old H2O ice depth
168real(dp),                 intent(in)    :: h2oice_depth_new ! New H2O ice depth
169real(dp),                 intent(in)    :: tsurf            ! Surface temperature
170real(dp), dimension(:,:), intent(in)    :: tsoil_ts_old     ! Old soil temperature time series
171real(dp), dimension(:,:), intent(in)    :: tsoil_ts_new     ! New soil temperature time series
172real(dp),                 intent(inout) :: d_h2oice         ! Evolution of perennial ice
173
174! LOCAL VARIABLES
175! ---------------
176real(dp)            :: Rz_old, Rz_new, R_dec, hum_dec, psv_max_old, psv_max_new
177integer(di)         :: iday
178real(dp), parameter :: coef_diff = 4.e-4 ! Diffusion coefficient
179real(dp), parameter :: zcdv = 0.0325     ! Drag coefficient
180
181! CODE
182! ----
183! Higher resistance due to growing lag layer (higher depth)
184Rz_old = h2oice_depth_old*zcdv/coef_diff ! Old resistance from PCM
185Rz_new = h2oice_depth_new*zcdv/coef_diff ! New resistance based on new depth
186R_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
189psv_max_old = 0._dp
190psv_max_new = 0._dp
191do 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)))
194end do
195
196! Lower humidity due to growing lag layer (higher depth)
197if (abs(psv_max_old) < tol) then
198    hum_dec = 1._dp
199else
200    hum_dec = psv_max_new/psv_max_old ! Decrease because of lower water vapor pressure at the new depth
201end if
202
203! Flux correction (decrease)
204d_h2oice = d_h2oice*R_dec*hum_dec
205
206END SUBROUTINE evolve_tend_h2o
207!=======================================================================
208
209END MODULE tendencies
Note: See TracBrowser for help on using the repository browser.