Functions

BlockTensorFactorization.Core.circle_curvatureMethod
circle_curvature(y::AbstractVector{<:Real}; h=1, estimate_endpoints=true)

Inverse radius of a the circle passing through each 3 adjacent points on y, (0,y[i-1]), (h,y[i]), and (2h,y[i+1]).

If estimate_endpoints=true, assumes the function that y comes from is 1 to the left of the given values, and 0 to the right. This is typical of relative error decay as a function of rank. If false, pads the boundary with the adjacent curvature.

See three_point_circle.

source
BlockTensorFactorization.Core.coarsenMethod
coarsen(Y::AbstractArray, scale::Integer; dims=1:ndims(Y))

Coarsens or downsamples Y by scale. Only keeps every scale entries along the dimensions specified.

Example

Y = randn(12, 12, 12)

coarsen(Y, 2) == Y[begin:2:end, begin:2:end, begin:2:end]

coarsen(Y, 4; dims=(1, 3)) == Y[begin:4:end, :, begin:4:end]

coarsen(Y, 3; dims=2) == Y[:, begin:3:end, :]
source
BlockTensorFactorization.Core.contractionsMethod
contractions(D::AbstractDecomposition)

A tuple of functions defining a recipe for reconstructing a full array from the factors of the decomposition.

Example

(op1, op2) = contractions(D)
(A, B, C) = factors(D)

array(D) == op2(op1(A, B), C)
source
BlockTensorFactorization.Core.cpproductMethod
cpproduct((A, B, C, ...))
cpproduct(A, B, C, ...)

Multiplies the inputs by treating them as matrices in a CP decomposition.

Example

cpproduct(A, B, C) == @einsum T[i, j, k] := A[i, r] * B[j, r] * C[k, r]

source
BlockTensorFactorization.Core.cubic_spline_coefficientsMethod
cubic_spline_coefficients(y::AbstractVector{<:Real}; h=1)

Calculates the list of coefficients a, b, c, d for an interpolating spline of y[i]=y(x[i]).

The spline is defined as $f(x) = g_i(x)$ on $x[i] \leq x \leq x[i+1]$ where

$g_i(x) = a[i](x-x[i])^3 + b[i](x-x[i])^2 + c[i](x-x[i]) + d[i]$

Uses the following boundary conditions

  • $g_1(x[1]-h) = 1$ (i.e. the $y$-intercept is $(0,1)$ for uniform spaced x=1:I)
  • $g_I(x[I]+h) = y[I]$ (i.e. repeated right end-point)
  • $g_I''(x[I]+h) = 0$ (i.e. flat/no-curvature one spacing after end-point)
source
BlockTensorFactorization.Core.curvatureMethod
curvature(y::AbstractVector{<:Real}; method=:finite_differences)

Approximates the signed curvature of a function given evenly spaced samples.

Possible methods

  • :finite_differences: Approximates first and second derivative with 3rd order finite differences. See d_dx and d2_dx2.
  • :splines: Curvature of a third order spline. See d_dx_and_d2_dx2_spline.
  • :circles: Inverse radius of a circle through rolling three points. See circle_curvature.
  • :breakpoints: WARNING does not compute a value that approximates the curvature of a continuous function. Computes the inverse least-squares error of f.(eachindex(y); z) and y for all z in eachindex(y) where f(x; z) = a + b(min(x, z) - z) + c(max(x, z) - z). Useful if y looks like two lines. See breakpoint_curvature.
source
BlockTensorFactorization.Core.d2_dx2Method
d2_dx2(y::AbstractVector{<:Real})

Approximate second derivative with finite elements.

$\frac{d^2}{dx^2}y[i] \approx \frac{1}{\Delta x^2}(y[i-1] - 2y[i] + y[i+1])$

Assumes y[i] = y(x[i]) are samples with unit spaced inputs Δx = x[i+1] - x[i] = 1.

source
BlockTensorFactorization.Core.d_dxMethod
d_dx(y::AbstractVector{<:Real})

