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 overloaded
  • source_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

Optional Arguments: Tuning Parameters

  • expansion_order::Int: order of multipole expansions; default is 5
  • multipole_acceptance::Float64: acceptance criterion for multipole expansions; default is 0.4
  • leaf_size_target::Union{Nothing,Int}: leaf size for target systems; if not provided, the minimum of the source leaf sizes is used
  • leaf_size_source::Union{Nothing,Int}: leaf size for source systems; if not provided, the default leaf size is used
  • error_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 is true
  • recenter::Bool: whether to recenter branches around their bodies, accounting for finite body radius; default is false
  • interaction_list_method::InteractionListMethod: method for building interaction lists; default is SelfTuningTargetStop()

Optional Arguments: Additional Options

  • farfield::Bool: whether to compute farfield interactions; default is true
  • nearfield::Bool: whether to compute nearfield interactions; default is true
  • self_induced::Bool: whether to compute self-induced interactions; default is true
  • upward_pass::Bool: whether to perform the upward pass; default is true
  • horizontal_pass::Bool: whether to perform the horizontal pass; default is true
  • downward_pass::Bool: whether to perform the downward pass; default is true
  • scalar_potential::Union{Bool,AbstractVector{Bool}}: whether to compute the scalar potential; default is false
  • gradient::Union{Bool,AbstractVector{Bool}}: whether to compute the vector field; default is true
  • hessian::Union{Bool,AbstractVector{Bool}}: whether to compute the vector gradient; default is false
  • extra_farfield::Bool: whether to compute extra farfield interactions; default is false
source
FastMultipole.tune_fmmFunction
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 defined
  • source_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; if nothing, the FMM will simply use the expansion_order keyword argument to fix the expansion order
  • expansion_order::Int: the max expansion order for the FMM; defaults to 4
  • leaf_size_source::Int: the leaf size for the source systems; defaults to default_leaf_size(source_systems)
  • max_expansion_order::Int: the maximum allowable expansion order if an error tolerance is requested; defaults to 20
  • multipole_acceptances::AbstractRange{Float64}: a range of multipole acceptance critia to test; defaults to range(0.3, stop=0.8, step=0.1)
  • lamb_helmholtz::Bool: whether to use the Lamb-Hellmholtz decomposition; defaults to false
  • verbose::Bool: whether to print progress information; defaults to true
  • kwargs...: additional keyword arguments to pass to the fmm! function

