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