source: branches/GRISLIv3/SOURCES/conserv-mass-adv-diff_sept2009_mod.f90 @ 409

Last change on this file since 409 was 389, checked in by dumas, 16 months ago

use only in module equat_adv_diff_2D_vect

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