Approximate first derivative with finite elements.

$\frac{d}{dx}y[i] \approx \frac{1}{2\Delta x}(y[i+1] - y[i-1])$

Assumes y[i] = y(x[i]) are samples with unit spaced inputs Δx = x[i+1] - x[i] = 1.

source
BlockTensorFactorization.Core.default_kwargsMethod
default_kwargs(Y; kwargs...)

Handles all keywords and options, and sets defaults if not provided.

Keywords & Defaults

Initialization

  • decomposition: nothing. Can provide a custom initialized AbstractDecomposition. Note this exact decomposition is mutated.
  • model: Tucker1, but overridden by the type of AbstractDecomposition if given decomposition
  • rank: nothing, but overridden by the rank of AbstractDecomposition if given decomposition. Automatically calls rank_detect_factorize if both rank & decomposition are not provided.
  • init: abs_randn for nonnegative inputs Y, randn otherwise
  • constrain_init: true. Ensures the initialization satisfies all given constraints. Defaults to false if given decomposition
  • freeze: the default frozen factors of the model
  • continuous_dims: missing. Dimensions of Y that come from discretizations of continuous data. If provided, multiscale_factorize is called and can speed up factorization. If continuous_dims==nothing, factorization will only happen at one scale. In the future, if continuous_dims==missing,factorize may guess if there are continuous dimensions.

Updates

  • objective: L2(). Objective to minimize
  • norm: l2norm. Norm to use for statistics, can be unrelated to the objective
  • random_order: false. Perform the updates in a random order each iteration, Overrides to true when recursive_random_order=true
  • group_updates_by_factor: false. Groups updates on the same factor together. Overrides to true when random_order=true. Useful when randomizing order of updates but you want to keep matching momentum-gradientstep-constraint together
  • recursive_random_order: false. Performs inner blocked updates (grouped updates) in a random order (recursively) each iteration. Note the outer most list of updates can be performed in order if random_order=false
  • do_subblock_updates: false. Performs gradient descent on subblocks within a factor separately. May result in smaller Lipschitz constants and hence larger step sizes being used.

Momentum

  • momentum: true
  • δ: 0.9999. Amount of momentum, between [0,1)
  • previous_iterates: 1. Number of pervious iterates to save and use between iterations

Constraints

  • constraints: nothing. Can be a list of ConstraintUpdate, or just one
  • final_constraints: nothing. Constraints to apply after the final iteration. Will apply constraints if constrain_output is true and none are given
  • constrain_output: false. Apply the final_constraints, will override to true if final_constraints are given

Stats

  • stats: [Iteration, ObjectiveValue, GradientNorm] or in the case of nonnegative Y, GradientNNCone in place of GradientNorm
  • converged: GradientNorm or in the case of nonnegative Y, GradientNNCone. What stat(s) to use for convergence. Will converge is any one of the provided stats is below their respective tolerance
  • tolerance: 1. A list the same length as converged
  • maxiter: 1000. Additional stopping criterion if the number of iterations exceeds this number
source
BlockTensorFactorization.Core.eachfactorindexMethod
eachfactorindex(D::AbstractDecomposition)

An iterable to index the factors of D; not necessarily 1:ndims(D). Does not include non-data factors like the core of a CPDecomposition.

For example, eachfactorindex(D::Tucker) == 0:ndims(D) since core(D) is the zeroth factor. eachfactorindex(D::CPDecomposition) == 1:ndims(D). core(D) exists, but it is frozen as the identity tensor. And eachfactorindex(D::Tucker1) == 0:1 representing the core and matrix factor of D.

source
BlockTensorFactorization.Core.eachfibreMethod
eachfibre(A::AbstractArray; n::Integer, kwargs...)

Creates views of A that are that n-fibres of A.

Shorthand for eachslice(A; dims=(1,...,n-1,n+1,...,ndims(A))).

See eachslice.

source
BlockTensorFactorization.Core.finalconstrain!Method
finalconstrain!(decomposition; constraints, final_constraints, kwargs...)

