MODULE stopping_crit_mod implicit none !======================================================================= contains !======================================================================= !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !!! !!! Purpose: Criterions to check if the PEM needs to call the PCM !!! Author: RV & LL, 02/2023 !!! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! SUBROUTINE stopping_crit_h2o_ice(cell_area,h2oice_ini_surf,ini_h2oice_sublim,h2o_ice,stopPEM,ngrid) use time_evol_mod, only: h2o_ice_crit use comslope_mod, only: subslope_dist, nslope implicit none !======================================================================= ! ! Routine to check if the h2o ice criterion to stop the PEM is reached ! !======================================================================= ! Inputs !------- integer, intent(in) :: ngrid ! # of physical grid points real, dimension(ngrid), intent(in) :: cell_area ! Area of the cells real, dimension(ngrid,nslope), intent(in) :: h2o_ice ! Actual density of h2o ice real, intent(in) :: h2oice_ini_surf ! Initial surface of sublimating h2o ice logical, dimension(ngrid,nslope), intent(in) :: ini_h2oice_sublim ! Grid points where h2o ice was initially sublimating ! Outputs !-------- integer, intent(inout) :: stopPEM ! Stopping criterion code ! Locals ! ------ integer :: i, islope ! Loop real :: h2oice_now_surf ! Current surface of h2o ice !======================================================================= ! Computation of the present surface of h2o ice still sublimating h2oice_now_surf = 0. do i = 1,ngrid do islope = 1,nslope 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) enddo enddo ! Check of the criterion if (h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)) then stopPEM = 1 write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold" write(*,*) "h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit)", h2oice_now_surf < h2oice_ini_surf*(1. - h2o_ice_crit) write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf write(*,*) "Percentage of change accepted =", h2o_ice_crit*100 else if (h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)) then stopPEM = 1 write(*,*) "Reason of stopping: the surface of h2o ice sublimating reaches the threshold" write(*,*) "h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit)", h2oice_now_surf > h2oice_ini_surf*(1. + h2o_ice_crit) write(*,*) "Initial surface of h2o ice sublimating =", h2oice_ini_surf write(*,*) "Current surface of h2o ice sublimating =", h2oice_now_surf write(*,*) "Percentage of change accepted =", h2o_ice_crit*100 endif if (abs(h2oice_ini_surf) < 1.e-5) stopPEM = 0 END SUBROUTINE stopping_crit_h2o_ice !======================================================================= SUBROUTINE stopping_crit_co2(cell_area,co2ice_ini_surf,ini_co2ice_sublim,co2_ice,stopPEM,ngrid,global_avg_press_PCM,global_avg_press_new,nslope) use time_evol_mod, only: co2_ice_crit, ps_criterion use comslope_mod, only: subslope_dist implicit none !======================================================================= ! ! Routine to check if the co2 and pressure criteria to stop the PEM are reached ! !======================================================================= ! Inputs !------- integer, intent(in) :: ngrid, nslope ! # of grid physical grid points real, dimension(ngrid), intent(in) :: cell_area ! Area of the cells real, dimension(ngrid,nslope), intent(in) :: co2_ice ! Actual density of co2 ice real, intent(in) :: co2ice_ini_surf ! Initial surface of sublimatingco2 ice logical, dimension(ngrid,nslope), intent(in) :: ini_co2ice_sublim ! Grid points where co2 ice was initially sublimating real, intent(in) :: global_avg_press_PCM ! Planet average pressure from the PCM start files real, intent(in) :: global_avg_press_new ! Planet average pressure from the PEM computations ! Outputs !-------- integer, intent(inout) :: stopPEM ! Stopping criterion code ! Locals ! ------ integer :: i, islope ! Loop real :: co2ice_now_surf ! Current surface of co2 ice !======================================================================= ! Computation of the present surface of co2 ice still sublimating co2ice_now_surf = 0. do i = 1,ngrid do islope = 1,nslope 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) enddo enddo ! Check of the criterion if (co2ice_now_surf < co2ice_ini_surf*(1. - co2_ice_crit)) then stopPEM = 3 write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold" write(*,*) "co2ice_now_surf < co2ice_ini_surf*(1. - co2_ice_crit)", co2ice_now_surf < co2ice_ini_surf*(1. - co2_ice_crit) write(*,*) "Initial surface of co2 ice sublimating =", co2ice_ini_surf write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf write(*,*) "Percentage of change accepted =", co2_ice_crit*100. else if (co2ice_now_surf > co2ice_ini_surf*(1. + co2_ice_crit)) then stopPEM = 3 write(*,*) "Reason of stopping: the surface of co2 ice sublimating reaches the threshold" write(*,*) "co2ice_now_surf > co2ice_ini_surf*(1. + co2_ice_crit)", co2ice_now_surf > co2ice_ini_surf*(1. + co2_ice_crit) write(*,*) "Current surface of co2 ice sublimating =", co2ice_now_surf write(*,*) "Initial surface of co2 ice sublimating =", co2ice_ini_surf write(*,*) "Percentage of change accepted =", co2_ice_crit*100. endif if (abs(co2ice_ini_surf) < 1.e-5) stopPEM = 0 if (global_avg_press_new < global_avg_press_PCM*(1. - ps_criterion)) then stopPEM = 4 write(*,*) "Reason of stopping: the global pressure reaches the threshold" write(*,*) "global_avg_press_new < global_avg_press_PCM*(1. - ps_criterion)", global_avg_press_new < global_avg_press_PCM*(1. - ps_criterion) write(*,*) "Initial global pressure =", global_avg_press_PCM write(*,*) "Current global pressure =", global_avg_press_new write(*,*) "Percentage of change accepted =", ps_criterion*100. else if (global_avg_press_new > global_avg_press_PCM*(1. + ps_criterion)) then stopPEM = 4 write(*,*) "Reason of stopping: the global pressure reaches the threshold" write(*,*) "global_avg_press_new > global_avg_press_PCM*(1. + ps_criterion)", global_avg_press_new > global_avg_press_PCM*(1. + ps_criterion) write(*,*) "Initial global pressure =", global_avg_press_PCM write(*,*) "Current global pressure =", global_avg_press_new write(*,*) "Percentage of change accepted =", ps_criterion*100. endif END SUBROUTINE stopping_crit_co2 END MODULE