source: trunk/LMDZ.COMMON/libf/evolution/stopping_crit_mod.F90 @ 3977

Last change on this file since 3977 was 3977, checked in by jbclement, 4 days ago

PEM:

  • All operations computed by the PEM on the PCM data (averages, minima) are now performed by XIOS with two dedicated file. One is for daily operation outputs ("Xoutdaily4pem*.nc") and the other for yearly operation outputs ("Xoutyearly4pem*.nc").
  • Deletion of the reshaping tool "reshape_XIOS_output" used to convert XIOS outputs onto the PCM grid. Thus, the PEM is now able to read directly the format of XIOS outputs.
  • Addition of subroutines to convert data between a lon x lat array and a vector.

JBC

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