source: trunk/SOURCES/Hudson_files/dragging_hudson_jorge_mod_old.f90 @ 334

Last change on this file since 334 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

File size: 6.6 KB
Line 
1!> \file dragging_hudson_jorge_mod_old.f90
2!! Module qui definie les zones de stream
3!<
4
5!> \namespace  dragging_hudson_jorge
6!! Definie les zones de stream avec uniquement:
7!! @note - un critere sur la hauteur d'eau
8!! @note jalv: module inspire par le dragging_hwat_contmaj mais simplifie
9!! \author ...
10!! \date ...
11!! @note Used module
12!! @note   - use module3D_phy
13!! @note   - use sedim_declar
14!<
15
16
17module dragging_hudson_jorge
18
19!jalv: module inspire par le dragging_hwat_contmaj mais simplifie
20! Definie les zones de stream avec uniquement:
21! * un critere sur la hauteur d'eau
22 
23
24  use module3d_phy
25  use sedim_declar
26
27implicit none
28logical,dimension(nx,ny) :: fleuvemx
29logical,dimension(nx,ny) :: fleuvemy
30logical,dimension(nx,ny) :: fleuve
31logical,dimension(nx,ny) :: cote
32
33real :: seuil_sedim              !< seuil sur hwater pour avoir du glissement
34real :: valmax
35integer :: imax,jmax
36integer :: i_moins1,i_plus1,j_moins1,j_plus1
37integer :: lmax=20
38integer :: idep,jdep,iloc,jloc
39
40integer :: itestd
41
42
43real :: tostick   !< pour la glace posee
44real :: tob_ile    !< pour les iles
45real :: cry_lim=50.  !< courbure limite pour le suivi des fleuves
46contains
47!-------------------------------------------------------------------------------
48
49!> SUBROUTINE: init_dragging
50!! Cette routine fait l'initialisation du dragging.
51!<
52subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
53
54implicit none
55 
56namelist/drag_hudson_jorge/hwatstream,cf,betamax,toblim,tostick,seuil_sedim   
57
58if (itracebug.eq.1)  call tracebug(' Dragging avec hwatermax')
59
60
61! formats pour les ecritures dans 42
62428 format(A)
63
64! lecture des parametres du run                      block dra_hudson_jorge
65!--------------------------------------------------------------------
66
67rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
68read(num_param,drag_hudson_jorge)
69
70write(num_rep_42,428)'!___________________________________________________________' 
71write(num_rep_42,428) '&drag_hudson_jorge              ! nom du bloc drag_hudson_jorge '
72write(num_rep_42,*)
73write(num_rep_42,*) 'hwatstream      = ',hwatstream 
74write(num_rep_42,*) 'cf             = ',cf
75write(num_rep_42,*) 'betamax         = ', betamax
76write(num_rep_42,*) 'toblim         = ', toblim
77write(num_rep_42,*) 'seuil_sedim         = ', seuil_sedim
78write(num_rep_42,*)'/'                           
79write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream en partant de la cote'
80write(num_rep_42,428) '!  si hwater > hwatstream '
81write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
82write(num_rep_42,428) '! seulement pour les points cotiers'
83write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
84write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
85write(num_rep_42,*)
86
87tostick=1.e5   ! valeurs pour les points non flgzmx
88tob_ile=betamax/2.
89moteurmax=toblim
90
91!-------------------------------------------------------------------
92! masque stream
93
94      mstream_mx(:,:)=1
95      mstream_my(:,:)=1
96
97! coefficient permettant de modifier le basal drag.
98      drag_mx(:,:)=1.
99      drag_my(:,:)=1.
100
101return
102end subroutine init_dragging
103!________________________________________________________________________________
104
105
106
107!-------------------------------------------------------------------------
108
109!> SUBROUTINE: dragging
110!! Defini la localisation des streams et le frottement basal
111!<
112subroutine dragging   ! defini la localisation des streams et le frottement basal
113
114
115fleuvemx(:,:)=.false.
116fleuvemy(:,:)=.false.
117fleuve(:,:)=.false.
118ilemx(:,:)=.false.
119ilemy(:,:)=.false.
120
121
122do j=1,ny
123   do i=1,nx
124      if ((hwater(i,j).gt.hwatstream).or.((mksedim(i,j).eq.2.).and.(hwater(i,j).ge.seuil_sedim))) then
125         fleuve(i,j)=.true.
126         fleuvemx(i,j)=.true.                  !ajoute par jalv
127         fleuvemy(i,j)=.true.                  !ajoute par jalv
128      else
129         fleuve(i,j)=.false.
130         fleuvemx(i,j)=.false.                 !ajoute par jalv
131         fleuvemy(i,j)=.false.                 !ajoute par jalv
132      end if
133   end do
134end do
135
136!jalv
137!call detect_assym(nx,ny,0,76,1,0,1,0,betamx,itestd)
138!if (itestd.gt.0) then
139!   write(6,*) 'avant calcul betax asymetrie sur betamx  pour time=',time
140!   stop
141!else
142!   write(6,*) 'avant calcul betax pas d asymetrie sur betamx  pour time=',time
143!end if
144
145!call detect_assym(nx,ny,0,76,1,0,1,0,betamy,itestd)
146!if (itestd.gt.0) then
147!    write(6,*) 'avant calcul betay asymetrie sur betamy  pour time=',time
148!    stop
149!else
150!    write(6,*) 'avant calcul betay pas d asymetrie sur betamy  pour time=',time
151!end if
152
153! pas de fleuve dans les endroits trop en pente
154!test jalv: j'autorise les fleuves dans les endroits trop en pente
155
156!fleuvemx(:,:)=fleuvemx(:,:).and.(abs(rog*Hmx(:,:)*sdx(:,:)).lt.toblim)
157!fleuvemy(:,:)=fleuvemy(:,:).and.(abs(rog*Hmy(:,:)*sdy(:,:)).lt.toblim)
158
159! calcul du frottement basal (ce bloc etait avant dans neffect)
160
161do j=1,ny         
162   do i=1,nx
163
164      if (fleuvemy(i,j)) then
165         betamy(i,j)=betamax
166!         betamy(i,j)=cf*neffmy(i,j)*20.
167         betamy(i,j)=15. 
168         betamy(i,j)=min(betamy(i,j),betamax)
169         betamy(i,j)=max(betamy(i,j),20.)
170      endif
171
172 !     if (cf*neffmy(i,j).gt.1500.) then
173 !        fleuvemy(i,j)=.false.
174!      endif
175
176      if (ilemy(i,j)) then
177!         betamy(i,j)=cf*neffmy(i,j)
178         betamy(i,j)=10.
179         betamy(i,j)=min(betamy(i,j),tob_ile)       
180      endif
181
182      if (flotmy(i,j)) then     ! flottant
183         betamy(i,j)=0.
184      endif
185
186   end do
187end do
188
189
190do j=1,ny           ! le noeud 1 n'est pas attribue
191   do i=1,nx
192      if (fleuvemx(i,j)) then
193         betamx(i,j)=betamax
194!         betamx(i,j)=cf*neffmx(i,j)*20.
195         betamx(i,j)=15.         
196         betamx(i,j)=min(betamx(i,j),betamax)
197         betamx(i,j)=max(betamx(i,j),20.)
198      endif
199
200!      if (cf*neffmx(i,j).gt.1500.) then
201!         fleuvemx(i,j)=.false.
202!      endif
203
204      if (ilemx(i,j)) then
205!         betamx(i,j)=cf*neffmx(i,j)
206         betamx(i,j)=10.
207         betamx(i,j)=min(betamx(i,j),tob_ile)       
208      endif
209
210      if (flotmx(i,j)) then     ! flottant
211         betamx(i,j)=0.
212      endif
213
214   end do
215end do
216
217
218!jalv
219!call detect_assym(nx,ny,0,76,1,0,1,0,betamx,itestd)
220!if (itestd.gt.0) then
221!   write(6,*) 'apres dragging asymetrie sur betamx  pour time=',time
222!   stop
223!else
224!   write(6,*) 'apres dragging pas d asymetrie sur betamx  pour time=',time
225!end if
226
227!call detect_assym(nx,ny,0,76,1,0,1,0,betamy,itestd)
228!if (itestd.gt.0) then
229!    write(6,*) 'apres dragging asymetrie sur betamy  pour time=',time
230!    stop
231!else
232!    write(6,*) 'apres dragging pas d asymetrie sur betamy  pour time=',time
233!end if
234
235return
236
237end subroutine dragging
238
239end module dragging_hudson_jorge
240
Note: See TracBrowser for help on using the repository browser.