References

WavePropBase.ambient_dimensionFunction
ambient_dimension(x)

Dimension of the ambient space where x lives. For geometrical objects this can differ from its geometric_dimension; for example a triangle in ℝ³ has ambient dimension 3 but geometric dimension 2, while a curve in ℝ³ has ambient dimension 3 but geometric dimension 1.

source
WavePropBase.boundaryFunction
boundary(ω)

Return the boundary of ω. For a mesh element gives the d-1 dimensional elements composing its boundary, while for an entity gives the corresponding d-1 dimensional entities.

source
WavePropBase.coordsFunction
coords(x)

Return an SVector giving a cartesian coordinate used to sort x. For points, this is simply the coordinates of the point, while for elements this may be the center of the element. You must overload this function for your own type if you want e.g. the clustering algorithms to work for x.

source
WavePropBase.dimensionFunction
dimension(space)

The length of a basis for space; i.e. the number of linearly independent elements required to span space.

source
WavePropBase.geometric_dimensionFunction
geometric_dimension(x)
geometric_dimension(Ω::Domain)

Number of degrees of freedom necessary to locally represent the geometrical object. For example, lines have geometric dimension of 1 (whether in ℝ² or in ℝ³), while surfaces have geometric dimension of 2.

When the argument is a Domain, return the largest geometric dimension encoutered.

source
WavePropBase.jacobianFunction
jacobian(f,x)

Given a (possibly vector-valued) function f : 𝐑ᵐ → 𝐅ᵐ, return the m × n matrix Aᵢⱼ = ∂fᵢ/∂xⱼ.

source
WavePropBase.normalFunction
normal(el,x̂)
normal(jac::SMatrix)

The unit normal vector at coordinate , guaranteed to be orthogonal to all columns of jacobian(el,x).

source
WavePropBase.Utils.STypeType
const SType{T} = Union{T,Type{T}}

Union type of T and its data type Type{T}. Used to simplify methods defined on singleton types where both foo(::T) and foo(::Type{T}) are required.

source
WavePropBase.Utils.cart2sphMethod
cart2sph(x,y,z)

Map cartesian coordinates x,y,z to spherical ones r, θ, φ representing the radius, elevation, and azimuthal angle respectively. The convention followed is that 0 ≤ θ ≤ π and -π < φ ≤ π.

source
WavePropBase.Utils.diagonalblockmatrix_to_matrixMethod
diagonalblockmatrix_to_matrix(A::Matrix{B}) where {B<:SMatrix}

Convert a diagonal block matrix A::AbstractVector{B}, where A is the list of diagonal blocks and B<:SMatrix, to the equivalent SparseMatrixCSC{T}, where T = eltype(B).

source
WavePropBase.Utils.matrix_to_blockmatrixMethod
matrix_to_blockmatrix(A::Matrix,B)

Convert a Matrix{T} to a Matrix{B}, where B<:Type{SMatrix}. The element type of B must match that of A, and the size of A must be divisible by the size of B along each dimension.

source
WavePropBase.Utils.sort_by_typeMethod
sort_by_type(v)

Sort the elements of v into vectors vi according to their type. Return a Dict{DataType,Vector} mapping each type to a vector of that type.

Examples

v = [1,"a",3,"b"]
dict = sort_by_type(v)
source
WavePropBase.Utils.sph2cartMethod
sph2cart(x,y,z)

Map spherical coordinates r,θ,φ representing the radius, elevation, and azimuthal angle respectively, to cartesian ones x, y, z .

source
WavePropBase.Utils.vector_to_blockvectorMethod
vector_to_blockvector(A::Vector,B)

Convert a Vector{T} to a Vector{B}, where B<:Type{SVector}. The element type of B must match that of A, and the size of A must be divisible by the size of B along each dimension.

source
WavePropBase.Utils.@interfaceMacro
@interface f [n=1]

Declare that the function f is an interface function. The call f(args...) resolves to M.f(args...) where M is parent module of the args[n] object.

The somewhat contrived example below illustrates how this can be used to have a generic method defined in module A applied to a type defined on module B which is independent of A but which implements the interface function f:

