source: trunk/SOURCES/Antarctique_general_files/dragging_map_Ant_mod.f90

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

initial import GRISLI trunk

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