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

Last change on this file since 3525 was 3498, checked in by jbclement, 3 weeks ago

PEM:

  • Correction of the variable name for the ice table depth in "pemetat0.F90". So it is now got as intended from the "startpem.nc" file;
  • Renaming of the tendencies in the PEM with the prefix 'd_' instead of 'tend_';
  • Modification of the PEM time step type from integer to real. As a consequence, all time variables are now of real type. This change adds the possibility to consider fractions of year as time step.

JBC

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