module A
    using WavePropBase.Utils
    Utils.@interface foo
    # a method which works on any type `x` implementing the `foo` function
    do_work(x) = 2*foo(x)
end

module B
    struct Foo end
    foo(x::Foo) = 1
end

using .A
using .B
foo = B.Foo()
A.do_work(foo)

# output

2

Note that if in the example above module A implements a generic version of foo, the call A.do_work(foo) would use that method instead based on the dispatch rules.

source
WavePropBase.Geometry.TAGSConstant
const TAGS::Dict{Int,Vector{Int}}

Global dictionary storing the used entity tags (the value) for a given dimension (the key).

source
WavePropBase.Geometry.ComplexPoint2DType
const ComplexPoint2D
const ComplexPoint2D(x1, x2)
const ComplexPoint2D(x::NTuple{2, ComplexF64})

A complex 2D point, stored in a StaticArray. ComplexPoint2D = SVector{2, ComplexF64}.

source
WavePropBase.Geometry.ComplexPoint3DType
const ComplexPoint3D
const ComplexPoint3D(x1, x2, x3)
const ComplexPoint3D(x::NTuple{3, ComplexF64})

A complex 3D point, stored in a StaticArray. ComplexPoint3D = SVector{3, ComplexF64}.

source
WavePropBase.Geometry.ElementaryEntityType
struct ElementaryEntity <: AbstractEntity

The most basic representation of an AbstractEntity.

Fields:

  • dim::UInt8: the geometrical dimension of the entity (e.g. line has dim=1, surface has dim=2, etc)
  • tag::Int64: an integer tag associated to the entity
  • boundary::Vector{AbstractEntity}: the entities of dimension dim-1 forming the entity's boundary
source
WavePropBase.Geometry.Point2DType
const Point2D
const Point2D(x1, x2)
const Point2D(x::NTuple{2, Float64})

A point in 2D space, stored in a StaticArray. Point2D = SVector{2, Float64}.

source
WavePropBase.Geometry.Point3DType
const Point3D
const Point3D(x1, x2, x3)
const Point3D(x::NTuple{3, Float64})

A point in 3D space, stored in a StaticArray. Point3D = SVector{3, Float64}.

source
WavePropBase.Geometry.PointEntityType
PointEntity{N,T} <: AbstractEntity

Zero-dimension geometrical entity. As a subtype of [AbstractEntity],(@ref) the (dim,tag) of all created point entities get added to the global ENTITIES. Intended usage is to build higher dimensionsional entities, and not to represent regular points such as grid points.

source
Base.:==Method
==(Ω1::AbstractEntity,Ω2::AbstractEntity)

Two entities are considered equal geometric_dimension(Ω1)==geometric_dimension(Ω2) and abs(tag(Ω1))=abs(tag(Ω2)). For entities of co-dimension one, the sign of tag(Ω) is used to determine the orientation of the normal vector.

Notice that this implies dim and tag of an elementary entity should uniquely define it (up to the sign of tag), and therefore global variables like TAGS are needed to make sure newly created AbstractEntity have a new (dim,tag) identifier.

source
Base.:==Method
===(Ω1::Domain,Ω2::Domain)

Two Domains are equal if all their entities are equal (regardless of order).

source
Base.inMethod
in(ω::ElementaryEntity,Ω::Domain)

Check whether an ElementaryEntity belongs to a Domain by recursively checking whether it belongs to its boundary.

source
Base.iterateFunction
iterate(Ω::Domain)

Iterating over a domain means iterating over its entities.

source
Base.keysMethod

Return all tags of the elementary entities in the domain Ω corresponding to the dimension d.

source
Base.keysMethod

Return all tags of the elementary entities in the domain Ω corresponding to the dimensions contained in dims.

source
Base.lengthMethod
length(Ω:::Domain)

The length of a domain corresponds to the number of elementary entities that make it.

source
WavePropBase.Geometry.assertequaldimMethod
assertequaldim(Ω1::Domain,Ω2::Domain)

Check that two domains have same dimension.

