Transition Path Theory
Transition Path Theory (TPT) provides a framework for analyzing rare transition events between metastable states in complex systems. This module implements TPT calculations for Langevin dynamics systems.
Theory Overview
TPT characterizes transition pathways between two metastable states A and B by computing:
- Forward and backward committor functions
- Reactive probability density
- Reactive current
- Transition rates
The calculations are performed on a triangulated mesh representing the system's state space.
CriticalTransitions.Langevin
— Typestruct Langevin{H, D, KE, T}
A struct representing a Langevin dynamical system with damping rate gamma
and temperature
beta``.
Fields
Hamiltonian::Any
: Function that computes the total energy (kinetic + potential) of the system.driftfree::Any
: Function representing the divergence-free part of the drift.kinetic::Any
: Function computing the kinetic energy of the system.gamma::Any
: Damping coefficient that determines the strength of friction.beta::Any
: Inverse temperature parameter (β = 1/kT) that sets the noise intensity.
Constructors
Langevin(Hamiltonian, driftfree, kinetic, gamma, beta)
defined at /home/runner/work/CriticalTransitions.jl/CriticalTransitions.jl/src/transition_path_theory/langevin.jl:14
.
Committor Functions
CriticalTransitions.committor
— Functioncommittor(
sys::Langevin,
mesh::CriticalTransitions.Mesh,
Aind,
Bind
) -> Vector{Float64}
Solve the committor equation for a two-dimensional Langevin system using finite elements.
Arguments
sys::Langevin
: The Langevin system containing kinetic energy, drift-free function, beta, and gamma.mesh::Mesh
: The mesh containing points and triangles.Aind::Vector{Int}
: Indices of mesh points corresponding to set A (Dirichlet boundary condition = 0).Bind::Vector{Int}
: Indices of mesh points corresponding to set B (Dirichlet boundary condition = 1).
Returns
- A vector of committor values of length
N
, where each entry corresponds to a mesh node.
Implementation Details
This function assembles a global matrix A
and right-hand-side vector b
for the finite element discretization of the committor equation in Langevin dynamics. The code imposes Dirichlet boundary conditions on specified nodes (Aind
set to 0, Bind
set to 1). It computes elementwise contributions with stima_Langevin
and stimavbdv
, applying exponential factors involving beta
and gamma
to incorporate potential and damping effects. Finally, it solves the resulting linear system for the committor values on the free (non-boundary) nodes.
Invariant Probability Density
CriticalTransitions.invariant_pdf
— Functioninvariant_pdf(
sys::Langevin,
mesh::CriticalTransitions.Mesh,
Amesh::CriticalTransitions.Mesh,
Bmesh::CriticalTransitions.Mesh
) -> Any
Compute the invariant probability density function (PDF) for a Langevin system over a given mesh.
Arguments
sys::Langevin
: The Langevin system containing the Hamiltonian and beta parameters.mesh::Mesh
: The main mesh containing points and triangles.Amesh::Mesh
: The mesh corresponding to region A.Bmesh::Mesh
: The mesh corresponding to region B.
Returns
- A scalar value representing the normalization constant
Z
of the invariant PDF.
Implementation Details
This function calculates the invariant PDF by integrating the exponential of the negative Hamiltonian over the areas of the triangles in the provided meshes. The normalization constant Z
is computed by summing the contributions from the main mesh, region A mesh, and region B mesh.
Reactive Current
reactive_current
Probability Calculations
CriticalTransitions.probability_reactive
— Functionprobability_reactive(
sys::Langevin,
mesh::CriticalTransitions.Mesh,
q,
qminus,
Z
) -> Any
Calculate the probability that a trajectory is reactive (transitions from A to B).
Arguments
sys::Langevin
: System with Hamiltonian and inverse temperature (beta)mesh::Mesh
: Mesh structure containing points and triangulationq
: Forward committor function valuesqminus
: Backward committor function valuesZ
: Partition function value
Returns
- Probability (float) of reactive trajectories normalized by partition function
CriticalTransitions.probability_last_A
— Functionprobability_last_A(
sys::Langevin,
mesh::CriticalTransitions.Mesh,
Ames::CriticalTransitions.Mesh,
qminus,
Z
) -> Any
Calculate the probability that the last visited metastable state was A.
Arguments
sys::Langevin
: System with Hamiltonian and inverse temperature (beta)mesh::Mesh
: Main mesh structure containing points and triangulationAmes::Mesh
: Mesh structure for region Aqminus
: Backward committor function valuesZ
: Partition function value
Returns
- Probability (float) that the system was last in state A, normalized by partition function
References
- W. E and E. Vanden-Eijnden, "Towards a Theory of Transition Paths", Journal of Statistical Physics, 2006