Reference

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
  • source_radius::TF: if this branch is a leaf, distance from center to the outer edge of farthest body contained in this branch; otherwise, this is the distance from center to the corner of its source_box
  • target_radius::TF: distance from center to the farthest body center contained in this branch
  • source_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, respectively
  • target_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, respectively
  • max_influence::TF: maximum influence of any body in this branch on any body in its child branches; used to enforce a relative error tolerance
source
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.DerivativesSwitchMethod
DerivativesSwitch(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 system
  • gradient::Vector{Bool}: a vector of ::Bool indicating whether the vector field should be computed for each target system
  • hessian::Vector{Bool}: a vector of ::Bool indicating whether the vector gradient should be computed for each target system
source
FastMultipole.DerivativesSwitchMethod
DerivativesSwitch(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 system
  • gradient::Bool: a ::Bool indicating whether the vector field should be computed for each target system
  • hessian::Bool: a ::Bool indicating whether the vector gradient should be computed for each target system
source
FastMultipole.DerivativesSwitchMethod
DerivativesSwitch(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 system
  • gradient::Bool: a ::Bool indicating whether the vector field should be computed for the target system
  • hessian::Bool: a ::Bool indicating whether the vector gradient should be computed for the target system
source
FastMultipole.ExpansionSwitchType
ExpansionSwitch

Switch indicating which expansions should be used:

  1. scalar potential (SP)
  2. vector potential via Lamb-Helmholtz decomposition (VP)
source
FastMultipole.ProbeSystemType
ProbeSystem

Convenience system for defining locations at which the potential, vector field, or vector gradient may be desired.

source
FastMultipole.TreeType

bodies[indexlist] is the same sort operation as performed by the tree sortedbodies[inverseindexlist] undoes the sort operation performed by the tree

source
FastMultipole.EmptyTreeMethod
EmptyTree(system)

Returns an empty tree. Used if system is empty.

Arguments

  • system: the system from which a tree is to be created

Returns

  • tree: if typeof(system)<:Tuple, a ::MultiTree is returned; otherwise, a ::SingleTree is returned
source
FastMultipole.allocate_buffersMethod
allocate_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!.

source
FastMultipole.allocate_small_buffersMethod
allocate_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.

source
FastMultipole.back_rotate_z!Method

Assumes eimϕs have already been computed. DOES NOT overwrite rotated weights (unlike other rotate functions); rather, accumulates on top of it.

source
FastMultipole.body_to_multipole!Method
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.buffer_to_system_strength!Method
buffer_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_bodyth 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 object
  • i_body::Int: the index of the body in source_system whose strength is to be set
  • source_buffer::Matrix{Float64}: the source buffer containing the body information
  • i_buffer::Int: the index of the body in source_buffer whose strength is to be set
source
FastMultipole.buffer_to_target_system!Method
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.calculate_ib!Method

assumes 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)

source
FastMultipole.calculate_pj!Method

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

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

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
source
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

  • target_buffers::Vector{<:Any}: buffers for target systems; if not provided, buffers are allocated using allocate_buffers
  • target_small_buffers::Vector{<:Any}: small buffers for target systems; if not provided, small buffers are allocated using allocate_small_buffers
  • source_buffers::Vector{<:Any}: buffers for source systems; if not provided, buffers are allocated using allocate_buffers
  • source_small_buffers::Vector{<:Any}: small buffers for source systems; if not provided, small buffers are allocated using allocate_small_buffers

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_recenter::Bool: whether to shrink and recenter branches around their bodies, accounting for finite body radius; default is true
  • interaction_list_method::InteractionListMethod: method for building interaction lists; default is SelfTuningTreeStop()

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
source
FastMultipole.get_n_bodiesMethod
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.get_normalMethod
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_positionMethod
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.index_by_sourceMethod
index_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.

source
FastMultipole.influence!Method
influence!(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 buffer
  • target_buffer::Matrix{TF}: target buffer used to compute the influence
  • source_system::{UserDefinedSystem}: system object used solely for dispatch
  • source_buffer::Matrix{TF}: source buffer used to compute the influence
source
FastMultipole.influence!Method
influence!(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 list
  • influences_per_system::Vector{Vector{Float64}}: vector of vectors containing the influence for each target system, sorted the same way as the buffers
  • target_buffers::NTuple{N,Matrix{Float64}}: target buffers used to compute the influence
  • source_systems::NTuple{N,<:{UserDefinedSystem}}: system objects used for dispatch
  • source_buffers::NTuple{N,Matrix{Float64}}: source buffers used to compute the influence
source
FastMultipole.local_errorMethod
Multipole coefficients up to expansion order `P+1` must be added to `multipole_branch` before calling this function.
source
FastMultipole.make_direct_assignments!Method
make_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.

source
FastMultipole.map_by_branchMethod
map_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.

source
FastMultipole.map_by_leafMethod
map_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.

source
FastMultipole.nearfield_device!Method
nearfield_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 which source_system acts
  • derivatives_switches::Union{DerivativesSwitch, NTuple{N,DerivativesSwitch}}: determines whether the scalar potential, vector field, and or vector gradient should be calculated
  • source_systems: user-defined system acting on target_system
source
FastMultipole.nearfield_device!Method
nearfield_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 which source_system acts
  • target_tree::Tree: octree object used to sort target_systems
  • derivatives_switches::Union{DerivativesSwitch, NTuple{N,DerivativesSwitch}}: determines whether the scalar potential, vector field, and or vector gradient should be calculated
  • source_systems: user-defined system acting on target_system
  • source_tree::Tree: octree object used to sort target_systems
  • direct_list::Vector{SVector{2,Int32}}: each element [i,j] maps nearfield interaction from source_tree.branches[j] on target_tree.branches[i]
source
FastMultipole.predict_errorMethod

performs 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

source
FastMultipole.rotate_multipole_y!Method

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

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

source
FastMultipole.shrink_branch!Method

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

source
FastMultipole.source_system_to_buffer!Method
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.strength_dimsMethod
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.strength_to_valueMethod
strength_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, where dim 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
source
FastMultipole.target_influence_to_buffer!Method
target_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_targetth 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 position
  • target_buffer[5:7, i_buffer] should be set to the vector field at the body position
  • target_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 the scalar_potential to the i_buffer body in target_buffer
  • set_gradient!(target_buffer, i_buffer, gradient): accumulates gradient to the i_buffer body in target_buffer
  • set_hessian!(target_buffer, i_buffer, hessian): accumulates hessian to the i_buffer body in target_buffer
source
FastMultipole.tune_fmmMethod
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.update_nonself_influence!Method
update_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.

source
FastMultipole.value_to_strength!Method
value_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 information
  • source_system::{UserDefinedSystem}: the user-defined system object, used solely for dispatch
  • i_body::Int: the index of the body in source_buffer to set the strength for
  • value::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 the i_body body in source_buffer, formatted as a vector (e.g. strength::SVector{dim,Float64})
source