Expansions#

class sumpy.expansion.ExpansionBase(kernel: Kernel, order: int, use_rscale: bool | None = None)[source]#
kernel#
order#
use_rscale#
abstract get_coefficient_identifiers() List[Hashable][source]#
Returns:

the identifiers of the coefficients that actually get stored.

abstract coefficients_from_source(kernel, avec, bvec, rscale, sac=None)[source]#

Form an expansion from a source point.

Parameters:
  • avec – vector from source to center.

  • bvec – vector from center to target. Not usually necessary, except for line-Taylor expansion.

  • sac – a symbolic assignment collection where temporary expressions are stored.

Returns:

a list of sympy expressions representing the coefficients of the expansion.

coefficients_from_source_vec(kernels, avec, bvec, rscale, weights, sac=None)[source]#

Form an expansion with a linear combination of kernels and weights.

Parameters:
  • avec – vector from source to center.

  • bvec – vector from center to target. Not usually necessary, except for line-Taylor expansion.

  • sac – a symbolic assignment collection where temporary expressions are stored.

Returns:

a list of sympy expressions representing the coefficients of the expansion.

loopy_expansion_formation(kernels: Sequence[Kernel], strength_usage: Sequence[int], nstrengths: int) TranslationUnit[source]#
Returns:

a loopy kernel that returns the coefficients for the expansion given by kernels with each kernel using the strength given by strength_usage.

abstract evaluate(kernel, coeffs, bvec, rscale, sac=None)[source]#
Returns:

a sympy expression corresponding to the evaluated expansion with the coefficients in coeffs.

loopy_evaluator(kernels: Sequence[Kernel]) TranslationUnit[source]#
Returns:

a loopy kernel that returns the evaluated target transforms of the potential given by kernels.

with_kernel(kernel: Kernel) ExpansionBase[source]#
copy(**kwargs) ExpansionBase[source]#
__len__()[source]#
__eq__(other)[source]#

Return self==value.

__ne__(other)[source]#

Return self!=value.

Expansion Wranglers#

class sumpy.expansion.ExpansionTermsWrangler(order: int, dim: int, max_mi: int | None = None)[source]#
order#
dim#
max_mi#
abstract get_coefficient_identifiers() List[Hashable][source]#
abstract get_full_kernel_derivatives_from_stored(stored_kernel_derivatives, rscale, sac=None)[source]#
abstract get_stored_mpole_coefficients_from_full(full_mpole_coefficients, rscale, sac=None)[source]#
get_full_coefficient_identifiers() List[Hashable][source]#

Returns identifiers for every coefficient in the complete expansion.

class sumpy.expansion.FullExpansionTermsWrangler(order: int, dim: int, max_mi: int | None = None)[source]#
class sumpy.expansion.LinearPDEBasedExpansionTermsWrangler(order: int, dim: int, knl: Kernel, max_mi: int | None = None)[source]#
__init__(order: int, dim: int, knl: Kernel, max_mi: int | None = None) None[source]#
Parameters:
  • order – order of the expansion

  • dim – number of dimensions

  • knl – kernel for the PDE

Expansion Factories#

class sumpy.expansion.ExpansionFactoryBase[source]#
abstract get_local_expansion_class(base_kernel)[source]#
Returns:

a subclass of ExpansionBase suitable for base_kernel.

abstract get_multipole_expansion_class(base_kernel)[source]#
Returns:

a subclass of ExpansionBase suitable for base_kernel.

class sumpy.expansion.DefaultExpansionFactory[source]#

An implementation of ExpansionFactoryBase that gives the β€˜best known’ expansion for each kernel.

class sumpy.expansion.VolumeTaylorExpansionFactory[source]#

An implementation of ExpansionFactoryBase that uses Volume Taylor expansions for each kernel.

Differential Operators#

Differential operator interface#

class sumpy.expansion.diff_op.LinearPDESystemOperator(dim, *eqs)[source]#

Represents a constant-coefficient linear differential operator of a vector-valued function with dim spatial variables. It is represented by a tuple of immutable dictionaries. The dictionary maps a DerivativeIdentifier to the coefficient. This object is immutable. Optionally supports a time variable as the last variable in the multi-index of the DerivativeIdentifier.

class sumpy.expansion.diff_op.DerivativeIdentifier(mi, vec_idx)#
sumpy.expansion.diff_op.make_identity_diff_op(ninput, noutput=1, time_dependent=False)[source]#

Returns the identity as a linear PDE system operator. if include_time is true, then the last dimension of the multi-index is time.

Parameters:
  • ninput – number of spatial variables of the function

  • noutput – number of output values of function

  • time_dependent – include time as a dimension

sumpy.expansion.diff_op.as_scalar_pde(pde: LinearPDESystemOperator, comp_idx: int) LinearPDESystemOperator[source]#

