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

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

PEM:
Massive structural refactor of the PEM codebase for improved readability, consistency and maintainability. The goal is to modernize, standardize and consolidate the code while removing legacy complexity. In detail, this change:

  • Performs large-scale cleanup removing unused code, obsolete routines, duplicated functionality and deprecated initialization logic;
  • Removes "*_mod" wrappers;
  • Replaces mixed naming conventions with clearer variable names, domain-based file/module names and purpose-based routine names;
  • Adds native reading/writing capabilities to the PEM removing the cumbersome dependency on Mars PCM subroutines;
  • Centralizes soil/ice/orbital/PEM configuration logic into dedicated modules;
  • Simplifies routines and updates calls/interfaces to match the new structure.

This refactor significantly clarifies the codebase and provides a cleaner foundation for forthcoming developments.
JBC

File size: 14.5 KB
Line 
1MODULE stopping_crit
2
3implicit none
4
5real :: h2oice_crit ! Percentage of change of the surface of h2o ice sublimating before stopping the PEM
6real :: co2ice_crit ! Percentage of change of the surface of co2 ice sublimating before stopping the PEM
7real :: ps_crit ! Percentage of change of averaged surface pressure before stopping the PEM
8
9type :: stopFlags
10    logical :: surface_h2oice_change = .false.
11    logical :: zero_dh2oice          = .false.
12    logical :: surface_co2ice_change = .false.
13    logical :: pressure_change       = .false.
14    logical :: nyears_max_reached    = .false.
15    logical :: nyears_maxorb_reached = .false.
16    logical :: nyears_tot_reached    = .false.
17    logical :: time_limit_reached    = .false.
18contains
19    procedure :: is_any_set
20    procedure :: stop_code
21    procedure :: stop_message
22end type
23
24!=======================================================================
25contains
26!=======================================================================
27
28!=======================================================================
29FUNCTION is_any_set(this) RESULT(stopPEM)
30! To detect if any flag is true
31
32class(stopFlags), intent(in) :: this
33logical :: stopPEM
34
35stopPEM = this%surface_h2oice_change .or. &
36          this%zero_dh2oice          .or. &
37          this%surface_co2ice_change .or. &
38          this%pressure_change       .or. &
39          this%nyears_max_reached    .or. &
40          this%nyears_maxorb_reached .or. &
41          this%nyears_tot_reached    .or. &
42          this%time_limit_reached
43
44END FUNCTION is_any_set
45!=======================================================================
46
47!=======================================================================
48! To give the digit code corresponding to the stopping flag
49FUNCTION stop_code(this) result(code)
50   
51class(StopFlags), intent(in) :: this
52integer :: code
53
54code = 0
55if     (this%surface_h2oice_change) then; code = 1
56elseif (this%zero_dh2oice)          then; code = 2
57elseif (this%surface_co2ice_change) then; code = 3
58elseif (this%pressure_change)       then; code = 4
59elseif (this%nyears_max_reached)    then; code = 5
60elseif (this%nyears_maxorb_reached) then; code = 6
61elseif (this%nyears_tot_reached)    then; code = 7
62elseif (this%time_limit_reached)    then; code = 8
63endif
64
65END FUNCTION stop_code
66!=======================================================================
67
68!=======================================================================
69! To give the message corresponding to the stopping flag
70FUNCTION stop_message(this) result(msg)
71   
72class(StopFlags), intent(in) :: this
73character(:), allocatable :: msg
74integer :: stopCode
75
76if     (this%surface_h2oice_change) then; msg = " **** STOPPING because surface of h2o ice has changed too much (see message above)."
77elseif (this%zero_dh2oice)          then; msg = " **** STOPPING because tendencies on h2o ice = 0 (see message above)."
78elseif (this%surface_co2ice_change) then; msg = " **** STOPPING because surface of co2 ice has changed too much (see message above)."
79elseif (this%pressure_change)       then; msg = " **** STOPPING because surface global pressure has changed too much (see message above)."
80elseif (this%nyears_max_reached)    then; msg = " **** STOPPING because maximum number of years is reached."
81elseif (this%nyears_maxorb_reached) then; msg = " **** STOPPING because maximum number of years due to orbital parameters is reached."
82elseif (this%nyears_tot_reached)    then; msg = " **** STOPPING because the number of years to be simulated is reached."
83elseif (this%time_limit_reached)    then; msg = " **** STOPPING because the job time limit is going to be reached."
84endif
85
86END FUNCTION stop_message
87!=======================================================================
88
89!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
90!!!
91!!! Purpose: Criterions to check if the PEM needs to call the PCM
92!!! Author: RV & LL, 02/2023
93!!!
94!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
95
96SUBROUTINE stopping_crit_h2o_ice(cell_area,h2oice_ini_surf,is_h2oice_sublim_ini,h2o_ice,stopCrit,ngrid)
97
98use comslope_mod,  only: subslope_dist, nslope
99
100implicit none
101
102!=======================================================================
103!
104! Routine to check if the h2o ice criterion to stop the PEM is reached
105!
106!=======================================================================
107! Inputs
108!-------
109integer,                          intent(in) :: ngrid                ! # of physical grid points
110real,    dimension(ngrid),        intent(in) :: cell_area            ! Area of the cells
111real,    dimension(ngrid,nslope), intent(in) :: h2o_ice              ! Actual density of h2o ice
112real,                             intent(in) :: h2oice_ini_surf      ! Initial surface of sublimating h2o ice
113logical, dimension(ngrid,nslope), intent(in) :: is_h2oice_sublim_ini ! Grid points where h2o ice was initially sublimating
114! Outputs
115!--------
116type(stopFlags), intent(inout) :: stopCrit ! Stopping criterion
117! Local variables
118! ---------------
119integer :: i, islope       ! Loop
120real    :: h2oice_now_surf ! Current surface of h2o ice
121
122!=======================================================================
123if (stopCrit%stop_code() > 0) return
124
125! Computation of the present surface of h2o ice still sublimating
126h2oice_now_surf = 0.
127do i = 1,ngrid
128    do islope = 1,nslope
129        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)
130    enddo
131enddo
132
133! Check of the criterion
134if (h2oice_now_surf < h2oice_ini_surf*(1. - h2oice_crit)) then
135    stopCrit%surface_h2oice_change = .true.
136    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
137    write(*,*) "h2oice_now_surf < h2oice_ini_surf*(1. - h2oice_crit)", h2oice_now_surf < h2oice_ini_surf*(1. - h2oice_crit)
138    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
139    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
140    write(*,*) "Percentage of change accepted =", h2oice_crit*100
141else if (h2oice_now_surf > h2oice_ini_surf*(1. + h2oice_crit)) then
142    stopCrit%surface_h2oice_change = .true.
143    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
144    write(*,*) "h2oice_now_surf > h2oice_ini_surf*(1. + h2oice_crit)", h2oice_now_surf > h2oice_ini_surf*(1. + h2oice_crit)
145    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
146    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
147    write(*,*) "Percentage of change accepted =", h2oice_crit*100
148endif
149
150END SUBROUTINE stopping_crit_h2o_ice
151!=======================================================================
152
153!=======================================================================
154SUBROUTINE 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)
155
156use comslope_mod,  only: subslope_dist
157
158implicit none
159
160!=======================================================================
161!
162! Routine to check if the co2 and pressure criteria to stop the PEM are reached
163!
164!=======================================================================
165! Inputs
166!-------
167integer,                          intent(in) :: ngrid, nslope          ! # of grid physical grid points
168real,    dimension(ngrid),        intent(in) :: cell_area              ! Area of the cells
169real,    dimension(ngrid,nslope), intent(in) :: co2_ice                ! Actual density of co2 ice
170real,                             intent(in) :: co2ice_sublim_surf_ini ! Initial surface of sublimatingco2 ice
171logical, dimension(ngrid,nslope), intent(in) :: is_co2ice_sublim_ini   ! Grid points where co2 ice was initially sublimating
172real,                             intent(in) :: ps_avg_global_ini      ! Planet average pressure from the PCM start files
173real,                             intent(in) :: ps_avg_global          ! Planet average pressure from the PEM computations
174! Outputs
175!--------
176type(stopFlags), intent(inout) :: stopCrit ! Stopping criterion
177
178! Local variables
179! ---------------
180integer :: i, islope       ! Loop
181real    :: co2ice_now_surf ! Current surface of co2 ice
182
183!=======================================================================
184if (stopCrit%stop_code() > 0) return
185
186! Computation of the present surface of co2 ice still sublimating
187co2ice_now_surf = 0.
188do i = 1,ngrid
189    do islope = 1,nslope
190        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)
191    enddo
192enddo
193
194! Check of the criterion
195if (co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2ice_crit)) then
196    stopCrit%surface_co2ice_change = .true.
197    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
198    write(*,*) "co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2ice_crit)", co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2ice_crit)
199    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_sublim_surf_ini
200    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
201    write(*,*) "Percentage of change accepted =", co2ice_crit*100.
202else if (co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2ice_crit)) then
203    stopCrit%surface_co2ice_change = .true.
204    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
205    write(*,*) "co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2ice_crit)", co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2ice_crit)
206    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
207    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_sublim_surf_ini
208    write(*,*) "Percentage of change accepted =", co2ice_crit*100.
209endif
210
211if (ps_avg_global < ps_avg_global_ini*(1. - ps_crit)) then
212    stopCrit%pressure_change = .true.
213    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
214    write(*,*) "ps_avg_global < ps_avg_global_ini*(1. - ps_crit)", ps_avg_global < ps_avg_global_ini*(1. - ps_crit)
215    write(*,*) "Initial global pressure =", ps_avg_global_ini
216    write(*,*) "Current global pressure =", ps_avg_global
217    write(*,*) "Percentage of change accepted =", ps_crit*100.
218else if (ps_avg_global > ps_avg_global_ini*(1. + ps_crit)) then
219    stopCrit%pressure_change = .true.
220    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
221    write(*,*) "ps_avg_global > ps_avg_global_ini*(1. + ps_crit)", ps_avg_global > ps_avg_global_ini*(1. + ps_crit)
222    write(*,*) "Initial global pressure =", ps_avg_global_ini
223    write(*,*) "Current global pressure =", ps_avg_global
224    write(*,*) "Percentage of change accepted =", ps_crit*100.
225endif
226
227END SUBROUTINE stopping_crit_co2
228!=======================================================================
229
230!=======================================================================
231SUBROUTINE 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)
232
233use evolution,      only: dt
234use comslope_mod,   only: subslope_dist, def_slope_mean
235use phys_constants, only: pi
236
237implicit none
238
239!=======================================================================
240!
241! Routine to check if the h2o is only exchanged between grid points
242!
243!=======================================================================
244! Inputs
245!-------
246integer,                       intent(in) :: ngrid, nslope            ! # of grid points, # of subslopes
247real, dimension(ngrid),        intent(in) :: cell_area                ! Area of each mesh grid (m^2)
248real, dimension(ngrid),        intent(in) :: delta_h2o_adsorbed       ! Mass of H2O adsorbed/desorbded in the soil (kg/m^2)
249real, dimension(ngrid),        intent(in) :: delta_h2o_icetablesublim ! Mass of H2O that condensed/sublimated at the ice table
250real, dimension(ngrid,nslope), intent(in) :: h2o_ice                  ! H2O ice (kg/m^2)
251real, dimension(ngrid,nslope), intent(in) :: d_h2oice                 ! Tendency of H2O ice (kg/m^2/year)
252! Outputs
253!--------
254type(stopFlags), intent(inout) :: stopCrit ! Stopping criterion
255real, intent(out) :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! Variables to conserve H2O
256! Local variables
257! ---------------
258integer :: i, islope ! Loop
259!=======================================================================
260! Checks to escape the subroutine (adjusment not relevant in 1D or if the PEM should stop already)
261if (ngrid == 1 .or. stopCrit%stop_code() > 0) then
262    S_atm_2_h2o = 1.
263    S_h2o_2_atm = S_atm_2_h2o
264    S_atm_2_h2oice = 1.
265    S_h2oice_2_atm = S_atm_2_h2oice
266    return
267endif
268
269! We compute the amount of water going out from the atmosphere (S_atm_2_h2o) and going into the atmophere (S_h2o_2_atm)
270S_atm_2_h2o = 0.
271S_h2o_2_atm = 0.
272S_atm_2_h2oice = 0.
273S_h2oice_2_atm = 0.
274do i = 1,ngrid
275    if (delta_h2o_adsorbed(i) > 0.) then
276        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_adsorbed(i)*cell_area(i)
277    else
278        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_adsorbed(i)*cell_area(i)
279    endif
280    if (delta_h2o_icetablesublim(i) > 0.) then
281        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_icetablesublim(i)*cell_area(i)
282    else
283        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_icetablesublim(i)*cell_area(i)
284    endif
285    do islope = 1,nslope
286        if (d_h2oice(i,islope) > 0.) then
287            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.)
288            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.)
289        else if (d_h2oice(i,islope) < 0. .and. h2o_ice(i,islope) > 0.) then
290            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.)
291            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.)
292        endif
293    enddo
294enddo
295
296! Since relative atmospheric water is kept constant, we need to equate condensing reservoirs to the sublimating ones.
297! It is not possible if one of them is missing.
298if ((abs(S_atm_2_h2oice) < 1.e-10 .or. abs(S_h2oice_2_atm) < 1.e-10)) then
299    write(*,*) "Reason of stopping: there is no sublimating or condensing h2o ice!"
300    write(*,*) "This can be due to the absence of h2o ice in the PCM run."
301    write(*,*) "Amount of condensing ice  =", S_atm_2_h2oice
302    write(*,*) "Amount of sublimating ice =", S_h2oice_2_atm
303    stopCrit%zero_dh2oice = .true.
304endif
305
306END SUBROUTINE stopping_crit_h2o
307!=======================================================================
308
309END MODULE
Note: See TracBrowser for help on using the repository browser.