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

Last change on this file since 142 was 142, checked in by dumas, 7 years ago

Debug output suppressed and bug correction in calving

File size: 46.1 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    if (iter_beta.eq.0) then
499       Call dragging         
500    endif
501    if (itracebug.eq.1)  call tracebug(' routine flottab apres call dragging')
502!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
503!!$if (itestf.gt.0) then
504!!$   write(6,*) 'dans flottab apres dragging asymetrie sur H  pour time=',time
505!!$   stop
506!!$else
507!!$   write(6,*) 'dans flottab aapres dragging pas d asymetrie sur H  pour time=',time
508!!$
509!!$end if
510
511    !   6- calcule les vitesses des points qui sont devenus gzm
512    !$OMP PARALLEL
513    !$OMP DO         
514    do j=1,ny
515       do i=2,nx-1
516          !      si le point etait posé (non gz) et devient gzmx
517          !      definir la direction de la vitesse (moyenne des points)
518
519          if ((.not.gz1mx(i,j)).and.(.not.fl1mx(i,j)).and.gzmx(i,j).and. &
520               (i.gt.2).and.(i.lt.nx))  then
521             uxs1(i,j)=(uxbar(i+1,j)+uxbar(i-1,j))/2.
522          endif
523
524       end do
525    end do
526    !$OMP END DO
527
528    !$OMP DO
529    do j=2,ny-1
530       do i=1,nx
531          !      si le point etait posé (non gz) et devient gzmy
532          !      definir la direction de la vitesse (moyenne des points)
533          if ((.not.gz1my(i,j)).and.(.not.fl1my(i,j)).and.gzmy(i,j).and. &
534               (j.gt.2).and.(j.lt.ny))  then
535             uys1(i,j)=(uybar(i,j+1)+uybar(i,j-1))/2.
536          endif
537
538       end do
539    end do
540    !$OMP END DO
541
542
543    !  7-On determine finalement la position des noeuds stream ou shelf
544    !  -------------------------------------------------------------
545
546    if (nt.ge.2) then   ! pour ne pas faire ce calcul lors du premier passage
547       !$OMP WORKSHARE
548       uxbar(:,:)=uxs1(:,:)
549       uybar(:,:)=uys1(:,:)
550       !$OMP END WORKSHARE
551    endif
552
553    !$OMP WORKSHARE
554    flgzmx(:,:)=(marine.and.(flotmx(:,:).or.gzmx(:,:).or.ilemx(:,:)))   &
555         .or.(.not.marine.and.flotmx(:,:))
556    flgzmy(:,:)=(marine.and.(flotmy(:,:).or.gzmy(:,:).or.ilemy(:,:)))   &
557         .or.(.not.marine.and.flotmy(:,:))
558    !$OMP END WORKSHARE
559
560
561    ! 8- Pour la fusion basale sous les ice shelves- region proche de la grounding line
562    !---------------------------------------------------------------------------------
563    !       fbm est vrai si le point est flottant mais un des voisins est pose
564    !_________________________________________________________________________
565    !$OMP DO
566    do j=2,ny-1
567       do i=2,nx-1
568
569          ! if (i.gt.2.AND.i.lt.nx) then   
570          fbm(i,j)=flot(i,j).and.      & 
571               ((.not.flot(i+1,j)).or.(.not.flot(i,j+1)) &
572               .or.(.not.flot(i-1,j)).or.(.not.flot(i,j-1)))
573          !     endif
574       end do
575    end do
576    !$OMP END DO
577
578
579    !  9-On determine maintenant la position du front de glace
580    !  -------------------------------------------------------
581    !  C'est ici que l'on determine la position et le type de front
582    !       print*,'on est dans flottab pour definir les fronts'
583
584
585
586!!$do i=3,nx-2
587!!$   do j=3,ny-2
588!!$      if (h(i,j).gt.1.1) ice(i,j)=1       
589!!$   end do
590!!$end do
591!                       print*, 'flolottab debug', H(71,25),flot(71,25),bm(71,25),ice(71,25),time
592!print*,'H(90,179) flottab 2',H(90,179),ice(90,179),flot(90,179)
593    !$OMP WORKSHARE
594    where (flot(:,:))
595!cdc 1m       where (H(:,:).gt.max(Hmin,Hmin+BM(:,:)-Bmelt(:,:)))
596      where (H(:,:).gt.max(BM(:,:)-Bmelt(:,:)+petit_H,0.)*dt)
597!cdc                     where (H(:,:).gt.0.)
598          ice(:,:)=1
599       elsewhere 
600          ice(:,:)=0
601       end where
602    elsewhere
603       where (H(:,:).gt.0.) 
604          ice(:,:)=1
605       elsewhere 
606          ice(:,:)=0
607       end where
608    end where
609    !$OMP END WORKSHARE
610    !$OMP END PARALLEL
611!print*,'H(90,179) flottab 3',H(90,179),ice(90,179),flot(90,179),Bsoc(90,179)+H(90,179)*ro/row -sealevel
612    call DETERMIN_TACHE 
613    !~
614    !~ synchro :    if (isynchro.eq.1) then
615    !~ !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
616    !~ !!$if (itestf.gt.0) then
617    !~ !!$   write(6,*) 'dans flottab avant DETERMIN_TACHE  asymetrie sur H  pour time=',time
618    !~ !!$   stop
619    !~ !!$else
620    !~ !!$   write(6,*) 'dans flottab apres  DETERMIN_TACHE pas d asymetrie sur H  pour time=',time
621    !~ !!$
622    !~ !!$end if
623    !~
624    !~
625    !~ !----------------------------------------------!
626    !~ !On determine les differents ice strean/shelf  !
627    !~ !      call DETERMIN_TACHE                      !
628    !~ !----------------------------------------------!
629    !~
630    !~
631    !~ !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
632    !~ !!$if (itestf.gt.0) then
633    !~ !!$   write(6,*) 'dans flottab apres DETERMIN_TACHE  asymetrie sur H  pour time=',time
634    !~ !!$   stop
635    !~ !!$else
636    !~ !!$   write(6,*) 'dans flottab apres DETERMIN_TACHE pas d asymetrie sur H  pour time=',time
637    !~ !!$
638    !~ !!$end if
639   
640    !~ !ice(:,:)=0  ! Attention, voir si ca marche toujours pour l'Antarctique et heminord !
641   
642    !On compte comme englacé uniquement les calottes dont une partie est posée     
643    !$OMP PARALLEL PRIVATE(smax_,smax_coord,smax_i,smax_j,mask_tache_ij)
644    !$OMP DO
645    do j=3,ny-2
646       do i=3,nx-2
647    test1:  if (.not.iceberg1D(table_out(i,j))) then ! on est pas sur un iceberg
648             if (nb_pts_tache(table_out(i,j)).ge.1) then
649                 ice(i,j)=1
650    !~             if (nb_pts_tache(table_out(i,j)).le.10) then  ! les petites iles sont en sia
651    !~ !    write(6,*) 'petite ile ',i,j
652    !~               flgzmx(i,j)=.false.
653    !~                flgzmx(i+1,j)=.false.
654    !~                flgzmy(i,j)=.false.
655    !~                flgzmy(i,j+1)=.false.
656    !~                gzmx(i:i+1,j)=.false.
657    !~                gzmy(i,j:j+1)=.false.
658    !~             endif
659    !~
660    !~
661    ! ici on est sur une tache non iceberg >= 5 points                       
662    ! on teste si la tache n'est pas completement ice stream
663   
664    test2:       if (icetrim(table_out(i,j))) then   ! on a une ile d'ice stream
665   
666       mask_tache_ij(:,:)=.false.
667       mask_tache_ij(:,:)=(table_out(:,:).eq.table_out(i,j))       ! pour toute la tache
668   
669         smax_=maxval(S(:,:),MASK=mask_tache_ij(:,:))
670         smax_coord(:)=maxloc(S(:,:),MASK=mask_tache_ij(:,:))
671         smax_i=smax_coord(1)
672         smax_j=smax_coord(2)
673   
674    !~ !!$               smax_i=0 ; smax_j=0 ; smax_=sealevel
675    !~ !!$               do ii=3,nx-2
676    !~ !!$                  do jj=3,ny-2
677    !~ !!$                     if (table_out(ii,jj).eq.table_out(i,j)) then
678    !~ !!$                        if (s(ii,jj).gt.smax_) then                           
679    !~ !!$                           smax_ =s(ii,jj)
680    !~ !!$                           smax_i=ii
681    !~ !!$                           smax_j=jj
682    !~ !!$                        endif
683    !~ !!$                     endif
684    !~ !!$                  end do
685    !~ !!$               end do
686   
687   
688                   gzmx(smax_i,smax_j)=.false. ; gzmx(smax_i+1,smax_j)=.false.
689                   gzmy(smax_i,smax_j)=.false. ; gzmx(smax_i,smax_j+1)=.false.
690                   flgzmx(smax_i,smax_j)=.false. ; flgzmx(smax_i+1,smax_j)=.false.
691                   flgzmy(smax_i,smax_j)=.false. ; flgzmx(smax_i,smax_j+1)=.false.
692     
693                   if (Smax_.le.sealevel) then
694                      write(num_tracebug,*)'Attention, une ile avec la surface sous l eau'
695                      write(num_tracebug,*)'time=',time,'   coord:',smax_i,smax_j
696                   end if
697   
698                endif test2 
699             end if                                             ! endif deplace
700!cdc transfere dans calving :             
701           else ! on est sur un iceberg                          !   test1
702                iceberg(i,j)=iceberg1D(table_out(i,j))
703!~              ice(i,j)=0
704!~              h(i,j)=0. !1. afq, we should put everything in calving!
705!~              surnet=H(i,j)*(1.-ro/row)
706!~              S(i,j)=surnet+sealevel
707!~              B(i,j)=S(i,j)-H(i,j)
708             
709          endif test1
710   
711   
712       end do
713    end do
714   
715    !$OMP END DO
716    !$OMP END PARALLEL
717
718    !~
719    !~ !----------------------------------------------
720    !~ ! On caracterise le front des ice shelfs/streams
721    !~
722    !~ !      call DETERMIN_FRONT   
723    !~
724    !~ !----------------------------------------------     
725    !~ !!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
726    !~ !!$if (itestf.gt.0) then
727    !~ !!$   write(6,*) 'dans flottab apres DETERMIN_front  asymetrie sur H  pour time=',time
728    !~ !!$   stop
729    !~ !!$else
730    !~ !!$   write(6,*) 'dans flottab apres DETERMIN_front pas d asymetrie sur H  pour time=',time
731    !~ !!$
732    !~ !!$end if
733    !~
734    !~ endif synchro
735
736    ! correction momentanée pour symetrie Heino
737    !where ((.not.flot(:,:)).and.(ice(:,:).eq.0)) H(:,:)=0.     
738
739    !fin de routine flottab2
740    !print*, 'front',front(50,30),ice(50,30),flotmx(i,j),uxbar(i,j)
741    !~ print*,'fin flottab',S(132,183),H(132,183),BSOC(132,183),B(132,183),sealevel
742    !~ print*,'fin flottab',flot(132,183),ice(132,183),gzmx(132,183),gzmy(132,183),ilemx(132,183),ilemy(132,183)
743    !read(5,*)
744
745    !call determin_front ! cette version ne conserve pas la masse !!!
746    call determin_front_tof ! version simplifiee
747
748    call determin_marais
749   
750  end subroutine flottab
751!-------------------------------------------------------------------- 
752
753!> SUBROUTINE: determin_tache
754!! Routine pour la dtermination du numero de tache a effectuer
755!>
756subroutine determin_tache
757
758!$ USE OMP_LIB
759
760implicit none
761integer :: indice
762integer :: label         ! no des taches rencontrées dans le mask
763integer :: label_max     ! no temporaire maxi de tache rencontrées
764!      integer :: mask_nb = 4
765integer,parameter :: mask_nb = 2   ! version ou on ne compte pas les diagonales
766integer :: vartemp       ! variable temporaire pour reorganiser compt
767!     integer,dimension(mask_nb) :: mask
768integer,dimension(mask_nb) :: mask
769
770
771! 1-initialisation
772!-----------------
773label_max=1      ! numero de la tache, la premiere tache est notée 1
774label=1
775do i=1,n_ta_max
776   compt(i)=i
777enddo
778!      table_in  = .false.
779!$OMP PARALLEL
780!$OMP WORKSHARE
781table_out(:,:) = 0
782iceberg1D(:)  = .true.
783icetrim (:) = .true.
784nb_pts_tache(:) = 0
785!$OMP END WORKSHARE
786!$OMP END PARALLEL
787!    open(unit=100,file="tache.data",status='replace')
788
789! 2-reperage des taches
790!----------------------
791
792do j=2,ny-1
793  do i=2,nx-1
794
795
796
797     IF (ice(i,j).ge.1) THEN ! on est sur la glace-----------------------------!
798                       
799        if ((ice(i-1,j).ge.1).or.(ice(i,j-1).ge.1)) then   !masque de 2 cases adjacentes
800           !      un des voisins est deja en glace
801           mask(1) = table_out(i-1,j)
802           mask(2) = table_out(i,j-1) 
803           label = label_max
804
805           !on determine la valeur de la tache minimun (>0) presente ds le masque
806           do indice=1,mask_nb
807              if (mask(indice).gt.0) label=min(label,mask(indice))
808           enddo
809!cdc       label=min(label,minval(mask(:), mask=mask > 0))
810
811           !on fixe la valeur de la tache voisine minimun au point etudie (via label)
812           table_out(i,j)=label
813           !si ce noeud est posé, alors la tache n'est pas un iceberg et iceberg=.F.
814           if (.not.FLOT(I,J)) then
815              iceberg1D(label)=.false.
816           endif
817
818           !si ce noeud est posé, alors la tache n'est pas un ice stream et icestrim=.F.
819           if ((.not.gzmx(i,j).and..not.gzmx(i+1,j)).and.  &
820               (.not.gzmy(i,j).and..not.gzmy(i,j+1))) then
821              icetrim(label)=.false.
822           endif
823
824           ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt
825           ! i.e. les plus grands numeros correspondent au plus petit
826           ! on lui affecte le numero de la tache fondamentale avec un signe -
827           ! pour indiquer le changement
828
829           do indice=1,mask_nb
830              if(mask(indice).gt.label) then
831                 compt(mask(indice))=-label
832              endif
833           enddo
834
835        else !aucun des voisins est une tache
836           table_out(i,j)= label_max
837           compt(label_max)=label_max
838           !si ce noeud est posé, alors la ache n'est pas un iceberg et iceberg=.F.
839           if (.not.FLOT(I,J)) then
840              iceberg1D(label_max)=.false.
841           endif
842
843           !si ce noeud est posé, alors le tache n'est pas un ice stream et icestrim=.F.
844           if ((.not.gzmx(i,j).and..not.gzmx(i+1,j)).and.  &
845               (.not.gzmy(i,j).and..not.gzmy(i,j+1))) then
846              icetrim(label)=.false.
847           endif
848             
849           label_max  = label_max+1
850           if (label_max.gt.n_ta_max) print*,'trop de taches=',label_max 
851        endif
852
853
854     else           !on est pas sur une tache----------------------------------------------
855        table_out(i,j)=0    ! Pas necessaire (reecrit 0 sur 0)             
856     endif          !---------------------------------------------------------------------
857
858
859  enddo
860enddo
861
862
863
864! On reorganise compt en ecrivant le numero de la tache fondamentale
865! i.e. du plus petit numero present sur la tache (Sans utiliser de recursivité)
866! On indique aussi le nb de point que contient chaque taches (nb_pts_tache)
867
868do indice=1,label_max
869   vartemp = compt(indice)
870   if (compt(indice).lt.0) then
871      compt(indice)= compt(-vartemp)
872      if (.not.iceberg1D(indice)) iceberg1D(-vartemp)=.false.
873      if (.not.icetrim(indice)) icetrim(-vartemp)=.false.
874   endif   
875enddo
876
877!$OMP PARALLEL
878!$OMP DO REDUCTION(+:nb_pts_tache)
879do j=1,ny
880   do i=1,nx
881      if (table_out(i,j).ne.0) then
882         table_out(i,j)=compt(table_out(i,j)) 
883         nb_pts_tache(compt(table_out(i,j)))= nb_pts_tache(compt(table_out(i,j)))+1   
884      endif
885   enddo
886enddo
887!$OMP END DO
888!$OMP END PARALLEL
889
890
891!!$tablebis(:,:)=table_out(:,:)
892!!$do j=1,ny
893!!$   do i=1,nx
894!!$      if (tablebis(i,j).ne.0) then     ! tache de glace
895!!$         numtache=table_out(i,j)
896!!$         nb_pt=count(table_out(:,:).eq.numtache)  ! compte tous les points de la tache
897!!$         nb_pts_tache(table_out(i,j))=nb_pt       !
898!!$
899!!$         where (table_out(:,:).eq.numtache)         
900!!$            tablebis(:,:)=0                       ! la table de tache est remise a 0 pour eviter de repasser
901!!$         end where
902!!$          write(6,*) i,j,nb_pt,table_out(i,j)
903!!$      endif
904!!$   enddo
905!!$enddo
906
907!!$do j=1,ny
908!!$   do i=1,nx
909!!$      debug_3d(i,j,56)=nb_pts_tache(table_out(i,j))
910!!$   end do
911!!$end do
912!!$     
913end subroutine determin_tache
914!----------------------------------------------------------------------
915!> SUBROUTINE: determin_front
916!!Routine pour la determination du front
917!>
918subroutine determin_front
919!!$ USE OMP_LIB
920integer :: i_moins1,i_plus1,i_plus2
921integer :: j_moins1,j_plus1,j_plus2
922     
923      !$OMP PARALLEL
924      !$OMP DO
925      do j=3,ny-2
926        do i=3,nx-2
927
928 surice:if  (ice(i,j).eq.0) then
929         do pmx=-1,1,2
930          do pmy=-1,1,2
931
932 diagice :  if (ice(i+pmx,j+pmy).eq.1) then
933
934              if ((ice(i+pmx,j)+ice(i,j+pmy).eq.2)) then    ! test (i) pour eviter les langues
935                                                            ! de glaces diagonales en coin(26dec04)
936                   if ((ice(i+2*pmx,j).eq.1.and.ice(i+2*pmx,j+pmy).eq.0).or.&
937                       (ice(i,j+2*pmy).eq.1.and.ice(i+pmx,j+2*pmy).eq.0)) then
938                           ice(i,j)=1
939                           h(i,j)=max(1.,h(i,j))
940                   endif         
941
942! test (i) pour eviter les langues de glaces diagonales :
943!                        mouvement du cheval aux echecs 
944
945               if ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.1)) then
946                   if (ice(i+2*pmx,j+pmy).eq.1.and. &
947                      (ice(i+2*pmx,j+2*pmy)+ice(i,j+2*pmy)).ge.1) then
948                           ice(i,j)=1
949                           h(i,j)=max(1.,h(i,j))
950                   endif       
951                   if (ice(i+pmx,j+2*pmy).eq.1.and. &
952                      (ice(i+2*pmx,j+2*pmy)+ice(i+2*pmx,j)).ge.1) then
953                           ice(i,j)=1
954                           h(i,j)=max(1.,h(i,j))
955                   endif       
956
957! test (ii) pour eviter les langues de glaces diagonales :
958!            le point glace ice(i+pmx,j+pmy) a :
959!                          - ses 4 voisins frontaux en glace
960!                          - mais 2 voisins vides diagonalement opposes   
961           
962               elseif ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.2)  &
963                                .and.ice(i+2*pmx,j+2*pmy).eq.0) then
964
965! test (iii) pour faire les tests (i) et (ii)
966                    ice(i,j)=1
967                    h(i,j)=max(1.,h(i,j))
968                    ice(i+2*pmx,j+2*pmy)=1
969                    h(i+2*pmx,j+2*pmy)=max(1.,h(i+2*pmx,j+2*pmy))
970               endif
971            endif
972           endif diagice   
973          enddo
974         enddo
975         endif surice       
976       end do
977      end do
978      !$OMP END DO
979      !$OMP ENd PARALLEL
980
981!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
982!!$if (itestf.gt.0) then
983!!$   write(6,*) 'dans front avant remplissage baies asymetrie sur H  pour time=',time
984!!$   stop
985!!$else
986!!$   write(6,*) 'dans front avant remplissage baies pas d asymetrie sur H  pour time=',time
987!!$
988!!$end if
989
990
991!     print*,'dans remplissage baies',time
992       
993baies: do k=1,2
994!$OMP PARALLEL
995!$OMP DO PRIVATE(i_moins1,j_moins1,i_plus1,j_plus1,i_plus2,j_plus2)
996do j=1,ny
997   do i=1,nx
998
999surice_xy: if  (ice(i,j).eq.0) then
1000   i_moins1=max(i-1,1)
1001   j_moins1=max(j-1,1)
1002   i_plus1=min(i+1,nx)
1003   j_plus1=min(j+1,ny)
1004   i_plus2=min(i+2,nx)
1005   j_plus2=min(j+2,ny)
1006 
1007! test (iii) pour trouver les baies de largeur 1 ou 2 cases
1008! et combler les trous si ce sont des baies
1009! si ce ne sont pas des baies, ne pas combler et creer des langues de glaces artificielles             
1010! baies horizontales 
1011 
1012         if (ice(i_moins1,j).eq.1.and.(ice(i_plus1,j).eq.1.or.ice(i_plus2,j).eq.1)) then
1013            if (ice(i,j_moins1).eq.1.or.ice(i,j_plus1).eq.1)  then    ! ice(i,j)=1
1014               ice(i,j)=1
1015               H(i,j)=max(1.,H(i,j))
1016            endif
1017         endif
1018
1019
1020         if (ice(i,j_moins1).eq.1.and.(ice(i,j_plus1).eq.1.or.ice(i,j_plus2).eq.1)) then
1021            if (ice(i_moins1,j).eq.1.or.ice(i_plus1,j).eq.1) then !ice(i,j)=1
1022               ice(i,j)=1
1023               H(i,j)=max(1.,H(i,j))
1024            endif
1025         endif
1026
1027      endif surice_xy
1028   end do
1029end do
1030!$OMP END DO
1031!$OMP END PARALLEL
1032end do baies
1033
1034!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
1035!!$if (itestf.gt.0) then
1036!!$   write(6,*) 'dans front apres remplissage baies asymetrie sur H  pour time=',time
1037!!$   stop
1038!!$else
1039!!$   write(6,*) 'dans front apres remplissage baies pas d asymetrie sur H  pour time=',time
1040!!$
1041!!$end if
1042
1043!$OMP PARALLEL
1044!$OMP DO
1045do j=2,ny-1
1046   do i=2,nx-1
1047
1048      if (ice(i,j).eq.1) then     !         test si ice=1
1049
1050! if ice, on determine front...
1051! ainsi, front=0 sur les zones = 0
1052
1053         front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1054         !front= le nb de faces en contact avec un voisin englacé
1055      endif
1056   end do
1057end do
1058!$OMP END DO
1059
1060! traitement des bords. On considere que l'exterieur n'a pas de glace
1061! attention ce n'est vrai que sur la grande grille
1062
1063!$OMP DO PRIVATE(i)
1064do j=2,ny-1
1065   i=1
1066   front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1067   i=nx
1068   front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1))
1069end do
1070!$OMP END DO
1071
1072!$OMP DO PRIVATE(j)
1073do i=2,nx-1
1074   j=1 
1075   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1))
1076   j=ny
1077   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1))
1078end do
1079!$OMP END DO
1080
1081! traitement des coins
1082
1083front(1,1)=ice(2,1)+ice(2,1)
1084front(1,ny)=ice(2,ny)+ice(1,ny-1)
1085front(nx,1)=ice(nx,2)+ice(nx-1,1)
1086front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny)
1087
1088!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
1089!!$if (itestf.gt.0) then
1090!!$   write(6,*) 'dans front apres front asymetrie sur H  pour time=',time
1091!!$   stop
1092!!$else
1093!!$   write(6,*) 'dans front apres front pas d asymetrie sur H  pour time=',time
1094!!$
1095!!$end if
1096
1097!   on ne compte pas les taches de glace de 2 cases (horizontales ou verticales)
1098!   en fait, si ces deux cases sont flottantes, il faut enlever les icebergs
1099!   de n'importe quelle taille).
1100!   si ces deux taches sont posées (ou une des deux), il n'y a pas assez de conditions aux limites
1101
1102!$OMP DO
1103do j=1,ny
1104   do i=1,nx-1
1105      if (front(i,j).eq.1) then
1106         if (front(i+1,j).eq.1) then
1107            ice(i,j)=0
1108            ice(i+1,j)=0
1109            front(i,j)=0
1110            front(i+1,j)=0
1111         endif
1112      endif
1113   end do
1114end do
1115!$OMP END DO
1116
1117!$OMP DO
1118do j=1,ny-1
1119   do i=1,nx
1120      if (front(i,j).eq.1) then
1121         if (front(i,j+1).eq.1) then
1122            ice(i,j)=0
1123            ice(i,j+1)=0
1124            front(i,j)=0
1125            front(i,j+1)=0
1126         endif
1127      end if
1128   end do
1129end do
1130!$OMP END DO
1131
1132!isolx signifie pas de voisins en x
1133!isoly signifie pas de voisins en y
1134!remarque :
1135!si isolx/y=.true. alors frontfacex/y=0 (a la fois +1 & -1 or +1-1=0)
1136
1137! calcul de frontfacex et isolx
1138!$OMP DO
1139do j=1,ny
1140   do i=2,nx-1
1141
1142      if (front(i,j).ge.1.and.front(i,j).le.3) then    !front(entre 1 et 3)
1143         
1144         if ((ice(i-1,j)+ice(i+1,j)).lt.2) then        ! il y a un front // a x
1145         
1146            if ((ice(i-1,j)+ice(i+1,j)).eq.0) then
1147               isolx(i,j)=.true.
1148            elseif (ice(i-1,j).eq.0) then
1149               frontfacex(i,j)=-1                      ! front  i-1 |i  i+1
1150            else
1151               frontfacex(i,j)=+1                      ! front  i-1  i| i+1
1152            endif
1153         endif
1154      end if !fin du test il y a un front
1155         
1156   end do
1157end do
1158!$OMP END DO
1159
1160! calcul de frontfacey et isoly
1161!$OMP DO
1162do j=2,ny-1
1163   do i=1,nx
1164
1165      if (front(i,j).ge.1.and.front(i,j).le.3) then   !front(entre 1 et 3)
1166         
1167         if ((ice(i,j-1)+ice(i,j+1)).lt.2) then       ! il y a un front // a y
1168             
1169            if ((ice(i,j-1)+ice(i,j+1)).eq.0) then
1170               isoly(i,j)=.true.                      !front   j-1 |j| j+1
1171            elseif (ice(i,j-1).eq.0) then
1172               frontfacey(i,j)=-1                     !front   j-1 |j j+1
1173            else
1174               frontfacey(i,j)=+1                     !front   j-1  j| j+1
1175            endif
1176         endif
1177     end if !fin du test il y a un front
1178         
1179   end do
1180end do
1181!$OMP END DO
1182
1183
1184! traitement des bords. On considere que l'exterieur n'a pas de glace
1185! attention ce n'est vrai que sur la grande grille
1186
1187!$OMP DO PRIVATE(i)
1188do j=2,ny-1
1189   i=1
1190   if (front(i,j).ge.1)  then
1191      if (ice(i+1,j).eq.0) then
1192         isolx(i,j)=.true.
1193      else
1194         frontfacex(i,j)=-1
1195      endif
1196   end if
1197   i=nx
1198   if (front(i,j).ge.1)  then
1199      if (ice(i-1,j).eq.0) then
1200         isolx(i,j)=.true.
1201      else
1202         frontfacex(i,j)=1
1203      endif
1204   end if
1205end do
1206!$OMP END DO
1207
1208!$OMP DO PRIVATE(j)
1209do i=2,nx-1
1210   j=1 
1211   if (front(i,j).ge.1)  then
1212      if (ice(i,j+1).eq.0) then
1213         isoly(i,j)=.true.
1214      else
1215         frontfacey(i,j)=-1
1216      endif
1217   end if
1218   j=ny
1219   if (front(i,j).ge.1)  then
1220      if (ice(i,j-1).eq.0) then
1221         isoly(i,j)=.true.
1222      else
1223         frontfacey(i,j)=1
1224      endif
1225   end if
1226end do
1227!$OMP END DO
1228!$OMP END PARALLEL
1229
1230return
1231end subroutine determin_front
1232!------------------------------------------------------------------------------
1233
1234subroutine determin_front_tof
1235
1236integer,dimension(nx,ny) :: nofront ! tableau de travail (points au dela du front)
1237!$OMP PARALLEL
1238!$OMP DO
1239do j=2,ny-1
1240   do i=2,nx-1
1241
1242      if (ice(i,j).eq.1) then     !         test si ice=1
1243
1244! if ice, on determine front...
1245! ainsi, front=0 sur les zones = 0
1246
1247         front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1248         !front= le nb de faces en contact avec un voisin englacé
1249      endif
1250   end do
1251end do
1252!$OMP END DO
1253
1254! traitement des bords. On considere que l'exterieur n'a pas de glace
1255! attention ce n'est vrai que sur la grande grille
1256
1257!$OMP DO PRIVATE(i)
1258do j=2,ny-1
1259   i=1
1260   front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1261   i=nx
1262   front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1))
1263end do
1264!$OMP END DO
1265
1266!$OMP DO PRIVATE(j)
1267do i=2,nx-1
1268   j=1 
1269   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1))
1270   j=ny
1271   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1))
1272end do
1273!$OMP END DO
1274!$OMP BARRIER
1275! traitement des coins
1276
1277front(1,1)=ice(2,1)+ice(2,1)
1278front(1,ny)=ice(2,ny)+ice(1,ny-1)
1279front(nx,1)=ice(nx,2)+ice(nx-1,1)
1280front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny)
1281
1282! Les points à plus d'un point de grille du bord sont front=-1
1283!$OMP WORKSHARE
1284nofront(:,:)=0
1285!$OMP END WORKSHARE
1286!$OMP BARRIER
1287
1288!$OMP DO
1289do j=2,ny-1
1290        do i=2,nx-1
1291                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
1292                        nofront(i,j)=-1
1293                endif
1294        enddo
1295enddo
1296!$OMP END DO
1297!$OMP BARRIER
1298
1299!$OMP WORKSHARE
1300where (nofront(:,:).eq.-1)
1301        front(:,:)=-1
1302endwhere
1303!$OMP END WORKSHARE
1304!$OMP END PARALLEL             
1305
1306end subroutine determin_front_tof
1307
1308
1309!> SUBROUTINE: determin_marais
1310!! afq -- Routine pour l'identification des "marais"
1311!! Un marais est une tache "shelf" entouré de points grounded
1312!! Copie sauvage de determin_tache, adapte au probleme du marais
1313!>
1314subroutine determin_marais
1315
1316!$ USE OMP_LIB
1317
1318implicit none
1319integer :: indice
1320integer :: label         ! no des taches rencontrées dans le mask
1321integer :: label_max     ! no temporaire maxi de tache rencontrées
1322!      integer :: mask_nb = 4
1323integer,parameter :: mask_nb = 2   ! version ou on ne compte pas les diagonales
1324integer :: vartemp       ! variable temporaire pour reorganiser compt
1325!     integer,dimension(mask_nb) :: mask
1326integer,dimension(mask_nb) :: mask
1327integer pm ! loop integer
1328
1329integer,dimension(nx,ny)        :: table_out_marais      !< pour les numeros des taches
1330integer,dimension(0:n_ta_max)     :: compt_marais        !< contient les equivalence entre les taches
1331integer,dimension(0:n_ta_max)     :: nb_pts_marais !< indique le nombre de points par tache
1332logical,dimension(0:n_ta_max)     :: marais      !< T si iceberg, F si calotte posee
1333 
1334! 1-initialisation
1335!-----------------
1336label_max=1      ! numero de la tache, la premiere tache est notée 1
1337label=1
1338do i=1,n_ta_max
1339   compt_marais(i)=i
1340enddo
1341!$OMP PARALLEL
1342!$OMP WORKSHARE
1343table_out_marais(:,:) = 0
1344marais(:)  = .true.
1345nb_pts_marais(:) = 0
1346!$OMP END WORKSHARE
1347!$OMP END PARALLEL
1348
1349! 2-reperage des taches
1350!----------------------
1351
1352do j=2,ny-1
1353  do i=2,nx-1
1354
1355
1356
1357     IF ((ice(i,j).ge.1).and.flot(i,j)) THEN ! on est sur la glace qui flotte-------------------!
1358     !IF ((ice(i,j).ge.1).and.flot(i,j).and.(H(i,j).gt.1.)) THEN ! on est sur la glace qui flotte-------------------!
1359               
1360        !write (*,*) "un point qui nous interesse!",i,j
1361       
1362        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
1363        !if (((ice(i-1,j).ge.1).and.flot(i-1,j).and.(H(i-1,j).gt.1.)).or.((ice(i,j-1).ge.1).and.flot(i,j-1).and.(H(i,j-1).gt.1.))) then   !masque de 2 cases adjacentes
1364           !      un des voisins est deja en glace
1365           mask(1) = table_out_marais(i-1,j)
1366           mask(2) = table_out_marais(i,j-1) 
1367           label = label_max
1368
1369           !on determine la valeur de la tache minimun (>0) presente ds le masque
1370           do indice=1,mask_nb
1371              if (mask(indice).gt.0) label=min(label,mask(indice))
1372           enddo
1373
1374           !on fixe la valeur de la tache voisine minimun au point etudie (via label)
1375           table_out_marais(i,j)=label
1376           
1377           !si un des voisins n'est pas glace alors la tache n'est pas un marais
1378           do pm=-1,1,2
1379              if ( (ice(i+pm,j).eq.0) .or. (ice(i,j+pm).eq.0) ) then
1380                 marais(label)=.false.
1381              endif
1382           enddo
1383
1384           ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt
1385           ! i.e. les plus grands numeros correspondent au plus petit
1386           ! on lui affecte le numero de la tache fondamentale avec un signe -
1387           ! pour indiquer le changement
1388
1389           do indice=1,mask_nb
1390              if(mask(indice).gt.label) then
1391                 compt_marais(mask(indice))=-label
1392              endif
1393           enddo
1394           
1395           
1396        else !aucun des voisins est une tache
1397           table_out_marais(i,j)= label_max
1398           compt_marais(label_max)=label_max
1399           
1400           !si un des voisins n'est pas glace alors la tache n'est pas un marais
1401           do pm=-1,1,2
1402              if ( (ice(i+pm,j).eq.0) .or. (ice(i,j+pm).eq.0) ) then
1403                 marais(label_max)=.false.
1404              endif
1405           enddo
1406
1407           
1408           label_max  = label_max+1
1409           if (label_max.gt.n_ta_max) print*,'trop de taches=',label_max 
1410        endif
1411
1412
1413     else           !on est pas sur une tache----------------------------------------------
1414        table_out_marais(i,j)=0    ! Pas necessaire (reecrit 0 sur 0)             
1415     endif          !---------------------------------------------------------------------
1416
1417
1418  enddo
1419enddo
1420
1421
1422
1423! On reorganise compt en ecrivant le numero de la tache fondamentale
1424! i.e. du plus petit numero present sur la tache (Sans utiliser de recursivité)
1425! On indique aussi le nb de point que contient chaque taches (nb_pts_tache)
1426
1427do indice=1,label_max
1428   vartemp = compt_marais(indice)
1429   if (compt_marais(indice).lt.0) then
1430      compt_marais(indice)= compt_marais(-vartemp)
1431      if (.not.marais(indice)) marais(-vartemp)=.false.
1432   endif   
1433enddo
1434
1435!$OMP PARALLEL
1436!$OMP DO REDUCTION(+:nb_pts_tache)
1437do j=1,ny
1438   do i=1,nx
1439      if (table_out_marais(i,j).ne.0) then
1440         table_out_marais(i,j)=compt_marais(table_out_marais(i,j)) 
1441         nb_pts_marais(compt_marais(table_out_marais(i,j)))= nb_pts_marais(compt_marais(table_out_marais(i,j)))+1   
1442      endif
1443   enddo
1444enddo
1445!$OMP END DO
1446!$OMP END PARALLEL
1447
1448do j=1,ny
1449   do i=1,nx
1450      flot_marais(i,j) = marais(table_out_marais(i,j))
1451   enddo
1452enddo
1453
1454end subroutine determin_marais
1455
1456end module flottab_mod
1457
1458
Note: See TracBrowser for help on using the repository browser.