Docstrings

This page lists all docstrings in Polynomials4ML.jl including for functions that are not part of the public API. Please check with the Public API which functionality for which we aim to guarantee semver-stability.

Polynomials4ML.AbstractP4MLBasisType

abstract type AbstractP4MLBasis end

Annotates types that map a low-dimensional input, scalar or SVector, to a vector of scalars (feature vector, embedding, basis...).

source
Polynomials4ML.AbstractP4MLTensorType

abstract type AbstractP4MLTensor end

Annotates layers that map a vector to a vector. Each of the vectors may represent a tensor (hence the name). Future interfaces may generalize the allowed dimensionality of inputs to allow tensorial shapes.

source
Polynomials4ML.CTrigBasisType

Complex trigonometric polynomials up to degree N (inclusive). The basis is constructed in the order

[1, exp(im*θ), exp(-im*θ), exp(2im*θ), exp(-2im*θ), ..., 
                                exp(N*im*θ), exp(-N*im*θ) ]

where θ is input variable.

source
Polynomials4ML.ChebBasisType

ChebBasis(N):

Chebyshev polynomials up to degree N-1 (inclusive). i.e basis with length N. The basis is ordered as

\[[1, x, 2x^2-1, 4x^3-3x, ..., 2xT_{N-1}(x)-T_{N-2}(x)]\]

where x is input variable.

The differences between ChebBasis and chebyshev_basis is that ChebBasis computes the basis on the fly when it is compiled and it does not store the recursion coefficients as in chebyshev_basis. There might be a small performance benefit from this.

Warning: ChebBasis and chebyshev_basis have different normalization.

source
Polynomials4ML.LinearLayerType

