Mesh

Base.eltypeMethod
eltype(Ωₕ::AbstractMeshType)

Returns the element type of the points of Ωₕ.

This function is a required part of the AbstractMeshType interface.

source
Bramble.boundary_indicesMethod
boundary_indices(Ωₕ::AbstractMeshType)

Returns an iterator or collection of all CartesianIndex points on the boundary of mesh Ωₕ.

source
Bramble.boundary_indicesMethod
boundary_indices(idxs::CartesianIndices)

Returns all boundary facets of a CartesianIndices domain as a tuple of CartesianIndices. Each element of the returned tuple represents a distinct boundary section, such as a face or edge of the domain.

Example

julia> domain = CartesianIndices((2, 2));
	   boundary_indices(domain)
(CartesianIndices((1:1, 1:2)), CartesianIndices((2:2, 1:2)), CartesianIndices((1:2, 1:1)), CartesianIndices((1:2, 2:2)))
source
Bramble.cell_measureMethod
cell_measure(Ωₕ::AbstractMeshType, idx)

Returns the measure of the cell $\square_{idx}$ centered at the index idx (can be a CartesianIndex or a Tuple):

  • 1D mesh, with idx = $(i,)$

\[\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}$

  • 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)$.

source
Bramble.change_points!Method
change_points!(Ωₕ::AbstractMeshType, [Ω::Domain], pts)

Changes the coordinates of the internal points of the mesh Ωₕ to the new coordinates specified in pts. This function assumes the points in pts are ordered and that the first and last of them coincide with the bounds of the Ω. if the domain Ω is passed as an argument, the markers of the mesh are also recalculated after this change.

source
Bramble.forward_spacingMethod
forward_spacing(Ωₕ::AbstractMeshType, idx)

Returns a tuple with the forward_spacing, for each submesh, at index idx.

  • 1D mesh, with idx = $(i,)$ or $i$

\[h_{x,i} \vcentcolon = x_{i+1} - x_i, \, i=1,\dots,N-1\]

and $h_{x,N} \vcentcolon = x_N - x_{N-1}$

  • 2D mesh, with idx = $(i,j)$

\[(h_{x,i}, h_{y,j}) \vcentcolon = (x_{i+1} - x_i, y_{j+1} - y_j)\]

  • 3D mesh, with idx = $(i,j,l)$

\[(h_{x,i}, h_{y,j}, h_{z,l}) \vcentcolon = (x_{i+1} - x_i, y_{j+1} - y_j, z_{l+1} - z_l)\]

source
Bramble.generate_indicesMethod
generate_indices([::Dimension], pts)

Returns the CartesianIndices of a mesh with pts[i] in each direction or just pts, if the argument is an Int.

source
Bramble.half_pointsMethod
half_points(Ωₕ::AbstractMeshType)

Returns the half points, for each submesh.

  • 1D mesh (with idx=$(i,)$ or $i$)

\[x_{i+1/2} \vcentcolon = x_i + \frac{h_{i+1}}{2}, \, i=1,\dots,N-1,\]

\[x_{N+1/2} \vcentcolon = x_{N}, \, x_{1/2} \vcentcolon = x_1.\]

  • 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}).\]

source
Bramble.half_spacingsMethod
half_spacings(Ωₕ::AbstractMeshType)

Returns the indexwise average of the space stepsize, for each submesh.

  • 1D mesh, with idx = $(i,)$ or $i$

\[h_{x,i+1/2} \vcentcolon = \frac{h_{x,i} + h_{x,i+1}}{2}, \, i=1,\dots,N-1,\]

\[h_{x,N+1/2} \vcentcolon = \frac{h_{N}}{2}`` and ``h_{x,1/2} \vcentcolon = \frac{h_1}{2}``.\]

  • 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})\]

source
Bramble.index_in_markerMethod
index_in_marker(Ωₕ::AbstractMeshType, label::Symbol)

Returns the BitVector associated with the marker with label of mesh Ωₕ.

source
Bramble.indicesMethod
indices(Ωₕ::AbstractMeshType)

Returns the CartesianIndices associated with the points of mesh Ωₕ. This function is a required part of the AbstractMeshType interface. Any concrete subtype of AbstractMeshType must implement this method and have a field called indices of type CartesianIndices.

source
Bramble.interior_indicesMethod
interior_indices(Ωₕ::AbstractMeshType)

Returns a CartesianIndices object representing the interior region of the mesh Ωₕ.

source
Bramble.interior_indicesMethod
interior_indices(indices::CartesianIndices)

