numinoisy: Generates time series of deterministic-behavior with...

View source: R/numinoisy.R

numinoisyR Documentation

Generates time series of deterministic-behavior with stochatic perturbations (measurement and/or dynamical noise)

Description

Generates time series from Ordinary Differential Equations perturbed by dynamical and/or measurement noises

Usage

numinoisy(
  x0,
  t,
  K,
  varData = NULL,
  txVarBruitA = NULL,
  txVarBruitM = NULL,
  varBruitA = NULL,
  varBruitM = NULL,
  taux = NULL,
  freq = NULL,
  variables = NULL,
  method = NULL
)

Arguments

x0

The initial conditions. Should be a vector which size must be equal to the model dimension dim(K)[2] (the number of variables of the model defined by matrix K).

t

A vector providing all the dates for which the output are expected.

K

The Ordinary Differential Equations used to model the dynamics. The number of column should correspond to the number of variables, the number of lines to the number of parameters following the convention defined by poLabs(nVar,dMax).

varData

A vector of size nVar providing the caracteristic variances of each variable of the dynamical systems in ODE defined by matrix K. If not provided, this variance is automatically estimated.

txVarBruitA

A vector defining the ratio of ADDITIVE noise for each variable of the dynamical system in ODE. The additive noise is added at the end of the numerical integration process. The ratio is defined relatively to the signal variance of each variable.

txVarBruitM

A vector defining the ratio of DYNAMICAL noise for each variable of the dynamical system in ODE. This noise is a perturbation added at each numerical integration step. The ratio is defined relatively to the signal variance of each variable.

varBruitA

A vector defining the variance of ADDITIVE noise for each variable of the dynamical system in ODE. The additive noise is added at the end of the numerical integration process.

varBruitM

A vector defining the variance of DYNAMICAL noise for each variable of the dynamical system in ODE. This noise is a perturbation added at each numerical integration step.

taux

Generates random gaps in time series. Parameter taux defines the ratio of data to be kept (e.g. for taux=0.75, 75 percents of the data are kept).

freq

Subsamples the time series. Parameter freq defines the periodicity of data kept (e.g. for freq=3, 1 data out of 3 is kept).

variables

Defines which variables must be generated.

method

Defines the numerical integration method to be used. The fourth-order Runge-Kutta method is used by default (method = 'rk4'). Other method may be used (such as 'ode45' or 'lsoda'), see function ode from package deSolve for details.

Value

A list of two variables:

$donnees The integrated trajectory (first column is the time, next columns are the model variables)

$bruitM The level of dynamical noise

$bruitA The level of additive noise

$vectBruitM The vector of the dynamical noise used to produce the time series

$vectBruitA The vector of the additive noise used to produce the time series

$ecart_type The level standard deviation

Author(s)

Sylvain Mangiarotti, Malika Chassan

Examples

#############
# Example 1 #
#############
# Rossler Model formulation
# The model dimension
nVar = 3
 # maximal polynomial degree
dMax = 2
a = 0.520
b = 2
c = 4
Eq1 <- c(0,-1, 0,-1, 0, 0, 0, 0, 0, 0)
Eq2 <- c(0, 0, 0, a, 0, 0, 1, 0, 0, 0)
Eq3 <- c(b,-c, 0, 0, 0, 0, 0, 1, 0, 0)
K <- cbind(Eq1, Eq2, Eq3)
# Edit the equations
visuEq(K, nVar, dMax)
# initial conditions
v0 <- c(-0.6, 0.6, 0.4)
# output time required
timeOut = (0:800)/50
# variance of additive noise
varBruitA = c(0,0,0)^2
# variance of multiplitive noise
varBruitM = c(2E-2, 0, 2E-2)^2
# numerical integration with noise
intgr <- numinoisy(v0, timeOut, K, varBruitA = varBruitA, varBruitM = varBruitM, freq = 1)
# Plot of the simulated time series obtained
dev.new()
plot(intgr$donnees[,2], intgr$donnees[,3], type='l',
      main='phase portrait', xlab='x(t)', ylab = 'y(t)')
dev.new()
oldpar <- par(no.readonly = TRUE)    
on.exit(par(oldpar))  
par(mfrow = c(3, 1))
plot(intgr$donnees[,1], intgr$donnees[,2], type='l',
      main='phase portrait', xlab='x(t)', ylab = 'y(t)')
lines(intgr$donnees[,1], intgr$vectBruitM[,2]*10, type='l',
      main='phase portrait', xlab='x(t)', ylab = 'e(t)*10', col='red')
plot(intgr$donnees[,1], intgr$donnees[,3], type='l',
      main='phase portrait', xlab='x(t)', ylab = 'y(t)')
lines(intgr$donnees[,1], intgr$vectBruitM[,3]*10, type='l',
      main='phase portrait', xlab='x(t)', ylab = 'e(t)*10', col='red')
plot(intgr$donnees[,1], intgr$donnees[,4], type='l',
      main='phase portrait', xlab='x(t)', ylab = 'y(t)')
lines(intgr$donnees[,1], intgr$vectBruitM[,4]*10, type='l',
      main='phase portrait', xlab='x(t)', ylab = 'e(t)*10', col='red')


GPoM documentation built on July 9, 2023, 6:23 p.m.