Applies final_constraints (or if its nothing, applies constraints) to the decomposition.

Any RescaleUpdate are applied (the factor is scaled), but the rescaling of other factors is skipped.

source
BlockTensorFactorization.Core.frozenMethod
frozen(D::AbstractDecomposition)

A tuple of Bools the same length as factors(D) showing which factors are "frozen" in the sense that a block decent algorithm should skip these factors when decomposing a tensor.

See isfrozen.

source
BlockTensorFactorization.Core.geomeanMethod
geomean(v)
geomean(v...)

Geometric mean of a collection: prod(v)^(1/length(v)).

If prod(v) is detected to be 0 or Inf, the safer (but slower) implementation exp(mean(log.(v))) is used.

source
BlockTensorFactorization.Core.group_by_factorMethod
group_by_factor(blockedupdate::BlockedUpdate)

Groups updates according to the factor they operate on.

If blockedupdate contains other BlockedUpdates, the inner updates are grouped when they all operate on the same factor.

Updates which do not have an assigned factor are grouped together.

The order which these groups appear in the output follows the same order as the first appearence of each unique factor that is operated on.

source
BlockTensorFactorization.Core.interlaceMethod
interlace(u, v)

Takes two iterables, u and v, and alternates elements from u and v into a vector. If u and v are not the same length, extra elements are put on the end of the vector.

source
BlockTensorFactorization.Core.interpolateMethod
interpolate(Y, scale; dims=1:ndims(Y), degree=0, kwargs...)

Interpolates Y to a larger array with repeated values.

Keywords

scale. How much to scale up the size of Y. A dimension with size k will be scaled to scale*k - (scale - 1) = scale*(k-1) + 1

dims:1:ndims(Y). Which dimensions to interpolate.

degree:0. What degree of interpolation to use. 0 is constant interpolation, 1 is linear.

Like the opposite of coarsen.

Example

julia> Y = collect(reshape(1:6, 2, 3))
2×3 Matrix{Int64}:
 1  3  5
 2  4  6

julia> interpolate(Y, 2)
3×5 Matrix{Int64}:
 1  1  3  3  5
 1  1  3  3  5
 2  2  4  4  6

julia> interpolate(Y, 3; dims=2)
2×7 Matrix{Int64}:
 1  1  1  3  3  3  5
 2  2  2  4  4  4  6

julia> interpolate(Y, 1) == Y
true
source
BlockTensorFactorization.Core.isnonnegativeMethod
isnonnegative(X::AbstractArray{<:Real})
isnonnegative(x::Real)

Checks if all entries of X are bigger or equal to zero. Will be a standard function in Base but using this for now: https://github.com/JuliaLang/julia/pull/53677

source
BlockTensorFactorization.Core.make_splineMethod
make_spline(y::AbstractVector{<:Real}; h=1)

Returns a function f(x) that is an interpolating/extrapolating spline for y, with uniform stepsize h between the x-values of the knots.

source
BlockTensorFactorization.Core.make_update!Method
make_update!(decomposition, Y; momentum, constraints, constrain_init, group_updates_by_factor, do_subblock_updates, kwargs...)

What one iteration of the algorithm looks like. One iteration is likely a full cycle through each block or factor of the model.

source
BlockTensorFactorization.Core.max_possible_rankMethod
max_possible_rank(Y, model)

Returns the maximum rank possible Y could have under the model.

For matrices I × J this is 1:min(I, J). This is can be extended to tensors for different type of decompositions.

Tucker-1 rank is ≤ min(I, prod(J1,...,JN)) for tensors I × J1 × … × JN.

The CP-rank is ≤ minimum_{n} (prod(I1,...,IN) / In) for tensors I1 × … × IN in general. Although some shapes have have tighter upper bounds. For example, 2 × I × I tensors over ℝ have a maximum rank of floor(3I/2).

source
BlockTensorFactorization.Core.multiscale_factorizeMethod
multiscale_factorize(Y; continuous_dims=1:ndims(Y), rank=1, model=Tucker1, kwargs...)