struct LinearLayer : This lux layer returns W * x if feature_first is true, otherwise it returns x * transpose(W), where W is the weight matrix`.

  • x::AbstractMatrix of size (in_dim, N) or (N, in_dim), where in_dim = feature dimension, N = batch size
  • W::AbstractMatrix of size (out_dim, in_dim)

Constructor

LinearLayer(in_dim, out_dim; feature_first = false)

Example

in_d, out_d = 4, 3 # feature dimensions
N = 10 # batch size

# feature_first = true
l = P4ML.LinearLayer(in_d, out_d; feature_first = true)
ps, st = LuxCore.setup(MersenneTwister(1234), l)
x = randn(in_d, N) # feature-first
out, st = l(x, ps, st)
println(out == W * x) # true

# feature_first = false
l2 = P4ML.LinearLayer(in_d, out_d; feature_first = true)
ps2, st2 = LuxCore.setup(MersenneTwister(1234), l2)
x = randn(N, in_d) # batch-first
out, st = l(x, ps, st)
println(out == x * transpose(W))) # true
source
Polynomials4ML.MonoBasisType

Standard Monomials basis. This should very rarely be used. Possibly useful in combination with a transformation of the inputs, e.g. exponential.

source
Polynomials4ML.OrthPolyBasis1D3TType

OrthPolyBasis1D3T: defines a basis of polynomials in terms of a 3-term recursion,

\[\begin{aligned} P_1(x) &= A_1 \\ P_2 &= A_2 x + B_2 \\ P_{n} &= (A_n x + B_n) P_{n-1}(x) + C_n P_{n-2}(x) \end{aligned}\]

Typically (but not necessarily) such bases are obtained by orthogonalizing the monomials with respect to a user-specified distribution, which can be either continuous or discrete but must have a density function. See also

  • legendre_basis
  • chebyshev_basis
  • jacobi_basis
source
Polynomials4ML.PooledSparseProductType

struct PooledSparseProduct : This implements a fused (tensor) product and pooling operation. Suppose we are given $N$ embeddings $\phi^{(i)}_{k_i}$ then the pooled sparse product generates feature vectors of the form

\[A_{k_1, \dots, k_N} = \sum_{j} \prod_{t = 1}^N \phi^{(t)}_{k_t}(x_j)\]

where $x_j$ are an list of inputs (multi-set).

Constructor

PooledSparseProduct(spec)

where spec is a list of $(k_1, \dots, k_N)$ tuples or vectors, or AbstractMatrix where each column specifies such a tuple.

source
Polynomials4ML.RTrigBasisType

RTrigBasis(N):

Real trigonometric polynomials up to degree N (inclusive). The basis is ordered as

[1, cos(θ), sin(θ), cos(2θ), sin(2θ), ..., cos(Nθ), sin(Nθ) ]

where θ is input variable.

source
Polynomials4ML.SparseSymmProdType

SparseSymmProd : sparse symmetric product with entries stored as tuples. Input is a vector A; each entry of the output vector AA is of the form

\[ {\bm A}_{i_1, \dots, i_N} = \prod_{t = 1}^N A_{i_t}.\]

Constructor

SparseSymmProd(spec)

where spec is a list of tuples or vectors, each of which specifies an AA basis function as described above. For example,

spec = [ (1,), (2,), (1,1), (1,2), (2,2), 
         (1,1,1), (1,1,2), (1,2,2), (2,2,2) ]
basis = SparseSymmProd(spec)         

defines a basis of 9 functions,

\[[ A_1, A_2, A_1^2, A_1 A_2, A_2^2, A_1^3, A_1^2 A_2, A_1 A_2^2, A_2^3 ]\]

source
Polynomials4ML.SparseSymmProdDAGType

struct SparseSymmProdDAG : alternative (recursive) implementation of SparseSymmProd. This has better theoretical performance for high correlation orders.

The potential downside is that it inserts auxiliary basis functions into the basis. This means, that the specification of the output will be different from the specification that is used to construct the basis. To that end, the field projection can be used to reduce it back to the original spec. E.g.,

basis = SparseSymmProd(spec)
basis_dag = SparseSymmProdDAG(spec)
A = randn(nA)
basis(A) ≈ basis_dag(A)[basis_dag.projection]   # true

However, the field projection is used only for information, and not to actually reduce the output. One could of course use it to compose the output with a projection matrix.

source
Polynomials4ML.SparseSymmProdDAGMethod

Construct the DAG used to evaluate an AA basis and returns it as a SparseSymmProdDAG

Arguments

  • spec : AA basis specification, list of vectors of integers / indices pointing into A

Kwargs:

  • filter = _-> true :
  • verbose = false : print some information about the
source
Polynomials4ML.SphericalCoordsType

struct SphericalCoords : a simple datatype storing spherical coordinates of a point (x,y,z) in the format (r, cosφ, sinφ, cosθ, sinθ). Use spher2cart and cart2spher to convert between cartesian and spherical coordinates.

source
Polynomials4ML.chebyshev_basisMethod

chebyshev_basis(N::Integer):

Constructs an OrthPolyBasis1D3T object representing a possibly rescaled version of the basis of Chebyshev polynomials of the first kind. N is the length of the basis, not the degree.

Careful: the normalisation may be non-standard.

source
Polynomials4ML.jacobi_basisMethod

jacobi_basis(N::Integer, α::Real, β::Real):

Constructs an OrthPolyBasis1D3T object representing a possibly rescaled version of the basis of Jacobi polynomials Jαβ. N is the length of the basis, not the degree.

Careful: the normalisation may be non-standard.

source
Polynomials4ML.legendre_basisMethod

legendre_basis(N::Integer):

Constructs an OrthPolyBasis1D3T object representing a possibly rescaled version of the basis of Legendre polynomials (L2 orthonormal on [-1, 1]). N is the length of the basis, not the degree.

Careful: the normalisation may be non-standard.

source
Polynomials4ML.luxMethod

lux(basis) : convert a basis / embedding object into a lux layer. This assumes that the basis accepts a number or short vector as input and produces an output that is a vector. It also assumes that batched operations are implemented, as well as some other functionality.

source
Polynomials4ML.orthpolybasisMethod

function orthpolybasis(...) : construct a univariate orthogonal polynomial basis with respect to some specified inner product. For the standard 3-term recursion polynomials, use legendre_basis, jacobi_basis or chebyshev_basis.

The orthpolybasis currently implements orthogonal polynomials for discrete weights with the following constructors:

orthpolybasis(N::Integer, W::DiscreteWeights{TW}; TX = Float64)
orthpolybasis(N::Integer, X::AbstractVector{<: Real}, 
              W::AbstractVector{<: Real}, normalizeW=false; kwargs...)
source
Polynomials4ML.Utils.gensparseMethod

gensparse(...) : utility function to generate high-dimensional sparse grids which are downsets. All arguments are keyword arguments (with defaults):

  • NU : maximum correlation order
  • minvv = 0 : minvv[i] gives the minimum value forvv[i]`
  • maxvv = Inf : maxvv[i] gives the minimum value forvv[i]`
  • tup2b = vv -> vv :
  • admissible = _ -> false : determines whether a tuple belongs to the downset
  • filter = _ -> true : a callable object that returns true of tuple is to be kept and

false otherwise (whether or not it is part of the downset!) This is used, e.g. to enfore conditions such as ∑ lₐ = even or |∑ mₐ| ≦ M

  • INT = Int : integer type to be used
  • ordered = false : whether only ordered tuples are produced; ordered tuples

correspond to permutation-invariant basis functions

source