Returns a scalar PDE that is satisfied by the comp_idx component of pde.

To do this, we first convert a system of PDEs into a matrix where each row represents one PDE of the system of PDEs and each column represents a component. We convert a derivative to a polynomial expression and multiply that by the coefficient and enter that value into the matrix.

eg:

\[\begin{split}\frac{\partial^2 u}{\partial x^2} + \ 2 \frac{\partial^2 v}{\partial x y} = 0 \\ 3 \frac{\partial^2 u}{\partial y^2} + \ \frac{\partial^2 v}{\partial x^2} = 0\end{split}\]

is converted into,

\[\begin{split}\begin{bmatrix} x^2 & 2xy \\ 2y^2 & x^2 \end{bmatrix}.\end{split}\]

Let \(r_i\) be the columns of the above matrix. In order find a scalar PDE for the \(i\)-th component, we need to find some polynomials, \(a_1, a_2, \ldots, a_n\) such that the vector \(\sum_i a_i r_i\) has zeros except for the \(i\)-th component. In other words, we need to find a vector of polynomials such that the inner product of it with each of the columns except for the \(i\)-th column is zero. i.e. \(a_1, a_2, \ldots, a_n\) is a syzygy of all columns except for the \(i\)-th column.

To calculate a module that annihilates all but the \(i\)-th column, we first calculate the syzygy module of each column. A syzygy of a column vector is a row vector which has inner product zero with the column vector. A syzygy module is the module of all such syzygies and is generated by a finite set of syzygies. After calculating the syzygy module of each column, we intersect them except for the \(i\)-th module resulting in a syzygy module that when multiplied by the matrix gives a matrix with zeros EXCEPT for the \(i\)-th column. Therefore the vector \(a_1, a_2, \ldots, a_n\) is in this module.

When there are multiple such scalar PDEs, we want to get a combination of them that has the smallest positive degree. To do this, we calculate a groebner basis of the polynomials. We use the Groebner basis property that the largest monomial of each polynomial generated by set of polynomials is divisible by the largest monomial of some polynomial in the Groebner basis. When a graded monomial ordering is used this implies that the degree of any polynomial generated is greater than or equal to the degree of a polynomial in the Groebner basis. We choose that polynomial as our scalar PDE.

Parameters:
  • pde – An instance of LinearPDESystemOperator

  • comp_idx – the index of the component of the PDE solution vector for which a scalar PDE is requested.

Local Expansions#

class sumpy.expansion.local.LocalExpansionBase(kernel, order, use_rscale=None, m2l_translation=None)[source]#

Base class for local expansions.

abstract translate_from(src_expansion, src_coeff_exprs, src_rscale, dvec, tgt_rscale, sac=None, m2l_translation_classes_dependent_data=None)[source]#

Translate from a multipole or local expansion to a local expansion

Parameters:
  • src_expansion – The source expansion to translate from.

  • src_coeff_exprs – An iterable of symbolic expressions representing the coefficients of the source expansion.

  • src_rscale – scaling factor for the source expansion.

  • dvec – symbolic expression for the distance between target and source centers.

  • tgt_rscale – scaling factor for the target expansion.

  • sac – An object of type sumpy.assignment_collection.SymbolicAssignmentCollection to collect common subexpressions or None.

  • m2l_translation_classes_dependent_data – An iterable of symbolic expressions representing the expressions returned by translation_classes_dependent_data().

class sumpy.expansion.local.VolumeTaylorLocalExpansion(kernel, order, use_rscale=None, m2l_translation=None)[source]#
class sumpy.expansion.local.H2DLocalExpansion(kernel, order, use_rscale=None, m2l_translation=None)[source]#
class sumpy.expansion.local.Y2DLocalExpansion(kernel, order, use_rscale=None, m2l_translation=None)[source]#
class sumpy.expansion.local.LineTaylorLocalExpansion(kernel, order, use_rscale=None, m2l_translation=None)[source]#

Multipole Expansions#

class sumpy.expansion.multipole.VolumeTaylorMultipoleExpansion(kernel, order, use_rscale=None)[source]#
class sumpy.expansion.multipole.H2DMultipoleExpansion(kernel, order, use_rscale=None)[source]#
class sumpy.expansion.multipole.Y2DMultipoleExpansion(kernel, order, use_rscale=None)[source]#

Multipole to Local Translations#

class sumpy.expansion.m2l.M2LTranslationClassFactoryBase[source]#
abstract get_m2l_translation_class(base_kernel, local_expansion_class)[source]#

Returns a subclass of M2LTranslationBase suitable for base_kernel and local_expansion_class.

class sumpy.expansion.m2l.NonFFTM2LTranslationClassFactory[source]#

An implementation of M2LTranslationClassFactoryBase that uses non FFT M2L translation class.

class sumpy.expansion.m2l.FFTM2LTranslationClassFactory[source]#

