test_stochasticIntegrators
Tests the Euler-Mayurama and SRK integrator described in:
Tang, X., Xiao, A. Efficient weak second-order stochastic Runge-Kutta methods for Itô stochastic differential equations. Bit Numer Math 57, 241-260 (2017). https://doi.org/10.1007/s10543-016-0618-9
- class test_stochasticIntegrators.ComplexOrnsteinUhlenbeckSystem(x0: float = 0.1, y0: float = -0.1, theta_x: float = 0.1, theta_y: float = 0.073, sigma_x1: float = 0.015, sigma_x2: float = 0.011, sigma_y1: float = 0, sigma_y2: float = 0.029)[source]
Bases:
objectA process defined by
dx = -theta*x*dt + sigma_x1*dW_1 + sigma_x2*dW_2 dy = -theta*y*dt + sigma_y1*dW_1 + sigma_y2*dW_2
- f(t: float, x: ndarray[tuple[Any, ...], dtype[float64]])[source]
Drift function for the complex OU process.
- property g
Return list of diffusion functions.
- class test_stochasticIntegrators.Example1System(y1_0: float = 1, y2_0: float = 1)[source]
Bases:
objectExample 1 dynamical system in:
Tang, X., Xiao, A. Efficient weak second-order stochastic Runge-Kutta methods for Itô stochastic differential equations. Bit Numer Math 57, 241-260 (2017). https://doi.org/10.1007/s10543-016-0618-9
- f(t: float, x: ndarray[tuple[Any, ...], dtype[float64]])[source]
Drift function for Example 1 system.
- property g
Return list of diffusion functions.
- class test_stochasticIntegrators.ExponentialSystem(x0: float = 1)[source]
Bases:
objectA simple deterministic system with one state: dx/dt = x*t.
- class test_stochasticIntegrators.OrnsteinUhlenbeckSystem(x0: float = 0.1, mu: float = 0, theta: float = 0.1, sigma: float = 0.01)[source]
Bases:
objectA process defined by
dx = theta*(mu - x)*dt + sigma*dW
- f(t: float, x: ndarray[tuple[Any, ...], dtype[float64]])[source]
Drift function for the OU process.
- property g
Return list of diffusion functions.
- test_stochasticIntegrators.estimateErrorAndEmpiricalVariance(computeTrajectory: Callable[[], ndarray[tuple[Any, ...], dtype[float64]]], G: Callable[[ndarray[tuple[Any, ...], dtype[float64]]], float], estimateGOnTrajectory: float, M1: int, M2: int)[source]
Computes the error and empirical variance according to the equations described in Section 4 of Tang & Xiao.
- Parameters:
computeTrajectory (Callable[[], npt.NDArray[np.float64]]) – A function that,
called (when)
from (realizes one simulation by forward propagating the system)
paper (In the) –
y_N.G (Callable[[npt.NDArray[np.float64]], float]) – A differentiable function
number. (that takes in a state vector and outputs a single)
estimateGOnTrajectory (float) – The estimate of the application of the function
system. (G on the random variable that is the last state of the propagated)
paper –
E[G(y(t_N))]M1 (int) – Number of batches.
M2 (int) – Number of trajectories per batch.
- Returns:
The integrator error (
\hat{\mu}in the paper), and the empirical variance (\hat{\sigma}_\mu^2in the paper).- Return type:
tuple[float, float]
- test_stochasticIntegrators.eulerMayuramaIntegrate(f: Callable[[float, ndarray[tuple[Any, ...], dtype[float64]]], ndarray[tuple[Any, ...], dtype[float64]]], g_list: List[Callable[[float, ndarray[tuple[Any, ...], dtype[float64]]], ndarray[tuple[Any, ...], dtype[float64]]]], x0: ndarray[tuple[Any, ...], dtype[float64]], dt: float, tf: float, rng_seed: int)[source]
- Euler-Mayurama integrator for the vector SDE:
dX = f(t,X) dt + sum_k g_list[k](t,X) dW_k
- Parameters:
f – Drift function.
g_list – List of diffusion functions.
x0 – Initial state.
dt – Time step.
tf – Final time.
rng_seed – Random seed.
- Returns:
Array of state trajectories, including time as the first column.
- test_stochasticIntegrators.getBasiliskSim(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], dt: float, x0: ndarray[tuple[Any, ...], dtype[float64]], f: Callable[[float, ndarray[tuple[Any, ...], dtype[float64]]], ndarray[tuple[Any, ...], dtype[float64]]], g: List[Callable[[float, ndarray[tuple[Any, ...], dtype[float64]]], ndarray[tuple[Any, ...], dtype[float64]]]], seed: int | None)[source]
Set up and return a Basilisk simulation for a given SDE and integrator method.
- Parameters:
method – Integration method (“W2Ito1”, “W2Ito2”, or “EulerMayurama”).
dt – Time step.
x0 – Initial state.
f – Drift function.
g – List of diffusion functions.
seed – RNG seed (or None for random).
- Returns:
Tuple of (scSim, stateModel, integratorObject, stateLogger).
- test_stochasticIntegrators.srk2Integrate(f: Callable[[float, ndarray[tuple[Any, ...], dtype[float64]]], ndarray[tuple[Any, ...], dtype[float64]]], g_list: List[Callable[[float, ndarray[tuple[Any, ...], dtype[float64]]], ndarray[tuple[Any, ...], dtype[float64]]]], x0: ndarray[tuple[Any, ...], dtype[float64]], dt: float, tf: float, rng_seed: int, alpha: ndarray[tuple[Any, ...], dtype[float64]], beta0: ndarray[tuple[Any, ...], dtype[float64]], beta1: ndarray[tuple[Any, ...], dtype[float64]], A0: ndarray[tuple[Any, ...], dtype[float64]], B0: ndarray[tuple[Any, ...], dtype[float64]], A1: ndarray[tuple[Any, ...], dtype[float64]], B1: ndarray[tuple[Any, ...], dtype[float64]], B2: ndarray[tuple[Any, ...], dtype[float64]])[source]
- Generic s-stage SRK integrator for vector SDE:
dX = f(t,X) dt + sum_k g_list[k](t,X) dW_k
Method described in Tang & Xiao. :param f: Drift function. :param g_list: List of diffusion functions. :param x0: Initial state. :param dt: Time step. :param tf: Final time. :param rng_seed: Random seed. :param alpha: SRK coefficients. :param beta0: SRK coefficients. :param beta1: SRK coefficients. :param A0: SRK coefficients. :param B0: SRK coefficients. :param A1: SRK coefficients. :param B1: SRK coefficients. :param B2: SRK coefficients.
- Returns:
Array of state trajectories, including time as the first column.
- test_stochasticIntegrators.test_deterministic(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], plot: bool = False)[source]
Test deterministic integration (no diffusion) for all integrator methods. Compares Basilisk and Python implementations against the analytical solution for the exponential system dx/dt = x*t.
- Parameters:
method – Integration method (“W2Ito1”, “W2Ito2”, or “EulerMayurama”).
plot – If True, plot the relative error.
- test_stochasticIntegrators.test_example1(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], plot: bool = False)[source]
Test Example 1 system from Tang & Xiao (2017) for all integrator methods. Compares Basilisk and Python implementations for a single realization.
- Parameters:
method – Integration method.
plot – If True, plot the state trajectories.
- test_stochasticIntegrators.test_ou(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], plot: bool = False)[source]
Test Ornstein-Uhlenbeck process integration for all integrator methods. Compares Basilisk and Python implementations for a single realization.
- Parameters:
method – Integration method.
plot – If True, plot the state trajectories.
- test_stochasticIntegrators.test_ouComplex(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], plot: bool = False)[source]
Test integration of a two-dimensional coupled Ornstein-Uhlenbeck process. Compares Basilisk and Python implementations for a single realization.
- Parameters:
method – Integration method.
plot – If True, plot the state trajectories.
- test_stochasticIntegrators.test_validateExample1(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], tf: float = 0.1, skipAssert: bool = True)[source]
Validate the weak accuracy of the integrators for Example 1 from Tang & Xiao (2017). Compares the empirical variance of the final state to the analytical value using multiple Monte Carlo batches.
- Parameters:
method – Integration method.
- test_stochasticIntegrators.test_validateOu(method: Literal['W2Ito1', 'W2Ito2', 'EulerMayurama'], figureOfMerit: Literal['mean', 'variance'], tf: float = 0.1)[source]
Validate the weak accuracy of the integrators for the Ornstein-Uhlenbeck process. Compares the empirical mean or variance of the final state to the analytical value using multiple Monte Carlo batches.
- Parameters:
method – Integration method.
figureOfMerit – “mean” or “variance” to validate.