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

Last change on this file since 3991 was 3991, checked in by jbclement, 6 days ago

PEM:
Apply documentation template everywhere: standardized headers format with short description, separators between functions/subroutines, normalized code sections, aligned dependencies/arguments/variables declaration.
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! DECLARATION
17! -----------
18implicit none
19
20! MODULE VARIABLES
21! ----------------
22real :: h2oice_crit ! Percentage of change of the surface of h2o ice sublimating before stopping the PEM
23real :: co2ice_crit ! Percentage of change of the surface of co2 ice sublimating before stopping the PEM
24real :: ps_crit     ! Percentage of change of averaged surface pressure before stopping the PEM
25
26type :: stopFlags
27    logical :: surface_h2oice_change = .false.
28    logical :: zero_dh2oice          = .false.
29    logical :: surface_co2ice_change = .false.
30    logical :: pressure_change       = .false.
31    logical :: nyears_max_reached    = .false.
32    logical :: nyears_maxorb_reached = .false.
33    logical :: nyears_tot_reached    = .false.
34    logical :: time_limit_reached    = .false.
35contains
36    procedure :: is_any_set
37    procedure :: stop_code
38    procedure :: stop_message
39end type
40
41contains
42!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
43
44!=======================================================================
45FUNCTION is_any_set(this) RESULT(stopPEM)
46!-----------------------------------------------------------------------
47! NAME
48!     is_any_set
49!
50! DESCRIPTION
51!     Detect if any stopping flag is set to true.
52!
53! AUTHORS & DATE
54!     JB Clement, 2025
55!
56! NOTES
57!
58!-----------------------------------------------------------------------
59
60! DECLARATION
61! -----------
62implicit none
63
64! ARGUMENTS
65! ---------
66class(stopFlags), intent(in) :: this    ! Stopping flags object
67logical                      :: stopPEM ! True if any flag is set
68
69! CODE
70! ----
71stopPEM = this%surface_h2oice_change .or. &
72          this%zero_dh2oice          .or. &
73          this%surface_co2ice_change .or. &
74          this%pressure_change       .or. &
75          this%nyears_max_reached    .or. &
76          this%nyears_maxorb_reached .or. &
77          this%nyears_tot_reached    .or. &
78          this%time_limit_reached
79
80END FUNCTION is_any_set
81!=======================================================================
82
83!=======================================================================
84FUNCTION stop_code(this) result(code)
85!-----------------------------------------------------------------------
86! NAME
87!     stop_code
88!
89! DESCRIPTION
90!     Return digit code corresponding to the active stopping flag.
91!
92! AUTHORS & DATE
93!     JB Clement, 2025
94!
95! NOTES
96!
97!-----------------------------------------------------------------------
98
99! DECLARATION
100! -----------
101implicit none
102
103! ARGUMENTS
104! ---------
105class(StopFlags), intent(in) :: this ! Stopping flags object
106integer                      :: code ! Code corresponding to active flag (0 if none)
107
108! CODE
109! ----
110code = 0
111if     (this%surface_h2oice_change) then; code = 1
112elseif (this%zero_dh2oice)          then; code = 2
113elseif (this%surface_co2ice_change) then; code = 3
114elseif (this%pressure_change)       then; code = 4
115elseif (this%nyears_max_reached)    then; code = 5
116elseif (this%nyears_maxorb_reached) then; code = 6
117elseif (this%nyears_tot_reached)    then; code = 7
118elseif (this%time_limit_reached)    then; code = 8
119endif
120
121END FUNCTION stop_code
122!=======================================================================
123
124!=======================================================================
125FUNCTION stop_message(this) result(msg)
126!-----------------------------------------------------------------------
127! NAME
128!     stop_message
129!
130! DESCRIPTION
131!     Return descriptive message corresponding to the active stopping flag.
132!
133! AUTHORS & DATE
134!     JB Clement, 2025
135!
136! NOTES
137!
138!-----------------------------------------------------------------------
139
140! DECLARATION
141! -----------
142implicit none
143
144! ARGUMENTS
145! ---------
146class(StopFlags), intent(in) :: this ! Stopping flags object
147character(:), allocatable    :: msg  ! Descriptive message for active flag
148
149! CODE
150! ----
151if     (this%surface_h2oice_change) then; msg = " **** STOPPING because surface of h2o ice has changed too much (see message above)."
152elseif (this%zero_dh2oice)          then; msg = " **** STOPPING because tendencies on h2o ice = 0 (see message above)."
153elseif (this%surface_co2ice_change) then; msg = " **** STOPPING because surface of co2 ice has changed too much (see message above)."
154elseif (this%pressure_change)       then; msg = " **** STOPPING because surface global pressure has changed too much (see message above)."
155elseif (this%nyears_max_reached)    then; msg = " **** STOPPING because maximum number of years is reached."
156elseif (this%nyears_maxorb_reached) then; msg = " **** STOPPING because maximum number of years due to orbital parameters is reached."
157elseif (this%nyears_tot_reached)    then; msg = " **** STOPPING because the number of years to be simulated is reached."
158elseif (this%time_limit_reached)    then; msg = " **** STOPPING because the job time limit is going to be reached."
159endif
160
161END FUNCTION stop_message
162!=======================================================================
163
164!=======================================================================
165SUBROUTINE stopping_crit_h2o_ice(cell_area,h2oice_ini_surf,is_h2oice_sublim_ini,h2o_ice,stopCrit,ngrid)
166
167!-----------------------------------------------------------------------
168! NAME
169!     stopping_crit_h2o_ice
170!
171! DESCRIPTION
172!     Check if H2O ice surface change criterion to stop PEM is reached.
173!
174! AUTHORS & DATE
175!     R. Vandemeulebrouck
176!     L. Lange, 2023
177!     JB Clement, 2023-2025
178!
179! NOTES
180!
181!-----------------------------------------------------------------------
182
183! DEPENDENCIES
184! ------------
185use comslope_mod,  only: subslope_dist, nslope
186
187! DECLARATION
188! -----------
189implicit none
190
191! ARGUMENTS
192! ---------
193integer,                          intent(in)    :: ngrid                ! # of physical grid points
194real,    dimension(ngrid),        intent(in)    :: cell_area            ! Area of the cells
195real,    dimension(ngrid,nslope), intent(in)    :: h2o_ice              ! Actual density of h2o ice
196real,                             intent(in)    :: h2oice_ini_surf      ! Initial surface of sublimating h2o ice
197logical, dimension(ngrid,nslope), intent(in)    :: is_h2oice_sublim_ini ! Grid points where h2o ice was initially sublimating
198type(stopFlags),                  intent(inout) :: stopCrit             ! Stopping criterion
199
200! LOCAL VARIABLES
201! ---------------
202integer :: i, islope       ! Loop
203real    :: h2oice_now_surf ! Current surface of h2o ice
204
205! CODE
206! ----
207if (stopCrit%stop_code() > 0) return
208
209! Computation of the present surface of h2o ice still sublimating
210h2oice_now_surf = 0.
211do i = 1,ngrid
212    do islope = 1,nslope
213        if (is_h2oice_sublim_ini(i,islope) .and. h2o_ice(i,islope) > 0.) h2oice_now_surf = h2oice_now_surf + cell_area(i)*subslope_dist(i,islope)
214    enddo
215enddo
216
217! Check of the criterion
218if (h2oice_now_surf < h2oice_ini_surf*(1. - h2oice_crit)) then
219    stopCrit%surface_h2oice_change = .true.
220    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
221    write(*,*) "h2oice_now_surf < h2oice_ini_surf*(1. - h2oice_crit)", h2oice_now_surf < h2oice_ini_surf*(1. - h2oice_crit)
222    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
223    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
224    write(*,*) "Percentage of change accepted =", h2oice_crit*100
225else if (h2oice_now_surf > h2oice_ini_surf*(1. + h2oice_crit)) then
226    stopCrit%surface_h2oice_change = .true.
227    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
228    write(*,*) "h2oice_now_surf > h2oice_ini_surf*(1. + h2oice_crit)", h2oice_now_surf > h2oice_ini_surf*(1. + h2oice_crit)
229    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
230    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
231    write(*,*) "Percentage of change accepted =", h2oice_crit*100
232endif
233
234END SUBROUTINE stopping_crit_h2o_ice
235!=======================================================================
236
237!=======================================================================
238SUBROUTINE stopping_crit_co2(cell_area,co2ice_sublim_surf_ini,is_co2ice_sublim_ini,co2_ice,stopCrit,ngrid,ps_avg_global_ini,ps_avg_global,nslope)
239!-----------------------------------------------------------------------
240! NAME
241!     stopping_crit_co2
242!
243! DESCRIPTION
244!     Check if CO2 ice surface and pressure change criteria are reached.
245!
246! AUTHORS & DATE
247!     R. Vandemeulebrouck
248!     L. Lange, 2023
249!     JB Clement, 2023-2025
250!
251! NOTES
252!
253!-----------------------------------------------------------------------
254
255! DEPENDENCIES
256! ------------
257use comslope_mod,  only: subslope_dist
258
259! DECLARATION
260! -----------
261implicit none
262
263! ARGUMENTS
264! ---------
265integer,                          intent(in)    :: ngrid, nslope          ! # of grid physical grid points
266real,    dimension(ngrid),        intent(in)    :: cell_area              ! Area of the cells
267real,    dimension(ngrid,nslope), intent(in)    :: co2_ice                ! Actual density of co2 ice
268real,                             intent(in)    :: co2ice_sublim_surf_ini ! Initial surface of sublimatingco2 ice
269logical, dimension(ngrid,nslope), intent(in)    :: is_co2ice_sublim_ini   ! Grid points where co2 ice was initially sublimating
270real,                             intent(in)    :: ps_avg_global_ini      ! Planet average pressure from the PCM start files
271real,                             intent(in)    :: ps_avg_global          ! Planet average pressure from the PEM computations
272type(stopFlags),                  intent(inout) :: stopCrit               ! Stopping criterion
273
274! LOCAL VARIABLES
275! ---------------
276integer :: i, islope       ! Loop
277real    :: co2ice_now_surf ! Current surface of co2 ice
278
279! CODE
280! ----
281if (stopCrit%stop_code() > 0) return
282
283! Computation of the present surface of co2 ice still sublimating
284co2ice_now_surf = 0.
285do i = 1,ngrid
286    do islope = 1,nslope
287        if (is_co2ice_sublim_ini(i,islope) .and. co2_ice(i,islope) > 0.) co2ice_now_surf = co2ice_now_surf + cell_area(i)*subslope_dist(i,islope)
288    enddo
289enddo
290
291! Check of the criterion
292if (co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2ice_crit)) then
293    stopCrit%surface_co2ice_change = .true.
294    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
295    write(*,*) "co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2ice_crit)", co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2ice_crit)
296    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_sublim_surf_ini
297    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
298    write(*,*) "Percentage of change accepted =", co2ice_crit*100.
299else if (co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2ice_crit)) then
300    stopCrit%surface_co2ice_change = .true.
301    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
302    write(*,*) "co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2ice_crit)", co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2ice_crit)
303    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
304    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_sublim_surf_ini
305    write(*,*) "Percentage of change accepted =", co2ice_crit*100.
306endif
307
308if (ps_avg_global < ps_avg_global_ini*(1. - ps_crit)) then
309    stopCrit%pressure_change = .true.
310    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
311    write(*,*) "ps_avg_global < ps_avg_global_ini*(1. - ps_crit)", ps_avg_global < ps_avg_global_ini*(1. - ps_crit)
312    write(*,*) "Initial global pressure =", ps_avg_global_ini
313    write(*,*) "Current global pressure =", ps_avg_global
314    write(*,*) "Percentage of change accepted =", ps_crit*100.
315else if (ps_avg_global > ps_avg_global_ini*(1. + ps_crit)) then
316    stopCrit%pressure_change = .true.
317    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
318    write(*,*) "ps_avg_global > ps_avg_global_ini*(1. + ps_crit)", ps_avg_global > ps_avg_global_ini*(1. + ps_crit)
319    write(*,*) "Initial global pressure =", ps_avg_global_ini
320    write(*,*) "Current global pressure =", ps_avg_global
321    write(*,*) "Percentage of change accepted =", ps_crit*100.
322endif
323
324END SUBROUTINE stopping_crit_co2
325!=======================================================================
326
327!=======================================================================
328SUBROUTINE stopping_crit_h2o(ngrid,nslope,cell_area,delta_h2o_adsorbed,delta_h2o_icetablesublim,h2o_ice,d_h2oice,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,stopCrit)
329!-----------------------------------------------------------------------
330! NAME
331!     stopping_crit_h2o
332!
333! DESCRIPTION
334!     Check if PEM oscillates infinitely with H2O (no net exchange).
335!
336! AUTHORS & DATE
337!     L. Lange, 2024
338!     JB Clement, 2025
339!
340! NOTES
341!
342!-----------------------------------------------------------------------
343
344! DEPENDENCIES
345! ------------
346use evolution,      only: dt
347use comslope_mod,   only: subslope_dist, def_slope_mean
348use phys_constants, only: pi
349
350! DECLARATION
351! -----------
352implicit none
353
354! ARGUMENTS
355! ---------
356integer,                       intent(in)    :: ngrid, nslope            ! # of grid points, # of subslopes
357real, dimension(ngrid),        intent(in)    :: cell_area                ! Area of each mesh grid (m^2)
358real, dimension(ngrid),        intent(in)    :: delta_h2o_adsorbed       ! Mass of H2O adsorbed/desorbded in the soil (kg/m^2)
359real, dimension(ngrid),        intent(in)    :: delta_h2o_icetablesublim ! Mass of H2O that condensed/sublimated at the ice table
360real, dimension(ngrid,nslope), intent(in)    :: h2o_ice                  ! H2O ice (kg/m^2)
361real, dimension(ngrid,nslope), intent(in)    :: d_h2oice                 ! Tendency of H2O ice (kg/m^2/year)
362real,                          intent(out)   :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! Variables to conserve H2O
363type(stopFlags),               intent(inout) :: stopCrit                 ! Stopping criterion
364
365! LOCAL VARIABLES
366! ---------------
367integer :: i, islope ! Loop
368
369! CODE
370! ----
371! Checks to escape the subroutine (adjusment not relevant in 1D or if the PEM should stop already)
372if (ngrid == 1 .or. stopCrit%stop_code() > 0) then
373    S_atm_2_h2o = 1.
374    S_h2o_2_atm = S_atm_2_h2o
375    S_atm_2_h2oice = 1.
376    S_h2oice_2_atm = S_atm_2_h2oice
377    return
378endif
379
380! We compute the amount of water going out from the atmosphere (S_atm_2_h2o) and going into the atmophere (S_h2o_2_atm)
381S_atm_2_h2o = 0.
382S_h2o_2_atm = 0.
383S_atm_2_h2oice = 0.
384S_h2oice_2_atm = 0.
385do i = 1,ngrid
386    if (delta_h2o_adsorbed(i) > 0.) then
387        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_adsorbed(i)*cell_area(i)
388    else
389        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_adsorbed(i)*cell_area(i)
390    endif
391    if (delta_h2o_icetablesublim(i) > 0.) then
392        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_icetablesublim(i)*cell_area(i)
393    else
394        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_icetablesublim(i)*cell_area(i)
395    endif
396    do islope = 1,nslope
397        if (d_h2oice(i,islope) > 0.) then
398            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.)
399            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.)
400        else if (d_h2oice(i,islope) < 0. .and. h2o_ice(i,islope) > 0.) then
401            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.)
402            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.)
403        endif
404    enddo
405enddo
406
407! Since relative atmospheric water is kept constant, we need to equate condensing reservoirs to the sublimating ones.
408! It is not possible if one of them is missing.
409if ((abs(S_atm_2_h2oice) < 1.e-10 .or. abs(S_h2oice_2_atm) < 1.e-10)) then
410    write(*,*) "Reason of stopping: there is no sublimating or condensing h2o ice!"
411    write(*,*) "This can be due to the absence of h2o ice in the PCM run."
412    write(*,*) "Amount of condensing ice  =", S_atm_2_h2oice
413    write(*,*) "Amount of sublimating ice =", S_h2oice_2_atm
414    stopCrit%zero_dh2oice = .true.
415endif
416
417END SUBROUTINE stopping_crit_h2o
418!=======================================================================
419
420END MODULE stopping_crit
Note: See TracBrowser for help on using the repository browser.