Computes the CartesianIndices representing the interior of a given domain, excluding all boundary points. This is achieved by shrinking the index range in each dimension by one from both ends. Dimensions with a length of one or less are returned unchanged.

Examples

julia> domain = CartesianIndices((3, 3)); interior_indices(domain)
CartesianIndices((2:2, 2:2))
julia> domain_2d_line = CartesianIndices((1, 5));
	   interior_indices(domain_2d_line)
CartesianIndices((1:1, 2:4))
source
Bramble.is_boundary_indexMethod
is_boundary_index(Ωₕ::AbstractMeshType, idx)

Checks if the CartesianIndex idx lies on the boundary of the mesh Ωₕ. This delegates to the implementation for CartesianIndices.

source
Bramble.is_boundary_indexMethod
is_boundary_index(idxs::CartesianIndices, idx)

Checks if a given index idx lies on the boundary of a CartesianIndices domain.

Example

julia> domain = CartesianIndices((3, 4));
	   is_boundary_index(domain, (1, 2));
true

julia> is_boundary_index(domain, (2, 2))
false
source
Bramble.pointMethod
point(Ωₕ::AbstractMeshType, idx)

Returns the tuple with the point from the mesh corresponding to index idx. See points.

source
Bramble.setMethod
set(Ωₕ::AbstractMeshType)

Returns the set of the domain over which the mesh Ωₕ is defined.

source
Bramble.spacingMethod
spacing(Ωₕ::AbstractMeshType, idx)

Returns a tuple with the spacing, for each submesh, at index idx.

  • 1D mesh, with idx = $(i,)$ or $i$

\[h_{x,i} \vcentcolon = x_i - x_{i-1}, \, i=2,\dots,N\]

and $h_{x,1} \vcentcolon = x_2 - x_1$

  • 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})\]

source
Bramble.MeshMarkersType
mutable struct Dict{Symbol, BitVector} <: AbstractDict{Symbol, BitVector}

Efficient storage type for mesh markers as a Dict of Symbols. For each label, a BitVector is assigned that determines, for a given index, if the corresponding geometric point is identified by the marker.

source
Bramble.__process_condition!Method
__process_condition!(mesh_marker, identifier, Ωₕ)

Core logic for evaluating a function-based (level-set) marker.

It iterates through every point in the mesh, evaluates the identifier function at that point's coordinates, and sets the marker to true if the function returns true.

source
Bramble._init_mesh_markersMethod
_init_mesh_markers(Ωₕ, domain_markers)

Internal helper function to create and initialize the MeshMarkers dictionary.

This function extracts all unique labels from the provided DomainMarkers object, which can come from symbol-, tuple-, or condition-based markers. It then prepares a MeshMarkers dictionary where each label is a key associated with a BitVector of falses, ready to be populated.

source
Bramble._mark_indices!Method
_mark_indices!(marker_set, linear_indices, indices_to_mark)

A utility function to efficiently update a boolean marker vector.

It sets the value to true at the linear positions corresponding to the CartesianIndex or collection of CartesianIndices provided in indices_to_mark.

source
Bramble._set_markers_symbols!Method
_set_markers_symbols!(mesh_markers, symbols, Ωₕ)

Processes markers that are identified by predefined symbols (e.g., :left, :top) or sets of those symbols.

source
Bramble.boundary_symbol_to_cartesianMethod
boundary_symbol_to_cartesian(indices::CartesianIndices)

Returns a named tuple connecting the facet labels of a set to the corresponding CartesianIndices.

Examples

julia> boundary_symbol_to_cartesian(CartesianIndices((1:3, 1:4)))
(left = CartesianIndices((1:1, 1:4)),
 right = CartesianIndices((3:3, 1:4)),
 top = CartesianIndices((1:3, 4:4)),
 bottom = CartesianIndices((1:3, 1:1)))
julia> boundary_symbol_to_cartesian(CartesianIndices((1:10, 1:20, 1:15)))
(left = CartesianIndices((1:10, 1:1, 1:15)),
 right = CartesianIndices((1:10, 20:20, 1:15)),
 top = CartesianIndices((1:10, 1:20, 15:15)),
 bottom = CartesianIndices((1:10, 1:20, 1:1)),
 front = CartesianIndices((10:10, 1:20, 1:15)),
 back = CartesianIndices((1:1, 1:20, 1:15)))
source
Bramble.boundary_symbol_to_dictMethod
boundary_symbol_to_dict(indices::CartesianIndices)

