Skip to content

Transport Maps

Transport maps for density transformation and variational inference.

Index

Types

UncertaintyQuantification.TransportMap Type
julia
TransportMap(map, target, transform_density, names)

A transport map used to transform between standard normal space Z and physical space X. The map is the optimized triangular transport map, target is the target density, and transform_density optionally specifies random variables for density transformation, and names is a vector of variable names.

source
UncertaintyQuantification.TransportMapFromSamples Type
julia
TransportMapFromSamples(map, samples, names)

A transport map constructed from samples in the physical space X which are mapped to the standard normal space Z. The map is the optimized composed map, samples is a DataFrame of samples in the physical space X used to fit the map, and names is a vector of variable names.

source

Functions

UncertaintyQuantification.mapfromdensity Function
julia
mapfromdensity(transportmap, target, quadrature, names, transform_density, optimizer, options)

Optimize a transport map from the given target density by minimizing the KL-divergence. The transportmap is the transport map to be optimized, target is the target density in the physical space, and quadrature specifies quadrature points in the standard normal space. The KL-divergence is evaluated at the given quadrature points. The names is a vector of variable names, and transform_density optionally specifies random variables for transformation. The optimizer specifies the optimization method from Optim.jl (default: LBFGS()), and options allows passing options to the optimizer (default: Optim.Options()). Returns a JointDistribution with the optimized TransportMap.

Alternative calls

julia
mapfromdensity(transportmap, target, quadrature, names, transform_density)  # optimizer = LBFGS(), options = Optim.Options()
source
UncertaintyQuantification.mapfromsamples Function
julia
mapfromsamples(transportmap, X, optimizer, options)

Fit a transport map from samples. The transportmap is the transport map to be optimized, and X is a DataFrame with samples in the physical space. The optimizer specifies the optimization method from Optim.jl (default: LBFGS()), and options allows passing options to the optimizer (default: Optim.Options()). Returns a JointDistribution with the optimized TransportMapFromSamples.

Alternative calls

julia
mapfromsamples(transportmap, X)  # optimizer = LBFGS(), options = Optim.Options()
source
UncertaintyQuantification.to_physical_space! Method
julia
to_physical_space!(tm::TransportMap, Z::DataFrame)

Transforms samples in Z from standard normal space to physical space using transport map tm.

source
UncertaintyQuantification.to_physical_space! Method
julia
to_physical_space!(tm::TransportMapFromSamples, Z::DataFrame)

Transforms samples in Z from standard normal space to physical space using transport map tm.

source
UncertaintyQuantification.to_standard_normal_space! Method
julia
to_standard_normal_space!(tm::TransportMap, X::DataFrame)

Transforms samples in X from physical space to standard normal space using transport map tm.

source
UncertaintyQuantification.to_standard_normal_space! Method
julia
to_standard_normal_space!(tm::TransportMapFromSamples, X::DataFrame)

Transforms samples in X from physical space to standard normal space using transport map tm.

source
UncertaintyQuantification.sample Function
julia
sample(tm::AbstractTransportMap, n::Integer=1)

Generate n samples in the physical space X using the transport map tm.

source
Distributions.pdf Method
julia
pdf(tm::TransportMap, x::AbstractVector{<:Real})

Evaluate the probability density function (pdf) of the transport map in the physical space, i.e., the pushforward density. x is a vector, representing a point in the M-dimensional target space. Returns a Float64.

source
Distributions.pdf Method
julia
pdf(tm::TransportMap, x::AbstractMatrix{<:Real})

Evaluate the probability density function (pdf) of the transport map in the physical space, i.e., the pushforward density. X is a matrix of points in the physical space for which the pdf is evaluated. Returns a Vector{Float64} of evaluated pdf values for each row of the matrix.

source
Distributions.pdf Method
julia
pdf(tm::TransportMapFromSamples, x)

Evaluate the probability density function (pdf) of the transport map in the physical space, i.e., the pullback density. The argument x can be either a vector or matrix of values where the pdf is evaluated. Returns a Float64 for a vector input or Vector{Float64} for a matrix input.

source
Distributions.logpdf Method
julia
logpdf(tm::AbstractTransportMap, x)

Log-probability density function of the transport map in the physical space X. The argument x can be either a vector or matrix of values where the log-pdf is evaluated. Returns a Float64 for a vector input or Vector{Float64} for a matrix input.

source
Statistics.mean Function
julia
mean(tm::TransportMap, quad::AbstractQuadratureWeights=SparseSmolyakWeights)

Get the mean value of the density approximated by the transport map. The mean value is computed using numerically in the reference space. The scheme can be set using the optional argument quad which uses SparseSmolyakWeights by default.

source
Statistics.mean Function
julia
mean(tm::TransportMapFromSamples, quad::AbstractQuadratureWeights=SparseSmolyakWeights)

Get the mean value of the density approximated by the transport map. The mean value is computed using numerically in the reference space. The scheme can be set using the optional argument quad which uses SparseSmolyakWeights by default.

source
Statistics.var Function
julia
var(tm::TransportMap, quad::AbstractQuadratureWeights=SparseSmolyakWeights)

Get the variance of the density approximated by the transport map. The variance is computed numerically in the reference space. The scheme can be set using the optional argument quad which uses SparseSmolyakWeights by default.

source
Statistics.var Function
julia
var(tm::TransportMapFromSamples, quad::AbstractQuadratureWeights=SparseSmolyakWeights)

Get the variance of the density approximated by the transport map. The variance is computed numerically in the reference space. The scheme can be set using the optional argument quad which uses SparseSmolyakWeights by default.

source
Statistics.std Function
julia
std(tm::TransportMap, quad::AbstractQuadratureWeights=SparseSmolyakWeights)

Get the standard deviation of the density approximated by the transport map. The standard deviation is computed numerically in the reference space. The scheme can be set using the optional argument quad which uses SparseSmolyakWeights by default.

source
Statistics.std Function
julia
std(tm::TransportMapFromSamples, quad::AbstractQuadratureWeights=SparseSmolyakWeights)

Get the standard deviation of the density approximated by the transport map. The standard deviation is computed numerically in the reference space. The scheme can be set using the optional argument quad which uses SparseSmolyakWeights by default.

source
Statistics.median Method
julia
median(tm::TransportMap)

Get the median of the density approximated by the transport map.

source
Statistics.median Method
julia
median(tm::TransportMapFromSamples)

Get the median of the density approximated by the transport map.

source
UncertaintyQuantification.variancediagnostic Function
julia
variancediagnostic(tm::TransportMap, Z::DataFrame)

Evaluate the 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. The argument Z is a DataFrame of samples in the standard normal space. Returns the evaluated variance diagnostic.

source