API
Documentation for Bramble.jl
's public API.
Geometries and meshes
Bramble.interval
— Methodinterval(x, y)
interval(x::CartesianProduct{1})
Returns a 1
-dimensional CartesianProduct representing the interval [x
,y
].
Example
julia> interval(0, 1)
CartesianProduct{1,Float64}((0.0,1.0))
Bramble.:×
— Function×(X::CartesianProduct, Y::CartesianProduct)
Returns the cartesian product of two CartesianProducts X
and Y
as a new CartesianProduct.
Example
julia> X = interval(0, 1); Y = interval(2, 3); X × Y;
Type: Float64
Dim: 2
Set: [0.0, 1.0] × [2.0, 3.0]
Bramble.dim
— Functiondim(X::CartesianProduct)
dim(::Type{<:CartesianProduct})
Returns the topological dimension of a CartesianProduct.
Example
julia> X = cartesianproduct(0, 1); dim(X)
1
julia> Y = cartesianproduct(((0,1), (4,5))); dim(Y)
2
dim(Ω::Domain)
dim(::Type{<:Domain})
Returns the topological dimension of the Domain Ω
.
Example
julia> I = interval(0.0, 1.0); dim(domain(I × I))
2
dim(Ωₕ::MeshType)
dim(::Type{<:MeshType})
Returns the tolopogical dimension of Ωₕ
.
Bramble.domain
— Methoddomain(X::CartesianProduct)
domain(Ω::CartesianProduct, markers::MarkersType)
Returns a Domain from a CartesianProduct, assuming the single Marker ":dirichlet" => x -> 0
. Alternatively, a set of Marker can be passed as an argument.
Example
julia> domain(interval(0,1))
Type: Float64
Dim: 1
Set: [0.0, 1.0]
Boundary markers: :dirichlet
julia> I = interval(0,1); m = markers( :dirichlet => @embed(I, x->x[1]-1)), :neumann => @embed(I, x->x[1]-0)); domain(interval(0,1), m)
Type: Float64
Dim: 1
Set: [0.0, 1.0]
Boundary markers: :dirichlet, :neumann
Bramble.create_markers
— Functioncreate_markers(m::MarkerType...)
Converts several Pair{Symbol,F}
(:symbol => func) to Markers. These are to be passed in the construction of a Domain. The functions need to be defined as BrambleBareFunctions.
Example
julia> create_markers( :dirichlet => @embed(X, x -> x[1]-1), :neumann => @embed(X, x -> x[2]-0) )
Bramble.markers
— MethodBramble.labels
— MethodBramble.@embed
— Macro@embed(X:CartesianProduct, f)
@embed(X:Domain, f)
@embed(X:MeshType, f)
@embed(X:SpaceType, f)
Returns a new wrapped version of function f
. If
X
is not a SpaceType, the topological dimension ofX
is used to caracterize the types in the returning BrambleBareFunction; in this case, the new function can be applied to tuples (with point coordinates) or cartesian indices (ifX
is a MeshType)X
is a SpaceType, the returning function type, BrambleGridSpaceFunction, is caracterized by the dimension of the elements of the space.
Example
julia> Ω = domain(interval(0,1) × interval(0,1)); f = @embed Ω x -> x[1]*x[2]+1; # or f = @embed(Ω, x -> x[1]*x[2]+1)
Bramble.mesh
— Methodmesh(Ω::Domain, npts::Int, unif::Bool)
mesh(Ω::Domain, npts::NTuple{D}, unif::NTuple{D})
Returns a Mesh1D or a MeshnD ($D=2,3$) defined on the Domain Ω
. The number of points for each coordinate projection mesh are given in the tuple npts
. The distribution of points on the submeshes are encoded in the tuple unif
.
For future reference, the mesh points are denoted as
- 1D mesh, with `npts` = ``N_x``
\[x_i, \, i=1,\dots,N.\]
- 2D mesh, with
npts
= ($N_x$, $N_y$)
\[(x_i,y_j), \, i=1,\dots,N_x, \, j=1,\dots,N_y\]
- 3D mesh, with
npts
= ($N_x$, $N_y$, $N_z$)
\[(x_i,y_j,z_l), \, i=1,\dots,N_x, \, j=1,\dots,N_y, \, l=1,\dots,N_z.\]
Examples
julia> I = interval(0,1); Ωₕ = mesh(domain(I), 10, true)
1D mesh
nPoints: 10
Markers: Dirichlet
julia> X = domain(interval(0,1) × interval(4,5)); Ωₕ = mesh(X, (10, 15), (true, false))
2D mesh
nPoints: 150
Markers: ["Dirichlet"]
Submeshes:
x direction | nPoints: 10
y direction | nPoints: 15
Bramble.points
— Functionpoints(Ωₕ::Mesh1D)
points(Ωₕ::Mesh1D, i)
points(Ωₕ::Mesh1D, Iterator)
Returns a vector with all the points $x_i, \, i=1,\dots,N$ in Ωₕ
. A second argument can be passed. If it is an Int
or a CartesianIndex{1}
, it returns the i
-th point of Ωₕ
, $x_i$. If the second argument is Iterator
then the function returns a generator iterating over the points.
points(Ωₕ::MeshnD)
points(Ωₕ::MeshnD{D}, idx)
points(Ωₕ::MeshnD{D}, Iterator)
Returns a tuple with the points of Ωₕ
. If the Tuple
idx
is passed as the second argument is passed, it returns the tuple with the point corresponding to that index. Alternatively, if Iterator
is passed as the second argument, a generator iterating over all points of the mesh is returned.
- 2D mesh, with
npts
= ($N_x$, $N_y$)
\[([x_i]_{i=1}^{N_x}, [y_j]_{j=1}^{N_y})\]
- 3D mesh, with
npts
= ($N_x$, $N_y$, $N_z$)
\[([x_i]_{i=1}^{N_x}, [y_j]_{j=1}^{N_y}, [z_l]_{l=1}^{N_z}).\]
Bramble.hₘₐₓ
— Functionhₘₐₓ(Ωₕ::Mesh1D)
Returns the maximum over the space stepsize $h_i$of mesh Ωₕ
\[h_{max} \vcentcolon = \max_{i=1,\dots,N} x_i - x_{i-1}.\]
hₘₐₓ(Ωₕ::MeshnD)
Returns the maximum diagonal of mesh Ωₕ
.
- 2D mesh
\[\max_{i,j} \Vert (h_{x,i}, h_{y,j}) \Vert_2\]
- 3D mesh
\[\max_{i,j,l} \Vert (h_{x,i}, h_{y,j}, h_{z,l}) \Vert_2\]
Bramble.npoints
— Functionnpoints(Ωₕ::Mesh1D)
npoints(Ωₕ::Mesh1D, Tuple)
Returns the number of points $x_i$ in Ωₕ
. If the second argument is passed, it returns the same information as a 1
-tuple.
Example
julia> Ωₕ = mesh(domain(interval(0,1)), 10, true); npoints(Ωₕ)
10
npoints(Ωₕ::MeshnD)
npoints(Ωₕ::MeshnD, Tuple)
Returns the number of points of mesh Ωₕ
. If Tuple
is passed as the second argument, it returns a tuple with the number of points of each submesh composing Ωₕ
.
Spaces
Bramble.gridspace
— Functiongridspace(Ωₕ::MeshType)
Constructor for a GridSpace defined on the mesh Ωₕ
. This builds the weights for the inner products mentioned in GridSpace as well as the differentiation matrices associated with the grid points of mesh Ωₕ
.
Bramble.element
— Functionelement(Wₕ::SpaceType)
Returns a VectorElement for GridSpace Wₕ
with uninitialized components.
element(Wₕ::SpaceType, α::Number)
Returns a VectorElement for GridSpace Wₕ
with all components equal to α
.
element(Wₕ::SpaceType, v::AbstractVector)
Returns a VectorElement for GridSpace Wₕ
with the same coefficients of v
.
Bramble.mesh
— Methodmesh(Wₕ::SpaceType)
mesh(::Type{<:SpaceType{MType}})
Returns the mesh
on which the GridSpace Wₕ
is defined. If the input argument is a type derived from SpaceType then the function returns the MeshType associated with it.
Bramble.avgₕ
— Functionavgₕ(Wₕ::SpaceType, f)
Returns a VectorElement with the average of function f
with respect to the cell_measure of mesh(Wₕ)
around each grid point. It can accept any function (like x->x[2]+x[1])
) or a BrambleBareFunction. The latter is preferred. It is defined as follows
- 1D case
\[\textrm{avg}ₕ(x_i) = \frac{1}{|\square_i|} \int_{\square_i} f(x) dx, \, i = 1,\dots,N\]
- 2D case
\[\textrm{avg}ₕ(x_i, y_j) = \frac{1}{|\square_{i,j}|} \iint_{\square_{i,j}} f(x,y) dA, \, i = 1,\dots,N_x, j = 1,\dots,N_y\]
- 3D case
\[\textrm{avg}ₕ(x_i, y_j, z_l) = \frac{1}{|\square_{i,j,l}|} \iiint_{\square_{i,j,l}} f(x,y,z) dV, \, i = 1,\dots,N_x, j = 1,\dots,N_y, l = 1,\dots,N_z\]
Please check the implementations of functions cell_measure (for the 1
-dimensional case) and cell_measure (for the n
-dimensional cases).
Bramble.avgₕ!
— Functionavgₕ!(uₕ::VectorElement, f)
In-place version of averaging operator avgₕ.
Bramble.Rₕ
— FunctionRₕ(Wₕ::SpaceType, f)
Standard nodal restriction operator. It returns a VectorElement with the result of evaluating the function f
at the points of mesh(Wₕ)
. It can accept any function (like x->x[2]+x[1])
) or a BrambleBareFunction. The latter is preferred.
- 1D case
\[\textrm{R}ₕ(x_i) = f(x_i), \, i = 1,\dots,N\]
- 2D case
\[\textrm{R}ₕ (x_i, y_j)= f(x_i, y_j), \, i = 1,\dots,N_x, j = 1,\dots,N_y\]
- 3D case
\[\textrm{R}ₕ (x_i, y_j, z_l)= f(x_i, y_j, z_l), \, i = 1,\dots,N_x, j = 1,\dots,N_y, l = 1,\dots,N_z\]
Bramble.Rₕ!
— FunctionRₕ!(uₕ::VectorElement, f)
In-place version of the restriction operator Rₕ.
Bramble.diff₋ₓ
— Functiondiff₋ₓ(Wₕ::SpaceType)
diff₋ₓ(Ωₕ::MeshType)
Returns a MatrixElement implementing the backward difference matrix for the mesh grid of Wₕ
, in the x
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by diff₋ₓ. It can also be calculated passing a mesh as argument.
diff₋ₓ(uₕ::VectorElement)
Returns the backward difference, in the x
direction, of the element uₕ
.
- 1D case
\[\textrm{diff}_{-x} \textrm{u}_h(x_i) \vcentcolon = \textrm{u}_h(x_i) - \textrm{u}_h(x_{i-1})\]
- 2D and 3D case
\[\textrm{diff}_{-x} \textrm{u}_h(x_i, \dots) \vcentcolon = \textrm{u}_h(x_i, \dots)-\textrm{u}_h(x_{i-1}, \dots)\]
diff₋ₓ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the backward difference matrix, in the x
direction, by Uₕ
.
Bramble.diff₋ᵧ
— Functiondiff₋ᵧ(Wₕ::SpaceType)
diff₋ᵧ(Ωₕ::MeshType)
Returns a MatrixElement implementing the backward difference matrix for the mesh grid of Wₕ
, in the y
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by diff₋ᵧ. It can also be calculated passing a mesh as argument.
diff₋ᵧ(uₕ::VectorElement)
Returns the backward difference, in the y
direction, of the element uₕ
.
- 2D and 3D case
\[\textrm{diff}_{-y} \textrm{u}_h(x_i, y_j,\dots) \vcentcolon = \textrm{u}_h(x_i, y_j,\dots)-\textrm{u}_h(x_i, y_{j-1}, \dots)\]
diff₋ᵧ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the backward difference matrix, in the y
direction, by Uₕ
.
Bramble.diff₋₂
— Functiondiff₋₂(Wₕ::SpaceType)
diff₋₂(Ωₕ::MeshType)
Returns a MatrixElement implementing the backward difference matrix for the mesh grid of Wₕ
, in the z
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by diff₋₂. It can also be calculated passing a mesh as argument.
diff₋₂(uₕ::VectorElement)
Returns the backward difference, in the z
direction, of the element uₕ
.
\[\textrm{diff}_{-z} \textrm{u}_h(x_i, y_j,z_l) \vcentcolon = \textrm{u}_h(x_i, y_j,z_l)-\textrm{u}_h(x_i, y_j, z_{l-1})\]
diff₋₂(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the backward difference matrix, in the z
direction, by Uₕ
.
Bramble.diff₋ₕ
— Functiondiff₋ₕ(Wₕ::SpaceType)
Returns a tuple of MatrixElements implementing the backward difference operators in the x
, y
, and z
directions. If the problem is 1
-dimensikonal, it returns a single MatrixElement.
diff₋ₕ(uₕ::VectorElement)
Returns a tuple of VectorElements implementing the backward difference operators in the x
, y
, and z
directions applied to uₕ
. If the problem is 1
-dimensional, it returns a single VectorElement.
diff₋ₕ(Uₕ::MatrixElement)
Returns a tuple of MatrixElements implementing the forward difference operators in the x
, y
, and z
directions applied to Uₕ
. If the problem is 1
-dimensional, it returns a single MatrixElement.
Bramble.diffₓ
— Functiondiffₓ(Wₕ::SpaceType)
diffₓ(Ωₕ::MeshType)
Returns a MatrixElement implementing the forward difference matrix for the mesh grid of Wₕ
, in the x
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by diffₓ. It can also accept a mesh as an argument.
diffₓ(uₕ::VectorElement)
Returns the forward difference, in the x
direction, of the element uₕ
.
- 1D case
\[\textrm{diff}_{x} \textrm{u}_h(x_i) \vcentcolon = \textrm{u}_h(x_{i+1}) - \textrm{u}_h(x_i)\]
- 2D and 3D case
\[\textrm{diff}_{x} \textrm{u}_h(x_i, \dots) \vcentcolon = \textrm{u}_h(x_{i+1}, \dots)-\textrm{u}_h(x_i, \dots)\]
diffₓ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the forward difference matrix, in the x
direction, by Uₕ
.
Bramble.diffᵧ
— Functiondiffᵧ(Wₕ::SpaceType)
diffᵧ(Ωₕ::MeshType)
Returns a MatrixElement implementing the forward difference matrix for the mesh grid of Wₕ
, in the y
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by diffᵧ. It can also accept a mesh as an argument.
diffᵧ(uₕ::VectorElement)
Returns the forward difference, in the y
direction, of the element uₕ
.
- 2D and 3D case
\[\textrm{diff}_{y} \textrm{u}_h(x_i, y_j,\dots) \vcentcolon = \textrm{u}_h(x_i, y_{j+1},\dots)-\textrm{u}_h(x_i, y_j, \dots)\]
diffᵧ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the forward difference matrix, in the y
direction, by Uₕ
.
Bramble.diff₂
— Functiondiff₂(Wₕ::SpaceType)
diff₂(Ωₕ::MeshType)
Returns a MatrixElement implementing the forward difference matrix for the mesh grid of Wₕ
, in the z
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by diff₂. It can also accept a mesh as an argument.
diff₂(uₕ::VectorElement)
Returns the forward difference, in the z
direction, of the element uₕ
.
\[\textrm{diff}_{-z} \textrm{u}_h(x_i, y_j,z_l) \vcentcolon = \textrm{u}_h(x_i, y_j,z_l)-\textrm{u}_h(x_i, y_j, z_{l-1})\]
diff₂(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the forward difference matrix, in the z
direction, by Uₕ
.
Bramble.diffₕ
— Methoddiffₕ(Wₕ::SpaceType)
Returns a tuple of MatrixElements implementing the forward difference operators in the x
, y
and z
directions. If the problem is 1
-dimensional, it returns a single MatrixElement.
Bramble.D₋ₓ
— FunctionD₋ₓ(Wₕ::SpaceType)
D₋ₓ(Ωₕ::MeshType)
Returns a MatrixElement implementing the backward finite difference matrix for the mesh grid of Wₕ
, in the x
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by D₋ₓ. It can also accept a mesh as an argument.
D₋ₓ(uₕ::VectorElement)
Returns the backward finite difference, in the x
direction, of the element uₕ
.
- 1D case
\[\textrm{D}_{-x} \textrm{u}_h (x_i) \vcentcolon = \frac{\textrm{u}_h(x_i) - \textrm{u}_h(x_{i-1})}{h_i}\]
- 2D and 3D case
\[\textrm{D}_{-x} \textrm{u}_h (x_i, \dots) \vcentcolon = \frac{\textrm{u}_h(x_i, \dots)-\textrm{u}_h(x_{i-1}, \dots)}{h_{x,i}}\]
D₋ₓ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the backward finite difference matrix D₋ₓ
with the MatrixElement Uₕ
.
Bramble.D₋ᵧ
— FunctionD₋ᵧ(Wₕ::SpaceType)
D₋ᵧ(Ωₕ::MeshType)
Returns a MatrixElement implementing the backward finite difference matrix for the mesh grid of Wₕ
, in the y
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by D₋ᵧ. It can also accept a mesh as an argument.
D₋ᵧ(uₕ::VectorElement)
Returns the backward finite difference, in the y
direction, of the element uₕ
.
\[\textrm{D}_{-y} \textrm{u}_h(x_i, y_j, \dots) \vcentcolon = \frac{\textrm{u}_h(x_i, y_j, \dots)-\textrm{u}_h(x_i, y_{j-1}, \dots)}{h_{y,j}}\]
D₋ᵧ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the backward finite difference matrix D₋ᵧ
with the MatrixElement Uₕ
.
Bramble.D₋₂
— FunctionD₋₂(Wₕ::SpaceType)
D₋₂(Ωₕ::MeshType)
Returns a MatrixElement implementing the backward finite difference matrix for the mesh grid of Wₕ
, in the z
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by D₋₂. It can also accept a mesh as an argument.
D₋₂(uₕ::VectorElement)
Returns the backward finite difference, in the z
direction, of the element uₕ
.
\[\textrm{D}_{-z} \textrm{u}_h(x_i, y_j, z_l) \vcentcolon = \frac{\textrm{u}_h(x_i, y_j, z_l)-\textrm{u}_h(x_i, y_j, z_)}{h_{z,l}}\]
D₋₂(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the backward finite difference matrix D₋₂
with the MatrixElement Uₕ
.
Bramble.∇₋ₕ
— Method∇₋ₕ(Wₕ::SpaceType)
Returns a tuple of MatrixElements implementing the backward finite difference operators in the x
, y
, and z
directions. If the problem is 1
-dimensional, it returns a single MatrixElement.
Bramble.jumpₓ
— Functionjumpₓ(Wₕ::SpaceType)
jumpₓ(Ωₕ::MeshType)
Returns a MatrixElement implementing the jump matrix for the mesh grid of Wₕ
, in the x
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by jumpₓ. It also accepts a mesh as an argument.
jumpₓ(uₕ::VectorElement)
Returns the jump, in the x
direction, of the element uₕ
.
- 1D case
\[\textrm{jump}_{x} \textrm{u}_h(x_i) \vcentcolon = \textrm{u}_h(x_i) - \textrm{u}_h(x_{i+1})\]
- 2D and 3D case
\[\textrm{jump}_{x} \textrm{u}_h(x_i, \dots) \vcentcolon = \textrm{u}_h(x_i, \dots)-\textrm{u}_h(x_{i+1}, \dots)\]
jumpₓ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the jump matrix jumpₓ with the MatrixElement Uₕ
.
Bramble.jumpᵧ
— Functionjumpᵧ(Wₕ::SpaceType)
jumpᵧ(Ωₕ::MeshType)
Returns a MatrixElement implementing the jump matrix for the mesh grid of Wₕ
, in the y
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by jumpᵧ. It also accepts a mesh as an argument.
jumpᵧ(uₕ::VectorElement)
Returns the jump, in the y
direction, of the element uₕ
.
\[\textrm{jump}_{y} \textrm{u}_h(x_i, y_j,\dots) \vcentcolon = \textrm{u}_h(x_i, y_j\dots)-\textrm{u}_h(x_i, y_{j+1}, \dots)\]
jumpᵧ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the jump matrix jumpᵧ with the MatrixElement Uₕ
.
Bramble.jumpₕ
— Methodjumpₕ(Wₕ::SpaceType)
Returns a tuple of MatrixElements implementing the jump operators in the x
, y
, and z
directions. If the problem is 1
-dimensional, it returns a single MatrixElement.
Bramble.jump₂
— Functionjump₂(Wₕ::SpaceType)
jump₂(Ωₕ::MeshType)
Returns a MatrixElement implementing the jump matrix for the mesh grid of Wₕ
, in the z
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by jump₂. It also accepts a mesh as an argument.
jump₂(uₕ::VectorElement)
Returns the jump, in the z
direction, of the element uₕ
.
\[\textrm{jump}_{z} \textrm{u}_h(x_i, y_j,z_l) \vcentcolon = \textrm{u}_h(x_i, y_j, z_l)-\textrm{u}_h(x_i, y_j, z_{l+1})\]
jump₂(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the jump matrix jump₂ with the MatrixElement Uₕ
.
Bramble.M₋ₕₓ
— FunctionM₋ₕₓ(Wₕ::SpaceType)
M₋ₕₓ(Ωₕ::MeshType)
Returns a MatrixElement implementing the average matrix for the mesh grid of Wₕ
, in the x
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by M₋ₕₓ. It also accepts a mesh as argument.
M₋ₕₓ(uₕ::VectorElement)
Returns the average, in the x
direction, of the element uₕ
.
- 1D case
\[\textrm{M}_{-hx} \textrm{u}_h(x_i) \vcentcolon = \frac{\textrm{u}_h(x_i) + \textrm{u}_h(x_{i-1})}{2}\]
- 2D and 3D case
\[\textrm{M}_{-hx} \textrm{u}_h(x_i, \dots) \vcentcolon = \frac{\textrm{u}_h(x_i, \dots)+\textrm{u}_h(x_{i-1}, \dots)}{2}\]
M₋ₕₓ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the average matrix M₋ₕₓ with the MatrixElement Uₕ
.
Bramble.M₋ₕᵧ
— FunctionM₋ₕᵧ(Wₕ::SpaceType)
M₋ₕᵧ(Ωₕ::MeshType)
Returns a MatrixElement implementing the average matrix for the mesh grid of Wₕ
, in the y
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by M₋ₕᵧ. It also accepts a mesh as argument.
M₋ₕᵧ(uₕ::VectorElement)
Returns the average, in the y
direction, of the element uₕ
.
\[\textrm{M}_{-hy} \textrm{u}_h(x_i, y_j,\dots) \vcentcolon = \textrm{u}_h(x_i, y_j\dots)-\textrm{u}_h(x_i, y_{j+1}, \dots)\]
M₋ₕᵧ(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the average matrix M₋ₕᵧ with the MatrixElement Uₕ
.
Bramble.M₋ₕ₂
— FunctionM₋ₕ₂(Wₕ::SpaceType)
M₋ₕ₂(Ωₕ::MeshType)
Returns a MatrixElement implementing the average matrix for the mesh grid of Wₕ
, in the z
direction. It is defined as being the (sparse) matrix representation of the linear operator defined by M₋ₕ₂. It also accepts a mesh as argument.
M₋ₕ₂(uₕ::VectorElement)
Returns the average, in the z
direction, of the element uₕ
.
\[\textrm{M}_{-hz} \textrm{u}_h(x_i, y_j,z_l) \vcentcolon = \frac{\textrm{u}_h(x_i, y_j, z_l)+\textrm{u}_h(x_i, y_j, z_{l-1})}{2}\]
M₋ₕ₂(Uₕ::MatrixElement)
Returns a MatrixElement resulting of the multiplication of the average matrix M₋ₕ₂ with the MatrixElement Uₕ
.
Bramble.M₋ₕ
— MethodM₋ₕ(Wₕ::SpaceType)
Returns a tuple of MatrixElements implementing the average operators in the x
, y
, and z
directions.
Bramble.innerₕ
— Functioninnerₕ(uₕ::VectorElement, vₕ::VectorElement)
innerₕ(Uₕ::VecOrMatElem, Vₕ::VecOrMatElem)
Returns the discrete $L^2$ inner product of the grid functions uₕ
and vₕ
- 1D case
\[(\textrm{u}_h, \textrm{v}_h)_h \vcentcolon = \sum_{i=1}^N |\square_{i}| \textrm{u}_h(x_i) \textrm{v}_h(x_i)\]
- 2D case
\[(\textrm{u}_h, \textrm{v}_h)_h \vcentcolon = \sum_{i=1}^{N_x} \sum_{j=1}^{N_y} |\square_{i,j}| \textrm{u}_h(x_i,y_j) \textrm{v}_h(x_i,y_j)\]
- 3D case
\[(\textrm{u}_h, \textrm{v}_h)_h \vcentcolon = \sum_{i=1}^{N_x} \sum_{j=1}^{N_y} \sum_{l=1}^{N_z} |\square_{i,j,l}| \textrm{u}_h(x_i,y_j) \textrm{v}_h(x_i,y_j)\]
Bramble.inner₊ₓ
— Functioninner₊ₓ(uₕ::VecOrMatElem, vₕ::VecOrMatElem)
Returns the discrete modified $L^2$ inner product of the grid functions uₕ
and vₕ
associated with the first variable. It accepts arguments of type VectorElement or MatrixElement, in any order.
For VectorElements, it is defined as
- 1D case
\[(\textrm{u}_h, \textrm{v}_h)_+ \vcentcolon = \sum_{i=1}^{N_x} h_{i} \textrm{u}_h(x_i) \textrm{v}_h(x_i)\]
- 2D case
\[(\textrm{u}_h, \textrm{v}_h)_{+x} \vcentcolon = \sum_{i=1}^{N_x}\sum_{j=1}^{N_y} h_{x,i} h_{y,j+1/2} \textrm{u}_h(x_i,y_j) \textrm{v}_h(x_i,y_j)\]
- 3D case
\[(\textrm{u}_h, \textrm{v}_h)_{+x} \vcentcolon = \sum_{i=1}^{N_x}\sum_{j=1}^{N_y}\sum_{l=1}^{N_z} h_{x,i} h_{y,j+1/2} h_{z,l+1/2} \textrm{u}_h(x_i,y_j,z_l) \textrm{v}_h(x_i,y_j,z_l).\]
Bramble.inner₊ᵧ
— Functioninner₊ᵧ(uₕ::VecOrMatElem, vₕ::VecOrMatElem)
Returns the discrete modified $L^2$ inner product of the grid functions uₕ
and vₕ
associated with the second variable. It accepts
- 2D case
\[(\textrm{u}_h, \textrm{v}_h)_{+y} \vcentcolon = \sum_{i=1}^{N_x}\sum_{j=1}^{N_y} h_{x,i} h_{y,j+1/2} \textrm{u}_h(x_i,y_j) \textrm{v}_h(x_i,y_j)\]
- 3D case
\[(\textrm{u}_h, \textrm{v}_h)_{+y} \vcentcolon = \sum_{i=1}^{N_x}\sum_{j=1}^{N_y}\sum_{l=1}^{N_z} h_{x,i+1/2} h_{y,j} h_{z,l+1/2} \textrm{u}_h(x_i,y_j,z_l) \textrm{v}_h(x_i,y_j,z_l).\]
Bramble.inner₊₂
— Functioninner₊₂(uₕ::VecOrMatElem, vₕ::VecOrMatElem)
Returns the discrete modified $L^2$ inner product of the grid functions uₕ
and vₕ
associated with the z
variable
\[(\textrm{u}_h, \textrm{v}_h)_{+z} \vcentcolon = \sum_{i=1}^{N_x}\sum_{j=1}^{N_y}\sum_{l=1}^{N_z} h_{x,i+1/2} h_{y,j+1/2} h_{z,l} \textrm{u}_h(x_i,y_j,z_l) \textrm{v}_h(x_i,y_j,z_l).\]
Bramble.inner₊
— Functioninner₊(uₕ::VecOrMatElem, vₕ::VecOrMatElem)
inner₊(uₕ::VecOrMatElem, vₕ::VecOrMatElem, Tuple)
inner₊(uₕ::NTuple, vₕ::NTuple)
Returns the discrete modified $L^2$ inner product of the grid functions uₕ
and vₕ
. It accepts arguments of type VectorElement or MatrixElement, in any order.
If the Tuple
argument is given, it returns D
-tuple of all $\textrm{inner}_{x_i,+}$ applied to its input arguments, where D
is the topological dimension of the mesh associated with the elements.
If NTuple
s of VectorElement or MatrixElement are passed as input arguments, it returns the sum of all inner products $(\textrm{u}_h[i],\textrm{v}_h[i])_{+x_i}$.
For VectorElements, the definition is given by
- 1D case
\[(\textrm{u}_h, \textrm{v}_h)_+ \vcentcolon = \sum_{i=1}^{N_x} h_{i} \textrm{u}_h(x_i) \textrm{v}_h(x_i)\]
- 2D case
\[(\textrm{u}_h, \textrm{v}_h)_+ \vcentcolon = (\textrm{u}_h, \textrm{v}_h)_{+x} + (\textrm{u}_h, \textrm{v}_h)_{+y}\]
- 3D case
\[(\textrm{u}_h, \textrm{v}_h)_+ \vcentcolon = (\textrm{u}_h, \textrm{v}_h)_{+x} + (\textrm{u}_h, \textrm{v}_h)_{+y} + (\textrm{u}_h, \textrm{v}_h)_{+z}.\]
See the definitions of inner₊ₓ, inner₊ᵧ and inner₊₂ for more details.
Bramble.normₕ
— Functionnormₕ(uₕ::VectorElement)
Returns the discrete $L^2$ norm of the grid function uₕ
, defined as
\[\Vert \textrm{u}_h \Vert_h \vcentcolon = \sqrt{(\textrm{u}_h, \textrm{u}_h)_h}\]
Bramble.norm₊
— Functionnorm₊(uₕ::VectorElement)
norm₊(uₕ::NTuple{D,VectorElement})
Returns the discrete modified $L^2$ norm of the grid function uₕ
. It also accepts a NTuple
of VectorElements.
For VectorElements uₕ
, it is defined as
\[\Vert \textrm{u}_h \Vert_+ = \sqrt{(\textrm{u}_h,\textrm{u}_h)_+}.\]
and for NTuple
s of VectorElements it returns
\[\Vert \textrm{u}_h \Vert_+ \vcentcolon = \sqrt{ \sum_{i=1}^D(\textrm{u}_h[i],\textrm{u}_h[i])_{+,x_i}}.\]
Bramble.snorm₁ₕ
— Functionsnorm₁ₕ(uₕ::VectorElement)
Returns the discrete version of the standard $H^1$ seminorm of VectorElement uₕ
.
\[|\textrm{u}_h|_{1h} \vcentcolon = \Vert \nabla_h \textrm{u}_h \Vert_h\]
Bramble.norm₁ₕ
— Functionnorm₁ₕ(uₕ::VectorElement)
Returns the discrete version of the standard $H^1$ norm of VectorElement uₕ
.
\[\Vert \textrm{u}_h \Vert_{1h} \vcentcolon = \sqrt{ \Vert \textrm{u}_h \Vert_h^2 + \Vert \nabla_h \textrm{u}_h \Vert_h^2 }\]
Linear and bilinear forms
Bramble.form
— Functionform(Wₕ::SpaceType, Vₕ::SpaceType, f)
Returns a bilinear form from a given expression and trial and test spaces.
form(Wₕ::SType, f::F)
Returns a linear form from a given expression and a test space.
Bramble.assemble
— Functionassemble(a::BilinearForm)
Returns the assembled matrix of a bilinear form.
assemble(a::BilinearForm, bcs::Constraints)
Returns the assembled matrix of a bilinear form with imposed constraints.
assemble(l::LinearForm)
Returns the assembled linear form as a vector.
assemble(l::LinearForm, bcs::Constraints)
Returns the assembled linear form with imposed constraints as a vector of numbers.
Bramble.assemble!
— Functionassemble!(A::AbstractMatrix, a::BilinearForm)
Copies the assembled matrix of a bilinear form to a given matrix.
assemble!(A::AbstractMatrix, a::BilinearFormType, bcs::Constraints)
Copies the assembled matrix of a bilinear form with imposed constraints to a given matrix.
assemble!(x::AbstractVector, l::LinearForm)
In-place assemble of a linear form into a given vector.
assemble!(x::VectorElement, l::LinearForm)
In-place assemble of a linear form into a given VectorElement.
assemble!(vec::AbstractVector, l::LinearForm, bcs::Constraints)
In-place assemble of a linear form with imposed constraints into a given vector.
Bramble.constraints
— Functionconstraints(pairs::NTuple{D,MarkerType}, type::Symbol = :dirichlet)
Returns a Constraints object from a tuple of Markers and a symbol defining the type of boundary condition. Currently, the only supported type is for Dirichlet boundary conditions. The default type is :dirichlet
.
constraints(f::BrambleBareFunction; type::Symbol = :dirichlet)
Returns a Constraints object from a single BrambleBareFunction and a symbol defining the type of boundary condition.
Bramble.symmetrize!
— Functionsymmetrize!(A, F, bcs::Constraints, Ωₕ::MeshType)
After Dirichlet boundary conditions are applied to matrix A
and vector F
using the Constraints object bcs
, this function allows to make A
symmetric, if the original matrix (before applying boundary conditions was symmetric). The algorithm goes as follows: for any given row i
where Dirichlet boundary conditions have been applied
- calculate `dᵢ = cᵢ .* F`, where `cᵢ` is the `i`-th column of `A`;
- replace `F` by substracting `dᵢ` to `F` (except for the `i`-th component)
- replace all elements in the `i`-th column of `A` (except the `i`-th by zero).