source: trunk/SOURCES/Antarctique_general_files/dragging_map_allow_mod.f90 @ 4

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

initial import GRISLI trunk

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