source: trunk/SOURCES/conserv-mass-adv-diff_sept2009_mod.f90 @ 225

Last change on this file since 225 was 224, checked in by dumas, 6 years ago

format of the parameter file created by the code, possibility to use it in reading

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