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

Last change on this file since 4065 was 4065, checked in by jbclement, 8 days 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: 15.7 KB
Line 
1MODULE stopping_crit
2!-----------------------------------------------------------------------
3! NAME
4!     stopping_crit
5!
6! DESCRIPTION
7!     Stopping criteria for PEM.
8!
9! AUTHORS & DATE
10!     JB Clement, 2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPENDENCIES
17! ------------
18use numerics, only: dp, qp, di, k4, minieps_qp
19
20! DECLARATION
21! -----------
22implicit none
23
24! PARAMETERS
25! ----------
26real(dp), protected :: h2oice_crit ! Change of the surface of h2o ice sublimating before stopping the PEM
27real(dp), protected :: co2ice_crit ! Change of the surface of co2 ice sublimating before stopping the PEM
28real(dp), protected :: ps_crit     ! Change of averaged surface pressure before stopping the PEM
29
30! VARIABLES
31! ---------
32type :: stopFlags
33    logical(k4) :: surface_h2oice_change  = .false.
34    logical(k4) :: zero_dh2oice           = .false.
35    logical(k4) :: surface_co2ice_change  = .false.
36    logical(k4) :: pressure_change        = .false.
37    logical(k4) :: nmax_yr_run_reached    = .false.
38    logical(k4) :: nmax_yr_runorb_reached = .false.
39    logical(k4) :: nmax_yr_sim_reached    = .false.
40    logical(k4) :: time_limit_reached     = .false.
41contains
42    procedure :: is_any_set
43    procedure :: stop_code
44    procedure :: stop_message
45end type
46
47contains
48!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
49
50!=======================================================================
51SUBROUTINE set_stopping_crit_config(h2oice_crit_in,co2ice_crit_in,ps_crit_in)
52!-----------------------------------------------------------------------
53! NAME
54!     set_stopping_crit_config
55!
56! DESCRIPTION
57!     Setter for 'stopping_crit' configuration parameters.
58!
59! AUTHORS & DATE
60!     JB Clement, 02/2026
61!
62! NOTES
63!
64!-----------------------------------------------------------------------
65
66! DEPENDENCIES
67! ------------
68use utility, only: real2str
69use display, only: print_msg
70
71! DECLARATION
72! -----------
73implicit none
74
75! ARGUMENTS
76! ---------
77real(dp), intent(in) :: h2oice_crit_in, co2ice_crit_in, ps_crit_in
78
79! CODE
80! ----
81h2oice_crit = h2oice_crit_in
82co2ice_crit = co2ice_crit_in
83ps_crit = ps_crit_in
84call print_msg('h2oice_crit = '//real2str(h2oice_crit))
85call print_msg('co2ice_crit = '//real2str(co2ice_crit))
86call print_msg('ps_crit     = '//real2str(ps_crit))
87
88END SUBROUTINE set_stopping_crit_config
89!=======================================================================
90
91!=======================================================================
92FUNCTION is_any_set(this) RESULT(stopPEM)
93!-----------------------------------------------------------------------
94! NAME
95!     is_any_set
96!
97! DESCRIPTION
98!     Detect if any stopping flag is set to true.
99!
100! AUTHORS & DATE
101!     JB Clement, 2025
102!
103! NOTES
104!
105!-----------------------------------------------------------------------
106
107! DECLARATION
108! -----------
109implicit none
110
111! ARGUMENTS
112! ---------
113class(stopFlags), intent(in) :: this    ! Stopping flags object
114logical(k4)                  :: stopPEM ! True if any flag is set
115
116! CODE
117! ----
118stopPEM = this%surface_h2oice_change  .or. &
119          this%zero_dh2oice           .or. &
120          this%surface_co2ice_change  .or. &
121          this%pressure_change        .or. &
122          this%nmax_yr_run_reached    .or. &
123          this%nmax_yr_runorb_reached .or. &
124          this%nmax_yr_sim_reached    .or. &
125          this%time_limit_reached
126
127END FUNCTION is_any_set
128!=======================================================================
129
130!=======================================================================
131FUNCTION stop_code(this) result(code)
132!-----------------------------------------------------------------------
133! NAME
134!     stop_code
135!
136! DESCRIPTION
137!     Return digit code corresponding to the active stopping flag.
138!
139! AUTHORS & DATE
140!     JB Clement, 2025
141!
142! NOTES
143!
144!-----------------------------------------------------------------------
145
146! DECLARATION
147! -----------
148implicit none
149
150! ARGUMENTS
151! ---------
152class(StopFlags), intent(in) :: this ! Stopping flags object
153integer(di)                  :: code ! Code corresponding to active flag (0 if none)
154
155! CODE
156! ----
157code = 0
158if     (this%surface_h2oice_change)  then; code = 1
159elseif (this%zero_dh2oice)           then; code = 2
160elseif (this%surface_co2ice_change)  then; code = 3
161elseif (this%pressure_change)        then; code = 4
162elseif (this%nmax_yr_run_reached)    then; code = 5
163elseif (this%nmax_yr_runorb_reached) then; code = 6
164elseif (this%nmax_yr_sim_reached)    then; code = 7
165elseif (this%time_limit_reached)     then; code = 8
166end if
167
168END FUNCTION stop_code
169!=======================================================================
170
171!=======================================================================
172FUNCTION stop_message(this) result(msg)
173!-----------------------------------------------------------------------
174! NAME
175!     stop_message
176!
177! DESCRIPTION
178!     Return descriptive message corresponding to the active stopping flag.
179!
180! AUTHORS & DATE
181!     JB Clement, 2025
182!
183! NOTES
184!
185!-----------------------------------------------------------------------
186
187! DECLARATION
188! -----------
189implicit none
190
191! ARGUMENTS
192! ---------
193class(StopFlags), intent(in) :: this ! Stopping flags object
194character(:), allocatable    :: msg  ! Descriptive message for active flag
195
196! CODE
197! ----
198if     (this%surface_h2oice_change)  then; msg = "**** STOPPING because surface of h2o ice has changed too much (see message above)."
199elseif (this%zero_dh2oice)           then; msg = "**** STOPPING because tendencies on h2o ice = 0 (see message above)."
200elseif (this%surface_co2ice_change)  then; msg = "**** STOPPING because surface of co2 ice has changed too much (see message above)."
201elseif (this%pressure_change)        then; msg = "**** STOPPING because surface global pressure has changed too much (see message above)."
202elseif (this%nmax_yr_run_reached)    then; msg = "**** STOPPING because maximum number of years is reached."
203elseif (this%nmax_yr_runorb_reached) then; msg = "**** STOPPING because maximum number of years due to orbital parameters is reached."
204elseif (this%nmax_yr_sim_reached)    then; msg = "**** STOPPING because the number of years to be simulated is reached."
205elseif (this%time_limit_reached)     then; msg = "**** STOPPING because the job time limit is going to be reached."
206end if
207
208END FUNCTION stop_message
209!=======================================================================
210
211!=======================================================================
212SUBROUTINE stopping_crit_h2oice(h2oice_sublim_coverage_ini,h2o_ice,d_h2oice,stopcrit)
213
214!-----------------------------------------------------------------------
215! NAME
216!     stopping_crit_h2oice
217!
218! DESCRIPTION
219!     Check if H2O ice surface has changed too much.
220!
221! AUTHORS & DATE
222!     R. Vandemeulebrouck
223!     L. Lange, 2023
224!     JB Clement, 2023-2025
225!
226! NOTES
227!
228!-----------------------------------------------------------------------
229
230! DEPENDENCIES
231! ------------
232use slopes,   only: subslope_dist, def_slope_mean
233use geometry, only: ngrid, nslope, cell_area
234use maths,    only: pi
235use display,  only: print_msg
236use utility,  only: real2str
237
238! DECLARATION
239! -----------
240implicit none
241
242! ARGUMENTS
243! ---------
244real(dp), dimension(:,:), intent(in)    :: h2o_ice, d_h2oice
245real(dp),                 intent(in)    :: h2oice_sublim_coverage_ini
246type(stopFlags),          intent(inout) :: stopcrit
247
248! LOCAL VARIABLES
249! ---------------
250integer(di) :: i, islope
251real(dp)    :: h2oice_now_surf
252
253! CODE
254! ----
255! Check to escape the subroutine (not relevant in 1D or if the PEM should stop already)
256if (ngrid == 1 .or. stopcrit%stop_code() > 0) return
257
258! Computation of the present surface of sublimating H2O ice
259h2oice_now_surf = 0._dp
260do i = 1,ngrid
261    do islope = 1,nslope
262        if (h2o_ice(i,islope) > 0._dp .and. d_h2oice(i,islope) < 0._dp) h2oice_now_surf = h2oice_now_surf + cell_area(i)*subslope_dist(i,islope)/cos(pi*def_slope_mean(islope)/180._dp)
263    end do
264end do
265
266! Check of the criterion
267if (h2oice_now_surf < h2oice_sublim_coverage_ini*(1._dp - h2oice_crit) .or. h2oice_now_surf > h2oice_sublim_coverage_ini*(1. + h2oice_crit)) then
268    stopcrit%surface_h2oice_change = .true.
269    call print_msg("Surface of H2O ice sublimating (ini|now) [m2] = "//real2str(h2oice_sublim_coverage_ini)//' | '//real2str(h2oice_now_surf))
270    call print_msg("Change (accepted|current) [%] = +/- "//real2str(h2oice_crit*100._dp)//' | '//real2str(100._dp*(h2oice_now_surf - h2oice_sublim_coverage_ini)/h2oice_sublim_coverage_ini))
271end if
272
273END SUBROUTINE stopping_crit_h2oice
274!=======================================================================
275
276!=======================================================================
277SUBROUTINE stopping_crit_co2ice(co2ice_sublim_coverage_ini,co2_ice,d_co2ice,stopcrit)
278!-----------------------------------------------------------------------
279! NAME
280!     stopping_crit_co2ice
281!
282! DESCRIPTION
283!     Check if CO2 ice surface has changed too much.
284!
285! AUTHORS & DATE
286!     R. Vandemeulebrouck
287!     L. Lange, 2023
288!     JB Clement, 2023-2025
289!
290! NOTES
291!
292!-----------------------------------------------------------------------
293
294! DEPENDENCIES
295! ------------
296use slopes,   only: subslope_dist, def_slope_mean
297use geometry, only: ngrid, nslope, cell_area
298use maths,    only: pi
299use display,  only: print_msg
300use utility,  only: real2str
301
302! DECLARATION
303! -----------
304implicit none
305
306! ARGUMENTS
307! ---------
308real(dp), dimension(:,:), intent(in)    :: co2_ice, d_co2ice
309real(dp),                 intent(in)    :: co2ice_sublim_coverage_ini
310type(stopFlags),          intent(inout) :: stopcrit
311
312! LOCAL VARIABLES
313! ---------------
314integer(di) :: i, islope
315real(dp)    :: co2ice_now_surf ! Current surface of co2 ice
316
317! CODE
318! ----
319! Check to escape the subroutine (not relevant in 1D or if the PEM should stop already)
320if (ngrid == 1 .or. stopcrit%stop_code() > 0) return
321
322! Computation of the present surface of co2 ice still sublimating
323co2ice_now_surf = 0.
324do i = 1,ngrid
325    do islope = 1,nslope
326        if (co2_ice(i,islope) > 0._dp .and. d_co2ice(i,islope) < 0._dp) co2ice_now_surf = co2ice_now_surf + cell_area(i)*subslope_dist(i,islope)/cos(pi*def_slope_mean(islope)/180._dp)
327    end do
328end do
329
330! Check of the criterion
331if (co2ice_now_surf < co2ice_sublim_coverage_ini*(1. - co2ice_crit) .or. co2ice_now_surf > co2ice_sublim_coverage_ini*(1._dp + co2ice_crit)) then
332    stopcrit%surface_co2ice_change = .true.
333    call print_msg("Initial surface of co2 ice sublimating (ini|now) [m2] = "//real2str(co2ice_sublim_coverage_ini)//' | '//real2str(co2ice_now_surf))
334    call print_msg("Change (accepted|current) [%] = +/- "//real2str(co2ice_crit*100._dp)//' | '//real2str(100._dp*(co2ice_now_surf - co2ice_sublim_coverage_ini)/co2ice_sublim_coverage_ini))
335end if
336
337END SUBROUTINE stopping_crit_co2ice
338!=======================================================================
339
340!=======================================================================
341SUBROUTINE stopping_crit_pressure(ps_avg_glob_ini,ps_avg_global,stopcrit)
342!-----------------------------------------------------------------------
343! NAME
344!     stopping_crit_pressure
345!
346! DESCRIPTION
347!     Check if pressure has changed too much.
348!
349! AUTHORS & DATE
350!     R. Vandemeulebrouck
351!     L. Lange, 2023
352!     JB Clement, 2023-2025
353!
354! NOTES
355!
356!-----------------------------------------------------------------------
357
358! DEPENDENCIES
359! ------------
360use display,  only: print_msg
361use utility,  only: real2str
362
363! DECLARATION
364! -----------
365implicit none
366
367! ARGUMENTS
368! ---------
369real(dp),        intent(in)    :: ps_avg_glob_ini, ps_avg_global
370type(stopFlags), intent(inout) :: stopcrit
371
372! CODE
373! ----
374! Check to escape the subroutine (not relevant if the PEM should stop already)
375if (stopcrit%stop_code() > 0) return
376
377if (ps_avg_global < ps_avg_glob_ini*(1._dp - ps_crit) .or. ps_avg_global > ps_avg_glob_ini*(1._dp + ps_crit)) then
378    stopcrit%pressure_change = .true.
379    call print_msg("Initial global pressure (ini|now) [Pa] = "//real2str(ps_avg_glob_ini)//' | '//real2str(ps_avg_global))
380    call print_msg("Change (accepted|current) [%] = +/- "//real2str(ps_crit*100._dp)//' | '//real2str(100._dp*(ps_avg_global - ps_avg_glob_ini)/ps_avg_glob_ini))
381end if
382
383END SUBROUTINE stopping_crit_pressure
384!=======================================================================
385
386!=======================================================================
387SUBROUTINE stopping_crit_h2o(delta_h2o_ads,delta_icetable,h2o_ice,d_h2oice,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,stopcrit)
388!-----------------------------------------------------------------------
389! NAME
390!     stopping_crit_h2o
391!
392! DESCRIPTION
393!     Check if PEM oscillates infinitely with H2O (no net exchange).
394!
395! AUTHORS & DATE
396!     L. Lange, 2024
397!     JB Clement, 2025
398!
399! NOTES
400!
401!-----------------------------------------------------------------------
402
403! DEPENDENCIES
404! ------------
405use geometry,  only: ngrid, nslope, cell_area
406use evolution, only: dt
407use slopes,    only: subslope_dist, def_slope_mean
408use maths,     only: pi
409use display,   only: print_msg
410use utility,   only: real2str
411
412! DECLARATION
413! -----------
414implicit none
415
416! ARGUMENTS
417! ---------
418real(dp), dimension(:),   intent(in)    :: delta_h2o_ads            ! Mass of H2O adsorbed/desorbded in the soil (kg/m^2)
419real(dp), dimension(:),   intent(in)    :: delta_icetable ! Mass of H2O that condensed/sublimated at the ice table
420real(dp), dimension(:,:), intent(in)    :: h2o_ice, d_h2oice
421real(qp),                 intent(out)   :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm
422type(stopFlags),          intent(inout) :: stopcrit
423
424! LOCAL VARIABLES
425! ---------------
426integer(di) :: i, islope
427
428! CODE
429! ----
430! Check to escape the subroutine (not relevant in 1D or if the PEM should stop already)
431if (ngrid == 1 .or. stopcrit%stop_code() > 0) then
432    S_atm_2_h2o = 1._qp
433    S_h2o_2_atm = S_atm_2_h2o
434    S_atm_2_h2oice = 1._qp
435    S_h2oice_2_atm = S_atm_2_h2oice
436    return
437end if
438
439! We compute the amount of water going out from the atmosphere (S_atm_2_h2o) and going into the atmophere (S_h2o_2_atm)
440S_atm_2_h2o = 0._qp
441S_h2o_2_atm = 0._qp
442S_atm_2_h2oice = 0._qp
443S_h2oice_2_atm = 0._qp
444do i = 1,ngrid
445    if (delta_h2o_ads(i) > 0._dp) then
446        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_ads(i)*cell_area(i)
447    else
448        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_ads(i)*cell_area(i)
449    end if
450    if (delta_icetable(i) > 0._dp) then
451        S_atm_2_h2o = S_atm_2_h2o + delta_icetable(i)*cell_area(i)
452    else
453        S_h2o_2_atm = S_h2o_2_atm + delta_icetable(i)*cell_area(i)
454    end if
455    do islope = 1,nslope
456        if (d_h2oice(i,islope) > 0._dp) then
457            S_atm_2_h2o = S_atm_2_h2o + d_h2oice(i,islope)*dt*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180._dp)
458            S_atm_2_h2oice = S_atm_2_h2oice + d_h2oice(i,islope)*dt*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180._dp)
459        else if (d_h2oice(i,islope) < 0._dp .and. h2o_ice(i,islope) > 0._dp) then
460            S_h2o_2_atm = S_h2o_2_atm - d_h2oice(i,islope)*dt*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180._dp)
461            S_h2oice_2_atm = S_h2oice_2_atm - d_h2oice(i,islope)*dt*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180._dp)
462        end if
463    end do
464end do
465
466! Since relative atmospheric water is kept constant, we need to equate condensing reservoirs to the sublimating ones.
467! It is not possible if one of them is missing.
468if ((abs(S_atm_2_h2oice) < minieps_qp .or. abs(S_h2oice_2_atm) < minieps_qp)) then
469    call print_msg("There is no sublimating or condensing H2O ice! This can be due to the absence of H2O ice in the PCM run.")
470    call print_msg("Amount of condensing ice  [kg] = "//real2str(S_atm_2_h2oice))
471    call print_msg("Amount of sublimating ice [kg] = "//real2str(S_h2oice_2_atm))
472    stopcrit%zero_dh2oice = .true.
473end if
474
475END SUBROUTINE stopping_crit_h2o
476!=======================================================================
477
478END MODULE stopping_crit
Note: See TracBrowser for help on using the repository browser.