source: branches/iLoveclim/SOURCES/dragging_hwat-contigu_mod.f90 @ 244

Last change on this file since 244 was 187, checked in by aquiquet, 6 years ago

Grisli-iloveclim branch merged to trunk at revision 185

File size: 7.3 KB
Line 
1!> \file dragging_hwat-contigu_mod.f90
2!! Defini les zones de stream avec 2 criteres
3!<
4
5!> \namespace  dragging_hwat_cont
6!! Defini les zones de stream avec 2 criteres
7!! \author ...
8!! \date ...
9!! @note Les 2 criteres sont:
10!! @note   * un critere sur la hauteur d'eau
11!! @note   * un critere de continuite
12!! depuis la cote
13!! @note Used module
14!! @note   - use module3D_phy
15!<
16module dragging_hwat_cont
17
18! Defini les zones de stream avec un critere sur la hauteur d'eau + un critere de continuite
19! depuis la cote
20
21use module3d_phy
22use fake_beta_iter_vitbil_mod
23
24implicit none
25!logical,dimension(nx,ny) :: fleuvemx
26!logical,dimension(nx,ny) :: fleuvemy
27real,dimension(nx,ny) :: hwatmx
28real,dimension(nx,ny) :: hwatmY
29
30real :: valmax
31integer :: imax,jmax
32integer :: i_moins1,i_plus1,j_moins1,j_plus1
33integer :: lmax=20
34integer :: idep,jdep,iloc,jloc
35
36real :: tostick
37
38contains
39!-------------------------------------------------------------------------------
40!> SUBROUTINE: init_dragging
41!! Cette routine fait l'initialisation du dragging.
42!>
43subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
44
45implicit none
46 
47namelist/drag_hwat_cont/hwatstream,cf,betamax,toblim,tostick   
48
49if (itracebug.eq.1)  call tracebug(' Dragging avec hwatermax')
50
51
52! formats pour les ecritures dans 42
53428 format(A)
54
55! lecture des parametres du run                      block draghwat
56!--------------------------------------------------------------------
57
58rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
59read(num_param,drag_hwat_cont)
60
61write(num_rep_42,428)'!___________________________________________________________' 
62write(num_rep_42,428) '&drag_hwat_cont              ! nom du bloc hwatermax '
63write(num_rep_42,*)
64write(num_rep_42,*) 'hwatsream      = ',hwatstream 
65write(num_rep_42,*) 'cf             = ',cf
66write(num_rep_42,*) 'betamax         = ', betamax
67write(num_rep_42,*) 'toblim         = ', toblim
68write(num_rep_42,*)'/'                           
69write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream en partant de la cote'
70write(num_rep_42,428) '!  si hwater > hwatstream '
71write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
72write(num_rep_42,428) '! betamax : (Pa) frottement maxi '
73write(num_rep_42,428) '! toblim : (Pa) tes pour les iles '
74write(num_rep_42,*)
75
76inv_beta=0
77tostick=1.e5   ! valeurs pour les points non flgzmx
78
79!-------------------------------------------------------------------
80! masque stream
81
82      mstream_mx(:,:)=1
83      mstream_my(:,:)=1
84
85
86! coefficient permettant de modifier le basal drag.
87      drag_mx(:,:)=1.
88      drag_my(:,:)=1.
89
90   
91      return
92      end subroutine init_dragging
93!________________________________________________________________________________
94
95!> SUBROUTINE: dragging
96!! Defini la localisation des streams et le frottement basal
97!>
98!-------------------------------------------------------------------------
99subroutine dragging   ! defini la localisation des streams et le frottement basal
100
101
102!         les iles n'ont pas de condition neff mais ont la condition toblim
103!         (ce bloc etait avant dans flottab)
104
105
106do j=2,ny
107   do i=2,nx
108      ilemx(i,j)=ilemx(i,j).and.(abs(tobmx(i,j)).lt.toblim)
109      ilemy(i,j)=ilemy(i,j).and.(abs(tobmy(i,j)).lt.toblim)
110   end do
111end do
112
113! calcul de la hauteur d'eau sur les mailles mineures
114do j=1,ny
115   do i=2,nx
116      hwatmx(i,j)=0.5*(hwater(i,j)+hwater(i-1,j))
117   end do
118end do
119
120do j=2,ny
121   do i=1,nx
122      hwatmy(i,j)=0.5*(hwater(i,j)+hwater(i,j-1))
123   end do
124end do
125
126
127debug_3D(:,:,1)=hwatmx(:,:)
128debug_3D(:,:,2)=hwatmy(:,:)
129
130
131! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
132
133! Partant de la cote, le fleuve est determine en remontant vers le point de hauteur d'eau max
134! tant que ce point est superieur a hwatermax.
135
136! Attention : ce systeme ne permet pas d'avoir des fleuves convergents
137!             et les fleuves n'ont une épaisseur que de 1 noeud.
138
139fleuvemx(:,:)=.false.
140fleuvemy(:,:)=.false.
141
142fleuve_mx: do j=1,ny 
143   do i=2,nx                      ! le noeud 1 n'est pas attribue
144
145cote_mx :  if (gzmx(i,j)) then    ! a ce niveau les gzmx sont sur la grounding line
146         idep=i
147         jdep=j
148
149suit_x : do l=1,lmax         ! debut de la boucle de suivi, lmax longueur maxi des fleuves
150           i_moins1=max(idep-1,2)
151           j_moins1=max(jdep-1,1)
152           i_plus1=min(idep+1,nx)
153           j_plus1=min(jdep+1,ny)
154           
155           valmax=0.
156
157! recherche du max en excluant les points flottants
158         do jloc=j_moins1,j_plus1
159            do iloc=i_moins1,i_plus1
160 
161               if ((hwatmx(iloc,jloc).gt.valmax).and.(.not.flotmx(iloc,jloc))) then
162                  imax=iloc
163                  jmax=jloc
164                  valmax=hwatmx(iloc,jloc)
165               endif
166            end do
167         end do
168
169         if (hwatmx(imax,jmax).gt.hwatstream) then
170            fleuvemx(imax,jmax)=.true.
171            idep=imax
172            jdep=jmax
173         end if
174      end do suit_x
175
176   end if cote_mx
177
178if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. 
179
180
181
182! calcul du frottement basal (ce bloc etait avant dans neffect)
183         if (gzmx(i,j)) then                          ! stream
184            betamx(i,j)=cf*neffmx(i,j) 
185         else if (flotmx(i,j).or.ilemx(i,j)) then     ! flottant ou ile
186            betamx(i,j)=0.
187         else                                         ! grounded, SIA
188            betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
189         endif     
190
191 end do
192end do fleuve_mx
193
194
195
196
197fleuve_my: do j=2,ny                  ! le noeud 1 n'est pas attribue
198   do i=1,nx   
199
200cote_my :  if (gzmy(i,j)) then    ! a ce niveau les gzmx sont sur la grounding line
201         idep=i
202         jdep=j
203
204suit_y : do l=1,lmax         ! debut de la boucle de suivi, lmax longueur maxi des fleuves
205           i_moins1=max(idep-1,1)
206           j_moins1=max(jdep-1,2)
207           i_plus1=min(idep+1,nx)
208           j_plus1=min(jdep+1,ny)
209           
210           valmax=0.
211
212! recherche du max en excluant les points flottants
213         do jloc=j_moins1,j_plus1
214            do iloc=i_moins1,i_plus1
215 
216               if ((hwatmy(iloc,jloc).gt.valmax).and.(.not.flotmy(iloc,jloc))) then
217                  imax=iloc
218                  jmax=jloc
219                  valmax=hwatmy(iloc,jloc)
220               endif
221            end do
222         end do
223
224         if (hwatmy(imax,jmax).gt.hwatstream) then
225            fleuvemy(imax,jmax)=.true.
226            idep=imax
227            jdep=jmax
228         end if
229      end do suit_y
230
231   end if cote_my
232
233if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. 
234
235
236 ! calcul du frottement basal (ce bloc etait avant dans neffect)
237   if (gzmy(i,j)) then                ! stream
238      betamy(i,j)=cf*neffmy(i,j) 
239   else if (flotmy(i,j).or.ilemy(i,j)) then     ! flottant ou ile
240      betamy(i,j)=0.
241   else                               ! grounded, SIA
242      betamy(i,j)=tostick             ! frottement glace posee
243   endif
244
245
246 end do
247end do fleuve_my
248
249
250
251where (fleuvemx(:,:))
252   debug_3D(:,:,3)=1
253elsewhere
254   debug_3D(:,:,3)=0
255endwhere
256
257where (flotmx(:,:))
258   debug_3D(:,:,3)=-1
259endwhere
260
261!_____________________
262
263
264where (fleuvemy(:,:))
265   debug_3D(:,:,4)=1
266elsewhere
267   debug_3D(:,:,4)=0
268endwhere
269
270where (flotmy(:,:))
271   debug_3D(:,:,4)=-1
272end where
273
274
275return
276end subroutine dragging
277
278end module dragging_hwat_cont
279
Note: See TracBrowser for help on using the repository browser.