source: trunk/SOURCES/dragging_hwat_contmaj_mod.f90 @ 111

Last change on this file since 111 was 66, checked in by aquiquet, 8 years ago

Dragging neff parameters in param_list && distinction between actual ice streams and SSA velocity used as Ub in SIA (slow SSA zones) in netcdf outputs

File size: 10.7 KB
Line 
1!> \file dragging_hwat_contmaj_mod.f90
2!! Defini les zones de stream avec 3 criteres
3!<
4
5!> \namespace  dragging_hwat_contmaj
6!! Defini les zones de stream avec 3 criteres
7!! \author ...
8!! \date ...
9!! @note Les trois criteres sont:
10!! @note   * un critere sur la hauteur d'eau
11!! @note   * un critere de continuite depuis la cote
12!! @note   * un critere sur la courbure du socle (si negatif, vallees)
13!! @note Used module
14!! @note   - use module3D_phy
15!! @note   - use param_phy_mod
16!<
17
18module dragging_hwat_contmaj
19
20! Defini les zones de stream avec :
21! * un critere sur la hauteur d'eau
22! * un critere de continuite depuis la cote
23! * un critere sur la courbure du socle (si negatif, vallees)
24
25use module3d_phy
26use param_phy_mod
27
28implicit none
29!logical,dimension(nx,ny) :: fleuvemx
30!logical,dimension(nx,ny) :: fleuvemy
31logical,dimension(nx,ny) :: fleuve
32logical,dimension(nx,ny) :: cote
33
34
35real :: valmax
36integer :: imax,jmax
37integer :: i_moins1,i_plus1,j_moins1,j_plus1
38integer :: lmax=20
39integer :: idep,jdep,iloc,jloc
40
41
42real :: tostick   ! pour la glace posee
43real :: tob_ile    ! pour les iles
44real :: cry_lim=50.  ! courbure limite pour le suivi des fleuves
45
46
47real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
48real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
49real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
50real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
51logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
52
53contains
54!-------------------------------------------------------------------------------
55
56!> SUBROUTINE: init_dragging
57!! Cette routine fait l'initialisation du dragging.
58!>
59subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
60
61implicit none
62 
63namelist/drag_hwat_cont/hwatstream,cf,betamax,toblim,tostick   
64
65if (itracebug.eq.1)  call tracebug(' dragging avec hwatermax')
66
67
68! formats pour les ecritures dans 42
69428 format(A)
70
71! lecture des parametres du run                      block draghwat
72!--------------------------------------------------------------------
73
74rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
75read(num_param,drag_hwat_cont)
76
77write(num_rep_42,428)'!___________________________________________________________' 
78write(num_rep_42,428) '&drag_hwat_cont              ! nom du bloc hwatermax '
79write(num_rep_42,*)
80write(num_rep_42,*) 'hwatstream      = ',hwatstream 
81write(num_rep_42,*) 'cf             = ',cf
82write(num_rep_42,*) 'betamax         = ', betamax
83write(num_rep_42,*) 'toblim         = ', toblim
84write(num_rep_42,*)'/'                           
85write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream en partant de la cote'
86write(num_rep_42,428) '!  si hwater > hwatstream '
87write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
88write(num_rep_42,428) '! seulement pour les points cotiers'
89write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
90write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
91write(num_rep_42,*)
92
93tostick=1.e5   ! valeurs pour les points non flgzmx
94tob_ile=betamax/2.
95moteurmax=toblim
96betamax_2d(:,:)=betamax
97
98!-------------------------------------------------------------------
99! masque stream
100
101      mstream_mx(:,:)=1
102      mstream_my(:,:)=1
103
104
105! coefficient permettant de modifier le basal drag.
106      drag_mx(:,:)=1.
107      drag_my(:,:)=1.
108
109where (socle_cry(:,:).lt.cry_lim)
110   debug_3D(:,:,25)=1
111elsewhere
112   debug_3D(:,:,25)=0
113endwhere
114
115   
116      return
117      end subroutine init_dragging
118!________________________________________________________________________________
119
120!> SUBROUTINE: dragging
121!! Defini la localisation des streams et le frottement basal
122!>
123!-------------------------------------------------------------------------
124subroutine dragging   ! defini la localisation des streams et le frottement basal
125
126
127!         les iles n'ont pas de condition neff mais ont la condition toblim
128!         (ce bloc etait avant dans flottab)
129
130
131do j=2,ny
132   do i=2,nx
133      ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
134      ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
135   end do
136end do
137
138! le calcul des fleuves se fait sur les mailles majeures en partant de la cote
139
140
141
142
143! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
144
145! Partant de la cote, le fleuve est determine en remontant vers le point de hauteur d'eau max
146! tant que ce point est superieur a hwatermax.
147
148! Attention : ce systeme ne permet pas d'avoir des fleuves convergents
149!             et les fleuves n'ont une épaisseur que de 1 noeud.
150
151fleuvemx(:,:)=.false.
152fleuvemy(:,:)=.false.
153fleuve(:,:)=.false.
154cote(:,:)=.false.
155
156! detection des cotes
157do  j=2,ny-1 
158   do i=2,nx-1 
159      if ((.not.flot(i,j)).and. & 
160      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
161         cote(i,j)=.true.
162      endif
163   end do
164end do
165
166fleuve_maj: do j=2,ny-1 
167ifleuve:   do i=2,nx-1                     
168
169cote_detect :  if (cote(i,j)) then
170         idep=i
171         jdep=j
172
173!         if (hwater(i,j).gt.hwatstream) then
174         if (socle_cry(i,j).lt.0.) then    ! dans une vallee
175            fleuve(i,j)=.true. 
176!           write(166,*)
177!           write(166,*)'cote',i,j,hwater(i,j),hwatstream
178         else
179            cote(i,j)=.false.
180            cycle ifleuve
181         endif
182
183suit : do l=1,lmax         ! debut de la boucle de suivi, lmax longueur maxi des fleuves
184           i_moins1=max(idep-1,2)
185           j_moins1=max(jdep-1,1)
186           i_plus1=min(idep+1,nx)
187           j_plus1=min(jdep+1,ny)
188           
189
190
191!!$! recherche du max en suivant la hauteur d'eau maxi chez les voisins
192!!$! * en excluant les points flottants
193!!$! * et ceux qui sont deja tagges fleuves   !version 18 juillet 2007
194
195!!$           valmax=0.
196!!$
197!!$         do jloc=j_moins1,j_plus1
198!!$            do iloc=i_moins1,i_plus1
199!!$
200!!$               if ((hwater(iloc,jloc).gt.valmax) &
201!!$                    .and.(.not.flot(iloc,jloc))  &
202!!$                    .and.(.not.fleuve(iloc,jloc))) then
203!!$                  imax=iloc
204!!$                  jmax=jloc
205!!$                  valmax=hwater(iloc,jloc)
206!!$               endif
207!!$            end do
208!!$         end do
209
210! recherche du max en suivant le socle le plus profond
211! * en excluant les points flottants
212! * et ceux qui sont deja tagges fleuves
213
214           valmax=1000.
215
216        do jloc=j_moins1,j_plus1
217            do iloc=i_moins1,i_plus1
218 
219               if ((B(iloc,jloc).lt.valmax)      & 
220                    .and.(.not.flot(iloc,jloc))  &
221                    .and.(.not.fleuve(iloc,jloc)).and.(socle_cry(iloc,jloc).lt.cry_lim)) then
222                  imax=iloc
223                  jmax=jloc
224                  valmax=B(iloc,jloc)
225               endif
226            end do
227         end do
228
229
230
231!         if (hwater(imax,jmax).gt.max(hwatstream-L*40.,200.)) then  ! version 18 juil
232         if ((hwater(imax,jmax).gt.hwatstream).and.(socle_cry(i,j).lt.cry_lim)) then
233            fleuve(imax,jmax)=.true.
234            idep=imax
235            jdep=jmax
236!            write(166,*) imax,jmax,hwater(imax,jmax)
237         else
238            fleuve(imax,jmax)=.false.
239            exit suit
240         end if
241 
242      end do suit
243
244   end if cote_detect
245
246end do ifleuve
247end do fleuve_maj
248
249
250! detection de problemes
251!write(166,*)
252!write(166,*) 'detection de pb dans fleuve   time=',time
253
254do j=2,ny
255   do i=2,nx
256      if (fleuve(i,j).and.(socle_cry(i,j).gt.cry_lim)) then
257         write(166,*) 'pb fleuve',time,i,j,hwater(i,j),socle_cry(i,j)
258      endif
259   end do
260end do
261
262do j=1,ny-1
263   do i=1,nx-1
264      if (fleuve(i,j)) then
265         fleuvemx(i,j)=.true.
266         fleuvemx(i+1,j)=.true.
267         fleuvemy(i,j)=.true.
268         fleuvemy(i,j+1)=.true.
269      end if
270   end do
271end do
272
273! pas de fleuve dans les endroits trop en pente
274
275fleuvemx(:,:)=fleuvemx(:,:).and.(abs(rog*Hmx(:,:)*sdx(:,:)).lt.toblim)
276fleuvemy(:,:)=fleuvemy(:,:).and.(abs(rog*Hmy(:,:)*sdy(:,:)).lt.toblim)
277
278do j=1,ny
279   do i=1,nx
280      if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. 
281
282! calcul du frottement basal (ce bloc etait avant dans neffect)
283
284      if (cotemx(i,j)) then                        ! point cotier
285         betamx(i,j)=cf*neffmx(i,j) 
286         betamx(i,j)=min(betamx(i,j),betamax)
287      else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then  ! stream
288         betamx(i,j)=betamax
289         betamx(i,j)=cf*neffmx(i,j)*20.                  ! attention a cette modif juillet 2007
290         betamx(i,j)=min(betamx(i,j),betamax)
291         betamx(i,j)=max(betamx(i,j),10.)
292
293         if (cf*neffmx(i,j).gt.1500.) then
294            gzmx(i,j)=.false.
295            betamx(i,j)=tostick
296         endif
297
298
299      else if (ilemx(i,j)) then
300         betamx(i,j)=cf*neffmx(i,j) 
301         betamx(i,j)=min(betamx(i,j),tob_ile)       
302      else if (flotmx(i,j)) then     ! flottant ou ile
303         betamx(i,j)=0.
304      else                                         ! grounded, SIA
305         betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
306      endif
307
308   end do
309end do
310
311
312do j=1,ny                  ! le noeud 1 n'est pas attribue
313   do i=1,nx   
314
315      if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. 
316
317
318 ! calcul du frottement basal (ce bloc etait avant dans neffect)
319
320      if (cotemy(i,j)) then                        ! point cotier
321         betamy(i,j)=cf*neffmy(i,j) 
322         betamy(i,j)=min(betamy(i,j),betamax)
323      else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then  ! stream
324         betamy(i,j)=betamax
325         betamy(i,j)=cf*neffmy(i,j)*20. 
326         betamy(i,j)=min(betamy(i,j),betamax)
327         betamy(i,j)=max(betamy(i,j),10.)
328
329         if (cf*neffmy(i,j).gt.1500.) then
330            gzmy(i,j)=.false.
331            betamy(i,j)=tostick
332         endif
333
334      else if (ilemy(i,j)) then
335         betamy(i,j)=cf*neffmy(i,j) 
336         betamy(i,j)=min(betamy(i,j),tob_ile)     
337      else if (flotmy(i,j)) then     ! flottant ou ile
338         betamy(i,j)=0.
339      else                               ! grounded, SIA
340         betamy(i,j)=tostick             ! frottement glace posee
341      endif
342   end do
343end do
344
345
346where (fleuve(:,:))
347   debug_3D(:,:,1)=1
348elsewhere (flot(:,:))
349   debug_3D(:,:,1)=2
350elsewhere
351   debug_3D(:,:,1)=0
352endwhere
353
354where (cote(:,:))
355   debug_3D(:,:,2)=1
356elsewhere
357   debug_3D(:,:,2)=0
358endwhere
359
360where (fleuvemx(:,:))
361   debug_3D(:,:,3)=1
362elsewhere
363   debug_3D(:,:,3)=0
364endwhere
365
366where (flotmx(:,:))
367   debug_3D(:,:,3)=-1
368endwhere
369
370!_____________________
371
372
373where (fleuvemy(:,:))
374   debug_3D(:,:,4)=1
375elsewhere
376   debug_3D(:,:,4)=0
377endwhere
378
379where (flotmy(:,:))
380   debug_3D(:,:,4)=-1
381end where
382
383debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) 
384debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) 
385
386return
387end subroutine dragging
388
389end module dragging_hwat_contmaj
390
Note: See TracBrowser for help on using the repository browser.