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

Last change on this file since 4092 was 4074, checked in by jbclement, 3 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: 15.7 KB
RevLine 
[3989]1MODULE stopping_crit
[3991]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!-----------------------------------------------------------------------
[2888]15
[4065]16! DEPENDENCIES
17! ------------
18use numerics, only: dp, qp, di, k4, minieps_qp
19
[3991]20! DECLARATION
21! -----------
[3130]22implicit none
[2888]23
[4065]24! PARAMETERS
25! ----------
[4074]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
[4065]28real(dp), protected :: ps_crit     ! Change of averaged surface pressure before stopping the PEM
[3989]29
[4065]30! VARIABLES
31! ---------
[3967]32type :: stopFlags
[4065]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.
[3967]41contains
42    procedure :: is_any_set
43    procedure :: stop_code
44    procedure :: stop_message
45end type
46
[3130]47contains
[3991]48!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
[3989]49
50!=======================================================================
[4065]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!=======================================================================
[3989]92FUNCTION is_any_set(this) RESULT(stopPEM)
[3991]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!-----------------------------------------------------------------------
[2888]106
[3991]107! DECLARATION
108! -----------
109implicit none
[3967]110
[3991]111! ARGUMENTS
112! ---------
113class(stopFlags), intent(in) :: this    ! Stopping flags object
[4065]114logical(k4)                  :: stopPEM ! True if any flag is set
[3991]115
116! CODE
117! ----
[4065]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. &
[3967]125          this%time_limit_reached
126
127END FUNCTION is_any_set
[3989]128!=======================================================================
[3967]129
130!=======================================================================
131FUNCTION stop_code(this) result(code)
[3991]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!-----------------------------------------------------------------------
[3967]145
[3991]146! DECLARATION
147! -----------
148implicit none
149
150! ARGUMENTS
151! ---------
152class(StopFlags), intent(in) :: this ! Stopping flags object
[4065]153integer(di)                  :: code ! Code corresponding to active flag (0 if none)
[3991]154
155! CODE
156! ----
[3967]157code = 0
[4065]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
[3967]167
168END FUNCTION stop_code
[3989]169!=======================================================================
[3967]170
171!=======================================================================
172FUNCTION stop_message(this) result(msg)
[3991]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!-----------------------------------------------------------------------
[3967]186
[3991]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! ----
[4074]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)."
[4065]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
[3967]207
208END FUNCTION stop_message
[3989]209!=======================================================================
[3967]210
[3991]211!=======================================================================
[4065]212SUBROUTINE stopping_crit_h2oice(h2oice_sublim_coverage_ini,h2o_ice,d_h2oice,stopcrit)
[2888]213
[3991]214!-----------------------------------------------------------------------
215! NAME
[4065]216!     stopping_crit_h2oice
[3991]217!
218! DESCRIPTION
[4065]219!     Check if H2O ice surface has changed too much.
[3991]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! ------------
[4065]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
[2888]237
[3991]238! DECLARATION
239! -----------
[3130]240implicit none
[2888]241
[3991]242! ARGUMENTS
243! ---------
[4065]244real(dp), dimension(:,:), intent(in)    :: h2o_ice, d_h2oice
245real(dp),                 intent(in)    :: h2oice_sublim_coverage_ini
246type(stopFlags),          intent(inout) :: stopcrit
[3991]247
248! LOCAL VARIABLES
[3989]249! ---------------
[4065]250integer(di) :: i, islope
251real(dp)    :: h2oice_now_surf
[2888]252
[3991]253! CODE
254! ----
[4074]255! Check to escape the subroutine (if the PEM should stop already)
256if (stopcrit%stop_code() > 0) return
[3430]257
[4065]258! Computation of the present surface of sublimating H2O ice
259h2oice_now_surf = 0._dp
[3130]260do i = 1,ngrid
261    do islope = 1,nslope
[4065]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
[2888]265
[3130]266! Check of the criterion
[4065]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
[3130]272
[4065]273END SUBROUTINE stopping_crit_h2oice
[3989]274!=======================================================================
[2888]275
[3149]276!=======================================================================
[4065]277SUBROUTINE stopping_crit_co2ice(co2ice_sublim_coverage_ini,co2_ice,d_co2ice,stopcrit)
[3991]278!-----------------------------------------------------------------------
279! NAME
[4065]280!     stopping_crit_co2ice
[3991]281!
282! DESCRIPTION
[4065]283!     Check if CO2 ice surface has changed too much.
[3991]284!
285! AUTHORS & DATE
286!     R. Vandemeulebrouck
287!     L. Lange, 2023
288!     JB Clement, 2023-2025
289!
290! NOTES
291!
292!-----------------------------------------------------------------------
[2888]293
[3991]294! DEPENDENCIES
295! ------------
[4065]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
[2888]301
[3991]302! DECLARATION
303! -----------
[3130]304implicit none
[2888]305
[3991]306! ARGUMENTS
307! ---------
[4065]308real(dp), dimension(:,:), intent(in)    :: co2_ice, d_co2ice
309real(dp),                 intent(in)    :: co2ice_sublim_coverage_ini
310type(stopFlags),          intent(inout) :: stopcrit
[3149]311
[3991]312! LOCAL VARIABLES
[3989]313! ---------------
[4065]314integer(di) :: i, islope
[4074]315real(dp)    :: co2ice_now_surf ! Current surface of CO2 ice
[2888]316
[3991]317! CODE
318! ----
[4065]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
[3430]321
[4074]322! Computation of the present surface of CO2 ice still sublimating
[3327]323co2ice_now_surf = 0.
[3130]324do i = 1,ngrid
325    do islope = 1,nslope
[4065]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
[2888]329
[3130]330! Check of the criterion
[4065]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.
[4074]333    call print_msg("Initial surface of CO2 ice sublimating (ini|now) [m2] = "//real2str(co2ice_sublim_coverage_ini)//' | '//real2str(co2ice_now_surf))
[4065]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
[2888]336
[4065]337END SUBROUTINE stopping_crit_co2ice
338!=======================================================================
[3130]339
[3989]340!=======================================================================
[4065]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!-----------------------------------------------------------------------
[2888]357
[4065]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
[3938]384!=======================================================================
[4065]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)
[3991]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!-----------------------------------------------------------------------
[3938]402
[3991]403! DEPENDENCIES
404! ------------
[4065]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
[3938]411
[3991]412! DECLARATION
413! -----------
[3938]414implicit none
415
[3991]416! ARGUMENTS
417! ---------
[4065]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
[3991]423
424! LOCAL VARIABLES
[3989]425! ---------------
[4065]426integer(di) :: i, islope
[3991]427
428! CODE
429! ----
[4065]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
[3977]433    S_h2o_2_atm = S_atm_2_h2o
[4065]434    S_atm_2_h2oice = 1._qp
[3977]435    S_h2oice_2_atm = S_atm_2_h2oice
436    return
[4065]437end if
[3938]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)
[4065]440S_atm_2_h2o = 0._qp
441S_h2o_2_atm = 0._qp
442S_atm_2_h2oice = 0._qp
443S_h2oice_2_atm = 0._qp
[3938]444do i = 1,ngrid
[4065]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)
[3938]447    else
[4065]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)
[3938]452    else
[4065]453        S_h2o_2_atm = S_h2o_2_atm + delta_icetable(i)*cell_area(i)
454    end if
[3938]455    do islope = 1,nslope
[4065]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
[3938]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.
[4065]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
[3938]474
475END SUBROUTINE stopping_crit_h2o
[3989]476!=======================================================================
[3938]477
[3991]478END MODULE stopping_crit
Note: See TracBrowser for help on using the repository browser.