Mesh
Bramble.MeshMarkers
— TypeMeshMarkers{D}
Type of dictionary to store the CartesianIndices
associated with a Marker.
Bramble.MeshType
— Type MeshType{D}
Abstract type for meshes. Meshes are only parametrized by their topological dimension D
`.
Base.eltype
— Methodeltype(Ωₕ::MeshType)
eltype(::Type{<:MeshType})
Returns the type of element of the points of the mesh.
Bramble._i2p
— Method_i2p(pts::NTuple{D, Vector{T}}, index::CartesianIndex{D})
Returns a D
tuple with the coordinates of the point in pts
associated with the CartesianIndex
given by ìndex`.
Bramble.dim
— Methoddim(Ωₕ::MeshType)
dim(::Type{<:MeshType})
Returns the tolopogical dimension of Ωₕ
.
Bramble.indices
— Methodindices(Ωₕ::MeshType)
Returns the CartesianIndices
associated with the points of mesh Ωₕ
.
Bramble.marker
— Methodmarker(Ωₕ::MeshType, str::Symbol)
Returns the Marker function with label str
.
Bramble.Mesh1D
— Typestruct Mesh1D{T} <: MeshType{1}
markers::MeshMarkers{1}
indices::CartesianIndices{1}
pts::Vector{T}
npts::Int
end
Structure to create a 1D mesh with npts
points of type T
. The points that define the mesh are stored in pts
and are identified, following the same order, with the indices in indices
. The variable markers
stores, for each Domain marker, the indices satisfying $f(x_i)=0$, where f
is the marker's function.
For future reference, the npts
entries of vector pts
are
\[x_i, \, i=1,\dots,N.\]
Bramble.boundary_indices
— Methodboundary_indices(Ωₕ::Mesh1D)
Returns the indices of the boundary points of mesh Ωₕ
.
Bramble.cell_measure
— Methodcell_measure(Ωₕ::Mesh1D, i)
cell_measure(Ωₕ::Mesh1D, Iterator)
Returns the measure of the cell
\[\square_{i} \vcentcolon = \left[x_i - \frac{h_{i}}{2}, x_i + \frac{h_{i+1}}{2} \right]\]
at CartesianIndex
i
in mesh Ωₕ
, which is given by $h_{i+1/2}$. If the second argument Iterator
is supplied, the function returns a generator iterating over all cell measures.
Bramble.generate_indices
— Methodgenerate_indices(npts::Int)
Returns a CartesianIndices
object for a vector of length npts
.
Bramble.half_points
— Methodhalf_points(Ωₕ::Mesh1D, i)
half_points(Ωₕ::Mesh1D, Iterator)
Returns the average of two neighboring, $x_{i+1/2}$, points in mesh Ωₕ
, at index i
. If the second argument Iterator
is supplied, the function returns a generator iterating over all half points.
\[x_{i+1/2} \vcentcolon = x_i + \frac{h_{i+1}}{2}, \, i=1,\dots,N-1,\]
$x_{N+1/2} \vcentcolon = x_{N}$ and $x_{1/2} \vcentcolon = x_1$.
Bramble.half_spacing
— Methodhalf_spacing(Ωₕ::Mesh1D, i)
half_spacing(Ωₕ::Mesh1D, Iterator)
Returns the indexwise average of the space stepsize, $h_{i+1/2}$, at index i
in mesh Ωₕ
. If the second argument Iterator
is supplied, the function returns a generator iterating over all half spacings.
\[h_{i+1/2} \vcentcolon = \frac{h_i + h_{i+1}}{2}, \, i=1,\dots,N-1,\]
$h_{N+1/2} \vcentcolon = \frac{h_{N}}{2}$ and $h_{1/2} \vcentcolon = \frac{h_1}{2}$.
Bramble.interior_indices
— Methodinterior_indices(Ωₕ::Mesh1D)
Returns the indices of the interior points of mesh Ωₕ
.
Bramble.npoints
— Methodnpoints(Ωₕ::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
Bramble.spacing
— Methodspacing(Ωₕ::Mesh1D, i)
spacing(Ωₕ::Mesh1D, Iterator)
Returns the space stepsize, $h_i$ at index i
in mesh Ωₕ
. If the second argument Iterator
is supplied, the function returns a generator iterating over all spacings.
\[h_i \vcentcolon = x_i - x_{i-1}, \, i=2,\dots,N\]
and $h_1 \vcentcolon = x_2 - x_1$.
Bramble.MeshnD
— Typestruct MeshnD{n,T} <: MeshType{n}
markers::MeshMarkers{n}
indices::CartesianIndices{n,NTuple{n,UnitRange{Int}}}
submeshes::NTuple{n,Mesh1D{T}}
end
Structure to store a cartesian nD-mesh ($2 \leq n \leq 3$) with points of type T
. For efficiency, the mesh points are not stored. Instead, we store the points of the 1D meshes that make up the nD mesh. To connect both nD and 1D meshes, we use the indices in indices
. The Domain markers are translated to markers
as for Mesh1D.
Bramble.MeshnD
— Method(Ωₕ::MeshnD)(i)
Returns the i
-th submesh of Ωₕ
.
Bramble.boundary_indices
— Methodboundary_indices(Ωₕ::MeshnD)
Returns the boundary indices of mesh Ωₕ
.
Bramble.cell_measure
— Methodcell_measure(Ωₕ::MeshnD, idx)
Returns the measure of the cell $\square_{idx}$ centered at the index idx
(can be a CartesianIndex
or a Tuple
):
- 2D mesh
\[ \square_{i,j} \vcentcolon = \left[x_i - \frac{h_{x,i}}{2}, x_i + \frac{h_{x,i+1}}{2} \right] \times \left[y_j - \frac{h_{y,j}}{2}, y_j + \frac{h_{y,j+1}}{2} \right]\]
with area $h_{x,i+1/2} h_{y,j+1/2}$, where idx
= $(i,j)$,
- 3D mesh
\[\square_{i,j,l} \vcentcolon = \left[x_i - \frac{h_{x,i}}{2}, x_i + \frac{h_{x,i+1}}{2}\right] \times \left[y_j - \frac{h_{y,j}}{2}, y_j + \frac{h_{y,j+1}}{2}\right] \times \left[z_l - \frac{h_{z,l}}{2}, z_l + \frac{h_{z,l+1}}{2}\right]\]
with volume $h_{x,i+1/2} h_{y,j+1/2} h_{z,l+1/2}$, where idx
= $(i,j,l)$.
Bramble.generate_indices
— Methodgenerate_indices(nPoints::NTuple)
Returns the CartesianIndices
of a mesh with nPoints[i]
in each direction.
Bramble.half_points
— Methodhalf_points(Ωₕ::MeshnD{D}, idx)
half_points(Ωₕ::MeshnD{D}, Iterator)
Returns a tuple with the half_points, for each submesh, of the points at index idx
. If Iterator
is passed as the second argument, a generator iterating over all half_points
of the mesh is returned.
- 2D mesh, with
idx
= $(i,j)$
\[(x_{i+1/2}, y_{j+1/2})\]
- 3D mesh, with
idx
= $(i,j,l)$
\[(x_{i+1/2}, y_{j+1/2}, z_{l+1/2}).\]
Bramble.half_spacing
— Methodhalf_spacing(Ωₕ::MeshnD, idx)
half_spacing(Ωₕ::MeshnD{D}, Iterator)
Returns a tuple with the half_spacing, for each submesh, at index idx
. If Iterator
is passed as the second argument, a generator iterating over all half_spacing
s of the mesh is returned.
- 2D mesh, with
idx
= $(i,j)$
\[(h_{x,i+1/2}, h_{y,j+1/2})\]
- 3D mesh, with
idx
= $(i,j,l)$
\[(h_{x,i+1/2}, h_{y,j+1/2}, h_{z,l+1/2})\]
Bramble.interior_indices
— Methodinterior_indices(Ωₕ::MeshnD)
Returns the interior indices of mesh Ωₕ
.
Bramble.is_boundary_index
— Methodis_boundary_index(idx, R::CartesianIndices)
Returns true if the index idx
is a boundary index of the mesh with indices stored in R
.
Bramble.npoints
— Methodnpoints(Ωₕ::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 Ωₕ
.
Bramble.spacing
— Methodspacing(Ωₕ::MeshnD, idx::NTuple)
spacing(Ωₕ::MeshnD{D}, Iterator)
Returns a tuple with the spacing, for each submesh, at index idx
. If Iterator
is passed as the second argument, a generator iterating over all spacing
s of the mesh is returned.
- 2D mesh, with
idx
= $(i,j)$
\[(h_{x,i}, h_{y,j}) \vcentcolon = (x_i - x_{i-1}, y_j - y_{j-1})\]
- 3D mesh, with
idx
= $(i,j,l)$
\[(h_{x,i}, h_{y,j}, h_{z,l}) \vcentcolon = (x_i - x_{i-1}, y_j - y_{j-1}, z_l - z_{l-1})\]