source: trunk/SOURCES/flottab2-0.7.f90 @ 182

Last change on this file since 182 was 180, checked in by aquiquet, 6 years ago

Cleaning-up: unused variables removed

File size: 42.4 KB
Line 
1!> \file flottab2-0.7.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!! @todo "nasty Island". If the bedrock is above sealevel force grounded (mk_init)
15!<
16module flottab_mod
17 
18  !$ USE OMP_LIB
19USE module3D_phy
20use module_choix
21
22
23
24 IMPLICIT NONE
25
26 real ::  surnet !< surnet hauteur de glace au dessus de la mer
27 real ::  archim !< test de flottaison
28! real, parameter :: Hmin=1.001   !< Hmin pour être considere comme point ice
29
30 integer:: itestf
31
32 logical,dimension(nx,ny) ::  gz1mx,gz1my
33 logical,dimension(nx,ny) ::  fl1mx,fl1my
34
35
36 real,dimension(nx,ny) ::  uxs1   !< uxbar a l'entree de flottab
37 real,dimension(nx,ny) ::  uys1   !< uybar a l'entree de flottab
38
39
40 integer pmx,pmy !pm=plus-moins -1 ou 1 pour x et y
41
42
43! Variables pour la determination des differents shelfs/stream
44! (representés comme des taches ou l'on resoud l'eq elliptique)
45!________________________________________________________________
46 integer,parameter               :: n_ta_max=2000!< nombre de tache max
47 integer,dimension(nx,ny)        :: table_out    !< pour les numeros des taches
48 integer,dimension(nx,ny)        :: tablebis     !< pour les numeros des taches
49 integer,dimension(0:n_ta_max)     :: compt        !< contient les equivalence entre les taches
50 integer,dimension(0:n_ta_max)     :: nb_pts_tache !< indique le nombre de points par tache
51 logical,dimension(0:n_ta_max)     :: iceberg1D      !< T si iceberg, F si calotte posee
52 
53 logical,dimension(nx,ny)        :: mask_tache_ij !< masque de travail
54                                                  !< vrai pour toute la tache de i,j
55 integer,dimension(2)            :: smax_coord    !< pour le maxloc des iles
56
57!  Variables pour determiner le point le plus haut (surf)
58!  d'une ile completement stream
59!_________________________________________________________
60
61! icetrim : T si ice stream, F si calotte posee(vertical shear)
62 logical,dimension(n_ta_max) :: icetrim 
63
64 integer ::  ii,jj
65 integer ::  smax_i
66 integer ::  smax_j
67 real    ::  smax_
68 integer ::  numtache
69 integer ::  nb_pt
70 real    ::  petit_H=0.001 ! pour test ice sur zone flottante
71contains
72! -----------------------------------------------------------------------------------
73!> SUBROUTINE: flottab()
74!! Cette routine determine les endroits ou la glace
75!! flotte , les iles, et la position du front de glace
76!! @note Il y a 4 sortes de zone
77!! @note    - Pose
78!  @note    - Grounding zone et streams    gzmx et gzmy
79!  @note    - Iles             ilemx, ilemy
80!  @note    - flottant  flot sur le noeud majeur, flotmx sur le noeud mineur
81!>
82  subroutine flottab()
83    !
84    !                                   Vince 5 Jan 95
85    !                                       Modifie 20 Jan 95
86    !                                       Modifie le 30 Novembre 98
87    !    Passage f90 + determination des fronts Vincent dec 2003
88    !    nettoyage et nouvelle détermination de gzmx et gzmy Cat le 10 juillet 2005
89    !    Re-nettoyage par Cat en aout 2006.
90    !    Le calcul de gzmx et gzmy pour les points intérieurs passe dans la subroutine dragging
91    !    pour les points cotiers, toujours fait dans flottab
92    !
93    !                                       -----------------
94    !
95    !      Cette routine determine les endroits ou la glace
96    !      flotte , les iles, et la position du front de glace
97    !
98    !      Il y a 4 sortes de zone
99    !      Pose
100    !      Grounding zone et streams    gzmx et gzmy
101    !      Iles             ilemx, ilemy
102    !      flottant  flot sur le noeud majeur, flotmx sur le noeud mineur
103
104    !   passage dans flottab tous les pas de temps dt  !
105    !
106    ! _________________________________________________________
107
108    !~ print*,'debut flottab',S(132,183),H(132,183),BSOC(132,183),B(132,183),sealevel
109    !~ print*,'debut flottab',flot(132,183),ice(132,183)
110!print*,'H(90,179) flottab 1',H(90,179),ice(90,179), flot(90,179)
111
112    if (itracebug.eq.1)  call tracebug(' Entree dans routine flottab')
113
114    SHELFY = .FALSE.
115
116
117    ! cas particulier des runs paleo ou on impose un masque grounded
118
119    !$OMP PARALLEL PRIVATE(archim,surnet)
120    if (igrdline.eq.2) then
121       !$OMP WORKSHARE
122       where ( mk_init(:,:).eq.1)                 ! pose
123          flot(:,:) = .False.
124          H(:,:)=max(H(:,:),(10.+sealevel-Bsoc(:,:))*row/ro)  ! pour avoir archim=10
125          S(:,:) = Bsoc(:,:) + H(:,:) 
126       elsewhere
127          flot(:,:) = .True.
128       end where
129       !$OMP END WORKSHARE
130    end if
131
132
133    ! 1-INITIALISATION
134    ! ----------------
135    ! initialisation des variables pour detecter les points qui se mettent
136    ! a flotter entre 2 dtt
137
138    appel_new_flot=.false.
139    !$OMP DO
140    do j=1,ny
141       do i=1,nx
142          new_flot_point(i,j)=.false.
143          new_flotmx(i,j)=.false.
144          new_flotmy(i,j)=.false.
145       enddo
146    enddo
147    !$OMP END DO
148
149    ! ICE(:,:)=(H(:,:).gt.1) ! ice=.true. si epaisseur > 1m
150
151    !$OMP WORKSHARE
152    ICE(:,:)=0
153    front(:,:)=0
154    frontfacex(:,:)=0
155    frontfacey(:,:)=0
156    isolx(:,:)=.false.
157    isoly(:,:)=.false.
158    cotemx(:,:)=.false.
159    cotemy(:,:)=.false.
160    boost=.false.
161    iceberg(:,:)=.false.
162    !$OMP END WORKSHARE
163
164    ! fin de l'initialisation
165    !_____________________________________________________________________
166
167    ! 2-TESTE LES NOUVEAUX POINTS FLOTTANTS
168    ! -------------------------------------
169
170    !$OMP DO
171    do j=1,ny
172       do i=1,nx
173
174          uxs1(i,j)=uxbar(i,j)
175          uys1(i,j)=uybar(i,j)
176
177          archim = Bsoc(i,j)+H(i,j)*ro/row -sealevel
178          !      if ((i.eq.132).and.(j.eq.183)) print*,'archim=',archim
179
180
181          arch: if ((ARCHIM.LT.0.).and.(H(I,J).gt.1.e-3)) then    ! le point flotte
182             mk(i,j)=1
183
184
185             ex_pose: if ((.not.FLOT(I,J)).and.(isynchro.eq.1)) then  !  il ne flottait pas avant
186                FLOT(I,J)=.true.
187                BOOST=.false.
188
189                if (igrdline.eq.1) then   ! en cas de grounding line prescrite
190                   flot(i,j)=.false.
191                   H(i,j)=(10.+sealevel-Bsoc(i,j))*row/ro  ! pour avoir archim=10
192                   new_flot_point(i,j)=.false.
193                endif
194
195             else                                       ! isynchro=0 ou il flottait déja
196
197                if (.not.FLOT(I,J)) then               ! il ne flottait pas (isynchro=0)
198                   new_flot_point(i,j)=.true.          ! signale un point qui ne flottait pas
199                   ! au pas de temps precedent
200                   flot(i,j)=.true.
201
202                   if (igrdline.eq.1) then   ! en cas de grounding line prescrite
203                      flot(i,j)=.false.
204                      H(i,j)=(10.+sealevel-Bsoc(i,j))*row/ro  ! pour avoir archim=10
205                      new_flot_point(i,j)=.false.
206                   endif
207
208                endif
209             endif ex_pose
210
211
212          else if ((H(i,j).gt.0.).and.(archim.GE.0.)) then    !    le point ne flotte pas et est englace
213             mk(i,j)=0
214
215             if(FLOT(I,J)) then         !  mais il flottait avant
216                FLOT(I,J)=.false.
217                BOOST=.false.
218             endif
219             !cdc correction topo pour suivre  variations sealevel
220             !cdd           S(i,j)=Bsoc(i,j)+H(i,j)
221             S(i,j)=Bsoc(i,j)+H(i,j)
222             B(i,j)=Bsoc(i,j)
223
224          else if  ((H(i,j).LE.0.).and.(archim.LT.0.)) then    !    terre deglace qui devient ocean
225!cdc             ice(i,j)=0
226!cdc         H(i,j)=1.
227!cdc  1m           H(i,j)=min(1.,max(0.,(sealevel - Bsoc(i,j))*row/ro-0.01))
228             flot(i,j)=.true.  !cdc points ocean sont flot meme sans glace
229             H(i,j)=0.
230             S(i,j)=sealevel
231             B(i,j)=S(i,j)-H(i,j)
232          endif arch
233
234          !   Si la glace flotte -> PRUDENCE !!!
235          !   S et B sont alors determines avec la condition de flottabilite
236
237          if (flot(i,j)) then
238             shelfy = .true.
239
240             surnet=H(i,j)*(1.-ro/row)
241             S(i,j)=surnet+sealevel
242             B(i,j)=S(i,j)-H(i,j)
243          end if
244
245       end do
246    end do
247    !$OMP END DO
248    !~ print*,'flottab 2',S(132,183),H(132,183),BSOC(132,183),B(132,183),sealevel
249    !~ print*,'flottab 2',flot(132,183),ice(132,183)
250    !~ if (S(132,183).LT.sealevel) print*,'BUGGGGGGGGGGGGG !!!!!!!!!!!!!!!!'
251       
252!!$ do i=1,nx
253!!$    do j=1,ny
254!!$       if (flot(i,j)) then
255!!$          mk(i,j)=1
256!!$       else
257!!$          mk(i,j)=0
258!!$       endif
259!!$    end do
260!!$ end do
261
262
263
264
265    !-----------------------------------------------------------------------
266    !$OMP DO
267    domain_x: do j=1,ny       
268       do i=2,nx
269
270          !    3_x    A- NOUVELLE DEFINITION DE FLOTMX, LES POINTS
271          !            AYANT UN DES VOISINS FLOTTANTS SONT FLOTMX ET GZMX
272          !         -------------------------------------------
273
274          !            fl1 est l'ancienne valeur de flotmx
275          gz1mx(i,j)=gzmx(i,j)         !     gz1 est l'ancienne valeur de gzmx
276          fl1mx(i,j)=flotmx(i,j)       !     fl1 est l'ancienne valeur de flotmx
277
278          flotmx(i,j)=flot(i,j).and.flot(i-1,j)
279
280          ! test pour detecter les nouveaux flotmx entre 2 dtt :
281
282          if (flotmx(i,j).and.(new_flot_point(i,j).or. &
283               new_flot_point(i-1,j))) then
284             appel_new_flot=.true.
285             new_flotmx(i,j)=.true.
286          endif
287
288
289          !   premiere determination de gzmx
290          !__________________________________________________________________________
291
292          ! gzmx si un des deux voisins est flottant et l'autre posé
293          !                                                                   i-1   i
294          gzmx(i,j)=((flot(i,j).and..not.flot(i-1,j)) & !         F    P
295               .or.(.not.flot(i,j).and.flot(i-1,j)))         !         P    F
296
297          !                 A condition d'etre assez proche de la flottaison
298          !                 sur le demi noeud condition archim < 100 m
299
300          archim=(Bsoc(i,j)+Bsoc(i-1,j))*0.5-sealevel+ro/row*Hmx(i,j)
301          gzmx(i,j)=gzmx(i,j).and.(archim.le.100.) 
302          cotemx(i,j)=gzmx(i,j)
303
304       end do
305    end do domain_x
306    !$OMP END DO
307    !if (itracebug.eq.1)  call tracebug('  routine flottab apres domain_x')
308
309    !     3_y    B- NOUVELLE DEFINITION DE FLOTMY
310    !         --------------------------------
311    !$OMP DO
312    domain_y: do j=2,ny       
313       do i=1,nx
314
315          gz1my(i,j)=gzmy(i,j)           !       gz1 est l'ancienne valeur de gzmy
316          fl1my(i,j)=flotmy(i,j)         !       fl1 est l'ancienne valeur de flotmy
317
318          flotmy(i,j)=flot(i,j).and.flot(i,j-1)
319
320          ! test pour detecter les nouveaux flotmy entre 2 dtt :
321
322          if (flotmy(i,j).and.(new_flot_point(i,j).or.     &
323               new_flot_point(i,j-1))) then
324             appel_new_flot=.true.
325             new_flotmy(i,j)=.true.
326          endif
327
328          !   premiere determination de gzmy
329          !__________________________________________________________________________
330
331          ! gzmy si un des deux voisins est flottant et l'autre posé
332
333          gzmy(i,j)=((flot(i,j).and..not.flot(i,j-1))          &
334               .or.(.not.flot(i,j).and.flot(i,j-1)))                 
335
336          !                 A condition d'etre assez proche de la flottaison
337          !                 sur le demi noeud condition archim > 100 m
338
339          archim=(Bsoc(i,j)+Bsoc(i,j-1))*0.5-sealevel+ro/row*Hmy(i,j)
340          gzmy(i,j)=gzmy(i,j).and.(archim.le.100.) 
341          cotemy(i,j)=gzmy(i,j)
342
343       end do
344    end do domain_y
345    !$OMP END DO
346
347
348    !-------------------------------------------------------------------------------------
349    ! attention : pour expériences Heino
350    !             gzmy(i,j)=gzmy_heino(i,j)
351    !             shelfy=.true.
352    !             shelfy=.false.
353    !____________________________________________________________________________________
354
355
356!!$
357!!$
358!!$! 4- Condition sur les bords
359!!$
360!!$     
361!!$       do i=2,nx
362!!$         flotmx(i,1) = (flot(i,1).or.flot(i-1,1))
363!!$         flotmy(i,1) = .false.
364!!$       end do
365!!$       
366!!$       do j=2,ny
367!!$          flotmy(1,j) = (flot(1,j).or.flot(1,j-1))
368!!$          flotmx(1,j) = .false.
369!!$       end do
370!!$
371!!$       flotmx(1,1) = .false.
372!!$       flotmy(1,1) = .false.
373!!$
374
375
376    !      4- determination des iles
377    !      -------------------------
378    !$OMP WORKSHARE
379    ilemx(:,:)=.false.
380    ilemy(:,:)=.false.
381    !$OMP END WORKSHARE
382
383    !       selon x
384    !$OMP DO
385    ilesx:  do j=2,ny-1
386       do i=3,nx-2
387          !                       F   G   F   
388          !                           x
389          ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
390          if ((flot(i-1,j).and..not.flot(i,j).and.flot(i+1,j)).and. &
391               (sdx(i,j).LT.1.E-02)) then
392             ilemx(i,j)=.true.
393             ilemx(i+1,j)=.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-1,j).and..not.flot(i,j)                  &
399               .and..not.flot(i+1,j)).and.flot(i+2,j).and.           &
400               (sdx(i,j).LT.1.E-02.and.sdx(i+1,j).LT.1.E-02)) then
401             ilemx(i,j)=.true.
402             ilemx(i+1,j)=.true.
403             ilemx(i+2,j)=.true.
404
405             !                      F   G   G   F
406             !                              x
407             ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
408          else if ((flot(i-2,j).and..not.flot(i-1,j)                &
409               .and..not.flot(i,j)).and.flot(i+1,j).and.             &
410               (sdx(i,j).LT.1.E-02.and.sdx(i-1,j).LT.1.E-02)) then
411             ilemx(i-1,j)=.true.
412             ilemx(i,j)=.true.
413             ilemx(i+1,j)=.true.
414
415             !                      F   G   G   G    F
416             !                              x
417             ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
418          else if ((i.lt.nx-2)                                      &
419               .and.(flot(i-2,j).and..not.flot(i-1,j)                &
420               .and..not.flot(i,j)).and..not.flot(i+1,j)             &
421               .and.flot(i+2,j).and.                                &
422               (sdx(i,j).LT.1.E-02.and.sdx(i-1,j).LT.1.E-02         &
423               .and.sdx(i+1,j).LT.1.E-02)) then
424             ilemx(i-1,j)=.true.
425             ilemx(i,j)=.true.
426             ilemx(i+1,j)=.true.
427             ilemx(i+2,j)=.true.
428
429          endif
430
431       end do
432    end do ilesx
433    !$OMP END DO
434
435    !       selon y
436    !$OMP DO
437    ilesy: do j=3,ny-2
438       do i=2,nx-1
439          !                       F   G   F   
440          !                           x
441          ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
442          if ((flot(i,j-1).and..not.flot(i,j).and.flot(i,j+1)).and. &
443               (sdy(i,j).LT.1.E-02)) then
444             ilemy(i,j)=.true.
445             ilemy(i,j+1)=.true.
446
447             !                      F  G   G   F
448             !                         x
449             ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
450          else if ((flot(i,j-1).and..not.flot(i,j)                  &
451               .and..not.flot(i,j+1)).and.flot(i,j+2).and.           &
452               (sdy(i,j).LT.1.E-02.and.sdy(i,j+1).LT.1.E-02)) then
453             ilemy(i,j)=.true.
454             ilemy(i,j+1)=.true.
455             ilemy(i,j+2)=.true.
456
457             !                      F   G   G   F
458             !                              x
459             ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
460          else if ((flot(i,j-2).and..not.flot(i,j-1)                 &
461               .and..not.flot(i,j)).and.flot(i,j+1).and.              &
462               (sdy(i,j).LT.1.E-02.and.sdy(i,j-1).LT.1.E-02)) then
463             ilemy(i,j-1)=.true.
464             ilemy(i,j)=.true.
465             ilemy(i,j+1)=.true.
466
467             !                      F   G   G   G    F
468             !                              x
469             ! modif tof 26/08/02 limite sur la pente (si diff S > 400 m)
470          else if ((j.lt.ny-2)                                      &
471               .and.(flot(i,j-2).and..not.flot(i,j-1)                &
472               .and..not.flot(i,j)).and..not.flot(i,j+1)             &
473               .and.flot(i,j+2).and.                                &
474               (sdy(i,j).LT.1.E-02.and.sdy(i,j-1).LT.1.E-02         &
475               .and.sdy(i,j+1).LT.1.E-02)) then
476             ilemy(i,j-1)=.true.
477             ilemy(i,j)=.true.
478             ilemy(i,j+1)=.true.
479             ilemy(i,j+2)=.true.
480          endif
481       end do
482    end do ilesy
483    !$OMP END DO
484    !$OMP END PARALLEL
485    ! fin des iles
486
487!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
488!!$if (itestf.gt.0) then
489!!$   write(6,*) 'dans flottab avant dragging asymetrie sur H  pour time=',time
490!!$   stop
491!!$else
492!!$   write(6,*) 'dans flottab avant dragging pas d asymetrie sur H  pour time=',time
493!!$
494!!$end if
495
496    ! 5- calcule les noeuds qui sont streams a l'interieur et donne le betamx et betamy
497    !----------------------------------------------------------------------------------
498    call dragging         
499   
500    if (itracebug.eq.1)  call tracebug(' routine flottab apres call dragging')
501!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
502!!$if (itestf.gt.0) then
503!!$   write(6,*) 'dans flottab apres dragging asymetrie sur H  pour time=',time
504!!$   stop
505!!$else
506!!$   write(6,*) 'dans flottab aapres dragging pas d asymetrie sur H  pour time=',time
507!!$
508!!$end if
509
510    !   6- calcule les vitesses des points qui sont devenus gzm
511    !$OMP PARALLEL
512    !$OMP DO         
513    do j=1,ny
514       do i=2,nx-1
515          !      si le point etait posé (non gz) et devient gzmx
516          !      definir la direction de la vitesse (moyenne des points)
517
518          if ((.not.gz1mx(i,j)).and.(.not.fl1mx(i,j)).and.gzmx(i,j).and. &
519               (i.gt.2).and.(i.lt.nx))  then
520             uxs1(i,j)=(uxbar(i+1,j)+uxbar(i-1,j))/2.
521          endif
522
523       end do
524    end do
525    !$OMP END DO
526
527    !$OMP DO
528    do j=2,ny-1
529       do i=1,nx
530          !      si le point etait posé (non gz) et devient gzmy
531          !      definir la direction de la vitesse (moyenne des points)
532          if ((.not.gz1my(i,j)).and.(.not.fl1my(i,j)).and.gzmy(i,j).and. &
533               (j.gt.2).and.(j.lt.ny))  then
534             uys1(i,j)=(uybar(i,j+1)+uybar(i,j-1))/2.
535          endif
536
537       end do
538    end do
539    !$OMP END DO
540
541
542    !  7-On determine finalement la position des noeuds stream ou shelf
543    !  -------------------------------------------------------------
544
545    if (nt.ge.2) then   ! pour ne pas faire ce calcul lors du premier passage
546       !$OMP WORKSHARE
547       uxbar(:,:)=uxs1(:,:)
548       uybar(:,:)=uys1(:,:)
549       !$OMP END WORKSHARE
550    endif
551
552    !$OMP WORKSHARE
553    flgzmx(:,:)=(marine.and.(flotmx(:,:).or.gzmx(:,:).or.ilemx(:,:)))   &
554         .or.(.not.marine.and.flotmx(:,:))
555    where (hmx(:,:).eq.0.)
556      flgzmx(:,:) = .false.
557    endwhere
558    flgzmy(:,:)=(marine.and.(flotmy(:,:).or.gzmy(:,:).or.ilemy(:,:)))   &
559         .or.(.not.marine.and.flotmy(:,:))
560    where (hmy(:,:).eq.0.)
561      flgzmy(:,:) = .false.
562    endwhere
563    !$OMP END WORKSHARE
564
565
566    ! 8- Pour la fusion basale sous les ice shelves- region proche de la grounding line
567    !---------------------------------------------------------------------------------
568    !       fbm est vrai si le point est flottant mais un des voisins est pose
569    !_________________________________________________________________________
570    !$OMP DO
571    do j=2,ny-1
572       do i=2,nx-1
573
574          ! if (i.gt.2.AND.i.lt.nx) then   
575          fbm(i,j)=flot(i,j).and.      & 
576               ((.not.flot(i+1,j)).or.(.not.flot(i,j+1)) &
577               .or.(.not.flot(i-1,j)).or.(.not.flot(i,j-1)))
578          !     endif
579       end do
580    end do
581    !$OMP END DO
582
583
584    !  9-On determine maintenant la position du front de glace
585    !  -------------------------------------------------------
586    !  C'est ici que l'on determine la position et le type de front
587    !       print*,'on est dans flottab pour definir les fronts'
588
589
590
591!!$do i=3,nx-2
592!!$   do j=3,ny-2
593!!$      if (h(i,j).gt.1.1) ice(i,j)=1       
594!!$   end do
595!!$end do
596!                       print*, 'flolottab debug', H(71,25),flot(71,25),bm(71,25),ice(71,25),time
597!print*,'H(90,179) flottab 2',H(90,179),ice(90,179),flot(90,179)
598    !$OMP WORKSHARE
599    where (flot(:,:))
600!cdc 1m       where (H(:,:).gt.max(Hmin,Hmin+BM(:,:)-Bmelt(:,:)))
601      where (H(:,:).gt.max(BM(:,:)-Bmelt(:,:)+petit_H,0.)*dt)
602!cdc                     where (H(:,:).gt.0.)
603          ice(:,:)=1
604       elsewhere 
605          ice(:,:)=0
606       end where
607    elsewhere
608       where (H(:,:).gt.0.) 
609          ice(:,:)=1
610       elsewhere 
611          ice(:,:)=0
612       end where
613    end where
614    !$OMP END WORKSHARE
615    !$OMP END PARALLEL
616!print*,'H(90,179) flottab 3',H(90,179),ice(90,179),flot(90,179),Bsoc(90,179)+H(90,179)*ro/row -sealevel
617 
618    !call determin_front ! cette version ne conserve pas la masse !!!
619    call determin_front_tof ! version simplifiee
620
621    ! pour sorties initMIP:
622    debug_3D(:,:,118) = ice(:,:)*(1-mk(:,:))
623    debug_3D(:,:,119) = ice(:,:)*mk(:,:)
624
625   
626  end subroutine flottab
627!-------------------------------------------------------------------- 
628
629!> SUBROUTINE: determin_tache
630!! Routine pour la dtermination du numero de tache a effectuer
631!>
632subroutine determin_tache
633
634!$ USE OMP_LIB
635
636implicit none
637integer :: indice
638integer :: label         ! no des taches rencontrées dans le mask
639integer :: label_max     ! no temporaire maxi de tache rencontrées
640!      integer :: mask_nb = 4
641integer,parameter :: mask_nb = 2   ! version ou on ne compte pas les diagonales
642!     integer,dimension(mask_nb) :: mask
643integer,dimension(mask_nb) :: mask
644
645
646! 1-initialisation
647!-----------------
648label_max=1      ! numero de la tache, la premiere tache est notée 1
649label=1
650do i=1,n_ta_max
651   compt(i)=i
652enddo
653!      table_in  = .false.
654!$OMP PARALLEL
655!$OMP WORKSHARE
656table_out(:,:) = 0
657iceberg1D(:)  = .true.
658icetrim (:) = .true.
659nb_pts_tache(:) = 0
660!$OMP END WORKSHARE
661!$OMP END PARALLEL
662!    open(unit=100,file="tache.data",status='replace')
663
664! 2-reperage des taches
665!----------------------
666
667do j=2,ny-1
668  do i=2,nx-1
669
670     IF (ice(i,j).ge.1) THEN ! on est sur la glace-----------------------------!
671                       
672        if ((ice(i-1,j).ge.1).or.(ice(i,j-1).ge.1)) then   !masque de 2 cases adjacentes
673           !      un des voisins est deja en glace
674           mask(1) = table_out(i-1,j)
675           mask(2) = table_out(i,j-1) 
676           label = label_max
677
678           !on determine la valeur de la tache minimun (>0) presente ds le masque
679           do indice=1,mask_nb
680              if (mask(indice).gt.0) label=min(label,mask(indice))
681           enddo
682!cdc       label=min(label,minval(mask(:), mask=mask > 0))
683
684           !on fixe la valeur de la tache voisine minimun au point etudie (via label)
685           table_out(i,j)=label
686           !si ce noeud est posé, alors la tache n'est pas un iceberg et iceberg=.F.
687           if (.not.FLOT(I,J)) then
688              iceberg1D(label)=.false.
689           endif
690
691           !si ce noeud est posé, alors la tache n'est pas un ice stream et icestrim=.F.
692           if ((.not.gzmx(i,j).and..not.gzmx(i+1,j)).and.  &
693               (.not.gzmy(i,j).and..not.gzmy(i,j+1))) then
694              icetrim(label)=.false.
695           endif
696
697           ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt
698           ! on lui affecte le numero de la tache fondamentale
699
700           do indice=1,mask_nb
701              if(mask(indice).gt.label) then
702                 compt(mask(indice))=-label
703              endif
704           enddo
705                                                                   ! exemple on est sur le point X : 5  X
706          do indice=1,mask_nb                                      !                                   20
707            if(mask(indice).gt.label) then                         ! mask(2)=20 > 5             
708              compt(mask(indice))=label                            ! compt(20)=5
709              if (.not.iceberg1D(mask(indice))) iceberg1D(label)=.false. ! si la tache n'etais pas un iceberg => iceberg =.false.  iceberg(5)=.false.               
710              if (.not.icetrim(mask(indice))) icetrim(label)=.false.
711              where (table_out(:,:).eq.mask(indice))               ! where table_out(:,:)=mask(2)=20
712                table_out(:,:)=label                        ! table_out(:,:)=label=5                 
713              endwhere
714            endif
715          enddo   
716
717        else !aucun des voisins est une tache
718           table_out(i,j)= label_max
719           compt(label_max)=label_max
720           !si ce noeud est posé, alors la ache n'est pas un iceberg et iceberg=.F.
721           if (.not.FLOT(I,J)) then
722              iceberg1D(label_max)=.false.
723           endif
724
725           !si ce noeud est posé, alors le tache n'est pas un ice stream et icestrim=.F.
726           if ((.not.gzmx(i,j).and..not.gzmx(i+1,j)).and.  &
727               (.not.gzmy(i,j).and..not.gzmy(i,j+1))) then
728              icetrim(label)=.false.
729           endif
730             
731           label_max  = label_max+1
732           if (label_max.gt.n_ta_max) print*,'ATTENTION trop de taches icebergs=',label_max 
733        endif
734
735
736     else           !on est pas sur une tache----------------------------------------------
737        table_out(i,j)=0    ! Pas necessaire (reecrit 0 sur 0)             
738     endif          !---------------------------------------------------------------------
739
740
741  enddo
742enddo
743
744
745
746! On reorganise compt en ecrivant le numero de la tache fondamentale
747! i.e. du plus petit numero present sur la tache (Sans utiliser de recursivité)
748! On indique aussi le nb de point que contient chaque taches (nb_pts_tache)
749
750  !$OMP PARALLEL
751  !$OMP DO
752  do j=1,ny
753    do i=1,nx
754      if (table_out(i,j).ne.0) then
755        nb_pts_tache(compt(table_out(i,j)))= nb_pts_tache(compt(table_out(i,j)))+1
756      endif
757    enddo
758  enddo
759  !$OMP END DO
760  !$OMP END PARALLEL
761
762  !On compte comme englacé uniquement les calottes dont une partie est posée     
763  !$OMP PARALLEL PRIVATE(smax_,smax_coord,smax_i,smax_j,mask_tache_ij)
764  !$OMP DO
765  do j=3,ny-2
766    do i=3,nx-2
767      test1: if (.not.iceberg1D(table_out(i,j))) then ! on est pas sur un iceberg
768        if (nb_pts_tache(table_out(i,j)).ge.1) then
769          ice(i,j)=1
770          ! ici on est sur une tache non iceberg >= 5 points                       
771          ! on teste si la tache n'est pas completement ice stream
772   
773          test2: if (icetrim(table_out(i,j))) then   ! on a une ile d'ice stream
774   
775            mask_tache_ij(:,:)=.false.
776            mask_tache_ij(:,:)=(table_out(:,:).eq.table_out(i,j))       ! pour toute la tache
777   
778            smax_=maxval(S(:,:),MASK=mask_tache_ij(:,:))
779            smax_coord(:)=maxloc(S(:,:),MASK=mask_tache_ij(:,:))
780            smax_i=smax_coord(1)
781            smax_j=smax_coord(2)
782 
783            gzmx(smax_i,smax_j)=.false. ; gzmx(smax_i+1,smax_j)=.false.
784            gzmy(smax_i,smax_j)=.false. ; gzmx(smax_i,smax_j+1)=.false.
785            flgzmx(smax_i,smax_j)=.false. ; flgzmx(smax_i+1,smax_j)=.false.
786            flgzmy(smax_i,smax_j)=.false. ; flgzmx(smax_i,smax_j+1)=.false.
787     
788            if (Smax_.le.sealevel) then
789              write(num_tracebug,*)'Attention, une ile avec la surface sous l eau'
790              write(num_tracebug,*)'time=',time,'   coord:',smax_i,smax_j
791            end if   
792          endif test2 
793        end if                                             ! endif deplace
794!cdc transfere dans calving :             
795      else ! on est sur un iceberg                          !   test1
796        iceberg(i,j)=iceberg1D(table_out(i,j))
797!~        ice(i,j)=0
798!~        h(i,j)=0. !1. afq, we should put everything in calving!
799!~        surnet=H(i,j)*(1.-ro/row)
800!~        S(i,j)=surnet+sealevel
801!~        B(i,j)=S(i,j)-H(i,j)             
802      endif test1   
803    end do
804  end do   
805  !$OMP END DO
806  !$OMP END PARALLEL
807
808debug_3D(:,:,121)=real(table_out(:,:))
809
810end subroutine determin_tache
811!----------------------------------------------------------------------
812!> SUBROUTINE: determin_front
813!!Routine pour la determination du front
814!>
815subroutine determin_front
816!!$ USE OMP_LIB
817integer :: i_moins1,i_plus1,i_plus2
818integer :: j_moins1,j_plus1,j_plus2
819     
820      !$OMP PARALLEL
821      !$OMP DO
822      do j=3,ny-2
823        do i=3,nx-2
824
825 surice:if  (ice(i,j).eq.0) then
826         do pmx=-1,1,2
827          do pmy=-1,1,2
828
829 diagice :  if (ice(i+pmx,j+pmy).eq.1) then
830
831              if ((ice(i+pmx,j)+ice(i,j+pmy).eq.2)) then    ! test (i) pour eviter les langues
832                                                            ! de glaces diagonales en coin(26dec04)
833                   if ((ice(i+2*pmx,j).eq.1.and.ice(i+2*pmx,j+pmy).eq.0).or.&
834                       (ice(i,j+2*pmy).eq.1.and.ice(i+pmx,j+2*pmy).eq.0)) then
835                           ice(i,j)=1
836                           h(i,j)=max(1.,h(i,j))
837                   endif         
838
839! test (i) pour eviter les langues de glaces diagonales :
840!                        mouvement du cheval aux echecs 
841
842               if ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.1)) then
843                   if (ice(i+2*pmx,j+pmy).eq.1.and. &
844                      (ice(i+2*pmx,j+2*pmy)+ice(i,j+2*pmy)).ge.1) then
845                           ice(i,j)=1
846                           h(i,j)=max(1.,h(i,j))
847                   endif       
848                   if (ice(i+pmx,j+2*pmy).eq.1.and. &
849                      (ice(i+2*pmx,j+2*pmy)+ice(i+2*pmx,j)).ge.1) then
850                           ice(i,j)=1
851                           h(i,j)=max(1.,h(i,j))
852                   endif       
853
854! test (ii) pour eviter les langues de glaces diagonales :
855!            le point glace ice(i+pmx,j+pmy) a :
856!                          - ses 4 voisins frontaux en glace
857!                          - mais 2 voisins vides diagonalement opposes   
858           
859               elseif ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.2)  &
860                                .and.ice(i+2*pmx,j+2*pmy).eq.0) then
861
862! test (iii) pour faire les tests (i) et (ii)
863                    ice(i,j)=1
864                    h(i,j)=max(1.,h(i,j))
865                    ice(i+2*pmx,j+2*pmy)=1
866                    h(i+2*pmx,j+2*pmy)=max(1.,h(i+2*pmx,j+2*pmy))
867               endif
868            endif
869           endif diagice   
870          enddo
871         enddo
872         endif surice       
873       end do
874      end do
875      !$OMP END DO
876      !$OMP ENd PARALLEL
877
878!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
879!!$if (itestf.gt.0) then
880!!$   write(6,*) 'dans front avant remplissage baies asymetrie sur H  pour time=',time
881!!$   stop
882!!$else
883!!$   write(6,*) 'dans front avant remplissage baies pas d asymetrie sur H  pour time=',time
884!!$
885!!$end if
886
887
888!     print*,'dans remplissage baies',time
889       
890baies: do k=1,2
891!$OMP PARALLEL
892!$OMP DO PRIVATE(i_moins1,j_moins1,i_plus1,j_plus1,i_plus2,j_plus2)
893do j=1,ny
894   do i=1,nx
895
896surice_xy: if  (ice(i,j).eq.0) then
897   i_moins1=max(i-1,1)
898   j_moins1=max(j-1,1)
899   i_plus1=min(i+1,nx)
900   j_plus1=min(j+1,ny)
901   i_plus2=min(i+2,nx)
902   j_plus2=min(j+2,ny)
903 
904! test (iii) pour trouver les baies de largeur 1 ou 2 cases
905! et combler les trous si ce sont des baies
906! si ce ne sont pas des baies, ne pas combler et creer des langues de glaces artificielles             
907! baies horizontales 
908 
909         if (ice(i_moins1,j).eq.1.and.(ice(i_plus1,j).eq.1.or.ice(i_plus2,j).eq.1)) then
910            if (ice(i,j_moins1).eq.1.or.ice(i,j_plus1).eq.1)  then    ! ice(i,j)=1
911               ice(i,j)=1
912               H(i,j)=max(1.,H(i,j))
913            endif
914         endif
915
916
917         if (ice(i,j_moins1).eq.1.and.(ice(i,j_plus1).eq.1.or.ice(i,j_plus2).eq.1)) then
918            if (ice(i_moins1,j).eq.1.or.ice(i_plus1,j).eq.1) then !ice(i,j)=1
919               ice(i,j)=1
920               H(i,j)=max(1.,H(i,j))
921            endif
922         endif
923
924      endif surice_xy
925   end do
926end do
927!$OMP END DO
928!$OMP END PARALLEL
929end do baies
930
931!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
932!!$if (itestf.gt.0) then
933!!$   write(6,*) 'dans front apres remplissage baies asymetrie sur H  pour time=',time
934!!$   stop
935!!$else
936!!$   write(6,*) 'dans front apres remplissage baies pas d asymetrie sur H  pour time=',time
937!!$
938!!$end if
939
940!$OMP PARALLEL
941!$OMP DO
942do j=2,ny-1
943   do i=2,nx-1
944
945      if (ice(i,j).eq.1) then     !         test si ice=1
946
947! if ice, on determine front...
948! ainsi, front=0 sur les zones = 0
949
950         front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
951         !front= le nb de faces en contact avec un voisin englacé
952      endif
953   end do
954end do
955!$OMP END DO
956
957! traitement des bords. On considere que l'exterieur n'a pas de glace
958! attention ce n'est vrai que sur la grande grille
959
960!$OMP DO PRIVATE(i)
961do j=2,ny-1
962   i=1
963   front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
964   i=nx
965   front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1))
966end do
967!$OMP END DO
968
969!$OMP DO PRIVATE(j)
970do i=2,nx-1
971   j=1 
972   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1))
973   j=ny
974   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1))
975end do
976!$OMP END DO
977
978! traitement des coins
979
980front(1,1)=ice(2,1)+ice(2,1)
981front(1,ny)=ice(2,ny)+ice(1,ny-1)
982front(nx,1)=ice(nx,2)+ice(nx-1,1)
983front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny)
984
985!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
986!!$if (itestf.gt.0) then
987!!$   write(6,*) 'dans front apres front asymetrie sur H  pour time=',time
988!!$   stop
989!!$else
990!!$   write(6,*) 'dans front apres front pas d asymetrie sur H  pour time=',time
991!!$
992!!$end if
993
994!   on ne compte pas les taches de glace de 2 cases (horizontales ou verticales)
995!   en fait, si ces deux cases sont flottantes, il faut enlever les icebergs
996!   de n'importe quelle taille).
997!   si ces deux taches sont posées (ou une des deux), il n'y a pas assez de conditions aux limites
998
999!$OMP DO
1000do j=1,ny
1001   do i=1,nx-1
1002      if (front(i,j).eq.1) then
1003         if (front(i+1,j).eq.1) then
1004            ice(i,j)=0
1005            ice(i+1,j)=0
1006            front(i,j)=0
1007            front(i+1,j)=0
1008         endif
1009      endif
1010   end do
1011end do
1012!$OMP END DO
1013
1014!$OMP DO
1015do j=1,ny-1
1016   do i=1,nx
1017      if (front(i,j).eq.1) then
1018         if (front(i,j+1).eq.1) then
1019            ice(i,j)=0
1020            ice(i,j+1)=0
1021            front(i,j)=0
1022            front(i,j+1)=0
1023         endif
1024      end if
1025   end do
1026end do
1027!$OMP END DO
1028
1029!isolx signifie pas de voisins en x
1030!isoly signifie pas de voisins en y
1031!remarque :
1032!si isolx/y=.true. alors frontfacex/y=0 (a la fois +1 & -1 or +1-1=0)
1033
1034! calcul de frontfacex et isolx
1035!$OMP DO
1036do j=1,ny
1037   do i=2,nx-1
1038
1039      if (front(i,j).ge.1.and.front(i,j).le.3) then    !front(entre 1 et 3)
1040         
1041         if ((ice(i-1,j)+ice(i+1,j)).lt.2) then        ! il y a un front // a x
1042         
1043            if ((ice(i-1,j)+ice(i+1,j)).eq.0) then
1044               isolx(i,j)=.true.
1045            elseif (ice(i-1,j).eq.0) then
1046               frontfacex(i,j)=-1                      ! front  i-1 |i  i+1
1047            else
1048               frontfacex(i,j)=+1                      ! front  i-1  i| i+1
1049            endif
1050         endif
1051      end if !fin du test il y a un front
1052         
1053   end do
1054end do
1055!$OMP END DO
1056
1057! calcul de frontfacey et isoly
1058!$OMP DO
1059do j=2,ny-1
1060   do i=1,nx
1061
1062      if (front(i,j).ge.1.and.front(i,j).le.3) then   !front(entre 1 et 3)
1063         
1064         if ((ice(i,j-1)+ice(i,j+1)).lt.2) then       ! il y a un front // a y
1065             
1066            if ((ice(i,j-1)+ice(i,j+1)).eq.0) then
1067               isoly(i,j)=.true.                      !front   j-1 |j| j+1
1068            elseif (ice(i,j-1).eq.0) then
1069               frontfacey(i,j)=-1                     !front   j-1 |j j+1
1070            else
1071               frontfacey(i,j)=+1                     !front   j-1  j| j+1
1072            endif
1073         endif
1074     end if !fin du test il y a un front
1075         
1076   end do
1077end do
1078!$OMP END DO
1079
1080
1081! traitement des bords. On considere que l'exterieur n'a pas de glace
1082! attention ce n'est vrai que sur la grande grille
1083
1084!$OMP DO PRIVATE(i)
1085do j=2,ny-1
1086   i=1
1087   if (front(i,j).ge.1)  then
1088      if (ice(i+1,j).eq.0) then
1089         isolx(i,j)=.true.
1090      else
1091         frontfacex(i,j)=-1
1092      endif
1093   end if
1094   i=nx
1095   if (front(i,j).ge.1)  then
1096      if (ice(i-1,j).eq.0) then
1097         isolx(i,j)=.true.
1098      else
1099         frontfacex(i,j)=1
1100      endif
1101   end if
1102end do
1103!$OMP END DO
1104
1105!$OMP DO PRIVATE(j)
1106do i=2,nx-1
1107   j=1 
1108   if (front(i,j).ge.1)  then
1109      if (ice(i,j+1).eq.0) then
1110         isoly(i,j)=.true.
1111      else
1112         frontfacey(i,j)=-1
1113      endif
1114   end if
1115   j=ny
1116   if (front(i,j).ge.1)  then
1117      if (ice(i,j-1).eq.0) then
1118         isoly(i,j)=.true.
1119      else
1120         frontfacey(i,j)=1
1121      endif
1122   end if
1123end do
1124!$OMP END DO
1125!$OMP END PARALLEL
1126
1127return
1128end subroutine determin_front
1129!------------------------------------------------------------------------------
1130
1131subroutine determin_front_tof
1132
1133integer,dimension(nx,ny) :: nofront ! tableau de travail (points au dela du front)
1134!$OMP PARALLEL
1135!$OMP DO
1136do j=2,ny-1
1137   do i=2,nx-1
1138
1139      if (ice(i,j).eq.1) then     !         test si ice=1
1140
1141! if ice, on determine front...
1142! ainsi, front=0 sur les zones = 0
1143
1144         front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1145         !front= le nb de faces en contact avec un voisin englacé
1146      endif
1147   end do
1148end do
1149!$OMP END DO
1150
1151! traitement des bords. On considere que l'exterieur n'a pas de glace
1152! attention ce n'est vrai que sur la grande grille
1153
1154!$OMP DO PRIVATE(i)
1155do j=2,ny-1
1156   i=1
1157   front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1158   i=nx
1159   front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1))
1160end do
1161!$OMP END DO
1162
1163!$OMP DO PRIVATE(j)
1164do i=2,nx-1
1165   j=1 
1166   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1))
1167   j=ny
1168   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1))
1169end do
1170!$OMP END DO
1171!$OMP BARRIER
1172! traitement des coins
1173
1174front(1,1)=ice(2,1)+ice(2,1)
1175front(1,ny)=ice(2,ny)+ice(1,ny-1)
1176front(nx,1)=ice(nx,2)+ice(nx-1,1)
1177front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny)
1178
1179! Les points à plus d'un point de grille du bord sont front=-1
1180!$OMP WORKSHARE
1181nofront(:,:)=0
1182!$OMP END WORKSHARE
1183!$OMP BARRIER
1184
1185!$OMP DO
1186do j=2,ny-1
1187        do i=2,nx-1
1188                if ((ice(i,j).eq.0).and.(front(i-1,j)+front(i+1,j)+front(i,j-1)+front(i,j+1)).eq.0) then
1189                        nofront(i,j)=-1
1190                endif
1191        enddo
1192enddo
1193!$OMP END DO
1194!$OMP BARRIER
1195
1196!$OMP WORKSHARE
1197where (nofront(:,:).eq.-1)
1198        front(:,:)=-1
1199endwhere
1200!$OMP END WORKSHARE
1201!$OMP END PARALLEL             
1202
1203end subroutine determin_front_tof
1204
1205
1206!> SUBROUTINE: determin_marais
1207!! afq -- Routine pour l'identification des "marais"
1208!! Un marais est une tache "shelf" entouré de points grounded
1209!! Copie sauvage de determin_tache, adapte au probleme du marais
1210!>
1211subroutine determin_marais
1212
1213  !$ USE OMP_LIB
1214
1215  implicit none
1216 
1217  integer :: indice
1218  integer :: label         ! no des taches rencontrées dans le mask
1219  integer :: label_max     ! no temporaire maxi de tache rencontrées
1220  integer,parameter :: mask_nb = 2   ! version ou on ne compte pas les diagonales
1221  integer,dimension(mask_nb) :: mask ! numero de tache des points adjacents
1222
1223  integer,dimension(nx,ny)          :: table_out_marais    !< numeros de tache d'un point ij
1224  integer,dimension(0:n_ta_max)     :: compt_marais        !< contient les equivalence entre les taches
1225  integer,dimension(0:n_ta_max)     :: nb_pts_marais       !< indique le nombre de points par tache
1226  logical,dimension(0:n_ta_max)     :: marais              !< T si flottants entoure de poses, F sinon
1227
1228
1229! 1-initialisation
1230!-----------------
1231  label_max=1      ! numero de la tache, la premiere tache est notée 1
1232  label=1
1233  do i=1,n_ta_max
1234    compt_marais(i)=i
1235  enddo
1236  !$OMP PARALLEL
1237  !$OMP WORKSHARE
1238  table_out_marais(:,:) = 0
1239  marais(:)  = .true.
1240  nb_pts_marais(:) = 0
1241  !$OMP END WORKSHARE
1242  !$OMP END PARALLEL
1243
1244! 2-reperage des taches
1245!----------------------
1246
1247  do j=2,ny-1
1248    do i=2,nx-1
1249      if ((ice(i,j).ge.1).and.flot(i,j)) then ! on est sur la glace qui flotte-------------------!
1250
1251        if (((ice(i-1,j).ge.1).and.flot(i-1,j)).or.((ice(i,j-1).ge.1).and.flot(i,j-1))) then   !masque de 2 cases adjacentes
1252          !  un des voisins est deja en glace
1253          mask(1) = table_out_marais(i-1,j)
1254          mask(2) = table_out_marais(i,j-1) 
1255          label = label_max
1256
1257          ! on determine la valeur de la tache minimun (>0) presente ds le masque
1258          do indice=1,mask_nb
1259            if (mask(indice).gt.0) label=min(label,mask(indice))
1260          enddo
1261
1262          ! on fixe la valeur de la tache voisine minimun au point etudie (via label)
1263          table_out_marais(i,j)=label
1264           
1265          !si un des voisins n'est pas glace alors la tache n'est pas un marais
1266          if ( (ice(i+1,j).eq.0) .or. (ice(i,j+1).eq.0) .or. (ice(i-1,j).eq.0) .or. (ice(i,j-1).eq.0) ) then
1267            marais(label)=.false.
1268          endif
1269
1270          ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt_marais
1271          ! on lui affecte le numero de la tache fondamentale
1272
1273                                                                   ! exemple on est sur le point X : 5  X
1274          do indice=1,mask_nb                                      !                                   20
1275            if(mask(indice).gt.label) then                         ! mask(2)=20 > 5             
1276              compt_marais(mask(indice))=label                     ! compt_marais(20)=5
1277              if (.not.marais(mask(indice))) marais(label)=.false. ! si la tache n'etais pas un marais => marais =.false.  marais(-(-5))=.false.               
1278              where (table_out_marais(:,:).eq.mask(indice))        ! where table_out_marais(:,:)=mask(2)=20
1279                table_out_marais(:,:)=label                        ! table_out_marais(:,:)=label=5                 
1280              endwhere                                               
1281            endif
1282          enddo
1283           
1284        else !aucun des voisins est une tache
1285          table_out_marais(i,j)= label_max
1286          compt_marais(label_max)=label_max
1287           
1288          ! si un des voisins n'est pas glace alors la tache n'est pas un marais
1289          if ( (ice(i+1,j).eq.0) .or. (ice(i,j+1).eq.0) .or. (ice(i-1,j).eq.0) .or. (ice(i,j-1).eq.0) ) then
1290            marais(label_max)=.false.
1291          endif
1292          label_max  = label_max+1
1293          if (label_max.gt.n_ta_max) print*,'ATTENTION trop de taches=',label_max 
1294        endif
1295      else          ! on est pas sur une tache--------------------------------------------
1296        table_out_marais(i,j)=0    ! Pas necessaire (reecrit 0 sur 0)             
1297      endif         !---------------------------------------------------------------------
1298    enddo
1299  enddo
1300
1301  marais(0)=.false.
1302
1303  !$OMP PARALLEL
1304  !$OMP DO
1305  do j=1,ny
1306    do i=1,nx
1307      if (table_out_marais(i,j).ne.0) then
1308        nb_pts_marais(compt_marais(table_out_marais(i,j)))= nb_pts_marais(compt_marais(table_out_marais(i,j)))+1   
1309      endif
1310      flot_marais(i,j) = marais(table_out_marais(i,j))
1311    enddo
1312  enddo
1313  !$OMP END DO
1314  !$OMP END PARALLEL
1315
1316  debug_3D(:,:,122)=real(table_out_marais(:,:))
1317
1318end subroutine determin_marais
1319
1320end module flottab_mod
1321
1322
Note: See TracBrowser for help on using the repository browser.