Like factorize but uses progressively finer sub-grids of Y to speed up convergence. This is only effective when the dimensions given by dims come from discretizations of continuous data.

For example, if Y has 3 dimensions where Y[i, j, k] are samples from a continuous 2D function fi(xj, y_k) on a grid, use

multiscale_factorize(Y; continuous_dims=(2,3))

since second and third dimensions are continuous.

source
BlockTensorFactorization.Core.nmode_productMethod
nmode_product(A::AbstractArray, B::AbstractMatrix, n::Integer)

Contracts the nth mode of A with the first mode of B. Equivalent to A ×ₙ B where

(A ×ₙ B)[i₁, …, i_N] = ∑ⱼ A[i₁, …, iₙ₋₁, j, iₙ₊₁, …, i_N] B[iₙ, j].

source
BlockTensorFactorization.Core.nmode_productMethod
nmode_product(A::AbstractArray, b::AbstractVector, n::Integer)

Contracts the nth mode of A with b. Equivalent to A ×ₙ b where

(A ×ₙ b)[i₁, …, iₙ₋₁, iₙ₊₁, …, i_N] = ∑_iₙ A[i₁, …, i_N] b[iₙ].

source
BlockTensorFactorization.Core.outer_productMethod
outer_product(vectors)
outer_product(vectors...)

Outer product of a collection of vectors.

For example,

outer_product(u, v) == u * v'

and

outer_product(u, v, w)[i, j, k] == u[i] * v[j] * w[k].

Returned array will have same number dimensions as the length of the collection.

source
BlockTensorFactorization.Core.parse_constraintsMethod
parse_constraints(constraints, decomposition; kwargs...)

Parses the constraints to make sure we have a valid list of ConstraintUpdate.

If only one AbstractConstraint is given, assume we want this constraint to apply to every factor in the decomposition, and make a ConstraintUpdate for each factor.

If we are given a list of AbstractConstraint, assume we want them to apply to each factor of the decomposition in order.

source
BlockTensorFactorization.Core.proj_one_hotMethod
proj_one_hot(x)

Projects an array x to the closest one hot array: an array with all 0's except for a single 1. Does not mutate; see proj_one_hot! for a mutating version.

This is not a unique projection if there are multiple largest entries. In this case, will pick one of the largest entries to be 1 and set the rest to zero.

source
BlockTensorFactorization.Core.projsplx!Method
projsplx!(y; sum=one(eltype(y)))

Projects (in Euclidean distance) the array y into the simplex.

See projsplx for a non-mutating version.

[1] Yunmei Chen and Xiaojing Ye, "Projection Onto A Simplex", 2011

source
BlockTensorFactorization.Core.rank_detect_factorizeMethod
rank_detect_factorize(Y; kwargs...)

Wraps factorize() with rank detection.

Selects the rank that maximizes the standard curvature of the Relative Error (as a function of rank).

Keywords

  • online_rank_estimation: false. Set to true to stop testing larger ranks after the first peak in curvature
  • curvature_method: :splines. Can also pick :finite_differences (faster but less accurate) or circles (fastest and smallest memory but more sensitive to results from factorize). Set to :breakpoints to pick the rank R that minimizes least-squares error in the model f(r) = a + b(min(r, R) - R) + c(max(r, R) - R) and the errors.
  • model: Tucker1. Only rank detection with Tucker1 and CPDecomposition is currently implemented
  • max_rank: max_possible_rank(Y, model). Test ranks from 1 up to max_rank. Defaults to largest possible rank under the model
  • rank: nothing. If a rank is passed, rank detection is ignored and factorize(Y; kwargs...) is called

Any other keywords from factorize, full list given by default_kwargs.

source
BlockTensorFactorization.Core.rankofMethod
rankof(D::AbstractDecomposition)

Internal dimension sizes for a decomposition. Returns the sizes of all factors if not defined for a concrete subtype of AbstractDecomposition.

