-
Notifications
You must be signed in to change notification settings - Fork 197
/
Copy pathexample17.py
214 lines (172 loc) · 5.19 KB
/
example17.py
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
"""
2D flow in a channel around a circular cylinder
Partially taken from cylinder problem, Chapter 21, FEniCS book
Uses Chorin & Temam's solution method
"""
from dolfin import *
# Constants related to the geometry
bmarg = 1.e-3 + DOLFIN_EPS
xmin = 0.0
xmax = 2.2
ymin = 0.0
ymax = 0.41
xcenter = 0.2
ycenter = 0.2
radius = 0.05
# generate coarse mesh, refine only near cylinder
domain = Rectangle(xmin,ymin,xmax,ymax) \
- Circle(xcenter,ycenter,radius,12)
mesh = Mesh(domain, 22)
# refine mesh twice around the cylinder
cylCenter = Point(xcenter, ycenter)
for refinements in [0,1]:
cell_markers = CellFunction("bool", mesh)
cell_markers.set_all(False)
for cell in cells(mesh):
p = cell.midpoint()
if (xcenter/2. < p[0] < (xmax - xcenter)/2.) and \
(ycenter/4. < p[1] < ymax - ycenter/4.) :
cell_markers[cell] = True
mesh = refine(mesh, cell_markers)
plot(mesh)
# timestepping
dt = .00125
endTime = 3. - 1.e-5
# kinematic viscosity
nu = .001
# boundary conditions using a mesh function
boundaries = MeshFunction("size_t", mesh, mesh.topology().dim()-1)
# Inflow boundary
class InflowBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and x[0] < xmin + bmarg
# No-slip boundary
class NoslipBoundary(SubDomain):
def inside(self, x, on_boundary):
dx = x[0] - xcenter
dy = x[1] - ycenter
r = sqrt(dx*dx + dy*dy)
return on_boundary and \
(x[1] < ymin + bmarg or x[1] > ymax - bmarg or \
r < radius + bmarg)
# Outflow boundary
class OutflowBoundary(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and x[0] > xmax - bmarg
# uncomment for more detail, or use PROGRESS
#set_log_level(DEBUG)
# Define function spaces
V = VectorFunctionSpace(mesh, "CG", 2)
Q = FunctionSpace(mesh, "CG", 1)
W = V * Q
# no-slip velocity b.c.
noslipBoundary = NoslipBoundary()
g0 = Constant( (0.,0.) )
bc0 = DirichletBC(W.sub(0), g0, noslipBoundary)
bc0u = DirichletBC(V, g0, noslipBoundary)
# inlet velocity b.c.
inflowBoundary = InflowBoundary()
g1 = Expression( ("4.*Um*(x[1]*(ymax-x[1]))/(ymax*ymax)" , "0.0"), \
Um=1.5, ymax=ymax)
bc1 = DirichletBC(W.sub(0), g1, inflowBoundary)
bc1u = DirichletBC(V, g1, inflowBoundary)
# outflow pressure b.c.
outflowBoundary = OutflowBoundary()
g2 = Constant(0.)
bc2 = DirichletBC(W.sub(1), g2, outflowBoundary)
bc2p = DirichletBC(Q, g2, outflowBoundary)
# outflow velocity b.c., same as inlet
bc3 = DirichletBC(W.sub(0), g1, outflowBoundary)
bc3u = DirichletBC(V, g1, outflowBoundary)
# collect b.c.
bcs = [bc0, bc1, bc2, bc3]
# functions
(us, ps) = TrialFunctions(W)
(vs, qs) = TestFunctions(W)
# weak form Stokes equation
Stokes = (inner(grad(us), grad(vs)) - div(vs)*ps + qs*div(us))*dx
f = Constant((0., 0.))
LStokes = inner(f, vs)*dx
# initial condition comes from Stokes equation
w = Function(W)
solve(Stokes == LStokes, w, bcs,solver_parameters=dict(linear_solver="lu"))
# Split the mixed solution using deepcopy
(uinit, pinit) = w.split(True)
umax = max(abs(uinit.vector().array()))
print "Worst possible Courant number, initial velocity=",dt*umax/mesh.hmin()
# Timestepping method starts here
step = 0
t=0
u0 = Function(V)
u0.assign(uinit)
utilde = Function(V)
u1 = Function(V)
p1 = Function(Q)
# Define test and trial functions
v = TestFunction(V)
q = TestFunction(Q)
u = TrialFunction(V)
p = TrialFunction(Q)
# STEP 1 (u will be utilde)
F1 = (1./dt)*inner(v, u - u0)*dx + inner(v, grad(u0)*u0)*dx \
+ nu*inner(grad(u), grad(v))*dx - inner(f, v)*dx
a1 = lhs(F1)
L1 = rhs(F1)
# STEP 2 (p)
a2 = inner(grad(q), grad(p))*dx
L2 = -(1./dt)*q*div(utilde)*dx
# STEP 3 (Velocity update)
a3 = inner(v, u)*dx
L3 = inner(v, utilde)*dx - dt*inner(v, grad(p1))*dx
# Assemble matrices
# cannot use symmetric b.c.!
A1 = assemble(a1)
bc0u.apply(A1)
bc1u.apply(A1)
bc3u.apply(A1)
A2 = assemble(a2)
bc2p.apply(A2)
A3 = assemble(a3)
bc0u.apply(A3)
bc1u.apply(A3)
bc3u.apply(A3)
# set up 3 quiet solvers for fixed matrices
solver1 = LUSolver(A1)
solver1.parameters["reuse_factorization"] = True
solver1.parameters["report"] = False
solver2 = LUSolver(A2)
solver2.parameters["reuse_factorization"] = True
solver2.parameters["report"] = False
solver3 = LUSolver(A3)
solver3.parameters["reuse_factorization"] = True
solver3.parameters["report"] = False
# timestepping
while t < endTime * (1. + 1.e-10):
t += dt
step += 1
# STEP 1 solution: solve for utilde
b = assemble(L1)
bc0u.apply(b)
bc1u.apply(b)
bc3u.apply(b)
solver1.solve( utilde.vector(), b )
# STEP 2 solution: solve for p1
b = assemble(L2)
bc2p.apply(b)
solver2.solve( p1.vector(), b )
# STEP 3 solution: solve for u1
b = assemble(L3)
bc0u.apply(b)
bc1u.apply(b)
bc3u.apply(b)
solver3.solve( u1.vector(), b )
# prepare for next time step
u0.assign(u1)
plot(u0,title="Velocity",rescale=False)
# print a little information each step
print "t=%f, max abs u=%e, max p=%e, v(.5,0)=%e" \
% (t, max(abs(u1.vector().array())), \
max(p1.vector().array()), u1([.35,0.2])[1] )
umax = max(abs(u1.vector().array()))
print "Worst possible Courant number=",dt*umax/mesh.hmin()
interactive()