source: branches/iLoveclim/SOURCES/Antarctique_general_files/dragging_LGM_mod.f90 @ 90

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

initial import GRISLI trunk

File size: 8.9 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_LGM
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
31use interface_input
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)
42real, dimension(nx,ny) :: work_tab         ! pour lire les masques
43
44logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta (aurel)
45
46character(len=100) :: masque_stream        !<  file masque streams autorises
47character(len=100) :: file_ncdf            !< fichier netcdf issue des fichiers .dat
48
49real    :: Tstream_lim                     !< limit for streaming (lower than melting point to account for spatial variability)
50real    :: beta_prescrit                    !< valeur pour les ice streams
51real    :: valmax
52integer :: imax,jmax
53integer :: i_moins1,i_plus1,j_moins1,j_plus1
54integer :: lmax=20
55integer :: idep,jdep,iloc,jloc
56
57
58real :: tostick      ! pour la glace posee
59real :: tob_ile      ! pour les iles
60
61contains
62!-------------------------------------------------------------------------------
63
64
65!> SUBROUTINE: init_dragging
66!! Cette routine fait l'initialisation du dragging.
67!<   
68
69subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
70
71implicit none
72
73integer :: bidon
74
75 
76namelist/drag_LGM/Tstream_lim,cf,beta_prescrit,betamax,masque_stream
77
78if (itracebug.eq.1) write(num_tracebug,*)' init_dragging avec dragging_LGM'
79
80
81! formats pour les ecritures dans 42
82428 format(A)
83
84! lecture des parametres du run                   
85!--------------------------------------------------------------------
86
87rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
88read(num_param,drag_LGM)
89
90write(num_rep_42,428) '!___________________________________________________________' 
91write(num_rep_42,428) '! module dragging_LGM '
92write(num_rep_42,drag_LGM)                   
93write(num_rep_42,428) '! critere de temperature pour passage en stream'
94write(num_rep_42,428) '! si T > Tlim'
95write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
96write(num_rep_42,428) '! seulement pour les points cotiers'
97write(num_rep_42,428) '! beta_prescrit : parametre de frottement maxi sous les streams '
98write(num_rep_42,428) '! masque_stream : nom du fichier qui contient le masque de streams autorises'
99write(num_rep_42,*)
100
101tostick=betamax   ! valeurs pour les points non flgzmx
102tob_ile=beta_prescrit/2.
103!moteurmax=1.e5
104!moteurmax=0.
105
106!-------------------------------------------------------------------
107! masque stream : lecture du masque
108
109masque_stream  = trim(dirnameinp)//trim(masque_stream)
110
111!call lect_input(1,'z',1,work_tab,masque_stream,file_ncdf)
112
113call lect_input(1,'work_tab1',1,work_tab,masque_stream,file_ncdf)
114mstream(:,:) = int(work_tab(:,:))
115
116
117! elargit les streams autorisés aux noeuds mineurs adjacents
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         gzmx(i:i+1,j)=.true.
175         gzmy(i,j:j+1)=.true.
176      endif
177   end do
178end do
179
180! le calcul des fleuves se fait sur les mailles majeures
181! normalement, les points cotiers sont deja tagges gzmx dans la routine flottab
182
183do j=2,ny-1
184   do i=2,nx-1
185      if ( (T(i,j,nz)-Tpmp(i,j,nz).gt.Tstream_lim).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
198dragmx : 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),beta_prescrit)
206      else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then  ! stream -> beta_prescrit
207         betamx(i,j)=beta_prescrit
208         betamx(i,j)=min(betamx(i,j),beta_prescrit)
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=beta_prescrit/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 dragmx
227
228dragmy : 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),beta_prescrit)
236      else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then  ! stream -> beta_prescrit
237         betamy(i,j)=beta_prescrit
238         betamy(i,j)=min(betamy(i,j),beta_prescrit)
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=beta_prescrit/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 dragmy
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
299drag_mx(:,:) = betamx(:,:)
300drag_my(:,:) = betamy(:,:)
301
302return
303end subroutine dragging
304
305end module dragging_LGM
306
307
308
Note: See TracBrowser for help on using the repository browser.