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

Last change on this file since 3956 was 3938, checked in by jbclement, 3 months ago

PEM:

  • Correction of the process to balance the H2O flux from and into the atmosphere accross reservoirs: (i) computation of H2O amount going from/in the atmosphere is corrected (missing dt and wrong condition); (ii) computation of the scaling factor is corrected because we need to balance all the icetable/adsorbed/surface ice flux but it affects only sublimating/condensing surface ice flux; (iii) process of balancing is modified to be applied to every flux by scaling instead of applying it only to positive or negative flux by trimming; (iv) balance is now fully processed by the tendency before evolving the ice. This is done through dedicated subroutines.
  • Correction of ice conversion between the layering algorithm and PEM ice variables (density factor was missing).
  • Addition of H2O flux balance and related stopping criterion for the layering algorithm.
  • Few updates for files in the deftank to be in line with the wiki Planeto pages.

JBC

File size: 10.5 KB
Line 
1MODULE stopping_crit_mod
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8
9!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
10!!!
11!!! Purpose: Criterions to check if the PEM needs to call the PCM
12!!! Author: RV & LL, 02/2023
13!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15
16SUBROUTINE stopping_crit_h2o_ice(cell_area,h2oice_ini_surf,is_h2oice_sublim_ini,h2o_ice,stopPEM,ngrid)
17
18use time_evol_mod, only: h2o_ice_crit
19use comslope_mod,  only: subslope_dist, nslope
20
21implicit none
22
23!=======================================================================
24!
25! Routine to check if the h2o ice criterion to stop the PEM is reached
26!
27!=======================================================================
28! Inputs
29!-------
30integer,                          intent(in) :: ngrid                ! # of physical grid points
31real,    dimension(ngrid),        intent(in) :: cell_area            ! Area of the cells
32real,    dimension(ngrid,nslope), intent(in) :: h2o_ice              ! Actual density of h2o ice
33real,                             intent(in) :: h2oice_ini_surf      ! Initial surface of sublimating h2o ice
34logical, dimension(ngrid,nslope), intent(in) :: is_h2oice_sublim_ini ! Grid points where h2o ice was initially sublimating
35! Outputs
36!--------
37integer, intent(inout) :: stopPEM ! Stopping criterion code
38! Locals
39! ------
40integer :: i, islope       ! Loop
41real    :: h2oice_now_surf ! Current surface of h2o ice
42
43!=======================================================================
44if (stopPEM > 0) return
45
46! Computation of the present surface of h2o ice still sublimating
47h2oice_now_surf = 0.
48do i = 1,ngrid
49    do islope = 1,nslope
50        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)
51    enddo
52enddo
53
54! Check of the criterion
55if (h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)) then
56    stopPEM = 1
57    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
58    write(*,*) "h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)", h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)
59    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
60    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
61    write(*,*) "Percentage of change accepted =", h2o_ice_crit*100
62else if (h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)) then
63    stopPEM = 1
64    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
65    write(*,*) "h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)", h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)
66    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
67    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
68    write(*,*) "Percentage of change accepted =", h2o_ice_crit*100
69endif
70
71END SUBROUTINE stopping_crit_h2o_ice
72
73!=======================================================================
74
75SUBROUTINE stopping_crit_co2(cell_area,co2ice_sublim_surf_ini,is_co2ice_sublim_ini,co2_ice,stopPEM,ngrid,ps_avg_global_ini,ps_avg_global,nslope)
76
77use time_evol_mod, only: co2_ice_crit, ps_criterion
78use comslope_mod,  only: subslope_dist
79
80implicit none
81
82!=======================================================================
83!
84! Routine to check if the co2 and pressure criteria to stop the PEM are reached
85!
86!=======================================================================
87! Inputs
88!-------
89integer,                          intent(in) :: ngrid, nslope          ! # of grid physical grid points
90real,    dimension(ngrid),        intent(in) :: cell_area              ! Area of the cells
91real,    dimension(ngrid,nslope), intent(in) :: co2_ice                ! Actual density of co2 ice
92real,                             intent(in) :: co2ice_sublim_surf_ini ! Initial surface of sublimatingco2 ice
93logical, dimension(ngrid,nslope), intent(in) :: is_co2ice_sublim_ini   ! Grid points where co2 ice was initially sublimating
94real,                             intent(in) :: ps_avg_global_ini      ! Planet average pressure from the PCM start files
95real,                             intent(in) :: ps_avg_global          ! Planet average pressure from the PEM computations
96! Outputs
97!--------
98integer, intent(inout) :: stopPEM ! Stopping criterion code
99
100! Locals
101! ------
102integer :: i, islope       ! Loop
103real    :: co2ice_now_surf ! Current surface of co2 ice
104
105!=======================================================================
106if (stopPEM > 0) return
107
108! Computation of the present surface of co2 ice still sublimating
109co2ice_now_surf = 0.
110do i = 1,ngrid
111    do islope = 1,nslope
112        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)
113    enddo
114enddo
115
116! Check of the criterion
117if (co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2_ice_crit)) then
118    stopPEM = 3
119    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
120    write(*,*) "co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2_ice_crit)", co2ice_now_surf < co2ice_sublim_surf_ini*(1. - co2_ice_crit)
121    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_sublim_surf_ini
122    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
123    write(*,*) "Percentage of change accepted =", co2_ice_crit*100.
124else if (co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2_ice_crit)) then
125    stopPEM = 3
126    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
127    write(*,*) "co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2_ice_crit)", co2ice_now_surf > co2ice_sublim_surf_ini*(1. + co2_ice_crit)
128    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
129    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_sublim_surf_ini
130    write(*,*) "Percentage of change accepted =", co2_ice_crit*100.
131endif
132
133if (ps_avg_global < ps_avg_global_ini*(1. - ps_criterion)) then
134    stopPEM = 4
135    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
136    write(*,*) "ps_avg_global < ps_avg_global_ini*(1. - ps_criterion)", ps_avg_global < ps_avg_global_ini*(1. - ps_criterion)
137    write(*,*) "Initial global pressure =", ps_avg_global_ini
138    write(*,*) "Current global pressure =", ps_avg_global
139    write(*,*) "Percentage of change accepted =", ps_criterion*100.
140else if (ps_avg_global > ps_avg_global_ini*(1. + ps_criterion)) then
141    stopPEM = 4
142    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
143    write(*,*) "ps_avg_global > ps_avg_global_ini*(1. + ps_criterion)", ps_avg_global > ps_avg_global_ini*(1. + ps_criterion)
144    write(*,*) "Initial global pressure =", ps_avg_global_ini
145    write(*,*) "Current global pressure =", ps_avg_global
146    write(*,*) "Percentage of change accepted =", ps_criterion*100.
147endif
148
149END SUBROUTINE stopping_crit_co2
150
151!=======================================================================
152
153SUBROUTINE 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,stopPEM)
154
155use time_evol_mod, only: dt
156use comslope_mod,  only: subslope_dist, def_slope_mean
157#ifndef CPP_STD
158    use comcstfi_h,   only: pi
159#else
160    use comcstfi_mod, only: pi
161#endif
162
163implicit none
164
165!=======================================================================
166!
167! Routine to check if the h2o is only exchanged between grid points
168!
169!=======================================================================
170! Inputs
171!-------
172integer,                       intent(in) :: ngrid, nslope            ! # of grid points, # of subslopes
173real, dimension(ngrid),        intent(in) :: cell_area                ! Area of each mesh grid (m^2)
174real, dimension(ngrid),        intent(in) :: delta_h2o_adsorbed       ! Mass of H2O adsorbed/desorbded in the soil (kg/m^2)
175real, dimension(ngrid),        intent(in) :: delta_h2o_icetablesublim ! Mass of H2O that condensed/sublimated at the ice table
176real, dimension(ngrid,nslope), intent(in) :: h2o_ice                  ! H2O ice (kg/m^2)
177real, dimension(ngrid,nslope), intent(in) :: d_h2oice                 ! Tendency of H2O ice (kg/m^2/year)
178! Outputs
179!--------
180integer, intent(inout) :: stopPEM ! Stopping criterion code
181real,    intent(out) :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! Variables to conserve H2O
182
183! Locals
184! ------
185integer :: i, islope ! Loop
186
187!=======================================================================
188if (stopPEM > 0) return
189
190! We compute the amount of water going out from the atmosphere (S_atm_2_h2o) and going into the atmophere (S_h2o_2_atm)
191S_atm_2_h2o = 0.
192S_h2o_2_atm = 0.
193S_atm_2_h2oice = 0.
194S_h2oice_2_atm = 0.
195do i = 1,ngrid
196    if (delta_h2o_adsorbed(i) > 0.) then
197        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_adsorbed(i)*cell_area(i)
198    else
199        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_adsorbed(i)*cell_area(i)
200    endif
201    if (delta_h2o_icetablesublim(i) > 0.) then
202        S_atm_2_h2o = S_atm_2_h2o + delta_h2o_icetablesublim(i)*cell_area(i)
203    else
204        S_h2o_2_atm = S_h2o_2_atm + delta_h2o_icetablesublim(i)*cell_area(i)
205    endif
206    do islope = 1,nslope
207        if (d_h2oice(i,islope) > 0.) then
208            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.)
209            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.)
210        else if (d_h2oice(i,islope) < 0. .and. h2o_ice(i,islope) > 0.) then
211            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.)
212            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.)
213        endif
214    enddo
215enddo
216
217! Since relative atmospheric water is kept constant, we need to equate condensing reservoirs to the sublimating ones.
218! It is not possible if one of them is missing.
219if (abs(S_atm_2_h2oice) < 1.e-10 .or. abs(S_h2oice_2_atm) < 1.e-10) then
220    write(*,*) "Reason of stopping: there is no sublimating or condensing h2o ice!"
221    write(*,*) "This can be due to the absence of h2o ice in the PCM run."
222    write(*,*) "Amount of condensing ice  =", S_atm_2_h2oice
223    write(*,*) "Amount of sublimating ice =", S_h2oice_2_atm
224    stopPEM = 2
225endif
226
227END SUBROUTINE stopping_crit_h2o
228
229END MODULE
Note: See TracBrowser for help on using the repository browser.