source: trunk/SOURCES/dragging_hwat-contigu_mod.f90 @ 23

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

initial import GRISLI trunk

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