This repository has been archived by the owner on Jan 13, 2024. It is now read-only.
/
infsup_stokes.py
214 lines (184 loc) · 8.28 KB
/
infsup_stokes.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
# Copyright (C) 2016-2022 by the multiphenics authors
#
# This file is part of multiphenics.
#
# multiphenics is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# multiphenics is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with multiphenics. If not, see <http://www.gnu.org/licenses/>.
#
from numpy import isclose
from dolfin import *
# import matplotlib.pyplot as plt
from multiphenics import *
"""
In this tutorial we compare the computation of the inf-sup constant
of a Stokes by standard FEniCS code and multiphenics code.
"""
# -------------------------------------------------- #
# MESH GENERATION #
# Create mesh
mesh = UnitSquareMesh(32, 32)
# Create boundaries
class Wall(SubDomain):
def inside(self, x, on_boundary):
return on_boundary and (x[1] < 0 + DOLFIN_EPS or x[1] > 1 - DOLFIN_EPS)
boundaries = MeshFunction("size_t", mesh, mesh.topology().dim() - 1)
boundaries.set_all(0)
wall = Wall()
wall.mark(boundaries, 1)
# Function spaces
V_element = VectorElement("Lagrange", mesh.ufl_cell(), 2)
Q_element = FiniteElement("Lagrange", mesh.ufl_cell(), 1)
def normalize(u1, u2, p):
u1.vector()[:] /= assemble(inner(grad(u1), grad(u1))*dx)
u1.vector().apply("insert")
u2.vector()[:] /= assemble(inner(grad(u2), grad(u2))*dx)
u2.vector().apply("insert")
p.vector()[:] /= assemble(p*p*dx)
p.vector().apply("insert")
# -------------------------------------------------- #
# STANDARD FEniCS FORMULATION BY FEniCS MixedElement #
def run_monolithic():
# Function spaces
W_element = MixedElement(V_element, Q_element)
W = FunctionSpace(mesh, W_element)
# Test and trial functions: monolithic
vq = TestFunction(W)
(v, q) = split(vq)
up = TrialFunction(W)
(u, p) = split(up)
# Variational forms
lhs = inner(grad(u), grad(v))*dx - div(v)*p*dx - div(u)*q*dx
rhs = - inner(p, q)*dx
# Boundary conditions
bc = DirichletBC(W.sub(0), Constant((0., 0.)), boundaries, 1)
# Assemble lhs and rhs matrices
LHS = assemble(lhs)
RHS = assemble(rhs)
# Solve
LHS = as_backend_type(LHS)
RHS = as_backend_type(RHS)
eigensolver = SLEPcEigenSolver(LHS, RHS, bc)
eigensolver.parameters["problem_type"] = "gen_non_hermitian"
eigensolver.parameters["spectrum"] = "target real"
eigensolver.parameters["spectral_transform"] = "shift-and-invert"
eigensolver.parameters["spectral_shift"] = 1.e-5
eigensolver.solve(1)
r, c = eigensolver.get_eigenvalue(0)
assert abs(c) < 1.e-10
assert r > 0., "r = " + str(r) + " is not positive"
print("Inf-sup constant (monolithic): ", sqrt(r))
# Extract eigenfunctions
r_fun, c_fun = Function(W), Function(W)
eigensolver.get_eigenpair(r_fun, c_fun, 0)
(u_fun, p_fun) = r_fun.split(deepcopy=True)
(u_fun_1, u_fun_2) = u_fun.split(deepcopy=True)
normalize(u_fun_1, u_fun_2, p_fun)
# plt.figure()
# plot(u_fun_1, title="Velocity 1 monolithic", mode="color")
# plt.figure()
# plot(u_fun_2, title="Velocity 2 monolithic", mode="color")
# plt.figure()
# plot(p_fun, title="Pressure monolithic", mode="color")
# plt.show()
return (r, u_fun_1, u_fun_2, p_fun)
(eig_m, u_fun_1_m, u_fun_2_m, p_fun_m) = run_monolithic()
# -------------------------------------------------- #
# multiphenics FORMULATION #
def run_block():
# Function spaces
W_element = BlockElement(V_element, Q_element)
W = BlockFunctionSpace(mesh, W_element)
# Test and trial functions
vq = BlockTestFunction(W)
(v, q) = block_split(vq)
up = BlockTrialFunction(W)
u, p = block_split(up)
# Variational forms
lhs = [[inner(grad(u), grad(v))*dx, - div(v)*p*dx],
[- div(u)*q*dx , 0 ]]
rhs = [[0 , 0 ],
[0 , - p*q*dx ]]
# Boundary conditions
wallc = DirichletBC(W.sub(0), Constant((0., 0.)), boundaries, 1)
bc = BlockDirichletBC([[wallc], []])
# Assemble lhs and rhs matrices
LHS = block_assemble(lhs)
RHS = block_assemble(rhs)
# Solve
eigensolver = BlockSLEPcEigenSolver(LHS, RHS, bc)
eigensolver.parameters["problem_type"] = "gen_non_hermitian"
eigensolver.parameters["spectrum"] = "target real"
eigensolver.parameters["spectral_transform"] = "shift-and-invert"
eigensolver.parameters["spectral_shift"] = 1.e-5
eigensolver.solve(1)
r, c = eigensolver.get_eigenvalue(0)
assert abs(c) < 1.e-10
assert r > 0., "r = " + str(r) + " is not positive"
print("Inf-sup constant (block): ", sqrt(r))
# Extract eigenfunctions
r_fun, c_fun = BlockFunction(W), BlockFunction(W)
eigensolver.get_eigenpair(r_fun, c_fun, 0)
(u_fun, p_fun) = r_fun.block_split()
(u_fun_1, u_fun_2) = u_fun.split(deepcopy=True)
normalize(u_fun_1, u_fun_2, p_fun)
# plt.figure()
# plot(u_fun_1, title="Velocity 1 block", mode="color")
# plt.figure()
# plot(u_fun_2, title="Velocity 2 block", mode="color")
# plt.figure()
# plot(p_fun, title="Pressure block", mode="color")
# plt.show()
return (r, u_fun_1, u_fun_2, p_fun)
(eig_b, u_fun_1_b, u_fun_2_b, p_fun_b) = run_block()
# -------------------------------------------------- #
# ERROR COMPUTATION #
def run_error(eig_m, eig_b, u_fun_1_m, u_fun_1_b, u_fun_2_m, u_fun_2_b, p_fun_m, p_fun_b):
err_inf_sup = abs(sqrt(eig_b) - sqrt(eig_m))/sqrt(eig_m)
print("Relative error for inf-sup constant equal to", err_inf_sup)
assert isclose(err_inf_sup, 0., atol=1.e-8)
# Even after normalization, eigenfunctions may have different signs. Try both and assume that the correct
# error computation is the one for which the error is minimum
err_1_plus = u_fun_1_b + u_fun_1_m
err_2_plus = u_fun_2_b + u_fun_2_m
err_p_plus = p_fun_b + p_fun_m
err_1_minus = u_fun_1_b - u_fun_1_m
err_2_minus = u_fun_2_b - u_fun_2_m
err_p_minus = p_fun_b - p_fun_m
err_1_plus_norm = assemble(inner(grad(err_1_plus), grad(err_1_plus))*dx)
err_2_plus_norm = assemble(inner(grad(err_2_plus), grad(err_2_plus))*dx)
err_p_plus_norm = assemble(err_p_plus*err_p_plus*dx)
err_1_minus_norm = assemble(inner(grad(err_1_minus), grad(err_1_minus))*dx)
err_2_minus_norm = assemble(inner(grad(err_2_minus), grad(err_2_minus))*dx)
err_p_minus_norm = assemble(err_p_minus*err_p_minus*dx)
u_fun_1_norm = assemble(inner(grad(u_fun_1_m), grad(u_fun_1_m))*dx)
u_fun_2_norm = assemble(inner(grad(u_fun_2_m), grad(u_fun_2_m))*dx)
p_fun_norm = assemble(p_fun_m*p_fun_m*dx)
def select_error(err_plus, err_plus_norm, err_minus, err_minus_norm, vec_norm, component_name):
ratio_plus = sqrt(err_plus_norm/vec_norm)
ratio_minus = sqrt(err_minus_norm/vec_norm)
if ratio_minus < ratio_plus:
print("Relative error for ", component_name, "component of eigenvector equal to", ratio_minus, "(the one with opposite sign was", ratio_plus, ")")
assert isclose(ratio_minus, 0., atol=1.e-6)
# plt.figure()
# plot(err_minus, title="Error " + component_name, mode="color")
# plt.show()
else:
print("Relative error for", component_name, "component of eigenvector equal to", ratio_plus, "(the one with opposite sign was", ratio_minus, ")")
assert isclose(ratio_plus, 0., atol=1.e-6)
# plt.figure()
# plot(err_plus, title="Error " + component_name, mode="color")
# plt.show()
select_error(err_1_plus, err_1_plus_norm, err_1_minus, err_1_minus_norm, u_fun_1_norm, "velocity 1")
select_error(err_2_plus, err_2_plus_norm, err_2_minus, err_2_minus_norm, u_fun_2_norm, "velocity 2")
select_error(err_p_plus, err_p_plus_norm, err_p_minus, err_p_minus_norm, p_fun_norm, "pressure")
run_error(eig_m, eig_b, u_fun_1_m, u_fun_1_b, u_fun_2_m, u_fun_2_b, p_fun_m, p_fun_b)