source: trunk/LMDZ.COMMON/libf/evolution/evol_ice_mod.F90 @ 3949

Last change on this file since 3949 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: 6.4 KB
RevLine 
[3149]1MODULE evol_ice_mod
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8
[3553]9SUBROUTINE evol_co2_ice(ngrid,nslope,co2_ice,d_co2ice,zshift_surf)
[3149]10
[3498]11use time_evol_mod, only: dt
[3553]12use glaciers_mod,  only: rho_co2ice
[3149]13
14implicit none
15
16!=======================================================================
17!
18! Routine to compute the evolution of CO2 ice
19!
20!=======================================================================
[3553]21! arguments:
22! ----------
23! INPUT
[3149]24integer, intent(in) :: ngrid, nslope ! # of grid points, # of subslopes
25
[3553]26! OUTPUT
27real, dimension(ngrid,nslope), intent(inout) :: co2_ice     ! Previous and actual density of CO2 ice
28real, dimension(ngrid,nslope), intent(inout) :: d_co2ice    ! Evolution of perennial ice over one year
29real, dimension(ngrid,nslope), intent(out) :: zshift_surf ! Elevation shift for the surface [m]
[3149]30
[3553]31! local:
32! ------
[3368]33real, dimension(ngrid,nslope) :: co2_ice_old ! Old density of CO2 ice
[3149]34!=======================================================================
35! Evolution of CO2 ice for each physical point
[3603]36write(*,*) '> Evolution of CO2 ice'
[3367]37
38co2_ice_old = co2_ice
[3498]39co2_ice = co2_ice + d_co2ice*dt
[3367]40where (co2_ice < 0.)
41    co2_ice = 0.
[3498]42    d_co2ice = -co2_ice_old/dt
[3149]43end where
[3367]44
[3553]45zshift_surf = d_co2ice*dt/rho_co2ice
46
[3149]47END SUBROUTINE evol_co2_ice
48
49!=======================================================================
50
[3554]51SUBROUTINE evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbed,delta_h2o_icetablesublim,h2o_ice,d_h2oice,zshift_surf,stopPEM)
[3149]52
[3938]53use time_evol_mod,     only: dt
54use glaciers_mod,      only: rho_h2oice
55use stopping_crit_mod, only: stopping_crit_h2o
[3149]56
57implicit none
58
59!=======================================================================
60!
[3368]61! Routine to compute the evolution of h2o ice
[3149]62!
63!=======================================================================
[3553]64! arguments:
65! ----------
66! INPUT
[3149]67integer,                intent(in) :: ngrid, nslope            ! # of grid points, # of subslopes
68real, dimension(ngrid), intent(in) :: cell_area                ! Area of each mesh grid (m^2)
[3556]69real, dimension(ngrid), intent(in) :: delta_h2o_adsorbed       ! Mass of H2O adsorbed/desorbded in the soil (kg/m^2)
[3149]70real, dimension(ngrid), intent(in) :: delta_h2o_icetablesublim ! Mass of H2O that have condensed/sublimated at the ice table (kg/m^2)
71
[3553]72! OUTPUT
[3938]73real, dimension(ngrid,nslope), intent(inout) :: h2o_ice  ! H2O ice (kg/m^2)
74real, dimension(ngrid,nslope), intent(inout) :: d_h2oice ! Tendency of H2O ice (kg/m^2/year)
[3571]75integer,                       intent(inout) :: stopPEM  ! Stopping criterion code
[3553]76real, dimension(ngrid,nslope), intent(out) :: zshift_surf ! Elevation shift for the surface [m]
[3149]77
[3553]78! local:
79! ------
[3938]80integer                       :: i, islope                                                ! Loop variables
81real                          :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! Variables to balance H2O
82real, dimension(ngrid,nslope) :: d_h2oice_new                                             ! Tendencies computed to keep balance
[3149]83!=======================================================================
[3603]84write(*,*) '> Evolution of H2O ice'
[3308]85
[3938]86if (ngrid == 1) then ! In 1D
87    h2o_ice = h2o_ice + d_h2oice*dt
88    zshift_surf = d_h2oice*dt/rho_h2oice
89else ! In 3D
90    call 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)
91    if (stopPEM > 0) return
[3149]92
[3938]93    call balance_h2oice_reservoirs(ngrid,nslope,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new)
94    h2o_ice = h2o_ice + d_h2oice_new*dt
95    zshift_surf = d_h2oice_new*dt/rho_h2oice
96endif
[3149]97
[3938]98END SUBROUTINE evol_h2o_ice
99
100!=======================================================================
101
102SUBROUTINE balance_h2oice_reservoirs(ngrid,nslope,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new)
103
104use time_evol_mod, only: dt
105
106implicit none
107
108!=======================================================================
109!
110! Routine to balance the H2O flux from and into the atmosphere accross reservoirs
111!
112!=======================================================================
113! arguments:
114! ----------
115! INPUT
116integer,                       intent(in) :: ngrid, nslope ! # of grid points, # of subslopes
117real,                          intent(in) :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! ! Variables to conserve H2O
118real, dimension(ngrid,nslope), intent(in) :: h2o_ice ! H2O ice (kg/m^2)
119
120! OUTPUT
121real, dimension(ngrid,nslope), intent(inout) :: d_h2oice ! Tendency of H2O ice (kg/m^2/year)
122real, dimension(ngrid,nslope), intent(out) :: d_h2oice_new ! Adjusted tendencies to keep balance between donor and recipient reservoirs
123
124! local:
125! ------
126integer :: i, islope
127real    :: S_target, S_target_subl_h2oice, S_target_cond_h2oice, S_ghostice, d_target
128!=======================================================================
129S_target = (S_atm_2_h2o + S_h2o_2_atm)/2.
130S_target_cond_h2oice = S_atm_2_h2oice + S_target - S_atm_2_h2o
131S_target_subl_h2oice = S_h2oice_2_atm + S_target - S_h2o_2_atm
132
133d_h2oice_new = 0.
134S_ghostice = 0.
135do i = 1,ngrid
136    do islope = 1,nslope
137        if (d_h2oice(i,islope) > 0.) then ! Condensing
138            d_h2oice_new(i,islope) = d_h2oice_new(i,islope)*S_target_cond_h2oice/S_atm_2_h2oice
139        else if (d_h2oice(i,islope) < 0.) then ! Sublimating
140            d_target = d_h2oice(i,islope)*S_target_subl_h2oice/S_h2oice_2_atm
141            if (abs(d_target*dt) < h2o_ice(i,islope)) then ! Enough ice to sublimate everything
142                d_h2oice_new(i,islope) = d_target
143            else ! Not enough ice to sublimate everything
144                ! We sublimate what we can
145                d_h2oice_new(i,islope) = h2o_ice(i,islope)/dt
146                ! It means the tendency is zero next time
[3498]147                d_h2oice(i,islope) = 0.
[3938]148                ! We compute the amount of H2O ice that we could not make sublimate
149                S_ghostice = S_ghostice + abs(d_target*dt) - h2o_ice(i,islope)
[3149]150            endif
[3938]151        endif
[3149]152    enddo
[3938]153enddo
[3149]154
[3938]155! We need to remove this ice unable to sublimate from places where ice condenseds in order to keep balance
156where (d_h2oice_new > 0.) d_h2oice_new = d_h2oice_new*(S_target_cond_h2oice - S_ghostice)/S_target_cond_h2oice
[3149]157
[3938]158END SUBROUTINE balance_h2oice_reservoirs
[3553]159
[3149]160END MODULE evol_ice_mod
Note: See TracBrowser for help on using the repository browser.