source: trunk/SOURCES/dragging_neff_slope_mod.f90 @ 63

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

Function of slope for betamax in dragging_neff_slope

File size: 13.7 KB
Line 
1!> \file dragging_neff_contmaj_mod.f90
2!! Defini les zones de stream avec 3 criteres
3!<
4
5!> \namespace  dragging_neff_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 pression effective
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_neff_slope
19
20! Defini les zones de stream avec :
21! * un critere sur la pression effective
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
27use interface_input
28use io_netcdf_grisli
29
30implicit none
31logical,dimension(nx,ny) :: fleuvemx
32logical,dimension(nx,ny) :: fleuvemy
33logical,dimension(nx,ny) :: fleuve
34logical,dimension(nx,ny) :: cote
35logical,dimension(nx,ny) :: slowssa ! not actual stream, but ssa used as Ub
36logical,dimension(nx,ny) :: slowssamx
37logical,dimension(nx,ny) :: slowssamy
38
39real,dimension(nx,ny) :: hires_slope ! slope comupted on a high resolution topography
40character(len=100) :: slope_fich     ! fichier grille
41character(len=100) :: file_ncdf      !< fichier netcdf issue des fichiers .dat
42
43real :: valmax
44integer :: imax,jmax
45integer :: i_moins1,i_plus1,j_moins1,j_plus1
46integer :: lmax=20
47integer :: idep,jdep,iloc,jloc
48
49
50real :: tostick   ! pour la glace posee
51real :: betamin   ! betamin : (Pa) frottement mini sous les streams
52real :: tob_ile    ! pour les iles
53real :: cry_lim=50.  ! courbure limite pour le suivi des fleuves
54
55real,dimension(nx,ny) :: neff ! pression effective noeuds majeurs
56real :: seuil_neff    ! seuil sur la pression effective pour avoir glissement en zone sediment
57real :: coef_gz       ! coef frottement zones stream std  (10)
58real :: coef_ile      ! coef frottement zones iles   (0.1)
59
60real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
61real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
62real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
63real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
64logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
65
66contains
67!-------------------------------------------------------------------------------
68
69!> SUBROUTINE: init_dragging
70!! Cette routine fait l'initialisation du dragging.
71!>
72subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
73
74implicit none
75 
76namelist/drag_neff_slope/hwatstream,cf,betamax,betamin,toblim,tostick,seuil_neff,coef_gz,coef_ile,slope_fich
77
78if (itracebug.eq.1)  call tracebug(' dragging avec hwatermax')
79
80
81! formats pour les ecritures dans 42
82428 format(A)
83
84! lecture des parametres du run                      block drag neff
85!--------------------------------------------------------------------
86
87rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
88read(num_param,drag_neff_slope)
89
90write(num_rep_42,428)'!___________________________________________________________' 
91write(num_rep_42,428) '&drag_neff_cont              ! nom du bloc dragging neff cont'
92write(num_rep_42,*)
93write(num_rep_42,*) 'hwatstream      = ',hwatstream 
94write(num_rep_42,*) 'cf              = ',cf
95write(num_rep_42,*) 'betamax         = ', betamax
96write(num_rep_42,*) 'betamin         = ', betamin
97write(num_rep_42,*) 'toblim          = ', toblim
98write(num_rep_42,*) 'tostick          = ', tostick
99write(num_rep_42,*) 'seuil_neff      = ', seuil_neff
100write(num_rep_42,*) 'coef_gz         = ', coef_gz
101write(num_rep_42,*) 'coef_ile        = ', coef_ile
102write(num_rep_42,'(A,A)') 'slope_fich = ', slope_fich 
103write(num_rep_42,*)'/'                           
104write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream en partant de la cote'
105write(num_rep_42,428) '!  si hwater > hwatstream '
106write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
107write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
108write(num_rep_42,428) '! betamin : (Pa) frottement mini sous les streams '
109write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
110write(num_rep_42,428) '! tostick : (Pa)  pour les points non flgzmx '
111write(num_rep_42,428) '! seuil_neff (Pa) seuil sur la pression effective pour avoir glissement'
112write(num_rep_42,428) '! coef_gz : coef frottement zones stream std'
113write(num_rep_42,428) '! coef_ile : coef frottement zones iles'
114write(num_rep_42,*)
115
116!tostick=1.e5   ! valeurs pour les points non flgzmx
117tob_ile=betamax/2.
118moteurmax=toblim
119
120! betamax_2d depends on sub-grid slopes.
121slope_fich=trim(dirnameinp)//trim(slope_fich)
122call lect_input(1,'slope',1,hires_slope(:,:),slope_fich,file_ncdf)
123
124write (*,*) slope_fich
125write (*,*) hires_slope(:,:)
126write (*,*) "fin test"
127! from slopes, we create an index between 0 & 1
128! 1 is very mountainous, 0 is flat
129hires_slope(:,:)=(max(min(hires_slope(:,:),5000.),500.)-500.)/4500.
130! now we compute the actual betamax used by the remplimat routines
131! /!\ the slope is used to modify the beta where we have a temperate base,
132! but NO ice stream... -> Slow SSA zone (SSA used to compute Ub)
133betamax_2d(:,:)=max ( tostick * (1. - (1 - betamax / tostick ) * hires_slope(:,:)**(1/10.)) , betamax )
134do j=1,ny
135   do i=1,nx
136      write(18745,*) betamax_2d (i,j)
137   enddo
138enddo
139
140!-------------------------------------------------------------------
141! masque stream
142      mstream_mx(:,:)=1
143      mstream_my(:,:)=1
144
145
146! coefficient permettant de modifier le basal drag.
147      drag_mx(:,:)=1.
148      drag_my(:,:)=1.
149
150where (socle_cry(:,:).lt.cry_lim)
151   debug_3D(:,:,25)=1
152elsewhere
153   debug_3D(:,:,25)=0
154endwhere
155
156   
157      return
158      end subroutine init_dragging
159!________________________________________________________________________________
160
161!> SUBROUTINE: dragging
162!! Defini la localisation des streams et le frottement basal
163!>
164!-------------------------------------------------------------------------
165subroutine dragging   ! defini la localisation des streams et le frottement basal
166
167
168!         les iles n'ont pas de condition neff mais ont la condition toblim
169!         (ce bloc etait avant dans flottab)
170
171
172do j=2,ny
173   do i=2,nx
174      ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
175      ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
176   end do
177end do
178
179! le calcul des fleuves se fait sur les mailles majeures en partant de la cote
180
181
182
183
184! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
185
186! Partant de la cote, le fleuve est determine en remontant vers le point de hauteur d'eau max
187! tant que ce point est superieur a hwatermax.
188
189! Attention : ce systeme ne permet pas d'avoir des fleuves convergents
190!             et les fleuves n'ont une épaisseur que de 1 noeud.
191
192fleuvemx(:,:)=.false.
193fleuvemy(:,:)=.false.
194fleuve(:,:)=.false.
195cote(:,:)=.false.
196
197! detection des cotes
198do  j=2,ny-1 
199   do i=2,nx-1 
200      if ((.not.flot(i,j)).and. & 
201      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
202         cote(i,j)=.true.
203      endif
204   end do
205end do
206
207! calcul de neff (pression effective sur noeuds majeurs)
208do j=1,ny-1
209   do i=1,nx-1
210        neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4
211   enddo
212enddo
213!aurel, for the borders:
214neff(:,ny)=1.e5
215neff(nx,:)=1.e5
216
217!!$
218!!$fleuve_maj: do j=2,ny-1 
219!!$ifleuve:   do i=2,nx-1                     
220!!$
221!!$cote_detect :  if (cote(i,j)) then
222!!$         idep=i
223!!$         jdep=j
224!!$
225!!$         if (socle_cry(i,j).lt.0.) then    ! dans une vallee
226!!$            fleuve(i,j)=.true.
227!!$         else
228!!$            cote(i,j)=.false.
229!!$            cycle ifleuve
230!!$         endif
231!!$
232!!$suit : do l=1,lmax         ! debut de la boucle de suivi, lmax longueur maxi des fleuves
233!!$           i_moins1=max(idep-1,2)
234!!$           j_moins1=max(jdep-1,1)
235!!$           i_plus1=min(idep+1,nx)
236!!$           j_plus1=min(jdep+1,ny)
237!!$           
238!!$! recherche du max en suivant le socle le plus profond
239!!$! * en excluant les points flottants
240!!$! * et ceux qui sont deja tagges fleuves
241!!$
242!!$           valmax=1000.
243!!$
244!!$           do jloc=j_moins1,j_plus1
245!!$              do iloc=i_moins1,i_plus1
246!!$
247!!$                 if ((B(iloc,jloc).lt.valmax)      &
248!!$                      .and.(.not.flot(iloc,jloc))  &
249!!$                      .and.(.not.fleuve(iloc,jloc)).and.(socle_cry(iloc,jloc).lt.cry_lim)) then
250!!$                    imax=iloc
251!!$                    jmax=jloc
252!!$                    valmax=B(iloc,jloc)
253!!$                 endif
254!!$              end do
255!!$           end do
256!!$
257!!$         if ((hwater(imax,jmax).gt.hwatstream).and.(socle_cry(i,j).lt.cry_lim)) then
258!!$            fleuve(imax,jmax)=.true.
259!!$            idep=imax
260!!$            jdep=jmax
261!!$         else
262!!$            fleuve(imax,jmax)=.false.
263!!$            exit suit
264!!$         end if
265!!$
266!!$      end do suit
267!!$
268!!$   end if cote_detect
269!!$
270!!$end do ifleuve
271!!$end do fleuve_maj
272
273! aurel, we add the neff threshold:
274where ((neff(:,:).le.seuil_neff).and.(.not.flot(:,:)).and.(H(:,:).gt.1.)) fleuve(:,:)=.true.
275
276do j=1,ny-1
277   do i=1,nx-1
278      if (fleuve(i,j)) then
279         fleuvemx(i,j)=.true.
280         fleuvemx(i+1,j)=.true.
281         fleuvemy(i,j)=.true.
282         fleuvemy(i,j+1)=.true.
283      end if
284   end do
285end do
286
287! we look for the warm based points that will not be treated as stream (ub from SSA):
288slowssa(:,:)=.false.
289slowssamx(:,:)=.false.
290slowssamy(:,:)=.false.
291do  j=1,ny
292   do i=1,nx
293      !if ((not(flot(i,j))).and.(hwater(i,j).gt.0.1)) slowssa(i,j)=.true.
294      if ((not(flot(i,j))).and.(ibase(i,j).ne.1).and.(H(i,j).gt.1.)) slowssa(i,j)=.true.
295   end do
296end do
297do j=1,ny-1
298   do i=1,nx-1
299      if (slowssa(i,j)) then
300         slowssamx(i,j)=.true.
301         slowssamx(i+1,j)=.true.
302         slowssamy(i,j)=.true.
303         slowssamy(i,j+1)=.true.
304      end if
305   end do
306end do
307
308do j=1,ny
309   do i=1,nx
310
311      ! the actual streams and the warm based points are gzmx now:
312      if ( ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) .or. ((.not.ilemx(i,j)).and.(slowssamx(i,j))) ) gzmx(i,j)=.true. 
313
314     
315! calcul du frottement basal (ce bloc etait avant dans neffect)
316
317      if (cotemx(i,j)) then                        ! point cotier
318         betamx(i,j)=cf*neffmx(i,j) 
319         betamx(i,j)=min(betamx(i,j),betamax)
320
321      else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then ! tous les points SSA
322
323         if (fleuvemx(i,j)) then                 ! the actual streams
324            betamx(i,j)=cf*neffmx(i,j)*coef_gz
325            betamx(i,j)=min(betamx(i,j),betamax)
326            betamx(i,j)=max(betamx(i,j),betamin)
327
328            if (cf*neffmx(i,j).gt.betamax*2.) then  ! a stream that's becoming grounded...
329               if (slowssamx(i,j)) then
330                  betamx(i,j)=betamax_2d(i,j)
331               else
332                  gzmx(i,j)=.false.
333                  betamx(i,j)=tostick
334               endif
335            endif
336
337         else                    ! not an actual stream, SSA is used to compute Ub
338            betamx(i,j)=betamax_2d(i,j)
339         endif
340         
341      else if (ilemx(i,j)) then
342         betamx(i,j)=cf*neffmx(i,j)*coef_ile 
343         betamx(i,j)=min(betamx(i,j),tob_ile)       
344      else if (flotmx(i,j)) then     ! flottant ou ile
345         betamx(i,j)=0.
346      else                                         ! grounded, SIA
347         betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
348      endif
349
350   end do
351end do
352
353do j=1,ny
354   do i=1,nx
355
356      ! the actual streams and the warm based points are gzmx now:
357      if ( ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) .or. ((.not.ilemy(i,j)).and.(slowssamy(i,j))) ) gzmy(i,j)=.true. 
358
359     
360! calcul du frottement basal (ce bloc etait avant dans neffect)
361
362      if (cotemy(i,j)) then                        ! point cotier
363         betamy(i,j)=cf*neffmy(i,j) 
364         betamy(i,j)=min(betamy(i,j),betamax)
365
366      else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then ! tous les points SSA
367
368         if (fleuvemy(i,j)) then                 ! the actual streams
369            betamy(i,j)=cf*neffmy(i,j)*coef_gz
370            betamy(i,j)=min(betamy(i,j),betamax)
371            betamy(i,j)=max(betamy(i,j),betamin)
372
373            if (cf*neffmy(i,j).gt.betamax*2.) then  ! a stream that's becoming grounded...
374               if (slowssamy(i,j)) then
375                  betamy(i,j)=betamax_2d(i,j)
376               else
377                  gzmy(i,j)=.false.
378                  betamy(i,j)=tostick
379               endif
380            endif
381
382         else                    ! not an actual stream, SSA is used to compute Ub
383            betamy(i,j)=betamax_2d(i,j)
384         endif
385         
386      else if (ilemy(i,j)) then
387         betamy(i,j)=cf*neffmy(i,j)*coef_ile 
388         betamy(i,j)=min(betamy(i,j),tob_ile)       
389      else if (flotmy(i,j)) then     ! flottant ou ile
390         betamy(i,j)=0.
391      else                                         ! grounded, SIA
392         betamy(i,j)=tostick                       ! frottement glace posee (1 bar)
393      endif
394
395   end do
396end do
397
398
399do j=2,ny-1
400   do i=2,nx-1
401      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
402          + (betamy(i,j)+betamy(i,j+1)))*0.25
403   end do
404end do
405
406
407where (fleuve(:,:))
408   debug_3D(:,:,1)=1
409elsewhere (flot(:,:))
410   debug_3D(:,:,1)=2
411elsewhere
412   debug_3D(:,:,1)=0
413endwhere
414
415where (cote(:,:))
416   debug_3D(:,:,2)=1
417elsewhere
418   debug_3D(:,:,2)=0
419endwhere
420
421where (fleuvemx(:,:))
422   debug_3D(:,:,3)=1
423elsewhere
424   debug_3D(:,:,3)=0
425endwhere
426
427where (flotmx(:,:))
428   debug_3D(:,:,3)=-1
429endwhere
430
431!_____________________
432
433
434where (fleuvemy(:,:))
435   debug_3D(:,:,4)=1
436elsewhere
437   debug_3D(:,:,4)=0
438endwhere
439
440where (flotmy(:,:))
441   debug_3D(:,:,4)=-1
442end where
443
444debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) 
445debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) 
446
447return
448end subroutine dragging
449
450end module dragging_neff_slope
451
Note: See TracBrowser for help on using the repository browser.