1 | !> \file flottab2-0.6.F90 |
---|
2 | !! Module pour determiner les endroits ou la glace |
---|
3 | !! flotte , les iles, et la position du front de glace |
---|
4 | !< |
---|
5 | |
---|
6 | !> \namespace flottab_mod |
---|
7 | !! Determine les endroits ou la glace |
---|
8 | !! flotte , les iles, et la position du front de glace |
---|
9 | !! \author Vincent & Cat |
---|
10 | !! \date 10 juillet 2005 |
---|
11 | !! @note Used module |
---|
12 | !! @note - use module3D_phy |
---|
13 | !! @note - use module_choix |
---|
14 | !< |
---|
15 | module flottab_mod |
---|
16 | |
---|
17 | USE module3D_phy |
---|
18 | use module_choix |
---|
19 | |
---|
20 | |
---|
21 | |
---|
22 | IMPLICIT NONE |
---|
23 | |
---|
24 | real :: surnet !surnet hauteur de glace au dessus de la mer |
---|
25 | real :: archim ! test de flottaison |
---|
26 | |
---|
27 | integer:: itestf |
---|
28 | |
---|
29 | logical,dimension(nx,ny) :: gz1mx,gz1my |
---|
30 | logical,dimension(nx,ny) :: fl1mx,fl1my |
---|
31 | |
---|
32 | |
---|
33 | real,dimension(nx,ny) :: uxs1 ! uxbar a l'entree de flottab |
---|
34 | real,dimension(nx,ny) :: uys1 ! uybar a l'entree de flottab |
---|
35 | |
---|
36 | |
---|
37 | integer pmx,pmy !pm=plus-moins -1 ou 1 pour x et y |
---|
38 | |
---|
39 | |
---|
40 | ! Variables pour la determination des differents shelfs/stream |
---|
41 | ! (representés comme des taches ou l'on resoud l'eq elliptique) |
---|
42 | !________________________________________________________________ |
---|
43 | integer,parameter :: n_ta_max=2000! nombre de tache max |
---|
44 | ! integer,dimension(0:nx,0:ny) :: table_out ! pour les numeros des taches |
---|
45 | ! integer,dimension(0:nx,0:ny) :: tablebis ! pour les numeros des taches |
---|
46 | integer,dimension(0:nx,0:ny) :: table_out ! pour les numeros des taches |
---|
47 | integer,dimension(0:nx,0:ny) :: tablebis ! pour les numeros des taches |
---|
48 | integer,dimension(n_ta_max) :: compt !contient les equivalence entre les taches |
---|
49 | integer,dimension(0:n_ta_max) :: nb_pts_tache ! indique le nombre de points par tache |
---|
50 | logical,dimension(0:n_ta_max) :: iceberg ! T si iceberg, F si calotte posée |
---|
51 | |
---|
52 | ! Variables pour determiner le point le plus haut (surf) |
---|
53 | ! d'une ile completement stream |
---|
54 | !_________________________________________________________ |
---|
55 | |
---|
56 | ! icetrim : T si ice stream, F si calotte posée(vertical shear) |
---|
57 | logical,dimension(0: n_ta_max) :: icetrim |
---|
58 | |
---|
59 | integer :: ii,jj |
---|
60 | integer :: smax_i |
---|
61 | integer :: smax_j |
---|
62 | integer :: nb_pt ! numero de la tache (var travail) |
---|
63 | real :: smax_ |
---|
64 | |
---|
65 | contains |
---|
66 | ! ----------------------------------------------------------------------------------- |
---|
67 | !> SUBROUTINE: flottab() |
---|
68 | !! Cette routine determine les endroits ou la glace |
---|
69 | !! flotte , les iles, et la position du front de glace |
---|
70 | !! @note Il y a 4 sortes de zone |
---|
71 | !! @note - Pose |
---|
72 | !! @note - Grounding zone et streams gzmx et gzmy |
---|
73 | !! @note - Iles ilemx, ilemy |
---|
74 | !! @note - flottant flot sur le noeud majeur, flotmx sur le noeud mineur |
---|
75 | !> |
---|
76 | subroutine flottab() |
---|
77 | ! |
---|
78 | ! Vince 5 Jan 95 |
---|
79 | ! Modifie 20 Jan 95 |
---|
80 | ! Modifie le 30 Novembre 98 |
---|
81 | ! Passage f90 + determination des fronts Vincent dec 2003 |
---|
82 | ! nettoyage et nouvelle détermination de gzmx et gzmy Cat le 10 juillet 2005 |
---|
83 | ! Re-nettoyage par Cat en aout 2006. |
---|
84 | ! Le calcul de gzmx et gzmy pour les points intérieurs passe dans la subroutine dragging |
---|
85 | ! pour les points cotiers, toujours fait dans flottab |
---|
86 | ! |
---|
87 | ! ----------------- |
---|
88 | ! |
---|
89 | ! Cette routine determine les endroits ou la glace |
---|
90 | ! flotte , les iles, et la position du front de glace |
---|
91 | ! |
---|
92 | ! Il y a 4 sortes de zone |
---|
93 | ! Pose |
---|
94 | ! Grounding zone et streams gzmx et gzmy |
---|
95 | ! Iles ilemx, ilemy |
---|
96 | ! flottant flot sur le noeud majeur, flotmx sur le noeud mineur |
---|
97 | |
---|
98 | ! passage dans flottab tous les pas de temps dt ! |
---|
99 | ! |
---|
100 | ! _________________________________________________________ |
---|
101 | |
---|
102 | |
---|
103 | |
---|
104 | if (itracebug.eq.1) call tracebug(' Entree dans routine flottab-0.6') |
---|
105 | |
---|
106 | |
---|
107 | |
---|
108 | CALL write_trace('flottab') |
---|
109 | SHELFY = .FALSE. |
---|
110 | |
---|
111 | |
---|
112 | ! 1-INITIALISATION |
---|
113 | ! ---------------- |
---|
114 | ! initialisation des variables pour detecter les points qui se mettent |
---|
115 | ! a flotter entre 2 dtt |
---|
116 | |
---|
117 | appel_new_flot=.false. |
---|
118 | do j=1,ny |
---|
119 | do i=1,nx |
---|
120 | new_flot_point(i,j)=.false. |
---|
121 | new_flotmx(i,j)=.false. |
---|
122 | new_flotmy(i,j)=.false. |
---|
123 | enddo |
---|
124 | enddo |
---|
125 | |
---|
126 | ! ICE(:,:)=(H(:,:).gt.1) ! ice=.true. si epaisseur > 1m |
---|
127 | |
---|
128 | ICE(:,:)=0 |
---|
129 | front(:,:)=0 |
---|
130 | frontfacex(:,:)=0 |
---|
131 | frontfacey(:,:)=0 |
---|
132 | isolx(:,:)=.FALSE. |
---|
133 | isoly(:,:)=.FALSE. |
---|
134 | cotemx(:,:)=.false. |
---|
135 | cotemy(:,:)=.false. |
---|
136 | boost=.false. |
---|
137 | |
---|
138 | ! fin de l'initialisation |
---|
139 | !_____________________________________________________________________ |
---|
140 | |
---|
141 | ! 2-TESTE LES NOUVEAUX POINTS FLOTTANTS |
---|
142 | ! ------------------------------------- |
---|
143 | |
---|
144 | |
---|
145 | do j=1,ny |
---|
146 | do i=1,nx |
---|
147 | |
---|
148 | uxs1(i,j)=uxbar(i,j) |
---|
149 | uys1(i,j)=uybar(i,j) |
---|
150 | |
---|
151 | archim = Bsoc(i,j)+H(i,j)*ro/row -sealevel |
---|
152 | |
---|
153 | |
---|
154 | arch: if ((ARCHIM.LT.0.).and.(H(I,J).gt.1.E-3)) then ! le point flotte |
---|
155 | |
---|
156 | |
---|
157 | ex_pose: if ((.not.FLOT(I,J)).and.(isynchro.eq.1)) then ! il ne flottait pas avant |
---|
158 | FLOT(I,J)=.true. |
---|
159 | BOOST=.false. |
---|
160 | |
---|
161 | ! Attention : avec le bloc dessous il faut faire un calcul de flottaison a la lecture |
---|
162 | |
---|
163 | if (igrdline.eq.1) then ! en cas de grounding line prescrite |
---|
164 | flot(i,j)=.false. |
---|
165 | H(i,j)=(10.+sealevel-Bsoc(i,j))*row/ro ! pour avoir archim=10 |
---|
166 | new_flot_point(i,j)=.false. |
---|
167 | endif |
---|
168 | |
---|
169 | |
---|
170 | |
---|
171 | |
---|
172 | else ! isynchro=0 ou il flottait déja |
---|
173 | |
---|
174 | if (.not.FLOT(I,J)) then ! il ne flottait pas (isynchro=0) |
---|
175 | new_flot_point(i,j)=.true. ! signale un point qui ne flottait pas |
---|
176 | ! au pas de temps precedent |
---|
177 | flot(i,j)=.true. |
---|
178 | |
---|
179 | if (igrdline.eq.1) then ! en cas de grounding line prescrite |
---|
180 | flot(i,j)=.false. |
---|
181 | H(i,j)=(10.+sealevel-Bsoc(i,j))*row/ro ! pour avoir archim=10 |
---|
182 | new_flot_point(i,j)=.false. |
---|
183 | endif |
---|
184 | |
---|
185 | endif |
---|
186 | endif ex_pose |
---|
187 | |
---|
188 | |
---|
189 | else ! le point ne flotte pas |
---|
190 | |
---|
191 | |
---|
192 | if(FLOT(I,J)) then ! mais il flottait avant |
---|
193 | FLOT(I,J)=.false. |
---|
194 | BOOST=.false. |
---|
195 | endif |
---|
196 | endif arch |
---|
197 | |
---|
198 | ! Si la glace flotte -> PRUDENCE !!! |
---|
199 | ! S et B sont alors determines avec la condition de flottabilite |
---|
200 | |
---|
201 | if (flot(i,j)) then |
---|
202 | shelfy = .true. |
---|
203 | |
---|
204 | surnet=H(i,j)*(1.-ro/row) |
---|
205 | S(i,j)=surnet+sealevel |
---|
206 | B(i,j)=S(i,j)-H(i,j) |
---|
207 | end if |
---|
208 | |
---|
209 | end do |
---|
210 | end do |
---|
211 | |
---|
212 | !----------------------------------------------------------------------- |
---|
213 | |
---|
214 | domain: do j=2,ny |
---|
215 | do i=2,nx |
---|
216 | |
---|
217 | ! 3_x A- NOUVELLE DEFINITION DE FLOTMX, LES POINTS |
---|
218 | ! AYANT UN DES VOISINS FLOTTANTS SONT FLOTMX ET GZMX |
---|
219 | ! ------------------------------------------- |
---|
220 | |
---|
221 | ! fl1 est l'ancienne valeur de flotmx |
---|
222 | gz1mx(i,j)=gzmx(i,j) ! gz1 est l'ancienne valeur de gzmx |
---|
223 | fl1mx(i,j)=flotmx(i,j) ! fl1 est l'ancienne valeur de flotmx |
---|
224 | |
---|
225 | flotmx(i,j)=flot(i,j).and.flot(i-1,j) |
---|
226 | |
---|
227 | ! test pour detecter les nouveaux flotmx entre 2 dtt : |
---|
228 | |
---|
229 | if (flotmx(i,j).and.(new_flot_point(i,j).or. & |
---|
230 | new_flot_point(i-1,j))) then |
---|
231 | appel_new_flot=.true. |
---|
232 | new_flotmx(i,j)=.true. |
---|
233 | endif |
---|
234 | |
---|
235 | |
---|
236 | ! premiere determination de gzmx |
---|
237 | !__________________________________________________________________________ |
---|
238 | |
---|
239 | ! gzmx si un des deux voisins est flottant et l'autre posé |
---|
240 | ! i-1 i |
---|
241 | gzmx(i,j)=((flot(i,j).and..not.flot(i-1,j)) & ! F P |
---|
242 | .or.(.not.flot(i,j).and.flot(i-1,j))) ! P F |
---|
243 | |
---|
244 | ! A condition d'etre assez proche de la flottaison |
---|
245 | ! sur le demi noeud condition archim > 100 m |
---|
246 | |
---|
247 | archim=(Bsoc(i,j)+Bsoc(i-1,j))*0.5-sealevel+ro/row*Hmx(i,j) |
---|
248 | gzmx(i,j)=gzmx(i,j).and.(archim.le.100.) |
---|
249 | cotemx(i,j)=gzmx(i,j) |
---|
250 | |
---|
251 | |
---|
252 | ! 3_y B- NOUVELLE DEFINITION DE FLOTMY |
---|
253 | ! -------------------------------- |
---|
254 | |
---|
255 | |
---|
256 | gz1my(i,j)=gzmy(i,j) ! gz1 est l'ancienne valeur de gzmy |
---|
257 | fl1my(i,j)=flotmy(i,j) ! fl1 est l'ancienne valeur de flotmy |
---|
258 | |
---|
259 | flotmy(i,j)=flot(i,j).and.flot(i,j-1) |
---|
260 | |
---|
261 | ! test pour detecter les nouveaux flotmy entre 2 dtt : |
---|
262 | |
---|
263 | if (flotmy(i,j).and.(new_flot_point(i,j).or. & |
---|
264 | new_flot_point(i,j-1))) then |
---|
265 | appel_new_flot=.true. |
---|
266 | new_flotmy(i,j)=.true. |
---|
267 | endif |
---|
268 | |
---|
269 | ! premiere determination de gzmy |
---|
270 | !__________________________________________________________________________ |
---|
271 | |
---|
272 | ! gzmy si un des deux voisins est flottant et l'autre posé |
---|
273 | |
---|
274 | gzmy(i,j)=((flot(i,j).and..not.flot(i,j-1)) & |
---|
275 | .or.(.not.flot(i,j).and.flot(i,j-1))) |
---|
276 | |
---|
277 | ! A condition d'etre assez proche de la flottaison |
---|
278 | ! sur le demi noeud condition archim > 100 m |
---|
279 | |
---|
280 | archim=(Bsoc(i,j)+Bsoc(i,j-1))*0.5-sealevel+ro/row*Hmy(i,j) |
---|
281 | gzmy(i,j)=gzmy(i,j).and.(archim.le.100.) |
---|
282 | cotemy(i,j)=gzmy(i,j) |
---|
283 | |
---|
284 | end do |
---|
285 | end do domain |
---|
286 | |
---|
287 | |
---|
288 | |
---|
289 | |
---|
290 | |
---|
291 | |
---|
292 | !------------------------------------------------------------------------------------- |
---|
293 | ! attention : pour expériences Heino |
---|
294 | ! gzmy(i,j)=gzmy_heino(i,j) |
---|
295 | ! shelfy=.true. |
---|
296 | ! shelfy=.false. |
---|
297 | !____________________________________________________________________________________ |
---|
298 | |
---|
299 | |
---|
300 | |
---|
301 | |
---|
302 | ! 4- Condition sur les bords |
---|
303 | |
---|
304 | do i=2,nx |
---|
305 | flotmx(i,1) = (flot(i,1).or.flot(i-1,1)) |
---|
306 | flotmy(i,1) = .false. |
---|
307 | end do |
---|
308 | |
---|
309 | do j=2,ny |
---|
310 | flotmy(1,j) = (flot(1,j).or.flot(1,j-1)) |
---|
311 | flotmx(1,j) = .false. |
---|
312 | end do |
---|
313 | |
---|
314 | flotmx(1,1) = .false. |
---|
315 | flotmy(1,1) = .false. |
---|
316 | |
---|
317 | |
---|
318 | ! 4- determination des iles |
---|
319 | ! ------------------------- |
---|
320 | |
---|
321 | ilemx(:,:)=.false. |
---|
322 | ilemy(:,:)=.false. |
---|
323 | |
---|
324 | ! selon x |
---|
325 | ilesx: do j=2,ny-1 |
---|
326 | do i=3,nx-2 |
---|
327 | ! F G F |
---|
328 | ! x |
---|
329 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
330 | if ((flot(i-1,j).and..not.flot(i,j).and.flot(i+1,j)).and. & |
---|
331 | (sdx(i,j).LT.1.E-02)) then |
---|
332 | ilemx(i,j)=.true. |
---|
333 | ilemx(i+1,j)=.true. |
---|
334 | |
---|
335 | ! F G G F |
---|
336 | ! x |
---|
337 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
338 | else if ((flot(i-1,j).and..not.flot(i,j) & |
---|
339 | .and..not.flot(i+1,j)).and.flot(i+2,j).and. & |
---|
340 | (sdx(i,j).LT.1.E-02.and.sdx(i+1,j).LT.1.E-02)) then |
---|
341 | ilemx(i,j)=.true. |
---|
342 | ilemx(i+1,j)=.true. |
---|
343 | ilemx(i+2,j)=.true. |
---|
344 | |
---|
345 | ! F G G F |
---|
346 | ! x |
---|
347 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
348 | else if ((flot(i-2,j).and..not.flot(i-1,j) & |
---|
349 | .and..not.flot(i,j)).and.flot(i+1,j).and. & |
---|
350 | (sdx(i,j).LT.1.E-02.and.sdx(i-1,j).LT.1.E-02)) then |
---|
351 | ilemx(i-1,j)=.true. |
---|
352 | ilemx(i,j)=.true. |
---|
353 | ilemx(i+1,j)=.true. |
---|
354 | |
---|
355 | ! F G G G F |
---|
356 | ! x |
---|
357 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
358 | else if ((i.lt.nx-2) & |
---|
359 | .and.(flot(i-2,j).and..not.flot(i-1,j) & |
---|
360 | .and..not.flot(i,j)).and..not.flot(i+1,j) & |
---|
361 | .and.flot(i+2,j).and. & |
---|
362 | (sdx(i,j).LT.1.E-02.and.sdx(i-1,j).LT.1.E-02 & |
---|
363 | .and.sdx(i+1,j).LT.1.E-02)) then |
---|
364 | ilemx(i-1,j)=.true. |
---|
365 | ilemx(i,j)=.true. |
---|
366 | ilemx(i+1,j)=.true. |
---|
367 | ilemx(i+2,j)=.true. |
---|
368 | |
---|
369 | endif |
---|
370 | |
---|
371 | end do |
---|
372 | end do ilesx |
---|
373 | |
---|
374 | ! selon y |
---|
375 | ilesy: do j=3,ny-2 |
---|
376 | do i=2,nx-1 |
---|
377 | ! F G F |
---|
378 | ! x |
---|
379 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
380 | if ((flot(i,j-1).and..not.flot(i,j).and.flot(i,j+1)).and. & |
---|
381 | (sdy(i,j).LT.1.E-02)) then |
---|
382 | ilemy(i,j)=.true. |
---|
383 | ilemy(i,j+1)=.true. |
---|
384 | |
---|
385 | ! F G G F |
---|
386 | ! x |
---|
387 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
388 | else if ((flot(i,j-1).and..not.flot(i,j) & |
---|
389 | .and..not.flot(i,j+1)).and.flot(i,j+2).and. & |
---|
390 | (sdy(i,j).LT.1.E-02.and.sdy(i,j+1).LT.1.E-02)) then |
---|
391 | ilemy(i,j)=.true. |
---|
392 | ilemy(i,j+1)=.true. |
---|
393 | ilemy(i,j+2)=.true. |
---|
394 | |
---|
395 | ! F G G F |
---|
396 | ! x |
---|
397 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
398 | else if ((flot(i,j-2).and..not.flot(i,j-1) & |
---|
399 | .and..not.flot(i,j)).and.flot(i,j+1).and. & |
---|
400 | (sdy(i,j).LT.1.E-02.and.sdy(i,j-1).LT.1.E-02)) then |
---|
401 | ilemy(i,j-1)=.true. |
---|
402 | ilemy(i,j)=.true. |
---|
403 | ilemy(i,j+1)=.true. |
---|
404 | |
---|
405 | ! F G G G F |
---|
406 | ! x |
---|
407 | ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m) |
---|
408 | else if ((j.lt.ny-2) & |
---|
409 | .and.(flot(i,j-2).and..not.flot(i,j-1) & |
---|
410 | .and..not.flot(i,j)).and..not.flot(i,j+1) & |
---|
411 | .and.flot(i,j+2).and. & |
---|
412 | (sdy(i,j).LT.1.E-02.and.sdy(i,j-1).LT.1.E-02 & |
---|
413 | .and.sdy(i,j+1).LT.1.E-02)) then |
---|
414 | ilemy(i,j-1)=.true. |
---|
415 | ilemy(i,j)=.true. |
---|
416 | ilemy(i,j+1)=.true. |
---|
417 | ilemy(i,j+2)=.true. |
---|
418 | endif |
---|
419 | end do |
---|
420 | end do ilesy |
---|
421 | ! fin des iles |
---|
422 | |
---|
423 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
424 | !!$if (itestf.gt.0) then |
---|
425 | !!$ write(6,*) 'dans flottab avant dragging asymetrie sur H pour time=',time |
---|
426 | !!$ stop |
---|
427 | !!$else |
---|
428 | !!$ write(6,*) 'dans flottab avant dragging pas d asymetrie sur H pour time=',time |
---|
429 | !!$ |
---|
430 | !!$end if |
---|
431 | |
---|
432 | ! 5- calcule les noeuds qui sont streams a l'interieur et donne le betamx et betamy |
---|
433 | !---------------------------------------------------------------------------------- |
---|
434 | Call dragging |
---|
435 | |
---|
436 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
437 | !!$if (itestf.gt.0) then |
---|
438 | !!$ write(6,*) 'dans flottab apres dragging asymetrie sur H pour time=',time |
---|
439 | !!$ stop |
---|
440 | !!$else |
---|
441 | !!$ write(6,*) 'dans flottab aapres dragging pas d asymetrie sur H pour time=',time |
---|
442 | !!$ |
---|
443 | !!$end if |
---|
444 | |
---|
445 | ! 6- calcule les vitesses des points qui sont devenus gzm |
---|
446 | |
---|
447 | do j=1,ny |
---|
448 | do i=2,nx-1 |
---|
449 | ! si le point etait posé (non gz) et devient gzmx |
---|
450 | ! definir la direction de la vitesse (moyenne des points) |
---|
451 | |
---|
452 | if ((.not.gz1mx(i,j)).and.(.not.fl1mx(i,j)).and.gzmx(i,j).and. & |
---|
453 | (i.gt.2).and.(i.lt.nx)) then |
---|
454 | uxs1(i,j)=(uxbar(i+1,j)+uxbar(i-1,j))/2. |
---|
455 | endif |
---|
456 | |
---|
457 | end do |
---|
458 | end do |
---|
459 | |
---|
460 | do j=2,ny-1 |
---|
461 | do i=1,nx |
---|
462 | ! si le point etait posé (non gz) et devient gzmy |
---|
463 | ! definir la direction de la vitesse (moyenne des points) |
---|
464 | if ((.not.gz1my(i,j)).and.(.not.fl1my(i,j)).and.gzmy(i,j).and. & |
---|
465 | (j.gt.2).and.(j.lt.ny)) then |
---|
466 | uys1(i,j)=(uybar(i,j+1)+uybar(i,j-1))/2. |
---|
467 | endif |
---|
468 | |
---|
469 | end do |
---|
470 | end do |
---|
471 | |
---|
472 | |
---|
473 | |
---|
474 | ! 7-On determine finalement la position des noeuds stream ou shelf |
---|
475 | ! ------------------------------------------------------------- |
---|
476 | |
---|
477 | if (nt.ge.2) then ! pour ne pas faire ce calcul lors du premier passage |
---|
478 | uxbar(:,:)=uxs1(:,:) |
---|
479 | uybar(:,:)=uys1(:,:) |
---|
480 | endif |
---|
481 | |
---|
482 | |
---|
483 | flgzmx(:,:)=(marine.and.(flotmx(:,:).or.gzmx(:,:).or.ilemx(:,:))) & |
---|
484 | .or.(.not.marine.and.flotmx(:,:)) |
---|
485 | flgzmy(:,:)=(marine.and.(flotmy(:,:).or.gzmy(:,:).or.ilemy(:,:))) & |
---|
486 | .or.(.not.marine.and.flotmy(:,:)) |
---|
487 | |
---|
488 | |
---|
489 | |
---|
490 | ! 8- Pour la fusion basale sous les ice shelves- region proche de la grounding line |
---|
491 | !--------------------------------------------------------------------------------- |
---|
492 | ! fbm est vrai si le point est flottant mais un des voisins est pose |
---|
493 | !_________________________________________________________________________ |
---|
494 | do j=2,ny-1 |
---|
495 | do i=2,nx-1 |
---|
496 | |
---|
497 | ! if (i.gt.2.AND.i.lt.nx) then |
---|
498 | fbm(i,j)=flot(i,j).and. & |
---|
499 | ((.not.flot(i+1,j)).or.(.not.flot(i,j+1)) & |
---|
500 | .or.(.not.flot(i-1,j)).or.(.not.flot(i,j-1))) |
---|
501 | ! endif |
---|
502 | end do |
---|
503 | end do |
---|
504 | |
---|
505 | |
---|
506 | |
---|
507 | |
---|
508 | |
---|
509 | |
---|
510 | ! 9-On determine maintenant la position du front de glace |
---|
511 | ! ------------------------------------------------------- |
---|
512 | ! C'est ici que l'on determine la position et le type de front |
---|
513 | ! print*,'on est dans flottab pour definir les fronts' |
---|
514 | |
---|
515 | |
---|
516 | |
---|
517 | !!$do i=3,nx-2 |
---|
518 | !!$ do j=3,ny-2 |
---|
519 | !!$ if (h(i,j).gt.1.1) ice(i,j)=1 |
---|
520 | !!$ end do |
---|
521 | !!$end do |
---|
522 | |
---|
523 | where (flot(:,:)) |
---|
524 | where (H(:,:).gt.(1.1)) |
---|
525 | ice(:,:)=1 |
---|
526 | elsewhere |
---|
527 | ice(:,:)=0 |
---|
528 | end where |
---|
529 | elsewhere |
---|
530 | where (H(:,:).gt.(.1)) |
---|
531 | ice(:,:)=1 |
---|
532 | elsewhere |
---|
533 | ice(:,:)=0 |
---|
534 | end where |
---|
535 | end where |
---|
536 | |
---|
537 | |
---|
538 | synchro : if (isynchro.eq.1) then |
---|
539 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
540 | !!$if (itestf.gt.0) then |
---|
541 | !!$ write(6,*) 'dans flottab avant DETERMIN_TACHE asymetrie sur H pour time=',time |
---|
542 | !!$ stop |
---|
543 | !!$else |
---|
544 | !!$ write(6,*) 'dans flottab apres DETERMIN_TACHE pas d asymetrie sur H pour time=',time |
---|
545 | !!$ |
---|
546 | !!$end if |
---|
547 | |
---|
548 | !----------------------------------------------! |
---|
549 | !On determine les differents ice strean/shelf ! |
---|
550 | call DETERMIN_TACHE ! |
---|
551 | !----------------------------------------------! |
---|
552 | |
---|
553 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
554 | !!$if (itestf.gt.0) then |
---|
555 | !!$ write(6,*) 'dans flottab apres DETERMIN_TACHE asymetrie sur H pour time=',time |
---|
556 | !!$ stop |
---|
557 | !!$else |
---|
558 | !!$ write(6,*) 'dans flottab apres DETERMIN_TACHE pas d asymetrie sur H pour time=',time |
---|
559 | !!$ |
---|
560 | !!$end if |
---|
561 | |
---|
562 | ! / test aurel |
---|
563 | do i= 1, nx |
---|
564 | do j=1,ny |
---|
565 | debug_3d(i,j,43)=table_out(i-1,j-1) |
---|
566 | end do |
---|
567 | end do |
---|
568 | ! fin test aurel / |
---|
569 | |
---|
570 | |
---|
571 | !ice(:,:)=0 Attention, voir si ca marche toujours pour l'Antarctique et heminord ! |
---|
572 | |
---|
573 | !On compte comme englacé uniquement les calottes dont une partie est posée |
---|
574 | |
---|
575 | do i=3,nx-2 |
---|
576 | do j=3,ny-2 |
---|
577 | test1: if (.not.iceberg(table_out(i,j))) then ! on est pas sur un iceberg |
---|
578 | if (nb_pts_tache(table_out(i,j)).ge.1) then |
---|
579 | ice(i,j)=1 |
---|
580 | if (nb_pts_tache(table_out(i,j)).le.5) then |
---|
581 | flgzmx(i,j)=.false. |
---|
582 | flgzmy(i,j)=.false. |
---|
583 | endif |
---|
584 | |
---|
585 | |
---|
586 | ! ici on est sur une tache non iceberg >= 5 points |
---|
587 | ! on teste si la tache n'est pas completement ice stream |
---|
588 | |
---|
589 | test2: if (icetrim(table_out(i,j))) then ! on a une ile d'ice stream !!!! test2 |
---|
590 | smax_i=0 ; smax_j=0 ; smax_=sealevel |
---|
591 | do ii=3,nx-2 |
---|
592 | do jj=3,ny-2 |
---|
593 | if (table_out(ii,jj).eq.table_out(i,j)) then |
---|
594 | if (s(ii,jj).gt.smax_) then |
---|
595 | smax_ =s(ii,jj) |
---|
596 | smax_i=ii |
---|
597 | smax_j=jj |
---|
598 | endif |
---|
599 | endif |
---|
600 | end do |
---|
601 | end do |
---|
602 | ! test aurel |
---|
603 | if((smax_i.eq.0).or.(smax_j.eq.0)) then |
---|
604 | print*, 'smax_i ou smax_j est nul', smax_i, smax_j,dt |
---|
605 | print*, 'je touche a rien du coup' |
---|
606 | ! aurel -> arbitrairement je prends le premier indice |
---|
607 | gzmx(3,3)=.false. ; gzmx(4,3)=.false. |
---|
608 | gzmy(3,3)=.false. ; gzmx(3,4)=.false. |
---|
609 | flgzmx(3,3)=.false. ; flgzmx(4,3)=.false. |
---|
610 | flgzmy(3,3)=.false. ; flgzmx(3,4)=.false. |
---|
611 | |
---|
612 | else |
---|
613 | gzmx(smax_i,smax_j)=.false. ; gzmx(smax_i+1,smax_j)=.false. |
---|
614 | gzmy(smax_i,smax_j)=.false. ; gzmx(smax_i,smax_j+1)=.false. |
---|
615 | flgzmx(smax_i,smax_j)=.false. ; flgzmx(smax_i+1,smax_j)=.false. |
---|
616 | flgzmy(smax_i,smax_j)=.false. ; flgzmx(smax_i,smax_j+1)=.false. |
---|
617 | endif |
---|
618 | ! fin test aurel |
---|
619 | |
---|
620 | endif test2 |
---|
621 | |
---|
622 | |
---|
623 | else ! on est sur un iceberg ! test1 |
---|
624 | ice(i,j)=0 |
---|
625 | h(i,j)=1. |
---|
626 | end if |
---|
627 | endif test1 |
---|
628 | |
---|
629 | |
---|
630 | end do |
---|
631 | end do |
---|
632 | |
---|
633 | |
---|
634 | !---------------------------------------------- |
---|
635 | ! On caracterise le front des ice shelfs/streams |
---|
636 | |
---|
637 | call DETERMIN_FRONT |
---|
638 | |
---|
639 | !---------------------------------------------- |
---|
640 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
641 | !!$if (itestf.gt.0) then |
---|
642 | !!$ write(6,*) 'dans flottab apres DETERMIN_front asymetrie sur H pour time=',time |
---|
643 | !!$ stop |
---|
644 | !!$else |
---|
645 | !!$ write(6,*) 'dans flottab apres DETERMIN_front pas d asymetrie sur H pour time=',time |
---|
646 | !!$ |
---|
647 | !!$end if |
---|
648 | |
---|
649 | endif synchro |
---|
650 | |
---|
651 | ! correction momentanée pour symetrie Heino |
---|
652 | !where ((.not.flot(:,:)).and.(ice(:,:).eq.0)) H(:,:)=0. |
---|
653 | |
---|
654 | !fin de routine flottab2 |
---|
655 | !print*, 'front',front(50,30),ice(50,30),flotmx(i,j),uxbar(i,j) |
---|
656 | end subroutine flottab |
---|
657 | !-------------------------------------------------------------------- |
---|
658 | !> SUBROUTINE: determin_tache |
---|
659 | !! Routine pour la dtermination du numero de tache a effectuer |
---|
660 | !> |
---|
661 | subroutine determin_tache |
---|
662 | |
---|
663 | implicit none |
---|
664 | integer :: indice |
---|
665 | integer :: label ! no des taches rencontrées dans le mask |
---|
666 | integer :: label_max ! no temporaire maxi de tache rencontrées |
---|
667 | ! integer :: mask_nb = 4 |
---|
668 | integer :: mask_nb = 2 ! version ou on ne compte pas les diagonales |
---|
669 | integer :: vartemp ! variable temporaire pour reorganiser compt |
---|
670 | ! integer,dimension(mask_nb) :: mask |
---|
671 | integer,dimension(2) :: mask |
---|
672 | |
---|
673 | ! 1-initialisation |
---|
674 | !----------------- |
---|
675 | label_max=1 ! numero de la tache, la premiere tache est notée 1 |
---|
676 | label=1 |
---|
677 | do i=1,n_ta_max |
---|
678 | compt(i)=i |
---|
679 | enddo |
---|
680 | ! table_in = .false. |
---|
681 | table_out = 0 |
---|
682 | iceberg = .true. |
---|
683 | icetrim = .true. |
---|
684 | nb_pts_tache = 0 |
---|
685 | ! open(unit=100,file="tache.data",status='replace') |
---|
686 | |
---|
687 | ! 2-reperage des taches |
---|
688 | !---------------------- |
---|
689 | do i=2,nx |
---|
690 | do j=2,ny |
---|
691 | |
---|
692 | if (label_max.eq.0) print*,'labmax=0',i,j |
---|
693 | |
---|
694 | IF (ice(i,j).ge.1) THEN ! on est sur la glace-----------------------------! |
---|
695 | ! if ((ice(i-1,j-1).ge.1).or.(ice(i-1,j).ge.1).or.& !masque de 4 cases adjacentes |
---|
696 | ! (ice(i-1,j+1).ge.1).or.(ice(i,j-1).ge.1)) then |
---|
697 | ! mask(1) = table_out(i-1,j-1) |
---|
698 | ! mask(2) = table_out(i-1,j) |
---|
699 | ! mask(3) = table_out(i-1,j+1) |
---|
700 | ! mask(4) = table_out(i,j-1) |
---|
701 | |
---|
702 | if ((ice(i-1,j).ge.1).or.(ice(i,j-1).ge.1)) then !masque de 2 cases adjacentes |
---|
703 | ! un des voisins est deja en glace |
---|
704 | mask(1) = table_out(i-1,j) |
---|
705 | mask(2) = table_out(i,j-1) |
---|
706 | label = label_max |
---|
707 | |
---|
708 | !on determine la valeur de la tache minimun (>0) presente ds le masque |
---|
709 | do indice=1,mask_nb |
---|
710 | if (mask(indice).gt.0) label=min(label,mask(indice)) |
---|
711 | enddo |
---|
712 | |
---|
713 | !on fixe la valeur de la tache voisine minimun au point etudié(via label) |
---|
714 | table_out(i,j)=label |
---|
715 | !si ce noeud est posé, alors le tache n'est pas un iceberg et iceberg=.F. |
---|
716 | if (.not.FLOT(I,J)) then |
---|
717 | iceberg(label)=.false. |
---|
718 | endif |
---|
719 | !si ce noeud est posé, alors la tache n'est pas un ice stream et icestrim=.F. |
---|
720 | if (.not.gzmx(I,J).AND..not.gzmy(I,J)) then |
---|
721 | icetrim(label)=.false. |
---|
722 | endif |
---|
723 | |
---|
724 | ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt |
---|
725 | ! i.e. les plus grands numeros correspondent au plus petit |
---|
726 | ! on lui affecte le numero de la tache fondamentale avec un signe - |
---|
727 | ! pour indiquer le changement |
---|
728 | |
---|
729 | do indice=1,mask_nb |
---|
730 | if(mask(indice).gt.label) then |
---|
731 | compt(mask(indice))=-label |
---|
732 | endif |
---|
733 | enddo |
---|
734 | |
---|
735 | else !aucun des voisins est une tache |
---|
736 | table_out(i,j)= label_max |
---|
737 | compt(label_max)=label_max |
---|
738 | !si ce noeud est posé, alors la ache n'est pas un iceberg et iceberg=.F. |
---|
739 | if (.not.FLOT(I,J)) then |
---|
740 | iceberg(label_max)=.false. |
---|
741 | endif |
---|
742 | !si ce noeud est posé, alors le tache n'est pas un ice stream et icestrim=.F. |
---|
743 | if (.not.gzmx(I,J).AND..not.gzmy(I,J)) then |
---|
744 | icetrim(label)=.false. |
---|
745 | endif |
---|
746 | |
---|
747 | label_max = label_max+1 |
---|
748 | if (label_max.gt.n_ta_max) print*,'trop de taches=',label_max |
---|
749 | endif |
---|
750 | |
---|
751 | else !on est pas sur une tache---------------------------------------------- |
---|
752 | table_out(i,j)=0 ! Pas necessaire (reecrit 0 sur 0) |
---|
753 | endif !--------------------------------------------------------------------- |
---|
754 | |
---|
755 | |
---|
756 | enddo |
---|
757 | enddo |
---|
758 | |
---|
759 | ! On reorganise compt en ecrivant le numero de la tache fondamentale |
---|
760 | ! i.e. du plus petit numero present sur la tache (Sans utiliser de recursivité) |
---|
761 | ! On indique aussi le nb de point que contient chaque taches (nb_pts_tache) |
---|
762 | do indice=1,label_max |
---|
763 | vartemp = compt(indice) |
---|
764 | if (compt(indice).lt.0) then |
---|
765 | compt(indice)= compt(-vartemp) |
---|
766 | if (.not.iceberg(indice)) iceberg(-vartemp)=.false. |
---|
767 | if (.not.icetrim(indice)) icetrim(-vartemp)=.false. |
---|
768 | endif |
---|
769 | enddo |
---|
770 | |
---|
771 | !!$do i=1,nx |
---|
772 | !!$ do j=1,ny |
---|
773 | !!$ if (table_out(i,j).ne.0) then |
---|
774 | !!$ table_out(i,j)=compt(table_out(i,j)) |
---|
775 | !!$ nb_pts_tache(compt(table_out(i,j)))= nb_pts_tache(compt(table_out(i,j)))+1 |
---|
776 | !!$ endif |
---|
777 | !!$ enddo |
---|
778 | !!$enddo |
---|
779 | |
---|
780 | tablebis(:,:)=table_out(:,:) |
---|
781 | do i=1,nx |
---|
782 | do j=1,ny |
---|
783 | if (tablebis(i,j).ne.0) then ! tache de glace |
---|
784 | numtache=table_out(i,j) |
---|
785 | nb_pt=count(table_out(:,:).eq.numtache) ! compte tous les points de la tache |
---|
786 | where (table_out(:,:).eq.numtache) |
---|
787 | nb_pts_tache(:,:)=nb_pt ! tous les points de la tache sont tagges nb_pt |
---|
788 | tablebis(:,:)=0 ! la table de tache est remise a 0 |
---|
789 | end where |
---|
790 | write(6,*) i,j,nb_pt,table_out(i,j) |
---|
791 | ! table_out(i,j)=compt(table_out(i,j)) |
---|
792 | ! nb_pts_tache(compt(table_out(i,j)))= nb_pts_tache(compt(table_out(i,j)))+1 |
---|
793 | endif |
---|
794 | enddo |
---|
795 | enddo |
---|
796 | |
---|
797 | |
---|
798 | ! 3-visualisation |
---|
799 | !---------------- |
---|
800 | !print*,'_pts_tache0:',nb_pts_tache(0:label_max+15) |
---|
801 | !print*,'lab tache', label_max |
---|
802 | !print*,'compt',compt |
---|
803 | !om_table='tache' |
---|
804 | !call printtable_i(table_out,nom_table) |
---|
805 | !nom_table='bergs' |
---|
806 | !call printtable_l(iceberg(table_out),nom_table) |
---|
807 | !nom_table='pt_ta' |
---|
808 | !call printtable_i(nb_pts_tache(table_out),nom_table) |
---|
809 | |
---|
810 | |
---|
811 | end subroutine determin_tache |
---|
812 | !---------------------------------------------------------------------- |
---|
813 | !> SUBROUTINE: determin_front |
---|
814 | !!Routine pour la determination du front |
---|
815 | !> |
---|
816 | subroutine determin_front |
---|
817 | |
---|
818 | integer :: i_moins1,i_plus1,i_plus2 |
---|
819 | integer :: j_moins1,j_plus1,j_plus2 |
---|
820 | |
---|
821 | do i=3,nx-2 |
---|
822 | do j=3,ny-2 |
---|
823 | |
---|
824 | surice:if (ice(i,j).eq.0) then |
---|
825 | do pmx=-1,1,2 |
---|
826 | do pmy=-1,1,2 |
---|
827 | |
---|
828 | diagice : if (ice(i+pmx,j+pmy).eq.1) then |
---|
829 | |
---|
830 | if ((ice(i+pmx,j)+ice(i,j+pmy).eq.2)) then ! test (i) pour eviter les langues |
---|
831 | ! de glaces diagonales en coin(26dec04) |
---|
832 | if ((ice(i+2*pmx,j).eq.1.and.ice(i+2*pmx,j+pmy).eq.0).or.& |
---|
833 | (ice(i,j+2*pmy).eq.1.and.ice(i+pmx,j+2*pmy).eq.0)) then |
---|
834 | ice(i,j)=1 |
---|
835 | h(i,j)=max(1.,h(i,j)) |
---|
836 | endif |
---|
837 | |
---|
838 | ! test (i) pour eviter les langues de glaces diagonales : |
---|
839 | ! mouvement du cheval aux echecs |
---|
840 | |
---|
841 | if ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.1)) then |
---|
842 | if (ice(i+2*pmx,j+pmy).eq.1.and. & |
---|
843 | (ice(i+2*pmx,j+2*pmy)+ice(i,j+2*pmy)).ge.1) then |
---|
844 | ice(i,j)=1 |
---|
845 | h(i,j)=max(1.,h(i,j)) |
---|
846 | endif |
---|
847 | if (ice(i+pmx,j+2*pmy).eq.1.and. & |
---|
848 | (ice(i+2*pmx,j+2*pmy)+ice(i+2*pmx,j)).ge.1) then |
---|
849 | ice(i,j)=1 |
---|
850 | h(i,j)=max(1.,h(i,j)) |
---|
851 | endif |
---|
852 | |
---|
853 | ! test (ii) pour eviter les langues de glaces diagonales : |
---|
854 | ! le point glace ice(i+pmx,j+pmy) a : |
---|
855 | ! - ses 4 voisins frontaux en glace |
---|
856 | ! - mais 2 voisins vides diagonalement opposes |
---|
857 | |
---|
858 | elseif ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.2) & |
---|
859 | .and.ice(i+2*pmx,j+2*pmy).eq.0) then |
---|
860 | |
---|
861 | ! test (iii) pour faire les tests (i) et (ii) |
---|
862 | ice(i,j)=1 |
---|
863 | h(i,j)=max(1.,h(i,j)) |
---|
864 | ice(i+2*pmx,j+2*pmy)=1 |
---|
865 | h(i+2*pmx,j+2*pmy)=max(1.,h(i+2*pmx,j+2*pmy)) |
---|
866 | endif |
---|
867 | endif |
---|
868 | endif diagice |
---|
869 | enddo |
---|
870 | enddo |
---|
871 | endif surice |
---|
872 | end do |
---|
873 | end do |
---|
874 | |
---|
875 | |
---|
876 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
877 | !!$if (itestf.gt.0) then |
---|
878 | !!$ write(6,*) 'dans front avant remplissage baies asymetrie sur H pour time=',time |
---|
879 | !!$ stop |
---|
880 | !!$else |
---|
881 | !!$ write(6,*) 'dans front avant remplissage baies pas d asymetrie sur H pour time=',time |
---|
882 | !!$ |
---|
883 | !!$end if |
---|
884 | |
---|
885 | |
---|
886 | ! print*,'dans remplissage baies',time |
---|
887 | baies: do k=1,2 |
---|
888 | do j=1,ny |
---|
889 | do i=1,nx |
---|
890 | |
---|
891 | surice_xy: if (ice(i,j).eq.0) then |
---|
892 | i_moins1=max(i-1,1) |
---|
893 | j_moins1=max(j-1,1) |
---|
894 | i_plus1=min(i+1,nx) |
---|
895 | j_plus1=min(j+1,ny) |
---|
896 | i_plus2=min(i+2,nx) |
---|
897 | j_plus2=min(j+2,ny) |
---|
898 | |
---|
899 | ! test (iii) pour trouver les baies de largeur 1 ou 2 cases |
---|
900 | ! et combler les trous si ce sont des baies |
---|
901 | ! si ce ne sont pas des baies, ne pas combler et creer des langues de glaces artificielles |
---|
902 | ! baies horizontales |
---|
903 | |
---|
904 | if (ice(i_moins1,j).eq.1.and.(ice(i_plus1,j).eq.1.or.ice(i_plus2,j).eq.1)) then |
---|
905 | if (ice(i,j_moins1).eq.1.or.ice(i,j_plus1).eq.1) then ! ice(i,j)=1 |
---|
906 | ice(i,j)=1 |
---|
907 | H(i,j)=max(1.,H(i,j)) |
---|
908 | endif |
---|
909 | endif |
---|
910 | |
---|
911 | |
---|
912 | if (ice(i,j_moins1).eq.1.and.(ice(i,j_plus1).eq.1.or.ice(i,j_plus2).eq.1)) then |
---|
913 | if (ice(i_moins1,j).eq.1.or.ice(i_plus1,j).eq.1) then !ice(i,j)=1 |
---|
914 | ice(i,j)=1 |
---|
915 | H(i,j)=max(1.,H(i,j)) |
---|
916 | endif |
---|
917 | endif |
---|
918 | |
---|
919 | endif surice_xy |
---|
920 | end do |
---|
921 | end do |
---|
922 | end do baies |
---|
923 | |
---|
924 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
925 | !!$if (itestf.gt.0) then |
---|
926 | !!$ write(6,*) 'dans front apres remplissage baies asymetrie sur H pour time=',time |
---|
927 | !!$ stop |
---|
928 | !!$else |
---|
929 | !!$ write(6,*) 'dans front apres remplissage baies pas d asymetrie sur H pour time=',time |
---|
930 | !!$ |
---|
931 | !!$end if |
---|
932 | |
---|
933 | |
---|
934 | do i=2,nx-1 |
---|
935 | do j=2,ny-1 |
---|
936 | |
---|
937 | if (ice(i,j).eq.1) then ! test si ice=1 |
---|
938 | |
---|
939 | ! if ice, on determine front... |
---|
940 | ! ainsi, front=0 sur les zones = 0 |
---|
941 | |
---|
942 | front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1)) |
---|
943 | !front= le nb de faces en contact avec un voisin englacé |
---|
944 | endif |
---|
945 | end do |
---|
946 | end do |
---|
947 | |
---|
948 | ! traitement des bords. On considere que l'exterieur n'a pas de glace |
---|
949 | ! attention ce n'est vrai que sur la grande grille |
---|
950 | |
---|
951 | |
---|
952 | do j=2,ny-1 |
---|
953 | i=1 |
---|
954 | front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1)) |
---|
955 | i=nx |
---|
956 | front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1)) |
---|
957 | end do |
---|
958 | |
---|
959 | do i=2,nx-1 |
---|
960 | j=1 |
---|
961 | front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)) |
---|
962 | j=ny |
---|
963 | front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1)) |
---|
964 | end do |
---|
965 | |
---|
966 | ! traitement des coins |
---|
967 | |
---|
968 | front(1,1)=ice(2,1)+ice(2,1) |
---|
969 | front(1,ny)=ice(2,ny)+ice(1,ny-1) |
---|
970 | front(nx,1)=ice(nx,2)+ice(nx-1,1) |
---|
971 | front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny) |
---|
972 | |
---|
973 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf) |
---|
974 | !!$if (itestf.gt.0) then |
---|
975 | !!$ write(6,*) 'dans front apres front asymetrie sur H pour time=',time |
---|
976 | !!$ stop |
---|
977 | !!$else |
---|
978 | !!$ write(6,*) 'dans front apres front pas d asymetrie sur H pour time=',time |
---|
979 | !!$ |
---|
980 | !!$end if |
---|
981 | |
---|
982 | ! on ne compte pas les taches de glace de 2 cases (horizontales ou verticales) |
---|
983 | ! en fait, si ces deux cases sont flottantes, il faut enlever les icebergs |
---|
984 | ! de n'importe quelle taille). |
---|
985 | ! si ces deux taches sont posées (ou une des deux), il n'y a pas assez de conditions aux limites |
---|
986 | |
---|
987 | |
---|
988 | do j=1,ny |
---|
989 | do i=1,nx-1 |
---|
990 | if (front(i,j).eq.1) then |
---|
991 | if (front(i+1,j).eq.1) then |
---|
992 | ice(i,j)=0 |
---|
993 | ice(i+1,j)=0 |
---|
994 | front(i,j)=0 |
---|
995 | front(i+1,j)=0 |
---|
996 | endif |
---|
997 | endif |
---|
998 | end do |
---|
999 | end do |
---|
1000 | |
---|
1001 | do j=1,ny-1 |
---|
1002 | do i=1,nx |
---|
1003 | if (front(i,j).eq.1) then |
---|
1004 | if (front(i,j+1).eq.1) then |
---|
1005 | ice(i,j)=0 |
---|
1006 | ice(i,j+1)=0 |
---|
1007 | front(i,j)=0 |
---|
1008 | front(i,j+1)=0 |
---|
1009 | endif |
---|
1010 | end if |
---|
1011 | end do |
---|
1012 | end do |
---|
1013 | |
---|
1014 | !isolx signifie pas de voisins en x |
---|
1015 | !isoly signifie pas de voisins en y |
---|
1016 | !remarque : |
---|
1017 | !si isolx/y=.true. alors frontfacex/y=0 (a la fois +1 & -1 or +1-1=0) |
---|
1018 | |
---|
1019 | ! calcul de frontfacex et isolx |
---|
1020 | do j=1,ny |
---|
1021 | do i=2,nx-1 |
---|
1022 | |
---|
1023 | if (front(i,j).ge.1.and.front(i,j).le.3) then !front(entre 1 et 3) |
---|
1024 | |
---|
1025 | if ((ice(i-1,j)+ice(i+1,j)).lt.2) then ! il y a un front // a x |
---|
1026 | |
---|
1027 | if ((ice(i-1,j)+ice(i+1,j)).eq.0) then |
---|
1028 | isolx(i,j)=.true. |
---|
1029 | elseif (ice(i-1,j).eq.0) then |
---|
1030 | frontfacex(i,j)=-1 ! front i-1 |i i+1 |
---|
1031 | else |
---|
1032 | frontfacex(i,j)=+1 ! front i-1 i| i+1 |
---|
1033 | endif |
---|
1034 | endif |
---|
1035 | end if !fin du test il y a un front |
---|
1036 | |
---|
1037 | end do |
---|
1038 | end do |
---|
1039 | |
---|
1040 | ! calcul de frontfacey et isoly |
---|
1041 | do j=2,ny-1 |
---|
1042 | do i=1,nx |
---|
1043 | |
---|
1044 | if (front(i,j).ge.1.and.front(i,j).le.3) then !front(entre 1 et 3) |
---|
1045 | |
---|
1046 | if ((ice(i,j-1)+ice(i,j+1)).lt.2) then ! il y a un front // a y |
---|
1047 | |
---|
1048 | if ((ice(i,j-1)+ice(i,j+1)).eq.0) then |
---|
1049 | isoly(i,j)=.true. !front j-1 |j| j+1 |
---|
1050 | elseif (ice(i,j-1).eq.0) then |
---|
1051 | frontfacey(i,j)=-1 !front j-1 |j j+1 |
---|
1052 | else |
---|
1053 | frontfacey(i,j)=+1 !front j-1 j| j+1 |
---|
1054 | endif |
---|
1055 | endif |
---|
1056 | end if !fin du test il y a un front |
---|
1057 | |
---|
1058 | end do |
---|
1059 | end do |
---|
1060 | |
---|
1061 | |
---|
1062 | |
---|
1063 | ! traitement des bords. On considere que l'exterieur n'a pas de glace |
---|
1064 | ! attention ce n'est vrai que sur la grande grille |
---|
1065 | |
---|
1066 | |
---|
1067 | do j=2,ny-1 |
---|
1068 | i=1 |
---|
1069 | if (front(i,j).ge.1) then |
---|
1070 | if (ice(i+1,j).eq.0) then |
---|
1071 | isolx(i,j)=.true. |
---|
1072 | else |
---|
1073 | frontfacex(i,j)=-1 |
---|
1074 | endif |
---|
1075 | end if |
---|
1076 | i=nx |
---|
1077 | if (front(i,j).ge.1) then |
---|
1078 | if (ice(i-1,j).eq.0) then |
---|
1079 | isolx(i,j)=.true. |
---|
1080 | else |
---|
1081 | frontfacex(i,j)=1 |
---|
1082 | endif |
---|
1083 | end if |
---|
1084 | end do |
---|
1085 | |
---|
1086 | do i=2,nx-1 |
---|
1087 | j=1 |
---|
1088 | if (front(i,j).ge.1) then |
---|
1089 | if (ice(i,j+1).eq.0) then |
---|
1090 | isoly(i,j)=.true. |
---|
1091 | else |
---|
1092 | frontfacey(i,j)=-1 |
---|
1093 | endif |
---|
1094 | end if |
---|
1095 | j=ny |
---|
1096 | if (front(i,j).ge.1) then |
---|
1097 | if (ice(i,j-1).eq.0) then |
---|
1098 | isoly(i,j)=.true. |
---|
1099 | else |
---|
1100 | frontfacey(i,j)=1 |
---|
1101 | endif |
---|
1102 | end if |
---|
1103 | end do |
---|
1104 | |
---|
1105 | |
---|
1106 | |
---|
1107 | return |
---|
1108 | end subroutine determin_front |
---|
1109 | !------------------------------------------------------------------------------ |
---|
1110 | |
---|
1111 | end module flottab_mod |
---|
1112 | |
---|
1113 | |
---|