source: trunk/LMDZ.COMMON/libf/evolution/glaciers_mod.F90 @ 3567

Last change on this file since 3567 was 3553, checked in by jbclement, 4 weeks ago

PEM:
Addition of the shifting of the soil temperature profile to follow the surface evolution due to ice condensation/sublimation + small cleanings/improvements.
JBC

File size: 11.8 KB
Line 
1MODULE glaciers_mod
2
3implicit none
4
5! Flags for ice management
6logical :: h2oice_flow  ! True by default, to compute H2O ice flow. Read in "run_PEM.def"
7logical :: co2ice_flow  ! True by default, to compute CO2 ice flow. Read in "run_PEM.def"
8logical :: metam_h2oice ! False by default, to compute H2O ice metamorphism. Read in "run_PEM.def"
9logical :: metam_co2ice ! False by default, to compute CO2 ice metamorphism. Read in "run_PEM.def"
10
11! Thresholds for ice management
12real :: inf_h2oice_threshold   ! To consider the amount of H2O ice as an infinite reservoir [kg.m-2]
13real :: metam_h2oice_threshold ! To consider frost is becoming perennial H2O ice [kg.m-2]
14real :: metam_co2ice_threshold ! To consider frost is becoming perennial CO2 ice [kg.m-2]
15
16real, parameter :: rho_co2ice = 1650. ! Density of CO2 ice [kg.m-3]
17real, parameter :: rho_h2oice = 920.  ! Density of H2O ice [kg.m-3]
18
19!=======================================================================
20contains
21!=======================================================================
22
23SUBROUTINE flow_co2glaciers(timelen,ngrid,nslope,iflat,subslope_dist,def_slope_mean,vmr_co2_PEM,ps_PCM,global_avg_ps_PCM,global_avg_ps_PEM,co2ice,flag_co2flow,flag_co2flow_mesh)
24
25!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
26!!!
27!!! Purpose: Main for CO2 glaciers evolution: compute maximum thickness, and do
28!!!          the ice transfer
29!!!
30!!!
31!!! Author: LL
32!!!
33!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
34
35implicit none
36
37! Inputs:
38!--------
39integer,                           intent(in) :: timelen, ngrid, nslope, iflat ! number of time sample, physical points, subslopes, index of the flat subslope
40real,    dimension(ngrid,nslope),  intent(in) :: subslope_dist                 ! Physical points x Slopes: Distribution of the subgrid slopes
41real,    dimension(ngrid),         intent(in) :: def_slope_mean                ! Physical points: values of the sub grid slope angles
42real,    dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM                   ! Physical x Time field : VMR of co2 in the first layer [mol/mol]
43real,    dimension(ngrid,timelen), intent(in) :: ps_PCM                        ! Physical x Time field: surface pressure given by the PCM [Pa]
44real,                              intent(in) :: global_avg_ps_PCM             ! Global averaged surface pressure from the PCM [Pa]
45real,                              intent(in) :: global_avg_ps_PEM             ! global averaged surface pressure during the PEM iteration [Pa]
46
47! Ouputs:
48!--------
49real, dimension(ngrid,nslope), intent(inout) :: co2ice            ! Physical x Slope field: co2 ice on the subgrid slopes [kg/m^2]
50real, dimension(ngrid,nslope), intent(inout) :: flag_co2flow      ! flag to see if there is flow on the subgrid slopes
51real, dimension(ngrid),        intent(inout) :: flag_co2flow_mesh ! same but within the mesh
52! Local
53!------
54real, dimension(ngrid,nslope) :: Tcond ! Physical field: CO2 condensation temperature [K]
55real, dimension(ngrid,nslope) :: hmax  ! Physical x Slope field: maximum thickness for co2  glacier before flow
56
57write(*,*) "Flow of CO2 glaciers"
58call computeTcondCO2(timelen,ngrid,nslope,vmr_co2_PEM,ps_PCM,global_avg_ps_PCM,global_avg_ps_PEM,Tcond)
59call compute_hmaxglaciers(ngrid,nslope,iflat,def_slope_mean,Tcond,"co2",hmax)
60call transfer_ice_duringflow(ngrid,nslope,iflat, subslope_dist,def_slope_mean,hmax,Tcond,"co2",co2ice,flag_co2flow,flag_co2flow_mesh)
61
62END SUBROUTINE flow_co2glaciers
63
64!=======================================================================
65
66SUBROUTINE flow_h2oglaciers(timelen,ngrid,nslope,iflat,subslope_dist,def_slope_mean,Tice,h2oice,flag_h2oflow,flag_h2oflow_mesh)
67
68!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
69!!!
70!!! Purpose: Main for H2O glaciers evolution: compute maximum thickness, and do
71!!!          the ice transfer
72!!!
73!!!
74!!! Author: LL
75!!!
76!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
77
78implicit none
79
80! arguments
81! ---------
82
83! Inputs:
84integer,                       intent(in) :: timelen, ngrid, nslope, iflat ! number of time sample, physical points, subslopes, index of the flat subslope
85real, dimension(ngrid,nslope), intent(in) :: subslope_dist  ! Physical points x Slopes : Distribution of the subgrid slopes
86real, dimension(ngrid),        intent(in) :: def_slope_mean ! Slopes: values of the sub grid slope angles
87real, dimension(ngrid,nslope), intent(in) :: Tice           ! Ice Temperature [K]
88! Ouputs:
89real, dimension(ngrid,nslope), intent(inout) :: h2oice            ! Physical x Slope field: co2 ice on the subgrid slopes [kg/m^2]
90real, dimension(ngrid,nslope), intent(inout) :: flag_h2oflow      ! flag to see if there is flow on the subgrid slopes
91real, dimension(ngrid),        intent(inout) :: flag_h2oflow_mesh ! same but within the mesh
92! Local
93real, dimension(ngrid,nslope) :: hmax ! Physical x Slope field: maximum thickness for co2  glacier before flow
94
95write(*,*) "Flow of H2O glaciers"
96call compute_hmaxglaciers(ngrid,nslope,iflat,def_slope_mean,Tice,"h2o",hmax)
97call transfer_ice_duringflow(ngrid,nslope,iflat, subslope_dist,def_slope_mean,hmax,Tice,"h2o",h2oice,flag_h2oflow,flag_h2oflow_mesh)
98
99END SUBROUTINE flow_h2oglaciers
100
101!=======================================================================
102
103SUBROUTINE compute_hmaxglaciers(ngrid,nslope,iflat,def_slope_mean,Tice,name_ice,hmax)
104
105use ice_table_mod, only: rho_ice
106use abort_pem_mod, only: abort_pem
107#ifndef CPP_STD
108    use comcstfi_h,   only: pi, g
109#else
110    use comcstfi_mod, only: pi, g
111#endif
112
113!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
114!!!
115!!! Purpose: Compute the maximum thickness of CO2 and H2O glaciers given a slope angle before initating flow
116!!!
117!!! Author: LL,based on  work by A.Grau Galofre (LPG) and Isaac Smith (JGR Planets 2022)
118!!!
119!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
120
121implicit none
122
123! Inputs
124! ------
125integer,                       intent(in) :: ngrid, nslope  ! # of grid points and subslopes
126integer,                       intent(in) :: iflat          ! index of the flat subslope
127real, dimension(nslope),       intent(in) :: def_slope_mean ! Slope field: Values of the subgrid slope angles [deg]
128real, dimension(ngrid,nslope), intent(in) :: Tice           ! Physical field:  ice temperature [K]
129character(3),                  intent(in) :: name_ice       ! Nature of ice
130! Outputs
131! -------
132real, dimension(ngrid,nslope), intent(out) :: hmax ! Physical grid x Slope field: maximum  thickness before flaw [m]
133! Local
134! -----
135real    :: tau_d      ! characteristic basal drag, understood as the stress that an ice mass flowing under its weight balanced by viscosity. Value obtained from I.Smith
136integer :: ig, islope ! loop variables
137real    :: slo_angle
138
139select case (trim(adjustl(name_ice)))
140    case('h2o')
141        tau_d = 1.e5
142    case('co2')
143        tau_d = 5.e3
144    case default
145        call abort_pem("compute_hmaxglaciers","Type of ice not known!",1)
146end select
147
148do ig = 1,ngrid
149    do islope = 1,nslope
150        if (islope == iflat) then
151            hmax(ig,islope) = 1.e8
152        else
153            slo_angle = abs(def_slope_mean(islope)*pi/180.)
154            hmax(ig,islope) = tau_d/(rho_ice(Tice(ig,islope),name_ice)*g*slo_angle)
155        endif
156    enddo
157enddo
158
159END SUBROUTINE compute_hmaxglaciers
160
161!=======================================================================
162
163SUBROUTINE transfer_ice_duringflow(ngrid,nslope,iflat,subslope_dist,def_slope_mean,hmax,Tice,name_ice,qice,flag_flow,flag_flowmesh)
164!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
165!!!
166!!! Purpose: Transfer the excess of ice from one subslope to another
167!!!          No transfer between mesh at the time
168!!! Author: LL
169!!!
170!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
171
172use ice_table_mod, only: rho_ice
173use abort_pem_mod, only: abort_pem
174#ifndef CPP_STD
175    use comcstfi_h,   only: pi
176#else
177    use comcstfi_mod, only: pi
178#endif
179
180implicit none
181
182! Inputs
183!-------
184integer,                       intent(in) :: ngrid, nslope  ! # of physical points and subslope
185integer,                       intent(in) :: iflat          ! index of the flat subslope
186real, dimension(ngrid,nslope), intent(in) :: subslope_dist  ! Distribution of the subgrid slopes within the mesh
187real, dimension(nslope),       intent(in) :: def_slope_mean ! values of the subgrid slopes
188real, dimension(ngrid,nslope), intent(in) :: hmax           ! maximum height of the  glaciers before initiating flow [m]
189real, dimension(ngrid,nslope), intent(in) :: Tice           ! Ice temperature[K]
190character(3),                  intent(in) :: name_ice       ! Nature of the ice
191! Outputs
192!--------
193real, dimension(ngrid,nslope), intent(inout) :: qice          ! CO2 in the subslope [kg/m^2]
194real, dimension(ngrid,nslope), intent(inout) :: flag_flow     ! boolean to check if there is flow on a subgrid slope
195real, dimension(ngrid),        intent(inout) :: flag_flowmesh ! boolean to check if there is flow in the mesh
196! Local
197!------
198integer :: ig, islope ! loop
199integer :: iaval      ! ice will be transfered here
200
201do ig = 1,ngrid
202    do islope = 1,nslope
203        if (islope /= iflat) then ! ice can be infinite on flat ground
204! First: check that CO2 ice must flow (excess of ice on the slope), ice can accumulate infinitely  on flat ground
205            if (qice(ig,islope) >= rho_ice(Tice(ig,islope),'h2o')*hmax(ig,islope)*cos(pi*def_slope_mean(islope)/180.)) then
206! Second: determine the flatest slopes possible:
207                if (islope > iflat) then
208                    iaval=islope-1
209                else
210                    iaval = islope + 1
211                endif
212                do while (iaval /= iflat .and. subslope_dist(ig,iaval) == 0)
213                    if (iaval > iflat) then
214                        iaval = iaval - 1
215                    else
216                        iaval = iaval + 1
217                    endif
218                enddo
219                qice(ig,iaval) = qice(ig,iaval) + (qice(ig,islope) - rho_ice(Tice(ig,islope),'h2o')*hmax(ig,islope)*cos(pi*def_slope_mean(islope)/180.)) &
220                               *subslope_dist(ig,islope)/subslope_dist(ig,iaval)*cos(pi*def_slope_mean(iaval)/180.)/cos(pi*def_slope_mean(islope)/180.)
221
222                qice(ig,islope) = rho_ice(Tice(ig,islope),'h2o')*hmax(ig,islope)*cos(pi*def_slope_mean(islope)/180.)
223                flag_flow(ig,islope) = 1.
224                flag_flowmesh(ig) = 1.
225            endif ! co2ice > hmax
226        endif ! iflat
227    enddo !islope
228enddo !ig
229
230END SUBROUTINE
231
232!=======================================================================
233
234SUBROUTINE computeTcondCO2(timelen,ngrid,nslope,vmr_co2_PEM,ps_PCM,global_avg_ps_PCM,global_avg_ps_PEM,Tcond)
235!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
236!!!
237!!! Purpose: Compute CO2 condensation temperature
238!!!
239!!! Author: LL
240!!!
241!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
242
243use constants_marspem_mod, only: alpha_clap_co2,beta_clap_co2
244
245implicit none
246
247! arguments:
248! ----------
249
250! INPUT
251integer,                        intent(in) :: timelen, ngrid, nslope ! # of timesample, physical points, subslopes
252real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM            ! Physical points x times field: VMR of CO2 in the first layer [mol/mol]
253real, dimension(ngrid,timelen), intent(in) :: ps_PCM                 ! Physical points x times field: surface pressure in the PCM [Pa]
254real,                           intent(in) :: global_avg_ps_PCM      ! Global averaged surfacepressure in the PCM [Pa]
255real,                           intent(in) :: global_avg_ps_PEM      ! Global averaged surface pressure computed during the PEM iteration
256
257! OUTPUT
258real, dimension(ngrid,nslope), intent(out) :: Tcond ! Physical points: condensation temperature of CO2, yearly averaged
259
260! LOCAL
261integer :: ig, it ! For loop
262real    :: ave    ! Intermediate to compute average
263
264do ig = 1,ngrid
265    ave = 0
266    do it = 1,timelen
267        ave = ave + beta_clap_co2/(alpha_clap_co2 - log(vmr_co2_PEM(ig,it)*ps_PCM(ig,it)*global_avg_ps_PCM/global_avg_ps_PEM/100))
268    enddo
269    Tcond(ig,:) = ave/timelen
270enddo
271
272END SUBROUTINE computeTcondCO2
273
274END MODULE glaciers_mod
Note: See TracBrowser for help on using the repository browser.