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

Last change on this file since 3566 was 3556, checked in by jbclement, 8 weeks ago

PEM:
As intended, years computed by the PEM runs are now the only ones to be counted for the duration of the PEM simulation. The possibility to count in addition years computed by the PCM runs is left as an option of "launchPEM.sh" with the variable 'counting' (0 = "only PEM runs count"; any other values = "PCM runs are taken into account") + several small corrections/improvements in the launching scripts.
JBC

File size: 6.9 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
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
[3553]53use time_evol_mod, only: dt
54use comslope_mod,  only: subslope_dist, def_slope_mean
55use glaciers_mod,  only: rho_h2oice
[3149]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!
[3368]67! Routine to compute the evolution of h2o ice
[3149]68!
69!=======================================================================
[3553]70! arguments:
71! ----------
72! INPUT
[3149]73integer,                intent(in) :: ngrid, nslope            ! # of grid points, # of subslopes
74real, dimension(ngrid), intent(in) :: cell_area                ! Area of each mesh grid (m^2)
[3556]75real, dimension(ngrid), intent(in) :: delta_h2o_adsorbed       ! Mass of H2O adsorbed/desorbded in the soil (kg/m^2)
[3149]76real, dimension(ngrid), intent(in) :: delta_h2o_icetablesublim ! Mass of H2O that have condensed/sublimated at the ice table (kg/m^2)
77
[3553]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]
[3149]83
[3553]84! local:
85! ------
[3206]86integer                       :: i, islope                                           ! Loop variables
[3149]87real                          :: pos_tend, neg_tend, real_coefficient, negative_part ! Variable to conserve h2o
[3553]88real, dimension(ngrid,nslope) :: new_tend                                            ! Tendencies computed in order to conserve h2o ice on the surface, only exchange between surface are done
[3149]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
[3554]95        if (delta_h2o_adsorbed(i) > 0.) then
96            pos_tend = pos_tend + delta_h2o_adsorbed(i)*cell_area(i)
[3149]97        else
[3554]98            neg_tend = neg_tend + delta_h2o_adsorbed(i)*cell_area(i)
[3149]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
[3498]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.)
[3149]109                else
[3498]110                    neg_tend = neg_tend - d_h2oice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(def_slope_mean(islope)*pi/180.)
[3149]111                endif
112            endif
113        enddo
114    enddo
[3308]115
[3498]116    new_tend = 0.
[3308]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!"
[3149]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
[3308]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
[3498]126            where (d_h2oice < 0.) ! We lower the sublimating rate by a coefficient
127                new_tend = d_h2oice*pos_tend/neg_tend
[3308]128            elsewhere ! We don't change the condensing rate
[3498]129                new_tend = d_h2oice
[3308]130            end where
131        else if (neg_tend < pos_tend .and. neg_tend > 0.) then ! More condensation on the planet than sublimation
[3498]132            where (d_h2oice < 0.) ! We don't change the sublimating rate
133                new_tend = d_h2oice
[3308]134            elsewhere ! We lower the condensing rate by a coefficient
[3498]135                new_tend = d_h2oice*neg_tend/pos_tend
[3308]136            end where
137        endif
[3149]138    endif
139
140    ! Evolution of the h2o ice for each physical point
[3498]141    h2o_ice = h2o_ice + new_tend*dt
[3149]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.
[3498]150                d_h2oice(i,islope) = 0.
[3149]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
[3308]163    do islope = 1,nslope
164        do i = 1,ngrid
[3498]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.)
[3308]166        enddo
167     enddo
[3149]168else ! ngrid == 1, i.e. in 1D
[3498]169    h2o_ice = h2o_ice + d_h2oice*dt
[3149]170endif
171
[3553]172zshift_surf = d_h2oice*dt/rho_h2oice
173
[3149]174END SUBROUTINE evol_h2o_ice
175
176END MODULE evol_ice_mod
Note: See TracBrowser for help on using the repository browser.