Skip to content

Optimization

Index

Types

TransportMaps.OptimizationHistory Type
julia
OptimizationHistory

A data structure to store the iteration history of adaptive transport map optimization from samples.

Fields

  • terms::Vector{Matrix{Int64}}: Multi-index sets at each iteration

  • train_objectives::Vector{Float64}: Normalized training objectives at each iteration

  • test_objectives::Vector{Float64}: Normalized test objectives at each iteration (empty if no test split)

  • gradients::Vector{Vector{Float64}}: Gradient values for candidate terms at each iteration

source
TransportMaps.OptimizationResult Type
julia
OptimizationResult

A data structure to store the optimization results for each component of a transport map when optimizing the map from samples.

Fields

  • train_objectives::Vector{Float64}: Train objective for each component

  • test_objectives::Vector{Float64}: Test objective for each component

  • optimization_results::Vector{Optim.MultivariateOptimizationResults}: Result from optimization with Optim.jl

source
TransportMaps.MapOptimizationResult Type
julia
MapOptimizationResult

A data struct to store the iteration history of adaptive transport map optimization from density.

Fields

  • maps::Vector{PolynomialMap}: Maps for each iteration

  • train_objectives::Vector{Float64}: Train objectives for each iteration

  • test_objectives::Vector{Float64}: Test objectives for each iteration

  • gradients::Vector{Vector{Float64}}: Gradients of train objectives

  • optimization_results::Vector{Optim.MultivariateOptimizationResults}: Results of optimization with Optim.jl

source

Functions

TransportMaps.optimize! Function
julia
optimize!(M::PolynomialMap, target::AbstractMapDensity, quadrature::AbstractQuadratureWeights;
          optimizer, options)

Optimize polynomial map coefficients to minimize KL divergence to a target density.

Arguments

  • M::PolynomialMap: The polynomial map to optimize.

  • target::AbstractMapDensity: Target map density object (provides the target density π(x) and any needed operations).

  • quadrature::AbstractQuadratureWeights: Quadrature points and weights used for numerical integration.

Keyword Arguments

  • optimizer: Optimizer from Optim.jl to use (default: LBFGS()).

  • options: Options passed to the optimizer (default: Optim.Options()).

Returns

  • Optimization result from Optim.jl. The optimized coefficients are written back into M.
source
julia
optimize!(M::PolynomialMap, samples::Matrix{Float64}, lm::AbstractLinearMap=LinearMap();
          optimizer::Optim.AbstractOptimizer = LBFGS(),
          options::Optim.Options = Optim.Options(),
          test_fraction::Float64 = 0.0)

Optimize polynomial map coefficients to minimize KL divergence to a target density.

Arguments

  • M::PolynomialMap: The polynomial map to optimize.

  • samples::Matrix{Float64}: A matrix of sample data used to initialize and fit the map. Columns are interpreted as components/dimensions and rows as individual sample points.

  • lm::AbstractLinearMap: A linear map used to standardize the samples before optimization (default: LinearMap(), identity map).

Keyword Arguments

  • optimizer::Optim.AbstractOptimizer: Optimizer from Optim.jl to use (default: LBFGS()).

  • options::Optim.Options: Options passed to the optimizer (default: Optim.Options()).

  • test_fraction::Float64: Fraction of samples to hold out for testing/validation (default: 0.0, i.e. no test split).

Returns

  • OptimizationResult: Optimization results containing training and test objectives for each component. The optimized coefficients are written back into M.
source
TransportMaps.optimize_adaptive_transportmap Function
julia
optimize_adaptive_transportmap(samples, maxterms; kwargs...)

Adaptively optimize a triangular transport map by greedily enriching the multi-index set.

Arguments

  • samples::Matrix{Float64}: Sample data where rows are samples and columns are dimensions

  • maxterms::Vector{Int64}: Maximum number of terms for each component

Keyword Arguments

  • lm::AbstractLinearMap=LinearMap(samples): Linear map to standardize samples

  • rectifier::AbstractRectifierFunction=Softplus(): Rectifier function to use

  • basis::AbstractPolynomialBasis=LinearizedHermiteBasis(): Polynomial basis

  • optimizer::Optim.AbstractOptimizer=LBFGS(): Optimization algorithm

  • options::Optim.Options=Optim.Options(): Optimizer options

  • test_fraction::Float64=0.0: Fraction of samples to use for testing (validation)

Returns

  • C::ComposedMap{LinearMap, PolynomialMap}: The optimized triangular transport map

  • iteration_histories::Vector{OptimizationHistory}: History of optimization for each component

source
julia
optimize_adaptive_transportmap(samples, maxterms, k_folds, lm, rectifier, basis; optimizer, options)

Adaptively optimize a triangular transport map using k-fold cross-validation to select the number of terms per component.

Arguments

  • samples::Matrix{Float64}: Sample data where rows are samples and columns are dimensions

  • maxterms::Vector{Int64}: Upper bound on the number of terms to consider for each component

  • k_folds::Int: Number of folds to use for cross-validation

  • lm::AbstractLinearMap: Linear map to standardize samples (default: LinearMap(samples))

  • rectifier::AbstractRectifierFunction: Rectifier function to use (default: Softplus())

  • basis::AbstractPolynomialBasis: Polynomial basis (default: LinearizedHermiteBasis())

Keyword Arguments

  • optimizer::Optim.AbstractOptimizer=LBFGS(): Optimization algorithm

  • options::Optim.Options=Optim.Options(): Optimizer options

Returns

  • ComposedMap{LinearMap, PolynomialMap}: The optimized triangular transport map trained on all samples

  • fold_histories::Vector{Vector{OptimizationHistory}}: Optimization history for each component and fold

  • selected_terms::Vector{Int}: Number of terms selected for each component

  • selected_fold::Vector{Int}: Index of the best fold for each component

source
julia
optimize_adaptive_transportmap(
    target::AbstractMapDensity,
    quadrature::AbstractQuadratureWeights,
    maxterms::Int;
    kwargs...
)

Adaptively optimize a triangular transport map from a target density by greedily enriching the multi-index set across all components simultaneously.

Arguments

  • target::AbstractMapDensity: Target density to approximate

  • quadrature::AbstractQuadratureWeights: Quadrature points and weights for integration

  • maxterms::Int: Maximum total number of terms to add across all components

Keyword Arguments

  • rectifier::AbstractRectifierFunction=Softplus(): Rectifier function to use

  • basis::AbstractPolynomialBasis=LinearizedHermiteBasis(): Polynomial basis

  • optimizer::Optim.AbstractOptimizer=LBFGS(): Optimization algorithm

  • options::Optim.Options=Optim.Options(): Optimizer options

  • validation_samples::Matrix{Float64}=Matrix{Float64}(undef,0,0): Samples for variance diagnostic validation

Returns

  • M::PolynomialMap: The optimized triangular transport map (with best validation variance diagnostic)

  • history::OptimizationHistory: History of optimization iterations

source
TransportMaps.variance_diagnostic Function
julia
variance_diagnostic(M::PolynomialMap, target::MapTargetDensity, Z::AbstractArray{<:Real})

Compute a variance-based diagnostic for assessing the quality of a transport map.

The diagnostic measures the variance of the log-ratio between the pushforward density and the reference density. A smaller variance indicates a better approximation of the target density by the transport map.

Arguments

  • M::PolynomialMap: The transport map to be evaluated

  • target::MapTargetDensity: The target density that the map should approximate

  • Z::AbstractArray{<:Real}: Sample points from the reference distribution, where each row is a sample and columns correspond to dimensions

Returns

  • Float64: The computed variance diagnostic
source