source: trunk/SOURCES/dragging_hwat_sedim_mod.f90 @ 170

Last change on this file since 170 was 168, checked in by dumas, 6 years ago

For compatibility update dragging modules

File size: 15.7 KB
Line 
1module dragging_hwat_sedim
2
3! Defini les zones de stream avec :
4! * un critere sur la hauteur d'eau
5! * un critere sur la courbure du socle (si negatif, vallees)
6! * un critere sur l'epaisseur de sediments (carte de Gabi & Laske)
7
8use module3d_phy
9use fake_beta_iter_vitbil_mod
10
11implicit none
12!logical,dimension(nx,ny) :: fleuvemx
13!logical,dimension(nx,ny) :: fleuvemy
14logical,dimension(nx,ny) :: fleuve
15logical,dimension(nx,ny) :: fleuve_sedim
16logical,dimension(nx,ny) :: fleuve_sedimx
17logical,dimension(nx,ny) :: fleuve_sedimy
18logical,dimension(nx,ny) :: cote
19
20
21real :: valmax
22integer :: imax,jmax
23integer :: i_moins1,i_plus1,j_moins1,j_plus1
24integer :: lmax=20
25integer :: idep,jdep,iloc,jloc
26
27
28real :: tostick   ! pour la glace posee
29real :: tob_ile    ! pour les iles
30real :: cry_lim=50.  ! courbure limite pour le suivi des fleuves
31
32real,dimension(nx,ny) :: h_sedim    ! epaisseur de sediment
33real,dimension(nx,ny) :: neff ! pression effective noeuds majeurs
34real :: seuil_hwater  ! seuil sur hwater pour avoir du glissement en zone sediment
35real :: seuil_sedim   ! seuil sur h_sedim pour avoir du glissement en zone sediment
36real :: seuil_neff    ! seuil sur la pression effective pour avoir glissement en zone sediment
37real :: coef_sedim    ! coef frottement zones sediments  (1)
38real :: coef_gz       ! coef frottement zones stream std  (10)
39real :: coef_ile      ! coef frottement zones iles   (0.1)
40real :: toto
41
42real :: betacotemax=100  ! valeur max pour betamx et betamy sur les cote
43real :: betacotemin=4
44real :: betasedmax=400    ! 400 dans c-HadCM3-evol-basal-melt-sed : cHadCM3s  et 250 dans c-HadCM3-evol-basal-melt-sed-betasedimin6
45real :: betasedmin=10   ! 10 dans c-HadCM3-evol-basal-melt-sed : cHadCM3s et 6 dans c-HadCM3-evol-basal-melt-sed-betasedimin6
46real :: betastreamax=500
47real :: betastreamin=10
48real :: betailemax=250
49real :: betailemin=6
50
51real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
52real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
53real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
54real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
55logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
56
57character(len=150) :: file_sedim   ! fichier epaisseur sediment
58
59contains
60!-------------------------------------------------------------------------------
61subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
62
63implicit none
64 
65namelist/drag_hwat_sedim/file_sedim,hwatstream,cf,betamax,toblim,seuil_hwater,seuil_sedim,seuil_neff,coef_sedim,coef_gz,coef_ile   
66
67if (itracebug.eq.1) write(num_tracebug,*)' dragging avec hwatermax'
68
69inv_beta=0
70
71! formats pour les ecritures dans 42
72428 format(A)
73
74! lecture des parametres du run                      block draghwat
75!--------------------------------------------------------------------
76
77rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
78read(num_param,drag_hwat_sedim)
79
80write(num_rep_42,428)'!___________________________________________________________' 
81write(num_rep_42,428) '&drag_hwat_sedim              ! nom du bloc hwat_sedim '
82write(num_rep_42,*)
83write(num_rep_42,*) 'file_sedim     = ', trim(file_sedim)
84write(num_rep_42,*) 'hwatstream     = ', hwatstream 
85write(num_rep_42,*) 'cf             = ', cf
86write(num_rep_42,*) 'betamax        = ', betamax
87write(num_rep_42,*) 'toblim         = ', toblim
88write(num_rep_42,*) 'seuil_hwater   = ', seuil_hwater
89write(num_rep_42,*) 'seuil_sedim    = ', seuil_sedim
90write(num_rep_42,*) 'seuil_neff     = ', seuil_neff
91write(num_rep_42,*) 'coef_sedim     = ', coef_sedim
92write(num_rep_42,*) 'coef_gz        = ', coef_gz
93write(num_rep_42,*) 'coef_ile       = ', coef_ile
94
95write(num_rep_42,*) '/'
96write(num_rep_42,428) '! file_sedim : fichier carte epaisseur sediment'                           
97write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream en partant de la cote si hwater > hwatstream '
98write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
99write(num_rep_42,428) '! seulement pour les points cotiers'
100write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
101write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
102write(num_rep_42,428) '! seuil_hwater (m) : seuil hwater pour avoir glissement sur zone sediment'
103write(num_rep_42,428) '! seuil_sedim (m) : seuil epaisseur sediment pour avoir glissement sur zone sediment'
104write(num_rep_42,428) '! seuil_neff (Pa) seuil sur la pression effective pour avoir glissement en zone sediment'
105write(num_rep_42,428) '! coef_sedim : coef frottement zones stream sediments'
106write(num_rep_42,428) '! coef_gz : coef frottement zones stream std'
107write(num_rep_42,428) '! coef_ile : coef frottement zones iles'
108write(num_rep_42,*)
109
110tostick=1.e5   ! valeurs pour les points non flgzmx
111tob_ile=betamax/2.
112moteurmax=toblim
113
114!-------------------------------------------------------------------
115! masque stream
116
117      mstream_mx(:,:)=1
118      mstream_my(:,:)=1
119
120
121! coefficient permettant de modifier le basal drag.
122      drag_mx(:,:)=1.
123      drag_my(:,:)=1.
124
125where (socle_cry(:,:).lt.cry_lim)
126   debug_3D(:,:,25)=1
127elsewhere
128   debug_3D(:,:,25)=0
129endwhere
130
131! lecture du fichier sediment :
132!open(20,file=TRIM(DIRNAMEINP)//trim(file_sedim)) !   'sediment_ij_hemin40.dat')
133!do j=1,ny
134!   do i=1,nx
135!      read(20,*),toto,toto,h_sedim(i,j)
136!   enddo
137!enddo
138!h_sedim(:,:)=max(h_sedim(:,:),0.)  ! pas d'epaisseurs de sediments negatives
139h_sedim(:,:)=10. 
140 
141      return
142      end subroutine init_dragging
143!________________________________________________________________________________
144
145
146!-------------------------------------------------------------------------
147subroutine dragging   ! defini la localisation des streams et le frottement basal
148
149
150!         les iles n'ont pas de condition neff mais ont la condition toblim
151!         (ce bloc etait avant dans flottab)
152
153
154do j=2,ny
155   do i=2,nx
156      ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
157      ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
158   end do
159end do
160
161! le calcul des fleuves se fait sur les mailles majeures en partant de la cote
162
163
164
165
166! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
167
168! Partant de la cote, le fleuve est determine en remontant vers le point de hauteur d'eau max
169! tant que ce point est superieur a hwatermax.
170
171! Attention : ce systeme ne permet pas d'avoir des fleuves convergents
172!             et les fleuves n'ont une épaisseur que de 1 noeud.
173
174fleuvemx(:,:)=.false.
175fleuvemy(:,:)=.false.
176fleuve(:,:)=.false.
177cote(:,:)=.false.
178
179! detection des cotes
180do  j=2,ny-1 
181   do i=2,nx-1 
182      if ((.not.flot(i,j)).and. & 
183      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
184         cote(i,j)=.true.
185      endif
186   end do
187end do
188
189! calcul de neff (pression effective sur noeuds majeurs)
190do j=1,ny-1
191   do i=1,nx-1
192        neff(i,j)=(neffmx(i,j)+neffmx(i,j+1)+neffmy(i,j)+neffmy(i+1,j))/4
193   enddo
194enddo
195
196! test sur hauteur d'eau et epaisseur sediment :
197!where (h_sedim(:,:).GT.seuil_sedim.and.hwater(:,:).GT.seuil_hwater.and.(.not.flot(:,:))) fleuve(:,:)=.true.
198where (h_sedim(:,:).GT.seuil_sedim.and.neff(:,:).LE.seuil_neff.and.hwater(:,:).GT.seuil_hwater.and.(.not.flot(:,:))) fleuve(:,:)=.true.
199
200fleuve_sedim(:,:)=fleuve(:,:)  ! points uniquement fleuves sediment
201
202
203fleuve_maj: do j=2,ny-1 
204ifleuve:   do i=2,nx-1                     
205
206cote_detect :  if (cote(i,j)) then
207   idep=i
208   jdep=j
209         
210      if (socle_cry(i,j).lt.0.) then    ! dans une vallee
211         fleuve(i,j)=.true. 
212!         write(166,*)
213!         write(166,*)'cote',i,j,hwater(i,j),hwatstream
214      else
215         cote(i,j)=.false.
216         cycle ifleuve
217      endif
218
219suit : do l=1,lmax         ! debut de la boucle de suivi, lmax longueur maxi des fleuves
220           i_moins1=max(idep-1,2)
221           j_moins1=max(jdep-1,1)
222           i_plus1=min(idep+1,nx)
223           j_plus1=min(jdep+1,ny)
224           
225
226
227!!$! recherche du max en suivant la hauteur d'eau maxi chez les voisins
228!!$! * en excluant les points flottants
229!!$! * et ceux qui sont deja tagges fleuves   !version 18 juillet 2007
230
231!!$           valmax=0.
232!!$
233!!$         do jloc=j_moins1,j_plus1
234!!$            do iloc=i_moins1,i_plus1
235!!$
236!!$               if ((hwater(iloc,jloc).gt.valmax) &
237!!$                    .and.(.not.flot(iloc,jloc))  &
238!!$                    .and.(.not.fleuve(iloc,jloc))) then
239!!$                  imax=iloc
240!!$                  jmax=jloc
241!!$                  valmax=hwater(iloc,jloc)
242!!$               endif
243!!$            end do
244!!$         end do
245
246! recherche du max en suivant le socle le plus profond
247! * en excluant les points flottants
248! * et ceux qui sont deja tagges fleuves
249
250           valmax=1000.
251
252        do jloc=j_moins1,j_plus1
253            do iloc=i_moins1,i_plus1
254 
255               if ((B(iloc,jloc).lt.valmax)      & 
256                    .and.(.not.flot(iloc,jloc))  &
257                    .and.(.not.fleuve(iloc,jloc)).and.(socle_cry(iloc,jloc).lt.cry_lim)) then
258                  imax=iloc
259                  jmax=jloc
260                 valmax=B(iloc,jloc)
261               endif
262            end do
263         end do
264
265
266
267
268         if ((hwater(imax,jmax).gt.hwatstream).and.(socle_cry(i,j).lt.cry_lim)) then
269            fleuve(imax,jmax)=.true.
270            idep=imax
271            jdep=jmax
272!            write(166,*) imax,jmax,hwater(imax,jmax)
273         else
274!            fleuve(imax,jmax)=.false.
275            exit suit
276         end if
277 
278      end do suit
279
280   end if cote_detect
281
282end do ifleuve
283end do fleuve_maj
284
285
286! detection de problemes
287!write(166,*)
288!write(166,*) 'detection de pb dans fleuve   time=',time
289!do j=2,ny
290!   do i=2,ny
291!      if (fleuve(i,j).and.(socle_cry(i,j).gt.cry_lim)) then
292!         write(166,*) 'pb fleuve',time,i,j,hwater(i,j),socle_cry(i,j)
293!      endif
294!   end do
295!end do
296
297do j=1,ny-1
298   do i=1,nx-1
299      if (fleuve(i,j)) then
300         fleuvemx(i,j)=.true.
301         fleuvemx(i+1,j)=.true.
302         fleuvemy(i,j)=.true.
303         fleuvemy(i,j+1)=.true.
304      end if
305      if (fleuve_sedim(i,j)) then
306         fleuve_sedimx(i,j)=.true.
307         fleuve_sedimx(i+1,j)=.true.
308         fleuve_sedimy(i,j)=.true.
309         fleuve_sedimy(i,j+1)=.true.
310      endif
311   end do
312end do
313
314! pas de fleuve dans les endroits trop en pente
315!fleuvemx(:,:)=fleuvemx(:,:).and.(abs(rog*Hmx(:,:)*sdx(:,:)).lt.toblim)
316!fleuvemy(:,:)=fleuvemy(:,:).and.(abs(rog*Hmy(:,:)*sdy(:,:)).lt.toblim)
317!fleuvemx(:,:)=fleuvemx(:,:).and.(abs(rog*Hmx(:,:)*sdx(:,:)).lt.toblim)
318!fleuvemy(:,:)=fleuvemy(:,:).and.(abs(rog*Hmy(:,:)*sdy(:,:)).lt.toblim)
319
320
321do j=1,ny
322   do i=1,nx
323      if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. 
324
325! calcul du frottement basal (ce bloc etait avant dans neffect)
326
327      if (cotemx(i,j)) then                        ! point cotier
328         betamx(i,j)=cf*neffmx(i,j) 
329! cat         betamx(i,j)=min(betamx(i,j),tobmax)
330! version Jorge
331! point cotier: frott maxi=40, mini=4
332!         betamx(i,j)=min(betamx(i,j),betacotemax)
333!         betamx(i,j)=max(betamx(i,j),betacotemin)
334! version Jorge stream  en zone sediment :
335      else if (fleuve_sedimx(i,j).and.gzmx(i,j).and.(.not.cotemx(i,j))) then  ! stream  en zone sediment
336! stream sediment: frott maxi=60, mini=6
337!cdc            betamx(i,j)=cf*neffmx(i,j)*3
338             betamx(i,j)=cf*neffmx(i,j)*coef_sedim
339!            betamx(i,j)=min(betamx(i,j),betasedmax)         
340!            betamx(i,j)=max(betamx(i,j),betasedmin)
341      else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then  ! stream normal
342! cat         betamx(i,j)=tobmax
343! cat         betamx(i,j)=cf*neffmx(i,j)*20.                  ! attention a cette modif juillet 2007
344! cat         betamx(i,j)=min(betamx(i,j),tobmax)
345! cat         betamx(i,j)=max(betamx(i,j),10.)
346! version Jorge
347! stream normal: frott maxi=80, mini=8
348         betamx(i,j)=cf*neffmx(i,j)*coef_gz         
349!         betamx(i,j)=min(betamx(i,j),betastreamax)           
350!         betamx(i,j)=max(betamx(i,j),betastreamin)
351
352! cat         if (cf*neffmx(i,j).gt.1500.) then
353! cat          gzmx(i,j)=.false.
354! cat           betamx(i,j)=tostick
355! cat        endif
356
357      else if (ilemx(i,j)) then
358         betamx(i,j)=cf*neffmx(i,j)*coef_ile
359! cat         betamx(i,j)=min(betamx(i,j),tob_ile)   
360! version Jorge
361! ile : frott maxi=20, mini=2   
362!         betamx(i,j)=min(betamx(i,j),betailemax)
363!         betamx(i,j)=max(betamx(i,j),betailemin)
364      else if (flotmx(i,j)) then     ! flottant ou ile
365         betamx(i,j)=0.
366      else                                         ! grounded, SIA
367         betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
368      endif
369   end do
370end do
371
372
373do j=1,ny                  ! le noeud 1 n'est pas attribue
374   do i=1,nx   
375
376      if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. 
377
378
379 ! calcul du frottement basal (ce bloc etait avant dans neffect)
380
381      if (cotemy(i,j)) then                        ! point cotier
382         betamy(i,j)=cf*neffmy(i,j) 
383! cat         betamy(i,j)=min(betamy(i,j),tobmax)
384! version Jorge
385! point cotier: frott maxi=40, mini=4
386!         betamy(i,j)=min(betamy(i,j),betacotemax)
387!         betamy(i,j)=max(betamy(i,j),betacotemin)
388! version Jorge stream  en zone sediment :
389      else if ((fleuve_sedimy(i,j).and.gzmx(i,j).and.(.not.cotemx(i,j)))) then  ! stream  en zone sediment
390! stream sediment: frott maxi=60, mini=6
391!cdc            betamy(i,j)=cf*neffmy(i,j)*3
392             betamy(i,j)=cf*neffmy(i,j)*coef_sedim
393!            betamy(i,j)=min(betamy(i,j),betasedmax)         
394!            betamy(i,j)=max(betamy(i,j),betasedmin)
395      else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then  ! stream normal
396! cat         betamy(i,j)=tobmax
397! cat         betamy(i,j)=cf*neffmy(i,j)*20.
398! cat         betamy(i,j)=min(betamy(i,j),tobmax)
399! cat         betamy(i,j)=max(betamy(i,j),10.)
400! version Jorge
401! stream normal: frott maxi=80, mini=8
402         betamy(i,j)=cf*neffmy(i,j)*coef_gz                     
403!         betamy(i,j)=min(betamy(i,j),betastreamax)           
404!         betamy(i,j)=max(betamy(i,j),betastreamin)
405
406! cat         if (cf*neffmy(i,j).gt.1500.) then
407! cat            gzmy(i,j)=.false.
408! cat            betamy(i,j)=tostick
409! cat         endif
410
411      else if (ilemy(i,j)) then
412         betamy(i,j)=cf*neffmy(i,j)*coef_ile
413! cat         betamy(i,j)=min(betamy(i,j),tob_ile)     
414! version Jorge
415! ile : frott maxi=20, mini=2   
416!         betamy(i,j)=min(betamy(i,j),betailemax)
417!         betamy(i,j)=max(betamy(i,j),betailemin)
418      else if (flotmy(i,j)) then     ! flottant ou ile
419         betamy(i,j)=0.
420      else                               ! grounded, SIA
421         betamy(i,j)=tostick             ! frottement glace posee
422      endif
423   end do
424end do
425
426!print*,'dragging debug h_sedim(78,144),hwater,seuil_sedim,seuil_hwater',h_sedim(78,144),seuil_sedim,hwater(78,144),seuil_hwater,fleuve(78,144),gzmx(78,144)
427!print*,'dragging debug h_sedim(90,157),hwater,seuil_sedim,seuil_hwater',h_sedim(90,157),seuil_sedim,hwater(90,157),seuil_hwater,fleuve(90,157),gzmx(90,157),betamx(90,157),betamx(78,144)
428
429
430where (fleuve(:,:))
431   debug_3D(:,:,1)=1
432elsewhere (flot(:,:))
433   debug_3D(:,:,1)=2
434elsewhere
435   debug_3D(:,:,1)=0
436endwhere
437
438where (cote(:,:))
439   debug_3D(:,:,2)=1
440elsewhere
441   debug_3D(:,:,2)=0
442endwhere
443
444where (fleuvemx(:,:))
445   debug_3D(:,:,3)=1
446elsewhere
447   debug_3D(:,:,3)=0
448endwhere
449
450where (flotmx(:,:))
451   debug_3D(:,:,3)=-1
452endwhere
453
454!_____________________
455
456
457where (fleuvemy(:,:))
458   debug_3D(:,:,4)=1
459elsewhere
460   debug_3D(:,:,4)=0
461endwhere
462
463where (flotmy(:,:))
464   debug_3D(:,:,4)=-1
465end where
466
467debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) 
468debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) 
469
470
471!print*,'point (79,145)',hwater(79,145),h_sedim(79,145),flot(79,145),fleuve_sedim(79,145),fleuve(79,145),betamx(79,145),fleuve_sedimx(79,145),gzmx(79,145)
472
473return
474end subroutine dragging
475
476end module dragging_hwat_sedim
477
Note: See TracBrowser for help on using the repository browser.