Skip to content

Maps

Index

Types

TransportMaps.PolynomialMap Type
julia
PolynomialMap <: AbstractTriangularMap

Triangular transport map with polynomial basis.

Fields

  • components::Vector{PolynomialMapComponent{<:AbstractPolynomialBasis}}: Vector of map components

  • reference::MapReferenceDensity: Reference density of the map

  • forwarddirection::Symbol: :target for map from density, :reference for map from samples

Constructors

  • PolynomialMap(dimension::Int,degree::Int, referencetype::Symbol=:normal, rectifier::AbstractRectifierFunction=Softplus(), basis::AbstractPolynomialBasis=LinearizedHermiteBasis(), map_type::Symbol=:total): Initialize polynomial map for map from density.

  • DiagonalMap(dimension::Int, degree::Int, referencetype::Symbol=:normal, rectifier::AbstractRectifierFunction=Softplus(), basis::AbstractPolynomialBasis=LinearizedHermiteBasis()): Initialize diagonal map with diagonal structure.

  • NoMixedMap(dimension::Int, degree::Int, referencetype::Symbol=:normal, rectifier::AbstractRectifierFunction=Softplus(), basis::AbstractPolynomialBasis=LinearizedHermiteBasis()): Initialize diagonal map without mixed terms.

source
TransportMaps.PolynomialMapComponent Type
julia
PolynomialMapComponent{T<:AbstractPolynomialBasis} <: AbstractMapComponent

A single component Mᵏ of a triangular transport map with polynomial basis.

Fields

  • basisfunctions::Vector{MultivariateBasis{T}}: Vector of multivariate basis functions

  • coefficients::Vector{Float64}: Coefficients for the basis functions

  • rectifier::AbstractRectifierFunction: Rectifier function applied to partial derivatives

  • index::Int: Index k indicating the component dimension

Constructors

  • PolynomialMapComponent(index::Int, degree::Int, rectifier::AbstractRectifierFunction=Softplus(), basis::AbstractPolynomialBasis=HermiteBasis(), map_type::Symbol=:total): Initialize component with specified degree and basis type.

  • PolynomialMapComponent(index::Int, degree::Int, rectifier::AbstractRectifierFunction, basis::AbstractPolynomialBasis, density::Distributions.UnivariateDistribution, map_type::Symbol=:total): Initialize component using analytical reference density.

  • PolynomialMapComponent(index::Int, degree::Int, rectifier::AbstractRectifierFunction, basis::AbstractPolynomialBasis, samples::AbstractMatrix{Float64}, map_type::Symbol=:total): Initialize component using sample data.

  • PolynomialMapComponent(multi_indices::Vector{Vector{Int}}, rectifier::AbstractRectifierFunction, basis::AbstractPolynomialBasis, samples::AbstractMatrix{Float64}): Initialize component with custom multi-index sets from samples.

  • PolynomialMapComponent(multi_indices::Vector{Vector{Int}}, rectifier::AbstractRectifierFunction, basis::AbstractPolynomialBasis, reference_density::Distributions.UnivariateDistribution): Initialize component with custom multi-index sets from density.

source
TransportMaps.LinearMap Type
julia
LinearMap <: AbstractLinearMap

A linear transformation map that standardizes data using mean and standard deviation.

Fields

  • μ::Vector{Float64}: Mean vector for each dimension

  • σ::Vector{Float64}: Standard deviation vector for each dimension

Constructors

  • LinearMap(samples::Matrix{Float64}): Compute empirical mean and standard deviation from samples.

  • LinearMap(μ::Vector{Float64}, σ::Vector{Float64}): Construct linear map with explicit mean and standard deviation.

source
TransportMaps.LaplaceMap Type
julia
LaplaceMap <: AbstractLinearMap

A linear transport map based on the Laplace approximation of a target density.

The Laplace approximation assumes the target density is approximately Gaussian around its mode. The map is defined by a location parameter (mode) and a scale parameter (Cholesky factor of the covariance).

