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

Last change on this file since 3786 was 3770, checked in by jbclement, 4 weeks ago

PEM:
Revision of the layering structure and algorithm:

  • 'stratum' components are now expressed in height
  • deletion of the (redundant) thickness feature of 'stratum'
  • 'stratif' in the PEM main program is renamed 'deposits'
  • addition of several procedures to get useful information about a stratum (major component or thickness)
  • all subsurface layers are now represented in the layering structure by strata with negative top elevation
  • simplification of the different situations arising from the input tendencies
  • porosity is determined for the entire stratum (and not anymore for each component) based on dominant component
  • sublimation of CO2 and H2O ice is now handled simultaneously (more realistic) in a stratum
  • linking the layering algorithm with the PEM initilization/finalization regarding PCM data and with the PEM stopping criteria
  • making separate cases for glaciers vs layering management
  • H2O sublimation flux correction is now handled with the layering when a dust lag layer layer is growing
  • update of 'h2o_ice_depth' and 'zdqsdif' accordingly at the PEM end for the PCM

JBC

File size: 6.9 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
151END MODULE
Note: See TracBrowser for help on using the repository browser.