source: branches/iLoveclim/SOURCES/conserv-mass-adv-diff_sept2009_mod.f90 @ 146

Last change on this file since 146 was 146, checked in by aquiquet, 7 years ago

Grisli-iLoveclim branch: merged to trunk at revision 145

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