In this tutorial, we will look at how to
CellFields
at arbitrary pointstriangulations. We will consider examples for
Let T1 and T2 be two triangulations of a domain Ω. Let Vi be the finite element space defined on the triangulation Ti for i=1,2. Let fh∈V1. The interpolation problem is to find gh∈V2 such that
dofV2k(gh)=dofV2k(fh),∀k∈{1,…,NV2dof}For the purpose of this tutorial we require Test
, Gridap
along with the
following submodules of Gridap
using Test
using Gridap
using Gridap.CellData
using Gridap.Visualization
We now create a computational domain on the unit square [0,1]2 consisting of 5 cells per direction
domain = (0,1,0,1)
partition = (5,5)
𝒯₁ = CartesianDiscreteModel(domain, partition)
Gridap
offers the feature to evaluate functions at arbitrary
points in the domain. This will be shown in the next
section. Interpolation then takes advantage of this feature to
obtain the FEFunction
in the new space from the old one by
evaluating the appropriate degrees of freedom. Interpolation works
using the composite type Interpolable
to tell Gridap
that the
argument can be interpolated between triangulations.
Let us define the infinite dimensional function
f(x) = x[1] + x[2]
This function will be interpolated to the source FESpace
V1. The space can be built using
reffe₁ = ReferenceFE(lagrangian, Float64, 1)
V₁ = FESpace(𝒯₁, reffe₁)
Finally to build the function fh, we do
fₕ = interpolate_everywhere(f,V₁)
To construct arbitrary points in the domain, we use Random
package:
using Random
pt = Point(rand(2))
pts = [Point(rand(2)) for i in 1:3]
The finite element function fh can be evaluated at arbitrary points (or array of points) by
fₕ(pt), fₕ.(pts)
We can also check our results using
@test fₕ(pt) ≈ f(pt)
@test fₕ.(pts) ≈ f.(pts)
Now let us define the new triangulation T2 of
Ω. We build the new triangulation using a partition of 20 cells per
direction. The map can be passed as an argument to
CartesianDiscreteModel
to define the position of the vertices in
the new mesh.
partition = (20,20)
𝒯₂ = CartesianDiscreteModel(domain,partition)
As before, we define the new FESpace
consisting of second order
elements
reffe₂ = ReferenceFE(lagrangian, Float64, 2)
V₂ = FESpace(𝒯₂, reffe₂)
Now we interpolate fh onto V2 to obtain the new function
gh. The first step is to create the Interpolable
version of
fh.
ifₕ = Interpolable(fₕ)
Then to obtain gh, we dispatch ifₕ
and the new FESpace
V2
to the interpolate_everywhere
method of Gridap
.
gₕ = interpolate_everywhere(ifₕ, V₂)
We can also use
interpolate
if interpolating only on the free dofs or
interpolate_dirichlet
if interpolating the Dirichlet dofs of the
FESpace
.
ḡₕ = interpolate(ifₕ, V₂)
The finite element function ˉgh is the same as gh in this example since all the dofs are free.
@test gₕ.cell_dof_values == ḡₕ.cell_dof_values
Now we obtain a finite element function using interpolate_dirichlet
g̃ₕ = interpolate_dirichlet(ifₕ, V₂)
Now ˜gh will be equal to 0 since there are
no Dirichlet nodes defined in the FESpace
. We can check by running
g̃ₕ.cell_dof_values
Like earlier we can check our results for gₕ
:
@test fₕ(pt) ≈ gₕ(pt) ≈ f(pt)
@test fₕ.(pts) ≈ gₕ.(pts) ≈ f.(pts)
We can visualize the results using Paraview
writevtk(get_triangulation(fₕ), "source", cellfields=["fₕ"=>fₕ])
writevtk(get_triangulation(gₕ), "target", cellfields=["gₕ"=>gₕ])
which produces the following output
The procedure is identical to Lagrangian finite element spaces, as discussed in the previous section. The extra thing here is that functions in Raviart-Thomas spaces are vector-valued. The degrees of freedom of the RT spaces are fluxes of the function across the edge of the element. Refer to the tutorial on Darcy equation with RT for more information on the RT elements.
Assuming a function
f(x) = VectorValue([x[1], x[2]])
on the domain, we build the associated finite dimensional version fh∈V1.
reffe₁ = ReferenceFE(raviart_thomas, Float64, 1) # RT space of order 1
V₁ = FESpace(𝒯₁, reffe₁)
fₕ = interpolate_everywhere(f, V₁)
As before, we can evaluate the RT function on any arbitrary point in the domain.
fₕ(pt), fₕ.(pts)
Constructing the target RT space and building the Interpolable
object,
reffe₂ = ReferenceFE(raviart_thomas, Float64, 1) # RT space of order 1
V₂ = FESpace(𝒯₂, reffe₂)
ifₕ = Interpolable(fₕ)
we can construct the new FEFunction
gh∈V2 from fh
gₕ = interpolate_everywhere(ifₕ, V₂)
Like earlier we can check our results
@test gₕ(pt) ≈ f(pt) ≈ fₕ(pt)
We can also interpolate vector-valued functions across triangulations. First, we define a vector-valued function on a two-dimensional mesh.
f(x) = VectorValue([x[1], x[1]+x[2]])
We then create a vector-valued reference element containing linear elements along with the source finite element space V1.
reffe₁ = ReferenceFE(lagrangian, VectorValue{2,Float64}, 1)
V₁ = FESpace(𝒯₁, reffe₁)
fₕ = interpolate_everywhere(f, V₁)
The target finite element space V2 can be defined in a similar manner.
reffe₂ = ReferenceFE(lagrangian, VectorValue{2,Float64}, 2)
V₂ = FESpace(𝒯₂, reffe₂)
The rest of the process is similar to the previous sections, i.e.,
define the Interpolable
version of fh and use
interpolate_everywhere
to find gh∈V₂.
ifₕ = Interpolable(fₕ)
gₕ = interpolate_everywhere(ifₕ, V₂)
We can then check the results
@test gₕ(pt) ≈ f(pt) ≈ fₕ(pt)
Similarly, it is possible to interpolate between multi-field finite element functions. First, we define the components h1(x),h2(x) of a multi-field function h(x) as follows.
h₁(x) = x[1]+x[2]
h₂(x) = x[1]
Next we create a Lagrangian finite element space containing linear elements.
reffe₁ = ReferenceFE(lagrangian, Float64, 1)
V₁ = FESpace(𝒯₁, reffe₁)
Next we create a MultiFieldFESpace
V1×V1 and
interpolate the function h(x) to the source space V1.
V₁xV₁ = MultiFieldFESpace([V₁,V₁])
fₕ = interpolate_everywhere([h₁, h₂], V₁xV₁)
Similarly, the target multi-field finite element space is created using Ω2.
reffe₂ = ReferenceFE(lagrangian, Float64, 2)
V₂ = FESpace(𝒯₂, reffe₂)
V₂xV₂ = MultiFieldFESpace([V₂,V₂])
Now, to find gh∈V2×V2, we first extract the components of
fh and obtain the Interpolable
version of the components.
fₕ¹, fₕ² = fₕ
ifₕ¹ = Interpolable(fₕ¹)
ifₕ² = Interpolable(fₕ²)
We can then use interpolate_everywhere
on the Interpolable
version of the components and obtain gh∈V2×V2 as
follows.
gₕ = interpolate_everywhere([ifₕ¹,ifₕ²], V₂xV₂)
We can then check the results of the interpolation, component-wise.
gₕ¹, gₕ² = gₕ
@test fₕ¹(pt) ≈ gₕ¹(pt)
@test fₕ²(pt) ≈ gₕ²(pt)
Gridap contributors acknowledge support received from Google, Inc. through the Google Summer of Code 2021 project A fast finite element interpolator in Gridap.jl.
This notebook was generated using Literate.jl.