Reference
FastMultipole.Branch
FastMultipole.DerivativesSwitch
FastMultipole.DerivativesSwitch
FastMultipole.DerivativesSwitch
FastMultipole.DerivativesSwitch
FastMultipole.ExpansionSwitch
FastMultipole.ProbeSystem
FastMultipole.Tree
FastMultipole.EmptyTree
FastMultipole.TreeByLevel
FastMultipole._rotate_multipole_y_n!
FastMultipole.allocate_buffers
FastMultipole.allocate_small_buffers
FastMultipole.back_rotate_local_y!
FastMultipole.back_rotate_multipole_y!
FastMultipole.back_rotate_z!
FastMultipole.body_to_multipole!
FastMultipole.buffer_to_system_strength!
FastMultipole.buffer_to_target_system!
FastMultipole.calculate_ib!
FastMultipole.calculate_pj!
FastMultipole.closest_corner
FastMultipole.d2rdx2
FastMultipole.data_per_body
FastMultipole.direct!
FastMultipole.direct!
FastMultipole.drdx
FastMultipole.fmm!
FastMultipole.get_n_bodies
FastMultipole.get_normal
FastMultipole.get_position
FastMultipole.index_by_source
FastMultipole.influence!
FastMultipole.influence!
FastMultipole.local_error
FastMultipole.local_power
FastMultipole.local_to_local!
FastMultipole.make_direct_assignments!
FastMultipole.map_by_branch
FastMultipole.map_by_leaf
FastMultipole.mirrored_source_to_vortex!
FastMultipole.multipole_error
FastMultipole.multipole_to_local!
FastMultipole.multipole_to_local!
FastMultipole.multipole_to_local!
FastMultipole.multipole_to_local!
FastMultipole.multipole_to_local!
FastMultipole.nearfield_device!
FastMultipole.nearfield_device!
FastMultipole.predict_error
FastMultipole.rotate_multipole_y!
FastMultipole.rotate_z!
FastMultipole.rotate_z_n!
FastMultipole.rotate_z_n_power!
FastMultipole.self_influence_matrices
FastMultipole.set_gradient!
FastMultipole.set_hessian!
FastMultipole.set_scalar_potential!
FastMultipole.shrink_branch!
FastMultipole.source_system_to_buffer!
FastMultipole.source_to_dipole!
FastMultipole.source_to_vortex_point!
FastMultipole.strength_dims
FastMultipole.strength_to_value
FastMultipole.target_influence_to_buffer!
FastMultipole.translate_local_z!
FastMultipole.translate_multipole_to_local_z!
FastMultipole.translate_multipole_to_local_z_m01_n
FastMultipole.translate_multipole_to_local_z_n!
FastMultipole.translate_multipole_z!
FastMultipole.tune_fmm
FastMultipole.update_Ts_0!
FastMultipole.update_eimϕs_0!
FastMultipole.update_eimϕs_n!
FastMultipole.update_nonself_influence!
FastMultipole.value_to_strength!
FastMultipole.Branch
— TypeBranch{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_index
fieldcenter::Vector{TF}
: center of this branch at which its multipole and local expansions are centeredsource_radius::TF
: if this branch is a leaf, distance fromcenter
to the outer edge of farthest body contained in this branch; otherwise, this is the distance fromcenter
to the corner of itssource_box
target_radius::TF
: distance fromcenter
to the farthest body center contained in this branchsource_box::Vector{TF}
: vector of length 6 containing the distances from the center to faces of a rectangular prism completely enclosing all bodies with their finite radius in the negative x, positive x, negative y, positive y, negative z, and positive z directions, respectivelytarget_box::Vector{TF}
: vector of length 3 containing the distances from the center to faces of a rectangular prism completely enclosing all body centers in the x, y, and z direction, respectivelymax_influence::TF
: maximum influence of any body in this branch on any body in its child branches; used to enforce a relative error tolerance
FastMultipole.DerivativesSwitch
— TypeDerivativesSwitch
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.
FastMultipole.DerivativesSwitch
— MethodDerivativesSwitch(scalar_potential, gradient, hessian)
Constructs a tuple of DerivativesSwitch
objects.
Arguments
scalar_potential::Vector{Bool}
: a vector of::Bool
indicating whether the scalar potential should be computed for each target systemgradient::Vector{Bool}
: a vector of::Bool
indicating whether the vector field should be computed for each target systemhessian::Vector{Bool}
: a vector of::Bool
indicating whether the vector gradient should be computed for each target system
FastMultipole.DerivativesSwitch
— MethodDerivativesSwitch(scalar_potential, gradient, hessian, target_systems)
Constructs a ::Tuple
of indentical DerivativesSwitch
objects of the same length as target_systems
(if it is a ::Tuple
), or a single DerivativesSwitch
(if target_system
is not a ::Tuple
)
Arguments
scalar_potential::Bool
: a::Bool
indicating whether the scalar potential should be computed for each target systemgradient::Bool
: a::Bool
indicating whether the vector field should be computed for each target systemhessian::Bool
: a::Bool
indicating whether the vector gradient should be computed for each target system
FastMultipole.DerivativesSwitch
— MethodDerivativesSwitch(scalar_potential, gradient, hessian)
Constructs a single DerivativesSwitch
object.
Arguments
scalar_potential::Bool
: a::Bool
indicating whether the scalar potential should be computed for the target systemgradient::Bool
: a::Bool
indicating whether the vector field should be computed for the target systemhessian::Bool
: a::Bool
indicating whether the vector gradient should be computed for the target system
FastMultipole.ExpansionSwitch
— TypeExpansionSwitch
Switch indicating which expansions should be used:
- scalar potential (
SP
) - vector potential via Lamb-Helmholtz decomposition (
VP
)
FastMultipole.ProbeSystem
— TypeProbeSystem
Convenience system for defining locations at which the potential, vector field, or vector gradient may be desired.
FastMultipole.Tree
— Typebodies[indexlist] is the same sort operation as performed by the tree sortedbodies[inverseindexlist] undoes the sort operation performed by the tree
FastMultipole.EmptyTree
— MethodEmptyTree(system)
Returns an empty tree. Used if system
is empty.
Arguments
system
: the system from which a tree is to be created
Returns
tree
: iftypeof(system)<:Tuple
, a::MultiTree
is returned; otherwise, a::SingleTree
is returned
FastMultipole.TreeByLevel
— FunctionDoesn't stop subdividing until ALL child branches have satisfied the leaf size.
FastMultipole._rotate_multipole_y_n!
— MethodThe first time this function is called, it should receive iζ=0, iT=0, and n=0. It will return iζ and iT appropriate for n+1.
FastMultipole.allocate_buffers
— Methodallocate_buffers(systems::Tuple, target::Bool)
Allocates buffers for the given systems. If target
is true
, it allocates space for position, scalar potential, gradient, and hessian matrix. Otherwise, it allocates enough memory for the user-defined source_system_to_buffer!
.
FastMultipole.allocate_small_buffers
— Methodallocate_small_buffers(systems::Tuple; target=false)
Allocates small buffers for the given systems. These buffers are used for temporary storage of body positions for octree sorting.
FastMultipole.back_rotate_local_y!
— MethodAssumes Ts, Hsπ2, and ηsmag have all been precomputed. Resets target_weights.
FastMultipole.back_rotate_multipole_y!
— MethodAssumes Ts, Hsπ2, and ζsmag have all been precomputed. Resets target_weights.
FastMultipole.back_rotate_z!
— MethodAssumes eimϕs have already been computed. DOES NOT overwrite rotated weights (unlike other rotate functions); rather, accumulates on top of it.
FastMultipole.body_to_multipole!
— Methodbody_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.buffer_to_system_strength!
— Methodbuffer_to_system_strength!(system::{UserDefinedSystem}, source_buffer::Matrix{Float64}, i_body::Int)
NOTE: this function is primarily used for the boundary element solver, and is not required for the FMM.
Updates the strength in system
for the i_body
th body using the strength information contained in source_buffer
. Should be overloaded for each user-defined system object used with the boundary element solver.
Arguments:
system::{UserDefinedSystem}
: the user-defined system objecti_body::Int
: the index of the body insource_system
whose strength is to be setsource_buffer::Matrix{Float64}
: the source buffer containing the body informationi_buffer::Int
: the index of the body insource_buffer
whose strength is to be set
FastMultipole.buffer_to_target_system!
— Methodbuffer_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_body
th body contained inside of target_system
,
target_buffer[4, i_buffer]
contains the scalar potential influence to be added to thei_target
body oftarget_system
target_buffer[5:7, i_buffer]
contains the vector field influence to be added to thei_target
body oftarget_system
target_buffer[8:16, i_buffer]
contains the vector field gradient to be added to thei_target
body 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_buffer
body intarget_buffer
get_gradient(target_buffer, i_buffer::Int)
: returns an SVector of length 3 containing the vector field induced at thei_buffer
body intarget_buffer
get_hessian(target_buffer, i_buffer::Int)
: returns an SMatrix of size 3x3 containing the vector gradient induced at thei_buffer
body 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.calculate_ib!
— Methodassumes j has already been calculated
Note: if X0real is replaced with X0real + Xw_real, etc., this becomes bnm (as used for volumes) instead of inm (as used for panels)
FastMultipole.calculate_pj!
— Methodassumes q has already been calculated
If X0real is replaced with X0real + Xv_real, etc., the result becomes jnm instead of pnm, as used for panels.
FastMultipole.closest_corner
— MethodVector from center 2 to the corner of box2 closest to center.
FastMultipole.d2rdx2
— Methoddrk/dxidx_j
FastMultipole.data_per_body
— Methoddata_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.direct!
— Methoddirect!(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 thei_body
bodyget_strength(source_buffer::Matrix, source_system::{UserDefinedSystem}, i_body::Int)
: returns an SVector containing the strength of thei_body
bodyget_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_vertex
vertex of thei_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.
FastMultipole.direct!
— Methoddirect!(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
Optional Arguments
scalar_potential::Bool
: either a::Bool
or a::AbstractVector{Bool}
of lengthlength(target_systems)
indicating whether each system should receive a scalar potential fromsource_systems
gradient::Bool
: either a::Bool
or a::AbstractVector{Bool}
of lengthlength(target_systems)
indicating whether each system should receive a vector field fromsource_systems
hessian::Bool
: either a::Bool
or a::AbstractVector{Bool}
of lengthlength(target_systems)
indicating whether each system should receive a vector gradient fromsource_systems
FastMultipole.drdx
— Methoddrj/dxi
FastMultipole.fmm!
— Functionfmm!(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 usingallocate_buffers
target_small_buffers::Vector{<:Any}
: small buffers for target systems; if not provided, small buffers are allocated usingallocate_small_buffers
source_buffers::Vector{<:Any}
: buffers for source systems; if not provided, buffers are allocated usingallocate_buffers
source_small_buffers::Vector{<:Any}
: small buffers for source systems; if not provided, small buffers are allocated usingallocate_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_recenter::Bool
: whether to shrink and recenter branches around their bodies, accounting for finite body radius; default istrue
interaction_list_method::InteractionListMethod
: method for building interaction lists; default isSelfTuningTreeStop()
Optional Arguments: Additional Options
farfield::Bool
: whether to compute farfield interactions; default istrue
nearfield::Bool
: whether to compute nearfield interactions; default istrue
self_induced::Bool
: whether to compute self-induced interactions; default istrue
upward_pass::Bool
: whether to perform the upward pass; default istrue
horizontal_pass::Bool
: whether to perform the horizontal pass; default istrue
downward_pass::Bool
: whether to perform the downward pass; default istrue
scalar_potential::Union{Bool,AbstractVector{Bool}}
: whether to compute the scalar potential; default isfalse
gradient::Union{Bool,AbstractVector{Bool}}
: whether to compute the vector field; default istrue
hessian::Union{Bool,AbstractVector{Bool}}
: whether to compute the vector gradient; default isfalse
FastMultipole.get_n_bodies
— Methodget_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.get_normal
— Methodget_normal(source_buffer, source_system::{UserDefinedSystem}, i)
OPTIONAL OVERLOAD:
Returns the unit normal vector for the i
th 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_position
— Methodget_position(system::{UserDefinedSystem}, i)
Returns a (static) vector of length 3 containing the x, y, and z coordinates of the position of the i
th body. Should be overloaded for each user-defined system object (where {UserDefinedSystem}
is replaced with the type of the user-defined system).
FastMultipole.index_by_source
— Methodindex_by_source(sorted_list::Vector{SVector{2,Int}}, leaf_index::Vector{Int})
Constructs an index map that maps each leaf to the range of indices in the sorted list where it acts as a source. It is possible for some leaves to never act as sources, in which case the range will be empty.
FastMultipole.influence!
— Methodinfluence!(influence, target_buffer, source_system, source_buffer)
NOTE: source_system
is provided solely for dispatch; it's member bodies will be out of order and should not be referenced.
NOTE: This function is primarily used for the boundary element solver, and is not required for the FMM.
Evaluate the influence as pertains to the boundary element influence matrix and overwrites it to influence
(which would need to be subtracted for it to act like the RHS of a linear system). Based on the current state of the target_buffer
and source_buffer
. Should be overloaded for each system type that is used in the boundary element solver.
Arguments:
influence::AbstractVector{TF}
: vector containing the influence for every body in the target buffertarget_buffer::Matrix{TF}
: target buffer used to compute the influencesource_system::{UserDefinedSystem}
: system object used solely for dispatchsource_buffer::Matrix{TF}
: source buffer used to compute the influence
FastMultipole.influence!
— Methodinfluence!(sorted_influences, influences_per_system, target_buffers, source_systems, source_buffers, source_tree)
Evaluate the influence as pertains to the boundary element influence matrix and subtracts it from sorted_influences
(which would act like the RHS of a linear system). Based on the current state of the target_buffers
and source_buffers
. Note that source_systems
is provided solely for dispatch. Note also that influences_per_system
is overwritten each time.
sorted_influences::Vector{Float64}
: single vector containing the influence for every body in the target buffers, sorted by source branch in the direct interaction listinfluences_per_system::Vector{Vector{Float64}}
: vector of vectors containing the influence for each target system, sorted the same way as the bufferstarget_buffers::NTuple{N,Matrix{Float64}}
: target buffers used to compute the influencesource_systems::NTuple{N,<:{UserDefinedSystem}}
: system objects used for dispatchsource_buffers::NTuple{N,Matrix{Float64}}
: source buffers used to compute the influence
FastMultipole.local_error
— MethodMultipole coefficients up to expansion order `P+1` must be added to `multipole_branch` before calling this function.
FastMultipole.local_power
— Methodperforms the lamb-helmholtz transformation assuming that the coordinate system is still aligned with the z axis
FastMultipole.local_to_local!
— MethodExpects ηsmag and Hsπ2 to be precomputed. Ts and eimϕs are computed here.
FastMultipole.make_direct_assignments!
— Methodmake_direct_assignments!(assignments, i_target_system, target_branches, i_source_system, source_branches, direct_list, n_threads, n_per_thread, interaction_list_method)
Assumes direct_list
is sorted by target branch index. Assigns ranges of interactions to each thread.
FastMultipole.map_by_branch
— Methodmap_by_branch(target_tree::Tree)
Constructs a mapping from each branch to the range of indices in the targets vector that correspond to the bodies in that branch.
This relies on the fact that tree.leaf_index
is sorted in order of increasing branch index.
FastMultipole.map_by_leaf
— Methodmap_by_leaf(source_tree::Tree)
Constructs a mapping from each leaf to the range of indices in the strengths vector that correspond to the bodies in that leaf.
FastMultipole.mirrored_source_to_vortex!
— Methodmy derivation
FastMultipole.multipole_error
— MethodMultipole coefficients up to expansion order `P+1` must be added to `multipole_branch` before calling this function.
FastMultipole.multipole_to_local!
— Methoddefaults to no error prediction
FastMultipole.multipole_to_local!
— MethodExpects ζsmag, ηsmag, and Hs_π2 to be computed a priori.
FastMultipole.multipole_to_local!
— MethodExpects ζsmag, ηsmag, and Hs_π2 to be computed a priori.
FastMultipole.multipole_to_local!
— MethodExpects ζsmag, ηsmag, and Hs_π2 to be computed a priori.
FastMultipole.multipole_to_local!
— MethodExpects ζsmag, ηsmag, and Hs_π2 to be computed a priori.
FastMultipole.nearfield_device!
— Methodnearfield_device!(target_systems, derivatives_switches, source_systems)
Dispatches nearfield_device!
without having to build a ::Tree
. Performs all interactions.
Arguments
target_systems
: user-defined system on whichsource_system
actsderivatives_switches::Union{DerivativesSwitch, NTuple{N,DerivativesSwitch}}
: determines whether the scalar potential, vector field, and or vector gradient should be calculatedsource_systems
: user-defined system acting ontarget_system
FastMultipole.nearfield_device!
— Methodnearfield_device!(target_systems, target_tree, derivatives_switches, source_systems, source_tree, direct_list)
User-defined function used to offload nearfield calculations to a device, such as GPU.
Arguments
target_systems
: user-defined system on whichsource_system
actstarget_tree::Tree
: octree object used to sorttarget_systems
derivatives_switches::Union{DerivativesSwitch, NTuple{N,DerivativesSwitch}}
: determines whether the scalar potential, vector field, and or vector gradient should be calculatedsource_systems
: user-defined system acting ontarget_system
source_tree::Tree
: octree object used to sorttarget_systems
direct_list::Vector{SVector{2,Int32}}
: each element[i,j]
maps nearfield interaction fromsource_tree.branches[j]
ontarget_tree.branches[i]
FastMultipole.predict_error
— Methodperforms a partial M2L for expansion_order+1, and predicts the error along the way;
expects that the upward pass has been performed up to expansion_order+1
FastMultipole.rotate_multipole_y!
— MethodRotate solid harmonic weights about the y axis by θ. Note that Hsπ2 and ζsmag must be updated a priori, but Ts is updated en situ. Resets rotated_weights before computing.
FastMultipole.rotate_z!
— MethodPerforms a z-axis rotation of the supplied solid harmonic coefficients. Computes e^{imϕ} as well.
FastMultipole.rotate_z_n!
— MethodPerforms a z-axis rotation of the supplied solid harmonic coefficients. Computes e^{imϕ} as well.
FastMultipole.rotate_z_n_power!
— MethodPerforms a z-axis rotation of the supplied solid harmonic coefficients. Computes e^{imϕ} and the multipole power at n as well.
FastMultipole.self_influence_matrices
— Methodself_influence_matrices(target_buffers, source_buffers, source_systems, target_tree, source_tree, derivatives_switches)
Constructs influence matrices for all leaves of the tree. (Assumes source tree and target trees are identical.)
FastMultipole.set_gradient!
— Methodset_gradient!(target_buffer, i_body, gradient)
Accumulates gradient
to target_buffer
.
FastMultipole.set_hessian!
— Methodset_hessian!(target_buffer, i_body, hessian)
Accumulates hessian
to target_buffer
.
FastMultipole.set_scalar_potential!
— Methodset_scalar_potential!(target_buffer, i_body, scalar_potential)
Accumulates scalar_potential
to target_buffer
.
FastMultipole.shrink_branch!
— MethodComputes the smallest bounding box to completely bound all child boxes.
Shrunk radii are merely the distance from the center to the corner of the box.
FastMultipole.source_system_to_buffer!
— Methodsource_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_body
th 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.source_to_dipole!
— Methodassumes source expansion coefficients have already been calculated
FastMultipole.source_to_vortex_point!
— Methodmight be slightly faster than mirroredsourceto_vortex
FastMultipole.strength_dims
— Methodstrength_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.strength_to_value
— Methodstrength_to_value(strength, source_system)
NOTE: this function is primarily used for the boundary element solver, and is not required for the FMM.
Converts the strength of a body in source_system
to a scalar value. Should be overloaded for each user-defined system object used with the boundary element solver.
Arguments:
strength::SVector{dim, Float64}
: the strength of the body, wheredim
is the number of components in the strength vector (e.g., 1 for a point source, 3 for a point vortex, etc.)source_system::{UserDefinedSystem}
: the user-defined system object, used solely for dispatch
FastMultipole.target_influence_to_buffer!
— Methodtarget_influence_to_buffer!(target_buffer, i_buffer, ::DerivativesSwitch{PS,GS,HS}, target_system::{UserDefinedSystem}, i_target) where {PS,GS,HS}
NOTE: this function is primarily used for the boundary element solver, and is not required for the FMM.
Updates the target_buffer
with influences from target_system
for the i_target
th body. Should be overloaded for each user-defined system object (where {UserDefinedSystem}
is replaced with the type of the user-defined system) to be used as a target, assuming the target buffer positions have already been set. It should behave as follows:
target_buffer[4, i_buffer]
should be set to the scalar potential at the body positiontarget_buffer[5:7, i_buffer]
should be set to the vector field at the body positiontarget_buffer[8:16, i_buffer]
should be set to the vector gradient at the body position
The following convenience functions can may be used to access the buffer:
set_scalar_potential!(target_buffer, i_buffer, scalar_potential)
: accumulates thescalar_potential
to thei_buffer
body intarget_buffer
set_gradient!(target_buffer, i_buffer, gradient)
: accumulatesgradient
to thei_buffer
body intarget_buffer
set_hessian!(target_buffer, i_buffer, hessian)
: accumulateshessian
to thei_buffer
body intarget_buffer
FastMultipole.translate_local_z!
— MethodOverwrites translated_weights
FastMultipole.translate_multipole_to_local_z!
— MethodOverwrites translated_weights
FastMultipole.translate_multipole_to_local_z_m01_n
— MethodThis function should receive n!tnp1 = 1/t when n=0
Calculates ϕn0, ϕn1, χn1 for P=n
FastMultipole.translate_multipole_to_local_z_n!
— MethodOverwrites translated_weights
FastMultipole.translate_multipole_z!
— MethodOverwrites translated_weights
FastMultipole.tune_fmm
— Methodtune_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_order
keyword 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 tofalse
verbose::Bool
: whether to print progress information; defaults totrue
kwargs...
: 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.update_Ts_0!
— MethodReturns sβ, cβ, and 1n appropriate for calling updateTsn! for n=1
FastMultipole.update_eimϕs_0!
— MethodReturns eiϕreal, eiϕimag, eimϕreal, and eimϕimag appropriate for inputs for n=1 in updateeimϕsn!
FastMultipole.update_eimϕs_n!
— MethodReturns eimϕreal and eimϕimag for n
FastMultipole.update_nonself_influence!
— Methodupdate_nonself_influence!(right_hand_side, nonself_matrices::Matrices, old_influence_storage::Vector, source_tree::Tree, target_tree::Tree, strengths_by_leaf::Vector{UnitRange{Int}}, index_map::Vector{UnitRange{Int}}, direct_list::Vector{SVector{2,Int32}}, targets_by_branch::Vector{UnitRange{Int}})
Updates the right-hand side vector based on the non-self influence matrices and the current strengths of the source systems. Does this by removing the old influence and adding the new.
FastMultipole.value_to_strength!
— Methodvalue_to_strength!(source_buffer, source_system, i_body, value)
NOTE: this function is primarily used for the boundary element solver, and is not required for the FMM.
Converts a scalar value to a vector strength of a body in source_system
. Should be overloaded for each user-defined system object used with the boundary element solver.
Arguments:
source_buffer::Matrix{Float64}
: the source buffer containing the body informationsource_system::{UserDefinedSystem}
: the user-defined system object, used solely for dispatchi_body::Int
: the index of the body insource_buffer
to set the strength forvalue::Float64
: the scalar value used to set the strength
The following convenience function may be helpful when accessing the buffer:
get_strength(source_buffer, source_system, i_body)
: returns the strength of thei_body
body insource_buffer
, formatted as a vector (e.g.strength::SVector{dim,Float64}
)