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

Last change on this file since 3149 was 3149, checked in by jbclement, 19 months ago

PEM:

  • Simplification of the algorithm managing the stopping criteria;
  • Complete rework of the ice management in the PEM (H2O & CO2);

    Subroutines to evolve the H2O and CO2 ice are now in the same module "evol_ice_mod.F90".
    Tendencies are computed from the variation of "ice + frost" between the 2 PCM runs.
    Evolving ice in the PEM is now called 'h2o_ice' or 'co2_ice' (not anymore in 'qsurf' and free of 'water_reservoir').
    Default value 'ini_h2o_bigreservoir' (= 10 m) initializes the H2O ice of the first PEM run where there is 'watercap'. For the next PEM runs, initialization is done with the value kept in "startpem.nc". CO2 ice is taken from 'perennial_co2ice' of the PCM (paleoclimate flag must be true).
    Simplification of the condition to compute the surface ice cover needed for the stopping criteria.
    Frost ('qsurf') is not evolved by the PEM and given back to the PCM.
    New default threshold value 'inf_h2oice_threshold' (= 2 m) to decide at the end of the PEM run if the H2O ice should be 'watercap' or not for the next PCM runs. If H2O ice cannot be 'watercap', then the remaining H2O ice is transferred to the frost ('qsurf').

  • Renaming of variables/subroutines for clarity;
  • Some cleanings throughout the code;
  • Small updates in files of the deftank.

JBC

