source: trunk/SOURCES/dragging_hwat_contmaj_mod.f90 @ 25

Last change on this file since 25 was 9, checked in by dumas, 9 years ago

Mise en place de Hemin-40 avec nouveaux module climat : climat_forcage_mois_mod.f90, ablation_mod.f90, pdd_declar_mod.f90. Suppression de l'appel à lect-clim-act-hemin40_mod.f90

File size: 10.7 KB
Line 
1!> \file dragging_hwat_contmaj_mod.f90
2!! Defini les zones de stream avec 3 criteres
3!<
4
5!> \namespace  dragging_hwat_contmaj
6!! Defini les zones de stream avec 3 criteres
7!! \author ...
8!! \date ...
9!! @note Les trois criteres sont:
10!! @note   * un critere sur la hauteur d'eau
11!! @note   * un critere de continuite depuis la cote
12!! @note   * un critere sur la courbure du socle (si negatif, vallees)
13!! @note Used module
14!! @note   - use module3D_phy
15!! @note   - use param_phy_mod
16!<
17
18module dragging_hwat_contmaj
19
20! Defini les zones de stream avec :
21! * un critere sur la hauteur d'eau
22! * un critere de continuite depuis la cote
23! * un critere sur la courbure du socle (si negatif, vallees)
24
25use module3d_phy
26use param_phy_mod
27
28implicit none
29logical,dimension(nx,ny) :: fleuvemx
30logical,dimension(nx,ny) :: fleuvemy
31logical,dimension(nx,ny) :: fleuve
32logical,dimension(nx,ny) :: cote
33
34
35real :: valmax
36integer :: imax,jmax
37integer :: i_moins1,i_plus1,j_moins1,j_plus1
38integer :: lmax=20
39integer :: idep,jdep,iloc,jloc
40
41
42real :: tostick   ! pour la glace posee
43real :: tob_ile    ! pour les iles
44real :: cry_lim=50.  ! courbure limite pour le suivi des fleuves
45
46
47real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
48real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
49real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
50real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
51logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
52
53contains
54!-------------------------------------------------------------------------------
55
56!> SUBROUTINE: init_dragging
57!! Cette routine fait l'initialisation du dragging.
58!>
59subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
60
61implicit none
62 
63namelist/drag_hwat_cont/hwatstream,cf,betamax,toblim,tostick   
64
65if (itracebug.eq.1)  call tracebug(' dragging avec hwatermax')
66
67
68! formats pour les ecritures dans 42
69428 format(A)
70
71! lecture des parametres du run                      block draghwat
72!--------------------------------------------------------------------
73
74rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
75read(num_param,drag_hwat_cont)
76
77write(num_rep_42,428)'!___________________________________________________________' 
78write(num_rep_42,428) '&drag_hwat_cont              ! nom du bloc hwatermax '
79write(num_rep_42,*)
80write(num_rep_42,*) 'hwatstream      = ',hwatstream 
81write(num_rep_42,*) 'cf             = ',cf
82write(num_rep_42,*) 'betamax         = ', betamax
83write(num_rep_42,*) 'toblim         = ', toblim
84write(num_rep_42,*)'/'                           
85write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream en partant de la cote'
86write(num_rep_42,428) '!  si hwater > hwatstream '
87write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
88write(num_rep_42,428) '! seulement pour les points cotiers'
89write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
90write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
91write(num_rep_42,*)
92
93tostick=1.e5   ! valeurs pour les points non flgzmx
94tob_ile=betamax/2.
95moteurmax=toblim
96
97!-------------------------------------------------------------------
98! masque stream
99
100      mstream_mx(:,:)=1
101      mstream_my(:,:)=1
102
103
104! coefficient permettant de modifier le basal drag.
105      drag_mx(:,:)=1.
106      drag_my(:,:)=1.
107
108where (socle_cry(:,:).lt.cry_lim)
109   debug_3D(:,:,25)=1
110elsewhere
111   debug_3D(:,:,25)=0
112endwhere
113
114   
115      return
116      end subroutine init_dragging
117!________________________________________________________________________________
118
119!> SUBROUTINE: dragging
120!! Defini la localisation des streams et le frottement basal
121!>
122!-------------------------------------------------------------------------
123subroutine dragging   ! defini la localisation des streams et le frottement basal
124
125
126!         les iles n'ont pas de condition neff mais ont la condition toblim
127!         (ce bloc etait avant dans flottab)
128
129
130do j=2,ny
131   do i=2,nx
132      ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
133      ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
134   end do
135end do
136
137! le calcul des fleuves se fait sur les mailles majeures en partant de la cote
138
139
140
141
142! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
143
144! Partant de la cote, le fleuve est determine en remontant vers le point de hauteur d'eau max
145! tant que ce point est superieur a hwatermax.
146
147! Attention : ce systeme ne permet pas d'avoir des fleuves convergents
148!             et les fleuves n'ont une épaisseur que de 1 noeud.
149
150fleuvemx(:,:)=.false.
151fleuvemy(:,:)=.false.
152fleuve(:,:)=.false.
153cote(:,:)=.false.
154
155! detection des cotes
156do  j=2,ny-1 
157   do i=2,nx-1 
158      if ((.not.flot(i,j)).and. & 
159      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
160         cote(i,j)=.true.
161      endif
162   end do
163end do
164
165fleuve_maj: do j=2,ny-1 
166ifleuve:   do i=2,nx-1                     
167
168cote_detect :  if (cote(i,j)) then
169         idep=i
170         jdep=j
171
172!         if (hwater(i,j).gt.hwatstream) then
173         if (socle_cry(i,j).lt.0.) then    ! dans une vallee
174            fleuve(i,j)=.true. 
175!           write(166,*)
176!           write(166,*)'cote',i,j,hwater(i,j),hwatstream
177         else
178            cote(i,j)=.false.
179            cycle ifleuve
180         endif
181
182suit : do l=1,lmax         ! debut de la boucle de suivi, lmax longueur maxi des fleuves
183           i_moins1=max(idep-1,2)
184           j_moins1=max(jdep-1,1)
185           i_plus1=min(idep+1,nx)
186           j_plus1=min(jdep+1,ny)
187           
188
189
190!!$! recherche du max en suivant la hauteur d'eau maxi chez les voisins
191!!$! * en excluant les points flottants
192!!$! * et ceux qui sont deja tagges fleuves   !version 18 juillet 2007
193
194!!$           valmax=0.
195!!$
196!!$         do jloc=j_moins1,j_plus1
197!!$            do iloc=i_moins1,i_plus1
198!!$
199!!$               if ((hwater(iloc,jloc).gt.valmax) &
200!!$                    .and.(.not.flot(iloc,jloc))  &
201!!$                    .and.(.not.fleuve(iloc,jloc))) then
202!!$                  imax=iloc
203!!$                  jmax=jloc
204!!$                  valmax=hwater(iloc,jloc)
205!!$               endif
206!!$            end do
207!!$         end do
208
209! recherche du max en suivant le socle le plus profond
210! * en excluant les points flottants
211! * et ceux qui sont deja tagges fleuves
212
213           valmax=1000.
214
215        do jloc=j_moins1,j_plus1
216            do iloc=i_moins1,i_plus1
217 
218               if ((B(iloc,jloc).lt.valmax)      & 
219                    .and.(.not.flot(iloc,jloc))  &
220                    .and.(.not.fleuve(iloc,jloc)).and.(socle_cry(iloc,jloc).lt.cry_lim)) then
221                  imax=iloc
222                  jmax=jloc
223                  valmax=B(iloc,jloc)
224               endif
225            end do
226         end do
227
228
229
230!         if (hwater(imax,jmax).gt.max(hwatstream-L*40.,200.)) then  ! version 18 juil
231         if ((hwater(imax,jmax).gt.hwatstream).and.(socle_cry(i,j).lt.cry_lim)) then
232            fleuve(imax,jmax)=.true.
233            idep=imax
234            jdep=jmax
235!            write(166,*) imax,jmax,hwater(imax,jmax)
236         else
237            fleuve(imax,jmax)=.false.
238            exit suit
239         end if
240 
241      end do suit
242
243   end if cote_detect
244
245end do ifleuve
246end do fleuve_maj
247
248
249! detection de problemes
250!write(166,*)
251!write(166,*) 'detection de pb dans fleuve   time=',time
252
253do j=2,ny
254   do i=2,nx
255      if (fleuve(i,j).and.(socle_cry(i,j).gt.cry_lim)) then
256         write(166,*) 'pb fleuve',time,i,j,hwater(i,j),socle_cry(i,j)
257      endif
258   end do
259end do
260
261do j=1,ny-1
262   do i=1,nx-1
263      if (fleuve(i,j)) then
264         fleuvemx(i,j)=.true.
265         fleuvemx(i+1,j)=.true.
266         fleuvemy(i,j)=.true.
267         fleuvemy(i,j+1)=.true.
268      end if
269   end do
270end do
271
272! pas de fleuve dans les endroits trop en pente
273
274fleuvemx(:,:)=fleuvemx(:,:).and.(abs(rog*Hmx(:,:)*sdx(:,:)).lt.toblim)
275fleuvemy(:,:)=fleuvemy(:,:).and.(abs(rog*Hmy(:,:)*sdy(:,:)).lt.toblim)
276
277do j=1,ny
278   do i=1,nx
279      if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. 
280
281! calcul du frottement basal (ce bloc etait avant dans neffect)
282
283      if (cotemx(i,j)) then                        ! point cotier
284         betamx(i,j)=cf*neffmx(i,j) 
285         betamx(i,j)=min(betamx(i,j),betamax)
286      else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then  ! stream
287         betamx(i,j)=betamax
288         betamx(i,j)=cf*neffmx(i,j)*20.                  ! attention a cette modif juillet 2007
289         betamx(i,j)=min(betamx(i,j),betamax)
290         betamx(i,j)=max(betamx(i,j),10.)
291
292         if (cf*neffmx(i,j).gt.1500.) then
293            gzmx(i,j)=.false.
294            betamx(i,j)=tostick
295         endif
296
297
298      else if (ilemx(i,j)) then
299         betamx(i,j)=cf*neffmx(i,j) 
300         betamx(i,j)=min(betamx(i,j),tob_ile)       
301      else if (flotmx(i,j)) then     ! flottant ou ile
302         betamx(i,j)=0.
303      else                                         ! grounded, SIA
304         betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
305      endif
306
307   end do
308end do
309
310
311do j=1,ny                  ! le noeud 1 n'est pas attribue
312   do i=1,nx   
313
314      if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. 
315
316
317 ! calcul du frottement basal (ce bloc etait avant dans neffect)
318
319      if (cotemy(i,j)) then                        ! point cotier
320         betamy(i,j)=cf*neffmy(i,j) 
321         betamy(i,j)=min(betamy(i,j),betamax)
322      else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then  ! stream
323         betamy(i,j)=betamax
324         betamy(i,j)=cf*neffmy(i,j)*20. 
325         betamy(i,j)=min(betamy(i,j),betamax)
326         betamy(i,j)=max(betamy(i,j),10.)
327
328         if (cf*neffmy(i,j).gt.1500.) then
329            gzmy(i,j)=.false.
330            betamy(i,j)=tostick
331         endif
332
333      else if (ilemy(i,j)) then
334         betamy(i,j)=cf*neffmy(i,j) 
335         betamy(i,j)=min(betamy(i,j),tob_ile)     
336      else if (flotmy(i,j)) then     ! flottant ou ile
337         betamy(i,j)=0.
338      else                               ! grounded, SIA
339         betamy(i,j)=tostick             ! frottement glace posee
340      endif
341   end do
342end do
343
344
345where (fleuve(:,:))
346   debug_3D(:,:,1)=1
347elsewhere (flot(:,:))
348   debug_3D(:,:,1)=2
349elsewhere
350   debug_3D(:,:,1)=0
351endwhere
352
353where (cote(:,:))
354   debug_3D(:,:,2)=1
355elsewhere
356   debug_3D(:,:,2)=0
357endwhere
358
359where (fleuvemx(:,:))
360   debug_3D(:,:,3)=1
361elsewhere
362   debug_3D(:,:,3)=0
363endwhere
364
365where (flotmx(:,:))
366   debug_3D(:,:,3)=-1
367endwhere
368
369!_____________________
370
371
372where (fleuvemy(:,:))
373   debug_3D(:,:,4)=1
374elsewhere
375   debug_3D(:,:,4)=0
376endwhere
377
378where (flotmy(:,:))
379   debug_3D(:,:,4)=-1
380end where
381
382debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) 
383debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) 
384
385return
386end subroutine dragging
387
388end module dragging_hwat_contmaj
389
Note: See TracBrowser for help on using the repository browser.