If one of the domain (or both) are empty, the assertion is assumed to be true.

source
WavePropBase.Geometry.new_tagMethod
new_tag(dim)

Generate a unique tag for an AbstractEntity of dimension dim.

The implementation consists of adding one to the maximum value of TAGS[dim]

See also: TAGS.

source
WavePropBase.keyMethod
key(e::AbstractEntity)

The (dim,tag) pair used as a key to identify geometrical entities.

source
WavePropBase.TreesModule
module Trees

Module defining some commonly used tree structures, as well as an interface for talking about trees inside the WaveProp organization.

source
WavePropBase.Trees.AbstractSplitterType
abstract type AbstractSplitter

An AbstractSplitter is used to split a ClusterTree. The interface requires the following methods:

  • should_split(clt,splitter) : return a Bool determining if the ClusterTree should be further divided
  • split!(clt,splitter) : perform the splitting of the ClusterTree handling the necessary data sorting.

See GeometricSplitter for an example of an implementation.

source
WavePropBase.Trees.CardinalitySplitterType
struct CardinalitySplitter <: AbstractSplitter

Used to split a ClusterTree along the largest dimension if length(tree)>nmax. The split is performed so the data is evenly distributed amongst all children.

source
WavePropBase.Trees.ClusterTreeType
mutable struct ClusterTree{T,S,D}

Tree structure used to cluster elements of type T into containers of type S. The method center(::T)::SVector is required for the clustering algorithms. An additional data field of type D can be associated with each node to store node-specific information (it defaults to D=Nothing).

Fields:

  • _elements::T : vector containing the sorted elements.
  • container::S : container for the elements in the current node.
  • index_range::UnitRange{Int} : indices of elements contained in the current node.
  • loc2glob::Vector{Int} : permutation from the local indexing system to the

original (global) indexing system used as input in the construction of the tree.

  • children::Vector{ClusterTree{N,T,D}}
  • parent::ClusterTree{N,T,D}
  • data::D : generic data field of type D.
source
WavePropBase.Trees.ClusterTreeMethod
ClusterTree(elements,splitter;[copy_elements=true])
ClusterTree{D}(points,splitter;[copy_elements=true])

Construct a ClusterTree from the given elements using the splitting strategy encoded in splitter. If copy_elements is set to false, the elements argument are directly stored in the ClusterTree and are permuted during the tree construction.

source
Base.parentMethod
parent(t::AbstractTree)

The node's parent. If t is a root, then parent(t)==t.

source
WavePropBase.Trees._binary_split!Method
_binary_split!(cluster::ClusterTree,dir,pos;parentcluster=cluster)
_binary_split!(cluster::ClusterTree,f;parentcluster=cluster)

Split a ClusterTree into two, sorting all elements in the process. For each resulting child assign child.parent=parentcluster.

Passing a dir and pos arguments splits the bounding_box box of node along direction dir at position pos, then sorts all points into the resulting left/right nodes.

If passed a predicate f, each point is sorted according to whether f(x) returns true (point sorted on the left node) or false (point sorted on the right node). At the end a minimal HyperRectangle containing all left/right points is created.

source
WavePropBase.Trees.depthFunction
depth(tree::AbstractTree,acc=0)

Recursive function to compute the depth of node in a a tree-like structure.

Overload this function if your structure has a more efficient way to compute depth (e.g. if it stores it in a field).

source
WavePropBase.Trees.filter_treeFunction
filter_tree(f,tree,isterminal=true)

Return a vector containing all the nodes of tree such that filter(node)==true. The argument isterminal can be used to control whether to continue the search on children of nodes for which f(node)==true.

source
WavePropBase.Trees.loc2globMethod
loc2glob(clt::ClusterTree)

The permutation from the (local) indexing system of the elements of the clt to the (global) indexes used upon the construction of the tree.

source
WavePropBase.InterpolationModule
module Interpolation

Module implementing various interpolation schemes commonly used in finite element and boundary element methods.

source
WavePropBase.Interpolation.AbstractElementType
abstract type AbstractElement{D,T}

