fmm.backend package¶
Subpackages¶
Submodules¶
fmm.backend.api module¶
Interface for compute backends.
fmm.backend.numba module¶
Compute operators, accelerated with Numba.
- fmm.backend.numba.l2l(key, local_expansions, l2l, key_to_index, nequivalent_points)¶
L2L operator applied to a parent key.
- keynp.int64
Operator applied to this key.
- local_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Local expansions, aligned by global index from key_to_index.
- l2lnp.array(shape=(8, ncheck_points, nequivalent_points), dtype=float)
Precomputed L2L operators.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- fmm.backend.numba.l2t(leaves, nleaves, key_to_index, key_to_leaf_index, targets, target_potentials, target_index_pointer, local_expansions, x0, r0, alpha_outer, equivalent_surface, nequivalent_points, p2p_parallel_function, dtype)¶
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_potentialsnp.array(shape=(ntargets, 4), dtype=float)
Target potentials (component 0), and x, y, and z components of potential gradient (components 1, 2, 3 resp.).
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- local_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Local expansions, aligned by global index from key_to_index.
- x0np.array(shape=(1, 3), dtype=np.float64)
Center of octree root node.
- r0np.float64
Half side length of octree root node.
- alpha_outerfloat
Relative size of outer surface.
- equivalent_surfacenp.array(shape=(nequivalent_points, 3), dtype=float)
Discretized equivalent surface.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- p2p_parallel_functionfunction
Parallel P2P function handle for this kernel.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
- fmm.backend.numba.m2l(keys, v_lists, u, s, vt, dc2e_inv_a, dc2e_inv_b, multipole_expansions, local_expansions, nequivalent_points, key_to_index, hash_to_index, scale)¶
M2L operator. Parallelized over all keys in a given level.
- keysnp.array(shape=(nkeys), dtype=np.int64)
All nodes at a given level of the octree.
- v_listsnp.array(shape=(ncomplete, nv_list), dtype=np.int64)
All V lists for nodes in octree.
- unp.array(shape=(ncheck_points, k), dtype=float)
Left singular vectors of M2L matrix for all transfer vectors at this level.
- snp.array(shape=(k, k), dtype=float)
Singular values of M2L matrix for all transfer vectors at this level.
- vtnp.array(shape=(k, nequivalent_points*316))
Right singular values of M2L matrix for all transfer vectors at this level.
- dc2e_inv_anp.array(dtype=float)
First component of inverse of downward check to equivalent Gram matrix.
- dc2e_inv_bnp.array(dtype=float)
Second component of inverse of downward check to equivalent Gram matrix.
- multipole_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Multipole expansions, aligned by global index from key_to_index.
- local_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Local expansions, aligned by global index from key_to_index.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- hash_to_indexnumba.typed.Dict.empty(key_type=numba.types.int64, value_type=numba.types.int64)
Map between transfer vector hash and index of right singular vector components at a given level.
- scalefloat
Kernel scale for keys at this level.
- fmm.backend.numba.m2l_core(key, v_list, u, s, vt, dc2e_inv_a, dc2e_inv_b, local_expansions, multipole_expansions, nequivalent_points, key_to_index, hash_to_index, scale)¶
Application of the M2L operator over the V list of a given node.
- keynp.int64
Operator applied to this key.
- v_listnp.array(shape=(nv_list), dtype=np.int64)
V list of this key.
- unp.array(shape=(ncheck_points, k), dtype=float)
Left singular vectors of M2L matrix for all transfer vectors at this level.
- snp.array(shape=(k, k), dtype=float)
Singular values of M2L matrix for all transfer vectors at this level.
- vtnp.array(shape=(k, nequivalent_points*316))
Right singular values of M2L matrix for all transfer vectors at this level.
- dc2e_inv_anp.array(dtype=float)
First component of inverse of downward check to equivalent Gram matrix.
- dc2e_inv_bnp.array(dtype=float)
Second component of inverse of downward check to equivalent Gram matrix.
- multipole_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Multipole expansions, aligned by global index from key_to_index.
- local_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Local expansions, aligned by global index from key_to_index.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- hash_to_indexnumba.typed.Dict.empty(key_type=numba.types.int64, value_type=numba.types.int64)
Map between transfer vector hash and index of right singular vector components at a given level.
- scalefloat
Kernel scale for keys at this level.
- fmm.backend.numba.m2m(keys, multipole_expansions, m2m, key_to_index, nequivalent_points)¶
M2M operator serially applied over keys in a given level. Parallelization doesn’t offer much benefit due to blocking writes to the parent multipole expansion.
As it’s not possible to know group siblings contained in the tree apriori due to the relative cost of searching the tree for siblings, it’s not possible to parallelize the operation over groups of siblings.
- keysnp.array(shape=(nkeys), dtype=np.int64)
All nodes at a given level of the octree.
- multipole_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Multipole expansions, aligned by global index from key_to_index.
- m2mnp.array(shape=(8, ncheck_points, nequivalent_points), dtype=float)
Precomputed M2M operators.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- fmm.backend.numba.m2t(leaves, nleaves, w_lists, targets, target_index_pointer, key_to_index, key_to_leaf_index, target_potentials, multipole_expansions, x0, r0, alpha_inner, equivalent_surface, nequivalent_points, p2p_function, gradient_function)¶
M2T operator parallelized over leaves.
The potential size of W lists make pre-allocation (c.f. L2T, near_field*) very expensive for the M2T operator, hence we opt for simple parallelization over leaves.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- w_listsnp.array(shape=(ncomplete, nw_list), dtype=np.int64)
All W lists for nodes in octree.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- target_potentialsnp.array(shape=(ntargets, 4), dtype=float)
Target potentials (component 0), and x, y, and z components of potential gradient (components 1, 2, 3 resp.).
- multipole_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Multipole expansions, aligned by global index from key_to_index.
- x0np.array(shape=(1, 3), dtype=np.float64)
Center of octree root node.
- r0np.float64
Half side length of octree root node.
- alpha_innerfloat
Relative size of inner surface.
- equivalent_surfacenp.array(shape=(nequivalent_points, 3), dtype=float)
Discretized equivalent surface.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- p2p_functionfunction
Serial P2P function handle for this kernel.
- grad_functionfunction
Serial gradient function handle for this kernel.
- fmm.backend.numba.near_field_node(leaves, nleaves, key_to_leaf_index, targets, target_index_pointer, sources, source_densities, source_index_pointer, max_points, target_potentials, p2p_parallel_function, dtype)¶
Compute near field interactions for sources and targets within each node.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
- source_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning sources by leaf.
- max_pointsint
Max points allowed in a node.
- target_potentialsnp.array(shape=(ntargets, 4), dtype=float)
Target potentials (component 0), and x, y, and z components of potential gradient (components 1, 2, 3 resp.).
- p2p_parallel_functionfunction
Parallel P2P function handle for this kernel.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
- fmm.backend.numba.near_field_u_list(leaves, nleaves, targets, target_index_pointer, sources, source_densities, source_index_pointer, key_to_index, key_to_leaf_index, u_lists, max_points, target_potentials, p2p_parallel_function, dtype)¶
Evaluate all near fields particles, corresponding to members of the U list, for all leaves.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
- source_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning sources by leaf.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- u_listsnp.array(shape=(ncomplete, nu_list), dtype=np.int64)
All U lists for nodes in octree.
- max_pointsint
Max points allowed in a node.
- target_potentialsnp.array(shape=(ntargets, 4), dtype=float)
Target potentials (component 0), and x, y, and z components of potential gradient (components 1, 2, 3 resp.).
- p2p_parallel_functionfunction
Parallel P2P function handle for this kernel.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
- fmm.backend.numba.p2m(leaves, nleaves, key_to_index, key_to_leaf_index, sources, source_densities, source_index_pointer, multipole_expansions, nequivalent_points, x0, r0, alpha_outer, check_surface, ncheck_points, uc2e_inv_a, uc2e_inv_b, p2p_function, scale_function, dtype)¶
P2M operator. Compute the multipole expansion from the sources at each leaf node.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
- source_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning sources by leaf.
- multipole_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Multipole expansions, aligned by global index from key_to_index.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- x0np.array(shape=(1, 3), dtype=np.float64)
Center of octree root node.
- r0np.float64
Half side length of octree root node.
- alpha_outerfloat
Relative size of outer surface.
- check_surfacenp.array(shape=(ncheck_points, 3), dtype=float)
Discretized check surface.
- ncheck_pointsint
Number of quadrature points on check_surface.
- uc2e_inv_anp.array(dtype=float)
First component of inverse of upward check to equivalent Gram matrix.
- uc2e_inv_bnp.array(dtype=float)
Second component of inverse of upward check to equivalent Gram matrix.
- p2p_functionfunction
Serial P2P function handle for this kernel.
- scale_functionfunction
Scale function handle for this kernel.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
- fmm.backend.numba.p2m_core(leaves, nleaves, key_to_index, nequivalent_points, ncheck_points, uc2e_inv_a, uc2e_inv_b, scales, multipole_expansions, check_potentials)¶
Parallelized loop applying P2M operator over all leaves.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- ncheck_pointsint
Number of quadrature points on check_surface.
- uc2e_inv_anp.array(dtype=float)
First component of inverse of upward check to equivalent Gram matrix.
- uc2e_inv_bnp.array(dtype=float)
Second component of inverse of upward check to equivalent Gram matrix.
- scalesnp.array(shape=(nleaves), dtype=float)
Scales calculated in prepare_p2m_data function
- multipole_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Multipole expansions, aligned by global index from key_to_index.
- check_potentialsnp.array(shape=(nleaves*ncheck_points), dtype=float)
Check potentials calculated in prepare_p2p_data function.
- fmm.backend.numba.prepare_l2t_data(leaves, nleaves, targets, target_index_pointer, equivalent_surface, nequivalent_points, x0, r0, alpha_outer, key_to_index, key_to_leaf_index, local_expansions, dtype)¶
Prepare L2T data to maximise cache re-use by allocating aligned arrays aligning source and target coordinates.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- equivalent_surfacenp.array(shape=(nequivalent_points, 3), dtype=float)
Discretized equivalent surface.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- x0np.array(shape=(1, 3), dtype=np.float64)
Center of octree root node.
- r0np.float64
Half side length of octree root node.
- alpha_outerfloat
Relative size of outer surface.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- local_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Local expansions, aligned by global index from key_to_index.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
- fmm.backend.numba.prepare_node_data(leaves, nleaves, key_to_leaf_index, targets, target_index_pointer, sources, source_densities, source_index_pointer, max_points, dtype)¶
Prepare source and target data within each node in order to maximize cache re-use.
Uses same strategy as over U lists.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
- source_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning sources by leaf.
- max_pointsint
Max points allowed in a node.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
5-tuple containing the re-batched sources, targets, source densities, source index pointers and target index pointers respectively.
- fmm.backend.numba.prepare_p2m_data(leaves, nleaves, sources, source_densities, source_index_pointer, key_to_leaf_index, x0, r0, alpha_outer, check_surface, ncheck_points, p2p_function, scale_function, dtype)¶
- Create aligned vector of scales, and check potentials indexed by leaves.
This maximizes cache re-use in the application of the P2M operator.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
source_index_pointer : np.array(shape=(nleaves+1), dtype=int) key_to_to_leaf_index : numba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- x0np.array(shape=(1, 3), dtype=np.float64)
Center of octree root node.
- r0np.float64
Half side length of octree root node.
- alpha_outerfloat
Relative size of outer surface.
- check_surfacenp.array(shape=(ncheck_points, 3), dtype=float)
Discretized check surface.
- ncheck_pointsint
Number of quadrature points on check_surface.
- p2p_functionfunction
Serial P2P function handle for this kernel.
- scale_functionfunction
Scale function handle for this kernel.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
- (np.array(float), np.array(float))
Tuple of scales and check potentials (ordered by leaf index) respectively.
- fmm.backend.numba.prepare_u_list_data(leaves, nleaves, targets, target_index_pointer, sources, source_densities, source_index_pointer, key_to_index, key_to_leaf_index, u_lists, max_points, dtype)¶
Prepare U list data to maximise cache re-use by allocating aligned arrays aligning source and target coordinates.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- targetsnp.array(shape=(ntargets, 3), dtype=float)
All target coordinates.
- target_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning targets by leaf.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
- source_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning sources by leaf.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- u_listsnp.array(shape=(ncomplete, nu_list), dtype=np.int64)
All U lists for nodes in octree.
- max_pointsint
Max points allowed in a node.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.
5-tuple containing the re-batched sources, targets, source densities, source index pointers and target index pointers respectively.
- fmm.backend.numba.s2l(leaves, nleaves, sources, source_densities, source_index_pointer, key_to_index, key_to_leaf_index, x_lists, local_expansions, x0, r0, alpha_inner, check_surface, nequivalent_points, dc2e_inv_a, dc2e_inv_b, scale_function, p2p_function, dtype)¶
S2L operator, parallelized simply over leaves.
Parallelization of the S2L operator doesn’t maximize cache re-use (c.f. P2M, L2T) This decision is taken as X lists are usually small, with relatively few nodes having X lists either. This makes the savings due to cache re-use competitive with the cost of array allocation.
- leaves: np.array(shape=(nleaves), dtype=np.int64)
Octree leaves.
- nleavesint
Number of leaves.
- sourcesnp.array(shape=(nsources, 3), dtype=float)
All source coordinates.
- source_densitiesnp.array(shape=(nsources), dtype=float)
Densities at source coordinates.
- source_index_pointernp.array(shape=(nleaves+1), dtype=int)
Index pointer aligning sources by leaf.
- key_to_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to global index.
- key_to_leaf_indexnumba.typed.Dict(key_type=np.int64, value_type=np.int64)
Map from key to leaf index.
- x_listsnp.array(shape=(ncomplete, nx_list), dtype=np.int64)
All X lists for nodes in octree.
- local_expansionsnp.array(shape=(nequivalent_points*ncomplete), dtype=float)
Local expansions, aligned by global index from key_to_index.
- x0np.array(shape=(1, 3), dtype=np.float64)
Center of octree root node.
- r0np.float64
Half side length of octree root node.
- alpha_innerfloat
Relative size of inner surface.`
- check_surfacenp.array(shape=(ncheck_points, 3), dtype=float)
Discretized check surface.
- nequivalent_pointsint
Number of quadrature points on equivalent_surface.
- dc2e_inv_anp.array(dtype=float)
First component of inverse of downward check to equivalent Gram matrix.
- dc2e_inv_bnp.array(dtype=float)
Second component of inverse of downward check to equivalent Gram matrix.
- scale_functionfunction
Scale function handle for this kernel.
- p2p_functionfunction
Serial P2P function handle for this kernel.
- dtypetype
Corresponds to precision of experiment ∈ {np.float32, np.float64}.