source: trunk/SOURCES/dragging_neff_slope_mod.f90 @ 84

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

Small bug correction in dragging modules for gfortran compatibility

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