source: trunk/SOURCES/calving_frange.f90 @ 150

Last change on this file since 150 was 148, checked in by aquiquet, 7 years ago

Minor bug correction in flottab + compatibility between ant40 and ant16 geometries + get rid of unused variable in eaubasale.

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