source: trunk/SOURCES/Greeneem_files/dragging_neem_mod.f90 @ 21

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

initial import GRISLI trunk

File size: 8.5 KB
Line 
1!> \file dragging-vit_bil_CISM_gen_neem_mod.f90
2!! Definition des zones de streams apres lecture d'une carte evaluee offline
3!<
4
5!> \namespace dragging_neem
6!! Definition zone stream. Input neem : critere sur vitesse de bilan et sur la 'courbure'
7!! \author ...
8!! \date ...
9!! @note Defini les zones de stream avec :
10!! @note   - un critere sur la hauteur d'eau
11!! @note   - un critere de fleuves actuels
12!! @note Attention il faut fonctionner en additionnant SIA et Stream pour éviter les pb
13!! de fleuves en pointille
14!! @note Used module
15!! @note   - use module3D_phy
16!! @note   - use param_phy_mod
17!<
18
19
20module dragging_neem
21
22! Defini les zones de stream avec :
23! * un critere sur la hauteur d'eau
24! * un critere de fleuves actuels
25! attention il faut fonctionner en additionnant SIA et Stream pour éviter les pb
26! de fleuves en pointille
27
28
29use module3d_phy
30use param_phy_mod
31
32implicit none
33logical,dimension(nx,ny) :: fleuvemx         !< fleuves sont les tableaux courants (dep. time)   
34logical,dimension(nx,ny) :: fleuvemy
35logical,dimension(nx,ny) :: fleuve     
36logical,dimension(nx,ny) :: cote   
37real, dimension(nx,ny) :: Vcol_x              !< uniquement pour compatibilite dans le spinup cat (aurel)
38real, dimension(nx,ny) :: Vcol_y              !< uniquement pour compatibilite dans le spinup cat (aurel)
39real, dimension(nx,ny) :: Vsl_x            !< pour compatibilite dragging beta (aurel)
40real, dimension(nx,ny) :: Vsl_y            !< pour compatibilite dragging beta (aurel)
41
42logical :: corr_def = .false.               !< for deformation correction, pour compatibilite beta (aurel)
43
44character(len=100) :: masque_stream        !<  file masque streams autorises
45
46real :: valmax
47integer :: imax,jmax
48integer :: i_moins1,i_plus1,j_moins1,j_plus1
49integer :: lmax=20
50integer :: idep,jdep,iloc,jloc
51
52
53real :: tostick      ! pour la glace posee
54real :: tob_ile      ! pour les iles
55
56contains
57!-------------------------------------------------------------------------------
58
59
60!> SUBROUTINE: init_dragging
61!! Cette routine fait l'initialisation du dragging.
62!<   
63
64subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
65
66implicit none
67
68integer :: bidon
69
70 
71namelist/drag_neem/hwatstream,cf,betamax,toblim,masque_stream
72
73if (itracebug.eq.1) write(num_tracebug,*)' init_dragging avec dragging_neem'
74
75
76! formats pour les ecritures dans 42
77428 format(A)
78
79! lecture des parametres du run                   
80!--------------------------------------------------------------------
81
82rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
83read(num_param,drag_neem)
84
85write(num_rep_42,428) '!___________________________________________________________' 
86write(num_rep_42,428) '! module dragging_neem '
87write(num_rep_42,drag_neem)                   
88write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream vitesses de bilan'
89write(num_rep_42,428) '!  si hwater > hwatstream '
90write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
91write(num_rep_42,428) '! seulement pour les points cotiers'
92write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
93write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
94write(num_rep_42,428) '! masque_stream : nom du fichier qui contient le masque de streams autorises'
95write(num_rep_42,*)
96
97tostick=1.e5   ! valeurs pour les points non flgzmx
98tob_ile=betamax/2.
99moteurmax=toblim
100!moteurmax=0.
101
102!-------------------------------------------------------------------
103! masque stream : lecture du masque
104! attention on prend la colonne 1 qui correspond aux noeuds majeurs
105masque_stream  = trim(dirnameinp)//trim(masque_stream)
106print*, 'le fichier des streams ', masque_stream
107call lect_datfile(nx,ny,mstream,1,masque_stream)     
108open (20,file=masque_stream)
109read(20,*)
110do j=1,ny
111   do i=1,nx
112      read(20,*) mstream(i,j),bidon, bidon! mstream_mx(i:i+1,j), mstream_my(i,j:j+1)
113   end do
114end do
115close(20)
116
117do j=1,ny
118   do i=1,nx
119      if (mstream(i,j).eq.1) then ! allowed ice streams
120         mstream_mx(i:i+1,j) = 1
121         mstream_my(i,j:j+1) = 1
122      end if
123   end do
124end do
125
126
127debug_3D(:,:,25) = mstream(:,:)
128
129
130! coefficient permettant de modifier le basal drag.
131      drag_mx(:,:)=1.
132      drag_my(:,:)=1.
133
134
135return
136end subroutine init_dragging
137!________________________________________________________________________________
138
139!> SUBROUTINE: dragging
140!! defini la localisation des streams et le frottement basal
141!<
142
143!-------------------------------------------------------------------------
144subroutine dragging   ! defini la localisation des streams et le frottement basal
145
146
147!         les iles n'ont pas de condition neff mais ont la condition toblim
148!         (ce bloc etait avant dans flottab)
149
150
151do j=2,ny
152   do i=2,nx
153      ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
154      ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
155   end do
156end do
157
158fleuvemx(:,:)=.false.
159fleuvemy(:,:)=.false.
160fleuve(:,:)=.false.
161cote(:,:)=.false.
162cotemx(:,:)=.false.
163cotemy(:,:)=.false.
164
165! detection des cotes
166do  j=2,ny-1 
167   do i=2,nx-1 
168      if ((.not.flot(i,j)).and. & 
169      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
170         cote(i,j)=.true.
171         cotemx(i:i+1,j) = .true.
172         cotemy(i,j:j+1) = .true.
173      endif
174   end do
175end do
176
177! le calcul des fleuves se fait sur les mailles majeures
178! normalement, les points cotiers sont deja tagges gzmx dans la routine flottab
179
180do j=2,ny-1
181   do i=2,nx-1
182     
183      if ((hwater(i,j).gt.hwatstream).and.(mstream(i,j).eq.1).and.(.not.flot(i,j))) then
184         fleuve(i,j)=.true.
185         fleuvemx(i:i+1,j)=.true.
186         fleuvemy(i,j:j+1)=.true.
187         gzmx(i:i+1,j)=.true.
188         gzmy(i,j:j+1)=.true.
189      endif
190   end do
191end do
192
193
194! calcul du frottement basal (ce bloc etait avant dans neffect)
195
196drag_mx : do j=2,ny-1
197   do i=2,nx-1
198
199      if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. 
200
201      if (cotemx(i,j)) then                        ! point cotier peut etre diminue si pression d'eau
202         betamx(i,j)=cf*neffmx(i,j) 
203         betamx(i,j)=min(betamx(i,j),betamax)
204      else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then  ! stream -> betamax
205         betamx(i,j)=betamax
206         betamx(i,j)=min(betamx(i,j),betamax)
207         betamx(i,j)=max(betamx(i,j),10.)
208
209!         if (cf*neffmx(i,j).gt.1500.) then
210!            gzmx(i,j)=.false.
211!            betamx(i,j)=tostick
212!         endif
213
214      else if (ilemx(i,j)) then                     ! ile : tob_ile=betamax/2 + peut etre diminue si pression d'eau
215         betamx(i,j)=cf*neffmx(i,j) 
216         betamx(i,j)=min(betamx(i,j),tob_ile)       
217      else if (flotmx(i,j)) then                    ! flottant ou ile
218         betamx(i,j)=0.
219      else                                         ! grounded, SIA
220         betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
221      endif
222
223   end do
224end do drag_mx
225
226drag_my : do j=2,ny-1
227   do i=2,nx-1
228
229      if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. 
230
231      if (cotemy(i,j)) then                        ! point cotier peut etre diminue si pression d'eau
232         betamy(i,j)=cf*neffmy(i,j) 
233         betamy(i,j)=min(betamy(i,j),betamax)
234      else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then  ! stream -> betamax
235         betamy(i,j)=betamax
236         betamy(i,j)=min(betamy(i,j),betamax)
237         betamy(i,j)=max(betamy(i,j),10.)
238
239!         if (cf*neffmy(i,j).gt.1500.) then
240!            gzmy(i,j)=.false.
241!            betamy(i,j)=tostick
242!         endif
243
244      else if (ilemy(i,j)) then                     ! ile : tob_ile=betamax/2 + peut etre diminue si pression d'eau
245         betamy(i,j)=cf*neffmy(i,j) 
246         betamy(i,j)=min(betamy(i,j),tob_ile)       
247      else if (flotmy(i,j)) then                    ! flottant ou ile
248         betamy(i,j)=0.
249      else                                         ! grounded, SIA
250         betamy(i,j)=tostick                       ! frottement glace posee (1 bar)
251      endif
252
253   end do
254end do drag_my
255
256
257where (fleuve(:,:))
258   debug_3D(:,:,1)=1
259elsewhere (flot(:,:))
260   debug_3D(:,:,1)=2
261elsewhere
262   debug_3D(:,:,1)=0
263endwhere
264
265where (cote(:,:))
266   debug_3D(:,:,2)=1
267elsewhere
268   debug_3D(:,:,2)=0
269endwhere
270
271where (fleuvemx(:,:))
272   debug_3D(:,:,3)=1
273elsewhere
274   debug_3D(:,:,3)=0
275endwhere
276
277where (flotmx(:,:))
278   debug_3D(:,:,3)=-1
279endwhere
280
281!_____________________
282
283
284where (fleuvemy(:,:))
285   debug_3D(:,:,4)=1
286elsewhere
287   debug_3D(:,:,4)=0
288endwhere
289
290where (flotmy(:,:))
291   debug_3D(:,:,4)=-1
292end where
293
294debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) 
295debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) 
296
297return
298end subroutine dragging
299
300end module dragging_neem
301
Note: See TracBrowser for help on using the repository browser.