Fields

  • mode::Vector{Float64}: Mode or mean vector of the approximation

  • chol::Matrix{Float64}: Lower Cholesky factor L of the covariance matrix (Σ = L * L')

Constructors

  • LaplaceMap(samples::Matrix{Float64}): Constructs a LaplaceMap from sample data by

computing the empirical mean and covariance.

  • LaplaceMap(density::MapTargetDensity, x0::Vector{Float64}; optimizer::Optim.AbstractOptimizer = LBFGS(), options::Optim.Options = Optim.Options()): Compute a Laplace approximation of a target density by finding the mode via optimization and computing the Hessian at the mode. If hessian_backend is nothing, uses the same backend as the density.
source
TransportMaps.ComposedMap Type
julia
ComposedMap{T<:AbstractLinearMap} <: AbstractComposedMap

A composed transport map consisting of a linear map followed by a polynomial map.

The composition is defined as S(x) = M(L(x)) where L is the linear map of type T and M is the polynomial map. The linear map can be a LinearMap or LaplaceMap.

Fields

  • linearmap<:T: The linear map component

  • polynomialmap::PolynomialMap: The polynomial map component

Constructors

  • ComposedMap(lm::AbstractLinearMap, pm::PolynomialMap)
source

Functions

TransportMaps.evaluate Function
julia
evaluate(mvb::MultivariateBasis{T}, z::Vector{<:Real}) where T<:AbstractPolynomialBasis

Evaluate the multivariate basis function at point z.

source
julia
evaluate(map_component::PolynomialMapComponent, z::AbstractVector{<:Real})

Evaluate the polynomial map component Mᵏ(z) at point z.

Computes Mᵏ(z) = f(z₁,...,zₖ₋₁,0) + ∫₀^{zₖ} g(∂f/∂xₖ) dxₖ, where g is the rectifier.

source
julia
evaluate(map_component::PolynomialMapComponent, Z::AbstractMatrix{<:Real})

Evaluate the map component at multiple points (row-wise) using multithreading.

Returns a vector where element i is Mᵏ(Z[i,:]).

source
julia
evaluate(M::PolynomialMap, z::AbstractVector{<:Real})

Evaluate the polynomial map at point z.

Returns a vector where each component i is Mⁱ(z₁,...,zᵢ).

source
julia
evaluate(M::PolynomialMap, Z::AbstractMatrix{<:Real})

Evaluate the polynomial map at multiple points (row-wise) using multithreading.

Each row of Z is a point to evaluate. Returns a matrix where row i contains M(Z[i,:]).

source
julia
evaluate(L::LinearMap, x::AbstractVector{<:Real})

Apply the linear transformation (x - μ) / σ to standardize the input.

source
julia
evaluate(L::LinearMap, X::AbstractMatrix{<:Real})

Apply the linear transformation to multiple points (row-wise).

source
julia
evaluate(L::LaplaceMap, x::AbstractVector{<:Real})

Apply the Laplace map transformation: L⁻¹(x - mode).

source
julia
evaluate(L::LaplaceMap, X::AbstractMatrix{<:Real})

Apply the Laplace map transformation to multiple points (row-wise).

source
julia
evaluate(C::ComposedMap, x::AbstractVector{<:Real})

Evaluate the composed map: S(x) = M(L(x)).

source
julia
evaluate(C::ComposedMap, X::AbstractMatrix{<:Real})

Evaluate the composed map for multiple points (row-wise).

source
TransportMaps.inverse Function
julia
inverse(map_component::PolynomialMapComponent, zₖ₋₁::AbstractVector{<:Real}, xₖ::Real)

Compute the inverse of the map component using one-dimensional root finding.

Given z[1:k-1] and the target value xₖ, finds zₖ such that Mᵏ(z) = xₖ.

source
julia
inverse(M::PolynomialMap, x::AbstractVector{<:Real}, k::Int=numberdimensions(M))

Compute the inverse of the first k components of the map at point x.

Returns z such that M(z)[1:k] = x[1:k].

source
julia
inverse(M::PolynomialMap, X::AbstractMatrix{<:Real}, k::Int=numberdimensions(M))

Compute the inverse at multiple points using multithreading.

Returns a matrix where row i contains M⁻¹(X[i,:])[1:k].

source
julia
inverse(L::LinearMap, y::AbstractVector{<:Real})

Invert the linear transformation: y * σ + μ to recover the original scale.

source
julia
inverse(L::LinearMap, Y::AbstractMatrix{<:Real})

Invert the transformation for multiple points (row-wise).

source
julia
inverse(L::LaplaceMap, y::AbstractVector{<:Real})

Invert the Laplace map: L * y + mode.

source
julia
inverse(L::LaplaceMap, Y::AbstractMatrix{<:Real})

Invert the transformation for multiple points (row-wise).

source
julia
inverse(C::ComposedMap, z::AbstractVector{<:Real})

Invert the composed map: S⁻¹(z) = L⁻¹(M⁻¹(z)).

source
julia
inverse(C::ComposedMap, Z::AbstractMatrix{<:Real})

Invert the composed map for multiple points (row-wise).

source
TransportMaps.pushforward Function
julia
pushforward(M::PolynomialMap, target::MapTargetDensity, z::AbstractVector{<:Real})

Compute the pushforward density ρ(z) = π(M(z)) |det J_M(z)|.

Maps the target density back through M to the reference space.

source
julia
pushforward(M::PolynomialMap, target::MapTargetDensity, Z::AbstractMatrix{<:Real})

Compute pushforward densities at multiple points using multithreading.

Returns a vector where element i is the pushforward density at Z[i,:].

source
TransportMaps.pullback Function
julia
pullback(M::PolynomialMap, x::AbstractVector{<:Real})

Compute the pullback density π̂(x) = ρ(M⁻¹(x)) |det J_{M⁻¹}(x)|.

Maps the reference density through M to approximate the target density.

source
julia
pullback(M::PolynomialMap, X::AbstractMatrix{<:Real})

Compute pullback densities at multiple points using multithreading.

Returns a vector where element i is the pullback density at X[i,:].

source
julia
pullback(C::ComposedMap, x::AbstractVector{<:Real})

Compute the pullback density: π(S(x)) * |det(∇S(x))|.

source
julia
pullback(C::ComposedMap, X::AbstractMatrix{<:Real})

Compute the pullback density for multiple points (row-wise).

source
TransportMaps.setcoefficients! Function
julia
setcoefficients!(map_component::PolynomialMapComponent, coefficients::AbstractVector{<:Real})

Set the coefficients of the map component.

source
julia
setcoefficients!(M::PolynomialMap, coefficients::AbstractVector{<:Real})

Set all coefficients of the polynomial map.

Coefficients are ordered by component: [c₁₁, c₁₂, ..., c₂₁, c₂₂, ..., cᵈₙ].

source
TransportMaps.getcoefficients Function
julia
getcoefficients(map_component::PolynomialMapComponent)

Get a copy of the coefficients from the map component.

source
julia
getcoefficients(M::PolynomialMap)

Get all coefficients from the polynomial map.

Returns a vector with coefficients ordered by component.

source
TransportMaps.setparameters! Function
julia
setparameters!(map::LinearMap, μ::AbstractVector{<:Real}, σ::AbstractVector{<:Real})

Set the mean and standard deviation parameters of the linear map.

source
TransportMaps.numbercoefficients Function
julia
numbercoefficients(map_component::PolynomialMapComponent)

Return the number of coefficients in the map component.

source
julia
numbercoefficients(M::PolynomialMap)

Return the total number of coefficients in the polynomial map.

source
TransportMaps.numberdimensions Function
julia
numberdimensions(M::PolynomialMap)

Return the number of dimensions (components) of the polynomial map.

source
julia
numberdimensions(L::LinearMap)

Return the number of dimensions of the linear map.

source
julia
numberdimensions(L::LaplaceMap)

Return the number of dimensions of the Laplace map.

source
julia
numberdimensions(C::ComposedMap)

Return the number of dimensions of the composed map.

source
TransportMaps.getmultiindexsets Function
julia
getmultiindexsets(map_component::PolynomialMapComponent)

Return the multi-indices as a matrix (one per row).

source
TransportMaps.conditional_density Function
julia
conditional_density(M::PolynomialMap, x_range::Real, x_given::AbstractVector{<:Real})

Compute the conditional density π(xₖ | x₁, ..., xₖ₋₁) at x_range given x_given.

source
julia
conditional_density(M::PolynomialMap, x_range::AbstractVector{<:Real}, x_given::AbstractVector{<:Real})

Compute conditional densities at multiple x_range values using multithreading.

source
TransportMaps.conditional_sample Function
julia
conditional_sample(M::PolynomialMap, x_given::AbstractVector{<:Real}, z_range::Real)

Generate a sample from π(xₖ | x₁, ..., xₖ₋₁) by pushing forward z_range ~ ρ(z_range).

source
julia
conditional_sample(M::PolynomialMap, x_given::AbstractVector{<:Real}, z_range::AbstractVector{<:Real})

Generate multiple samples from π(xₖ | x₁, ..., xₖ₋₁) using multithreading.

source
TransportMaps.multivariate_conditional_density Function
julia
multivariate_conditional_density(M::PolynomialMap, x::AbstractVector{<:Real})

Compute the multivariate conditional density π(xⱼ, xⱼ₊₁, ..., xₖ | x₁, ..., xⱼ₋₁) as ∏ᵢ₌ⱼᵏ π(xᵢ | x₁, ..., xᵢ₋₁).

source
julia
multivariate_conditional_density(M::PolynomialMap, x_range::AbstractVector{<:Real}, x_given::AbstractVector{<:Real})

Compute the conditional density π(xⱼ, ..., xₖ | x₁, ..., xⱼ₋₁) where x_range = [xⱼ, ..., xₖ].

source
TransportMaps.multivariate_conditional_sample Function
julia
multivariate_conditional_sample(M::PolynomialMap, x_given::AbstractVector{<:Real}, z_range::AbstractVector{<:Real})

Generate samples from π(xⱼ, ..., xₖ | x₁, ..., xⱼ₋₁) by sequentially pushing forward z_range values.

source