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

Last change on this file since 3554 was 3554, checked in by jbclement, 6 days ago

PEM:
Follow-up of previous commit (r3553).
JBC

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