1 | !> \file dragging_hwatstream_sedim_hudson.f90 |
---|
2 | !! Defini les zones de stream |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace dragging_hwatstream_sedim_hudson |
---|
6 | !! Defini les zones de stream |
---|
7 | !! @note On applique un critere sur la hauteur d'eau |
---|
8 | !! \author ... |
---|
9 | !! \date ... |
---|
10 | !! @note Used module |
---|
11 | !! @note - use module3D_phy |
---|
12 | !! @note - use sedim_declar |
---|
13 | !< |
---|
14 | |
---|
15 | module dragging_hwatstream_sedim_hudson |
---|
16 | ! Defini les zones de stream avec un critere sur la hauteur d'eau |
---|
17 | ! modification pour ajouter un frottement faible dans les zones sediments |
---|
18 | |
---|
19 | USE module3D_phy |
---|
20 | USE sedim_declar |
---|
21 | |
---|
22 | |
---|
23 | implicit none |
---|
24 | real :: tostick |
---|
25 | real :: Cfs !< coefficient du frottement des zones sediments |
---|
26 | real :: seuil_sedim !< seuil sur hw_mx pour avoir du glissement |
---|
27 | real :: coefmax |
---|
28 | real :: coefslid |
---|
29 | real :: coefdrag |
---|
30 | real, dimension(nx,ny) :: hw_mx |
---|
31 | real, dimension(nx,ny) :: hw_my |
---|
32 | |
---|
33 | contains |
---|
34 | !------------------------------------------------------------------------------- |
---|
35 | subroutine init_dragging !< Cette routine fait l'initialisation du dragging. |
---|
36 | |
---|
37 | implicit none |
---|
38 | |
---|
39 | namelist/drag_hwat_stream/hwatstream,cf,betamax,toblim,tostick,Cfs,seuil_sedim,coefmax |
---|
40 | |
---|
41 | if (itracebug.eq.1) call tracebug(' Dragging avec hwatermax') |
---|
42 | |
---|
43 | |
---|
44 | ! formats pour les ecritures dans 42 |
---|
45 | 428 format(A) |
---|
46 | |
---|
47 | ! lecture des parametres du run block drag_hwat_stream |
---|
48 | !-------------------------------------------------------------------- |
---|
49 | |
---|
50 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
51 | read(num_param,drag_hwat_stream) |
---|
52 | |
---|
53 | write(num_rep_42,428)'!___________________________________________________________' |
---|
54 | write(num_rep_42,428) '&drag_hwat_stream ! nom du bloc drag_hwat_stream ' |
---|
55 | write(num_rep_42,*) |
---|
56 | write(num_rep_42,*) 'hwatstream = ', hwatstream |
---|
57 | write(num_rep_42,*) 'cf = ', cf |
---|
58 | write(num_rep_42,*) 'betamax = ', betamax |
---|
59 | write(num_rep_42,*) 'toblim = ', toblim |
---|
60 | write(num_rep_42,*) 'Cfs = ', Cfs |
---|
61 | write(num_rep_42,*) 'seuil_sedim = ', seuil_sedim |
---|
62 | write(num_rep_42,*) 'coefmax = ', coefmax |
---|
63 | write(num_rep_42,*)'/' |
---|
64 | write(num_rep_42,428) '! hwatstream (m) : critere de passage en stream si hwater > hwatstream ' |
---|
65 | write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff' |
---|
66 | write(num_rep_42,428) '! betamax : (Pa) frottement maxi ' |
---|
67 | write(num_rep_42,428) '! toblim : (Pa) tes pour les iles ' |
---|
68 | write(num_rep_42,428) '! Cfs : coefficient de la loi de frottement sediment' |
---|
69 | write(num_rep_42,428) '! seuil_sedim : seuil pour passer en stream si zone sediment' |
---|
70 | write(num_rep_42,428) '! coefmax : hauteur d eau max dans le sediment (=hwatermax)' |
---|
71 | write(num_rep_42,*) |
---|
72 | |
---|
73 | tostick=1.e5 ! valeurs pour les points non flgzmx |
---|
74 | |
---|
75 | coefmax=hwatermax |
---|
76 | !------------------------------------------------------------------- |
---|
77 | ! masque stream |
---|
78 | |
---|
79 | mstream_mx(:,:)=1 |
---|
80 | mstream_my(:,:)=1 |
---|
81 | |
---|
82 | ! coefficient permettant de modifier le basal drag. |
---|
83 | drag_mx(:,:)=1. |
---|
84 | drag_my(:,:)=1. |
---|
85 | |
---|
86 | coefdrag=rog/Cfs |
---|
87 | Cfs=0. |
---|
88 | write(42,*) 'dragging_hwatermax_heino' |
---|
89 | write(42,*) 'coeff frottement zone sediment (Pa/m) Coefdrag=',coefdrag |
---|
90 | |
---|
91 | return |
---|
92 | end subroutine init_dragging |
---|
93 | |
---|
94 | |
---|
95 | ! defini la localisation des streams et le frottement basal |
---|
96 | subroutine dragging |
---|
97 | |
---|
98 | |
---|
99 | ! les iles n'ont pas de condition neff mais ont la condition toblim |
---|
100 | ! (ce bloc etait avant dans flottab) |
---|
101 | |
---|
102 | call moy_mxmy(nx,ny,Hwater,hw_mx,hw_my) |
---|
103 | |
---|
104 | do 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 |
---|
109 | end do |
---|
110 | |
---|
111 | |
---|
112 | gzm_x: do j=1,ny ! 2,ny |
---|
113 | do i=2,nx ! -1 |
---|
114 | |
---|
115 | ! gzmx pour les hauteurs d'eau superieures a un seuil |
---|
116 | ! attention, il semble qu'il y avait une erreur dans la precedente |
---|
117 | ! formulation car le test etait seulement sur hwater(i,j). |
---|
118 | ! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line |
---|
119 | |
---|
120 | ! mais gzmx peut être quand même vrai selon une conditon sur la pression effective |
---|
121 | gzmx(i,j)=gzmx(i,j).or.((Neffmx(i,j).lt.neffratio*hmx(i,j)) & |
---|
122 | .and..not.flotmx(i,j)) |
---|
123 | |
---|
124 | ! ce cas est un ilemx traité plus loin |
---|
125 | ! if (i.gt.2.AND.i.lt.nx) gzmx(i,j)=gzmx(i,j).or.(.not.flot(i,j).and. & |
---|
126 | ! .not.flot(i-1,j).and.flot(i+1,j).and.flot(i-2,j)) ! F P P F |
---|
127 | |
---|
128 | ! On rajoute une condition sur la contrainte basale |
---|
129 | ! si tobmx(i,j) > toblim (de l'ordre de 1 bar), passage a pose |
---|
130 | |
---|
131 | gzmx(i,j)=gzmx(i,j).and.(abs(tobmx(i,j)).lt.toblim) |
---|
132 | gzmx(i,j)=gzmx(i,j).or.((hwater(i,j).gt.hwatstream).and..not.flotmx(i,j)) |
---|
133 | |
---|
134 | ! calcul du frottement basal (ce bloc etait avant dans neffect) |
---|
135 | if (mkxsedim(i,j).eq.2) then ! loi sediment |
---|
136 | if (hw_mx(i,j).ge.seuil_sedim) then |
---|
137 | if (.not.flot(i,j)) then |
---|
138 | gzmx(i,j)=.true. |
---|
139 | coefslid=hw_mx(i,j) |
---|
140 | coefslid=max(coefslid,seuil_sedim) |
---|
141 | coefslid=min(coefslid/coefmax,1.) |
---|
142 | betamx(i,j)=min(coefdrag/coefslid,1.e5) |
---|
143 | flgzmx(i,j)=.true. |
---|
144 | mslid_mx(i,j)=5 |
---|
145 | shelfy=.true. |
---|
146 | endif |
---|
147 | endif |
---|
148 | else if (gzmx(i,j)) then ! stream |
---|
149 | betamx(i,j)=cf*neffmx(i,j) |
---|
150 | else if (flotmx(i,j).or.ilemx(i,j)) then ! flottant ou ile |
---|
151 | betamx(i,j)=0. |
---|
152 | else ! grounded, SIA |
---|
153 | betamx(i,j)=tostick ! frottement glace posee (1 bar) |
---|
154 | endif |
---|
155 | |
---|
156 | end do |
---|
157 | end do gzm_x |
---|
158 | |
---|
159 | gzm_y: do j=2,ny !-1 |
---|
160 | do i=1,nx ! 2,nx |
---|
161 | |
---|
162 | ! gzmy pour les hauteurs d'eau superieures a un seuil |
---|
163 | ! attention, il semble qu'il y avait une erreur dans la precedente |
---|
164 | ! formulation car le test etait seulement sur hwater(i,j). |
---|
165 | |
---|
166 | ! mais gzmx peut être quand même vrai selon une conditon sur la pression effective |
---|
167 | gzmy(i,j)=gzmy(i,j).or.((Neffmy(i,j).lt.neffratio*hmy(i,j)) & |
---|
168 | .and..not.flotmy(i,j)) |
---|
169 | |
---|
170 | ! ce cas est un ilemy traité plus loin |
---|
171 | ! if (j.gt.2.AND.j.lt.ny) gzmy(i,j)=gzmy(i,j).or.(flot(i,j).and. & |
---|
172 | ! .not.flot(i,j-1).and.flot(i,j+1).and.flot(i,j-2)) |
---|
173 | |
---|
174 | ! On rajoute une condition sur la contrainte basale |
---|
175 | ! si tobmy(i,j) > toblim (de l'ordre de 1 bar), passage a pose |
---|
176 | gzmy(i,j)=gzmy(i,j).and.(abs(tobmy(i,j)).lt.toblim) |
---|
177 | gzmy(i,j)=gzmy(i,j).or.((hwater(i,j).gt.hwatstream).and..not.flotmy(i,j)) |
---|
178 | |
---|
179 | ! calcul du frottement basal (ce bloc etait avant dans neffect) |
---|
180 | if (mkysedim(i,j).eq.2) then ! loi sediment |
---|
181 | if (hw_my(i,j).ge.seuil_sedim) then |
---|
182 | if (.not.flot(i,j)) then |
---|
183 | gzmy(i,j)=.true. |
---|
184 | coefslid=hw_my(i,j) |
---|
185 | coefslid=max(coefslid,seuil_sedim) |
---|
186 | coefslid=min(coefslid/coefmax,1.) |
---|
187 | betamy(i,j)=min(coefdrag/coefslid,1.e5) |
---|
188 | flgzmy(i,j)=.true. |
---|
189 | mslid_my(i,j)=5 |
---|
190 | shelfy=.true. |
---|
191 | endif |
---|
192 | endif |
---|
193 | else if (gzmy(i,j)) then ! stream |
---|
194 | betamy(i,j)=cf*neffmy(i,j) |
---|
195 | else if (flotmy(i,j).or.ilemy(i,j)) then ! flottant ou ile |
---|
196 | betamy(i,j)=0. |
---|
197 | else ! grounded, SIA |
---|
198 | betamy(i,j)=tostick ! frottement glace posee |
---|
199 | endif |
---|
200 | |
---|
201 | end do |
---|
202 | end do gzm_y |
---|
203 | |
---|
204 | return |
---|
205 | end subroutine dragging |
---|
206 | |
---|
207 | end module dragging_hwatstream_sedim_hudson |
---|
208 | |
---|