source: branches/GRISLIv3/SOURCES/calving_frange.f90 @ 409

Last change on this file since 409 was 383, checked in by dumas, 16 months ago

use only in module calving_frange

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