Keywords: periodic boundary conditions, multiple fields, mean value constraint
Figure 1: Left: Computational domain Ω with boundaries Γ1,
Γ3 (periodic boundary conditions) and Γ2, Γ4 (homogeneous
Dirichlet boundary conditions). Right: Magnitude of the resulting velocity field.
This example is a translation of the step-45 example from deal.ii which solves Stokes flow on a quarter circle. In particular it shows how to use periodic boundary conditions, how to solve a problem with multiple unknown fields, and how to enforce a specific mean value of the solution. For the mesh generation we use Gmsh.jl and then use FerriteGmsh.jl to import the mesh into Ferrite's format.
The strong form of Stokes flow with velocity u and pressure p can be written as follows: −Δu+∇p=(exp(−100||x−(0.75,0.1)||2),0)=:b∀x∈Ω,−∇⋅u=0∀x∈Ω,
The corresponding weak form reads as follows: Find (u,p)∈U×L2 s.t. ∫Ω[[δu⊗∇]:[u⊗∇]−(∇⋅δu) p ]dΩ=∫Ωδu⋅b dΩ∀δu∈U,∫Ω−(∇⋅u) δp dΩ=0∀δp∈L2,
With this formulation and boundary conditions for u the pressure will
only be determined up to a constant. We will therefore add an additional constraint which
fixes this constant (see deal.ii
step-11 for some more
discussion around this). In particular, we will enforce the mean value of the pressure on
the boundary to be 0, i.e. ∫Γp dΓ=0. One option is to
enforce this using a Lagrange multiplier. This would give a contribution λ∫Γδp dΓ to the second equation in the weak form above,
and a third equation δλ∫Γp dΓ=0 so that we
can solve for λ. However, since we in this case are not interested in computing
λ, and since the constraint is linear, we can directly embed this constraint
using an AffineConstraint
in Ferrite.
After FE discretization we obtain a linear system of the form K__ a_=f_, where K__=[K__uuK__TpuK__pu0__],a_=[a_ua_p],f_=[f_u0_],
The affine constraint to enforce zero mean pressure on the boundary is obtained from C__p a_p=0_, where (C__p)1j=∫Γϕpj dΓ.
Note
The constraint matrix C__p is the same matrix we would have obtained when assembling the system with the Lagrange multiplier. In that case the full system would be K__=[K__uuK__Tpu0__K__pu0__C__Tp0__C__p0],a_=[a_ua_pa_λ],f_=[f_u0_0_].
What follows is a program spliced with comments.
using Ferrite, FerriteGmsh, Gmsh, Tensors, LinearAlgebra, SparseArrays
Gmsh.jl
¶In the setup_grid
function below we use the
Gmsh.jl
package for setting up the geometry and
performing the meshing. We will not discuss this part in much detail but refer to the
Gmsh API documentation instead. The
most important thing to note is the mesh periodicity constraint that is applied between
the "inlet" and "outlet" parts using gmsh.model.set_periodic
. This is necessary to later
on apply a periodicity constraint for the approximated velocity field.
function setup_grid(h = 0.05)
# Initialize gmsh
Gmsh.initialize()
gmsh.option.set_number("General.Verbosity", 2)
# Add the points
o = gmsh.model.geo.add_point(0.0, 0.0, 0.0, h)
p1 = gmsh.model.geo.add_point(0.5, 0.0, 0.0, h)
p2 = gmsh.model.geo.add_point(1.0, 0.0, 0.0, h)
p3 = gmsh.model.geo.add_point(0.0, 1.0, 0.0, h)
p4 = gmsh.model.geo.add_point(0.0, 0.5, 0.0, h)
# Add the lines
l1 = gmsh.model.geo.add_line(p1, p2)
l2 = gmsh.model.geo.add_circle_arc(p2, o, p3)
l3 = gmsh.model.geo.add_line(p3, p4)
l4 = gmsh.model.geo.add_circle_arc(p4, o, p1)
# Create the closed curve loop and the surface
loop = gmsh.model.geo.add_curve_loop([l1, l2, l3, l4])
surf = gmsh.model.geo.add_plane_surface([loop])
# Synchronize the model
gmsh.model.geo.synchronize()
# Create the physical domains
gmsh.model.add_physical_group(1, [l1], -1, "Γ1")
gmsh.model.add_physical_group(1, [l2], -1, "Γ2")
gmsh.model.add_physical_group(1, [l3], -1, "Γ3")
gmsh.model.add_physical_group(1, [l4], -1, "Γ4")
gmsh.model.add_physical_group(2, [surf])
# Add the periodicity constraint using 4x4 affine transformation matrix,
# see https://en.wikipedia.org/wiki/Transformation_matrix#Affine_transformations
transformation_matrix = zeros(4, 4)
transformation_matrix[1, 2] = 1 # -sin(-pi/2)
transformation_matrix[2, 1] = -1 # cos(-pi/2)
transformation_matrix[3, 3] = 1
transformation_matrix[4, 4] = 1
transformation_matrix = vec(transformation_matrix')
gmsh.model.mesh.set_periodic(1, [l1], [l3], transformation_matrix)
# Generate a 2D mesh
gmsh.model.mesh.generate(2)
# Save the mesh, and read back in as a Ferrite Grid
grid = mktempdir() do dir
path = joinpath(dir, "mesh.msh")
gmsh.write(path)
togrid(path)
end
# Finalize the Gmsh library
Gmsh.finalize()
return grid
end
setup_grid (generic function with 2 methods)
As mentioned in the introduction we will use a quadratic approximation for the velocity
field and a linear approximation for the pressure to ensure that we fulfill the LBB
condition. We create the corresponding FE values with interpolations ipu
for the
velocity and ipp
for the pressure. Note that we specify linear geometric mapping
(ipg
) for both the velocity and pressure because our grid contains linear
triangles. However, since linear mapping is default this could have been skipped.
We also construct facet-values for the pressure since we need to integrate along
the boundary when assembling the constraint matrix C__.
function setup_fevalues(ipu, ipp, ipg)
qr = QuadratureRule{RefTriangle}(2)
cvu = CellValues(qr, ipu, ipg)
cvp = CellValues(qr, ipp, ipg)
qr_facet = FacetQuadratureRule{RefTriangle}(2)
fvp = FacetValues(qr_facet, ipp, ipg)
return cvu, cvp, fvp
end
setup_fevalues (generic function with 1 method)
The setup_dofs
function creates the DofHandler
, and adds the two fields: a
vector field :u
with interpolation ipu
, and a scalar field :p
with interpolation
ipp
.
function setup_dofs(grid, ipu, ipp)
dh = DofHandler(grid)
add!(dh, :u, ipu)
add!(dh, :p, ipp)
close!(dh)
return dh
end
setup_dofs (generic function with 1 method)
Now it is time to setup the ConstraintHandler
and add our boundary conditions and the
mean value constraint. This is perhaps the most interesting section in this example, and
deserves some attention.
Let's first discuss the assembly of the constraint matrix C__
and how to create an AffineConstraint
from it. This is done in the
setup_mean_constraint
function below. Assembling this is not so different from standard
assembly in Ferrite: we loop over all the facets, loop over the quadrature points, and loop
over the shape functions. Note that since there is only one constraint the matrix will
only have one row.
After assembling C
we construct an AffineConstraint
from it. We select the constrained
dof to be the one with the highest weight (just to avoid selecting one with 0 or a very
small weight), then move the remaining to the right hand side. As an example, consider the
case where the constraint equation C__p a_p is
w10p10+w23p23+w154p154=0
AffineConstraint
constructor expects.
Note
If all nodes along the boundary are equidistant all the weights would be the same. In this case we can construct the constraint without having to do any integration by simply finding all degrees of freedom that are located along the boundary (and using 1 as the weight). This is what is done in the deal.ii step-11 example.
function setup_mean_constraint(dh, fvp)
assembler = Ferrite.COOAssembler()
# All external boundaries
set = union(
getfacetset(dh.grid, "Γ1"),
getfacetset(dh.grid, "Γ2"),
getfacetset(dh.grid, "Γ3"),
getfacetset(dh.grid, "Γ4"),
)
# Allocate buffers
range_p = dof_range(dh, :p)
element_dofs = zeros(Int, ndofs_per_cell(dh))
element_dofs_p = view(element_dofs, range_p)
element_coords = zeros(Vec{2}, 3)
Ce = zeros(1, length(range_p)) # Local constraint matrix (only 1 row)
# Loop over all the boundaries
for (ci, fi) in set
Ce .= 0
getcoordinates!(element_coords, dh.grid, ci)
reinit!(fvp, element_coords, fi)
celldofs!(element_dofs, dh, ci)
for qp in 1:getnquadpoints(fvp)
dΓ = getdetJdV(fvp, qp)
for i in 1:getnbasefunctions(fvp)
Ce[1, i] += shape_value(fvp, qp, i) * dΓ
end
end
# Assemble to row 1
assemble!(assembler, [1], element_dofs_p, Ce)
end
C, _ = finish_assemble(assembler)
# Create an AffineConstraint from the C-matrix
_, J, V = findnz(C)
_, constrained_dof_idx = findmax(abs2, V)
constrained_dof = J[constrained_dof_idx]
V ./= V[constrained_dof_idx]
mean_value_constraint = AffineConstraint(
constrained_dof,
Pair{Int, Float64}[J[i] => -V[i] for i in 1:length(J) if J[i] != constrained_dof],
0.0,
)
return mean_value_constraint
end
setup_mean_constraint (generic function with 1 method)
We now setup all the boundary conditions in the setup_constraints
function below.
Since the periodicity constraint for this example is between two boundaries which are not
parallel to each other we need to i) compute the mapping between each mirror facet and the
corresponding image facet (on the element level) and ii) describe the dof relation between
dofs on these two facets. In Ferrite this is done by defining a transformation of entities
on the image boundary such that they line up with the matching entities on the mirror
boundary. In this example we consider the inlet Γ1 to be the image, and the
outlet Γ3 to be the mirror. The necessary transformation to apply then becomes a
rotation of π/2 radians around the out-of-plane axis. We set up the rotation matrix
R
, and then compute the mapping between mirror and image facets using
collect_periodic_facets
where the rotation is applied to the coordinates. In the
next step we construct the constraint using the PeriodicDirichlet
constructor.
We pass the constructor the computed mapping, and also the rotation matrix. This matrix is
used to rotate the dofs on the mirror surface such that we properly constrain
ux-dofs on the mirror to −uy-dofs on the image, and
uy-dofs on the mirror to ux-dofs on the image.
For the remaining part of the boundary we add a homogeneous Dirichlet boundary condition
on both components of the velocity field. This is done using the Dirichlet
constructor, which we have discussed in other tutorials.
function setup_constraints(dh, fvp)
ch = ConstraintHandler(dh)
# Periodic BC
R = rotation_tensor(π / 2)
periodic_faces = collect_periodic_facets(dh.grid, "Γ3", "Γ1", x -> R ⋅ x)
periodic = PeriodicDirichlet(:u, periodic_faces, R, [1, 2])
add!(ch, periodic)
# Dirichlet BC
Γ24 = union(getfacetset(dh.grid, "Γ2"), getfacetset(dh.grid, "Γ4"))
dbc = Dirichlet(:u, Γ24, (x, t) -> [0, 0], [1, 2])
add!(ch, dbc)
# Compute mean value constraint and add it
mean_value_constraint = setup_mean_constraint(dh, fvp)
add!(ch, mean_value_constraint)
# Finalize
close!(ch)
update!(ch, 0)
return ch
end
setup_constraints (generic function with 1 method)
Assembly of the global system is also something that we have seen in many previous
tutorials. One interesting thing to note here is that, since we have two unknown fields,
we use the dof_range
function to make sure we assemble the element contributions
to the correct block of the local stiffness matrix ke
.
function assemble_system!(K, f, dh, cvu, cvp)
assembler = start_assemble(K, f)
ke = zeros(ndofs_per_cell(dh), ndofs_per_cell(dh))
fe = zeros(ndofs_per_cell(dh))
range_u = dof_range(dh, :u)
ndofs_u = length(range_u)
range_p = dof_range(dh, :p)
ndofs_p = length(range_p)
ϕᵤ = Vector{Vec{2, Float64}}(undef, ndofs_u)
∇ϕᵤ = Vector{Tensor{2, 2, Float64, 4}}(undef, ndofs_u)
divϕᵤ = Vector{Float64}(undef, ndofs_u)
ϕₚ = Vector{Float64}(undef, ndofs_p)
for cell in CellIterator(dh)
reinit!(cvu, cell)
reinit!(cvp, cell)
ke .= 0
fe .= 0
for qp in 1:getnquadpoints(cvu)
dΩ = getdetJdV(cvu, qp)
for i in 1:ndofs_u
ϕᵤ[i] = shape_value(cvu, qp, i)
∇ϕᵤ[i] = shape_gradient(cvu, qp, i)
divϕᵤ[i] = shape_divergence(cvu, qp, i)
end
for i in 1:ndofs_p
ϕₚ[i] = shape_value(cvp, qp, i)
end
# u-u
for (i, I) in pairs(range_u), (j, J) in pairs(range_u)
ke[I, J] += (∇ϕᵤ[i] ⊡ ∇ϕᵤ[j]) * dΩ
end
# u-p
for (i, I) in pairs(range_u), (j, J) in pairs(range_p)
ke[I, J] += (-divϕᵤ[i] * ϕₚ[j]) * dΩ
end
# p-u
for (i, I) in pairs(range_p), (j, J) in pairs(range_u)
ke[I, J] += (-divϕᵤ[j] * ϕₚ[i]) * dΩ
end
# rhs
for (i, I) in pairs(range_u)
x = spatial_coordinate(cvu, qp, getcoordinates(cell))
b = exp(-100 * norm(x - Vec{2}((0.75, 0.1)))^2)
bv = Vec{2}((b, 0.0))
fe[I] += (ϕᵤ[i] ⋅ bv) * dΩ
end
end
assemble!(assembler, celldofs(cell), ke, fe)
end
return K, f
end
assemble_system! (generic function with 1 method)
We now have all the puzzle pieces, and just need to define the main function, which puts them all together.
function main()
# Grid
h = 0.05 # approximate element size
grid = setup_grid(h)
# Interpolations
ipu = Lagrange{RefTriangle, 2}()^2 # quadratic
ipp = Lagrange{RefTriangle, 1}() # linear
# Dofs
dh = setup_dofs(grid, ipu, ipp)
# FE values
ipg = Lagrange{RefTriangle, 1}() # linear geometric interpolation
cvu, cvp, fvp = setup_fevalues(ipu, ipp, ipg)
# Boundary conditions
ch = setup_constraints(dh, fvp)
# Global tangent matrix and rhs
coupling = [true true; true false] # no coupling between pressure test/trial functions
K = allocate_matrix(dh, ch; coupling = coupling)
f = zeros(ndofs(dh))
# Assemble system
assemble_system!(K, f, dh, cvu, cvp)
# Apply boundary conditions and solve
apply!(K, f, ch)
u = K \ f
apply!(u, ch)
# Export the solution
VTKGridFile("stokes-flow", grid) do vtk
write_solution(vtk, dh, u)
end
return
end
main (generic function with 1 method)
Run it!
main()
The resulting magnitude of the velocity field is visualized in Figure 1.
This notebook was generated using Literate.jl.