Title: | Langevin Analysis in One and Two Dimensions |
---|---|
Description: | Estimate drift and diffusion functions from time series and generate synthetic time series from given drift and diffusion coefficients. |
Authors: | Philip Rinn [aut, cre], Pedro G. Lind [aut], David Bastine [ctb] |
Maintainer: | Philip Rinn <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.3.2 |
Built: | 2024-10-30 06:44:40 UTC |
Source: | CRAN |
Langevin1D
calculates the Drift and Diffusion vectors (with errors)
for a given time series.
Langevin1D( data, bins, steps, sf = ifelse(is.ts(data), frequency(data), 1), bin_min = 100, reqThreads = -1, kernel = FALSE, h )
Langevin1D( data, bins, steps, sf = ifelse(is.ts(data), frequency(data), 1), bin_min = 100, reqThreads = -1, kernel = FALSE, h )
data |
a vector containing the time series or a time-series object. |
bins |
a scalar denoting the number of |
steps |
a vector giving the |
sf |
a scalar denoting the sampling frequency (optional if |
bin_min |
a scalar denoting the minimal number of events per |
reqThreads |
a scalar denoting how many threads to use. Defaults to
|
kernel |
a logical denoting if the kernel based Nadaraya-Watson estimator should be used to calculate drift and diffusion vectors. |
h |
a scalar denoting the bandwidth of the data. Defaults to Scott's
variation of Silverman's rule of thumb. Only used if |
Langevin1D
returns a list with thirteen (six if kernel
is TRUE
) components:
D1 |
a vector of the Drift coefficient for each |
eD1 |
a vector of the error of the Drift coefficient for each
|
D2 |
a vector of the Diffusion coefficient for each |
eD2 |
a vector of the error of the Diffusion coefficient for
each |
D4 |
a vector of the fourth Kramers-Moyal coefficient for each
|
mean_bin |
a vector of the mean value per |
density |
a vector of the number of events per |
M1 |
a matrix of the first conditional moment for each
|
eM1 |
a matrix of the error of the first conditional moment
for each |
M2 |
a matrix of the second conditional moment for each
|
eM2 |
a matrix of the error of the second conditional moment
for each |
M4 |
a matrix of the fourth conditional moment for each
|
U |
a vector of the |
Philip Rinn
# Set number of bins, steps and the sampling frequency bins <- 20 steps <- c(1:5) sf <- 1000 #### Linear drift, constant diffusion # Generate a time series with linear D^1 = -x and constant D^2 = 1 x <- timeseries1D(N = 1e6, d11 = -1, d20 = 1, sf = sf) # Do the analysis est <- Langevin1D(data = x, bins = bins, steps = steps, sf = sf) # Plot the result and add the theoretical expectation as red line plot(est$mean_bin, est$D1) lines(est$mean_bin, -est$mean_bin, col = "red") plot(est$mean_bin, est$D2) abline(h = 1, col = "red") #### Cubic drift, constant diffusion # Generate a time series with cubic D^1 = x - x^3 and constant D^2 = 1 x <- timeseries1D(N = 1e6, d13 = -1, d11 = 1, d20 = 1, sf = sf) # Do the analysis est <- Langevin1D(data = x, bins = bins, steps = steps, sf = sf) # Plot the result and add the theoretical expectation as red line plot(est$mean_bin, est$D1) lines(est$mean_bin, est$mean_bin - est$mean_bin^3, col = "red") plot(est$mean_bin, est$D2) abline(h = 1, col = "red")
# Set number of bins, steps and the sampling frequency bins <- 20 steps <- c(1:5) sf <- 1000 #### Linear drift, constant diffusion # Generate a time series with linear D^1 = -x and constant D^2 = 1 x <- timeseries1D(N = 1e6, d11 = -1, d20 = 1, sf = sf) # Do the analysis est <- Langevin1D(data = x, bins = bins, steps = steps, sf = sf) # Plot the result and add the theoretical expectation as red line plot(est$mean_bin, est$D1) lines(est$mean_bin, -est$mean_bin, col = "red") plot(est$mean_bin, est$D2) abline(h = 1, col = "red") #### Cubic drift, constant diffusion # Generate a time series with cubic D^1 = x - x^3 and constant D^2 = 1 x <- timeseries1D(N = 1e6, d13 = -1, d11 = 1, d20 = 1, sf = sf) # Do the analysis est <- Langevin1D(data = x, bins = bins, steps = steps, sf = sf) # Plot the result and add the theoretical expectation as red line plot(est$mean_bin, est$D1) lines(est$mean_bin, est$mean_bin - est$mean_bin^3, col = "red") plot(est$mean_bin, est$D2) abline(h = 1, col = "red")
Langevin2D
calculates the Drift (with error) and Diffusion matrices
for given time series.
Langevin2D( data, bins, steps, sf = ifelse(is.mts(data), frequency(data), 1), bin_min = 100, reqThreads = -1 )
Langevin2D( data, bins, steps, sf = ifelse(is.mts(data), frequency(data), 1), bin_min = 100, reqThreads = -1 )
data |
a matrix containing the time series as columns or a time-series object. |
bins |
a scalar denoting the number of |
steps |
a vector giving the |
sf |
a scalar denoting the sampling frequency (optional if |
bin_min |
a scalar denoting the minimal number of events per |
reqThreads |
a scalar denoting how many threads to use. Defaults to
|
Langevin2D
returns a list with nine components:
D1 |
a tensor with all values of the drift coefficient.
Dimension is |
eD1 |
a tensor with all estimated errors of the drift
coefficient. Dimension is |
D2 |
a tensor with all values of the diffusion coefficient.
Dimension is |
mean_bin |
a matrix of the mean value per |
density |
a matrix of the number of events per |
M1 |
a tensor of the first moment for each |
eM1 |
a tensor of the standard deviation of the first
moment for each bin (line label) and each |
M2 |
a tensor of the second moment for each bin (line
label) and each |
U |
a matrix of the |
Philip Rinn
plot method for class "Langevin". This method is only implemented for one-dimensional analysis for now.
## S3 method for class 'Langevin' plot(x, pch = 20, ...)
## S3 method for class 'Langevin' plot(x, pch = 20, ...)
x |
an object of class "Langevin". |
pch |
Either an integer specifying a symbol or a single character to be
used as the default in plotting points. See points for possible values
and their interpretation. Default is |
... |
Arguments to be passed to methods, such as |
Philip Rinn
print method for class "Langevin".
## S3 method for class 'Langevin' print(x, digits = max(3, getOption("digits") - 3), ...)
## S3 method for class 'Langevin' print(x, digits = max(3, getOption("digits") - 3), ...)
x |
an object of class "Langevin". |
digits |
integer, used for number formatting with |
... |
further arguments to be passed to or from other methods. They are ignored in this function. |
The function print.Langevin()
returns an overview of the
estimated drift and diffusion coefficients.
Philip Rinn
summary method for class "Langevin".
## S3 method for class 'Langevin' summary(object, ..., digits = max(3, getOption("digits") - 3))
## S3 method for class 'Langevin' summary(object, ..., digits = max(3, getOption("digits") - 3))
object |
an object of class "Langevin". |
... |
arguments to be passed to or from other methods. They are ignored in this function. |
digits |
integer, used for number formatting with |
The function summary.Langevin()
returns a summary of the
estimated drift and diffusion coefficients
Philip Rinn
timeseries1D
generates a one-dimensional Langevin process using a
simple Euler integration. The drift function is a cubic polynomial, the
diffusion function a quadratic.
timeseries1D( N, startpoint = 0, d13 = 0, d12 = 0, d11 = -1, d10 = 0, d22 = 0, d21 = 0, d20 = 1, sf = 1000, dt = 0 )
timeseries1D( N, startpoint = 0, d13 = 0, d12 = 0, d11 = -1, d10 = 0, d22 = 0, d21 = 0, d20 = 1, sf = 1000, dt = 0 )
N |
a scalar denoting the length of the time-series to generate. |
startpoint |
a scalar denoting the starting point of the time series. |
d13 , d12 , d11 , d10
|
scalars denoting the coefficients for the drift polynomial. |
d22 , d21 , d20
|
scalars denoting the coefficients for the diffusion polynomial. |
sf |
a scalar denoting the sampling frequency. |
dt |
a scalar denoting the maximal time step of integration. Default
|
timeseries1D
returns a time-series object of length
N
with the generated time-series.
Philip Rinn
# Generate standardized Ornstein-Uhlenbeck-Process (d11=-1, d20=1) # with integration time step 0.01 and sampling frequency 1 s <- timeseries1D(N=1e4, sf=1, dt=0.01); t <- 1:1e4; plot(t, s, t="l", main=paste("mean:", mean(s), " var:", var(s)));
# Generate standardized Ornstein-Uhlenbeck-Process (d11=-1, d20=1) # with integration time step 0.01 and sampling frequency 1 s <- timeseries1D(N=1e4, sf=1, dt=0.01); t <- 1:1e4; plot(t, s, t="l", main=paste("mean:", mean(s), " var:", var(s)));
timeseries2D
generates a two-dimensional Langevin process using a
simple Euler integration. The drift function is a cubic polynomial, the
diffusion function a quadratic.
timeseries2D( N, startpointx = 0, startpointy = 0, D1_1 = matrix(c(0, -1, rep(0, 14)), nrow = 4), D1_2 = matrix(c(0, 0, 0, 0, -1, rep(0, 11)), nrow = 4), g_11 = matrix(c(1, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), g_12 = matrix(c(0, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), g_21 = matrix(c(0, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), g_22 = matrix(c(1, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), sf = 1000, dt = 0 )
timeseries2D( N, startpointx = 0, startpointy = 0, D1_1 = matrix(c(0, -1, rep(0, 14)), nrow = 4), D1_2 = matrix(c(0, 0, 0, 0, -1, rep(0, 11)), nrow = 4), g_11 = matrix(c(1, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), g_12 = matrix(c(0, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), g_21 = matrix(c(0, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), g_22 = matrix(c(1, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 3), sf = 1000, dt = 0 )
N |
a scalar denoting the length of the time-series to generate. |
startpointx |
a scalar denoting the starting point of the time series x. |
startpointy |
a scalar denoting the starting point of the time series y. |
D1_1 |
a 4x4 matrix denoting the coefficients of D1 for x. |
D1_2 |
a 4x4 matrix denoting the coefficients of D1 for y. |
g_11 |
a 3x3 matrix denoting the coefficients of g11 for x. |
g_12 |
a 3x3 matrix denoting the coefficients of g12 for x. |
g_21 |
a 3x3 matrix denoting the coefficients of g21 for y. |
g_22 |
a 3x3 matrix denoting the coefficients of g22 for y. |
sf |
a scalar denoting the sampling frequency. |
dt |
a scalar denoting the maximal time step of integration. Default
|
The elements of the matrices are defined by the corresponding
equations for the drift and diffusion terms:
with for
.
with for
timeseries2D
returns a time-series object with the generated
time-series as columns.
Philip Rinn