Elements given by a fixed interpolation schemes mapping points on the the domain D<:AbstractReferenceShape (of singleton type) to a value of type T.

Instances el of AbstractElement are expected to implement:

  • el(x̂): evaluate the interpolation scheme at the (reference) coordinate x̂ ∈ D.
  • jacobian(el,x̂) : evaluate the jacobian matrix of the interpolation at the (reference) coordinate x ∈ D.
Note

For performance reasons, both el(x̂) and jacobian(el,x̂) should take as input a StaticVector and output a static vector or static array.

source
WavePropBase.Interpolation.LagrangeElementType
struct LagrangeElement{D,Np,T} <: AbstractElement{D,T}

Standard element over D <: AbstractReferenceShape commonly used in finite element methods. The underlying polynomial space is Pk{D,K}, and its interpolant maps the Np reference_nodes in D to Np values of type T stored in the field vals.

source
WavePropBase.Interpolation.PkType
struct Pk{D,K} <: AbstractPolynomialSpace{D}

The space of all polynomials over D of degree ≤K.

When D is a hypercube in d dimensions, the precise definition is Pk{D,K} = span{𝐱ᶿ : max(θ)≤ K}; when D is a d-dimensional simplex, the space is Pk{D,K} = span{𝐱ᶿ : sum(θ)≤ K}, where θ ∈ 𝐍ᵈ is a multi-index.

source
WavePropBase.Interpolation.TensorLagInterpType
struct TensorLagInterp{N,Td,T}

Generic Lagrange interpolation over an N-dimensional tensor grid. The implementation uses a multidimensional generalization of the barycentric formula.

The main constructor takes an SVector{N,Vector{Td}} containig the N one-dimensional nodes and an Array{N,T} of the function vals at the tensor product grid formed by the one-dimensional nodes.

Examples:

nx = 10
ny = 12
x   = [0.5+0.5cos((2k-1)*π/2nx) for k in 1:nx] # Chebyshev nodes
y   = [0.5+0.5cos((2k-1)*π/2ny) for k in 1:ny] # Chebyshev nodes
f   = (x) -> cos(x[1]*x[2])
vals = [f((x,y)) for x in x, y in y]
p   = TensorLagInterp(SVector(x,y),vals)
p((0.1,0.2)) ≈ f((0.1,0.2))
source
Base.splitMethod
split(rec::AbstractHyperRectangle,[axis]::Int,[place])

Split a hyperrectangle in two along the axis direction at the position place. Returns a tuple with the two resulting hyperrectangles.

When no place is given, defaults to splitting in the middle of the axis.

When no axis and no place is given, defaults to splitting along the largest axis.

source
WavePropBase.Interpolation.cheb2nodesMethod
cheb2nodes(n,a,b)

Return the n Chebyshev points of the second kind on the interval [a,b]. The nodes are nested in the following sense: cheb2nodes(n,a,b) == cheb2nodes(2n-1,a,b)[1:2:end].

source
WavePropBase.Interpolation.lagrange_basisMethod
lagrange_basis(nodes,[sp::AbstractPolynomialSpace])

Return the set of n polynomials in sp taking the value of 1 on node i and 0 on nodes j ≂̸ i for 1 ≤ i ≤ n. For N-dimensional tensor-product nodes represented in the form of an SVector{N,Vector{T}}, the argument sp may be ommited.

Danger

It is assumed that the value of a function on nodes uniquely determine a polynomial in sp.

source
WavePropBase.Interpolation.reference_nodesFunction
reference_nodes(el::LagrangeElement)

Return the reference nodes on domain(el) used for the polynomial interpolation. The function values on these nodes completely determines the interpolating polynomial.

We use the same convention as gmsh for defining the reference nodes and their order (see node ordering on gmsh documentation).

source
WavePropBase.distanceMethod
distance(x::SVector,r::HyperRectangle)

The (minimal) Euclidean distance between the point x and any point y ∈ r.

source
WavePropBase.distanceMethod
distance(r1::AbstractHyperRectangle,r2::AbstractHyperRectangle)

