1 | !> \file dragging_hudson_jorge_mod.f90 |
---|
2 | !! Module qui definie les zones de stream |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace dragging_calc_beta |
---|
6 | !! Definie les zones de stream avec uniquement: |
---|
7 | !! @note - un critere sur la hauteur d'eau |
---|
8 | !! @note - presence des sediments |
---|
9 | !! @note jalv: module inspire par le dragging_hwat_contmaj mais simplifie |
---|
10 | !! \author ... |
---|
11 | !! \date ... |
---|
12 | !! @note Used module |
---|
13 | !! @note - use module3D_phy |
---|
14 | !! @note - use sedim_declar |
---|
15 | !< |
---|
16 | |
---|
17 | module dragging_hudson_jorge_mod |
---|
18 | |
---|
19 | !jalv: module inspire par le dragging_hwat_contmaj mais simplifie |
---|
20 | ! Definie les zones de stream avec uniquement: |
---|
21 | ! * un critere sur la hauteur d'eau |
---|
22 | ! * presence des sediments |
---|
23 | |
---|
24 | |
---|
25 | use module3d_phy |
---|
26 | use sedim_declar |
---|
27 | |
---|
28 | implicit none |
---|
29 | logical,dimension(nx,ny) :: fleuvemx |
---|
30 | logical,dimension(nx,ny) :: fleuvemy |
---|
31 | logical,dimension(nx,ny) :: fleuve |
---|
32 | logical,dimension(nx,ny) :: cote |
---|
33 | |
---|
34 | real :: seuil_sedim !< seuil sur hwater pour avoir du glissement |
---|
35 | real :: valmax |
---|
36 | integer :: imax,jmax |
---|
37 | integer :: i_moins1,i_plus1,j_moins1,j_plus1 |
---|
38 | integer :: lmax=20 |
---|
39 | integer :: idep,jdep,iloc,jloc |
---|
40 | |
---|
41 | integer :: itestd |
---|
42 | |
---|
43 | real,dimension(nx,ny) :: betamx_temp !< pour lissage du betamx |
---|
44 | real,dimension(nx,ny) :: betamy_temp !< pour lissage du betamy |
---|
45 | |
---|
46 | |
---|
47 | real :: tostick !< pour la glace posee |
---|
48 | real :: tob_ile !< pour les iles |
---|
49 | real :: cry_lim=50. !< courbure limite pour le suivi des fleuves |
---|
50 | contains |
---|
51 | !------------------------------------------------------------------------------- |
---|
52 | |
---|
53 | !> SUBROUTINE: init_dragging |
---|
54 | !! Cette routine fait l'initialisation du dragging |
---|
55 | !< |
---|
56 | subroutine init_dragging ! Cette routine fait l'initialisation du dragging. |
---|
57 | |
---|
58 | implicit none |
---|
59 | |
---|
60 | namelist/drag_hudson_jorge/hwatstream,cf,betamax,toblim,seuil_sedim |
---|
61 | |
---|
62 | if (itracebug.eq.1) call tracebug(' dragging avec hwatermax') |
---|
63 | |
---|
64 | |
---|
65 | ! formats pour les ecritures dans 42 |
---|
66 | 428 format(A) |
---|
67 | |
---|
68 | ! lecture des parametres du run block dra_hudson_jorge |
---|
69 | !-------------------------------------------------------------------- |
---|
70 | |
---|
71 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
72 | read(num_param,drag_hudson_jorge) |
---|
73 | |
---|
74 | write(num_rep_42,428)'!___________________________________________________________' |
---|
75 | write(num_rep_42,428) '&drag_hudson_jorge ! nom du bloc drag_hudson_jorge ' |
---|
76 | write(num_rep_42,*) |
---|
77 | write(num_rep_42,*) 'hwatstream = ',hwatstream |
---|
78 | write(num_rep_42,*) 'cf = ',cf |
---|
79 | write(num_rep_42,*) 'betamax = ', betamax |
---|
80 | write(num_rep_42,*) 'toblim = ', toblim |
---|
81 | write(num_rep_42,*) 'seuil_sedim = ', seuil_sedim |
---|
82 | write(num_rep_42,*)'/' |
---|
83 | write(num_rep_42,428) '! hwatstream (m) : critere de passage en stream en partant de la cote' |
---|
84 | write(num_rep_42,428) '! si hwater > hwatstream ' |
---|
85 | write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff' |
---|
86 | write(num_rep_42,428) '! seulement pour les points cotiers' |
---|
87 | write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams ' |
---|
88 | write(num_rep_42,428) '! toblim : (Pa) pour les iles ' |
---|
89 | write(num_rep_42,*) |
---|
90 | |
---|
91 | tostick=1.e5 ! valeurs pour les points non flgzmx |
---|
92 | tob_ile=betamax/2. |
---|
93 | moteurmax=toblim |
---|
94 | |
---|
95 | !------------------------------------------------------------------- |
---|
96 | ! masque stream |
---|
97 | |
---|
98 | mstream_mx(:,:)=1 |
---|
99 | mstream_my(:,:)=1 |
---|
100 | |
---|
101 | ! coefficient permettant de modifier le basal drag. |
---|
102 | drag_mx(:,:)=1. |
---|
103 | drag_my(:,:)=1. |
---|
104 | |
---|
105 | return |
---|
106 | end subroutine init_dragging |
---|
107 | !________________________________________________________________________________ |
---|
108 | |
---|
109 | !> SUBROUTINE: dragging |
---|
110 | !! Defini la localisation des streams et le frottement basal |
---|
111 | !< |
---|
112 | subroutine dragging ! defini la localisation des streams et le frottement basal |
---|
113 | |
---|
114 | ! les iles n'ont pas de condition neff mais ont la condition toblim |
---|
115 | ! (ce bloc etait avant dans flottab) |
---|
116 | |
---|
117 | |
---|
118 | do j=2,ny |
---|
119 | do i=2,nx |
---|
120 | ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim) |
---|
121 | ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim) |
---|
122 | end do |
---|
123 | end do |
---|
124 | |
---|
125 | |
---|
126 | ! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line |
---|
127 | |
---|
128 | |
---|
129 | fleuvemx(:,:)=.false. |
---|
130 | fleuvemy(:,:)=.false. |
---|
131 | fleuve(:,:)=.false. |
---|
132 | cote(:,:)=.false. |
---|
133 | |
---|
134 | |
---|
135 | ! detection des cotes |
---|
136 | do j=2,ny-1 |
---|
137 | do i=2,nx-1 |
---|
138 | if ((.not.flot(i,j)).and. & |
---|
139 | ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then |
---|
140 | cote(i,j)=.true. |
---|
141 | endif |
---|
142 | end do |
---|
143 | end do |
---|
144 | |
---|
145 | |
---|
146 | ! critere de passage en stream |
---|
147 | do j=2,ny-1 |
---|
148 | do i=2,nx-1 |
---|
149 | if ((hwater(i,j).gt.hwatstream).or.((mksedim(i,j).eq.2).and.(hwater(i,j).ge.seuil_sedim))) then |
---|
150 | fleuve(i,j)=.true. |
---|
151 | else |
---|
152 | fleuve(i,j)=.false. |
---|
153 | endif |
---|
154 | enddo |
---|
155 | enddo |
---|
156 | |
---|
157 | do j=1,ny-1 |
---|
158 | do i=1,nx-1 |
---|
159 | if (fleuve(i,j)) then |
---|
160 | fleuvemx(i,j)=.true. |
---|
161 | fleuvemx(i+1,j)=.true. |
---|
162 | fleuvemy(i,j)=.true. |
---|
163 | fleuvemy(i,j+1)=.true. |
---|
164 | end if |
---|
165 | end do |
---|
166 | end do |
---|
167 | |
---|
168 | ! pas de fleuve dans les endroits trop en pente |
---|
169 | |
---|
170 | !fleuvemx(:,:)=fleuvemx(:,:).and.(abs(rog*Hmx(:,:)*sdx(:,:)).lt.toblim) |
---|
171 | !fleuvemy(:,:)=fleuvemy(:,:).and.(abs(rog*Hmy(:,:)*sdy(:,:)).lt.toblim) |
---|
172 | |
---|
173 | !..................... |
---|
174 | |
---|
175 | !jalv |
---|
176 | !call detect_assym(nx,ny,0,76,1,0,1,0,betamx,itestd) |
---|
177 | !if (itestd.gt.0) then |
---|
178 | ! write(6,*) 'avant calcul betax asymetrie sur betamx pour time=',time |
---|
179 | ! stop |
---|
180 | !else |
---|
181 | ! write(6,*) 'avant calcul betax pas d asymetrie sur betamx pour time=',time |
---|
182 | !end if |
---|
183 | |
---|
184 | !call detect_assym(nx,ny,0,76,1,0,1,0,betamy,itestd) |
---|
185 | !if (itestd.gt.0) then |
---|
186 | ! write(6,*) 'avant calcul betay asymetrie sur betamy pour time=',time |
---|
187 | ! stop |
---|
188 | !else |
---|
189 | ! write(6,*) 'avant calcul betay pas d asymetrie sur betamy pour time=',time |
---|
190 | !end if |
---|
191 | |
---|
192 | |
---|
193 | |
---|
194 | !...........................................! |
---|
195 | ! frottement selon x (avant dans neffect) ! ! test jalv: tous les cas ont la meme loi pour beta: beta=cf*neff |
---|
196 | !...........................................! ! mais le sueil maxi est different pour chaque cas |
---|
197 | |
---|
198 | |
---|
199 | |
---|
200 | do j=1,ny |
---|
201 | do i=1,nx |
---|
202 | |
---|
203 | if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. |
---|
204 | |
---|
205 | if (cotemx(i,j)) then ! point cotier: frott maxi=40, mini=4 |
---|
206 | betamx(i,j)=cf*neffmx(i,j) |
---|
207 | betamx(i,j)=min(betamx(i,j),40.) |
---|
208 | betamx(i,j)=max(betamx(i,j),4.) |
---|
209 | else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then ! stream normal: frott maxi=80, mini=8 |
---|
210 | betamx(i,j)=cf*neffmx(i,j) |
---|
211 | betamx(i,j)=min(betamx(i,j),80.) |
---|
212 | betamx(i,j)=max(betamx(i,j),8.) |
---|
213 | |
---|
214 | if (mksedim(i,j).eq.2) then ! jalv: stream en zone sediment |
---|
215 | betamx(i,j)=cf*neffmx(i,j) |
---|
216 | betamx(i,j)=min(betamx(i,j),60.) !stream sediment: frott maxi=60, mini=6 |
---|
217 | betamx(i,j)=max(betamx(i,j),6.) |
---|
218 | endif |
---|
219 | |
---|
220 | else if (ilemx(i,j)) then ! ile : frott maxi=20, mini=2 |
---|
221 | betamx(i,j)=cf*neffmx(i,j) |
---|
222 | betamx(i,j)=min(betamx(i,j),20.) |
---|
223 | betamx(i,j)=max(betamx(i,j),2.) |
---|
224 | else if (flotmx(i,j)) then ! flottant : frott 0. |
---|
225 | betamx(i,j)=0. |
---|
226 | else ! grounded, SIA |
---|
227 | ! betamx(i,j)=tostick ! frottement glace posee (1 bar) divise par mil !jalv |
---|
228 | betamx(i,j)=200. ! jalv: le cas restant glace posee a priori: frott=200 |
---|
229 | endif |
---|
230 | |
---|
231 | ! if (cf*neffmx(i,j).gt.1500.) then |
---|
232 | ! gzmx(i,j)=.false. |
---|
233 | ! betamx(i,j)=tostick |
---|
234 | ! endif |
---|
235 | |
---|
236 | end do |
---|
237 | end do |
---|
238 | |
---|
239 | |
---|
240 | !................................................! |
---|
241 | !...test..jalv..lissage..des..beta...............! |
---|
242 | !...on..calcule..la..moyenne..ponderee..des......! |
---|
243 | !...points..aux..alentours.......................! |
---|
244 | |
---|
245 | |
---|
246 | do j=1,ny |
---|
247 | do i=3,nx-2 |
---|
248 | |
---|
249 | betamx_temp(i,j)=(betamx(i,j)+(0.3*betamx(i-1,j))+(0.2*betamx(i-2,j))+(0.3*betamx(i+1,j))+(0.2*betamx(i+2,j)))/2 |
---|
250 | |
---|
251 | end do |
---|
252 | end do |
---|
253 | |
---|
254 | |
---|
255 | |
---|
256 | betamx(:,:)=betamx_temp(:,:) |
---|
257 | |
---|
258 | |
---|
259 | |
---|
260 | !...........................................! |
---|
261 | ! frottement selon y (avant dans neffect) ! ! test jalv: tous les cas ont la meme loi pour beta: beta=cf*neff |
---|
262 | !...........................................! ! mais le seuil maxi est different pour chaque cas |
---|
263 | |
---|
264 | |
---|
265 | do j=1,ny |
---|
266 | do i=1,nx |
---|
267 | |
---|
268 | if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. |
---|
269 | |
---|
270 | |
---|
271 | if (cotemy(i,j)) then ! point cotier: froot maxi=40, mini=4 |
---|
272 | betamy(i,j)=cf*neffmy(i,j) |
---|
273 | betamy(i,j)=min(betamy(i,j),40.) |
---|
274 | betamy(i,j)=max(betamy(i,j),4.) |
---|
275 | else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then ! stream normal: frott maxi=80, mini=8 |
---|
276 | betamy(i,j)=cf*neffmy(i,j) |
---|
277 | betamy(i,j)=min(betamy(i,j),80.) |
---|
278 | betamy(i,j)=max(betamy(i,j),8.) |
---|
279 | |
---|
280 | if (mksedim(i,j).eq.2) then ! jalv: stream en zone sediment |
---|
281 | betamy(i,j)=cf*neffmy(i,j) |
---|
282 | betamy(i,j)=min(betamy(i,j),60.) !stream sediment: frott maxi=60, mini=6. |
---|
283 | betamy(i,j)=max(betamy(i,j),6.) |
---|
284 | endif |
---|
285 | |
---|
286 | else if (ilemy(i,j)) then ! ile : frott maxi=20, mini=2 |
---|
287 | betamy(i,j)=cf*neffmy(i,j) |
---|
288 | betamy(i,j)=min(betamy(i,j),20.) |
---|
289 | betamy(i,j)=max(betamy(i,j),2.) |
---|
290 | else if (flotmy(i,j)) then ! flottant : frott 0 |
---|
291 | betamy(i,j)=0. |
---|
292 | else ! grounded, SIA |
---|
293 | ! betamy(i,j)=tostick/1000 ! frottement glace posee divise par mil !jalv |
---|
294 | betamy(i,j)=200. ! jalv: le cas restant glace posee a priori: frott maxi 200 |
---|
295 | endif |
---|
296 | |
---|
297 | ! if (cf*neffmy(i,j).gt.1500.) then |
---|
298 | ! gzmy(i,j)=.false. |
---|
299 | ! betamy(i,j)=tostick |
---|
300 | ! endif |
---|
301 | |
---|
302 | end do |
---|
303 | end do |
---|
304 | |
---|
305 | |
---|
306 | |
---|
307 | !................................................! |
---|
308 | !...test..jalv..lissage..des..beta...............! |
---|
309 | !...on..calcule..la..moyenne..ponderee..des......! |
---|
310 | !...points..aux..alentours.......................! |
---|
311 | |
---|
312 | |
---|
313 | do j=3,ny-2 |
---|
314 | do i=1,nx |
---|
315 | |
---|
316 | betamy_temp(i,j)=(betamy(i,j)+(0.3*betamy(i,j-1))+(0.2*betamy(i,j-2))+(0.3*betamy(i,j+1))+(0.2*betamy(i,j+2)))/2 |
---|
317 | |
---|
318 | end do |
---|
319 | end do |
---|
320 | |
---|
321 | betamy(:,:)=betamy_temp(:,:) |
---|
322 | |
---|
323 | |
---|
324 | |
---|
325 | !jalv |
---|
326 | !call detect_assym(nx,ny,0,76,1,0,1,0,betamx,itestd) |
---|
327 | !if (itestd.gt.0) then |
---|
328 | ! write(6,*) 'apres dragging asymetrie sur betamx pour time=',time |
---|
329 | ! stop |
---|
330 | !else |
---|
331 | ! write(6,*) 'apres dragging pas d asymetrie sur betamx pour time=',time |
---|
332 | !end if |
---|
333 | |
---|
334 | !call detect_assym(nx,ny,0,76,1,0,1,0,betamy,itestd) |
---|
335 | !if (itestd.gt.0) then |
---|
336 | ! write(6,*) 'apres dragging asymetrie sur betamy pour time=',time |
---|
337 | ! stop |
---|
338 | !else |
---|
339 | ! write(6,*) 'apres dragging pas d asymetrie sur betamy pour time=',time |
---|
340 | !end if |
---|
341 | |
---|
342 | return |
---|
343 | end subroutine dragging |
---|
344 | |
---|
345 | end module dragging_hudson_jorge |
---|