References
WavePropBase.ENTITIES
— Constantconst ENTITIES
Global dictionary storing the used entity tags (the value) for a given dimension (the key).
WavePropBase.SType
— Typeconst 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.
WavePropBase.TAGS
— Constantconst TAGS::Dict{Int,Vector{Int}}
Global dictionary storing the used entity tags (the value) for a given dimension (the key).
WavePropBase.etype_to_vtk_cell_type
— Constantconst etype_to_vtk_cell_type
Dictionary 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)
WavePropBase.AbstractElement
— Typeabstract 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) coordinatex̂ ∈ D
.jacobian(el,x̂)
: evaluate the jacobian matrix of the interpolation at the (reference) coordinatex ∈ D
.
For performance reasons, both el(x̂)
and jacobian(el,x̂)
should take as input a StaticVector
and output a static vector or static array.
WavePropBase.AbstractEntity
— Typeabstract type AbstractEntity
Entity of geometrical nature. Identifiable throught its (dim,tag)
key.
WavePropBase.AbstractHyperRectangle
— Typeabstract type AbstractHyperRectangle{N,T} <: AbstractElement{ReferenceHyperCube{N},SVector{N,T}}
Axis-aligned hyperrectangle in N
dimensions with coordinates of type SVector{N,T}
.
WavePropBase.AbstractKernel
— Typeabstract type AbstractKernel{T}
A kernel functions K
with the signature K(target,source)::T
.
See also: GenericKernel
, SingleLayerKernel
, DoubleLayerKernel
, AdjointDoubleLayerKernel
, HyperSingularKernel
WavePropBase.AbstractMesh
— Typeabstract type AbstractMesh{N,T}
An abstract mesh structure in dimension N
with primite data of type T
(e.g. Float64
for double precision representation).
Concrete subtypes of AbstractMesh
should implement ElementIterator
for accessing the mesh elements, and/or NodeIterator
for accesing the mesh nodes.
See also: GenericMesh
, UniformCartesianMesh
WavePropBase.AbstractPDE
— Typeabstract type AbstractPDE{N}
A partial differential equation in dimension N
. AbstractPDE
types are used to define AbstractPDEKernel
s.
WavePropBase.AbstractPDEKernel
— Typeabstract type AbstractPDEKernel{T,Op} <: AbstractKernel{T}
An AbstractKernel
with an associated pde::Op
field.
WavePropBase.AbstractPolynomialSpace
— Typeabstract type AbstractPolynomialSpace{D}
A polynomial space over D
. This is a vector space under polynomial addition and scalar multiplication.
WavePropBase.AbstractQuadratureRule
— Typeabstract type AbstractQuadratureRule{D}
A quadrature rule for integrating a function over the domain D
.
Calling q()
returns the nodes x
and weights w
for performing integration over domain(q)
.
Calling q(el)
returns the nodes and weights for integrating over el
.
WavePropBase.AbstractReferenceShape
— Typeabstract type AbstractReferenceShape{N}
A reference domain/shape in ℜᴺ
.
Used mostly for defining more complex shapes as transformations mapping an AbstractReferenceShape
into some region of ℜᴹ
.
See e.g. ReferenceLine
or ReferenceTriangle
for some examples of concrete subtypes.
WavePropBase.AbstractSingularityHandler
— Typeabstract type AbstractSingularityHandler{R}
Used for handling localized integrand singularities in R
.
WavePropBase.AbstractSplitter
— Typeabstract type AbstractSplitter
An AbstractSplitter
is used to split a ClusterTree
. The interface requires the following methods:
should_split(clt,splitter)
: return aBool
determining if theClusterTree
should be further dividedsplit!(clt,splitter)
: perform the splitting of theClusterTree
handling the necessary data sorting.
See GeometricSplitter
for an example of an implementation.
WavePropBase.AbstractTree
— Typeabstract type AbstracTree end
Supertype for tree-like objects.
WavePropBase.AdjointDoubleLayerKernel
— Typestruct AdjointDoubleLayerKernel{T,Op} <: AbstractPDEKernel{T,Op}
Given an operator Op
, construct its free-space adjoint double-layer kernel. This corresponds to the transpose(γ₁,ₓ[G])
, where G
is the SingleLayerKernel
. For operators such as Laplace
or Helmholtz
, this is simply the normal derivative of the fundamental solution respect to the target variable.
WavePropBase.CardinalitySplitter
— Typestruct 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.
See also: AbstractSplitter
WavePropBase.ClusterTree
— Typemutable struct ClusterTree{T,S,D}
Tree structure used to cluster elements of type V = eltype(T)
into containers of type S
. The method center(::V)::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 typeD
.
WavePropBase.ClusterTree
— MethodClusterTree(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.
WavePropBase.ComplexPoint2D
— Typeconst ComplexPoint2D
const ComplexPoint2D(x1, x2)
const ComplexPoint2D(x::NTuple{2, ComplexF64})
A complex 2D point, stored in a StaticArray. ComplexPoint2D = SVector{2, ComplexF64}.
WavePropBase.ComplexPoint3D
— Typeconst ComplexPoint3D
const ComplexPoint3D(x1, x2, x3)
const ComplexPoint3D(x::NTuple{3, ComplexF64})
A complex 3D point, stored in a StaticArray. ComplexPoint3D = SVector{3, ComplexF64}.
WavePropBase.CustomQuadratureRule
— Typestruct CustomQuadratureRule{D,N,T} <: AbstractQuadratureRule{D}
N
-point user-defined quadrature rule for integrating over D
.
WavePropBase.Domain
— Typestruct Domain
Represent a physical domain as a union of entities.
See also: AbstractEntity
, ElementaryEntity
.
WavePropBase.DoubleLayerKernel
— Typestruct DoubleLayerKernel{T,Op} <: AbstractPDEKernel{T,Op}
Given an operator Op
, construct its free-space double-layer kernel. This corresponds to the γ₁
trace of the SingleLayerKernel
. For operators such as Laplace
or Helmholtz
, this is simply the normal derivative of the fundamental solution respect to the source variable.
WavePropBase.Duffy
— Typestruct 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.
WavePropBase.DyadicMaxDepthSplitter
— Typestruct DyadicMaxDepthSplitter <: AbstractSplitter
Similar to DyadicSplitter
, but splits nodes until a maximum depth
is reached.
See also: AbstractSplitter
WavePropBase.DyadicMinimalSplitter
— Typestruct DyadicMinimalSplitter <: AbstractSplitter
Similar to DyadicSplitter
, but the boundin boxes are shrank to the minimal axis-aligned boxes at the end.
See also: AbstractSplitter
WavePropBase.DyadicSplitter
— Typestruct DyadicSplitter <: AbstractSplitter
Used to split an N
dimensional ClusterTree
into 2^N
children until at most nmax
points are contained in node.
See also: AbstractSplitter
WavePropBase.Elastostatic
— Typestruct Elastostatic{N,T} <: AbstractPDE{N}
Elastostatic equation in N
dimensions: μΔu + (μ+λ)∇(∇⋅u) = 0. Note that the displacement u is a vector of length N
since this is a vectorial problem.
WavePropBase.ElementIterator
— Typestruct 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.
WavePropBase.ElementaryEntity
— Typestruct ElementaryEntity <: AbstractEntity
The most basic representation of an AbstractEntity
.
Fields:
dim::UInt8
: the geometrical dimension of the entity (e.g. line hasdim=1
, surface hasdim=2
, etc)tag::Int64
: an integer tag associated to the entityboundary::Vector{AbstractEntity}
: the entities of dimensiondim-1
forming the entity's boundary
WavePropBase.ElementaryEntity
— MethodElementaryEntity(dim,tag)
Construct an ElementaryEntity
with an empty boundary .
WavePropBase.Fejer
— Typestruct Fejer{N}
N
-point Fejer's first quadrature rule for integrating a function over [0,1]
. Exactly integrates all polynomials of degree ≤ N-1
.
import WavePropBase as WPB
q = WPB.Fejer(;order=10)
WPB.integrate(cos,q) ≈ sin(1) - sin(0)
# output
true
WavePropBase.Gauss
— Typestruct Gauss{D,N} <: AbstractQuadratureRule{D}
Tabulated N
-point symmetric Gauss quadrature rule for integration over D
.
WavePropBase.GaussLegendre
— Typestruct GaussLegendre{N}
N
-point Gauss-Legendre quadrature rule for integrating a function over [0,1]
. Exactly integrates all polynomials of degree ≤ 2N-1
.
WavePropBase.GaussLegendre
— MethodGaussLegendre(;order)
Construct a GaussLegendre
of the desired order over the [0,1]
interval.
WavePropBase.GenericKernel
— Typestruct GenericKernel{T,F} <: AbstractKernel{T}
An AbstractKernel
with kernel
of type F
.
WavePropBase.GenericMesh
— Typestruct GenericMesh{N,T} <: AbstractMesh{N,T}
Data structure representing a generic mesh in an ambient space of dimension N
, with data of type T
.
WavePropBase.GeometricMinimalSplitter
— Typestruct GeometricMinimalSplitter <: AbstractSplitter
Like GeometricSplitter
, but shrinks the children's containters.
WavePropBase.GeometricSplitter
— Typestruct GeometricSplitter <: AbstractSplitter
Used to split a ClusterTree
in half along the largest axis.
WavePropBase.Helmholtz
— Typestruct Helmholtz{N,T}
Helmholtz equation in N
dimensions: Δu + k²u = 0.
WavePropBase.HyperCube
— Typestruct HyperCube{N,T}
Axis-aligned hypercube in N
dimensions given by low_corner::SVector{N,T}
and side::T
.
WavePropBase.HyperRectangle
— Typestruct HyperRectangle{N,T}
Axis-aligned hyperrectangle in N
dimensions given by low_corner::SVector{N,T}
and high_corner::SVector{N,T}
.
WavePropBase.HyperSingularKernel
— Typestruct HyperSingularKernel{T,Op} <: AbstractPDEKernel{T,Op}
Given an operator Op
, construct its free-space hypersingular kernel. This corresponds to the transpose(γ₁,ₓγ₁[G])
, where G
is the SingleLayerKernel
. For operators such as Laplace
or Helmholtz
, this is simply the normal derivative of the fundamental solution respect to the target variable of the DoubleLayerKernel
.
WavePropBase.IMT
— Typestruct IMT{A,P} <: AbstractSingularityHandler{ReferenceLine}
One-dimensional change of variables mapping [0,1] -> [0,1]
with the property that all derivatives vanish at the point x=0
.
See Davis and Rabinowitz.
WavePropBase.IntegralOperator
— Typestruct IntegralOperator{T} <: AbstractMatrix{T}
A discrete linear integral operator given by
\[I[u](x) = \int_{\Gamma\_s} K(x,y)u(y) ds_y, x \in \Gamma_{t}\]
where $\Gamma_s$ and $\Gamma_t$ are the source and target domains, respectively.
WavePropBase.IntegralPotential
— Typestruct IntegralPotential
Represent a potential given by a kernel
and a source_mesh
over which integration is performed.
IntegralPotential
s are created using IntegralPotential(kernel, source_mesh)
.
Evaluating an integral potential requires a density σ
(defined over the quadrature nodes of the source mesh) and a point x
at which to evaluate the integral
\[\int_{\Gamma} K(oldsymbol{x},oldsymbol{y})\sigma(y) ds_y, x \not \in \Gamma\]
WavePropBase.Kress
— Typestruct 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
.
WavePropBase.KressP
— Typestruct 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 a
1-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
.
WavePropBase.LagrangeElement
— Typestruct LagrangeElement{D,Np,T} <: AbstractElement{D,T}
Standard element over D <: AbstractReferenceShape
commonly used in finite element methods. The underlying polynomial space is PolynomialSpace{D,K}
, and its interpolant maps the Np
reference_nodes
in D
to Np
values of type T
stored in the field vals
.
WavePropBase.LagrangeLine
— Typeconst LagrangeLine = LagrangeElement{ReferenceLine}
WavePropBase.LagrangePoint
— Typeconst LagrangePoint{N,T} = LagrangeElement{ReferencePoint,1,SVector{N,T}}
WavePropBase.LagrangeSquare
— Typeconst LagrangeSquare = LagrangeElement{ReferenceSquare}
WavePropBase.LagrangeTetrahedron
— Typeconst LagrangeTetrahedron = LagrangeElement{ReferenceTetrahedron}
WavePropBase.LagrangeTriangle
— Typeconst LagrangeTriangle = LagrangeElement{ReferenceTriangle}
WavePropBase.Laplace
— Typestruct Laplace{N}
Laplace equation in N
dimension: Δu = 0.
WavePropBase.NodeIterator
— Typestruct NodeIterator{M}
Iterator for all the nodes in a mesh of type M
.
WavePropBase.NystromDensity
— Typestruct NystromDensity{T} <: AbstractVector{T}
Discrete density with values vals
on the qnodes
of a NystromMesh
WavePropBase.NystromDensity
— MethodNystromDensity(f::Function, X::NystromMesh)
Return a NystromDensity
with values f(q)
at the qnodes
of X
.
Note that the argument passsed to f
is a QuadratureNode
, so that f
may depend on quantities other than the coords
of the quadrature node (such as the normal
vector).
See also: QuadratureNode
WavePropBase.NystromMesh
— Typestruct NystromMesh{N,T} <: AbstractMesh{N,T}
A mesh data structure for solving boundary integral equation using Nyström methods.
A NystromMesh
can be constructed from a mesh::AbstractMesh
and a dictionary etype2qrule
mapping element types in mesh
(the keys) to an appropriate quadrature rule for integration over elements of that type (the value).
The degrees of freedom in a NystromMesh
are associated to nodal values at the quadrature nodes, and are represented using QuadratureNode
s.
WavePropBase.NystromMesh
— MethodNystromMesh(Ω;meshsize,qorder)
Create a NystromMesh
with elements of size meshsize
and quadrature order qorder
.
A mesh is first generated using meshgen
, and then a NystromMesh
is created on top of it with the quadrature information.
WavePropBase.NystromMesh
— MethodNystromMesh(msh::AbstractMesh,e2qrule::Dict)
NystromMesh(msh::AbstractMesh;qorder)
Construct a NystromMesh
with the quadrature q = e2qrule[E]
applied to each element type E
in msh. If an order
keyword is passed, a default quadrature of the desired order is used for each element type.
WavePropBase.ParametricElement
— TypeParametricElement{F,D,T} <: AbstractElement{D,T}
An element represented through a explicit function f
mapping D
into the element.
See also: AbstractElement
WavePropBase.ParametricEntity
— TypeParametricEntity <: AbstractEntity
A geometric entity given by an explicit parametrization
field. A ParametricEntity
ent
should be a callable through ent(x)
where x ∈ domain
.
See also: AbstractEntity
WavePropBase.PlotPoints
— Typestruct PlotPoints
Structure used for dispatching SVector
to plot recipes without type-piracy.
WavePropBase.PlotTree
— Typestruct PlotTree
Used to plot entire tree associated with a tree node, instead of just the node.
WavePropBase.Point1D
— Typeconst Point1D
const Point1D(x1)
A point in 1D space, stored in a StaticArray. Point1D = SVector{1, Float64}.
WavePropBase.Point2D
— Typeconst Point2D
const Point2D(x1, x2)
const Point2D(x::NTuple{2, Float64})
A point in 2D space, stored in a StaticArray. Point2D = SVector{2, Float64}.
WavePropBase.Point3D
— Typeconst 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}.
WavePropBase.PointEntity
— TypePointEntity{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.
WavePropBase.PolynomialSpace
— Typestruct PolynomialSpace{D,K} <: AbstractPolynomialSpace{D}
The space of all polynomials of degree ≤K
, commonly referred to as ℙₖ
.
The type parameter D
, of singleton type, is used to determine the reference domain of the polynomial basis. In particular, when D
is a hypercube in d
dimensions, the precise definition is ℙₖ = span{𝐱ᶿ : 0≤max(θ)≤ K}
; when D
is a d
-dimensional simplex, the space is ℙₖ = span{𝐱ᶿ : 0≤sum(θ)≤ K}
, where θ ∈ 𝐍ᵈ
is a multi-index.
See also: monomial_basis
, lagrange_basis
WavePropBase.PrincipalComponentSplitter
— Typestruct PrincipalComponentSplitter <: AbstractSplitter
WavePropBase.PseudoBlockMatrix
— Typestruct PseudoBlockMatrix{T<:SMatrix,S} <: AbstractMatrix{T}
A struct which behaves identically to a Matrix{T}
, but with the underlying data
stored as a Matrix{S}
, where S::Number = eltype(T)
. This allows for the use of blas
routines under-the-hood, while providing a convenient interface for handling matrices over tensors.
WavePropBase.QuadratureNode
— TypeQuadratureNode{N,T<:Real}
A point in ℝᴺ
with a weight
for performing numerical integration.
A QuadratureNode
can optionally store a normal
vector and a curvature
scalar, given by the divergence of the normal vector.
WavePropBase.ReferenceCube
— Typeconst ReferenceCube = ReferenceHyperCube{3}
Singleton type representing the unit cube [0,1]³
.
WavePropBase.ReferenceHyperCube
— Typestruct ReferenceHyperCube{N}
Singleton type representing the axis-aligned hypercube in N
dimensions with the lower corner at the origin and the upper-corner at (1,1,…,1)
.
WavePropBase.ReferenceLine
— Typeconst ReferenceLine = ReferenceHyperCube{1}
Singleton type representing the [0,1]
segment.
WavePropBase.ReferencePoint
— Typestruct ReferencePoint
Singleton type representing a reference zero-dimensional entitty (i.e. a point).
WavePropBase.ReferenceSquare
— Typeconst ReferenceSquare = ReferenceHyperCube{2}
Singleton type representing the square with vertices (0,0),(0,1),(1,1),(1,0)
WavePropBase.ReferenceTetrahedron
— Typestruct ReferenceTetrahedron
Singleton type representing the tetrahedron with vertices (0,0,0),(0,0,1),(0,1,0),(1,0,0)
WavePropBase.ReferenceTriangle
— Typestruct ReferenceTriangle
Singleton type representing the triangle with vertices (0,0),(1,0),(0,0)
WavePropBase.SingleLayerKernel
— Typestruct SingleLayerKernel{T,Op} <: AbstractPDEKernel{T,Op}
The free-space single-layer kernel (i.e. the fundamental solution) of an OP <: AbstractPDE
.
WavePropBase.SubMesh
— Typestruct 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.
WavePropBase.TensorLagInterp
— Typestruct 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:
using WavePropBase: TensorLagInterp
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(vals,(x,y))
p((0.1,0.2)) ≈ f((0.1,0.2))
# output
true
WavePropBase.TensorProductQuadrature
— TypeTensorProductQuadrature{N,Q}
A tensor-product of one-dimension quadrature rules. Integrates over [0,1]^N
.
Examples
qx = Fejer(10)
qy = TrapezoidalOpen(15)
q = TensorProductQuadrature(qx,qy)
WavePropBase.TensorProductSingularityHandler
— Typestruct TensorProductSingularityHandler{S} <: AbstractSingularityHandler{ReferenceSquare}
A tensor product of two one-dimensional AbstractSingularityHandler
s for performing integration over the ReferenceSquare
.
WavePropBase.Trapezoidal
— Typestruct 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)
WavePropBase.TrapezoidalOpen
— Typestruct TrapezoidalOpen{N} <: AbstractQuadratureRule{ReferenceLine}
Open trapezoidal rule.
WavePropBase.TriangularMesh
— Typestruct TriangularMesh{N,T} <: AbstractMesh{N,T}
A simple mesh structure containing only (flat) triangles in N
dimensions. The underlying data type is given by T
, which defaults to Float64
.
TriangularMesh
es are presented as a collection of nodes
and a connectivity
matrix of size 3 × Nt
, where Nt
is the number of triangles in the mesh.
WavePropBase.UniformCartesianMesh
— Typestruct 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.
WavePropBase.UniformCartesianMesh
— MethodUniformCartesianMesh(domain::HyperRectangle,sz::NTuple)
UniformCartesianMesh(domain::HyperRectangle;step::NTuple)
Construct a uniform UniformCartesianMesh
with sz[d]
elements along dimension d
. If the kwarg step
is passed, construct a UniformCartesianMesh
with elements of approximate size step
.
WavePropBase.Yukawa
— Typestruct Yukawa{N,T}
Yukawa equation in N
dimensions: Δu - λ²u = 0.
This is the same as the Helmholtz equation with k = iλ
, but more efficient functions are used to compute the kernels of the integral operators when λ
is real (the default).
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 ElementaryEntity
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.
Base.:==
— Method===(Ω1::Domain,Ω2::Domain)
Two Domain
s are equal if all their entities are equal (regardless of order).
Base.in
— Methodin(ω::ElementaryEntity,Ω::Domain)
Check whether an ElementaryEntity
belongs to a Domain
by recursively checking whether it belongs to its boundary.
Base.iterate
— Functioniterate(Ω::Domain)
Iterating over a domain means iterating over its entities.
Base.keys
— MethodReturn all tags of the elementary entities in the domain Ω
corresponding to the dimension d
.
Base.keys
— MethodReturn all tags of the elementary entities in the domain Ω
corresponding to the dimensions contained in dims
.
Base.length
— Methodlength(Ω:::Domain)
The length of a domain corresponds to the number of elementary entities that make it.
Base.parent
— Methodparent(t::AbstractTree)
The node's parent. If t
is a root, then parent(t)==t
.
Base.split
— Methodsplit(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.
WavePropBase._basis_vdim
— Method_basis_vdim(ℒ::AbstractPDE,ℙₖ::PolynomialSpace)
Return a set of polynomial pₙ
such that ℒ[pₙ]
gives the monomial basis for ℙₖ
.
WavePropBase._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 container
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.
WavePropBase._get_gauss_qcoords_and_qweights
— Method_get_gauss_and_qweights(R::Type{<:AbstractReferenceShape{D}}, N) where D
Returns the N
-point symmetric gaussian qnodes and qweights (x, w)
for integration over R
.
WavePropBase._meshgen
— Method_meshgen(f,d,sz)
Create a UniformCartesianMeshof
dpush-forward map. The cartesian mesh has size
sz`, and is uniform in parameter space.
WavePropBase._normal
— Method_normal(jac::SMatrix{M,N})
Given a an M
by N
matrix representing the jacobian of a codimension one object, compute the normal vector.
WavePropBase._vtk_cells
— Function_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.
WavePropBase._vtk_points
— Method_vtk_points(mesh::GenericMesh)
Creates the matrix of nodes in the format required by WriteVTK
.
WavePropBase.abstractmethod
— Methodabstractmethod
A method of an abstract type
for which concrete subtypes are expected to provide an implementation.
WavePropBase.ambient_dimension
— Methodambient_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.
WavePropBase.assert_extension
— Functionassert_extension(fname,ext,[msg])
Check that fname
is of extension ext
. Print the message msg
as an assertion error otherwise.
WavePropBase.assertequaldim
— Methodassertequaldim(Ω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.
WavePropBase.blockmatrix_to_matrix
— Methodblockmatrix_to_matrix(A::Matrix{B}) where {B<:SMatrix}
Convert a Matrix{B}
, where B<:SMatrix
, to the equivalent Matrix{T}
, where T = eltype(B)
WavePropBase.blockvector_to_vector
— Methodblockvector_to_vector(A::Vector{B}) where {B<:SVector}
Convert a Vector{B}
, where B<:SVector
, to the equivalent Vector{T}
, where T = eltype(B)
WavePropBase.boundary
— Methodboundary(Ω::Domain)
Return a domain comprising the external boundary of Ω.
See also: external_boundary
WavePropBase.cart2pol
— Methodcart2pol(x,y)
Map cartesian coordinates x,y
to polar coordinates r,θ
. The convention followed is that -π ≤ θ ≤ π
.
WavePropBase.cart2sph
— Methodcart2sph(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 -π < φ ≤ π
.
WavePropBase.center
— Methodcenter(Ω)
Center of the smallest ball containing Ω
.
WavePropBase.cheb1nodes
— Methodcheb1nodes(n,a,b)
Return the n
Chebyshev points of the first kind on the interval [a,b]
.
WavePropBase.cheb2nodes
— Methodcheb2nodes(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]
.
WavePropBase.children
— Methodchildren(t::AbstractTree)
Iterable collection of the node's children.
WavePropBase.clear_entities!
— Methodclear_entities!()
Empty the global variables used to keep track of the various entities created.
WavePropBase.container
— Methodcontainer(clt::ClusterTree)
Return the object enclosing all the elements of the clt
.
WavePropBase.container_type
— Methodcontainer_type(clt::ClusterTree)
Type used to enclose the elements of clt
.
WavePropBase.coords
— Methodcoords(x)
Return an SVector
with the cartesian coordinates associated to a geometrical object x
.
WavePropBase.cross_product_matrix
— Methodcross_product_matrix(v)
Returns the matrix Aᵥ
associated with the cross product v × ϕ
so that v × ϕ = Aᵥϕ
.
WavePropBase.curvature
— Methodcurvature(el, x̂)
The curvature of el
at the parametric coordinate x̂
, given by the following formula
\[\kappa = \nabla\cdot\frac{\nabla p}{|\nabla p|} = \frac{\Delta p}{|\nabla p|} + \frac{\nabla^\perp p\nabla^2 p\nabla p}{|\nabla p|^3}\]
WavePropBase.decompose
— Methoddecompose(s::AbstractReferenceShape,x,[target_shape])
Decompose an AbstractReferenceShape
into LagrangeElement
s so that x
is a vertex of the children elements.
Examples
s = ReferenceLine() el1, el2 = decompose(s,0.3) el1(1) == el2(0) == 0.3 # true
WavePropBase.decrement_index
— Functiondecrement_index(I::CartesianIndex,k[,n=1])
Equivalent to increment_index
(I,k,-n)
WavePropBase.degree
— Methoddegree(el::LagrangeElement)
The polynomial degree of the element. A LagrangeElement
of degree K
and domain D
belongs to the space PolynomialSpace{D,K}
.
WavePropBase.depth
— Functiondepth(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).
WavePropBase.diameter
— Methoddiameter(Ω)
Largest distance between x
and y
for x,y ∈ Ω
.
WavePropBase.dim_correction
— Methoddim_correction(pde,X,Y,S,D[;location,p,derivative,tol])
Given a pde
and a (possibly innacurate) discretizations of its single and double-layer operators S
and D
with domain Y
and range X
, compute corrections δS
and δD
such that S + δS
and D + δD
are more accurate approximations of the underlying integral operators.
The following optional keyword arguments are available:
p::DimParameters
: parameters associated with the density interpolation methodderivative
: if true, compute the correction to the adjoint double-layer and hypersingular operators instead. In this case,S
andD
should contain a discretization of adjoint double-layer and hypersingular operators, respectively.tol
: distance beyond which interactions are considered sufficiently far so that no correction is needed.
WavePropBase.dimension
— Methoddimension(space)
The length of a basis for space
; i.e. the number of linearly independent elements required to span space
.
WavePropBase.distance
— Methoddistance(X::ClusterTree, Y::ClusterTree)
Distance between the containers of X
and Y
.
WavePropBase.distance
— Methoddistance(x::SVector,r::HyperRectangle)
The (minimal) Euclidean distance between the point x
and any point y ∈ r
.
WavePropBase.distance
— Methoddistance(Ω₁,Ω₂)
Minimal Euclidean distance between a point x ∈ Ω₁
and y ∈ Ω₂
.
WavePropBase.dom2elt
— Methoddom2elt(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 Ω
.
WavePropBase.dom2elt
— Methoddom2elt(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.
WavePropBase.dom2elt
— Methoddom2elt(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.
WavePropBase.domain
— Methoddomain(f)
Given a function-like object f: Ω → R
, return Ω
.
WavePropBase.element_index_for_point
— Methodelement_index_for_point(p::SVector,m::UniformCartesianMesh)
Given a point p
, return the index I
of the element in m
containing p
.
WavePropBase.element_type
— Methodelement_type(clt::ClusterTree)
Type of elements sorted in clt
.
WavePropBase.elements
— Methodelements(clt::ClusterTree)
Iterable list of the elements inside clt
.
WavePropBase.enable_debug
— Functionenable_debug(mod)
Activate debugging messages in module mod
.
WavePropBase.entities
— Methodentities(Ω::Domain)
Return a vector of all elementary entities making up a domain.
WavePropBase.external_boundary
— MethodReturn the external boundaries inside a domain.
WavePropBase.fibonnaci_points_sphere
— Methodfibonnaci_points_sphere(N,r,c)
Return N
points distributed (roughly) in a uniform manner on the sphere of radius r
centered at c
.
WavePropBase.filter_tree
— Functionfilter_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
.
WavePropBase.filter_tree!
— Functionfilter_tree!(filter,nodes,tree,[isterminal=true])
Like filter_tree
, but appends results to nodes
.
WavePropBase.flip_normal
— Methodflip_normal(Γ::Domain)
Reverse the orientation of the normal vector in the entities of a domain.
WavePropBase.geometric_dimension
— Methodgeometric_dimension(x::AbstractEntity)
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.
WavePropBase.global_add_entity!
— Methodglobal_add_entity!(ent::AbstractEntity)
Add ent
to the global dictionary ENTITIES
and update TAGS
with its (dim,tag)
key. This function should be called by the inner constructor of every AbstractEntity
; see the constructor of ElementaryEntity
for an example.
WavePropBase.image
— Methodimage(f)
Given a function-like object f: Ω → R
, return R
.
WavePropBase.image
— Methodimage(f)
The a function-like object f: Ω → R
, return R
.
WavePropBase.increment_index
— Functionincrement_index(I::CartesianIndex,k[,n=1])
Increment I
by n
along the dimension k
. This is equivalent to I += n*eₖ
, where eₖ
is a vector with with 1
at the k
-th coordinate and zeros elsewhere.
WavePropBase.index_range
— Methodindex_range(clt::ClusterTree)
Indices of elements in root_elements(clt)
which lie inside clt
.
WavePropBase.integrate
— Methodintegrate(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.
WavePropBase.integrate
— Methodintegrate(f,msh::NystromMesh)
Compute ∑ᵢ f(qᵢ)wᵢ
, where the qᵢ
are the qnodes
s of msh
, and wᵢ
are the qweights
.
Note that you must define f(::QuadratureNode)
: use coords(q)
and normal(q)
if you need to access the coordinate or normal vector at que quadrature node.
WavePropBase.internal_boundary
— MethodReturn the internal boundaries inside a domain.
WavePropBase.isleaf
— Methodisleaf(t)
Return true
if t
has no children.
WavePropBase.isroot
— Methodisroot(t)
Return true
if t
is its own parent.
WavePropBase.jacobian
— Methodjacobian(f,x)
Given a (possibly vector-valued) functor f : 𝐑ᵐ → 𝐅ⁿ
, return the n × m
matrix Aᵢⱼ = ∂fᵢ/∂xⱼ
. By default a finite-difference approximation is performed, but you should overload this method for specific f
if better performance and/or precision is required.
Note: both x
and f(x)
are expected to be of SVector
type.
WavePropBase.key
— Methodkey(e::AbstractEntity)
The (dim,tag)
pair used as a key to identify geometrical entities.
WavePropBase.lagrange_basis
— Methodlagrange_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.
It is assumed that the value of a function on nodes
uniquely determine a polynomial in sp
.
WavePropBase.line
— Methodline(a,b)
Create a straight line connecting points a
and b
. Returns an instance of ParametricEntity
.
WavePropBase.loc2glob
— Methodloc2glob(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.
WavePropBase.matrix_to_blockmatrix
— Methodmatrix_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.
WavePropBase.meshgen!
— Methodmeshgen!(mesh,Ω,sz)
Similar to meshgen
, but append entries to mesh
.
WavePropBase.meshgen
— Methodmeshgen(Ω::Domain,num_elements)
meshgen(Ω::Domain;meshsize)
Generate a GenericMesh
for the domain Ω
with num_elements
per entity. To specify a different number of elements per entity, num_elements
should be a vector with as many elements as there are entities in Ω
. Alternatively, a meshsize
can be passed.
Requires the entities forming Ω
to be ParametricEntity
.
The quality of the generated mesh created usign meshgen
depends on the quality of the underlying parametrization. For complex surfaces, you are better off using a proper mesher such as gmsh
.
WavePropBase.monomial_basis
— Functionmonomial_basis(sp::PolynomialSpace)
Return an NTuple
containing a basis of monomials 𝐱ᶿ
spanning the polynomial space PolynomialSpace
.
WavePropBase.near_interaction_list
— Methodnear_interaction_list(X,Y::AbstractMesh; tol)
For each element el
of type E
in Y
, return the indices of the points in X
which are closer than tol
to the center
of el
.
This function returns a dictionary where e.g. dict[E][5] --> Vector{Int}
gives the indices of points in X
which are closer than tol
to the center of the fifth element of type E
.
WavePropBase.nearest_point_to_element
— Methodnearest_point_to_element(X::NystroMesh,Y::NystromMesh; tol)
For each element el
, return a list with the indices of all points in X
for which el
is the nearest element. Ignore indices for which the distance exceeds tol
.
WavePropBase.new_tag
— Methodnew_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
.
WavePropBase.normal
— Methodnormal(el,x̂)
The unit normal vector at coordinate x̂
, guaranteed to be orthogonal to all columns of jacobian(el,x)
.
WavePropBase.notimplemented
— Methodnotimplemented()
Things which should probably be implemented at some point.
WavePropBase.partition_by_depth
— Methodpartition_by_depth(tree)
Given a tree
, return a partition
vector whose i
-th entry stores all the nodes in tree
with depth=i-1
. Empty nodes are not added to the partition.
WavePropBase.pde
— Methodpde(K::AbstractPDEKernel)
Return the underlying AbstractPDE
when K
correspond to the kernel related to the underlying Greens function of a PDE.
WavePropBase.pol2cart
— Methodpol2cart(r,θ)
Map polar coordinates r,θ
to cartesian coordinates x,y
.
WavePropBase.print_threads_info
— Methodprint_threads_info()
Prints in console the total number of threads.
WavePropBase.qcoords
— Methodqcoords(Y)
Return the coordinate of the quadrature nodes associated with Y
.
WavePropBase.qnodes
— Methodqnodes(Q::NystromMesh)
Return a vector of QuadratureNode
s associated to the mesh Q
.
WavePropBase.qrule_for_reference_shape
— Methodqrule_for_reference_shape(ref,order)
Given a ref
erence shape and a desired quadrature order
, return an appropiate quadrature rule.
WavePropBase.quadgen
— Methodquadgen(el,q::AbstractQuadratureRule)
Return a set of points and weights for integrating over el
based on the quadrature rule q
.
WavePropBase.qweights
— Methodqweights(Y)
Return the quadrature weights associated with Y
.
WavePropBase.radius
— Methodradius(Ω)
Half the diameter
.
WavePropBase.reference_nodes
— Methodreference_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).
WavePropBase.return_type
— Methodreturn_type(f[,args...])
The type returned by f(args...)
, where args
is a tuple of types. Falls back to Base.promote_op
by default.
A functors of type T
with a knonw return type should extend return_type(::T,args...)
to avoid relying on promote_op
.
WavePropBase.root_elements
— Methodroot_elements(clt::ClusterTree)
The elements contained in the root of the tree to which clt
belongs.
WavePropBase.section
— Methodsection(rec::HyperRectangle{D},ax)
Return the dimension D-1
HyperRectangle
obtained by deleting the coordinates at the ax
dimension.
WavePropBase.should_split
— Methodshould_split(clt::ClusterTree, depth, splitter::AbstractSplitter)
Determine whether or not a ClusterTree
should be further divided.
WavePropBase.skeleton
— Methodskeleton(Ω::Domain)
Return all the boundaries of the domain, i.e. the domain's skeleton.
WavePropBase.sort_by_type
— Methodsort_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)
WavePropBase.sort_in_cartesian_mesh
— Methodsort_in_cartesian_mesh(X,msh)
Given a collection of points X
and a mesh msh
, return a dictionary mapping the keys of the mesh elements, given as a CartesianIndex
, to the indices of points inside that element.
WavePropBase.sph2cart
— Methodsph2cart(x,y,z)
Map spherical coordinates r,θ,φ
representing the radius, elevation, and azimuthal angle respectively, to cartesian ones x, y, z
.
WavePropBase.split!
— Methodsplit!(clt::ClusterTree,splitter::AbstractSplitter)
Divide clt
using the strategy implemented by splitter
. This function is reponsible of assigning the children
and parent
fields, as well as of permuting the data of clt
.
WavePropBase.svector
— Methodsvector(f,n)
Just like Base.ntuple
, but convert output to an SVector
.
WavePropBase.tag
— Methodtag(::AbstractEntity)
Integer tag used to idetify geometrical entities.
WavePropBase.uniform_points_circle
— Methoduniform_points_circle(N,r,c)
Return N
points uniformly distributed on a circle of radius r
centered at c
.
WavePropBase.vector_to_blockvector
— Methodvector_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.
WavePropBase.vtk_mesh_file
— Methodvtk_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.
WavePropBase.@interface
— Macro@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: @interface
@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.