The (minimal) Euclidean distance between a point x ∈ r1 and y ∈ r2.

source
WavePropBase.IntegrationModule
module Integration

Methods for integrating over instances of AbstractReferenceShape.

Besides some standard quadrature rules used to integrate smooth functions, it also defines singular integration routines useful for (weakly) singular integrands.

source
WavePropBase.Integration.AbstractQuadratureRuleType
abstract type AbstractQuadratureRule{D<:AbstractReferenceShape}

A quadrature rule for integrating a function over the domain D.

An instance q of AbstractQuadratureRule{D} is expected to implement the following methods:

  • q() : return the nodes x and weights w of the quadrature rule on the reference domain D. For performance reasons, the result should depend only on the type of q.
source
WavePropBase.Integration.DuffyType
struct Duffy <: AbstractSingularityHandler{RefereceTriangle}

Change of variables mapping the ReferenceSquare to the RefereceTriangle with the property that the jacobian vanishes at the (1,0) vertex of the triangle.

Useful for integrating functions with a singularity on the (1,0) edge of the reference triangle.

source
WavePropBase.Integration.FejerType
struct Fejer{N}

N-point Fejer's first quadrature rule for integrating a function over [0,1]. Exactly integrates all polynomials of degree ≤ N-1.

source
WavePropBase.Integration.KressType
struct Kress{P} <: AbstractSingularityHandler{ReferenceLine}

Change of variables mapping [0,1] to [0,1] with the property that the first P-1 derivatives of the transformation vanish at x=0.

source
WavePropBase.Integration.KressPType
struct KressP{P} <: AbstractSingularityHandler{ReferenceLine}

Like Kress, this change of variables maps the interval [0,1] onto itself, but the first P derivatives of the transformation vanish at both endpoints.

This change of variables can be used to periodize integrals in the following sense. Suppose we wish to compute the integral of f(x) from 0 to 1 where f is not a1-periodic function. If ϕ is an object of type KressP, then using it as a change of variables in the integration yields a similar integral from 0 to 1 (the interval 0≤0≤1 is mappend onto itself), but with integrand given by g(x) = f(ϕ(x))ϕ'(x). Since ϕ' vanishes (together with P of its derivatives), the function g(x) is now periodic (up to derivatives of order up to P) at the endpoints. Thus quadrature rules designed for periodic functions like the TrapezoidalOpen can be used to obtain high order convergence of g, which in turn yields a modified quadrature rule when viewed as a quadrature rule for f.

source
WavePropBase.Integration.SingularQuadratureRuleType
SingularQuadratureRule{D,Q,S} <: AbstractQuadratureRule{D}

A quadrature rule over D intended to integrate functions which are singular at a known point s ∈ D.

A singular quadrature is rule is composed of a regular quadrature rule (e.g. GaussLegendre) and a AbstractSingularityHandler to transform the regular quadrature. The regular quadrature rule generates nodes and weights on the domain(sing_handler), and those are mapped into an appropriate quadrature over D = range(sing_handler) using the singularity handler.

source
WavePropBase.Integration.TrapezoidalType
struct Trapezoidal{N} <: AbstractQuadratureRule{ReferenceLine}

Closed N-point trapezoidal rule for integrating a function over the interval [0,1].

Examples:

q    = Trapezoidal(10)
f(x) = exp(x)*cos(x)
integrate(f,q)
source
WavePropBase.Integration.integrateMethod
integrate(f,q::AbstractQuadrature)
integrate(f,x,w)

Integrate the function f using the quadrature rule q. This is simply sum(f.(x) .* w), where x and w are the quadrature nodes and weights, respectively.

source
WavePropBase.Mesh.ElementIteratorType
struct ElementIterator{E,M}

Return an iterator for all elements of type E on a mesh of type M.

Besides the methods listed in the iterator iterface of Julia, some functions also require the getindex(iter,i::Int) method for accessing the i-th element directly.

source
WavePropBase.Mesh.GenericMeshType
struct GenericMesh{N,T} <: AbstractMesh{N,T}

Data structure representing a generic mesh in an ambient space of dimension N, with data of type T.

