| Safe Haskell | None |
|---|---|
| Language | Haskell2010 |
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 qualified Data.Vector.Storable as V
import qualified Numeric.DDE as DDE
ikedaRhs beta ((DDE.State xs), (DDE.Hist hs), _) = DDE.State $ V.fromList [x']
where
-- Ikeda DDE definition
x' = (-x + beta * (sin x_tauD)) / tau
-- Constants
tau = 0.01
-- Dynamical variable x(t)
x = V.head xs
-- Delay term x(t - tau_D)
x_tauD = V.head hs
model beta hStep len1 totalIter = DDE.integ DDE.rk4 state0 hist0 len1 hStep (ikedaRhs beta) inp
where
-- Initial conditions:
-- dynamical state and delay history.
state0 = DDE.State $ V.fromList [pi/2]
hist0 = V.replicate len1 (pi/2)
-- Input is ignored in ikedaRhs
inp = DDE.Input $ V.replicate (totalIter + 1) 0
-- 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- integ :: Stepper1 -> State -> Vector Double -> Int -> Double -> RHS -> Input -> (State, Vector Double)
- integ' :: (State -> (Double, Double) -> (Double, Double) -> State) -> Int -> Int -> Int -> (State, Vector Double, Input) -> (State, Vector Double)
- integRk4_2D :: Int -> Double -> RHS -> Input -> (State, Vector Double)
- integHeun2_2D :: Int -> Double -> RHS -> Input -> (State, Vector Double)
- newtype Input = Input {}
- newtype InputSnapshot = Inp {}
- newtype State = State {}
- newtype HistorySnapshot = Hist {}
- newtype Stepper1 = Stepper1 {}
- rk4 :: Stepper1
- heun2 :: Stepper1
Integrators
Arguments
| :: Stepper1 | |
| -> State | Initial state x(t), y(t),... |
| -> Vector Double | Initial history for the delayed variable |
| -> Int | Delay length in samples |
| -> Double | Integration step |
| -> RHS | Derivative (DDE right-hand side) |
| -> Input | External forcing |
| -> (State, Vector Double) |
Generic integrator that records the whole time trace x(t)
Arguments
| :: (State -> (Double, Double) -> (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 Double, Input) | Initial state vector, initial history, and external forcing |
| -> (State, Vector Double) | Final state and recorded state of the first variable. The latter could be a Matrix if multiple variables are needed |
Generic integrator for DDEs with a single delay
Arguments
| :: Int | Delay length in samples |
| -> Double | Integration time step |
| -> RHS | DDE model |
| -> Input | External forcing |
| -> (State, Vector Double) |
RK4 integrator shortcut for 2D DDEs with zero initial conditions
Arguments
| :: Int | Delay length in samples |
| -> Double | Integration time step |
| -> RHS | DDE model |
| -> Input | External forcing |
| -> (State, Vector Double) |
Shortcut for Heun's 2nd order 2D DDEs with zero initial conditions
State of a dynamical system, it can be a vector of any length (x(t), y(t), ...).
newtype HistorySnapshot Source #
Contains only the required snapshot of history to make steppers (e.g. Heun) work. There could be several delay variables
Steppers
Stepper for DDEs with a single delay
>>>_stepper stepSize rhs xyState xTau1_2 u1_2