1 | !> \file sliding_Bindshadler_mod.f90 |
---|
2 | !! Module pour definir la loi de glissement Binshadler |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace sliding_Bindschadler |
---|
6 | !! Module pour definir la loi de glissement Binshadler |
---|
7 | !! \author CatRitz et al. |
---|
8 | !! \date 2001 |
---|
9 | !! @note Used module |
---|
10 | !! @note - use module3D_phy |
---|
11 | !! @note - use param_phy_mod |
---|
12 | !< |
---|
13 | |
---|
14 | module sliding_Bindschadler |
---|
15 | ! loi de glissement Binshadler modifie Ritz (Ritz et al. 2001) |
---|
16 | |
---|
17 | USE module3D_phy |
---|
18 | USE param_phy_mod |
---|
19 | |
---|
20 | implicit none |
---|
21 | real :: coefneff=1 |
---|
22 | real :: kweert ! coefficient de la loi de Weertman simple |
---|
23 | |
---|
24 | contains |
---|
25 | !------------------------------------------------------------------------------- |
---|
26 | subroutine init_sliding |
---|
27 | |
---|
28 | implicit none |
---|
29 | |
---|
30 | namelist/slid_bindsh/kweert,loigliss,coefbmax |
---|
31 | |
---|
32 | ! loigliss est declare dans 3D-... pourrait etre dans les divers modules de sliding |
---|
33 | ! coefbmax est declare dans 3D-.. facteur de normalisation pour l'influence de l'eau |
---|
34 | |
---|
35 | |
---|
36 | if (itracebug.eq.1) call tracebug(' Bindschadler entree dans routine init_sliding') |
---|
37 | |
---|
38 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
39 | read(num_param,slid_bindsh) |
---|
40 | |
---|
41 | write(num_rep_42,*)'!___________________________________________________________' |
---|
42 | write(num_rep_42,*)'! glissement module sliding_Bindschadler ' |
---|
43 | write(num_rep_42,*) |
---|
44 | write(num_rep_42,slid_bindsh) |
---|
45 | write(num_rep_42,*)'! kweert : coefficent, loigliss le type de loi' |
---|
46 | write(num_rep_42,*)'! coefbmax : facteur de normalisation pour influence eau' |
---|
47 | write(num_rep_42,*)'!___________________________________________________________' |
---|
48 | |
---|
49 | |
---|
50 | coefbmax=min(coefbmax,hwatermax) |
---|
51 | |
---|
52 | |
---|
53 | ! tableaux pour le glissement |
---|
54 | mslid_mx(:,:)=1 |
---|
55 | mslid_my(:,:)=1 |
---|
56 | slid_mx(:,:)=1. |
---|
57 | slid_my(:,:)=1. |
---|
58 | |
---|
59 | |
---|
60 | return |
---|
61 | end subroutine init_sliding |
---|
62 | |
---|
63 | |
---|
64 | !________________________________________________________________________________ |
---|
65 | subroutine sliding |
---|
66 | |
---|
67 | implicit none |
---|
68 | real :: neffmin |
---|
69 | |
---|
70 | ! cette routine calcule le terme ddbx et ddby dans le cas Bindshadler modifie |
---|
71 | ! loigliss=2 |
---|
72 | |
---|
73 | ! ------------------------------------------------- |
---|
74 | ! GLISSEMENT |
---|
75 | ! DDBX et DDBY termes dus au glissement |
---|
76 | ! relation avec la vitesse de glissement UXB et UYB |
---|
77 | ! UXB=-DDBX*SDX et UYB=-DDBY*SDY |
---|
78 | ! ------------------------------------------------- |
---|
79 | |
---|
80 | |
---|
81 | if (itracebug.eq.1) call tracebug(' Entree dans routine sliding Bindshadler') |
---|
82 | |
---|
83 | !coefmax=1. |
---|
84 | ddby(:,:) = 0. |
---|
85 | ddbx(:,:) = 0. |
---|
86 | |
---|
87 | slidy: do j=2,NY |
---|
88 | do I=2,NX-1 |
---|
89 | |
---|
90 | ! write(6,*) 'boucle y, i,j',i,j |
---|
91 | |
---|
92 | pose_y: if (.not.flgzmy(i,j)) then |
---|
93 | |
---|
94 | finey: if (HMY(I,J).lt.1.) then ! glace tres fine |
---|
95 | ddby(i,j)=0. |
---|
96 | |
---|
97 | else |
---|
98 | |
---|
99 | ! le glissement depend de l'eau basale. |
---|
100 | ! ------------------------------------- |
---|
101 | ! la pression effective a ete calculee dans Neffect |
---|
102 | ! Neffmy sur la demi maille y |
---|
103 | ! coefmybmelt est calculee dans neffect et correspond a |
---|
104 | ! la hauteur d'eau sur les demi mailles. |
---|
105 | |
---|
106 | |
---|
107 | sec_y: if ((COEFMYBMELT(i,j).le.0.).and.(.not.gzmy(i,j))) then |
---|
108 | DDBY(I,J)=0. |
---|
109 | ! write(6,*) '2',i,j,COEFMYBMELT(i,j) |
---|
110 | else |
---|
111 | neffmin=max(Neffmy(i,j),0.01) |
---|
112 | neffmin=max(Neffmy(i,j),1.e6) |
---|
113 | coefneff=min(abs(rog*hmy(i,j)/neffmin),1/0.0001) |
---|
114 | coefneff=min(coefneff,50.) |
---|
115 | coefneff=min(coefneff,5.) |
---|
116 | |
---|
117 | if ((Neffmy(i,j).gt.rog*2000.).and. & |
---|
118 | (COEFMYBMELT(i,j).lt.coefbmax)) then |
---|
119 | |
---|
120 | ! introduction d'un coefficient supplementaire |
---|
121 | ! pour avoir une limitation liee a la fusion basale |
---|
122 | |
---|
123 | COEFNeff=coefNeff*min(1.,COEFMYBMELT(i,j)/COEFBMAX) |
---|
124 | ! write(6,*) '3',i,j,coefneff |
---|
125 | endif |
---|
126 | |
---|
127 | ddby(i,j)=kweert*coefneff*(rog*hmy(i,j))**2.*sqrt(slope2my(i,j)) |
---|
128 | |
---|
129 | end if sec_y |
---|
130 | end if finey |
---|
131 | end if pose_y |
---|
132 | |
---|
133 | end do |
---|
134 | end do slidy |
---|
135 | |
---|
136 | |
---|
137 | |
---|
138 | slidx: do j=2,NY-1 |
---|
139 | do I=2,NX |
---|
140 | |
---|
141 | ! write(6,*) 'boucle x, i,j',i,j |
---|
142 | |
---|
143 | pose_x: if (.not.flgzmx(i,j)) then |
---|
144 | |
---|
145 | finex: if (HMX(I,J).lt.1.) then ! glace tres fine |
---|
146 | ddbx(i,j)=0. |
---|
147 | |
---|
148 | else |
---|
149 | |
---|
150 | ! le glissement depend de l'eau basale. |
---|
151 | ! ------------------------------------- |
---|
152 | ! la pression effective a ete calculee dans Neffect |
---|
153 | ! Neffmy sur la demi maille y |
---|
154 | ! coefmybmelt est calculee dans neffect et correspond a |
---|
155 | ! la hauteur d'eau sur les demi mailles. |
---|
156 | |
---|
157 | |
---|
158 | sec_x: if ((COEFMXBMELT(i,j).le.0.).and.(.not.gzmx(i,j))) then |
---|
159 | DDBX(I,J)=0. |
---|
160 | |
---|
161 | else |
---|
162 | |
---|
163 | neffmin=max(Neffmx(i,j),0.01) |
---|
164 | neffmin=max(Neffmx(i,j),1.e6) |
---|
165 | |
---|
166 | coefneff=min(abs(rog*hmx(i,j)/neffmin),1/0.0001) |
---|
167 | coefneff=min(coefneff,50.) |
---|
168 | coefneff=min(coefneff,5.) |
---|
169 | |
---|
170 | if ((Neffmx(i,j).gt.rog*2000.).and. & |
---|
171 | (COEFMXBMELT(i,j).lt.coefbmax)) then |
---|
172 | |
---|
173 | ! introduction d'un coefficient supplementaire |
---|
174 | ! pour avoir une limitation liee a la fusion basale |
---|
175 | |
---|
176 | COEFNeff=coefNeff*min(1.,COEFMXBMELT(i,j)/COEFBMAX) |
---|
177 | |
---|
178 | end if |
---|
179 | |
---|
180 | ddbx(i,j)=kweert*coefneff*(rog*hmx(i,j))**2.*sqrt(slope2mx(i,j)) |
---|
181 | |
---|
182 | end if sec_x |
---|
183 | end if finex |
---|
184 | |
---|
185 | end if pose_x |
---|
186 | |
---|
187 | end do |
---|
188 | end do slidx |
---|
189 | |
---|
190 | |
---|
191 | |
---|
192 | return |
---|
193 | end subroutine sliding |
---|
194 | |
---|
195 | end module sliding_bindschadler |
---|
196 | |
---|