dde-0.2.0: Delay differential equations

Safe HaskellNone
LanguageHaskell2010

Numeric.DDE

Contents

Description

Delay differential equations (DDE)

Example: Ikeda DDE

Below is a complete example simulating the Ikeda DDE defined as: tau * x(t)/dt = -x + beta * sin[x(t - tau_D)].

import           Linear ( V1 (..) )
import qualified Data.Vector.Storable as V
import qualified Numeric.DDE as DDE

ikedaRhs beta = DDE.RHS derivative
  where
    derivative ((V1 x), (DDE.Hist (V1 x_tauD)), _) = V1 x'
      where
        -- Ikeda DDE definition
        x' = (-x + beta * (sin x_tauD)) / tau

        -- Constants
        tau = 0.01

model beta hStep len1 totalIter = (state1, V.map (\(V1 x) -> x) trace)
  where
    -- Initial conditions:
    -- dynamical state and delay history.
    state0 = V1 (pi/2)
    hist0 = V.replicate len1 state0

    -- Input is ignored in ikedaRhs
    inp = DDE.Input $ V.replicate (totalIter + 1) 0

    (state1, trace) = DDE.integ DDE.rk4 state0 hist0 len1 hStep (ikedaRhs beta) inp

-- Control parameter
beta = 2.6

main = do
  let hStep = 0.001  -- Integration step
      tauD = 1.0  -- Delay time
      samplesPerDelay = round $ tauD / hStep
      delays = 8
      total = delays * samplesPerDelay

  let (state1, trace) = model beta hStep samplesPerDelay total

  mapM_ print $ V.toList trace

Synopsis

Integrators

integ Source #

Arguments

:: (Functor state, Storable (state Double), VectorSpace (state Double), Num (Scalar (state Double))) 
=> Stepper 
-> state Double

Initial state vector (x(t), y(t),...)

-> Vector (state Double)

Initial history for delayed variables

-> Int

Delay length in samples

-> Scalar (state Double)

Integration step

-> RHS (state Double)

Derivative (DDE right-hand side)

-> Input

External forcing

-> (state Double, Vector (state Double)) 

Generic integrator that records the whole time trace x(t) (single delay time).

integ' Source #

Arguments

:: Storable state 
=> (state -> (HistorySnapshot state, HistorySnapshot state) -> (Double, Double) -> state)

Iterator describing a DDE system

-> Int

Delay length in samples

-> Int

Number of last samples to record

-> Int

Total number of iterations

-> (state, Vector state, Input)

Initial state vector, initial history, and external forcing

-> (state, Vector state)

Final state and recorded state of the first variable. The latter is a vector of vectors (matrix) when multiple variables are involved.

Generic integrator for DDEs (single delay time). Records all dynamical variables.

integRk4 Source #

Arguments

:: Int

Delay length in samples

-> Double

Integration time step

-> RHS (V1 Double)

DDE model

-> Input

External forcing

-> (V1 Double, Vector (V1 Double)) 

RK4 integrator shortcut for 1D DDEs with zero initial conditions

integHeun2 Source #

Arguments

:: Int

Delay length in samples

-> Double

Integration time step

-> RHS (V1 Double)

DDE model

-> Input

External forcing

-> (V1 Double, Vector (V1 Double)) 

Shortcut for Heun's 2nd order 1D DDEs with zero initial conditions

integRk4_2D Source #

Arguments

:: Int

Delay length in samples

-> Double

Integration time step

-> RHS (V2 Double)

DDE model

-> Input

External forcing

-> (V2 Double, Vector (V2 Double)) 

RK4 integrator shortcut for 2D DDEs with zero initial conditions

integHeun2_2D Source #

Arguments

:: Int

Delay length in samples

-> Double

Integration time step

-> RHS (V2 Double)

DDE model

-> Input

External forcing

-> (V2 Double, Vector (V2 Double)) 

Shortcut for Heun's 2nd order 2D DDEs with zero initial conditions

newtype Input Source #

Vector of input data points

Constructors

Input 

Fields

newtype InputSnapshot Source #

Input u(t) is one-dimensional

Constructors

Inp 

Fields

newtype HistorySnapshot state Source #

Contains only the required snapshot of history to make steppers (e.g. Heun) work. There could be several delay variables

Constructors

Hist 

Fields

Steppers

newtype RHS state Source #

DDE right-hand side.

Parameter state is and abstraction of a dynamical system's state, i.e. it can be a vector of any length (x(t), y(t), ...).

Constructors

RHS 

Fields

type Stepper = forall state. (Functor state, VectorSpace (state Double), Num (Scalar (state Double))) => Scalar (state Double) -> RHS (state Double) -> state Double -> (HistorySnapshot (state Double), HistorySnapshot (state Double)) -> (Double, Double) -> state Double Source #

DDE stepper (all delays are equal).

Stepper is a function of the following arguments:

  • Integration step
  • DDE right-hand side
  • Current state vector (x(t), y(t), ...)
  • Two subsequent history snapshots
  • Two subsequent inputs

The result (step) is a new state vector.

rk4 :: Stepper Source #

Fourth order Runge-Kutta stepper

heun2 :: Stepper Source #

Second order Heun's stepper