Skip to content

Transport Maps#300

Merged
AnderGray merged 25 commits intomasterfrom
transport-maps
Mar 31, 2026
Merged

Transport Maps#300
AnderGray merged 25 commits intomasterfrom
transport-maps

Conversation

@lukasfritsch
Copy link
Copy Markdown
Member

This PR adds triangular transport maps (TMs) using TransportMaps.jl as the backend.

Transport Maps are added in inputs as a subtype of MultivariateDistribution, making them compatible with the JointDistribution interface. There are two different types added:

  • TransportMap which are constructed as a mapping from a target density in the physical space to a standard normal density
  • TransportMapFromSamples which are constructed starting with samples of the target density and map those to the standard normal space.

Also, TMs are used for model updating with the bayesianupdating function.
There, TransportMapBayesian wraps all the different objects and options of the map construction which is passed to the bayesianupdating function with a similar interface as for the sampling methods. The function returns the optimized map as a JointDistribution.

This PR also adds extensive documentation, tests and a literate example for Bayesian updating with transport maps.

@codecov
Copy link
Copy Markdown

codecov Bot commented Mar 10, 2026

Codecov Report

❌ Patch coverage is 95.80153% with 11 lines in your changes missing coverage. Please review.

Files with missing lines Patch % Lines
src/modelupdating/bayesianTM.jl 91.76% 7 Missing ⚠️
src/inputs/transportmaps.jl 97.60% 4 Missing ⚠️

📢 Thoughts on this report? Let us know!

Copy link
Copy Markdown
Contributor

@AnderGray AnderGray left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Very extensive and very well written, with great documentation and tests. This is awesome.

I may have found a couple of things:

  • Is the AutoMooncake() option working as intended? When I try it in your example, the code errors.
  • I don't know if it's only because it's in development mode, but some of docs/transportmaps hasn't rendered properly, in particular I'm seeing some unrendered math, and some loose bullet points.
  • You compute the mode of the transport map is the mapping of the zero in standard normal space. Is that always correct? It may be a good approximation in many cases, but I am finding examples where it is not the actual mode. Below is a simple example.
using UncertaintyQuantification
using Distributions
using ForwardDiff
using Plots

# Smooth target on all of R: skewed Gaussian mixture
dist = MixtureModel(
    [Normal(-1.0, 0.7), Normal(1.2, 0.9)],
    [0.75, 0.25],
)

# TransportMaps expects vector input even in 1D
logq(x::AbstractVector) = logpdf(dist, x[1])
dlogq(x::AbstractVector) = ForwardDiff.gradient(y -> logpdf(dist, y[1]), x)
target = MapTargetDensity(logq, dlogq)

# Learn a 1D transport map
M = PolynomialMap(1, 8, Normal(), Softplus())
quadrature = GaussHermiteWeights(35, 1)

tm_opt = mapfromdensity(M, target, quadrature, [:x1])
x_mode_tm = UncertaintyQuantification.mode(tm_opt.d)

# ------------------------------------------------------------
# Plot TM density and mode

xs = range(-4.0, 4.0, length=1000)
tm_pdf = [pdf(tm_opt, [x]) for x in xs]

plot(xs, tm_pdf, label="Transport map", lw=2)
vline!([x_mode_tm], label="TM mode", linestyle=:dashdot)
Image

In this example, the reported mode does not appear to coincide with the mode of the learned transport-map density.

@lukasfritsch
Copy link
Copy Markdown
Member Author

Thanks for the review @AnderGray!

  • I think the issue of AutoMooncake() not working is actually an issue of the transformation between the Beta and Normal distributions. By default, the prior is transformed to be standard normal such that the reference and target density are defined on the same support (-Inf to Inf). The mapping with auto-diff works well with a different prior, e.g. Uniform prior, however, it takes a while to setup the automatic differentiation.
  • You are right, I fixed that.
  • After some talk with @jgrashorn, we figured out that mapping zeros from the standard normal space through the map should give the median and not the mode. I have adjusted it, by renaming the function to median.

@lukasfritsch
Copy link
Copy Markdown
Member Author

This is the example again with median instead of mode:

using UncertaintyQuantification
using Distributions
using ForwardDiff
using Plots

# Smooth target on all of R: skewed Gaussian mixture
dist = MixtureModel(
    [Normal(-1.0, 0.7), Normal(1.2, 0.9)],
    [0.75, 0.25],
)

# TransportMaps expects vector input even in 1D
logq(x::AbstractVector) = logpdf(dist, x[1])
dlogq(x::AbstractVector) = ForwardDiff.gradient(y -> logpdf(dist, y[1]), x)
target = MapTargetDensity(logq, dlogq)

# Learn a 1D transport map
M = PolynomialMap(1, 8, Normal(), Softplus())
quadrature = GaussHermiteWeights(35, 1)

tm_opt = mapfromdensity(M, target, quadrature, [:x1])
x_mean_tm = median(tm_opt.d)
x_med_ana = median(dist)

# ------------------------------------------------------------
# Plot TM density and mode

xs = range(-4.0, 4.0, length=1000)
tm_pdf = [pdf(tm_opt, [x]) for x in xs]

mix_pdf = [pdf(dist, x) for x in xs]

plot(xs, tm_pdf, label="Transport map", lw=2)
plot!(xs, mix_pdf, label="Traget pdf", lw=2, ls=:dash)
vline!([x_mean_tm], label="Median TM", linestyle=:dashdot)
vline!([x_med_ana], label="Median", linestyle=:dashdot)

We obtain this plot:
medians

@lukasfritsch lukasfritsch requested a review from AnderGray March 31, 2026 12:51
Copy link
Copy Markdown
Contributor

@AnderGray AnderGray left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Great stuff ! Thanks for your changes

@AnderGray AnderGray merged commit 1311993 into master Mar 31, 2026
9 checks passed
@AnderGray AnderGray deleted the transport-maps branch March 31, 2026 14:04
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

feature New feature or feature request

Projects

None yet

Development

Successfully merging this pull request may close these issues.

2 participants