File size: 6.5 KB
Line 
1MODULE evol_ice_mod
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8
9SUBROUTINE evol_co2_ice(ngrid,nslope,co2_ice,tend_co2_ice)
10
11use time_evol_mod, only: dt_pem
12
13implicit none
14
15!=======================================================================
16!
17! Routine to compute the evolution of CO2 ice
18!
19!=======================================================================
20!   arguments:
21!   ----------
22!   INPUT
23integer, intent(in) :: ngrid, nslope ! # of grid points, # of subslopes
24
25!   OUTPUT
26real, dimension(ngrid,nslope), intent(inout) :: co2_ice      ! Previous and actual density of CO2 ice
27real, dimension(ngrid,nslope), intent(inout) :: tend_co2_ice ! Evolution of perennial ice over one year
28
29!   local:
30!   ----
31
32!=======================================================================
33! Evolution of CO2 ice for each physical point
34write(*,*) 'Evolution of co2 ice'
35co2_ice = co2_ice + tend_co2_ice*dt_pem
36where (co2_ice < 0.)
37    co2_ice = 0.
38    tend_co2_ice = 0.
39end where
40
41END SUBROUTINE evol_co2_ice
42
43!=======================================================================
44
45SUBROUTINE evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbded,delta_h2o_icetablesublim,h2o_ice,tend_h2o_ice,stopPEM)
46
47use time_evol_mod,          only: dt_pem
48use comslope_mod,           only: subslope_dist, def_slope_mean
49use criterion_pem_stop_mod, only: criterion_waterice_stop
50
51#ifndef CPP_STD
52    use comcstfi_h,   only: pi
53#else
54    use comcstfi_mod, only: pi
55#endif
56
57implicit none
58
59!=======================================================================
60!
61! Routine that compute the evolution of the h2o ice
62!
63!=======================================================================
64!   arguments:
65!   ----------
66!   INPUT
67integer,                intent(in) :: ngrid, nslope            ! # of grid points, # of subslopes
68real, dimension(ngrid), intent(in) :: cell_area                ! Area of each mesh grid (m^2)
69real, dimension(ngrid), intent(in) :: delta_h2o_adsorbded      ! Mass of H2O adsorbded/desorbded in the soil (kg/m^2)
70real, dimension(ngrid), intent(in) :: delta_h2o_icetablesublim ! Mass of H2O that have condensed/sublimated at the ice table (kg/m^2)
71
72!   OUTPUT
73real, dimension(ngrid,nslope), intent(inout) :: h2o_ice      ! Previous and actual density of h2o ice (kg/m^2)
74real, dimension(ngrid,nslope), intent(inout) :: tend_h2o_ice ! Evolution of perennial ice over one year (kg/m^2/year)
75integer,                       intent(inout) :: stopPEM      ! Stopping criterion code
76
77!   local:
78!   ------
79integer                       :: i, j, islope                                        ! Loop variables
80real                          :: pos_tend, neg_tend, real_coefficient, negative_part ! Variable to conserve h2o
81real, dimension(ngrid,nslope) :: new_tendencies                                      ! Tendencies computed in order to conserve h2o ice on the surface, only exchange between surface are done
82!=======================================================================
83if (ngrid /= 1) then ! Not in 1D
84    ! We compute the amount of condensing and sublimating h2o ice
85    pos_tend = 0.
86    neg_tend = 0.
87    do i = 1,ngrid
88        if (delta_h2o_adsorbded(i) > 0.) then
89            pos_tend = pos_tend + delta_h2o_adsorbded(i)*cell_area(i)
90        else
91            neg_tend = neg_tend + delta_h2o_adsorbded(i)*cell_area(i)
92        endif
93        if (delta_h2o_icetablesublim(i) > 0.) then
94            pos_tend = pos_tend + delta_h2o_icetablesublim(i)*cell_area(i)
95        else
96            neg_tend = neg_tend + delta_h2o_icetablesublim(i)*cell_area(i)
97        endif
98        do islope = 1,nslope
99            if (h2o_ice(i,islope) > 0.) then
100                if (tend_h2o_ice(i,islope) > 0.) then
101                    pos_tend = pos_tend + tend_h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
102                else
103                    neg_tend = neg_tend - tend_h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
104                endif
105            endif
106        enddo
107    enddo
108    ! We adapt the tendencies to conserve h2o and do only exchange between grid points
109    if (neg_tend > pos_tend .and. pos_tend > 0.) then ! More sublimation on the planet than condensation
110        where (tend_h2o_ice < 0.) ! We lower the sublimating rate by a coefficient
111            new_tendencies = tend_h2o_ice*pos_tend/neg_tend
112        elsewhere ! We don't change the condensing rate
113            new_tendencies = tend_h2o_ice
114        end where
115    else if (neg_tend < pos_tend .and. neg_tend > 0.) then ! More condensation on the planet than sublimation
116        where (tend_h2o_ice < 0.) ! We don't change the sublimating rate
117            new_tendencies = tend_h2o_ice
118        elsewhere ! We lower the condensing rate by a coefficient
119            new_tendencies = tend_h2o_ice*neg_tend/pos_tend
120        end where
121    else if (abs(pos_tend) < 1.e-10 .or. abs(neg_tend) == 1.e-10) then
122        write(*,*) "Reason of stopping: there is no sublimating or condensing h20 ice!"
123        write(*,*) "Tendencies on ice sublimating =", neg_tend
124        write(*,*) "Tendencies on ice increasing =", pos_tend
125        write(*,*) "This can be due to the absence of h2o ice in the PCM run!"
126        stopPEM = 2
127        new_tendencies = 0.
128    endif
129
130    ! Evolution of the h2o ice for each physical point
131    h2o_ice = h2o_ice + new_tendencies*dt_pem
132
133    ! We compute the amount of h2o that is sublimated in excess
134    negative_part = 0.
135    do i = 1,ngrid
136        do islope = 1,nslope
137            if (h2o_ice(i,islope) < 0.) then
138                negative_part = negative_part - h2o_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
139                h2o_ice(i,islope) = 0.
140                tend_h2o_ice(i,islope) = 0.
141            endif
142        enddo
143    enddo
144
145    ! We compute a coefficient by which we should remove the ice that has been added to places
146    ! even if this ice was contributing to an unphysical negative amount of ice at other places
147    if (abs(pos_tend) < 1.e-10) then
148        real_coefficient = 0.
149    else
150        real_coefficient = negative_part/pos_tend
151    endif
152    ! In the place of accumulation of ice, we remove a bit of ice in order to conserve h2o
153    where (new_tendencies > 0.) h2o_ice = h2o_ice - new_tendencies*real_coefficient*dt_pem*cos(def_slope_mean(islope)*pi/180.)
154else ! ngrid == 1, i.e. in 1D
155    h2o_ice = h2o_ice + tend_h2o_ice*dt_pem
156endif
157
158END SUBROUTINE evol_h2o_ice
159
160END MODULE evol_ice_mod
Note: See TracBrowser for help on using the repository browser.