Returns a dictionary connecting the facet labels of a set to the corresponding `CartesianIndices` (see [`boundary_symbol_to_cartesian`](@ref)).
source
Bramble.process_label_for_mesh!Method
process_label_for_mesh!(npts, markers_mesh, set_labels)

Initializes boolean vectors for a given set of labels within the mesh markers dictionary.

For each label in set_labels, this function creates a BitVector of length npts (the total number of points in the mesh), initializes it with all false values, and assigns it to the corresponding key in the markers_mesh dictionary. This prepares the storage for later marking which points belong to which labeled region.

Arguments

  • npts: The total number of points in the mesh.
  • markers_mesh: The MeshMarkers dictionary to be modified in-place.
  • set_labels: An iterator or collection of Symbol labels to initialize.
source
Bramble.set_markers!Method
set_markers!(Ωₕ::AbstractMeshType, domain_markers)

Populates the markers of a mesh Ωₕ based on a DomainMarkers object.

This is the main entry point for applying user-defined labels to the mesh points. It initializes the marker storage and then delegates to specialized functions for each type of marker identifier (symbols, tuples of symbols, and functions).

source
Bramble.Mesh1DType
mutable struct Mesh1D{BT<:Bramble.Backend, CI, VT<:(AbstractVector), T} <: Bramble.AbstractMeshType{1}

A mutable structure representing a 1D mesh.

This struct holds all the geometric and topological information for a one-dimensional grid. It includes the coordinates of the grid points (pts), the underlying geometric interval (set), and a dictionary of markers for labeling specific points or regions. Key geometric quantities like cell centers (half_pts) and cell widths (half_spacings) are pre-computed and stored for efficiency, which is particularly useful in numerical methods like the finite volume method.

The struct is mutable to allow for in-place modifications, such as mesh refinement.

Fields

  • set: the geometric domain, a 1D CartesianProduct (interval), over which the mesh is defined.

  • markers: a dictionary mapping Symbol labels to BitVectors, marking specific points on the mesh.

  • indices: the CartesianIndices of the grid, allowing for array-like iteration and indexing over the points.

  • backend: the computational backend used for linear algebra operations.

  • pts: a vector holding the coordinates of the grid points, $x_i$.

  • half_pts: a vector of pre-computed cell centers (midpoints), $x_{i+1/2}$.

  • half_spacings: a vector of pre-computed cell widths, $h_{i+1/2}$.

  • collapsed: a boolean flag indicating if the domain is degenerate (a single point).

For future reference, the entries of vector pts are

\[x_i, \, i=1,\dots,N.\]

source
Bramble.MeshnDType
mutable struct MeshnD{D, BT<:Bramble.Backend, CI, M1T<:Bramble.AbstractMeshType{1}, T} <: Bramble.AbstractMeshType{D}

Structure to store a cartesian nD-mesh ($2 \leq n \leq 3$). 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 DomainMarkers are translated to MeshMarkers as for Mesh1D.

Fields

  • set: the D-dimensional CartesianProduct (hyperrectangle) defining the geometric domain.

  • markers: a dictionary mapping Symbol labels to BitVectors, marking grid points.

  • indices: the CartesianIndices for the full D-dimensional grid, allowing for multi-dimensional indexing.

  • backend: the computational backend used for linear algebra operations.

  • submeshes: a tuple of D Mesh1D objects, representing the grid along each spatial dimension.

source
Bramble._meshMethod
_mesh(Ω::Domain, npts, unif, backend)

Internal constructor for a D-dimensional, tensor-product MeshnD.

This function orchestrates the creation of a structured multidimensional mesh. It first builds the underlying 1D submeshes for each dimension and then combines them into a single MeshnD object. It also handles the important edge case of "collapsed" dimensions (where an interval is just a point), forcing the number of grid points in that dimension to be 1.

Arguments

  • Ω: The D-dimensional Domain to be meshed.
  • npts: An NTuple{D, Int} specifying the number of points in each dimension.
  • unif: An NTuple{D, Bool} specifying if the grid is uniform in each dimension.
  • backend: The linear algebra backend.
source
Bramble.submeshesMethod
submeshes(Ω::Domain, npts, unif, backend)

Creates the component 1D submeshes for a tensor-product grid.

This function takes a D-dimensional Domain and generates a tuple of D independent Mesh1D objects. Each submesh corresponds to one of the spatial dimensions of the original domain.

Arguments

  • Ω: The D-dimensional Domain.
  • npts: A tuple containing the number of points for each dimension.
  • unif: A tuple of booleans indicating if the grid is uniform in each dimension.
  • backend: The linear algebra backend.
source