Examples

CPDecomposition: size of second dimension for the factors Tucker: size of the core factor Tucker1: size of the first dimension of the core factor

source
BlockTensorFactorization.Core.reshape_ndimsMethod
reshape_ndims(x, n)

Reshapes x to a higher order array with n dimensions.

When n > ndims(x), extra dimensions are prepended. Otherwise, trailing dimensions are collapsed.

Example

julia> x = [1, 2, 3]
3-element Vector{Int64}:
 1
 2
 3

julia> reshape_ndims(x, 1)
3-element Vector{Int64}:
 1
 2
 3

julia> reshape_ndims(x, 2)
1×3 Matrix{Int64}:
 1  2  3

julia> reshape_ndims(x, 3)
1×1×3 Array{Int64, 3}:
[:, :, 1] =
 1

[:, :, 2] =
 2

[:, :, 3] =
 3

julia> A = reshape(collect(1:12), 2, 2, 3)
2×2×3 Array{Int64, 3}:
[:, :, 1] =
 1  3
 2  4

[:, :, 2] =
 5  7
 6  8

[:, :, 3] =
  9  11
 10  12

julia> reshape_ndims(A, 2)
2×6 Matrix{Int64}:
 1  3  5  7   9  11
 2  4  6  8  10  12

Credit: https://discourse.julialang.org/t/outer-product-broadcast/103731/7

source
BlockTensorFactorization.Core.slicewise_dotMethod
slicewise_dot(A::AbstractArray, B::AbstractArray; dims=1, dimsA=dims, dimsB=dims)

Contracts all but the dimensions dimsA and dimsB of A and B. Entry-wise

(slicewise_dot(A,B; dimsA, dimsB))[i_dimsA, i_dimsB] = A[…, i_dimsA, ⋯] ⋅ B[…, i_dimsB, ⋯].

When dims==n, equivalent to A ⋅ₙ B. If A and B are both matrices, slicewise_dot(A, B) == A'B.

source
BlockTensorFactorization.Core.standard_curvatureMethod
standard_curvature(y::AbstractVector{<:Real}; method=:finite_differences)

Approximates the signed curvature of a function, scaled to the unit box $[0,1]^2$.

Assumes the function is 1 at 0 and (after x dimension is scaled) 0 at 1.

See curvature.

Possible methods

  • :finite_differences: Approximates first and second derivative with 3rd order finite differences. See d_dx and d2_dx2.
  • :splines: Curvature of a third order spline. See d_dx_and_d2_dx2_spline.
  • :circles: Inverse radius of a circle through rolling three points. See circle_curvature.
  • :breakpoints: WARNING does not compute a value that approximates the curvature of a continuous function. Computes the inverse least-squares error of f.(eachindex(y); z) and y for all z in eachindex(y) where f(x; z) = a + b(min(x, z) - z) + c(max(x, z) - z). Useful if y looks like two lines. See breakpoint_curvature.
source
BlockTensorFactorization.Core.three_point_circleMethod
three_point_circle((a,f),(b,g),(c,h))

Calculates radius r and center point (p, q) of the circle passing through the three points in the xy-plane.

Example

r, (p, q) = three_point_circle((1,2), (2,1), (5,2))
(r, (p, q)) == (√5, (3, 3))
source
BlockTensorFactorization.Core.tuckerproductMethod
tuckerproduct(G, (A, B, ...))
tuckerproduct(G, A, B, ...)

Multiplies the inputs by treating the first argument as the core and the rest of the arguments as matrices in a Tucker decomposition.

Example

tuckerproduct(G, (A, B, C)) == G ×₁ A ×₂ B ×₃ C
tuckerproduct(G, (A, B, C); exclude=2) == G ×₁ A ×₃ C
tuckerproduct(G, (A, B, C); exclude=2, excludes_missing=false) == G ×₁ A ×₃ C
tuckerproduct(G, (A, C); exclude=2, excludes_missing=true) == G ×₁ A ×₃ C
source