1 | !> \file dragging-vit_bil_LBq_gen_mod.f90 |
---|
2 | !! Module qui defini les zones de stream d'apres les vitesses de bilan |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace dragging_vit_bil_LBq_gen |
---|
6 | !! Module qui defini les zones de stream |
---|
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 | module dragging_vit_bil_LBq_gen |
---|
20 | |
---|
21 | ! Defini les zones de stream avec : |
---|
22 | ! * un critere sur la hauteur d'eau |
---|
23 | ! * un critere de fleuves actuels |
---|
24 | ! attention il faut fonctionner en additionnant SIA et Stream pour éviter les pb |
---|
25 | ! de fleuves en pointille |
---|
26 | |
---|
27 | |
---|
28 | use module3d_phy |
---|
29 | use param_phy_mod |
---|
30 | use interface_input |
---|
31 | implicit none |
---|
32 | logical,dimension(nx,ny) :: fleuvemx !< fleuves sont les tableaux courants (dep. time) |
---|
33 | logical,dimension(nx,ny) :: fleuvemy |
---|
34 | logical,dimension(nx,ny) :: fleuve |
---|
35 | logical,dimension(nx,ny) :: cote |
---|
36 | real, dimension(nx,ny) :: Vcol_x !< uniquement pour compatibilite dans le spinup |
---|
37 | real, dimension(nx,ny) :: Vcol_y !< uniquement pour compatibilite |
---|
38 | |
---|
39 | character(len=100) :: balance_vel_file !< vitesses de bilan |
---|
40 | real, dimension(nx,ny) :: vitbil_max !< vitesse de bilan max sur la maille |
---|
41 | real :: seuil_vel !< seuil sur les vitesses pour definir le masque stream |
---|
42 | |
---|
43 | real :: valmax |
---|
44 | integer :: imax,jmax |
---|
45 | integer :: i_moins1,i_plus1,j_moins1,j_plus1 |
---|
46 | integer :: lmax=20 |
---|
47 | integer :: idep,jdep,iloc,jloc |
---|
48 | |
---|
49 | |
---|
50 | real :: tostick ! pour la glace posee |
---|
51 | real :: tob_ile ! pour les iles |
---|
52 | |
---|
53 | |
---|
54 | contains |
---|
55 | !------------------------------------------------------------------------------- |
---|
56 | |
---|
57 | !> SUBROUTINE: init_dragging |
---|
58 | !! Cette routine fait l'initialisation du dragging. |
---|
59 | !< |
---|
60 | subroutine init_dragging ! Cette routine fait l'initialisation du dragging. |
---|
61 | |
---|
62 | implicit none |
---|
63 | |
---|
64 | namelist/drag_vit_bil_LBq_gen/hwatstream,cf,betamax,toblim,seuil_vel,balance_vel_file |
---|
65 | |
---|
66 | if (itracebug.eq.1) call tracebug(' Init_dragging avec dragging_vit_bil_LBq_15') |
---|
67 | |
---|
68 | |
---|
69 | ! formats pour les ecritures dans 42 |
---|
70 | 428 format(A) |
---|
71 | |
---|
72 | ! lecture des parametres du run |
---|
73 | !-------------------------------------------------------------------- |
---|
74 | |
---|
75 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
76 | read(num_param,drag_vit_bil_LBq_gen) |
---|
77 | |
---|
78 | write(num_rep_42,428) '!___________________________________________________________' |
---|
79 | write(num_rep_42,428) '! module dragging_vit_bil_LBq_gen ' |
---|
80 | write(num_rep_42,drag_vit_bil_LBq_gen) |
---|
81 | write(num_rep_42,428) '! hwatstream (m) : critere de passage en stream en partant de la cote' |
---|
82 | write(num_rep_42,428) '! si hwater > hwatstream ' |
---|
83 | write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff' |
---|
84 | write(num_rep_42,428) '! seulement pour les points cotiers' |
---|
85 | write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams ' |
---|
86 | write(num_rep_42,428) '! toblim : (Pa) pour les iles ' |
---|
87 | write(num_rep_42,428) '! seuil_vel (m/an) : seuil sur les vitesses pour definir le masque stream ' |
---|
88 | write(num_rep_42,428) '! balance_vel_file : nom du fichier qui contient les vitesse de bilan' |
---|
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 | !------------------------------------------------------------------- |
---|
97 | ! masque stream : lecture des vitesses de bilan |
---|
98 | ! on part des vitesses de Le Brocq qui sont sur les noeuds majeurs |
---|
99 | |
---|
100 | balance_vel_file = trim(dirnameinp)//trim(balance_vel_file) |
---|
101 | |
---|
102 | call lect_input(3,'vitbil',1,vitbil_max,balance_vel_file,trim(dirnameinp)//trim(runname)//'.nc') |
---|
103 | ! call lect_datfile(nx,ny,vitbil_max,1,balance_vel_file) |
---|
104 | |
---|
105 | mstream_mx(:,:)=0 |
---|
106 | mstream_my(:,:)=0 |
---|
107 | |
---|
108 | do j=1,ny |
---|
109 | do i=1,nx |
---|
110 | if (vitbil_max(i,j).gt.seuil_vel) then |
---|
111 | mstream_mx(i:i+1,j) = 1 |
---|
112 | mstream_my(i,j:j+1) = 1 |
---|
113 | mstream(i,j) = 1 ! allowed ice streams |
---|
114 | end if |
---|
115 | end do |
---|
116 | end do |
---|
117 | |
---|
118 | ! coefficient permettant de modifier le basal drag. |
---|
119 | drag_mx(:,:)=1. |
---|
120 | drag_my(:,:)=1. |
---|
121 | |
---|
122 | |
---|
123 | return |
---|
124 | end subroutine init_dragging |
---|
125 | !________________________________________________________________________________ |
---|
126 | |
---|
127 | !> SUBROUTINE: dragging |
---|
128 | !! Defini la localisation des streams et le frottement basal |
---|
129 | !< |
---|
130 | !------------------------------------------------------------------------- |
---|
131 | subroutine dragging ! defini la localisation des streams et le frottement basal |
---|
132 | |
---|
133 | |
---|
134 | ! les iles n'ont pas de condition neff mais ont la condition toblim |
---|
135 | ! (ce bloc etait avant dans flottab) |
---|
136 | |
---|
137 | |
---|
138 | do j=2,ny |
---|
139 | do i=2,nx |
---|
140 | ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim) |
---|
141 | ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim) |
---|
142 | end do |
---|
143 | end do |
---|
144 | |
---|
145 | fleuvemx(:,:)=.false. |
---|
146 | fleuvemy(:,:)=.false. |
---|
147 | fleuve(:,:)=.false. |
---|
148 | cote(:,:)=.false. |
---|
149 | cotemx(:,:)=.false. |
---|
150 | cotemy(:,:)=.false. |
---|
151 | |
---|
152 | ! detection des cotes |
---|
153 | do j=2,ny-1 |
---|
154 | do i=2,nx-1 |
---|
155 | if ((.not.flot(i,j)).and. & |
---|
156 | ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then |
---|
157 | cote(i,j)=.true. |
---|
158 | cotemx(i:i+1,j) = .true. |
---|
159 | cotemy(i,j:j+1) = .true. |
---|
160 | endif |
---|
161 | end do |
---|
162 | end do |
---|
163 | |
---|
164 | ! le calcul des fleuves se fait sur les mailles majeures |
---|
165 | ! normalement, les points cotiers sont deja tagges gzmx dans la routine flottab |
---|
166 | |
---|
167 | do j=2,ny-1 |
---|
168 | do i=2,nx-1 |
---|
169 | if ((hwater(i,j).gt.hwatstream).and.(mstream(i,j).eq.1).and.(.not.flot(i,j))) then |
---|
170 | fleuve(i,j)=.true. |
---|
171 | fleuvemx(i:i+1,j)=.true. |
---|
172 | fleuvemy(i,j:j+1)=.true. |
---|
173 | gzmx(i:i+1,j)=.true. |
---|
174 | gzmy(i,j:j+1)=.true. |
---|
175 | endif |
---|
176 | end do |
---|
177 | end do |
---|
178 | |
---|
179 | |
---|
180 | ! calcul du frottement basal (ce bloc etait avant dans neffect) |
---|
181 | |
---|
182 | drag_mx : do j=2,ny-1 |
---|
183 | do i=2,nx-1 |
---|
184 | |
---|
185 | if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. |
---|
186 | |
---|
187 | if (cotemx(i,j)) then ! point cotier peut etre diminue si pression d'eau |
---|
188 | betamx(i,j)=cf*neffmx(i,j) |
---|
189 | betamx(i,j)=min(betamx(i,j),betamax) |
---|
190 | else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then ! stream -> betamax |
---|
191 | betamx(i,j)=betamax |
---|
192 | betamx(i,j)=min(betamx(i,j),betamax) |
---|
193 | betamx(i,j)=max(betamx(i,j),10.) |
---|
194 | |
---|
195 | ! if (cf*neffmx(i,j).gt.1500.) then |
---|
196 | ! gzmx(i,j)=.false. |
---|
197 | ! betamx(i,j)=tostick |
---|
198 | ! endif |
---|
199 | |
---|
200 | else if (ilemx(i,j)) then ! ile : tob_ile=betamax/2 + peut etre diminue si pression d'eau |
---|
201 | betamx(i,j)=cf*neffmx(i,j) |
---|
202 | betamx(i,j)=min(betamx(i,j),tob_ile) |
---|
203 | else if (flotmx(i,j)) then ! flottant ou ile |
---|
204 | betamx(i,j)=0. |
---|
205 | else ! grounded, SIA |
---|
206 | betamx(i,j)=tostick ! frottement glace posee (1 bar) |
---|
207 | endif |
---|
208 | |
---|
209 | end do |
---|
210 | end do drag_mx |
---|
211 | |
---|
212 | drag_my : do j=2,ny-1 |
---|
213 | do i=2,nx-1 |
---|
214 | |
---|
215 | if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. |
---|
216 | |
---|
217 | if (cotemy(i,j)) then ! point cotier peut etre diminue si pression d'eau |
---|
218 | betamy(i,j)=cf*neffmy(i,j) |
---|
219 | betamy(i,j)=min(betamy(i,j),betamax) |
---|
220 | else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then ! stream -> betamax |
---|
221 | betamy(i,j)=betamax |
---|
222 | betamy(i,j)=min(betamy(i,j),betamax) |
---|
223 | betamy(i,j)=max(betamy(i,j),10.) |
---|
224 | |
---|
225 | ! if (cf*neffmy(i,j).gt.1500.) then |
---|
226 | ! gzmy(i,j)=.false. |
---|
227 | ! betamy(i,j)=tostick |
---|
228 | ! endif |
---|
229 | |
---|
230 | else if (ilemy(i,j)) then ! ile : tob_ile=betamax/2 + peut etre diminue si pression d'eau |
---|
231 | betamy(i,j)=cf*neffmy(i,j) |
---|
232 | betamy(i,j)=min(betamy(i,j),tob_ile) |
---|
233 | else if (flotmy(i,j)) then ! flottant ou ile |
---|
234 | betamy(i,j)=0. |
---|
235 | else ! grounded, SIA |
---|
236 | betamy(i,j)=tostick ! frottement glace posee (1 bar) |
---|
237 | endif |
---|
238 | |
---|
239 | end do |
---|
240 | end do drag_my |
---|
241 | |
---|
242 | |
---|
243 | where (fleuve(:,:)) |
---|
244 | debug_3D(:,:,1)=1 |
---|
245 | elsewhere (flot(:,:)) |
---|
246 | debug_3D(:,:,1)=2 |
---|
247 | elsewhere |
---|
248 | debug_3D(:,:,1)=0 |
---|
249 | endwhere |
---|
250 | |
---|
251 | where (cote(:,:)) |
---|
252 | debug_3D(:,:,2)=1 |
---|
253 | elsewhere |
---|
254 | debug_3D(:,:,2)=0 |
---|
255 | endwhere |
---|
256 | |
---|
257 | where (fleuvemx(:,:)) |
---|
258 | debug_3D(:,:,3)=1 |
---|
259 | elsewhere |
---|
260 | debug_3D(:,:,3)=0 |
---|
261 | endwhere |
---|
262 | |
---|
263 | where (flotmx(:,:)) |
---|
264 | debug_3D(:,:,3)=-1 |
---|
265 | endwhere |
---|
266 | |
---|
267 | !_____________________ |
---|
268 | |
---|
269 | |
---|
270 | where (fleuvemy(:,:)) |
---|
271 | debug_3D(:,:,4)=1 |
---|
272 | elsewhere |
---|
273 | debug_3D(:,:,4)=0 |
---|
274 | endwhere |
---|
275 | |
---|
276 | where (flotmy(:,:)) |
---|
277 | debug_3D(:,:,4)=-1 |
---|
278 | end where |
---|
279 | |
---|
280 | debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) |
---|
281 | debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) |
---|
282 | |
---|
283 | return |
---|
284 | end subroutine dragging |
---|
285 | |
---|
286 | end module dragging_vit_bil_LBq_gen |
---|
287 | |
---|