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