References
IFGF.IFGFOp — Typestruct IFGFOp{T} <: AbstractMatrix{T}High-level structure representing a linear operator with elements of type T with entry i,j given by K(X[i],Y[j]), where K is the underlying kernel function, and X and Y are the elements of the target_tree and source_tree, respectively.
Fields:
kernel: underlying kernel functiontarget_tree: tree structure used for clustering the target pointssource_tree: tree structure used for clustering the source pointsp::NTuple: number of chebyshev poins per dimension used in the interpolation of far interactionsds_func: function used to compute and appropriate mesh size in parameter space. See e.g.cone_domain_size_funcadm_func: function used to determine if the block with targets on boxAand sources on boxBadmit a low-rank approximation_source_tree_depth_partition: vector withi-th entry containing all nodes of the of depthi-1, where the root is given a depth of0.
IFGF.SourceTree — Typeconst SourceTree{N,T,V}Type alias for a ClusterTree of points (represented as SVector{N,T}) storing data of type SourceTreeData{N,T,V}. See the documentation of ClusterTree for more details.
IFGF.SourceTreeData — Typemutable struct SourceTreeData{N,T,V}Struct used to store data associated with a node of the SourceTree. The parametric type N represents the ambient dimension, T represents the primitive type used to represent the points (e.g. Float32 or Float64 for single/double precision), and V is the type of the density that will be interpolated.
Fields:
msh: uniform cartesian mesh of the interpolation space, indexable by aCartesianIndexconedatadict: dictionary mapping a cone tagI::CartesianIndexto anArraystoring its interpolation values/coefficient.farlist: nodes on the target tree that can be evaluated through interpolationnearlist: nodes on the target tree that must be evaluated by a direct sum.
IFGF.TargetTree — Typeconst TargetTree{N,T}Type alias for a ClusterTree of points (represented as SVector{N,T}) containing no additional data. See the documentation of ClusterTree in WavePropBase for more information.
IFGF._addnewcone! — Method_addnewcone!(node::SourceTree,I::CartesianIndex,p)If node does not contain cone I, add it to the conedict.
IFGF._density_type_from_kernel_type — Method_density_type_from_kernel_type(T)Helper function to compute the expected density V type associated with a kernel type T. For T<:Number, simply return T. For T<:SMatrix, return an SVector with elements of the type eltype(T) and of length size(T,2) so that multiplying an element of type T with an element of type V makes sense.
IFGF._init_msh! — Method_init_msh!(node,domain,ds)Create a uniform cartesian mesh of the interpolation domain of node with element spacing ds.
IFGF.allocate_interpolant_data! — Methodallocate_interpolant_data!(node::SourceTree,::Val{P},[I])Allocate the memory necessary to construct the I cone interpolant of node with P points, where P is a tuple of integers. If I is absent, allocate data for all cone interpolants of node.
IFGF.assemble_ifgf — Methodassemble_ifgf(K,X,Y;nmax=100,adm=modified_admissibility_condition,order=nothing,tol=1e-4,Δs=1.0,threads=true)Construct an approximation of the matrix [K(x,y) for x ∈ X, y ∈ Y]. If the kernel function K is oscillatory, you should implement the method wavenumber(::typeof(K)) --> Number to return its wavenumber.
The error in the approximation is controlled by two parameters: Δs and order. These parameters dictate the mesh size and interpolation order of the parametric space used to build the cone interpolants. While the order is held constant across all source cells, the Δs parameter controls only the size of the low-frequency cells; for cells of large acoustic size, the mesh size is adapted so as to keep the interpolation error constant across different levels of the source tree (see cone_domain_size_func for implementation details).
If the keyword order is ommited, a tolerance tol must be passed. It will then be used to estimate an appropriate interpolation order.
nmax can be used to control the maximum number of points contained in the target/source trees.
The adm function gives a criterion to determine if the interaction between a target cell and a source cell should be computed using an interpolation scheme.
IFGF.cart2interp — Methodcart2interp(x,c,h)
cart2interp(x,source::SourceTree)Map cartesian coordinates x into interpolation coordinates s centered at c and scaled by h. Interpolation coordinates are similar to hyperspherical coordinates (e.g. polar in two-dimensions and spherical in three-dimensions), except that the radial coordinates is replaced by s=h/r.
If passed a SourceTree, use its center and radius as c and h respectively.
IFGF.centered_factor — Methodcentered_factor(K,x,Y::SourceTree)Return a representative value of K(x,y) for y ∈ Y. This is used to help the interpolation of I(x) = ∑ᵢ K(x,yᵢ) where yᵢ ∈ Y and points x in the far field of Y. Defaults to K(x,yc) where yc is the center of Y.
IFGF.cheb2nodes — Methodcheb2nodes(n,[domain])
cheb2nodes(p::NTuple,[domain])
cheb2nodes(::Val{P},[domain])Return the n Chebyshev nodes of the second kind on the interval [-1,1].
If passed a tuple p, return the tensor product of one-dimensinal nodes of size p.
If passed a Val{P}, statically compute the tensor product nodes.
Finally, passing a domain::HyperRectangle shifts the interpolation nodes to it.
lb = (0.,0.)
ub = (3.,2.)
domain = HyperRectangle(lb,ub)
cheb2nodes((10,5),domain) # generate `10×5` Chebyshev nodes on the `[0,3]×[0,2]` rectangle.IFGF.cheb_error_estimate — Methodcheb_error_estimate(coefs::AbstractArray{T,N},dim)Given an N dimensional array of coefficients , estimate the relative error commited by the Chebyshev interpolant along dimension dim.
IFGF.compute_cone_list! — Methodcompute_cone_list!(ifgf::IFGFOperator)For each node in source_tree(ifgf), compute the cones which are necessary to cover both the target points on its farlist as well as the interpolation points on its parent's cones. The polynomial_order(ifgf) and meshsize_func(ifgf) are used to construct the appropriate interpolation domains and points.
IFGF.compute_interaction_list! — Functioncompute_interaction_list!(ifgf::IFGFOperator,adm)For each node in source_tree(ifgf), compute its nearlist and farlist. The nearlist contains all nodes in target_tree(ifgf) for which the interaction must be computed using the direct summation. The farlist includes all nodes in target_tree for which an accurate interpolant in available in the source node.
The adm function is called through adm(target,source) to determine if target is in the farlist of source. When that is not the case we recurse on their children unless both are a leaves, in which case target is added to the nearlist of source.
IFGF.compute_interpolation_domain — Methodcompute_interpolation_domain(source::SourceTree)Compute the minimal domain (as a HyperRectangle) in interpolation space for which source must provide a covering through its cone interpolants.
IFGF.cone_domain_size_func — Methodcone_domain_size_func(Δs₀::NTuple{N,T},k)Returns an anonymous function (node) -> Δs that computes an appropriate size Δs for the interpolation domain of node given an intial size Δs₀ and a wavenumber k. The function is constructed so as to scale Δs₀ by the inverse of the acoustic size of node.
IFGF.cone_index — Methodcone_index(x::SVector,tree::SourceTree)Given a point x in cartesian coordinates, return the index I of the cone in source to whcih x belongs, as well as its parametric coordinate s.
See also: cart2interp
IFGF.element_index — Methodelement_index(s::SVector,m::UniformCartesianMesh)Given a point s ∈ m, return the index I of the element in m containing s.
IFGF.free_interpolant_data! — Methodfree_interpolant_data!(node,[I])Free the memory used to store the I cone interpolant of node. If I is absent, free data associated to all cone interpolants of node
IFGF.initialize_source_tree — Methodinitialize_source_tree(;points,datatype,splitter,Xdatatype)Create the tree-structure for clustering the source points using the splitting strategy of splitter, returning a SourceTree built through the empty constructor (i.e. data = SourceTreeData()).
Arguments
points: vector of source points.splitter: splitting strategy.Xdatatype: type container of target points.Vdatatype: type of value stored at interpolation nodes
IFGF.initialize_target_tree — Methodinitialize_target_tree(;points,splitter)Create the tree-structure for clustering the target points using the splitting strategy of splitter. Returns a TargetTree.
IFGF.interp2cart — Methodinterp2cart(s,xc,h)
interp2cart(s,source::SourceTree)Map interpolation coordinate s into the corresponding cartesian coordinates.
Interpolation coordinates are similar to polar (in 2d) or spherical (in 3d) coordinates centered at the xc, except that the radial coordinate r is replaced by s = h/r.
If passed a SourceTree as second argument, call interp2cart(s,xc,h) with xc the center of the SourceTree and h its radius.
IFGF.inv_centered_factor — Methodinv_centered_factor(K,x,Y)Inverse of centered_factor.
IFGF.modified_admissibility_condition — Methodmodified_admissible_condition(target,source,[η])A target and source are admissible under the modified admissiblility condition (MAC) if the target box lies farther than r*η away, where r is the radius of the source box and η >= 1 is an adjustable parameter. By default, η = N / √N, where N is the ambient dimension.
IFGF.near_interaction! — Methodnear_interaction!(C,K,X,Y,σ,I,J)Compute C[i] <-- C[i] + ∑ⱼ K(X[i],Y[j])*σ[j] for i ∈ I, j ∈ J. The default implementation simply does a double loop over all i ∈ I and j ∈ J, but you override this method for the type of your own custom kernel K if you have a fast way of computing the full product in-place (e.g. using SIMD instructions).
IFGF.transfer_factor — Methodtransfer_factor(K,x,Y::SourceTree)Ratio between the centered factor of Y and the centered factor of its parent.
IFGF.wavenumber — Functionwavenumber(K::Function)For oscillatory kernels, return the characteristic wavenumber (i.e 2π divided by he wavelength). For non-oscillatory kernels, return 0.
LinearAlgebra.mul! — Methodmul!(y,A::IFGFOp,x,a,b;global_index=true,threads=true)Like the mul! function from LinearAlgebra. Setting global_index=false uses the internal indexing system of the IFGFOp, and therefore does not permute the vectors x and y to perform the multiplication.
IFGF.@hprofile — Macro@hprofileA macro which
- resets the default
TimerOutputs.get_defaulttimerto zero - execute the code block
- print the profiling details
This is useful as a coarse-grained profiling strategy to get a rough idea of where time is spent. Note that this relies on TimerOutputs annotations manually inserted in the code.
WavePropBase.Trees.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.Trees.AbstractSplitter — Typeabstract type AbstractSplitterAn AbstractSplitter is used to split a ClusterTree. The interface requires the following methods:
should_split(clt,splitter): return aBooldetermining if theClusterTreeshould be further dividedsplit!(clt,splitter): perform the splitting of theClusterTreehandling the necessary data sorting.
See GeometricSplitter for an example of an implementation.
WavePropBase.Trees.DyadicSplitter — Typestruct DyadicSplitter <: AbstractSplitterUsed to split an N dimensional ClusterTree into 2^N children until at most nmax points are contained in node.
See also: AbstractSplitter
WavePropBase.Trees.CardinalitySplitter — Typestruct CardinalitySplitter <: AbstractSplitterUsed 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.Trees.GeometricSplitter — Typestruct GeometricSplitter <: AbstractSplitterUsed to split a ClusterTree in half along the largest axis.
Missing docstring for mul!. Check Documenter's build log for details.