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

Last change on this file since 334 was 285, checked in by aquiquet, 5 years ago

Small modifications for gfortran enabled compilations

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