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

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

iter_beta flag is obsolete now

File size: 42.5 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
642integer :: vartemp       ! variable temporaire pour reorganiser compt
643!     integer,dimension(mask_nb) :: mask
644integer,dimension(mask_nb) :: mask
645
646
647! 1-initialisation
648!-----------------
649label_max=1      ! numero de la tache, la premiere tache est notée 1
650label=1
651do i=1,n_ta_max
652   compt(i)=i
653enddo
654!      table_in  = .false.
655!$OMP PARALLEL
656!$OMP WORKSHARE
657table_out(:,:) = 0
658iceberg1D(:)  = .true.
659icetrim (:) = .true.
660nb_pts_tache(:) = 0
661!$OMP END WORKSHARE
662!$OMP END PARALLEL
663!    open(unit=100,file="tache.data",status='replace')
664
665! 2-reperage des taches
666!----------------------
667
668do j=2,ny-1
669  do i=2,nx-1
670
671     IF (ice(i,j).ge.1) THEN ! on est sur la glace-----------------------------!
672                       
673        if ((ice(i-1,j).ge.1).or.(ice(i,j-1).ge.1)) then   !masque de 2 cases adjacentes
674           !      un des voisins est deja en glace
675           mask(1) = table_out(i-1,j)
676           mask(2) = table_out(i,j-1) 
677           label = label_max
678
679           !on determine la valeur de la tache minimun (>0) presente ds le masque
680           do indice=1,mask_nb
681              if (mask(indice).gt.0) label=min(label,mask(indice))
682           enddo
683!cdc       label=min(label,minval(mask(:), mask=mask > 0))
684
685           !on fixe la valeur de la tache voisine minimun au point etudie (via label)
686           table_out(i,j)=label
687           !si ce noeud est posé, alors la tache n'est pas un iceberg et iceberg=.F.
688           if (.not.FLOT(I,J)) then
689              iceberg1D(label)=.false.
690           endif
691
692           !si ce noeud est posé, alors la tache n'est pas un ice stream et icestrim=.F.
693           if ((.not.gzmx(i,j).and..not.gzmx(i+1,j)).and.  &
694               (.not.gzmy(i,j).and..not.gzmy(i,j+1))) then
695              icetrim(label)=.false.
696           endif
697
698           ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt
699           ! on lui affecte le numero de la tache fondamentale
700
701           do indice=1,mask_nb
702              if(mask(indice).gt.label) then
703                 compt(mask(indice))=-label
704              endif
705           enddo
706                                                                   ! exemple on est sur le point X : 5  X
707          do indice=1,mask_nb                                      !                                   20
708            if(mask(indice).gt.label) then                         ! mask(2)=20 > 5             
709              compt(mask(indice))=label                            ! compt(20)=5
710              if (.not.iceberg1D(mask(indice))) iceberg1D(label)=.false. ! si la tache n'etais pas un iceberg => iceberg =.false.  iceberg(5)=.false.               
711              if (.not.icetrim(mask(indice))) icetrim(label)=.false.
712              where (table_out(:,:).eq.mask(indice))               ! where table_out(:,:)=mask(2)=20
713                table_out(:,:)=label                        ! table_out(:,:)=label=5                 
714              endwhere
715            endif
716          enddo   
717
718        else !aucun des voisins est une tache
719           table_out(i,j)= label_max
720           compt(label_max)=label_max
721           !si ce noeud est posé, alors la ache n'est pas un iceberg et iceberg=.F.
722           if (.not.FLOT(I,J)) then
723              iceberg1D(label_max)=.false.
724           endif
725
726           !si ce noeud est posé, alors le tache n'est pas un ice stream et icestrim=.F.
727           if ((.not.gzmx(i,j).and..not.gzmx(i+1,j)).and.  &
728               (.not.gzmy(i,j).and..not.gzmy(i,j+1))) then
729              icetrim(label)=.false.
730           endif
731             
732           label_max  = label_max+1
733           if (label_max.gt.n_ta_max) print*,'ATTENTION trop de taches icebergs=',label_max 
734        endif
735
736
737     else           !on est pas sur une tache----------------------------------------------
738        table_out(i,j)=0    ! Pas necessaire (reecrit 0 sur 0)             
739     endif          !---------------------------------------------------------------------
740
741
742  enddo
743enddo
744
745
746
747! On reorganise compt en ecrivant le numero de la tache fondamentale
748! i.e. du plus petit numero present sur la tache (Sans utiliser de recursivité)
749! On indique aussi le nb de point que contient chaque taches (nb_pts_tache)
750
751  !$OMP PARALLEL
752  !$OMP DO
753  do j=1,ny
754    do i=1,nx
755      if (table_out(i,j).ne.0) then
756        nb_pts_tache(compt(table_out(i,j)))= nb_pts_tache(compt(table_out(i,j)))+1
757      endif
758    enddo
759  enddo
760  !$OMP END DO
761  !$OMP END PARALLEL
762
763  !On compte comme englacé uniquement les calottes dont une partie est posée     
764  !$OMP PARALLEL PRIVATE(smax_,smax_coord,smax_i,smax_j,mask_tache_ij)
765  !$OMP DO
766  do j=3,ny-2
767    do i=3,nx-2
768      test1: if (.not.iceberg1D(table_out(i,j))) then ! on est pas sur un iceberg
769        if (nb_pts_tache(table_out(i,j)).ge.1) then
770          ice(i,j)=1
771          ! ici on est sur une tache non iceberg >= 5 points                       
772          ! on teste si la tache n'est pas completement ice stream
773   
774          test2: if (icetrim(table_out(i,j))) then   ! on a une ile d'ice stream
775   
776            mask_tache_ij(:,:)=.false.
777            mask_tache_ij(:,:)=(table_out(:,:).eq.table_out(i,j))       ! pour toute la tache
778   
779            smax_=maxval(S(:,:),MASK=mask_tache_ij(:,:))
780            smax_coord(:)=maxloc(S(:,:),MASK=mask_tache_ij(:,:))
781            smax_i=smax_coord(1)
782            smax_j=smax_coord(2)
783 
784            gzmx(smax_i,smax_j)=.false. ; gzmx(smax_i+1,smax_j)=.false.
785            gzmy(smax_i,smax_j)=.false. ; gzmx(smax_i,smax_j+1)=.false.
786            flgzmx(smax_i,smax_j)=.false. ; flgzmx(smax_i+1,smax_j)=.false.
787            flgzmy(smax_i,smax_j)=.false. ; flgzmx(smax_i,smax_j+1)=.false.
788     
789            if (Smax_.le.sealevel) then
790              write(num_tracebug,*)'Attention, une ile avec la surface sous l eau'
791              write(num_tracebug,*)'time=',time,'   coord:',smax_i,smax_j
792            end if   
793          endif test2 
794        end if                                             ! endif deplace
795!cdc transfere dans calving :             
796      else ! on est sur un iceberg                          !   test1
797        iceberg(i,j)=iceberg1D(table_out(i,j))
798!~        ice(i,j)=0
799!~        h(i,j)=0. !1. afq, we should put everything in calving!
800!~        surnet=H(i,j)*(1.-ro/row)
801!~        S(i,j)=surnet+sealevel
802!~        B(i,j)=S(i,j)-H(i,j)             
803      endif test1   
804    end do
805  end do   
806  !$OMP END DO
807  !$OMP END PARALLEL
808
809debug_3D(:,:,121)=real(table_out(:,:))
810
811end subroutine determin_tache
812!----------------------------------------------------------------------
813!> SUBROUTINE: determin_front
814!!Routine pour la determination du front
815!>
816subroutine determin_front
817!!$ USE OMP_LIB
818integer :: i_moins1,i_plus1,i_plus2
819integer :: j_moins1,j_plus1,j_plus2
820     
821      !$OMP PARALLEL
822      !$OMP DO
823      do j=3,ny-2
824        do i=3,nx-2
825
826 surice:if  (ice(i,j).eq.0) then
827         do pmx=-1,1,2
828          do pmy=-1,1,2
829
830 diagice :  if (ice(i+pmx,j+pmy).eq.1) then
831
832              if ((ice(i+pmx,j)+ice(i,j+pmy).eq.2)) then    ! test (i) pour eviter les langues
833                                                            ! de glaces diagonales en coin(26dec04)
834                   if ((ice(i+2*pmx,j).eq.1.and.ice(i+2*pmx,j+pmy).eq.0).or.&
835                       (ice(i,j+2*pmy).eq.1.and.ice(i+pmx,j+2*pmy).eq.0)) then
836                           ice(i,j)=1
837                           h(i,j)=max(1.,h(i,j))
838                   endif         
839
840! test (i) pour eviter les langues de glaces diagonales :
841!                        mouvement du cheval aux echecs 
842
843               if ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.1)) then
844                   if (ice(i+2*pmx,j+pmy).eq.1.and. &
845                      (ice(i+2*pmx,j+2*pmy)+ice(i,j+2*pmy)).ge.1) then
846                           ice(i,j)=1
847                           h(i,j)=max(1.,h(i,j))
848                   endif       
849                   if (ice(i+pmx,j+2*pmy).eq.1.and. &
850                      (ice(i+2*pmx,j+2*pmy)+ice(i+2*pmx,j)).ge.1) then
851                           ice(i,j)=1
852                           h(i,j)=max(1.,h(i,j))
853                   endif       
854
855! test (ii) pour eviter les langues de glaces diagonales :
856!            le point glace ice(i+pmx,j+pmy) a :
857!                          - ses 4 voisins frontaux en glace
858!                          - mais 2 voisins vides diagonalement opposes   
859           
860               elseif ((ice(i+2*pmx,j+pmy)+ice(i+pmx,j+2*pmy).eq.2)  &
861                                .and.ice(i+2*pmx,j+2*pmy).eq.0) then
862
863! test (iii) pour faire les tests (i) et (ii)
864                    ice(i,j)=1
865                    h(i,j)=max(1.,h(i,j))
866                    ice(i+2*pmx,j+2*pmy)=1
867                    h(i+2*pmx,j+2*pmy)=max(1.,h(i+2*pmx,j+2*pmy))
868               endif
869            endif
870           endif diagice   
871          enddo
872         enddo
873         endif surice       
874       end do
875      end do
876      !$OMP END DO
877      !$OMP ENd PARALLEL
878
879!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
880!!$if (itestf.gt.0) then
881!!$   write(6,*) 'dans front avant remplissage baies asymetrie sur H  pour time=',time
882!!$   stop
883!!$else
884!!$   write(6,*) 'dans front avant remplissage baies pas d asymetrie sur H  pour time=',time
885!!$
886!!$end if
887
888
889!     print*,'dans remplissage baies',time
890       
891baies: do k=1,2
892!$OMP PARALLEL
893!$OMP DO PRIVATE(i_moins1,j_moins1,i_plus1,j_plus1,i_plus2,j_plus2)
894do j=1,ny
895   do i=1,nx
896
897surice_xy: if  (ice(i,j).eq.0) then
898   i_moins1=max(i-1,1)
899   j_moins1=max(j-1,1)
900   i_plus1=min(i+1,nx)
901   j_plus1=min(j+1,ny)
902   i_plus2=min(i+2,nx)
903   j_plus2=min(j+2,ny)
904 
905! test (iii) pour trouver les baies de largeur 1 ou 2 cases
906! et combler les trous si ce sont des baies
907! si ce ne sont pas des baies, ne pas combler et creer des langues de glaces artificielles             
908! baies horizontales 
909 
910         if (ice(i_moins1,j).eq.1.and.(ice(i_plus1,j).eq.1.or.ice(i_plus2,j).eq.1)) then
911            if (ice(i,j_moins1).eq.1.or.ice(i,j_plus1).eq.1)  then    ! ice(i,j)=1
912               ice(i,j)=1
913               H(i,j)=max(1.,H(i,j))
914            endif
915         endif
916
917
918         if (ice(i,j_moins1).eq.1.and.(ice(i,j_plus1).eq.1.or.ice(i,j_plus2).eq.1)) then
919            if (ice(i_moins1,j).eq.1.or.ice(i_plus1,j).eq.1) then !ice(i,j)=1
920               ice(i,j)=1
921               H(i,j)=max(1.,H(i,j))
922            endif
923         endif
924
925      endif surice_xy
926   end do
927end do
928!$OMP END DO
929!$OMP END PARALLEL
930end do baies
931
932!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
933!!$if (itestf.gt.0) then
934!!$   write(6,*) 'dans front apres remplissage baies asymetrie sur H  pour time=',time
935!!$   stop
936!!$else
937!!$   write(6,*) 'dans front apres remplissage baies pas d asymetrie sur H  pour time=',time
938!!$
939!!$end if
940
941!$OMP PARALLEL
942!$OMP DO
943do j=2,ny-1
944   do i=2,nx-1
945
946      if (ice(i,j).eq.1) then     !         test si ice=1
947
948! if ice, on determine front...
949! ainsi, front=0 sur les zones = 0
950
951         front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
952         !front= le nb de faces en contact avec un voisin englacé
953      endif
954   end do
955end do
956!$OMP END DO
957
958! traitement des bords. On considere que l'exterieur n'a pas de glace
959! attention ce n'est vrai que sur la grande grille
960
961!$OMP DO PRIVATE(i)
962do j=2,ny-1
963   i=1
964   front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
965   i=nx
966   front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1))
967end do
968!$OMP END DO
969
970!$OMP DO PRIVATE(j)
971do i=2,nx-1
972   j=1 
973   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1))
974   j=ny
975   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1))
976end do
977!$OMP END DO
978
979! traitement des coins
980
981front(1,1)=ice(2,1)+ice(2,1)
982front(1,ny)=ice(2,ny)+ice(1,ny-1)
983front(nx,1)=ice(nx,2)+ice(nx-1,1)
984front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny)
985
986!!$call detect_assym(nx,ny,0,41,1,0,1,0,H,itestf)
987!!$if (itestf.gt.0) then
988!!$   write(6,*) 'dans front apres front asymetrie sur H  pour time=',time
989!!$   stop
990!!$else
991!!$   write(6,*) 'dans front apres front pas d asymetrie sur H  pour time=',time
992!!$
993!!$end if
994
995!   on ne compte pas les taches de glace de 2 cases (horizontales ou verticales)
996!   en fait, si ces deux cases sont flottantes, il faut enlever les icebergs
997!   de n'importe quelle taille).
998!   si ces deux taches sont posées (ou une des deux), il n'y a pas assez de conditions aux limites
999
1000!$OMP DO
1001do j=1,ny
1002   do i=1,nx-1
1003      if (front(i,j).eq.1) then
1004         if (front(i+1,j).eq.1) then
1005            ice(i,j)=0
1006            ice(i+1,j)=0
1007            front(i,j)=0
1008            front(i+1,j)=0
1009         endif
1010      endif
1011   end do
1012end do
1013!$OMP END DO
1014
1015!$OMP DO
1016do j=1,ny-1
1017   do i=1,nx
1018      if (front(i,j).eq.1) then
1019         if (front(i,j+1).eq.1) then
1020            ice(i,j)=0
1021            ice(i,j+1)=0
1022            front(i,j)=0
1023            front(i,j+1)=0
1024         endif
1025      end if
1026   end do
1027end do
1028!$OMP END DO
1029
1030!isolx signifie pas de voisins en x
1031!isoly signifie pas de voisins en y
1032!remarque :
1033!si isolx/y=.true. alors frontfacex/y=0 (a la fois +1 & -1 or +1-1=0)
1034
1035! calcul de frontfacex et isolx
1036!$OMP DO
1037do j=1,ny
1038   do i=2,nx-1
1039
1040      if (front(i,j).ge.1.and.front(i,j).le.3) then    !front(entre 1 et 3)
1041         
1042         if ((ice(i-1,j)+ice(i+1,j)).lt.2) then        ! il y a un front // a x
1043         
1044            if ((ice(i-1,j)+ice(i+1,j)).eq.0) then
1045               isolx(i,j)=.true.
1046            elseif (ice(i-1,j).eq.0) then
1047               frontfacex(i,j)=-1                      ! front  i-1 |i  i+1
1048            else
1049               frontfacex(i,j)=+1                      ! front  i-1  i| i+1
1050            endif
1051         endif
1052      end if !fin du test il y a un front
1053         
1054   end do
1055end do
1056!$OMP END DO
1057
1058! calcul de frontfacey et isoly
1059!$OMP DO
1060do j=2,ny-1
1061   do i=1,nx
1062
1063      if (front(i,j).ge.1.and.front(i,j).le.3) then   !front(entre 1 et 3)
1064         
1065         if ((ice(i,j-1)+ice(i,j+1)).lt.2) then       ! il y a un front // a y
1066             
1067            if ((ice(i,j-1)+ice(i,j+1)).eq.0) then
1068               isoly(i,j)=.true.                      !front   j-1 |j| j+1
1069            elseif (ice(i,j-1).eq.0) then
1070               frontfacey(i,j)=-1                     !front   j-1 |j j+1
1071            else
1072               frontfacey(i,j)=+1                     !front   j-1  j| j+1
1073            endif
1074         endif
1075     end if !fin du test il y a un front
1076         
1077   end do
1078end do
1079!$OMP END DO
1080
1081
1082! traitement des bords. On considere que l'exterieur n'a pas de glace
1083! attention ce n'est vrai que sur la grande grille
1084
1085!$OMP DO PRIVATE(i)
1086do j=2,ny-1
1087   i=1
1088   if (front(i,j).ge.1)  then
1089      if (ice(i+1,j).eq.0) then
1090         isolx(i,j)=.true.
1091      else
1092         frontfacex(i,j)=-1
1093      endif
1094   end if
1095   i=nx
1096   if (front(i,j).ge.1)  then
1097      if (ice(i-1,j).eq.0) then
1098         isolx(i,j)=.true.
1099      else
1100         frontfacex(i,j)=1
1101      endif
1102   end if
1103end do
1104!$OMP END DO
1105
1106!$OMP DO PRIVATE(j)
1107do i=2,nx-1
1108   j=1 
1109   if (front(i,j).ge.1)  then
1110      if (ice(i,j+1).eq.0) then
1111         isoly(i,j)=.true.
1112      else
1113         frontfacey(i,j)=-1
1114      endif
1115   end if
1116   j=ny
1117   if (front(i,j).ge.1)  then
1118      if (ice(i,j-1).eq.0) then
1119         isoly(i,j)=.true.
1120      else
1121         frontfacey(i,j)=1
1122      endif
1123   end if
1124end do
1125!$OMP END DO
1126!$OMP END PARALLEL
1127
1128return
1129end subroutine determin_front
1130!------------------------------------------------------------------------------
1131
1132subroutine determin_front_tof
1133
1134integer,dimension(nx,ny) :: nofront ! tableau de travail (points au dela du front)
1135!$OMP PARALLEL
1136!$OMP DO
1137do j=2,ny-1
1138   do i=2,nx-1
1139
1140      if (ice(i,j).eq.1) then     !         test si ice=1
1141
1142! if ice, on determine front...
1143! ainsi, front=0 sur les zones = 0
1144
1145         front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1146         !front= le nb de faces en contact avec un voisin englacé
1147      endif
1148   end do
1149end do
1150!$OMP END DO
1151
1152! traitement des bords. On considere que l'exterieur n'a pas de glace
1153! attention ce n'est vrai que sur la grande grille
1154
1155!$OMP DO PRIVATE(i)
1156do j=2,ny-1
1157   i=1
1158   front(i,j)=(ice(i+1,j)+ice(i,j+1)+ice(i,j-1))
1159   i=nx
1160   front(i,j)=(ice(i-1,j)+ice(i,j+1)+ice(i,j-1))
1161end do
1162!$OMP END DO
1163
1164!$OMP DO PRIVATE(j)
1165do i=2,nx-1
1166   j=1 
1167   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j+1))
1168   j=ny
1169   front(i,j)=(ice(i-1,j)+ice(i+1,j)+ice(i,j-1))
1170end do
1171!$OMP END DO
1172!$OMP BARRIER
1173! traitement des coins
1174
1175front(1,1)=ice(2,1)+ice(2,1)
1176front(1,ny)=ice(2,ny)+ice(1,ny-1)
1177front(nx,1)=ice(nx,2)+ice(nx-1,1)
1178front(nx,ny)=ice(nx,ny-1)+ice(nx-1,ny)
1179
1180! Les points à plus d'un point de grille du bord sont front=-1
1181!$OMP WORKSHARE
1182nofront(:,:)=0
1183!$OMP END WORKSHARE
1184!$OMP BARRIER
1185
1186!$OMP DO
1187do j=2,ny-1
1188        do i=2,nx-1
1189                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
1190                        nofront(i,j)=-1
1191                endif
1192        enddo
1193enddo
1194!$OMP END DO
1195!$OMP BARRIER
1196
1197!$OMP WORKSHARE
1198where (nofront(:,:).eq.-1)
1199        front(:,:)=-1
1200endwhere
1201!$OMP END WORKSHARE
1202!$OMP END PARALLEL             
1203
1204end subroutine determin_front_tof
1205
1206
1207!> SUBROUTINE: determin_marais
1208!! afq -- Routine pour l'identification des "marais"
1209!! Un marais est une tache "shelf" entouré de points grounded
1210!! Copie sauvage de determin_tache, adapte au probleme du marais
1211!>
1212subroutine determin_marais
1213
1214  !$ USE OMP_LIB
1215
1216  implicit none
1217 
1218  integer :: indice
1219  integer :: label         ! no des taches rencontrées dans le mask
1220  integer :: label_max     ! no temporaire maxi de tache rencontrées
1221  integer,parameter :: mask_nb = 2   ! version ou on ne compte pas les diagonales
1222  integer,dimension(mask_nb) :: mask ! numero de tache des points adjacents
1223
1224  integer,dimension(nx,ny)          :: table_out_marais    !< numeros de tache d'un point ij
1225  integer,dimension(0:n_ta_max)     :: compt_marais        !< contient les equivalence entre les taches
1226  integer,dimension(0:n_ta_max)     :: nb_pts_marais       !< indique le nombre de points par tache
1227  logical,dimension(0:n_ta_max)     :: marais              !< T si flottants entoure de poses, F sinon
1228
1229
1230! 1-initialisation
1231!-----------------
1232  label_max=1      ! numero de la tache, la premiere tache est notée 1
1233  label=1
1234  do i=1,n_ta_max
1235    compt_marais(i)=i
1236  enddo
1237  !$OMP PARALLEL
1238  !$OMP WORKSHARE
1239  table_out_marais(:,:) = 0
1240  marais(:)  = .true.
1241  nb_pts_marais(:) = 0
1242  !$OMP END WORKSHARE
1243  !$OMP END PARALLEL
1244
1245! 2-reperage des taches
1246!----------------------
1247
1248  do j=2,ny-1
1249    do i=2,nx-1
1250      if ((ice(i,j).ge.1).and.flot(i,j)) then ! on est sur la glace qui flotte-------------------!
1251
1252        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
1253          !  un des voisins est deja en glace
1254          mask(1) = table_out_marais(i-1,j)
1255          mask(2) = table_out_marais(i,j-1) 
1256          label = label_max
1257
1258          ! on determine la valeur de la tache minimun (>0) presente ds le masque
1259          do indice=1,mask_nb
1260            if (mask(indice).gt.0) label=min(label,mask(indice))
1261          enddo
1262
1263          ! on fixe la valeur de la tache voisine minimun au point etudie (via label)
1264          table_out_marais(i,j)=label
1265           
1266          !si un des voisins n'est pas glace alors la tache n'est pas un marais
1267          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
1268            marais(label)=.false.
1269          endif
1270
1271          ! si 2 taches differentes sont dans le masque, il faut les identifier dans compt_marais
1272          ! on lui affecte le numero de la tache fondamentale
1273
1274                                                                   ! exemple on est sur le point X : 5  X
1275          do indice=1,mask_nb                                      !                                   20
1276            if(mask(indice).gt.label) then                         ! mask(2)=20 > 5             
1277              compt_marais(mask(indice))=label                     ! compt_marais(20)=5
1278              if (.not.marais(mask(indice))) marais(label)=.false. ! si la tache n'etais pas un marais => marais =.false.  marais(-(-5))=.false.               
1279              where (table_out_marais(:,:).eq.mask(indice))        ! where table_out_marais(:,:)=mask(2)=20
1280                table_out_marais(:,:)=label                        ! table_out_marais(:,:)=label=5                 
1281              endwhere                                               
1282            endif
1283          enddo
1284           
1285        else !aucun des voisins est une tache
1286          table_out_marais(i,j)= label_max
1287          compt_marais(label_max)=label_max
1288           
1289          ! si un des voisins n'est pas glace alors la tache n'est pas un marais
1290          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
1291            marais(label_max)=.false.
1292          endif
1293          label_max  = label_max+1
1294          if (label_max.gt.n_ta_max) print*,'ATTENTION trop de taches=',label_max 
1295        endif
1296      else          ! on est pas sur une tache--------------------------------------------
1297        table_out_marais(i,j)=0    ! Pas necessaire (reecrit 0 sur 0)             
1298      endif         !---------------------------------------------------------------------
1299    enddo
1300  enddo
1301
1302  marais(0)=.false.
1303
1304  !$OMP PARALLEL
1305  !$OMP DO
1306  do j=1,ny
1307    do i=1,nx
1308      if (table_out_marais(i,j).ne.0) then
1309        nb_pts_marais(compt_marais(table_out_marais(i,j)))= nb_pts_marais(compt_marais(table_out_marais(i,j)))+1   
1310      endif
1311      flot_marais(i,j) = marais(table_out_marais(i,j))
1312    enddo
1313  enddo
1314  !$OMP END DO
1315  !$OMP END PARALLEL
1316
1317  debug_3D(:,:,122)=real(table_out_marais(:,:))
1318
1319end subroutine determin_marais
1320
1321end module flottab_mod
1322
1323
Note: See TracBrowser for help on using the repository browser.