1 | !> \file prescribe-H-i2s_mod.f90 |
---|
2 | !! |
---|
3 | !! Determine mask and thickness for prescription of H |
---|
4 | !! include ice2sea prescribed retreat (Cat april 2012) |
---|
5 | !< |
---|
6 | |
---|
7 | !> \namespace prescribe_H |
---|
8 | !! |
---|
9 | !! Determine mask and thickness for prescription of H |
---|
10 | !! |
---|
11 | !! \author CatRitz |
---|
12 | !! \date june 2009 |
---|
13 | !! |
---|
14 | !! @note use module3D_phy |
---|
15 | !! @note use param_phy_mod |
---|
16 | !< |
---|
17 | |
---|
18 | module prescribe_H |
---|
19 | |
---|
20 | use module3D_phy |
---|
21 | use param_phy_mod |
---|
22 | use interface_input |
---|
23 | |
---|
24 | ! pour les reculs de type ice2sea |
---|
25 | use toy_retreat_mod |
---|
26 | |
---|
27 | implicit none |
---|
28 | ! real,dimension(nx,ny) :: Hp !< H value if prescribed |
---|
29 | ! real,dimension(nx,ny) :: Hp0 !< H value if prescribed (reference value) |
---|
30 | ! real,dimension(nx,ny) :: Delta_H !< Delta_H value if prescribed |
---|
31 | ! integer,dimension(nx,ny) :: i_delta_H !< 1 if Delta_H is prescribed on this node, else 0 |
---|
32 | ! integer,dimension(nx,ny) :: i_Hp !< 1 if H is prescribed on this node, else 0 |
---|
33 | ! integer,dimension(nx,ny) :: i_Hp0 !< i_hp mask reference value does not change with time |
---|
34 | ! integer,dimension(nx,ny) :: MK_gl0 !< mask grounding line initial |
---|
35 | ! integer,dimension(nx,ny) :: MK_flot0 !< mask float initial |
---|
36 | |
---|
37 | ! pour grounding line retreat, ice2sea |
---|
38 | ! nouvelle version interactive pour imposer des contraintes variees. |
---|
39 | ! sanity check, temps de depart en fonction des bassins, .... |
---|
40 | |
---|
41 | |
---|
42 | ! nouvelle version recul grounding line pour ice2sea |
---|
43 | ! ------------------------------------------------------ |
---|
44 | |
---|
45 | |
---|
46 | |
---|
47 | contains |
---|
48 | |
---|
49 | !______________________________________________________________________________________ |
---|
50 | !> initialise prescribe H |
---|
51 | !! calculate the initial grounding line position |
---|
52 | !! @return MK_flot0 and Mk_gl0 |
---|
53 | |
---|
54 | subroutine init_prescribe_H |
---|
55 | |
---|
56 | implicit none |
---|
57 | integer :: voisin !< pour test sur des masques |
---|
58 | |
---|
59 | if (itracebug.eq.1) call tracebug(' Entree dans routine init_prescribe_H') |
---|
60 | |
---|
61 | ! keep the mask of floating points |
---|
62 | where (flot(:,:)) |
---|
63 | MK_flot0(:,:)= 1 |
---|
64 | elsewhere |
---|
65 | MK_flot0(:,:)= 0 |
---|
66 | end where |
---|
67 | |
---|
68 | |
---|
69 | ! determine the initial grounding line |
---|
70 | |
---|
71 | MK_gl0(:,:)=0 |
---|
72 | HP0(:,:)=H0(:,:) |
---|
73 | |
---|
74 | |
---|
75 | do j=2,ny-1 |
---|
76 | do i=2,nx-1 |
---|
77 | voisin = MK_flot0(i-1,j) + MK_flot0(i+1,j) + MK_flot0(i,j-1) + MK_flot0(i,j+1) |
---|
78 | if ((MK_flot0(i,j).eq.0).and.(voisin.gt.0)) then ! grounded and at least one neighbour floating |
---|
79 | MK_gl0(i,j) = 1 ! tagged mask grounding line |
---|
80 | |
---|
81 | end if |
---|
82 | end do |
---|
83 | end do |
---|
84 | |
---|
85 | |
---|
86 | call prescribe_grid_boundary_0 |
---|
87 | |
---|
88 | |
---|
89 | ! Dans le cas du Groenland ice2sea, mk_init = 3 -> noeuds qui ne bougent pas |
---|
90 | |
---|
91 | if (geoplace(1:4).eq.'GI2S') then |
---|
92 | |
---|
93 | do j = 1,ny |
---|
94 | do i = 1,nx |
---|
95 | if (mk_init(i,j).eq.3) then |
---|
96 | i_hp0 (i,j) = 1 |
---|
97 | Hp0 (i,j) = H0(i,j) |
---|
98 | endif |
---|
99 | end do |
---|
100 | end do |
---|
101 | end if |
---|
102 | |
---|
103 | |
---|
104 | if (igrdline.eq.3) then |
---|
105 | where (flot(:,:)) |
---|
106 | mk_gr(:,:) = 0 |
---|
107 | elsewhere |
---|
108 | mk_gr(:,:) = 1 |
---|
109 | end where |
---|
110 | call init_proto_recul |
---|
111 | end if |
---|
112 | |
---|
113 | end subroutine init_prescribe_H |
---|
114 | !______________________________________________________________________________________ |
---|
115 | !> function prescribe_present_H_gl |
---|
116 | !! calculate the initial grounding line position |
---|
117 | !! @return i_hp(:,:) and hp(:,:) |
---|
118 | |
---|
119 | subroutine prescribe_present_H_gl |
---|
120 | |
---|
121 | implicit none |
---|
122 | |
---|
123 | if (itracebug.eq.1) call tracebug(' Entree dans routine prescribe_present_H_gl') |
---|
124 | |
---|
125 | where ((MK_flot0(:,:).eq. 1).or.(MK_gl0(:,:) .eq. 1)) ! floating or grounding line |
---|
126 | i_hp(:,:) = 1 ! thickness prescribed to present value |
---|
127 | hp(:,:) = Hp0(:,:) |
---|
128 | end where |
---|
129 | if (itracebug.eq.1) call tracebug(' fin prescribe_present_H_gl') |
---|
130 | |
---|
131 | |
---|
132 | end subroutine prescribe_present_H_gl |
---|
133 | |
---|
134 | !______________________________________________________________________________________ |
---|
135 | !> function prescribe_present_H_gl_bmelt |
---|
136 | !! calculate the initial grounding line position |
---|
137 | !! only grounding line points are prescribed to present value |
---|
138 | !! @return i_hp(:,:) and hp(:,:) |
---|
139 | |
---|
140 | subroutine prescribe_present_H_gl_bmelt |
---|
141 | |
---|
142 | implicit none |
---|
143 | |
---|
144 | if (itracebug.eq.1) call tracebug(' Entree dans routine prescribe_present_H_gl_bmelt') |
---|
145 | |
---|
146 | ! where (MK_gl0(:,:) .eq. 1) ! grounding line only !cdc pour calcule fonte basale |
---|
147 | ! i_hp(:,:) = 1 ! thickness prescribed to present value |
---|
148 | ! hp(:,:) = Hp0(:,:) |
---|
149 | ! end where |
---|
150 | i_hp(:,:) = 0 |
---|
151 | if (itracebug.eq.1) call tracebug(' fin prescribe_present_H_gl_bmelt') |
---|
152 | |
---|
153 | |
---|
154 | end subroutine prescribe_present_H_gl_bmelt |
---|
155 | |
---|
156 | !______________________________________________________________________________________ |
---|
157 | !> function prescribe_fixed_points |
---|
158 | !! fixes thickness for some points |
---|
159 | !! @return i_hp(:,:) and hp(:,:) |
---|
160 | |
---|
161 | subroutine prescribe_fixed_points |
---|
162 | |
---|
163 | if (itracebug.eq.1) call tracebug(' Entree dans routine prescribe_fixed_points') |
---|
164 | where (i_hp0(:,:).eq.1) ! les points i_hp0 le sont pour tout le run |
---|
165 | i_hp(:,:) = i_hp0(:,:) |
---|
166 | Hp(:,:) = Hp0(:,:) |
---|
167 | end where |
---|
168 | |
---|
169 | if (itracebug.eq.1) call tracebug(' fin prescribe_fixed_points') |
---|
170 | |
---|
171 | |
---|
172 | end subroutine prescribe_fixed_points |
---|
173 | |
---|
174 | !______________________________________________________________________________________ |
---|
175 | !> function prescribe_paleo_gl_shelf |
---|
176 | !! calculate the grounding line thickness |
---|
177 | !! no ice shelves |
---|
178 | !! @return i_hp(:,:) and hp(:,:) |
---|
179 | |
---|
180 | subroutine prescribe_paleo_gl_shelf |
---|
181 | |
---|
182 | implicit none |
---|
183 | |
---|
184 | if (itracebug.eq.1) call tracebug(' Entree dans routine prescribe_paleo_gl_shelf') |
---|
185 | |
---|
186 | |
---|
187 | ! noeuds qui doivent être imposés |
---|
188 | |
---|
189 | where ((MK_flot0(:,:).eq. 1).or.(MK_gl0(:,:) .eq. 1)) ! floating or grounding line |
---|
190 | i_hp(:,:) = 1 ! thickness prescribed |
---|
191 | end where |
---|
192 | |
---|
193 | ! valeur imposee |
---|
194 | where (MK_flot0(:,:).eq. 1) ! paleo shelf epaisseur a 1 |
---|
195 | hp(:,:) = 1. |
---|
196 | end where |
---|
197 | |
---|
198 | |
---|
199 | where (MK_gl0(:,:) .eq. 1) |
---|
200 | hp(:,:) = (sealevel-Bsoc(:,:))*row/ro +1. |
---|
201 | end where |
---|
202 | |
---|
203 | |
---|
204 | end subroutine prescribe_paleo_gl_shelf |
---|
205 | |
---|
206 | !______________________________________________________________________________________ |
---|
207 | !> function prescribe_grid_boundary |
---|
208 | !! prescribe the grid boundary thickness to 0 |
---|
209 | !! @return i_hp(:,:) and hp(:,:) |
---|
210 | |
---|
211 | subroutine prescribe_grid_boundary_0 |
---|
212 | |
---|
213 | implicit none |
---|
214 | if (itracebug.eq.1) call tracebug(' Entree dans routine prescribe_grid_boundary_0') |
---|
215 | |
---|
216 | ! prescribe thickness on the grid boundaries |
---|
217 | i_hp(:,:) = 0 |
---|
218 | |
---|
219 | i_hp(1,:) = 1 |
---|
220 | i_hp(nx,:)= 1 |
---|
221 | hp(1,:) = 0 |
---|
222 | hp(nx,:) = 0 |
---|
223 | |
---|
224 | i_hp(:,1) = 1 |
---|
225 | i_hp(:,ny)= 1 |
---|
226 | hp(:,1) = 0 |
---|
227 | hp(:,ny) = 0 |
---|
228 | |
---|
229 | ! valeurs de reference |
---|
230 | i_hp0(1,:) = 1 |
---|
231 | i_hp0(nx,:) = 1 |
---|
232 | hp0(1,:) = 0 |
---|
233 | hp0(nx,:) = 0 |
---|
234 | |
---|
235 | i_hp0(:,1) = 1 |
---|
236 | i_hp0(:,ny) = 1 |
---|
237 | hp0(:,1) = 0 |
---|
238 | hp0(:,ny) = 0 |
---|
239 | |
---|
240 | |
---|
241 | end subroutine prescribe_grid_boundary_0 |
---|
242 | !______________________________________________________________________________________ |
---|
243 | |
---|
244 | !______________________________________________________________________________________ |
---|
245 | !> break ice shelves |
---|
246 | !! |
---|
247 | !! |
---|
248 | subroutine break_all_ice_shelves |
---|
249 | |
---|
250 | implicit none |
---|
251 | if (itracebug.eq.1) call tracebug(' Entree dans routine break_all_ice_shelves ') |
---|
252 | |
---|
253 | |
---|
254 | debug_3D(:,:,89)=0. |
---|
255 | |
---|
256 | where (flot(:,:)) ! floating from the beginning and eventually from grounding line retreat |
---|
257 | i_hp(:,:) = 1 ! thickness prescribed to 1 m |
---|
258 | hp (:,:) = 1. |
---|
259 | hp0 (:,:) = 1. |
---|
260 | H(:,:) =1. |
---|
261 | debug_3D(:,:,89)=1. |
---|
262 | |
---|
263 | end where |
---|
264 | where (MK_flot0(:,:).eq. 1) |
---|
265 | i_hp(:,:) = 1 ! thickness prescribed to 1 m |
---|
266 | hp (:,:) = 1. |
---|
267 | hp0 (:,:) = 1. |
---|
268 | H(:,:) =1. |
---|
269 | end where |
---|
270 | |
---|
271 | |
---|
272 | end subroutine break_all_ice_shelves |
---|
273 | |
---|
274 | ! |
---|
275 | !______________________________________________________________________________________ |
---|
276 | !> break ice shelves |
---|
277 | !! |
---|
278 | !! |
---|
279 | subroutine melt_ice_shelves |
---|
280 | |
---|
281 | implicit none |
---|
282 | integer :: nbvoisins,nbdeglac,iv,jv |
---|
283 | real :: hmoy,hmin,hmax |
---|
284 | |
---|
285 | if (itracebug.eq.1) call tracebug(' Entree dans routine melt_ice_shelves') |
---|
286 | |
---|
287 | do j = 2,ny-1 |
---|
288 | do i=2, nx-1 |
---|
289 | if (flot(i,j)) then ! floating |
---|
290 | i_hp(i,j) = 1 ! thickness prescribed |
---|
291 | hmoy=0. |
---|
292 | hmin=0. |
---|
293 | hmax=0. |
---|
294 | nbvoisins=0 |
---|
295 | nbdeglac=0 |
---|
296 | do jv=-1,1 |
---|
297 | do iv=1,1 |
---|
298 | if (flot(i+iv,j+jv)) then |
---|
299 | nbvoisins=nbvoisins+1 |
---|
300 | hmoy=hmoy+h(i+iv,j+jv) |
---|
301 | hmin=min(hmin,h(i+iv,j+jv)) |
---|
302 | hmax=max(hmax,h(i+iv,j+jv)) |
---|
303 | if (h(i+jv,j+jv).lt.10.) nbdeglac=nbdeglac+1 |
---|
304 | end if |
---|
305 | end do |
---|
306 | end do |
---|
307 | hmoy=hmoy/(max(nbvoisins,1)) |
---|
308 | |
---|
309 | if(H(i,j).gt.hmoy) then ! melt only if smooth |
---|
310 | hp(i,j)=H(i,j)-1.*dt |
---|
311 | else if (hmax-hmin.lt.10.) then |
---|
312 | hp(i,j)=H(i,j)-1.*dt |
---|
313 | else if ((H(i,j).lt.50.).and.(Hmin.lt.30.)) then |
---|
314 | hp(i,j)=1. |
---|
315 | else if ((H(i,j).lt.50.).and.(nbdeglac.ge.4)) then |
---|
316 | hp(i,j)=1. |
---|
317 | endif |
---|
318 | hp(i,j) = max(hp(i,j),1.) |
---|
319 | end if |
---|
320 | end do |
---|
321 | end do |
---|
322 | |
---|
323 | end subroutine melt_ice_shelves |
---|
324 | |
---|
325 | |
---|
326 | !______________________________________________________________________________________ |
---|
327 | !> function prescribe_present_H_gl |
---|
328 | !! calculate the initial grounding line position |
---|
329 | !! @return i_hp(:,:) and hp(:,:) |
---|
330 | |
---|
331 | subroutine prescribe_present_H_gl_copy |
---|
332 | |
---|
333 | implicit none |
---|
334 | if (itracebug.eq.1) call tracebug(' Entree dans routine prescribe_present_H_gl-copy') |
---|
335 | |
---|
336 | do j=1,ny |
---|
337 | do i=1,nx |
---|
338 | if ((MK_flot0(i,j).eq. 1).or.(MK_gl0(i,j) .eq. 1)) then |
---|
339 | i_hp(i,j) = 1 ! thickness prescribed to present value |
---|
340 | hp(i,j) = Hp0(i,j) |
---|
341 | end if |
---|
342 | end do |
---|
343 | end do |
---|
344 | |
---|
345 | !!$where ((MK_flot0(:,:).eq. 1).or.(MK_gl0(:,:) .eq. 1)) ! floating or grounding line |
---|
346 | !!$ i_hp(:,:) = 1 ! thickness prescribed to present value |
---|
347 | !!$! hp(:,:) = 1 !Hp0(:,:) |
---|
348 | !!$end where |
---|
349 | |
---|
350 | |
---|
351 | end subroutine prescribe_present_H_gl_copy |
---|
352 | !______________________________________________________________________________________ |
---|
353 | |
---|
354 | |
---|
355 | end module prescribe_H |
---|