| Title: | Variational Flow-Based Inference for Rare Events and Large Deviations |
|---|---|
| Description: | Variational flow-based methods for modeling rare events using Kullback–Leibler (KL) divergence, normalizing flows, Girsanov change of measure, and Freidlin–Wentzell action functionals. The package provides tools for rare-event inference, minimum-action paths, and quasi-potential computation in stochastic dynamical systems. Methods are based on Rezende and Mohamed (2015) <doi:10.48550/arXiv.1505.05770>, Girsanov (1960) <doi:10.1137/1105027>, and Freidlin and Wentzell (2012, ISBN:978-0387955477). |
| Authors: | Pietro Piu [aut, cre] |
| Maintainer: | Pietro Piu <[email protected]> |
| License: | GPL (>= 3) |
| Version: | 0.1.0 |
| Built: | 2026-05-16 08:10:14 UTC |
| Source: | https://github.com/pietropiu-labstats/rareflow |
Computes the ELBO:
ELBOflow(flow, Qobs, pxgivenz, nmc = 256)ELBOflow(flow, Qobs, pxgivenz, nmc = 256)
flow |
A flow model created by |
Qobs |
Observed empirical distribution (probability vector). |
pxgivenz |
A function mapping a latent vector z to a categorical pmf. |
nmc |
Number of Monte Carlo samples. |
where:
q(z) is the flow-based variational posterior
p(z) is a standard Gaussian prior
p(x | z) is provided by the user as pxgivenz
The expectation is approximated via Monte Carlo sampling.
A numeric ELBO value.
f <- makeflow("planar", list(u = 0.1, w = 0.2, b = 0)) px <- function(z) c(0.3, 0.4, 0.3) ELBOflow(f, Qobs = c(0.2, 0.5, 0.3), pxgivenz = px, nmc = 100)f <- makeflow("planar", list(u = 0.1, w = 0.2, b = 0)) px <- function(z) c(0.3, 0.4, 0.3) ELBOflow(f, Qobs = c(0.2, 0.5, 0.3), pxgivenz = px, nmc = 100)
Computes the Freidlin-Wentzell quasi-potential between and ,
constructs a tilted likelihood proportional to ,
and fits a flow-based variational posterior.
fitflow_FW( observed, states = NULL, flowtype = "maf", flowspec = list(), inittheta = NULL, drift, x0, x1, T = 200, dt = 0.01, eps = 0.1, nmc = 256, control = list() )fitflow_FW( observed, states = NULL, flowtype = "maf", flowspec = list(), inittheta = NULL, drift, x0, x1, T = 200, dt = 0.01, eps = 0.1, nmc = 256, control = list() )
observed |
Empirical distribution Q. |
states |
Optional category names. |
flowtype |
Flow type. |
flowspec |
Structural parameters for the flow. |
inittheta |
Optional initial theta. |
drift |
Drift function |
x0 |
Starting point. |
x1 |
Target point. |
T |
Number of time steps. |
dt |
Time step. |
eps |
Noise strength (small parameter). |
nmc |
Monte Carlo samples. |
control |
Control list for |
This is useful for rare-event inference in small-noise diffusions, where the quasi-potential acts as an effective energy landscape.
Output of fitflowvariational().
Applies a Girsanov change of measure to tilt the likelihood and then
fits a flow-based variational posterior using fitflowvariational().
fitflow_girsanov( observed, states = NULL, flowtype = "maf", flowspec = list(), inittheta = NULL, base_pxgivenz, theta_path, Winc, dt, nmc = 256, control = list() )fitflow_girsanov( observed, states = NULL, flowtype = "maf", flowspec = list(), inittheta = NULL, base_pxgivenz, theta_path, Winc, dt, nmc = 256, control = list() )
observed |
Empirical distribution Q (probability vector). |
states |
Optional category names. |
flowtype |
Flow type ("maf", "splinepwlin", "planar", "radial"). |
flowspec |
Structural parameters for the flow. |
inittheta |
Optional initial theta for trainable flows. |
base_pxgivenz |
Likelihood |
theta_path |
Drift-tilting function or vector for Girsanov. |
Winc |
Brownian increments. |
dt |
Time step. |
nmc |
Monte Carlo samples. |
control |
Control list for |
This is useful when the target distribution arises from a drift-tilted diffusion process, where the Radon-Nikodym derivative is given by the Girsanov theorem.
Output of fitflowvariational().
This function fits a variational posterior q(z) using a chosen normalizing flow. The objective is the ELBO:
fitflowvariational( observed, states = NULL, flowtype = c("maf", "splinepwlin", "planar", "radial"), flowspec = list(), inittheta = NULL, pxgivenz, nmc = 256, control = list() )fitflowvariational( observed, states = NULL, flowtype = c("maf", "splinepwlin", "planar", "radial"), flowspec = list(), inittheta = NULL, pxgivenz, nmc = 256, control = list() )
observed |
Observed empirical distribution Q (probability vector). |
states |
Optional vector of category names. |
flowtype |
One of "maf", "splinepwlin", "planar", "radial". |
flowspec |
A list specifying structural parameters (d, K, etc.). |
inittheta |
Optional initial parameter vector for trainable flows. |
pxgivenz |
A function mapping latent z to a categorical pmf. |
nmc |
Number of Monte Carlo samples for ELBO estimation. |
control |
List of control parameters passed to |
The flow parameters (theta) are optimized via optim() when applicable
(MAF and spline flows). Planar and radial flows have no trainable parameters.
This function performs generic variational inference using a chosen
normalizing flow and a user-provided likelihood pxgivenz.
For specialized rare-event inference using:
Girsanov change of measure
Freidlin–Wentzell quasi-potential
see the wrapper functions:
fitflow_girsanov()
fitflow_FW()
These wrappers construct a tilted likelihood and then call
fitflowvariational() internally.
A list containing:
flow: the fitted flow model
elbo: final ELBO value
theta: optimized parameter vector (if applicable)
convergence: optim() convergence code
Computes the discrete Freidlin-Wentzell action for a path
represented as a matrix of size . The continuous action is:
Freidlin_Wentzell_action(phi, drift, dt)Freidlin_Wentzell_action(phi, drift, dt)
phi |
Matrix of path values of dimension |
drift |
Drift function |
dt |
Time step. |
and the discrete approximation is:
Numeric action value.
Computes an approximate Freidlin-Wentzell quasi-potential between two
points and by minimizing the FW action functional
over discretized paths.
FW_quasipotential( x0, x1, drift, T = 200, dt = 0.01, niter = 200, stepsize = 0.1 )FW_quasipotential( x0, x1, drift, T = 200, dt = 0.01, niter = 200, stepsize = 0.1 )
x0 |
Starting point (numeric vector). |
x1 |
Target point (numeric vector). |
drift |
Drift function |
T |
Number of time steps. |
dt |
Time step. |
niter |
Number of gradient descent iterations. |
stepsize |
Gradient descent step size. |
The algorithm:
Initializes a straight-line path between and .
Performs simple gradient descent on the FW action.
This is a naive but effective illustrative method for low-dimensional systems. More advanced solvers (string method, MAM, etc.) can be plugged in.
A list with:
path: matrix of size
action: FW action of the optimized path
Computes the Radon-Nikodym derivative (log form) associated with a Girsanov change of measure for an SDE:
girsanov_logratio(theta_path, Winc, dt)girsanov_logratio(theta_path, Winc, dt)
theta_path |
Numeric vector of drift tilts |
Winc |
Numeric vector of Brownian increments |
dt |
Time step size. |
tilted by an alternative drift:
.
The log-likelihood ratio is:
This function returns the log-ratio for a given path of drift tilts
theta_path, Brownian increments Winc, and time step dt.
A numeric log-likelihood ratio.
Computes the KL divergence D(Q || P) between two discrete probability distributions Q and P.
KLdiv(Q, P)KLdiv(Q, P)
Q |
Observed distribution (probability vector). |
P |
Baseline or reference distribution (probability vector). |
A numeric KL divergence value.
A readable and structured implementation of a Masked Autoregressive Flow. This version supports:
arbitrary dimension d
K sequential flow steps
a single parameter vector theta containing all weights
mafflowmodel(d = 3, K = 2, theta = NULL)mafflowmodel(d = 3, K = 2, theta = NULL)
d |
Dimension of the latent space. |
K |
Number of flow steps. |
theta |
Optional parameter vector. If NULL, random initialization. |
A flow model object with methods:
sampleq(n)
logq(z0)
applyflow(z0)
Same structure as makeneurolik(), but with a different default separation
parameter. This can be used to model simple biological switching systems
or coarse-grained gene expression states.
makebiolik(a = 0.2)makebiolik(a = 0.2)
a |
Separation parameter controlling the spacing between the two logits. |
A function mapping a latent vector z to a probability vector.
A unified constructor for all flow models in the rareflow package.
This function dispatches to the appropriate flow implementation
based on the flowtype argument.
makeflow(flowtype, params = list())makeflow(flowtype, params = list())
flowtype |
Character string specifying the flow type. |
params |
A named list of parameters required by the chosen flow. |
Supported flow types:
"planar" -> planarflowmodel()
"radial" -> radialflowmodel()
"maf" -> mafflowmodel()
"splinepwlin" -> splinepwlinflowmodel()
A flow model object with methods:
sampleq(n)
logq(z0)
applyflow(z0)
# Create a planar flow f <- makeflow("planar", list(u = 0.1, w = 0.2, b = 0)) s <- f$sampleq(10) # Create a 2D spline flow f2 <- makeflow("splinepwlin", list(d = 2, K = 8))# Create a planar flow f <- makeflow("planar", list(u = 0.1, w = 0.2, b = 0)) s <- f$sampleq(10) # Create a 2D spline flow f2 <- makeflow("splinepwlin", list(d = 2, K = 8))
Constructs a simple 3-category likelihood model based on a latent vector z.
The likelihood is defined by two logistic separators:
makeneurolik(a = 0.3)makeneurolik(a = 0.3)
a |
Separation parameter controlling the spacing between the two logits. |
p1 = sigmoid(mean(z) - a) p2 = sigmoid(mean(z) + a)
producing a 3-class probability vector:
(1 - p1, p1 - p2, p2)
This likelihood is useful for toy neural classification models or simple latent-to-categorical mappings.
A function mapping a latent vector z to a probability vector.
A simple and readable implementation of a 1-dimensional planar flow: z_K = z_0 + u * h(wz_0 + b)
planarflowmodel(u, w, b)planarflowmodel(u, w, b)
u |
Scalar parameter controlling the magnitude of the deformation. |
w |
Scalar parameter controlling the slope of the activation. |
b |
Scalar bias term. |
where:
h is a smooth activation (tanh)
the log-determinant is computed analytically
This flow is mainly useful for pedagogical purposes or as a lightweight building block in variational inference.
A flow model object with methods:
sampleq(n)
logq(z0)
applyflow(z0)
A readable implementation of a 1-dimensional radial flow:
radialflowmodel(z_ref, alpha, beta)radialflowmodel(z_ref, alpha, beta)
z_ref |
Reference point for the radial transformation. |
alpha |
Positive scalar controlling the denominator. |
beta |
Scalar controlling the magnitude of the deformation. |
z_K = z_0 + beta / (alpha + |z_0 - z_ref|) * (z_0 - z_ref)
where:
z_ref is a reference point
alpha > 0 ensures numerical stability
beta controls the strength of the radial deformation
The log-determinant is computed analytically for the 1D case.
A flow model object with methods:
sampleq(n)
logq(z0)
applyflow(z0)
Computes the classical Sanov upper bound:
sanovprob(Q, P, n)sanovprob(Q, P, n)
Q |
Observed empirical distribution. |
P |
True distribution. |
n |
Sample size. |
where Q is the empirical distribution and P is the true distribution.
A numeric upper bound.
A readable implementation of a monotone piecewise-linear spline flow. Each dimension is transformed independently using learned spline parameters.
splinepwlinflowmodel(d = 2, K = 8, theta = NULL)splinepwlinflowmodel(d = 2, K = 8, theta = NULL)
d |
Dimension of the latent space. |
K |
Number of spline bins. |
theta |
Optional parameter vector. If NULL, random initialization. |
The spline flow uses:
bins with learned widths and heights
a softmax transformation to ensure positivity and normalization
a sigmoid reparameterization for numerical stability
The flow is invertible and differentiable almost everywhere.
A flow model object with methods:
sampleq(n)
logq(z0)
applyflow(z0)