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

Last change on this file since 4110 was 4110, checked in by jbclement, 9 days ago

PEM:

  • Introduction of a configurable display/logging system with options 'out2term', 'out2log', 'verbosity_lvl'. All messages now use verbosity levels ('LVL_NFO', 'LVL_WRN', 'LVL_ERR' and 'LVL_DBG').
  • Code encapsulation improvements with systematic privacy/protection of module variables.
  • Addition of workflow safety checks for required executables, dependencies (e.g. 'ncdump'), input files and callphys keys.
  • Renaming of PEM starting and diagnostic files ("startevol.nc" into "startevo.nc", "diagevol.nc" into "diagevo.nc").

JBC

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