1 | !> \file conserv-mass-adv-diff_paral_mod.f90 |
---|
2 | !! Module qui calcule l evolution de l epaisseur en resolvant |
---|
3 | !! une equation d'advection diffusion |
---|
4 | !! version adaptee a un domaine parallele (MISMIP) |
---|
5 | !! en utilisant use reso_adv_diff_2D_paral |
---|
6 | !< |
---|
7 | |
---|
8 | !> \namespace equat_adv_diff_2D_paral |
---|
9 | !! Calcule l evolution de l epaisseur en resolvant |
---|
10 | !! une equation d'advection diffusion |
---|
11 | !! @note Version avec : call resol_adv_diff_2D_vect (Dmx,Dmy,advmx,advmy,H_p,i_Hp,bilmass,vieuxH,H) |
---|
12 | !! @note - le terme vecteur est passe en dummy parametres |
---|
13 | !! l epaisseur peut etre imposee |
---|
14 | !! \author CatRitz |
---|
15 | !! \date june 2009 |
---|
16 | !! @note Used modules: |
---|
17 | !! @note - use module3D_phy |
---|
18 | !! @note - use reso_adv_diff_2D |
---|
19 | !! @todo essayer dans le futur la methode de Richard dUX/dx = H dU/dx + U dH/dx |
---|
20 | !< |
---|
21 | module equat_adv_diff_2D_paral ! Cat mars 2012 |
---|
22 | use module3D_phy |
---|
23 | use reso_adv_diff_2D_paral |
---|
24 | use prescribe_H |
---|
25 | |
---|
26 | implicit none |
---|
27 | |
---|
28 | real, dimension(nx,ny) :: advmx !< partie advective en x |
---|
29 | real, dimension(nx,ny) :: advmy !< partie advective en y |
---|
30 | real, dimension(nx,ny) :: dmx !< partie advective en x |
---|
31 | real, dimension(nx,ny) :: dmy !< partie advective en y |
---|
32 | real, dimension(nx,ny) :: VieuxH !< ancienne valeur de H |
---|
33 | real, dimension(nx,ny) :: bilmass !< bilan de masse pour la colonne |
---|
34 | logical, dimension(nx,ny) :: zonemx !< pour separer advection-diffusion |
---|
35 | logical, dimension(nx,ny) :: zonemy !< pour separer advection-diffusion |
---|
36 | real :: adv_frac !< fraction du flux traitee en advection |
---|
37 | integer :: itesti |
---|
38 | integer :: itour |
---|
39 | |
---|
40 | |
---|
41 | contains |
---|
42 | |
---|
43 | !> SUBROUTINE: init_icethick |
---|
44 | !! Definition et initialisation des parametres qui gerent la repartition advection diffusion |
---|
45 | !! @return adv_frac |
---|
46 | |
---|
47 | subroutine init_icethick |
---|
48 | |
---|
49 | |
---|
50 | write(num_rep_42,*)'Conservation de la masse avec equation advection-diffusion ' |
---|
51 | write(num_rep_42,*)'la repartition depend de adv_frac' |
---|
52 | write(num_rep_42,*)'>1 -> advection seule' |
---|
53 | write(num_rep_42,*)'0 -> diffusion seule' |
---|
54 | write(num_rep_42,*)'0<*<1 -> fraction de l advection' |
---|
55 | write(num_rep_42,*) '-1 -> zones diffusion + zones advecttion' |
---|
56 | |
---|
57 | !varname='adv_frac' |
---|
58 | |
---|
59 | adv_frac=2. ! 0.8 |
---|
60 | |
---|
61 | !adv_frac=0 |
---|
62 | write(num_rep_42,*)'adv_frac=',adv_frac |
---|
63 | |
---|
64 | call init_reso_adv_diff_2D |
---|
65 | call init_prescribe_H ! initialize grounding line mask |
---|
66 | |
---|
67 | Dx1=1./Dx |
---|
68 | itour=0 |
---|
69 | |
---|
70 | |
---|
71 | end subroutine init_icethick |
---|
72 | |
---|
73 | !------------------------------------------------------------------ |
---|
74 | !> SUBROUTINE: icethick3 |
---|
75 | !! Routine pour le calcule de l'epaisseur de glace |
---|
76 | !< |
---|
77 | subroutine icethick3 |
---|
78 | |
---|
79 | implicit none |
---|
80 | real,dimension(nx,ny) :: Dminx,Dminy |
---|
81 | real,dimension(nx,ny) :: Uxdiff,Uydiff ! vitesse due a la diffusion |
---|
82 | integer :: it1,it2,it3 |
---|
83 | |
---|
84 | real aux ! pour le calcul du critere |
---|
85 | real maxdia ! sur le pas de temps |
---|
86 | integer imaxdia,jmaxdia |
---|
87 | |
---|
88 | |
---|
89 | if (itracebug.eq.1) call tracebug(" Entree dans routine icethick") |
---|
90 | |
---|
91 | ! Copie de H sur vieuxH |
---|
92 | ! -------------------- |
---|
93 | where (mk0(:,:).eq.0) |
---|
94 | vieuxH(:,:)=0. |
---|
95 | elsewhere (i_Hp(:,:).eq.1) ! previous prescribed thickness |
---|
96 | vieuxH(:,:)=Hp(:,:) ! H may have been changed by calving |
---|
97 | elsewhere |
---|
98 | vieuxH(:,:)=H(:,:) |
---|
99 | end where |
---|
100 | |
---|
101 | ! mk0 est le masque à ne jamais dépasser |
---|
102 | ! mk0=0 -> pas de glace |
---|
103 | ! mk0=-1 -> on garde la même epaisseur |
---|
104 | ! pour l'Antarctique masque mko=1 partout |
---|
105 | |
---|
106 | |
---|
107 | Maxdia = -1.0 |
---|
108 | imaxdia=0 |
---|
109 | jmaxdia=0 |
---|
110 | |
---|
111 | ! attention limitation ! |
---|
112 | uxbar(:,:)=max(uxbar(:,:),-3000.) |
---|
113 | uxbar(:,:)=min(uxbar(:,:),3000.) |
---|
114 | uybar(:,:)=max(uybar(:,:),-3000.) |
---|
115 | uybar(:,:)=min(uybar(:,:),3000.) |
---|
116 | |
---|
117 | ! le pas de temps est choisi pour rendre la matrice advective diagonale dominante |
---|
118 | do i=2,nx-1 |
---|
119 | do j=2,ny-1 |
---|
120 | aux = (abs(uxbar(i,j)) + abs(uxbar(i+1,j))) & |
---|
121 | + (abs(uybar(i,j)) + abs(uybar(i,j+1))) |
---|
122 | aux=aux*dx1*0.5 |
---|
123 | if (aux.gt.maxdia) then |
---|
124 | maxdia=aux |
---|
125 | imaxdia=i |
---|
126 | jmaxdia=j |
---|
127 | endif |
---|
128 | end do |
---|
129 | end do |
---|
130 | |
---|
131 | |
---|
132 | if (maxdia.ge.(testdiag/dtmax)) then |
---|
133 | dt = testdiag/Maxdia |
---|
134 | dt=max(dt,dtmin) |
---|
135 | else |
---|
136 | dt = dtmax |
---|
137 | end if |
---|
138 | |
---|
139 | ! next_time ajuste le pas de temps pour faciliter la synchronisation |
---|
140 | ! avec le pas de temps long (dtt) |
---|
141 | |
---|
142 | call next_time(time,dt,dtt,dtmax,dtmin,isynchro,itracebug,num_tracebug) |
---|
143 | |
---|
144 | |
---|
145 | ! calcul diffusivite - vitesse |
---|
146 | !----------------------------------------------------------------- |
---|
147 | |
---|
148 | Uxdiff(:,:)=diffmx(:,:)*(-sdx(:,:)) ! vitesse qui peut s'exprimer en terme diffusif |
---|
149 | Uydiff(:,:)=diffmy(:,:)*(-sdy(:,:)) ! dans les where qui suivent elle est limitee a uxbar |
---|
150 | ! pour eviter des problemes dans le cas 0< adv_frac <1 |
---|
151 | dmx(:,:)=diffmx(:,:) |
---|
152 | dmy(:,:)=diffmy(:,:) |
---|
153 | |
---|
154 | ! schema amont pour la diffusion en x |
---|
155 | where (Uxbar(:,:).ge.0.) |
---|
156 | dmx(:,:)=diffmx(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=1) |
---|
157 | uxdiff(:,:)=min(uxdiff(:,:),uxbar(:,:)) ! pour tenir compte limitation utot |
---|
158 | elsewhere |
---|
159 | dmx(:,:)=diffmx(:,:)*H(:,:) |
---|
160 | uxdiff(:,:)=max(uxdiff(:,:),uxbar(:,:)) ! pour tenir compte limitation utot |
---|
161 | end where |
---|
162 | |
---|
163 | |
---|
164 | |
---|
165 | ! schema amont pour la diffusion en y |
---|
166 | where (uybar(:,:).ge.0.) |
---|
167 | dmy(:,:)=diffmy(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=2) |
---|
168 | uydiff(:,:)=min(uydiff(:,:),uybar(:,:)) ! pour tenir compte limitation utot |
---|
169 | elsewhere |
---|
170 | dmy(:,:)=dmy(:,:)*H(:,:) |
---|
171 | uydiff(:,:)=max(uydiff(:,:),uybar(:,:)) ! pour tenir compte limitation utot |
---|
172 | end where |
---|
173 | |
---|
174 | |
---|
175 | |
---|
176 | ! tests pour la répartition advection - diffusion |
---|
177 | !------------------------------------------------ |
---|
178 | |
---|
179 | if (adv_frac.gt.1.) then ! advection seulement |
---|
180 | advmx(:,:)=Uxbar(:,:) |
---|
181 | advmy(:,:)=Uybar(:,:) |
---|
182 | Dmx(:,:)=0. |
---|
183 | Dmy(:,:)=0. |
---|
184 | |
---|
185 | else if (abs(adv_frac).lt.1.e-8) then ! diffusion seulement |
---|
186 | advmx(:,:)=0. |
---|
187 | advmy(:,:)=0. |
---|
188 | |
---|
189 | ! D reste la valeur au dessus) |
---|
190 | |
---|
191 | else if ((adv_frac.ge.1.e-8).and.(adv_frac.le.1.)) then ! advection-diffusion) |
---|
192 | |
---|
193 | ! selon x -------------------------------- |
---|
194 | |
---|
195 | ! advection = adv_frac* tout ce qui n'est pas diffusion |
---|
196 | ! diffusion = ce qui peut être exprime en diffusion |
---|
197 | ! + une partie (1-adv_frac) de l'advection |
---|
198 | |
---|
199 | where ((abs(sdx(:,:)).gt.1.e-8).and.(Hmx(:,:).gt.2.)) ! cas general |
---|
200 | advmx(:,:)=(Uxbar(:,:)-Uxdiff(:,:)) ! tout ce qui n'est pas diffusion |
---|
201 | advmx(:,:)=advmx(:,:)*adv_frac ! la fraction adv_frac du precedent |
---|
202 | Dminx(:,:)=-(Uxbar(:,:)-advmx(:,:))/sdx(:,:) ! ce qui reste exprime en diffusivite |
---|
203 | ! a multiplier par H |
---|
204 | |
---|
205 | ! schema amont pour la diffusivite |
---|
206 | where (uxbar(:,:).ge.0.) |
---|
207 | dmx(:,:)=dminx(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=1) |
---|
208 | elsewhere |
---|
209 | dmx(:,:)=dminx(:,:)*H(:,:) |
---|
210 | end where |
---|
211 | |
---|
212 | elsewhere ! zones trop plates ou trop fines |
---|
213 | Dmx(:,:)=0. |
---|
214 | advmx(:,:)=Uxbar(:,:) |
---|
215 | end where |
---|
216 | |
---|
217 | |
---|
218 | ! selon y -------------------------------- |
---|
219 | |
---|
220 | ! advection = adv_frac* tout ce qui n'est pas diffusion |
---|
221 | ! diffusion = ce qui peut être exprime en diffusion |
---|
222 | ! + une partie (1-adv_frac) de l'advection |
---|
223 | |
---|
224 | where ((abs(sdy(:,:)).gt.1.e-8).and.(Hmy(:,:).gt.2.)) ! cas general |
---|
225 | advmy(:,:)=(Uybar(:,:)-Uydiff(:,:)) ! tout ce qui n'est pas diffusion |
---|
226 | advmy(:,:)=advmy(:,:)*adv_frac ! la fraction adv_frac du precedent |
---|
227 | Dminy(:,:)=-(Uybar(:,:)-advmy(:,:))/sdy(:,:) ! ce qui reste exprime en diffusivite |
---|
228 | ! a multiplier par H |
---|
229 | |
---|
230 | ! schema amont pour la diffusivite |
---|
231 | where (uybar(:,:).ge.0.) |
---|
232 | dmy(:,:)=dminy(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=2) |
---|
233 | elsewhere |
---|
234 | dmy(:,:)=dminy(:,:)*H(:,:) |
---|
235 | end where |
---|
236 | |
---|
237 | |
---|
238 | |
---|
239 | elsewhere ! zones trop plates ou trop fines |
---|
240 | Dmy(:,:)=0. |
---|
241 | advmy(:,:)=Uybar(:,:) |
---|
242 | end where |
---|
243 | |
---|
244 | |
---|
245 | else if (adv_frac.lt.0) then ! diffusif dans zones SIA, advectif ailleurs |
---|
246 | |
---|
247 | zonemx(:,:)=flgzmx(:,:) |
---|
248 | zonemy(:,:)=flgzmy(:,:) |
---|
249 | |
---|
250 | ! selon x -------------- |
---|
251 | where (zonemx(:,:)) |
---|
252 | advmx(:,:)=Uxbar(:,:) |
---|
253 | Dmx(:,:)=0. |
---|
254 | elsewhere |
---|
255 | advmx(:,:)=0. |
---|
256 | end where |
---|
257 | |
---|
258 | ! selon y -------------- |
---|
259 | where (zonemy(:,:)) |
---|
260 | advmy(:,:)=Uybar(:,:) |
---|
261 | Dmy(:,:)=0. |
---|
262 | elsewhere |
---|
263 | advmy(:,:)=0. |
---|
264 | end where |
---|
265 | |
---|
266 | end if |
---|
267 | |
---|
268 | |
---|
269 | where (flot(:,:) ) |
---|
270 | tabtest(:,:)=1. |
---|
271 | elsewhere |
---|
272 | tabtest(:,:)=0. |
---|
273 | end where |
---|
274 | |
---|
275 | !------------------------------------------------------------------detect |
---|
276 | !!$ tabtest(:,:)=dmx(:,:) |
---|
277 | !!$ |
---|
278 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,tabtest,itesti) |
---|
279 | !!$ |
---|
280 | !!$if (itesti.gt.0) then |
---|
281 | !!$ write(6,*) 'asymetrie sur dmx avant resol pour time=',time,itesti |
---|
282 | !!$ stop |
---|
283 | !!$else |
---|
284 | !!$ write(6,*) ' pas d asymetrie sur dmx avant resol pour time=',time,itesti |
---|
285 | !!$end if |
---|
286 | !----------------------------------------------------------- fin detect |
---|
287 | |
---|
288 | |
---|
289 | !------------------------------------------------------------------detect |
---|
290 | !!$ tabtest(:,:)=dmy(:,:) |
---|
291 | !!$ |
---|
292 | !!$call detect_assym(nx,ny,0,41,1,0,1,1,tabtest,itesti) |
---|
293 | !!$ |
---|
294 | !!$if (itesti.gt.0) then |
---|
295 | !!$ write(6,*) 'asymetrie sur dmy avant resol pour time=',time,itesti |
---|
296 | !!$ stop |
---|
297 | !!$else |
---|
298 | !!$ write(6,*) ' pas d asymetrie sur dmy avant resol pour time=',time,itesti |
---|
299 | !!$end if |
---|
300 | !----------------------------------------------------------- fin detect |
---|
301 | |
---|
302 | !------------------------------------------------------------------detect |
---|
303 | !!$ tabtest(:,:)=advmx(:,:) |
---|
304 | !!$ |
---|
305 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,tabtest,itesti) |
---|
306 | !!$ |
---|
307 | !!$if (itesti.gt.0) then |
---|
308 | !!$ write(6,*) 'asymetrie sur advmx avant resol pour time=',time,itesti |
---|
309 | !!$ stop |
---|
310 | !!$else |
---|
311 | !!$ write(6,*) ' pas d asymetrie sur advdmx avant resol pour time=',time,itesti |
---|
312 | !!$end if |
---|
313 | !----------------------------------------------------------- fin detect |
---|
314 | |
---|
315 | !------------------------------------------------------------------detect |
---|
316 | !!$ tabtest(:,:)=advmy(:,:) |
---|
317 | !!$ |
---|
318 | !!$call detect_assym(nx,ny,0,41,1,0,-1,1,tabtest,itesti) |
---|
319 | !!$ |
---|
320 | !!$if (itesti.gt.0) then |
---|
321 | !!$ write(6,*) 'asymetrie sur advmy avant resol pour time=',time,itesti |
---|
322 | !!$ stop |
---|
323 | !!$else |
---|
324 | !!$ write(6,*) ' pas d asymetrie sur advmy avant resol pour time=',time,itesti |
---|
325 | !!$end if |
---|
326 | !----------------------------------------------------------- fin detect |
---|
327 | |
---|
328 | ! nouveau calcul de dt |
---|
329 | aux=maxval( (abs(advmx(:,:))+abs(advmy(:,:)))/dx+(abs(dmx(:,:))+abs(dmy(:,:)))/dx/dx) |
---|
330 | |
---|
331 | ! write(166,*) 'critere pour dt',time-dt,testdiag/aux,dt |
---|
332 | |
---|
333 | |
---|
334 | if (aux.gt.1.e-20) then |
---|
335 | if (testdiag/aux.lt.dt) then |
---|
336 | time=time-dt |
---|
337 | dt=testdiag/aux*4. |
---|
338 | if (itracebug.eq.1) call tracebug("aux tres petit,deuxieme appel a next_time dt*4") |
---|
339 | call next_time(time,dt,dtt,dtmax,dtmin,isynchro,itracebug,num_tracebug) |
---|
340 | end if |
---|
341 | end if |
---|
342 | |
---|
343 | |
---|
344 | timemax=time |
---|
345 | |
---|
346 | ! etait avant dans step |
---|
347 | if (time.lt.TGROUNDED) then |
---|
348 | MARINE=.false. |
---|
349 | else |
---|
350 | MARINE=.true. |
---|
351 | endif |
---|
352 | ! fin vient de step |
---|
353 | |
---|
354 | ! les variables dtdx et dtdx2 sont globales |
---|
355 | Dtdx2=Dt/(Dx**2) |
---|
356 | dtdx=dt/dx |
---|
357 | |
---|
358 | |
---|
359 | debug_3D(:,:,45)=dmx(:,:) |
---|
360 | debug_3D(:,:,46)=dmy(:,:) |
---|
361 | debug_3D(:,:,47)=advmx(:,:) |
---|
362 | debug_3D(:,:,48)=advmy(:,:) |
---|
363 | |
---|
364 | |
---|
365 | bilmass(:,:)=bm(:,:)-bmelt(:,:) ! surface and bottom mass balance |
---|
366 | |
---|
367 | ! diverses precription de l'epaisseur en fonction de l'objectif |
---|
368 | !------------------------------------------------------------------- |
---|
369 | call prescribe_fixed_points ! ceux qui sont fixes pour tout le run (bords de grille, regions) |
---|
370 | |
---|
371 | |
---|
372 | if ((igrdline.eq.1)) then ! present grounding line |
---|
373 | |
---|
374 | if ((time.lt.time_gl(1)).or.(nt.lt.2)) then |
---|
375 | call prescribe_present_H_gl ! prescribe present thickness at the grounding line |
---|
376 | else |
---|
377 | call prescribe_retreat |
---|
378 | endif |
---|
379 | end if |
---|
380 | |
---|
381 | if ((igrdline.eq.2)) then ! paleo grounding line |
---|
382 | call prescribe_paleo_gl_shelf |
---|
383 | end if |
---|
384 | |
---|
385 | if ((igrdline.eq.3)) then |
---|
386 | call prescr_ice2sea_retreat |
---|
387 | endif |
---|
388 | |
---|
389 | |
---|
390 | if (time.ge.t_break) then ! action sur les ice shelves |
---|
391 | call melt_ice_shelves |
---|
392 | ! call break_all_ice_shelves |
---|
393 | end if |
---|
394 | |
---|
395 | |
---|
396 | !!$! impose une variation d'épaisseur d'apres une carte a un temps donne (a affiner plus tard pour des runs transient) |
---|
397 | !!$if (time.eq.t_break) then |
---|
398 | !!$ call lect_prescribed_deltaH |
---|
399 | !!$else if (time.gt.t_break) then |
---|
400 | !!$ call prescribe_deltaH |
---|
401 | !!$end if |
---|
402 | !!$ |
---|
403 | !!$if (time.eq.t_break) then ! si appele apres lect_prescribed_deltaH, cumule les effets |
---|
404 | !!$ call break_all_ice_shelves |
---|
405 | !!$end if |
---|
406 | |
---|
407 | debug_3D(:,:,87)=S(:,:)-debug_3D(:,:,88) |
---|
408 | |
---|
409 | debug_3D(:,:,86)=Hp(:,:) |
---|
410 | debug_3D(:,:,85)=i_Hp(:,:) |
---|
411 | |
---|
412 | !!$where (i_hp(:,:).eq.1) |
---|
413 | !!$ vieuxh(:,:)=Hp(:,:) |
---|
414 | !!$endwhere |
---|
415 | |
---|
416 | |
---|
417 | ! Appel a la routine de resolution resol_adv_diff_2D_vect |
---|
418 | !------------------------------------------------------------ |
---|
419 | ! On ecrit les vitesses sous la forme |
---|
420 | ! Ux = -Dmx * sdx + advmx |
---|
421 | |
---|
422 | ! Dmx, Dmy sont les termes diffusifs de la vitesse |
---|
423 | ! advmx, advmy sont les termes advectifs. |
---|
424 | ! la repartition entre les deux depend de adv_frac. |
---|
425 | |
---|
426 | ! i_Hp est un tableau des points ou l'epaisseur est imposee a la valeur Hp |
---|
427 | ! bilmass est le bilan de masse (surface et fond) |
---|
428 | ! vieuxH est la precedente valeur de H |
---|
429 | ! H est le retour de la nouvelle valeur. |
---|
430 | |
---|
431 | |
---|
432 | call resol_adv_diff_2D_vect (Dmx,Dmy,advmx,advmy,Hp,i_Hp,bilmass,vieuxH,H) |
---|
433 | |
---|
434 | ! remise à 0 des epaisseurs negatives, on garde la difference dans ablbord |
---|
435 | where (H(:,:).lt.0.) |
---|
436 | ablbord(:,:)=H(:,:)/dt |
---|
437 | elsewhere |
---|
438 | ablbord(:,:)=0. |
---|
439 | endwhere |
---|
440 | |
---|
441 | H(:,:)=max(H(:,:),0.) ! pas d'epaisseur negative |
---|
442 | |
---|
443 | |
---|
444 | ! calcul du masque "ice" |
---|
445 | where (flot(:,:)) ! points flottants, sera éventuellement réévalué dans flottab |
---|
446 | H(:,:)=max(H(:,:),1.) ! dans la partie marine l'épaisseur mini est 1 m |
---|
447 | where(H(:,:).gt.1.) |
---|
448 | ice(:,:)=1 |
---|
449 | elsewhere |
---|
450 | ice(:,:)=0 |
---|
451 | endwhere |
---|
452 | elsewhere ! points posés |
---|
453 | where(H(:,:).gt.0.) |
---|
454 | ice(:,:)=1 |
---|
455 | elsewhere |
---|
456 | ice(:,:)=0 |
---|
457 | endwhere |
---|
458 | endwhere |
---|
459 | |
---|
460 | |
---|
461 | ! eventuellement retirer apres spinup ou avoir un cas serac |
---|
462 | Hdot(:,:)=(H(:,:)-vieuxH(:,:))/dt |
---|
463 | |
---|
464 | if (igrdline.ne.3) then |
---|
465 | Hdot(:,:)=min(Hdot(:,:),10.) |
---|
466 | Hdot(:,:)=max(Hdot(:,:),-10.) |
---|
467 | endif |
---|
468 | |
---|
469 | where (i_hp(:,:).ne.1) |
---|
470 | H(:,:)=vieuxH(:,:)+Hdot(:,:)*dt |
---|
471 | end where |
---|
472 | |
---|
473 | |
---|
474 | |
---|
475 | ! si mk0=-1, on garde l'epaisseur precedente |
---|
476 | ! en attendant un masque plus propre |
---|
477 | |
---|
478 | where(mk0(:,:).eq.-1) |
---|
479 | H(:,:)=vieuxH(:,:) |
---|
480 | Hdot(:,:)=0. |
---|
481 | end where |
---|
482 | |
---|
483 | |
---|
484 | |
---|
485 | !calul de l'ablation sur les bords (pourrait n'être appelé que pour les sorties courtes) |
---|
486 | if (isynchro.eq.1) call ablation_bord |
---|
487 | |
---|
488 | |
---|
489 | if (itracebug.eq.1) call tracebug(' Fin routine icethick') !, maxval(H) ',maxval(H(:,:)) |
---|
490 | |
---|
491 | |
---|
492 | end subroutine icethick3 |
---|
493 | |
---|
494 | end module equat_adv_diff_2D_paral |
---|