Reference
API
The following functions are the primary user-facing API of the FastMultipole package.
FastMultipole.fmm! — Function
fmm!(target_systems::Tuple, source_systems::Tuple; optargs...)Dispatches fmm! with automatic tree creation.
Arguments
target_systems::Union{Tuple, {UserDefinedSystem}}: either a system object for which compatibility functions have been overloaded, or a tuple of system objects for which compatibility functions have been overloadedsource_systems::Union{Tuple, {UserDefinedSystem}}: either a system object for which compatibility functions have been overloaded, or a tuple of system objects for which compatibility functions have been overloaded
Note: a convenience function fmm!(system) is provided, which is equivalent to fmm!(system, system). This is for situations where all system(s) act on all other systems, including themselves.
Optional Arguments: Allocation
target_buffers::Vector{<:Any}: buffers for target systems; if not provided, buffers are allocated usingFastMultipole.allocate_bufferstarget_small_buffers::Vector{<:Any}: small buffers for target systems; if not provided, small buffers are allocated usingFastMultipole.allocate_small_bufferssource_buffers::Vector{<:Any}: buffers for source systems; if not provided, buffers are allocated usingFastMultipole.allocate_bufferssource_small_buffers::Vector{<:Any}: small buffers for source systems; if not provided, small buffers are allocated usingFastMultipole.allocate_small_buffers
Optional Arguments: Tuning Parameters
expansion_order::Int: order of multipole expansions; default is 5multipole_acceptance::Float64: acceptance criterion for multipole expansions; default is 0.4leaf_size_target::Union{Nothing,Int}: leaf size for target systems; if not provided, the minimum of the source leaf sizes is usedleaf_size_source::Union{Nothing,Int}: leaf size for source systems; if not provided, the default leaf size is usederror_tolerance::Union{Nothing,ErrorMethod}: error tolerance for multipole to local translations; if not provided, no error treatment is performed
Optional Arguments: Tree Options
shrink::Bool: whether to shrink branches around their bodies, accounting for finite body radius; default istruerecenter::Bool: whether to recenter branches around their bodies, accounting for finite body radius; default isfalseinteraction_list_method::InteractionListMethod: method for building interaction lists; default isSelfTuningTargetStop()
Optional Arguments: Additional Options
farfield::Bool: whether to compute farfield interactions; default istruenearfield::Bool: whether to compute nearfield interactions; default istrueself_induced::Bool: whether to compute self-induced interactions; default istrueupward_pass::Bool: whether to perform the upward pass; default istruehorizontal_pass::Bool: whether to perform the horizontal pass; default istruedownward_pass::Bool: whether to perform the downward pass; default istruescalar_potential::Union{Bool,AbstractVector{Bool}}: whether to compute the scalar potential; default isfalsegradient::Union{Bool,AbstractVector{Bool}}: whether to compute the vector field; default istruehessian::Union{Bool,AbstractVector{Bool}}: whether to compute the vector gradient; default isfalseextra_farfield::Bool: whether to compute extra farfield interactions; default isfalse
FastMultipole.tune_fmm — Function
tune_fmm(target_systems, source_systems; optargs...)Tune the Fast Multipole Method (FMM) parameters for optimal performance on the given target and source systems, optionally subject to an error tolerance.
Arguments
target_systems::Union{Tuple,{UserDefinedSystem}}: a user-defined system object (or a tuple of them) for which the FMM interface functions have been definedsource_systems::Union{Tuple,{UserDefinedSystem}}: a user-defined system object (or a tuple of them) for which the FMM interface functions have been defined
Keyword Arguments
error_tolerance::Union{Nothing,Float64}: the error tolerance for the FMM; ifnothing, the FMM will simply use theexpansion_orderkeyword argument to fix the expansion orderexpansion_order::Int: the max expansion order for the FMM; defaults to 4leaf_size_source::Int: the leaf size for the source systems; defaults todefault_leaf_size(source_systems)max_expansion_order::Int: the maximum allowable expansion order if an error tolerance is requested; defaults to 20multipole_acceptances::AbstractRange{Float64}: a range of multipole acceptance critia to test; defaults torange(0.3, stop=0.8, step=0.1)lamb_helmholtz::Bool: whether to use the Lamb-Hellmholtz decomposition; defaults tofalseverbose::Bool: whether to print progress information; defaults totruekwargs...: additional keyword arguments to pass to thefmm!function
Returns
tuned_params::NamedTuple: a named tuple containing the best parameters found during tuning, which can be used in subsequentfmm!calls by splatting it as a keyword argument:`:leaf_size_source::Int: the optimal leaf size for the source systemsexpansion_order::Int: the optimal expansion order for the FMMmultipole_acceptance::Float64: the optimal multipole acceptance criterion
cache::Tuple: a tuple containing the cache used during tuning, which can be reused for subsequentfmm!calls by splatting it as a keyword argument
FastMultipole.direct! — Method
direct!(systems; derivatives_switches)Applies all interactions of systems acting on itself without multipole acceleration.
Arguments
systems: either- a system object for which compatibility functions have been overloaded, or
- a tuple of system objects for which compatibility functions have been overloaded
(target_systems, source_systems)wheretarget_systemsandsource_systemsare each either a system or tuple of systems for which compatibility functions have been overloaded
Optional Arguments
scalar_potential::Bool: either a::Boolor a::AbstractVector{Bool}of lengthlength(target_systems)indicating whether each system should receive a scalar potential fromsource_systemsgradient::Bool: either a::Boolor a::AbstractVector{Bool}of lengthlength(target_systems)indicating whether each system should receive a vector field fromsource_systemshessian::Bool: either a::Boolor a::AbstractVector{Bool}of lengthlength(target_systems)indicating whether each system should receive a vector gradient fromsource_systemsn_threads::Int: the number of threads to use for parallelization; defaults toThreads.nthreads()
Compatibility Functions
The following functions must be overloaded by the user to interface their code with the FastMultipole package.
FastMultipole.source_system_to_buffer! — Function
source_system_to_buffer!(buffer::Matrix, i_buffer, system::{UserDefinedSystem}, i_body)Compatibility function used to sort source systems. It should be overloaded for each system (where {UserDefinedSystem} is replaced with the type of the user-defined system) to be used as a source and should behave as follows. For the i_bodyth body contained inside of system,
buffer[1:3, i_buffer]should be set to the x, y, and z coordinates of the body position used for sorting into the octreebuffer[4, i_buffer]should be set to the radius beyond which a multipole expansion is allowed to be evaluated (e.g. for panels, or other bodies of finite area/volume)buffer[5:4+strength_dims, i_buffer]should be set to the body strength, which is a vector of lengthstrength_dims
Any additional information required for either forming multipole expansions or computing direct interactions should be stored in the rest of the column.
If a body contains vertices that are required for, e.g. computing multipole coefficients of dipole panels, these must be stored immediately following the body strength, and should be listed in a counter-clockwise order. For example, if I am using vortex tri-panels with strength_dims=3, I would set buffer[8:10,i] .= v1, buffer[11:13,i] .= v2, and bufer[14:16,i] .= v3, where v1, v2, and v3 are listed according to the right-hand-rule with thumb aligned with the panel normal vector.
Note that any system acting only as a target need not overload source_system_to_buffer!.
FastMultipole.data_per_body — Function
data_per_body(system::{UserDefinedSystem})Returns the number of values used to represent a single body in a source system. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
FastMultipole.get_position — Function
get_position(system::{UserDefinedSystem}, i)Returns a (static) vector of length 3 containing the x, y, and z coordinates of the position of the ith body. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
FastMultipole.get_previous_influence — Function
get_previous_influence(system::{UserDefinedSystem}, i)Returns the influence of the ith body in system from the previous FMM call. The relative error is predicted by dividing the absolute error by the result of this function. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
NOTE: If not overloaded, the default behavior is to return zero for both the scalar potential and vector field, effectively ignoring the relative tolerance in favor of an absolute tolerance.
Arguments:
system::{UserDefinedSystem}: the user-defined system objecti::Int: the index of the body within the system
Returns:
previous_potential::Float64: the previous scalar potential at theith bodyprevious_vector::SVector{3,Float64}: the previous vector field at theith body
FastMultipole.strength_dims — Function
strength_dims(system::{UserDefinedSystem})Returns the cardinality of the vector used to define the strength of each body inside system. E.g., a point mass would return 1, and a point dipole would return 3. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
FastMultipole.get_normal — Function
get_normal(source_buffer, source_system::{UserDefinedSystem}, i)OPTIONAL OVERLOAD:
Returns the unit normal vector for the ith body of source_buffer. May be (optionally) overloaded for a user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system); otherwise, the default behavior assumes counter-clockwise ordered vertices. Note that whatever method is used should match source_system_to_buffer! for each system.
FastMultipole.get_n_bodies — Function
get_n_bodies(system::{UserDefinedSystem})Returns the number of bodies contained inside system. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
FastMultipole.body_to_multipole! — Function
body_to_multipole!(system::{UserDefinedSystem}, multipole_coefficients, buffer, expansion_center, bodies_index, harmonics, expansion_order)Calculates the multipole coefficients due to the bodies contained in buffer[:,bodies_index] and accumulates them in multipole_coefficients. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
Typically, this is done using one of the convience functions contained within FastMultipole in one line, as
body_to_multipole!(system::MySystem, args...) = body_to_multipole!(Point{Vortex}, system, args...)FastMultipole.direct! — Method
direct!(target_buffer, target_index, derivatives_switch::DerivativesSwitch{PS,GS,HS}, ::{UserDefinedSystem}, source_buffer, source_index) where {PS,GS,HS}Calculates direct (nearfield) interactions of source_system on target_buffer. Should be overloaded or each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system), for all source bodies in source_index, at all target bodies in target_index, as follows:
# loop over source bodies
for i_source in source_index
# extract source body information here...
# loop over target bodies
for i_target in target_index
# get target position
target_position = get_position(target_buffer, i_target)
# evaluate influence here...
# update appropriate quantities
if PS
set_scalar_potential!(target_buffer, i_target, scalar_potential)
end
if GS
set_gradient!(target_buffer, i_target, gradient)
end
if HS
set_hessian!(target_buffer, i_target, hessian)
end
end
endNote that ::{UserDefinedSystem} is used purely for overloading the method for the appropriate system, and should NOT be accessed in this function, since it will NOT be indexed according to source_index. Rather, source_buffer, which is updated using source_system_to_buffer!, should be accessed.
The following convenience getter functions are available for accessing the source system:
get_position(source_system::{UserDefinedSystem}, i_body::Int): returns an SVector of length 3 containing the position of thei_bodybodyget_strength(source_buffer::Matrix, source_system::{UserDefinedSystem}, i_body::Int): returns an SVector containing the strength of thei_bodybodyget_vertex(source_buffer::Matrix, source_system::{UserDefinedSystem}, i_body::Int, i_vertex::Int): returns an SVector containing the x, y, and z coordinates of thei_vertexvertex of thei_bodybody
Note also that the compile time parameters PS, GS, and HS are used to determine whether the scalar potential and vector field should be computed, respectively. This allows us to skip unnecessary calculations and improve performance.
FastMultipole.buffer_to_target_system! — Function
buffer_to_target_system!(target_system::{UserDefinedSystem}, i_target, ::DerivativesSwitch{PS,GS,HS}, target_buffer, i_buffer) where {PS,GS,HS}Compatibility function used to update target systems. It should be overloaded for each system (where {UserDefinedSystem} is replaced with the type of the user-defined system) to be a target and should behave as follows. For the i_bodyth body contained inside of target_system,
target_buffer[4, i_buffer]contains the scalar potential influence to be added to thei_targetbody oftarget_systemtarget_buffer[5:7, i_buffer]contains the vector field influence to be added to thei_targetbody oftarget_systemtarget_buffer[8:16, i_buffer]contains the vector field gradient to be added to thei_targetbody oftarget_system
Note that any system acting only as a source (and not as a target) need not overload buffer_to_target_system!.
The following convenience functions can may be used to access the buffer:
get_scalar_potential(target_buffer, i_buffer::Int): returns the scalar potential induced at thei_bufferbody intarget_bufferget_gradient(target_buffer, i_buffer::Int): returns an SVector of length 3 containing the vector field induced at thei_bufferbody intarget_bufferget_hessian(target_buffer, i_buffer::Int): returns an SMatrix of size 3x3 containing the vector gradient induced at thei_bufferbody intarget_buffer
For some slight performance improvements, the booleans PS, GS, and HS can be used as a switch to indicate whether the scalar potential, vector field, and vector gradient are to be stored, respectively. Since they are compile-time parameters, if statements relying on them will not incur a runtime cost.
FastMultipole.has_vector_potential — Function
has_vector_potential(system::{UserDefinedSystem})Returns true if the system induces a vector potential, false otherwise. Should be overloaded for each user-defined system object (where {UserDefinedSystem} is replaced with the type of the user-defined system).
Data Structures
The following data structures are used by the FastMultipole package.
FastMultipole.DerivativesSwitch — Type
DerivativesSwitchSwitch indicating whether the scalar potential, vector potential, gradient, and/or hessian should be computed for a target system. Information is stored as type parameters, allowing the compiler to compile away if statements.
FastMultipole.Branch — Type
Branch{TF,N}Branch object used to sort more than one system into an octree. Type parameters represent:
TF: the floating point type (would be a dual number if using algorithmic differentiation)N: the number of systems represented
Fields
bodies_index::Vector{UnitRange}: vector of unit ranges indicating the index of bodies in each represented system, respectivelyn_branches::Int: number of child branches corresponding to this branchbranch_index::UnitRange: indices of this branch's child branchesi_parent::Int: index of this branch's parenti_leaf::Int: if this branch is a leaf, what is its index in its parent<:Tree'sleaf_indexfieldcenter::Vector{TF}: center of this branch at which its multipole and local expansions are centeredradius::TF: distance fromcenterto the farthest body contained in this branch (accounting for finite body radius if bodies are sources)box::Vector{TF}: vector of length 3 containing the distances from the center to faces of a rectangular prism completely enclosing all bodies in the x, y, and z direction, respectivelymin_potential::TF: maximum influence of any body in this branch on any body in its child branches; used to enforce a relative error tolerancemin_gradient::TF: maximum gradient magnitude of any body in this branch on any body in its child branches; used to enforce a relative error tolerance
FastMultipole.Tree — Type
Tree{TF,N}Tree object used to sort N systems into an octree.
Fields
branches::Vector{Branch{TF,N}}: a vector ofBranchobjects composing the treeexpansions::Array{TF,4}: 4-dimensional array whose(1,i,j,k)th element contains the real part of thejth expansion coefficient of thekth branch, and whose(2,i,j,k)th element contains the imaginary part. Ifi==1, the coefficient corresponds to the scalar potential; ifi==2, the coefficient corresponds to the ''\chi'' part of the Lamb-Helmholtz decomposition of the vector potential.levels_index::Vector{UnitRange{Int64}}: vector of unit ranges indicating the indices of branches at each level of the treeleaf_index::Vector{Int}: vector of indices of branches that are leavessort_index_list::NTuple{N,Vector{Int}}: tuple of vectors of indices used to sort the bodies in each system into the treeinverse_sort_index_list::NTuple{N,Vector{Int}}: tuple of vectors of indices used to undo the sort operation performed bysort_index_listbuffers::Vector{Matrix{TF}}: vector of buffers used to store the bodies computed influence of each system in the tree, as explained inFastMultipole.allocate_bufferssmall_buffers::Vector{Matrix{TF}}: vector of buffers used to pidgeon-hole sort bodies into the tree, as explained inFastMultipole.allocate_small_buffersexpansion_order::Int64: the maximum storable expansion orderleaf_size::SVector{N,Int64}: maximum number of bodies in a leaf for each system; if multiple systems are represented, the actual maximum depends on theInteractionListMethodused to create the tree
FastMultipole.ProbeSystem — Type
ProbeSystem{TF}Convenience system for defining locations at which the potential, vector field, or vector gradient may be desired. Interface functions are already defined and overloaded.
Fields
position::Vector{SVector{3,TF}}: vector of probe positionsscalar_potential::Vector{TF}: vector of scalar potential values at the positionsgradient::Vector{SVector{3,TF}}: vector of vector field values at the positionshessian::Vector{SMatrix{3,3,TF,9}}: vector of Hessian matrices at the positions
FastMultipole.Cache — Type
Cache{TF,NT,NS}Cache object used to store system buffers to avoid repeated allocations.
Fields
target_buffers::Vector{Matrix{TF}}: vector of buffers for target systemssource_buffers::Vector{Matrix{TF}}: vector of buffers for source systemstarget_small_buffers::Vector{Matrix{TF}}: vector of small buffers used for pidgeon-hole sorting target systems into the octreesource_small_buffers::Vector{Matrix{TF}}: vector of small buffers used for pidgeon-hole sorting source systems into the octree
Additional Functions
The following functions are used internally by the FastMultipole package, but may be useful to understand for advanced users.
FastMultipole.allocate_buffers — Function
allocate_buffers(systems::Tuple, target::Bool)Allocates buffers for the given systems.
Arguments
systems::Tuple: tuple of systems for which to allocate bufferstarget::Bool: iftrue, allocates space for body position, scalar potential, gradient, hessian matrix, and previous influence for relative error tolerance; otherwise, allocates space for information provided by the user-definedsource_system_to_buffer!TF: element type of the buffers
Returns
buffers::Vector{Matrix{TF}}: vector of matrices, one for each system, as follows:if
target==true, each matrix has size(18,N), whereNis the number of bodies in the system; indices correspond to:(1:3,i)- x, y, and z coordinates of theith body(4,i)- scalar potential induced at theith body(5:7,i)- x, y, and z components of the vector field induced at theith body(8:16,i)- hessian of the scalar potential induced at theith body(17,i)- estimate of the scalar potential induced at theith body used for relative error tolerance(18,i)- estimate of the vector field magnitude induced at theith body used for relative error tolerance
if
target==false, each matrix has size(M,N), whereNis the number of bodies in the system, andMis determined by the user-definedsource_system_to_buffer!function
FastMultipole.allocate_small_buffers — Function
allocate_small_buffers(systems::Tuple; target=false)Allocates small buffers for the given systems to be used in pidgeon-hole sorting bodies into the octree.
Arguments
systems::Tuple: tuple of systems for which to allocate buffersTF: element type of the buffers
Returns
small_buffers::Vector{Matrix{TF}}: vector of matrices of size(5,N)allocated for each system, where:Nis the number of bodies in the system- indices
(1:3,i)correspond to the x, y, and z coordinates of theith body - index
(4,i)corresponds to the estimated scalar potential of theith body - index
(5,i)corresponds to the estimated vector field magnitude of theith body