1 | !> \file dragging_coulomb_friction_mod.f90 |
---|
2 | !! Coulomb friction for a given value of the non-linearity |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace dragging_coulomb_friction |
---|
6 | !! Coulomb friction as defined in Pattyn,TC,2017 and reference therein |
---|
7 | !! \author Aurelien Quiquet |
---|
8 | !! \date 28/01/2021 |
---|
9 | !! @note Non-linearity is provided in the param file. |
---|
10 | !! @note * m should vary from 1 to infinity |
---|
11 | !! @note * m=1 (q=1) -> linear friction as in Quiquet et al. (2018) |
---|
12 | !! @note * m=infinity (q=0) -> -1 in param_list |
---|
13 | !! @note * m=3 (q=1/3) is a classic (Weertman, 1974; Schoof, 2007) |
---|
14 | !! @note Used module |
---|
15 | !! @note - use module3D_phy |
---|
16 | !! @note - use param_phy_mod |
---|
17 | !! @note - use fake_beta_iter_vitbil_mod |
---|
18 | !< |
---|
19 | |
---|
20 | module dragging_coulomb_friction |
---|
21 | |
---|
22 | !______________________________________________________________________________ |
---|
23 | ! |
---|
24 | ! I took the formulation of Pattyn in his The Cryos. 2017 paper. It looks like: |
---|
25 | ! |
---|
26 | ! taub = [ tauc / ( |ub|^(1-q) uO^q) ] ub |
---|
27 | ! |
---|
28 | ! with: -tauc is our cf x N (expressed differently in Pattyn). |
---|
29 | ! u0=100 m/yr |
---|
30 | ! |
---|
31 | ! What is inside [] is what we used as -Beta in Quiquet et al. (2018). |
---|
32 | ! |
---|
33 | ! For q=1, we get back our standard formulation (except that cf'=cf/100) |
---|
34 | ! |
---|
35 | ! In the param_list we will give: |
---|
36 | ! - cf as before |
---|
37 | ! - m=1/q (if m=-1 is provided -> q=0) |
---|
38 | ! - number of iterations (as it is non-linear) |
---|
39 | ! - betamin and betamax |
---|
40 | !______________________________________________________________________________ |
---|
41 | |
---|
42 | |
---|
43 | use module3d_phy, only: nx,ny |
---|
44 | use fake_beta_iter_vitbil_mod |
---|
45 | |
---|
46 | implicit none |
---|
47 | |
---|
48 | real :: betamin ! betamin : (Pa) frottement mini sous les streams |
---|
49 | |
---|
50 | real :: cf ! A friction coefficient |
---|
51 | real :: m_nolin ! Schoof non-linear exponent |
---|
52 | real :: q_nolin ! q=1/m |
---|
53 | |
---|
54 | real, dimension(nx,ny) :: Vcol_x !< uniquement pour compatibilite avec spinup cat |
---|
55 | real, dimension(nx,ny) :: Vcol_y !< uniquement pour compatibilite avec spinup cat |
---|
56 | real, dimension(nx,ny) :: Vsl_x !< uniquement pour compatibilite avec spinup cat |
---|
57 | real, dimension(nx,ny) :: Vsl_y !< uniquement pour compatibilite avec spinup cat |
---|
58 | logical :: corr_def = .false. !< for deformation correction, pour compatibilite beta |
---|
59 | |
---|
60 | contains |
---|
61 | !------------------------------------------------------------------------------- |
---|
62 | |
---|
63 | !> SUBROUTINE: init_dragging |
---|
64 | !! Cette routine fait l'initialisation du dragging. |
---|
65 | !> |
---|
66 | subroutine init_dragging ! Cette routine fait l'initialisation du dragging. |
---|
67 | |
---|
68 | use module3d_phy, only: niter_nolin,betamax,betamax_2d,inv_beta,mstream_mx,mstream_my,drag_mx,drag_my,num_rep_42 |
---|
69 | use runparam, only: itracebug |
---|
70 | |
---|
71 | implicit none |
---|
72 | |
---|
73 | namelist/drag_coulomb_friction/cf,m_nolin,niter_nolin,betamax,betamin |
---|
74 | |
---|
75 | if (itracebug.eq.1) call tracebug(' dragging avec neff et slope') |
---|
76 | |
---|
77 | |
---|
78 | ! formats pour les ecritures dans 42 |
---|
79 | 428 format(A) |
---|
80 | |
---|
81 | ! lecture des parametres du run block drag neff |
---|
82 | !-------------------------------------------------------------------- |
---|
83 | |
---|
84 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
85 | read(num_param,drag_coulomb_friction) |
---|
86 | |
---|
87 | write(num_rep_42,428)'!___________________________________________________________' |
---|
88 | write(num_rep_42,428) '&drag_coulomb_friction ! nom du bloc dragging coulomb friction' |
---|
89 | write(num_rep_42,*) |
---|
90 | write(num_rep_42,*) 'cf = ', cf |
---|
91 | write(num_rep_42,*) 'm_nolin = ', m_nolin |
---|
92 | write(num_rep_42,*) 'niter_nolin = ', niter_nolin |
---|
93 | write(num_rep_42,*) 'betamax = ', betamax |
---|
94 | write(num_rep_42,*) 'betamin = ', betamin |
---|
95 | write(num_rep_42,*)'/' |
---|
96 | write(num_rep_42,428) '! cf: a friction coefficient (to be tuned)' |
---|
97 | write(num_rep_42,428) '! m_nolin: non-linear exponent, from 1 to infinity (put -1 for infinity)' |
---|
98 | write(num_rep_42,428) '! m_nolin=1/q in Pattyn TC 2017, q in [0:1]' |
---|
99 | write(num_rep_42,428) '! niter_nolin: number of iterations to solve the non-linearity (expensive!)' |
---|
100 | |
---|
101 | inv_beta=0 |
---|
102 | !------------------------------------------------------------------- |
---|
103 | ! masque stream |
---|
104 | mstream_mx(:,:)=1 |
---|
105 | mstream_my(:,:)=1 |
---|
106 | |
---|
107 | |
---|
108 | ! coefficient permettant de modifier le basal drag. |
---|
109 | drag_mx(:,:)=1. |
---|
110 | drag_my(:,:)=1. |
---|
111 | |
---|
112 | betamax_2d(:,:) = betamax |
---|
113 | |
---|
114 | if (m_nolin.lt.0) then |
---|
115 | q_nolin = 0 |
---|
116 | else |
---|
117 | q_nolin = 1 / m_nolin |
---|
118 | endif |
---|
119 | |
---|
120 | return |
---|
121 | end subroutine init_dragging |
---|
122 | !________________________________________________________________________________ |
---|
123 | |
---|
124 | !> SUBROUTINE: dragging |
---|
125 | !! Calcul le frottement basal |
---|
126 | !> |
---|
127 | !------------------------------------------------------------------------- |
---|
128 | subroutine dragging ! defini la localisation des streams et le frottement basal |
---|
129 | !$ USE OMP_LIB |
---|
130 | |
---|
131 | use module3d_phy, only: nx,ny,betamax,beta_centre,betamx,betamy,neffmx,neffmy,hwater,flot,ux,uy |
---|
132 | |
---|
133 | implicit none |
---|
134 | |
---|
135 | real,dimension(nx,ny) :: tauc_mx ! Yield stress, x direction |
---|
136 | real,dimension(nx,ny) :: tauc_my ! Yield stress, y direction |
---|
137 | |
---|
138 | real,dimension(nx,ny) :: neff ! pression effective noeuds majeurs |
---|
139 | real,dimension(nx,ny) :: hwatmx ! hauteur d'eau staggered grille - afq |
---|
140 | real,dimension(nx,ny) :: hwatmy ! hauteur d'eau staggered grille - afq |
---|
141 | |
---|
142 | real,parameter :: u0 = 100d0 ! threshold sliding speed |
---|
143 | |
---|
144 | !$OMP PARALLEL |
---|
145 | |
---|
146 | ! calcul de neff (pression effective sur noeuds majeurs) |
---|
147 | if (sum(neffmx(:,:)).le.0.) neffmx(:,:) =1.e8 |
---|
148 | if (sum(neffmy(:,:)).le.0.) neffmy(:,:) =1.e8 |
---|
149 | |
---|
150 | !$OMP DO |
---|
151 | do j=1,ny-1 |
---|
152 | do i=1,nx-1 |
---|
153 | neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4 |
---|
154 | enddo |
---|
155 | enddo |
---|
156 | !$OMP END DO |
---|
157 | !aurel, for the borders: |
---|
158 | !$OMP WORKSHARE |
---|
159 | neff(:,ny)=1.e5 |
---|
160 | neff(nx,:)=1.e5 |
---|
161 | ! calcul de hwat (staggered) |
---|
162 | hwatmx(:,:)=0. |
---|
163 | hwatmy(:,:)=0. |
---|
164 | !$OMP END WORKSHARE |
---|
165 | !$OMP DO |
---|
166 | do i=2,nx |
---|
167 | hwatmx(i,:)=(hwater(i-1,:)+hwater(i,:))/2. |
---|
168 | enddo |
---|
169 | !$OMP END DO |
---|
170 | !$OMP DO |
---|
171 | do j=2,ny |
---|
172 | hwatmy(:,j)=(hwater(:,j-1)+hwater(:,j))/2. |
---|
173 | enddo |
---|
174 | !$OMP END DO |
---|
175 | |
---|
176 | !_____________________________ |
---|
177 | ! Coulomb friction computation: |
---|
178 | |
---|
179 | !$OMP WORKSHARE |
---|
180 | |
---|
181 | tauc_mx(:,:)= cf*neffmx(:,:) |
---|
182 | tauc_my(:,:)= cf*neffmy(:,:) |
---|
183 | |
---|
184 | where (abs(ux(:,:,nz)).gt.1) |
---|
185 | betamx(:,:) = tauc_mx(:,:) * ( abs(ux(:,:,nz))**(q_nolin-1.) ) / ( u0 ** q_nolin ) |
---|
186 | elsewhere |
---|
187 | betamx(:,:) = tauc_mx(:,:) / ( u0 ** q_nolin ) |
---|
188 | endwhere |
---|
189 | where (abs(uy(:,:,nz)).gt.1) |
---|
190 | betamy(:,:) = tauc_my(:,:) / ( abs(uy(:,:,nz))**(q_nolin-1.) ) / ( u0 ** q_nolin ) |
---|
191 | elsewhere |
---|
192 | betamy(:,:) = tauc_my(:,:) / ( u0 ** q_nolin ) |
---|
193 | endwhere |
---|
194 | |
---|
195 | betamx(:,:)=max(betamx(:,:),betamin) |
---|
196 | betamy(:,:)=max(betamy(:,:),betamin) |
---|
197 | betamx(:,:)=min(betamx(:,:),betamax) |
---|
198 | betamy(:,:)=min(betamy(:,:),betamax) |
---|
199 | |
---|
200 | where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax |
---|
201 | where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax |
---|
202 | |
---|
203 | !$OMP END WORKSHARE |
---|
204 | |
---|
205 | |
---|
206 | ! calcul de gzmx |
---|
207 | |
---|
208 | ! points flottants |
---|
209 | !$OMP DO |
---|
210 | do j=2,ny |
---|
211 | do i=2,nx |
---|
212 | if (flot(i,j).and.(flot(i-1,j))) then |
---|
213 | betamx(i,j)=0. |
---|
214 | end if |
---|
215 | if (flot(i,j).and.(flot(i,j-1))) then |
---|
216 | betamy(i,j)=0. |
---|
217 | end if |
---|
218 | end do |
---|
219 | end do |
---|
220 | !$OMP END DO |
---|
221 | |
---|
222 | !$OMP DO |
---|
223 | do j=2,ny-1 |
---|
224 | do i=2,nx-1 |
---|
225 | beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) & |
---|
226 | + (betamy(i,j)+betamy(i,j+1)))*0.25 |
---|
227 | end do |
---|
228 | end do |
---|
229 | !$OMP END DO |
---|
230 | !$OMP END PARALLEL |
---|
231 | |
---|
232 | return |
---|
233 | end subroutine dragging |
---|
234 | |
---|
235 | |
---|
236 | |
---|
237 | subroutine mstream_dragging ! defini la localisation des streams |
---|
238 | |
---|
239 | |
---|
240 | use module3d_phy, only: nx,ny,betamax,betamx,betamy,fleuvemx,fleuvemy,gzmx,gzmy,flgzmx,flgzmy,flot,flotmx,flotmy |
---|
241 | |
---|
242 | implicit none |
---|
243 | |
---|
244 | !$OMP PARALLEL |
---|
245 | |
---|
246 | !$OMP WORKSHARE |
---|
247 | fleuvemx(:,:)=.false. |
---|
248 | fleuvemy(:,:)=.false. |
---|
249 | gzmx(:,:)=.true. |
---|
250 | gzmy(:,:)=.true. |
---|
251 | flgzmx(:,:)=.false. |
---|
252 | flgzmy(:,:)=.false. |
---|
253 | !$OMP END WORKSHARE |
---|
254 | |
---|
255 | |
---|
256 | ! calcul de gzmx |
---|
257 | |
---|
258 | ! points flottants : flgzmx mais pas gzmx |
---|
259 | !$OMP DO |
---|
260 | do j=2,ny |
---|
261 | do i=2,nx |
---|
262 | if (flot(i,j).and.(flot(i-1,j))) then |
---|
263 | flotmx(i,j)=.true. |
---|
264 | flgzmx(i,j)=.true. |
---|
265 | gzmx(i,j)=.false. |
---|
266 | betamx(i,j)=0. |
---|
267 | end if |
---|
268 | if (flot(i,j).and.(flot(i,j-1))) then |
---|
269 | flotmy(i,j)=.true. |
---|
270 | flgzmy(i,j)=.true. |
---|
271 | gzmy(i,j)=.false. |
---|
272 | betamy(i,j)=0. |
---|
273 | end if |
---|
274 | end do |
---|
275 | end do |
---|
276 | !$OMP END DO |
---|
277 | |
---|
278 | !--------- autres criteres |
---|
279 | |
---|
280 | ! points poses |
---|
281 | !$OMP WORKSHARE |
---|
282 | gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax) ! Pas de calcul pour les points |
---|
283 | gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax) ! au fort beta |
---|
284 | |
---|
285 | flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:) |
---|
286 | flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:) |
---|
287 | |
---|
288 | fleuvemx(:,:)=gzmx(:,:) |
---|
289 | fleuvemy(:,:)=gzmy(:,:) |
---|
290 | !$OMP END WORKSHARE |
---|
291 | |
---|
292 | !$OMP END PARALLEL |
---|
293 | |
---|
294 | return |
---|
295 | end subroutine mstream_dragging |
---|
296 | |
---|
297 | |
---|
298 | |
---|
299 | |
---|
300 | |
---|
301 | |
---|
302 | end module dragging_coulomb_friction |
---|
303 | |
---|