An implementation of M2LTranslationClassFactoryBase that uses FFT M2L translation class.

class sumpy.expansion.m2l.DefaultM2LTranslationClassFactory[source]#

An implementation of M2LTranslationClassFactoryBase that gives the β€˜best known’ translation type for each kernel and local expansion class

class sumpy.expansion.m2l.M2LTranslationBase[source]#

Base class for Multipole to Local Translation

abstract translate(tgt_expansion, src_expansion, src_coeff_exprs, src_rscale, dvec, tgt_rscale, sac=None, translation_classes_dependent_data=None)[source]#
loopy_translate(tgt_expansion, src_expansion)[source]#
translation_classes_dependent_data(tgt_expansion, src_expansion, src_rscale, dvec, sac) Tuple[Any][source]#

Return an iterable of expressions that needs to be precomputed for multipole-to-local translations that depend only on the distance between the multipole center and the local center which is given as dvec.

Since there are only a finite number of different values for the distance between per level, these can be precomputed for the tree. In boxtree, these distances are referred to as translation classes.

When FFT is turned on, the output expressions are assumed to be transformed into Fourier space at the end by the caller.

translation_classes_dependent_ndata(tgt_expansion, src_expansion)[source]#

Return the number of expressions returned by translation_classes_dependent_data(). This method exists because calculating the number of expressions using the above method might be costly and translation_classes_dependent_data() cannot be memoized due to it having side effects through the argument sac.

abstract preprocess_multipole_exprs(tgt_expansion, src_expansion, src_coeff_exprs, sac, src_rscale)[source]#

Return the preprocessed multipole expansion for an optimized M2L. Preprocessing happens once per source box before M2L translation is done.

These expressions are used in a separate loopy kernel to avoid having to process for each target and source box pair. When FFT is turned on, the output expressions are assumed to be transformed into Fourier space at the end by the caller. When FFT is turned off, the output expressions are equal to the multipole expansion coefficients with zeros added to make the M2L computation a circulant matvec.

preprocess_multipole_nexprs(tgt_expansion, src_expansion)[source]#

Return the number of expressions returned by preprocess_multipole_exprs(). This method exists because calculating the number of expressions using the above method might be costly and it cannot be memoized due to it having side effects through the argument sac.

abstract postprocess_local_exprs(tgt_expansion, src_expansion, m2l_result, src_rscale, tgt_rscale, sac)[source]#

Return postprocessed local expansion for an optimized M2L. Postprocessing happens once per target box just after the M2L translation is done and before storing the expansion coefficients for the local expansion.

When FFT is turned on, the output expressions are assumed to have been transformed from Fourier space back to the original space by the caller.

postprocess_local_nexprs(tgt_expansion, src_expansion)[source]#

Return the number of expressions given as input to postprocess_local_exprs(). This method exists because calculating the number of expressions using the above method might be costly and it cannot be memoized due to it having side effects through the argument sac.

use_fft: ClassVar[bool] = False#
use_preprocessing: ClassVar[bool] = False#
class sumpy.expansion.m2l.VolumeTaylorM2LTranslation[source]#
class sumpy.expansion.m2l.VolumeTaylorM2LWithFFT[source]#
class sumpy.expansion.m2l.FourierBesselM2LTranslation[source]#

Estimating Expansion Orders#

class sumpy.expansion.level_to_order.FMMLibExpansionOrderFinder(tol, extra_order=0)[source]#

Return expansion orders that meet the tolerance for a given level using routines wrapped from pyfmmlib.

__init__(tol, extra_order=0)[source]#
Parameters:
  • tol – error tolerance

  • extra_order – order increase to accommodate, say, the taking of derivatives of the FMM expansions.

__call__(kernel, kernel_args, tree, level)[source]#

Call self as a function.

class sumpy.expansion.level_to_order.SimpleExpansionOrderFinder(tol, err_const_laplace=0.01, err_const_helmholtz=100, scaling_const_helmholtz=4, extra_order=1)[source]#

This models the Laplace truncation error as:

\[C_{\text{lap}} \left(\frac{\sqrt{d}}{3}\right)^{p+1}.\]

For the Helmholtz kernel, an additional term is added:

\[C_{\text{helm}} \frac 1{p!} \left(C_{\text{helmscale}} \cdot \frac{hk}{2\pi}\right)^{p+1},\]

where \(d\) is the number of dimensions, \(p\) is the expansion order, \(h\) is the box size, and \(k\) is the wave number.

__init__(tol, err_const_laplace=0.01, err_const_helmholtz=100, scaling_const_helmholtz=4, extra_order=1)[source]#
Parameters:

extra_order – order increase to accommodate, say, the taking of derivatives of the FMM expansions.

__call__(kernel, kernel_args, tree, level)[source]#

Call self as a function.