Returns

  • tuned_params::NamedTuple: a named tuple containing the best parameters found during tuning, which can be used in subsequent fmm! calls by splatting it as a keyword argument:`:

    • leaf_size_source::Int: the optimal leaf size for the source systems
    • expansion_order::Int: the optimal expansion order for the FMM
    • multipole_acceptance::Float64: the optimal multipole acceptance criterion
  • cache::Tuple: a tuple containing the cache used during tuning, which can be reused for subsequent fmm! calls by splatting it as a keyword argument

source
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) where target_systems and source_systems are each either a system or tuple of systems for which compatibility functions have been overloaded

Optional Arguments

  • scalar_potential::Bool: either a ::Bool or a ::AbstractVector{Bool} of length length(target_systems) indicating whether each system should receive a scalar potential from source_systems
  • gradient::Bool: either a ::Bool or a ::AbstractVector{Bool} of length length(target_systems) indicating whether each system should receive a vector field from source_systems
  • hessian::Bool: either a ::Bool or a ::AbstractVector{Bool} of length length(target_systems) indicating whether each system should receive a vector gradient from source_systems
  • n_threads::Int: the number of threads to use for parallelization; defaults to Threads.nthreads()
source

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 octree
  • buffer[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 length strength_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!.

source
FastMultipole.data_per_bodyFunction
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).

source
FastMultipole.get_positionFunction
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).

source
FastMultipole.get_previous_influenceFunction
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 object
  • i::Int: the index of the body within the system

Returns:

  • previous_potential::Float64: the previous scalar potential at the ith body
  • previous_vector::SVector{3,Float64}: the previous vector field at the ith body
source
FastMultipole.strength_dimsFunction
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).

source
FastMultipole.get_normalFunction
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.

source
FastMultipole.get_n_bodiesFunction
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).

source
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...)
source
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
end

Note 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 the i_body body
  • get_strength(source_buffer::Matrix, source_system::{UserDefinedSystem}, i_body::Int): returns an SVector containing the strength of the i_body body
  • get_vertex(source_buffer::Matrix, source_system::{UserDefinedSystem}, i_body::Int, i_vertex::Int): returns an SVector containing the x, y, and z coordinates of the i_vertex vertex of the i_body body

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.

source
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 the i_target body of target_system
  • target_buffer[5:7, i_buffer] contains the vector field influence to be added to the i_target body of target_system
  • target_buffer[8:16, i_buffer] contains the vector field gradient to be added to the i_target body of target_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 the i_buffer body in target_buffer
  • get_gradient(target_buffer, i_buffer::Int): returns an SVector of length 3 containing the vector field induced at the i_buffer body in target_buffer
  • get_hessian(target_buffer, i_buffer::Int): returns an SMatrix of size 3x3 containing the vector gradient induced at the i_buffer body in target_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.

source
FastMultipole.has_vector_potentialFunction
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).

source

Data Structures

The following data structures are used by the FastMultipole package.

FastMultipole.DerivativesSwitchType
DerivativesSwitch

Switch 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.

source
FastMultipole.BranchType
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, respectively
  • n_branches::Int: number of child branches corresponding to this branch
  • branch_index::UnitRange: indices of this branch's child branches
  • i_parent::Int: index of this branch's parent
  • i_leaf::Int: if this branch is a leaf, what is its index in its parent <:Tree's leaf_index field
  • center::Vector{TF}: center of this branch at which its multipole and local expansions are centered
  • radius::TF: distance from center to 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, respectively
  • min_potential::TF: maximum influence of any body in this branch on any body in its child branches; used to enforce a relative error tolerance
  • min_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
source
FastMultipole.TreeType
Tree{TF,N}

Tree object used to sort N systems into an octree.

Fields

  • branches::Vector{Branch{TF,N}}: a vector of Branch objects composing the tree
  • expansions::Array{TF,4}: 4-dimensional array whose (1,i,j,k)th element contains the real part of the jth expansion coefficient of the kth branch, and whose (2,i,j,k)th element contains the imaginary part. If i==1, the coefficient corresponds to the scalar potential; if i==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 tree
  • leaf_index::Vector{Int}: vector of indices of branches that are leaves
  • sort_index_list::NTuple{N,Vector{Int}}: tuple of vectors of indices used to sort the bodies in each system into the tree
  • inverse_sort_index_list::NTuple{N,Vector{Int}}: tuple of vectors of indices used to undo the sort operation performed by sort_index_list
  • buffers::Vector{Matrix{TF}}: vector of buffers used to store the bodies computed influence of each system in the tree, as explained in FastMultipole.allocate_buffers
  • small_buffers::Vector{Matrix{TF}}: vector of buffers used to pidgeon-hole sort bodies into the tree, as explained in FastMultipole.allocate_small_buffers
  • expansion_order::Int64: the maximum storable expansion order
  • leaf_size::SVector{N,Int64}: maximum number of bodies in a leaf for each system; if multiple systems are represented, the actual maximum depends on the InteractionListMethod used to create the tree
source
FastMultipole.ProbeSystemType
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 positions
  • scalar_potential::Vector{TF}: vector of scalar potential values at the positions
  • gradient::Vector{SVector{3,TF}}: vector of vector field values at the positions
  • hessian::Vector{SMatrix{3,3,TF,9}}: vector of Hessian matrices at the positions
source
FastMultipole.CacheType
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 systems
  • source_buffers::Vector{Matrix{TF}}: vector of buffers for source systems
  • target_small_buffers::Vector{Matrix{TF}}: vector of small buffers used for pidgeon-hole sorting target systems into the octree
  • source_small_buffers::Vector{Matrix{TF}}: vector of small buffers used for pidgeon-hole sorting source systems into the octree
source

Additional Functions

The following functions are used internally by the FastMultipole package, but may be useful to understand for advanced users.

FastMultipole.allocate_buffersFunction
allocate_buffers(systems::Tuple, target::Bool)

Allocates buffers for the given systems.

Arguments

  • systems::Tuple: tuple of systems for which to allocate buffers
  • target::Bool: if true, 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-defined source_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), where N is the number of bodies in the system; indices correspond to:

      • (1:3,i) - x, y, and z coordinates of the ith body
      • (4,i) - scalar potential induced at the ith body
      • (5:7,i) - x, y, and z components of the vector field induced at the ith body
      • (8:16,i) - hessian of the scalar potential induced at the ith body
      • (17,i) - estimate of the scalar potential induced at the ith body used for relative error tolerance
      • (18,i) - estimate of the vector field magnitude induced at the ith body used for relative error tolerance
    • if target==false, each matrix has size (M,N), where N is the number of bodies in the system, and M is determined by the user-defined source_system_to_buffer! function

source
FastMultipole.allocate_small_buffersFunction
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 buffers
  • TF: element type of the buffers

Returns

  • small_buffers::Vector{Matrix{TF}}: vector of matrices of size (5,N) allocated for each system, where:

    • N is the number of bodies in the system
    • indices (1:3,i) correspond to the x, y, and z coordinates of the ith body
    • index (4,i) corresponds to the estimated scalar potential of the ith body
    • index (5,i) corresponds to the estimated vector field magnitude of the ith body
source