-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathStokesCylindrical.jl
319 lines (298 loc) · 12.1 KB
/
StokesCylindrical.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
function form_stokes_cylindrical(grid::CartesianGrid,eta_s::Matrix,eta_n::Matrix,eta_vx::Matrix,eta_vy::Matrix,rhoX::Matrix,rhoY::Matrix,bc::BoundaryConditions,gx::Float64,gy::Float64)
# Form the Stokes system.
# For cylindrical axisymmetric coordinates.
# note that vx and x corresponds to the radial velocity and radius
# y and vy corresponds to the axial coordinate (up/down).
#
# Inputs:
# grid - the cartesian grid
# eta_s - viscosity at the basic nodes
# eta_n - viscosity at the cell centers
# eta_vx - viscosity at vx nodes
# eta_vy - viscosity at vy nodes
# rhoX - density at the vx nodes
# rhoY - density at the vy nodes
# bc - a vector describing the boundary conditions along the [left,right,top,bottom]
# gx,gy - gravitational body force in the x and y direction
# Outputs:
# L,R - the left hand side (matrix) and right hand side (vector) of the stokes system
k::Int64 = 1 # index into dof arrays
nx = grid.nx
ny = grid.ny
nn = nx*ny
nnz = 2*11*nn + 5*nn # total number of nonzeros in matrix (not including BCs)
row_index = zeros(Int64,nnz) # up to 5 nonzeros per row
col_index = zeros(Int64,nnz)
value = zeros(Float64, nnz)
dx = grid.W/(grid.nx-1)
dy = grid.H/(grid.ny-1)
kcont = 2*eta_s[1,1]/(dx+dy)# scaling factor for continuity equation
kcont = 1e20/(dx+dy)*2
kbond = 1.# scaling factor for dirichlet bc equations.
R=zeros(3*nn,1)
# loop over j
for j in 1:nx
# loop over i
for i in 1:ny
# dxp is the dx in the +x direction, dxm is dx in the -x direction, dxc is the spacing between cell centers
dxp = j<nx ? grid.x[j+1] - grid.x[j] : grid.x[j] - grid.x[j-1]
dxm = j>1 ? grid.x[j] - grid.x[j-1] : grid.x[j+1] - grid.x[j]
dxc = 0.5*(dxp+dxm)
# dyp and dym are spacing between vx nodes in the +y and -y directions
dyp = i<ny ? grid.yc[i+1]- grid.yc[i] : grid.yc[i] -grid.yc[i-1]
dym = i>1 ? grid.yc[i] - grid.yc[i-1] : grid.yc[i+1] -grid.yc[i]
dyc = 0.5*(dyp+dym)
# discretize the x-stokes - note that numbering in comments refers to Gerya Figure 7.18a
# and equation 7.22
this_row = vxdof(i,j,ny)
# Boundary cases first...
if j==1 || j == nx # left boundary or right boundary
# vx = 0
row_index[k] = this_row
col_index[k] = this_row
value[k] = kbond
k+=1
R[this_row] = 0.0 *kbond
elseif i==1 # top
if grid.yperiodic
# dvx/dy = 0 (free slip)
row_index[k] = this_row
col_index[k] = this_row
value[k] = -kbond
k+=1
row_index[k] = this_row
col_index[k] = vxdof(grid.ny,j,ny)
value[k] = kbond
k+=1
R[this_row] = 0.0*kbond
else
# dvx/dy = 0 (free slip)
row_index[k] = this_row
col_index[k] = this_row
value[k] = -kbond
k+=1
row_index[k] = this_row
col_index[k] = vxdof(i+1,j,ny)
value[k] = kbond
k+=1
R[this_row] = 0.0*kbond
end
else
# vx1 vx(i,j-1)
row_index[k] = this_row
col_index[k] = vxdof(i,j-1,ny)
value[k] = 2*eta_n[i,j]/dxm/dxc - 2.0/grid.x[j]*eta_vx[i,j]/(dxp+dxm)
k+=1
# vx2 vx(i-1,j)
row_index[k] = this_row
col_index[k] = vxdof(i-1,j,ny)
value[k] = eta_s[i-1,j]/dym/dyc
k+=1
# vx3 vx(i,j)
row_index[k] = this_row
col_index[k] = this_row
value[k] = -2*eta_n[i,j+1]/dxp/dxc -2*eta_n[i,j]/dxm/dxc - eta_s[i,j]/dyp/dyc - eta_s[i-1,j]/dym/dyc - 2.0/grid.x[j]^2*eta_vx[i,j]
if i == ny #vx4 vx(i+1,j)
# if i == nx, dvx/dy = 0 -> vx3 == vx4 (see Gerya fig 7.18a)
value[k] += eta_s[i,j]/dyp/dyc
k+=1
else
k+=1
# vx4
# enforce dvx/dy = 0 (free slip)
row_index[k] = this_row
col_index[k] = vxdof(i+1,j,ny)
value[k] = eta_s[i,j]/dyp/dyc
k+=1
end
# vx5 vx(i,j+1)
row_index[k] = this_row
col_index[k] = vxdof(i,j+1,ny)
value[k] = 2*eta_n[i,j+1]/dxp/dxc + 2.0/grid.x[j]*eta_vx[i,j]/(dxp+dxm)
k+=1
# vy1 vy(i-1,j)
row_index[k] = this_row
col_index[k] = vydof(i-1,j,ny)
value[k] = eta_s[i-1,j]/dxc/dyc
k+=1
# vy2
row_index[k] = this_row
col_index[k] = vydof(i,j,ny)
value[k] = -eta_s[i,j]/dxc/dyc
k+=1
# vy3
row_index[k] = this_row
col_index[k] = vydof(i-1,j+1,ny)
value[k] = -eta_s[i-1,j]/dxc/dyc
k+=1
# vy4
row_index[k] = this_row
col_index[k] = vydof(i,j+1,ny)
value[k] = eta_s[i,j]/dxc/dyc
k+=1
# P1
row_index[k] = this_row
col_index[k] = pdof(i,j,ny)
value[k] = kcont/dxc
k+=1
# P2
row_index[k] = this_row
col_index[k] = pdof(i,j+1,ny)
value[k] = -kcont/dxc
k+=1
R[this_row] = -gx*rhoX[i,j]
end
# END X-STOKES
# BEGIN Y-STOKES
dxp = j < nx ? grid.xc[j+1] - grid.xc[j] : grid.xc[j] -grid.xc[j-1]
dxm = j > 1 ? grid.xc[j] - grid.xc[j-1] : grid.xc[j+1]-grid.xc[j]
dxc = j > 1 ? grid.x[j] - grid.x[j-1] : grid.x[j+1] - grid.x[j]
dyp = i < ny ? grid.y[i+1] - grid.y[i] : grid.y[i] - grid.y[i-1]
dym = i > 1 ? grid.y[i] - grid.y[i-1] : grid.y[i+1] - grid.y[i]
dyc = i < ny ? grid.yc[i+1] - grid.yc[i] : grid.yc[i] - grid.yc[i-1]
this_row = vydof(i,j,ny)
if (!grid.yperiodic && i==1) || i == ny
# top row / bottom row
row_index[k] = this_row
col_index[k] = this_row
value[k] = kbond
k+=1
R[this_row] = 0.0*kbond
elseif j==1
# left boundary - free slip
row_index[k] = this_row
col_index[k] = this_row
value[k] = kbond
k+=1
row_index[k] = this_row
col_index[k] = vydof(i,j+1,ny)
value[k] = -kbond
k+=1
R[this_row] = 0.0*kbond
else
#vy1 vy(i,j-1)
row_index[k] = this_row
col_index[k] = vydof(i,j-1,ny)
value[k] = eta_s[i,j-1]/dxm/dxc - eta_vy[i,j]/grid.xc[j]/(dxp+dxm)
k+=1
#vy2 vy(i-1,j)
row_index[k] = this_row
col_index[k] = vydof(i-1,j,ny)
value[k] = 2*eta_n[i,j]/dym/dyc
k+=1
#vy3
row_index[k] = this_row
col_index[k] = this_row
value[k] = -2*eta_n[i+1,j]/dyp/dyc -2*eta_n[i,j]/dym/dyc - eta_s[i,j]/dxp/dxc - eta_s[i,j-1]/dxm/dxc
if j == nx
# free slip - vy5 = vx3.
value[k] += eta_s[i,j]/dxp/dxc + eta_vy[i,j]/grid.xc[j]/(dxp+dxm)
end
k+=1
#vy4 vy(i+1,j)
row_index[k] = this_row
col_index[k] = vydof(i+1,j,ny)
value[k] = 2*eta_n[i+1,j]/dyp/dyc
k+=1
#vy5 vy(i,j+1)
if j<nx
row_index[k] = this_row
col_index[k] = vydof(i,j+1,ny)
value[k] = eta_s[i,j]/dxp/dxc + eta_vy[i,j]/grid.xc[j]/(dxp+dxm)
k+=1
end
#vx1 vx(i,j-1)
row_index[k] = this_row
col_index[k] = vxdof(i,j-1,ny)
value[k] = eta_s[i,j-1]/dxc/dyc - eta_vy[i,j]/grid.xc[j]/dyc/2.0
k+=1
#vx2 vx(i+1,j-1)
row_index[k] = this_row
col_index[k] = vxdof(i+1,j-1,ny)
value[k] = -eta_s[i,j-1]/dxc/dyc + eta_vy[i,j]/grid.xc[j]/dyc/2.0
k+=1
#vx3 vx(i,j)
row_index[k] = this_row
col_index[k] = vxdof(i,j,ny)
value[k] = -eta_s[i,j]/dxc/dyc - eta_vy[i,j]/grid.xc[j]/dyc/2.0
k+=1
#vx4 vx(i+1,j)
row_index[k] = this_row
col_index[k] = vxdof(i+1,j,ny)
value[k] = eta_s[i,j]/dxc/dyc + eta_vy[i,j]/grid.xc[j]/dyc/2.0
k+=1
#P1
row_index[k] = this_row
col_index[k] = pdof(i,j,ny)
value[k] = kcont/dyc
k+=1
#P2
row_index[k] = this_row
col_index[k] = pdof(i+1,j,ny)
value[k] = -kcont/dyc
k+=1
R[this_row] = -gy*rhoY[i,j]
end
# END Y-STOKES
# discretize the continuity equation
# dvx/dx + dvy/dy = 0
this_row = pdof(i,j,ny)
if i==1 || j == 1 || (i==2 && j == 2)
row_index[k] = this_row
col_index[k] = this_row
value[k] = kbond
k+=1
R[this_row] = 0.0
else
dxm = grid.x[j] - grid.x[j-1]
dym = grid.y[i] - grid.y[i-1]
xc = grid.xc[j]
row_index[k] = this_row
col_index[k] = vxdof(i,j,ny)
value[k] = kcont/dxm*(grid.x[j]/grid.xc[j])
k+=1
row_index[k] = this_row
col_index[k] = vxdof(i,j-1,ny)
value[k] = -kcont/dxm*(grid.x[j-1]/grid.xc[j])
k+=1
row_index[k] = this_row
col_index[k] = vydof(i,j,ny)
value[k] = kcont/dym
k+=1
row_index[k] = this_row
col_index[k] = vydof(i-1,j,ny)
value[k] = -kcont/dym
k+=1
row_index[k] = this_row
col_index[k] = this_row
value[k] = 0.0
k+=1
R[this_row] = 0.0
end
# END CONTINUITY
end
end
@views row_index = row_index[1:(k-1)]
@views col_index = col_index[1:(k-1)]
@views value = value[1:(k-1)]
L = sparse(row_index,col_index,value)
return L,R
end
function compute_velocity_divergence(grid::CartesianGrid,vx::Matrix{Float64},vy::Matrix{Float64})
# This function computes the divergence of the velocity field.
# I wrote the function only to verify that the cylindrical version of the continuity equation is satisfied.
# inputs - grid, velocities at the velocity nodes.
nx = grid.nx
ny = grid.ny
divv = zeros(ny-1,nx-1)
# compute cell-by-cell velocity divergence
# 1/r d/dr(r vr)
for i in 2:ny
for j in 2:nx
dvxdx = 1.0/grid.xc[j]/(grid.x[j]-grid.x[j-1])*(grid.x[j]*vx[i,j] - grid.x[j-1]*vx[i,j-1])
dvydy = (vy[i,j]-vy[i-1,j])/(grid.y[i]-grid.y[i-1])
divv[i-1,j-1] = dvxdx + dvydy
end
end
return divv
end