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 | 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 | real :: mass_total |
---|
37 | integer :: testH ! somme des flux entrants en un point ij |
---|
38 | |
---|
39 | contains |
---|
40 | |
---|
41 | !--------------------------------------------------------------------------------------- |
---|
42 | subroutine init_calving |
---|
43 | |
---|
44 | namelist/calving/Hcoup_plateau,Hcoup_abysses,prof_plateau,prof_abysses,ifrange,meth_hcoup |
---|
45 | |
---|
46 | |
---|
47 | |
---|
48 | mass_total = 0d0 |
---|
49 | calv(:,:)=0. |
---|
50 | calv_dtt(:,:)=0. |
---|
51 | |
---|
52 | ! formats pour les ecritures dans 42 |
---|
53 | 428 format(A) |
---|
54 | |
---|
55 | ! lecture des parametres du run block eaubasale1 |
---|
56 | !-------------------------------------------------------------------- |
---|
57 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
58 | read(num_param,calving) |
---|
59 | |
---|
60 | |
---|
61 | write(num_rep_42,428)'!___________________________________________________________' |
---|
62 | write(num_rep_42,428) '&calving ! nom du bloc calving méthode Vincent ' |
---|
63 | write(num_rep_42,*) |
---|
64 | write(num_rep_42,*) 'Hcoup_plateau = ', Hcoup_plateau |
---|
65 | write(num_rep_42,*) 'Hcoup_abysses = ', Hcoup_abysses |
---|
66 | write(num_rep_42,*) 'prof_plateau = ', prof_plateau |
---|
67 | write(num_rep_42,*) 'prof_abysses = ', prof_abysses |
---|
68 | write(num_rep_42,*) 'ifrange = ', ifrange |
---|
69 | write(num_rep_42,*) 'meth_hcoup = ', meth_hcoup |
---|
70 | write(num_rep_42,*)'/' |
---|
71 | |
---|
72 | write(num_rep_42,428) '! Hcoup epaisseurs de coupure pour les zones peu prodondes et profondes' |
---|
73 | write(num_rep_42,428) '! Hcoup_plateau<Hcoup_abysses && prof_plateau<prof_abysses' |
---|
74 | write(num_rep_42,428) '! prof profondeur delimitant les zones peu prodondes et profondes' |
---|
75 | write(num_rep_42,428) '! ifrange=0 -> pas de traitement particulier sur les bords' |
---|
76 | write(num_rep_42,428) '! ifrange=1 -> traitement de Vincent avec ice shelves frangeants partout' |
---|
77 | write(num_rep_42,428) '! ifrange=2 -> ice shelves frangeants seulement si bm-bmelt positif' |
---|
78 | write(num_rep_42,*) '! meth_hcoup pour faire eventuellement varier Hcoup avec le climat' |
---|
79 | write(num_rep_42,*) '! Hcoup_plateau=',Hcoup_plateau |
---|
80 | write(num_rep_42,*) '! Hcoup_abysses=',Hcoup_abysses |
---|
81 | write(num_rep_42,*) '! prof_plateau=',prof_plateau |
---|
82 | write(num_rep_42,*) '! prof_abysses=',prof_abysses |
---|
83 | write(num_rep_42,*) |
---|
84 | |
---|
85 | |
---|
86 | ! afq -- coupure depend de la profondeur: |
---|
87 | ! |
---|
88 | ! hcoup prof_abysses |
---|
89 | ! ^ v |
---|
90 | ! | _______ hcoup_abysses |
---|
91 | ! | / |
---|
92 | ! | / |
---|
93 | ! | / |
---|
94 | ! | / |
---|
95 | ! | hcoup_plateau _______/ |
---|
96 | ! ^ |
---|
97 | ! prof_plateau |
---|
98 | |
---|
99 | Hcoup(:,:) = min ( max( & |
---|
100 | (-(Bsoc0(:,:)-sealevel) - prof_plateau)/(prof_abysses-prof_plateau) & |
---|
101 | *(hcoup_abysses-hcoup_plateau)+hcoup_plateau & |
---|
102 | , hcoup_plateau), hcoup_abysses ) |
---|
103 | |
---|
104 | if (meth_hcoup.eq.1) then |
---|
105 | Hcoup(:,:)=coefbmshelf*Hcoup(:,:) |
---|
106 | Hcoup(:,:)=min( max(Hcoup (:,:),Hcoup_plateau),Hcoup_abysses) |
---|
107 | else if (meth_hcoup.eq.2) then |
---|
108 | Hcoup(:,:)=coefbmshelf*Hcoup(:,:) |
---|
109 | Hcoup(:,:)=max(Hcoup(:,:),Hcoup_plateau) |
---|
110 | else if (meth_hcoup.eq.3) then |
---|
111 | Hcoup(:,:)=coefbmshelf*Hcoup(:,:) |
---|
112 | Hcoup(:,:)=min(max(Hcoup(:,:),0.),Hcoup_abysses) |
---|
113 | endif |
---|
114 | |
---|
115 | |
---|
116 | end subroutine init_calving |
---|
117 | !--------------------------------------------------------------------------------------- |
---|
118 | subroutine calving |
---|
119 | |
---|
120 | integer :: I_did_something ! pour la boucle sur le calving |
---|
121 | |
---|
122 | ! initialisation calving |
---|
123 | calv(:,:)=0. |
---|
124 | mass_total=0d0 |
---|
125 | |
---|
126 | ! calcul du dhdt lagrangien |
---|
127 | dhdt(:,:)=bm(:,:)-bmelt(:,:)-H(:,:)*(epsxx(:,:)+epsyy(:,:)) |
---|
128 | |
---|
129 | ! calcul du bilan surface et fond : divise par 2 car utilise dans des moyennes |
---|
130 | bil_tot(:,:)=0.5*(bm(:,:)-bmelt(:,:)) |
---|
131 | |
---|
132 | Hcoup(:,:) = min ( max( & |
---|
133 | (-(Bsoc(:,:)-sealevel) - prof_plateau)/(prof_abysses-prof_plateau) & |
---|
134 | *(hcoup_abysses-hcoup_plateau)+hcoup_plateau & |
---|
135 | , hcoup_plateau), hcoup_abysses ) |
---|
136 | |
---|
137 | if (meth_hcoup.eq.1) then |
---|
138 | Hcoup(:,:)=coefbmshelf*Hcoup(:,:) |
---|
139 | Hcoup(:,:)=min( max(Hcoup (:,:),Hcoup_plateau),Hcoup_abysses) |
---|
140 | else if (meth_hcoup.eq.2) then |
---|
141 | Hcoup(:,:)=coefbmshelf*Hcoup(:,:) |
---|
142 | Hcoup(:,:)=max(Hcoup(:,:),Hcoup_plateau) |
---|
143 | else if (meth_hcoup.eq.3) then |
---|
144 | Hcoup(:,:)=coefbmshelf*Hcoup(:,:) |
---|
145 | Hcoup(:,:)=min(max(Hcoup(:,:),0.),Hcoup_abysses) |
---|
146 | endif |
---|
147 | |
---|
148 | ! hauteur au dessus de la coupure |
---|
149 | hmhc(:,:)=H(:,:)-Hcoup(:,:) |
---|
150 | |
---|
151 | ! coupure de l'ice shelf |
---|
152 | !--------------------------- |
---|
153 | ! on divise le test en 2 'ifext' et 'ifint' pour reduire les tests redondants |
---|
154 | |
---|
155 | |
---|
156 | MULTI_CALV_LOOP : do l=1,10 |
---|
157 | ! calcul du dhdt lagrangien |
---|
158 | dhdt(:,:)=bm(:,:)-bmelt(:,:)-H(:,:)*(epsxx(:,:)+epsyy(:,:)) |
---|
159 | ! hauteur au dessus de la coupure |
---|
160 | hmhc(:,:)=H(:,:)-Hcoup(:,:) |
---|
161 | I_did_something = 0 |
---|
162 | do j=2,ny-1 |
---|
163 | do i=2,nx-1 |
---|
164 | |
---|
165 | |
---|
166 | ifext: if((flot(i,j)).and.(h(i,j).le.hcoup(i,j)).and.(h(i,j).gt.0.)) then |
---|
167 | ! ifext: pour les noeuds flottants avec h < hcoup |
---|
168 | |
---|
169 | !ifint: if((front(i,j).gt.0).and.(front(i,j).lt.4)) then |
---|
170 | !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 |
---|
171 | 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 |
---|
172 | ! ifint: le point doit avoir au - 1 voisin mais ne pas etre entouré de glace |
---|
173 | ! ce qui evite la formation des polynies dans les shelfs |
---|
174 | |
---|
175 | |
---|
176 | ! on regarde dhdt minimum pour voir si 1 des voisins peut alimenter le point |
---|
177 | |
---|
178 | ! hmhc est l'épaisseur en plus de l'épaisseur de coupure |
---|
179 | ! ux/dx=1/(deltat) ou deltat est le temps nécessaire |
---|
180 | ! à la glace pour passer d'un noeud à l'autre |
---|
181 | ! la deuxieme partie du test revient à dh/dt*deltat > hmhc |
---|
182 | |
---|
183 | ! Rajout vince : si le point amont est pose -> test vrai. |
---|
184 | |
---|
185 | ! rappel des differents ifrange |
---|
186 | !------------------------------- |
---|
187 | !Attention ifrange 1,2 ont sans doute un bug (il faudrait abs(uxbar)) |
---|
188 | |
---|
189 | ! ifrange=1 Rajout vince : si un point voisin est pose -> test vrai. |
---|
190 | ! comme dans article fenno2007 |
---|
191 | |
---|
192 | ! ifrange=2 pour les points ayant des voisins poses, garde le point |
---|
193 | ! seulement si le bilan surface fond est positif. |
---|
194 | ! bm(i,j)-bmelt(i,j) > 0. |
---|
195 | |
---|
196 | ! ifrange=3 pas de test sur l'epaisseur du point amont. |
---|
197 | ! test dhdt > -hmhc/deltat= -hmhc abs(U)/dx (correction du bug) |
---|
198 | ! + test sur bm(i,j)-bmelt(i,j) > 0. |
---|
199 | |
---|
200 | ! ifrange = 4 idem 3 mais avec test sur l'epaisseur du point amont |
---|
201 | |
---|
202 | ! ifrange = 0 pas de traitement specifique près du continent |
---|
203 | |
---|
204 | ! ifrange = -1 ancienne version MIS11 |
---|
205 | |
---|
206 | ! ifrange = 5 semble idem 0 |
---|
207 | |
---|
208 | ! ifrange = 6 ?? |
---|
209 | |
---|
210 | ! ifrange = 7 calving drastique sauf si coefbmshelf < 0.5 |
---|
211 | |
---|
212 | !cdc |
---|
213 | ! ifrange=44 idem 4 mais avec test sur i et j couple // afq: probably useless |
---|
214 | ! ifrange=8 C. Dumas : coupure si alimentation ne permet pas H > Hcoup |
---|
215 | |
---|
216 | if (ifrange.eq.1) then ! Rajout vince : si un point voisin est pose -> test vrai. |
---|
217 | ! comme dans article fenno2007 |
---|
218 | |
---|
219 | testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
220 | .and.(dhdt(i-1,j).gt.(hmhc(i-1,j)*uxbar(i-1,j)/dx)))& |
---|
221 | .or.(.not.flot(i-1,j)) ) ! |
---|
222 | |
---|
223 | testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
224 | .and.(dhdt(i+1,j).gt.(hmhc(i+1,j)*uxbar(i+1,j)/dx))) & |
---|
225 | .or.(.not.flot(i+1,j)) ) ! |
---|
226 | |
---|
227 | testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
228 | .and.(dhdt(i,j-1).gt.(hmhc(i,j-1)*uybar(i,j-1)/dx)))& |
---|
229 | .or.(.not.flot(i,j-1)) ) ! |
---|
230 | |
---|
231 | testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
232 | .and.(dhdt(i,j+1).gt.(hmhc(i,j+1)*uybar(i,j+1)/dx)))& |
---|
233 | .or.(.not.flot(i,j+1)) ) ! |
---|
234 | |
---|
235 | else if (ifrange.eq.2) then ! pour les points ayant des voisins posés, seulement si le bilan |
---|
236 | ! surface fond est positif. bm(i,j)-bmelt(i,j) > 0. |
---|
237 | |
---|
238 | bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.) |
---|
239 | |
---|
240 | testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
241 | .and. (dhdt(i-1,j).gt.(hmhc(i-1,j)*uxbar(i-1,j)/dx))) & |
---|
242 | .or.((.not.flot(i-1,j)).and.bilan_surf_fond )) ! |
---|
243 | |
---|
244 | testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
245 | .and.(dhdt(i+1,j).gt.(hmhc(i+1,j)*uxbar(i+1,j)/dx))) & |
---|
246 | .or.((.not.flot(i+1,j)).and.bilan_surf_fond) ) ! |
---|
247 | |
---|
248 | testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
249 | .and.(dhdt(i,j-1).gt.(hmhc(i,j-1)*uybar(i,j-1)/dx)))& |
---|
250 | .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ) ! |
---|
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)*uybar(i,j+1)/dx)))& |
---|
254 | .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ) ! |
---|
255 | |
---|
256 | else if (ifrange.eq.3) then ! nouvelle formulation Cat mars 08 |
---|
257 | |
---|
258 | ! pas de test sur l'epaisseur du point amont. |
---|
259 | ! test dhdt > -hmhc/deltat= -hmhc abs(U)/dx |
---|
260 | |
---|
261 | ! pour les points ayant des voisins posés, seulement si le bilan |
---|
262 | ! surface fond est positif. bm(i,j)-bmelt(i,j) > 0. |
---|
263 | |
---|
264 | |
---|
265 | bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.) |
---|
266 | |
---|
267 | testmij=( (uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
268 | .and.(dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j))/dx))) & |
---|
269 | .or.((.not.flot(i-1,j)).and.bilan_surf_fond ) ! |
---|
270 | |
---|
271 | testpij=( (uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
272 | .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx))) & |
---|
273 | .or.((.not.flot(i+1,j)).and.bilan_surf_fond ) ! |
---|
274 | |
---|
275 | testimj=(( uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
276 | .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx))) & |
---|
277 | .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ! |
---|
278 | |
---|
279 | testipj=((uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
280 | .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx))) & |
---|
281 | .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ! |
---|
282 | |
---|
283 | else if (ifrange.eq.4) then ! idem 3 mais avec test sur l'épaisseur du point amont |
---|
284 | |
---|
285 | bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.) |
---|
286 | |
---|
287 | testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
288 | .and. (dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) & |
---|
289 | .or.((.not.flot(i-1,j)).and.bilan_surf_fond )) ! |
---|
290 | |
---|
291 | testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
292 | .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j)/dx)))) & |
---|
293 | .or.((.not.flot(i+1,j)).and.bilan_surf_fond) ) ! |
---|
294 | |
---|
295 | testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
296 | .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1)/dx))))& |
---|
297 | .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ) ! |
---|
298 | |
---|
299 | testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
300 | .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1)/dx))))& |
---|
301 | .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ) ! |
---|
302 | |
---|
303 | else if (ifrange.eq.0) then ! pas de traitement special pres du continent |
---|
304 | |
---|
305 | testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
306 | .and.(dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) |
---|
307 | |
---|
308 | testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
309 | .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx))) |
---|
310 | |
---|
311 | testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
312 | .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx))) |
---|
313 | |
---|
314 | testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
315 | .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx))) |
---|
316 | |
---|
317 | else if (ifrange.eq.-1) then ! ancienne version MIS11 |
---|
318 | |
---|
319 | testmij=((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
320 | .and.(dhdt(i-1,j).gt.(hmhc(i-1,j)*uxbar(i-1,j)/dx))) |
---|
321 | |
---|
322 | testpij=((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
323 | .and.(dhdt(i+1,j).gt.(hmhc(i+1,j)*uxbar(i+1,j)/dx))) |
---|
324 | |
---|
325 | testimj=((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
326 | .and.(dhdt(i,j-1).gt.(hmhc(i,j-1)*uybar(i,j-1)/dx))) |
---|
327 | |
---|
328 | testipj=((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
329 | .and.(dhdt(i,j+1).gt.(hmhc(i,j+1)*uybar(i,j+1)/dx))) |
---|
330 | |
---|
331 | |
---|
332 | else if (ifrange.eq.5) then ! pas de traitement special pres du continent |
---|
333 | |
---|
334 | testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
335 | .and.((bil_tot(i-1,j)+bil_tot(i,j)) & |
---|
336 | .gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) |
---|
337 | |
---|
338 | testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
339 | .and.((bil_tot(i+1,j)+bil_tot(i,j)) & |
---|
340 | .gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx))) |
---|
341 | |
---|
342 | testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
343 | .and.((bil_tot(i,j-1)+bil_tot(i,j)) & |
---|
344 | .gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx))) |
---|
345 | |
---|
346 | testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
347 | .and.((bil_tot(i,j+1)+bil_tot(i,j)) & |
---|
348 | .gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx))) |
---|
349 | |
---|
350 | else if (ifrange.eq.6) then ! pas de traitement special pres du continent |
---|
351 | |
---|
352 | testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
353 | .and.(2.*bil_tot(i,j) & |
---|
354 | .gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) |
---|
355 | |
---|
356 | testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
357 | .and.(2.*bil_tot(i,j) & |
---|
358 | .gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j))/dx))) |
---|
359 | |
---|
360 | testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
361 | .and.(2.*bil_tot(i,j) & |
---|
362 | .gt.(-hmhc(i,j-1)*abs(uybar(i,j-1))/dx))) |
---|
363 | |
---|
364 | testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
365 | .and.(2.*bil_tot(i,j) & |
---|
366 | .gt.(-hmhc(i,j+1)*abs(uybar(i,j+1))/dx))) |
---|
367 | |
---|
368 | else if (ifrange.eq.7) then ! drastique calving sauf si coefbmelt.lt.0.5 |
---|
369 | |
---|
370 | testmij= ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
371 | .and.(coefbmshelf.lt.1.)) |
---|
372 | testpij= ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
373 | .and.(coefbmshelf.lt.1.)) |
---|
374 | testimj= ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
375 | .and.(coefbmshelf.lt.1.)) |
---|
376 | testipj= ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
377 | .and.(coefbmshelf.lt.1.)) |
---|
378 | |
---|
379 | else if (ifrange.eq.44) then ! idem 4 mais avec test sur i et j couple |
---|
380 | |
---|
381 | bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.) |
---|
382 | |
---|
383 | testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.) & ! voisin (i-1,j) amont et > hcoup |
---|
384 | .and. (dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) & |
---|
385 | .or.((.not.flot(i-1,j)).and.bilan_surf_fond )) ! |
---|
386 | |
---|
387 | testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.) & ! voisin (i+1,j) amont et > hcoup |
---|
388 | .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j)/dx)))) & |
---|
389 | .or.((.not.flot(i+1,j)).and.bilan_surf_fond) ) ! |
---|
390 | |
---|
391 | testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.) & ! voisin (i,j-1) amont et > hcoup |
---|
392 | .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1)/dx))))& |
---|
393 | .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ) ! |
---|
394 | |
---|
395 | testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.) & ! voisin (i,j+1) amont et > hcoup |
---|
396 | .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1)/dx))))& |
---|
397 | .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ) ! |
---|
398 | |
---|
399 | else if (ifrange.eq.8) then ! test C Dumas : si somme flux entrant permet H > Hcoup : vrai |
---|
400 | ! bilan_surf_fond=(bm(i,j)-bmelt(i,j).gt.0.) |
---|
401 | |
---|
402 | testmij=( ((hmhc(i-1,j).gt.0.).and.(uxbar(i,j).ge.0.)) & ! voisin (i-1,j) amont et > hcoup |
---|
403 | ! .and. (dhdt(i-1,j).gt.(-hmhc(i-1,j)*abs(uxbar(i-1,j)/dx)))) & |
---|
404 | .or.((.not.flot(i-1,j)).and.bilan_surf_fond )) ! |
---|
405 | |
---|
406 | testpij=( ((hmhc(i+1,j).gt.0.).and.(uxbar(i+1,j).le.0.)) & ! voisin (i+1,j) amont et > hcoup |
---|
407 | ! .and.(dhdt(i+1,j).gt.(-hmhc(i+1,j)*abs(uxbar(i+1,j)/dx)))) & |
---|
408 | .or.((.not.flot(i+1,j)).and.bilan_surf_fond) ) ! |
---|
409 | |
---|
410 | testimj=( ((hmhc(i,j-1).gt.0.).and.(uybar(i,j).ge.0.)) & ! voisin (i,j-1) amont et > hcoup |
---|
411 | ! .and.(dhdt(i,j-1).gt.(-hmhc(i,j-1)*abs(uybar(i,j-1)/dx))))& |
---|
412 | .or.((.not.flot(i,j-1)).and.bilan_surf_fond ) ) ! |
---|
413 | |
---|
414 | testipj=( ((hmhc(i,j+1).gt.0.).and.(uybar(i,j+1).le.0.)) & ! voisin (i,j+1) amont et > hcoup |
---|
415 | ! .and.(dhdt(i,j+1).gt.(-hmhc(i,j+1)*abs(uybar(i,j+1)/dx))))& |
---|
416 | .or.((.not.flot(i,j+1)).and.bilan_surf_fond ) ) ! |
---|
417 | |
---|
418 | testH=0. |
---|
419 | if (testmij) testH=testH + abs(uxbar(i,j))*H(i-1,j) |
---|
420 | if (testpij) testH=testH + abs(uxbar(i+1,j))*H(i+1,j) |
---|
421 | if (testimj) testH=testH + abs(uybar(i,j))*H(i,j-1) |
---|
422 | if (testipj) testH=testH + abs(uybar(i,j+1))*H(i,j+1) |
---|
423 | |
---|
424 | if (testH.gt.Hcoup(i,j)) then |
---|
425 | testmij=.true. |
---|
426 | else |
---|
427 | testmij=.false. |
---|
428 | testpij=.false. |
---|
429 | testimj=.false. |
---|
430 | testipj=.false. |
---|
431 | endif |
---|
432 | |
---|
433 | endif |
---|
434 | |
---|
435 | |
---|
436 | |
---|
437 | ! detection des polynies dans les shelfs. on regarde vers l'aval |
---|
438 | |
---|
439 | ! 1 des 3 voisins aval > hcoup du cote west (i-1) |
---|
440 | iin2=max(1,i-2) |
---|
441 | iin3=max(1,i-3) |
---|
442 | |
---|
443 | avalw=((hmhc(i-1,j).gt.0.).or.(hmhc(iin2,j).gt.0.) & |
---|
444 | .or.(hmhc(iin3,j).gt.0.)).and.(uxbar(i,j).le.0.) |
---|
445 | |
---|
446 | !avalw=((hmhc(i-1,j).gt.0..and.flot(i-1,j)).or.(hmhc(iin2,j).gt.0..and.flot(iin2,j)) & |
---|
447 | ! .or.(hmhc(iin3,j).gt.0..and.flot(iin3,j))).and.(uxbar(i,j).le.0.) |
---|
448 | |
---|
449 | ! 1 des 3 voisins aval > hcoup du cote est (i+1) |
---|
450 | iin2=min(i+2,nx) |
---|
451 | iin3=min(i+3,nx) |
---|
452 | |
---|
453 | avale=((hmhc(i+1,j).gt.0.).or.(hmhc(iin2,j).gt.0.) & |
---|
454 | .or.(hmhc(iin3,j).gt.0.)).and.(uxbar(i+1,j).ge.0) |
---|
455 | |
---|
456 | !avale=((hmhc(i+1,j).gt.0..and.flot(i+1,j)).or.(hmhc(iin2,j).gt.0..and.flot(iin2,j)) & |
---|
457 | ! .or.(hmhc(iin3,j).gt.0..and.flot(iin3,j))).and.(uxbar(i+1,j).ge.0) |
---|
458 | |
---|
459 | ! 1 des 3 voisins aval > hcoup du cote sud (j-1) |
---|
460 | jin2=max(1,j-2) |
---|
461 | jin3=max(1,j-3) |
---|
462 | |
---|
463 | avals=((hmhc(i,j-1).gt.0.).or.(hmhc(i,jin2).gt.0.) & |
---|
464 | .or.(hmhc(i,jin3).gt.0.)).and.(uybar(i,j).le.0.) |
---|
465 | |
---|
466 | !avals=((hmhc(i,j-1).gt.0..and.flot(i,j-1)).or.(hmhc(i,jin2).gt.0..and.flot(i,jin2)) & |
---|
467 | ! .or.(hmhc(i,jin3).gt.0..and.flot(i,jin3))).and.(uybar(i,j).le.0.) |
---|
468 | |
---|
469 | ! 1 des 3 voisins aval > hcoup du cote nord (j-1) |
---|
470 | jin2=min(j+2,ny) |
---|
471 | jin3=min(j+3,ny) |
---|
472 | |
---|
473 | avaln=((hmhc(i,j+1).gt.0.).or.(hmhc(i,jin2).gt.0.) & |
---|
474 | .or.(hmhc(i,jin2).gt.0.)) .and.(uybar(i,j+1).ge.0.) |
---|
475 | |
---|
476 | !avaln=((hmhc(i,j+1).gt.0..and.flot(i,j+1)).or.(hmhc(i,jin2).gt.0..and.flot(i,jin2)) & |
---|
477 | ! .or.(hmhc(i,jin2).gt.0..and.flot(i,jin2))) .and.(uybar(i,j+1).ge.0.) |
---|
478 | |
---|
479 | interieur=(avalw.or.avale).and.(avals.or.avaln) |
---|
480 | |
---|
481 | |
---|
482 | |
---|
483 | if ((.not.(testmij.or.testpij.or.testimj.or.testipj)) & ! pas suffisament alimente |
---|
484 | .and.(.not.interieur)) then ! et pas interieur |
---|
485 | |
---|
486 | !cdc pour ne pas faire de calving si H > 350m : if ((.not.(testmij.or.testpij.or.testimj.or.testipj)).and.(H(i,j).LT.350)) then |
---|
487 | ! mab: if none of the neighbours doesn't provide ice for the considered point, |
---|
488 | ! mab: then the ice is cut at this point and given to calv(i,j) |
---|
489 | !if ((.not.(testmij.or.testpij.or.testimj.or.testipj))) then |
---|
490 | ! mab: new version so that grid points with h=1m (all points were calving takes place are afterwards |
---|
491 | ! mab: set to a height of 1m) are not taken into account, provides a problem as soon as you start |
---|
492 | ! mab: coupling GRISLI to the iceberg modell |
---|
493 | !if(h(i,j).gt.1) then |
---|
494 | calv(i,j)=-h(i,j) |
---|
495 | !cdc H(i,j)=1. |
---|
496 | !cdc 1m H(i,j)=min(1.,max(0.,(sealevel - Bsoc(i,j))*row/ro-0.01)) |
---|
497 | H(i,j)=0. |
---|
498 | S(i,j)=H(i,j)*(1.-ro/row) + sealevel |
---|
499 | B(i,j)=S(i,j) - H(i,j) |
---|
500 | ! ATTENTION ne pas mettre ice=0 sinon degradation bilan d'eau (bm et bmelt non comptabilises dans ce cas) |
---|
501 | |
---|
502 | I_did_something = I_did_something + 1 |
---|
503 | endif |
---|
504 | |
---|
505 | end if ifint |
---|
506 | end if ifext |
---|
507 | end do |
---|
508 | end do |
---|
509 | if (I_did_something.eq.0) then |
---|
510 | ! print*,'stop MULTI_CALV_LOOP l=',l |
---|
511 | EXIT MULTI_CALV_LOOP |
---|
512 | else |
---|
513 | ! print*,'calving continue l',l,I_did_something |
---|
514 | endif |
---|
515 | enddo MULTI_CALV_LOOP |
---|
516 | |
---|
517 | |
---|
518 | ! on met en calving les points detectes iceberg : |
---|
519 | where (iceberg(:,:).and.(H(:,:).gt.0.)) |
---|
520 | calv(:,:)=-h(:,:) |
---|
521 | ice(:,:)=0 |
---|
522 | H(:,:)=0. |
---|
523 | S(:,:)=H(:,:)*(1.-ro/row) + sealevel |
---|
524 | B(:,:)=S(:,:) - H(:,:) |
---|
525 | endwhere |
---|
526 | |
---|
527 | calv_dtt(:,:) = calv_dtt(:,:) + calv(:,:) ! somme du calving sur dtt |
---|
528 | ! calv_dtt est remis a 0 dans steps_time_loop (tous les dtt) |
---|
529 | |
---|
530 | end subroutine calving |
---|
531 | !------------------------------------------------------------------------------------------ |
---|
532 | end module calving_frange |
---|
533 | |
---|