source: branches/iLoveclim/SOURCES/calving_frange.f90

Last change on this file was 244, checked in by aquiquet, 5 years ago

Grisli-iloveclim branch merged to trunk at revision 243

File size: 21.3 KB
Line 
1!> \file calving_frange.f90
2!! Module et routines qui calculent le calving avec la methode de Vincent
3!<
4
5
6!> \namespace module3D_phy
7!! Module et routines qui calculent le calving avec la methode de Vincent
8!! \author ...
9!! \date ...
10!! @note recopie direct de icethick5-ant. A revoir en particulier, accroche avec la glace posee
11!! @note Used modules
12!! @note   - use module3D_phy
13!! @todo il faudrait que le calving soit inactif dans les zones ou l epaisseur est imposee
14!<
15
16module calving_frange
17
18  use module3D_phy
19  use bilan_eau_mod
20  implicit none
21
22  real, dimension (nx,ny) :: hmhc  ! hauteur au dessus de la coupure
23  real, dimension (nx,ny) :: bil_tot ! bilan surface et fond (bm-bmelt)
24  real, dimension (nx,ny) :: hcoup ! epaiseur de coupure au temps time
25  real :: hcoup_plateau  ! coupure points peu profonds
26  real :: hcoup_abysses  ! coupure points ocean profond
27  real :: prof_plateau  ! profondeur max des points peu profonds
28  real :: prof_abysses  ! profondeur min des points ocean profond
29  integer :: meth_hcoup  ! pour avoir hcoup dépendant du coefbmshelf
30  integer :: ifrange        !  0 pas de traitement particulier pres du bord, 1 -> franges
31  integer ::  iin2,jin2,iin3,jin3             ! pour la detection polynies
32  logical :: testmij,testpij,testimj,testipj
33  logical :: bilan_surf_fond                  ! vrai si bm-bmelt est positif
34  logical :: avalw,avale,avals,avaln,interieur
35
36
37contains
38
39  !---------------------------------------------------------------------------------------
40  subroutine init_calving
41
42    namelist/calving/Hcoup_plateau,Hcoup_abysses,prof_plateau,prof_abysses,ifrange,meth_hcoup
43
44
45
46    calv(:,:)=0. 
47    calv_dtt(:,:)=0.
48
49    ! formats pour les ecritures dans 42
50428 format(A)
51
52    ! lecture des parametres du run                      block eaubasale1
53    !--------------------------------------------------------------------
54    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
55    read(num_param,calving)
56
57
58    write(num_rep_42,428)'!___________________________________________________________' 
59    write(num_rep_42,428) '&calving              ! nom du bloc calving méthode Vincent '
60    write(num_rep_42,*)
61    write(num_rep_42,*)   'Hcoup_plateau  = ', Hcoup_plateau
62    write(num_rep_42,*)   'Hcoup_abysses  = ', Hcoup_abysses
63    write(num_rep_42,*)   'prof_plateau  = ', prof_plateau
64    write(num_rep_42,*)   'prof_abysses  = ', prof_abysses
65    write(num_rep_42,*)   'ifrange    = ', ifrange
66    write(num_rep_42,*)   'meth_hcoup = ', meth_hcoup
67    write(num_rep_42,*)'/'
68
69    write(num_rep_42,428) '! Hcoup epaisseurs de coupure pour les zones peu prodondes et profondes'
70    write(num_rep_42,428) '! Hcoup_plateau<Hcoup_abysses && prof_plateau<prof_abysses'
71    write(num_rep_42,428) '! prof profondeur delimitant les zones peu prodondes et profondes'
72    write(num_rep_42,428) '! ifrange=0 -> pas de traitement particulier sur les bords'
73    write(num_rep_42,428) '! ifrange=1 -> traitement de Vincent avec ice shelves frangeants partout'
74    write(num_rep_42,428) '! ifrange=2 -> ice shelves frangeants seulement si bm-bmelt positif'
75    write(num_rep_42,*)   '! meth_hcoup pour faire eventuellement varier Hcoup avec le climat'
76
77
78    ! afq -- coupure depend de la profondeur:
79    !
80    !  hcoup                       prof_abysses
81    !   ^                               v
82    !   |                                _______ hcoup_abysses
83    !   |                               /
84    !   |                              /
85    !   |                             /
86    !   |                            /
87    !   |      hcoup_plateau _______/
88    !                               ^
89    !                          prof_plateau     
90
91    Hcoup(:,:) = min ( max(                                         &
92       (-(Bsoc0(:,:)-sealevel_2d(:,:)) - prof_plateau)/(prof_abysses-prof_plateau)    &
93       *(hcoup_abysses-hcoup_plateau)+hcoup_plateau               &
94       , hcoup_plateau), hcoup_abysses )
95
96    if (meth_hcoup.eq.1) then
97        Hcoup(:,:)=coefbmshelf*Hcoup(:,:)
98        Hcoup(:,:)=min( max(Hcoup (:,:),Hcoup_plateau),Hcoup_abysses)
99    else if (meth_hcoup.eq.2) then
100        Hcoup(:,:)=coefbmshelf*Hcoup(:,:)
101        Hcoup(:,:)=max(Hcoup(:,:),Hcoup_plateau)
102    else if (meth_hcoup.eq.3) then
103        Hcoup(:,:)=coefbmshelf*Hcoup(:,:)
104        Hcoup(:,:)=min(max(Hcoup(:,:),0.),Hcoup_abysses)   
105    endif
106
107
108  end subroutine init_calving
109  !---------------------------------------------------------------------------------------
110  subroutine calving
111
112    integer :: I_did_something  ! pour la boucle sur le calving
113
114    ! initialisation calving
115    calv(:,:)=0. 
116
117    ! calcul du dhdt lagrangien
118    dhdt(:,:)=bm(:,:)-bmelt(:,:)-H(:,:)*(epsxx(:,:)+epsyy(:,:))
119
120    ! calcul du bilan surface et fond : divise par 2 car utilise dans des moyennes
121    bil_tot(:,:)=0.5*(bm(:,:)-bmelt(:,:))
122
123    Hcoup(:,:) = min ( max(                                         &
124       (-(Bsoc(:,:)-sealevel_2d(:,:)) - prof_plateau)/(prof_abysses-prof_plateau)    &
125       *(hcoup_abysses-hcoup_plateau)+hcoup_plateau               &
126       , hcoup_plateau), hcoup_abysses )
127
128    if (meth_hcoup.eq.1) then
129        Hcoup(:,:)=coefbmshelf*Hcoup(:,:)
130        Hcoup(:,:)=min( max(Hcoup (:,:),Hcoup_plateau),Hcoup_abysses)
131    else if (meth_hcoup.eq.2) then
132        Hcoup(:,:)=coefbmshelf*Hcoup(:,:)
133        Hcoup(:,:)=max(Hcoup(:,:),Hcoup_plateau)
134    else if (meth_hcoup.eq.3) then
135        Hcoup(:,:)=coefbmshelf*Hcoup(:,:)
136        Hcoup(:,:)=min(max(Hcoup(:,:),0.),Hcoup_abysses)
137    endif
138
139    ! hauteur au dessus de la coupure
140    hmhc(:,:)=H(:,:)-Hcoup(:,:)
141
142    ! coupure de l'ice shelf
143    !---------------------------
144    ! on divise le test en 2 'ifext' et 'ifint' pour reduire les tests redondants 
145
146
147
148    MULTI_CALV_LOOP : do l=1,10
149      ! calcul du dhdt lagrangien
150      dhdt(:,:)=bm(:,:)-bmelt(:,:)-H(:,:)*(epsxx(:,:)+epsyy(:,:))
151      ! hauteur au dessus de la coupure
152      hmhc(:,:)=H(:,:)-Hcoup(:,:)
153      I_did_something = 0
154      do j=2,ny-1
155        do i=2,nx-1
156
157
158          ifext:  if((flot(i,j)).and.(h(i,j).le.hcoup(i,j)).and.(h(i,j).gt.0.))  then
159              ! ifext: pour les noeuds flottants englaces avec h < hcoup
160
161              !ifint:    if((front(i,j).gt.0).and.(front(i,j).lt.4)) then
162              !cdc pb avec front    ifint:    if((front(i,j).lt.4).or.((front(i-1,j)+front(i+1,j)+front(i,j-1)+front(i,j+1)).lt.16)) then
163              ifint:    if((H(i-1,j).lt.2.).or.(H(i+1,j).lt.2.).or.(H(i,j-1).lt.2.).or.(H(i,j+1).lt.2)) then ! si on est au bord avec test sur H
164                  ! ifint: le point doit avoir au - 1 voisin mais ne pas etre entouré de glace
165                  ! ce qui evite la formation des polynies dans les shelfs
166
167
168                  ! on regarde dhdt minimum pour voir si 1 des voisins peut alimenter le point
169
170                  !                     hmhc est l'épaisseur en plus de l'épaisseur de coupure
171                  !                     ux/dx=1/(deltat) ou deltat est le temps nécessaire
172                  !                     à la glace pour passer d'un noeud à l'autre
173                  !                     la deuxieme partie du test revient à dh/dt*deltat > hmhc 
174
175                  !                     Rajout vince : si le point amont est pose -> test vrai.
176
177                  ! rappel des differents ifrange
178                  !-------------------------------
179                  !Attention ifrange 1,2  ont sans doute un bug (il faudrait abs(uxbar))
180
181                  ! ifrange=1               Rajout vince : si un point voisin est pose -> test vrai.
182                  !                         comme dans article fenno2007
183
184                  ! ifrange=2               pour les points ayant des voisins poses, garde le point
185                  !                         seulement si le bilan  surface fond est positif.
186                  !                         bm(i,j)-bmelt(i,j) > 0.
187
188                  ! ifrange=3               pas de test sur l'epaisseur du point amont.
189                  !                         test dhdt > -hmhc/deltat= -hmhc abs(U)/dx (correction du bug)
190                  !                         + test sur bm(i,j)-bmelt(i,j) > 0.
191
192                  ! ifrange = 4             idem 3 mais avec test sur l'epaisseur du point amont
193
194                  ! ifrange = 0             pas de traitement specifique près du continent
195
196                  ! ifrange = -1            ancienne version MIS11
197
198                  ! ifrange = 5             semble idem 0
199
200                  ! ifrange = 6             ??
201
202                  ! ifrange = 7             calving drastique sauf si coefbmshelf < 0.5
203
204
205                  if (ifrange.eq.1) then           !   Rajout vince : si un point voisin est pose -> test vrai.
206                      !   comme dans article fenno2007
207
208                      testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
209                         .and.(dhdt(i-1,j).gt.(hmhc(i-1,j)*uxbar(i-1,j)/dx)))& 
210                         .or.(.not.flot(i-1,j)) ) !
211
212                      testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
213                         .and.(dhdt(i+1,j).gt.(hmhc(i+1,j)*uxbar(i+1,j)/dx))) &
214                         .or.(.not.flot(i+1,j)) ) !
215
216                      testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
217                         .and.(dhdt(i,j-1).gt.(hmhc(i,j-1)*uybar(i,j-1)/dx)))&
218                         .or.(.not.flot(i,j-1)) ) !
219
220                      testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
221                         .and.(dhdt(i,j+1).gt.(hmhc(i,j+1)*uybar(i,j+1)/dx)))&
222                         .or.(.not.flot(i,j+1)) ) !
223
224                  else if (ifrange.eq.2) then       ! pour les points ayant des voisins posés, seulement si le bilan
225                      ! surface fond est positif. bm(i,j)-bmelt(i,j) > 0.
226
227                      bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.)
228
229                      testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
230                         .and.  (dhdt(i-1,j).gt.(hmhc(i-1,j)*uxbar(i-1,j)/dx))) & 
231                         .or.((.not.flot(i-1,j)).and.bilan_surf_fond )) !
232
233                      testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
234                         .and.(dhdt(i+1,j).gt.(hmhc(i+1,j)*uxbar(i+1,j)/dx))) &
235                         .or.((.not.flot(i+1,j)).and.bilan_surf_fond) ) !
236
237                      testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
238                         .and.(dhdt(i,j-1).gt.(hmhc(i,j-1)*uybar(i,j-1)/dx)))&
239                         .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ) !
240
241                      testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
242                         .and.(dhdt(i,j+1).gt.(hmhc(i,j+1)*uybar(i,j+1)/dx)))&
243                         .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ) !
244
245                  else if (ifrange.eq.3) then       ! nouvelle formulation Cat mars 08
246
247                      ! pas de test sur l'epaisseur du point amont.
248                      ! test dhdt > -hmhc/deltat= -hmhc abs(U)/dx
249
250                      ! pour les points ayant des voisins posés, seulement si le bilan
251                      ! surface fond est positif. bm(i,j)-bmelt(i,j) > 0.
252
253
254                      bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.)
255
256                      testmij=( (uxbar(i,j).ge.0.) &  ! voisin (i-1,j) amont et > hcoup
257                         .and.(dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j))/dx))) & 
258                         .or.((.not.flot(i-1,j)).and.bilan_surf_fond ) !
259
260                      testpij=( (uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
261                         .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx))) &
262                         .or.((.not.flot(i+1,j)).and.bilan_surf_fond )  !
263
264                      testimj=(( uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
265                         .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx))) &
266                         .or.((.not.flot(i,j-1)).and.bilan_surf_fond )  !
267
268                      testipj=((uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
269                         .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx))) &
270                         .or.((.not.flot(i,j+1)).and.bilan_surf_fond )  !
271
272                  else if (ifrange.eq.4) then       ! idem 3 mais avec test sur l'épaisseur du point amont
273
274                      bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.)
275
276                      testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
277                         .and.  (dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) & 
278                         .or.((.not.flot(i-1,j)).and.bilan_surf_fond )) !
279
280                      testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
281                         .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j)/dx)))) &
282                         .or.((.not.flot(i+1,j)).and.bilan_surf_fond) ) !
283
284                      testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
285                         .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1)/dx))))&
286                         .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ) !
287
288                      testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
289                         .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1)/dx))))&
290                         .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ) !
291
292                  else if (ifrange.eq.0) then           ! pas de traitement special pres du continent
293
294                      testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
295                         .and.(dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx))))
296
297                      testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
298                         .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx)))
299
300                      testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
301                         .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx)))
302
303                      testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
304                         .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx)))
305
306                  else if (ifrange.eq.-1) then           ! ancienne version MIS11
307
308                      testmij=((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
309                         .and.(dhdt(i-1,j).gt.(hmhc(i-1,j)*uxbar(i-1,j)/dx)))   
310
311                      testpij=((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
312                         .and.(dhdt(i+1,j).gt.(hmhc(i+1,j)*uxbar(i+1,j)/dx)))
313
314                      testimj=((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
315                         .and.(dhdt(i,j-1).gt.(hmhc(i,j-1)*uybar(i,j-1)/dx)))
316
317                      testipj=((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
318                         .and.(dhdt(i,j+1).gt.(hmhc(i,j+1)*uybar(i,j+1)/dx)))
319
320
321                  else if (ifrange.eq.5) then           ! pas de traitement special pres du continent
322
323                      testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
324                         .and.((bil_tot(i-1,j)+bil_tot(i,j)) &
325                         .gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx))))
326
327                      testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
328                         .and.((bil_tot(i+1,j)+bil_tot(i,j)) &
329                         .gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx)))
330
331                      testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
332                         .and.((bil_tot(i,j-1)+bil_tot(i,j)) &
333                         .gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx)))
334
335                      testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
336                         .and.((bil_tot(i,j+1)+bil_tot(i,j)) &
337                         .gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx)))
338
339                  else if (ifrange.eq.6) then           ! pas de traitement special pres du continent
340
341                      testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
342                         .and.(2.*bil_tot(i,j) &
343                         .gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx))))
344
345                      testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
346                         .and.(2.*bil_tot(i,j) &
347                         .gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx)))
348
349                      testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
350                         .and.(2.*bil_tot(i,j) &
351                         .gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx)))
352
353                      testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
354                         .and.(2.*bil_tot(i,j) &
355                         .gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx)))
356
357                  else if (ifrange.eq.7) then           ! drastique calving sauf si coefbmelt.lt.0.5
358
359                      testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)  &  ! voisin (i-1,j) amont et > hcoup
360                         .and.(coefbmshelf.lt.1.))
361                      testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup
362                         .and.(coefbmshelf.lt.1.))
363                      testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)  &  ! voisin (i,j-1) amont et > hcoup
364                         .and.(coefbmshelf.lt.1.))
365                      testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup
366                         .and.(coefbmshelf.lt.1.))
367
368
369                  endif
370
371
372
373                  ! detection des polynies dans les shelfs. on regarde vers l'aval
374
375                  ! 1 des 3 voisins aval > hcoup du cote west (i-1)
376                  iin2=max(1,i-2)
377                  iin3=max(1,i-3)
378
379                  avalw=((hmhc(i-1,j).gt.0.).or.(hmhc(iin2,j).gt.0.) &
380                     .or.(hmhc(iin3,j).gt.0.)).and.(uxbar(i,j).le.0.)                         
381
382
383                  ! 1 des 3 voisins aval > hcoup du cote est (i+1)
384                  iin2=min(i+2,nx)
385                  iin3=min(i+3,nx)
386
387                  avale=((hmhc(i+1,j).gt.0.).or.(hmhc(iin2,j).gt.0.) &   
388                     .or.(hmhc(iin3,j).gt.0.)).and.(uxbar(i+1,j).ge.0)                         
389
390                  ! 1 des 3 voisins aval > hcoup du cote sud (j-1)
391                  jin2=max(1,j-2)
392                  jin3=max(1,j-3)
393
394                  avals=((hmhc(i,j-1).gt.0.).or.(hmhc(i,jin2).gt.0.) &   
395                     .or.(hmhc(i,jin3).gt.0.)).and.(uybar(i,j).le.0.) 
396
397                  ! 1 des 3 voisins aval > hcoup du cote nord (j-1)
398                  jin2=min(j+2,ny)
399                  jin3=min(j+3,ny)
400
401                  avaln=((hmhc(i,j+1).gt.0.).or.(hmhc(i,jin2).gt.0.) &   
402                     .or.(hmhc(i,jin2).gt.0.)) .and.(uybar(i,j+1).ge.0.)                         
403
404                  interieur=(avalw.or.avale).and.(avals.or.avaln)
405
406
407
408                  if ((.not.(testmij.or.testpij.or.testimj.or.testipj))  &  ! pas suffisament alimente
409                     .and.(.not.interieur)) then                ! et pas interieur
410                      calv(i,j)=-h(i,j) 
411                      !cdc                     H(i,j)=1.
412                      !cdc 1m                     H(i,j)=min(1.,max(0.,(sealevel_2d(i,j) - Bsoc(i,j))*row/ro-0.01))
413                      H(i,j)=0.
414                      S(i,j)=H(i,j)*(1.-ro/row) + sealevel_2d(i,j) !afq -- WARNING: est-ce qu'on veut vraiment mettre S a la valeur locale du niveau marin?
415                      B(i,j)=S(i,j) - H(i,j)
416                      ! ATTENTION ne pas mettre ice=0 sinon degradation bilan d'eau (bm et bmelt non comptabilises dans ce cas)
417
418                      I_did_something = I_did_something + 1
419                      !   if (l.ge.2) then
420                      !     print*,'calving l ij',l,i,j
421                      !   endif
422                  endif
423
424              end if ifint
425          end if ifext
426        end do
427      end do
428      if (I_did_something.eq.0) then
429          !        print*,'stop MULTI_CALV_LOOP l=',l
430          EXIT MULTI_CALV_LOOP
431      else
432          !        print*,'calving continue l',l,I_did_something
433      endif
434    enddo MULTI_CALV_LOOP
435
436
437    ! on met en calving les points detectes iceberg :
438    where (iceberg(:,:).and.(H(:,:).gt.0.))
439        calv(:,:)=-h(:,:)
440        ice(:,:)=0
441        H(:,:)=0.
442        S(:,:)=H(:,:)*(1.-ro/row) + sealevel_2d(:,:) !afq -- WARNING: est-ce qu'on veut vraiment mettre S a la valeur locale du niveau marin?
443        B(:,:)=S(:,:) - H(:,:)
444    endwhere
445
446    calv_dtt(:,:) = calv_dtt(:,:) + calv(:,:) ! somme du calving sur dtt 
447    ! calv_dtt est remis a 0 dans steps_time_loop (tous les dtt)
448  end subroutine calving
449  !------------------------------------------------------------------------------------------
450end module calving_frange
451
Note: See TracBrowser for help on using the repository browser.