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

Last change on this file since 3386 was 3339, checked in by jbclement, 13 months ago

PEM:
Correction of the way sublimating ice is identified at the beginning of the PEM + some updates for the default variables and the display of information.
JBC

File size: 7.0 KB
RevLine 
[3149]1MODULE stopping_crit_mod
[2888]2
[3130]3implicit none
[2888]4
[3149]5!=======================================================================
[3130]6contains
[3149]7!=======================================================================
[2888]8
9!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
10!!!
[3130]11!!! Purpose: Criterions to check if the PEM needs to call the PCM
[2888]12!!! Author: RV & LL, 02/2023
13!!!
14!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
15
[3327]16SUBROUTINE stopping_crit_h2o_ice(cell_area,h2oice_ini_surf,ini_h2oice_sublim,h2o_ice,stopPEM,ngrid)
[2888]17
[3159]18use time_evol_mod, only: h2o_ice_crit
[3130]19use comslope_mod,  only: subslope_dist, nslope
[2888]20
[3130]21implicit none
[2888]22
23!=======================================================================
24!
[3149]25! Routine to check if the h2o ice criterion to stop the PEM is reached
[2888]26!
27!=======================================================================
[3327]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) :: ini_h2oice_sublim ! Grid points where h2o ice was initially sublimating
35! Outputs
36!--------
[3149]37integer, intent(inout) :: stopPEM ! Stopping criterion code
[3327]38! Locals
39! ------
40integer :: i, islope       ! Loop
41real    :: h2oice_now_surf ! Current surface of h2o ice
[2888]42
43!=======================================================================
[3327]44! Computation of the present surface of h2o ice still sublimating
45h2oice_now_surf = 0.
[3130]46do i = 1,ngrid
47    do islope = 1,nslope
[3327]48        if (ini_h2oice_sublim(i,islope) .and. h2o_ice(i,islope) > 0.) h2oice_now_surf = h2oice_now_surf + cell_area(i)*subslope_dist(i,islope)
[2888]49    enddo
[3130]50enddo
[2888]51
[3130]52! Check of the criterion
[3327]53if (h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)) then
[3149]54    stopPEM = 1
[3159]55    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
[3327]56    write(*,*) "h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)", h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)
[3339]57    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
[3327]58    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
[3159]59    write(*,*) "Percentage of change accepted =", h2o_ice_crit*100
[3327]60else if (h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)) then
[3149]61    stopPEM = 1
[3159]62    write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold"
[3327]63    write(*,*) "h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)", h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)
[3339]64    write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf
[3327]65    write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf
[3159]66    write(*,*) "Percentage of change accepted =", h2o_ice_crit*100
[3143]67endif
[3130]68
[3327]69if (abs(h2oice_ini_surf) < 1.e-5) stopPEM = 0
[3130]70
[3149]71END SUBROUTINE stopping_crit_h2o_ice
[2888]72
[3149]73!=======================================================================
[2888]74
[3327]75SUBROUTINE stopping_crit_co2(cell_area,co2ice_ini_surf,ini_co2ice_sublim,co2_ice,stopPEM,ngrid,global_avg_press_PCM,global_avg_press_new,nslope)
[2888]76
[3159]77use time_evol_mod, only: co2_ice_crit, ps_criterion
[3130]78use comslope_mod,  only: subslope_dist
[2888]79
[3130]80implicit none
[2888]81
82!=======================================================================
83!
[3149]84! Routine to check if the co2 and pressure criteria to stop the PEM are reached
[2888]85!
86!=======================================================================
[3327]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_ini_surf      ! Initial surface of sublimatingco2 ice
93logical, dimension(ngrid,nslope), intent(in) :: ini_co2ice_sublim    ! Grid points where co2 ice was initially sublimating
94real,                             intent(in) :: global_avg_press_PCM ! Planet average pressure from the PCM start files
95real,                             intent(in) :: global_avg_press_new ! Planet average pressure from the PEM computations
96! Outputs
97!--------
[3149]98integer, intent(inout) :: stopPEM ! Stopping criterion code
99
[3327]100! Locals
101! ------
102integer :: i, islope       ! Loop
103real    :: co2ice_now_surf ! Current surface of co2 ice
[2888]104
105!=======================================================================
[3327]106! Computation of the present surface of co2 ice still sublimating
107co2ice_now_surf = 0.
[3130]108do i = 1,ngrid
109    do islope = 1,nslope
[3327]110        if (ini_co2ice_sublim(i,islope) .and. co2_ice(i,islope) > 0.) co2ice_now_surf = co2ice_now_surf + cell_area(i)*subslope_dist(i,islope)
[3130]111    enddo
112enddo
[2888]113
[3130]114! Check of the criterion
[3327]115if (co2ice_now_surf < co2ice_ini_surf*(1. - co2_ice_crit)) then
[3149]116    stopPEM = 3
117    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
[3327]118    write(*,*) "co2ice_now_surf < co2ice_ini_surf*(1. - co2_ice_crit)", co2ice_now_surf < co2ice_ini_surf*(1. - co2_ice_crit)
[3339]119    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_ini_surf
[3327]120    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
[3159]121    write(*,*) "Percentage of change accepted =", co2_ice_crit*100.
[3327]122else if (co2ice_now_surf > co2ice_ini_surf*(1. + co2_ice_crit)) then
[3149]123    stopPEM = 3
124    write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold"
[3327]125    write(*,*) "co2ice_now_surf > co2ice_ini_surf*(1. + co2_ice_crit)", co2ice_now_surf > co2ice_ini_surf*(1. + co2_ice_crit)
126    write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf
127    write(*,*) "Initial surface of co2 ice sublimating =", co2ice_ini_surf
[3159]128    write(*,*) "Percentage of change accepted =", co2_ice_crit*100.
[3130]129endif
[2888]130
[3327]131if (abs(co2ice_ini_surf) < 1.e-5) stopPEM = 0
[2888]132
[3149]133if (global_avg_press_new < global_avg_press_PCM*(1. - ps_criterion)) then
134    stopPEM = 4
135    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
136    write(*,*) "global_avg_press_new < global_avg_press_PCM*(1. - ps_criterion)", global_avg_press_new < global_avg_press_PCM*(1. - ps_criterion)
[3339]137    write(*,*) "Initial global pressure =", global_avg_press_PCM
[3149]138    write(*,*) "Current global pressure =", global_avg_press_new
[3130]139    write(*,*) "Percentage of change accepted =", ps_criterion*100.
[3149]140else if (global_avg_press_new > global_avg_press_PCM*(1. + ps_criterion)) then
141    stopPEM = 4
142    write(*,*) "Reason of stopping: the global pressure reaches the threshold"
143    write(*,*) "global_avg_press_new > global_avg_press_PCM*(1. + ps_criterion)", global_avg_press_new > global_avg_press_PCM*(1. + ps_criterion)
[3339]144    write(*,*) "Initial global pressure =", global_avg_press_PCM
[3149]145    write(*,*) "Current global pressure =", global_avg_press_new
[3143]146    write(*,*) "Percentage of change accepted =", ps_criterion*100.
147endif
[3130]148
[3149]149END SUBROUTINE stopping_crit_co2
[2888]150
151END MODULE
Note: See TracBrowser for help on using the repository browser.