source
WavePropBase.Mesh.SubMeshType
struct SubMesh{N,T} <: AbstractMesh{N,T}

Create a view of a parent mesh over a given domain.

A submesh implements the interface for AbstractMesh; therefore you can iterate over elements of the submesh just like you would with a mesh.

source
WavePropBase.Mesh.UniformCartesianMeshType
struct UniformCartesianMesh{N,T} <: AbstractMesh{N,T}

An N-dimensional cartesian grid given as the tensor-product of N one-dimensional LinRange{T} grids.

Iterating over a UniformCartesianMesh generates the elements which compose the mesh; i.e. the HyperRectangle cells.

source
WavePropBase.Mesh.dom2eltMethod
dom2elt(m::GenericMesh,Ω,E)

Compute the element indices idxs of the elements of type E composing Ω, so that m[E][idxs] gives all the elements of type E meshing Ω.

source
WavePropBase.Mesh.dom2eltMethod
dom2elt(m::GenericMesh,Ω)

Return a Dict with keys being the element types of m, and values being the indices of the elements in Ω of that type.

source
WavePropBase.Mesh.dom2eltMethod
dom2elt(m::SubMesh,[E])

A dictionary with keys being the element types of m, and values being the element indices in the parent mesh. If a type E is given, return the values associated with that key.

source
WavePropBase.IO.etype_to_vtk_cell_typeConstant
const etype_to_vtk_cell_type

OrderedDictionary mapping internal element types to a tuple containing: - the corresponding WriteVTK cell types (following the convention chosen by VTK, see below); - the indices in the elements column that defines the element. This is because we want to support the export of more than just the flat elements available in the VTK specification, hence which may require a conversion of some sort.

See VTK specification [Fig. 2] on http://www.vtk.org/VTK/img/file-formats.pdf

  • VTK_VERTEX (=1)
  • VTKPOLYVERTEX (=2)
  • VTK_LINE (=3)
  • VTKPOLYLINE (=4)
  • VTK_TRIANGLE (=5)
  • VTKTRIANGLESTRIP (=6)
  • VTK_POLYGON (=7)
  • VTK_PIXEL (=8)
  • VTK_QUAD (=9)
  • VTK_TETRA (=10)
  • VTK_VOXEL (=11)
  • VTK_HEXAHEDRON (=12)
  • VTK_WEDGE (=13)
  • VTK_PYRAMID (=14)
source
WavePropBase.IO._vtk_cellsFunction
_vtk_cells(mesh::GenericMesh, E::DataType)
_vtk_cells(mesh::GenericMesh, Ω::Domain)

Creates the vector of all elements contained in the mesh in the format required by WriteVTK for a particular element type E<:AbstractElement or a Domain instance.

source
WavePropBase.IO.vtk_mesh_fileMethod
vtk_mesh_file(mesh::GenericMesh[, Ω::Domain], name::String)

Creates a VTK file (.vtu) with name name containing the mesh information. It is possible to export only a Domain (i.e. only a part of the mesh).

Note that all the heavy lifting is done by the package WriteVTK. We refer to its documentation for more information.

Warning: the output is not an actual file (on disk). To save it, simply write:

vtk_mesh_file(mesh, "my_mesh") |> vtk_save

To add some data (scalar or vector values at the mesh nodes or mesh cells) for visualization purposes, you can do for instance:

vtkfile = vtk_mesh_file(mesh, name)
vtkfile["my_point_data", VTKPointData()] = pdata
vtkfile["my_cell_data", VTKCellData()] = cdata
vtk_save(vtkfile)

It is possible also to export a partition Ωs::Vector{Domain} using Multiblock files (.vtm), for instance like so

vtmfile = vtk_multiblock(name)
for (Ω, pdata) in zip(Ωs, pdatas)
    vtkfile = vtk_mesh_file(mesh, Ω, name)
    vtkfile["my_point_data", VTKPointData()] = pdata
end
vtk_save(vtmfile)

To save a sequence of solutions (time steps, iterations), simply append the number of the element to the file name. Paraview will recognize the sequence automatically.

source