-- Hoogle documentation, generated by Haddock -- See Hoogle, http://www.haskell.org/hoogle/ -- | Delay differential equations -- -- Please see the README on Github at -- https://github.com/masterdezign/dde#readme @package dde @version 0.0.1 module Numeric.DDE.Types -- | DDE right-hand side type RHS = (State, HistorySnapshot, InputSnapshot) -> State -- | Contains only the required snapshot of history to make steppers (e.g. -- Heun) work. There could be several delay variables newtype HistorySnapshot Hist :: Vector Double -> HistorySnapshot [_histsnap] :: HistorySnapshot -> Vector Double -- | Vector of input data points newtype Input Input :: Vector Double -> Input [_input] :: Input -> Vector Double -- | Input u(t) is one-dimensional newtype InputSnapshot Inp :: Double -> InputSnapshot [_insnap] :: InputSnapshot -> Double -- | State of a dynamical system, it can be a vector of any length (x(t), -- y(t), ...). newtype State State :: Vector Double -> State [_state] :: State -> Vector Double -- | Stepper for DDEs with a single delay -- --
-- >>> _stepper stepSize rhs xyState xTau1_2 u1_2 --newtype Stepper1 Stepper1 :: (Double -> RHS -> State -> (Double, Double) -> (Double, Double) -> State) -> Stepper1 [_stepper] :: Stepper1 -> Double -> RHS -> State -> (Double, Double) -> (Double, Double) -> State -- |
-- 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 --module Numeric.DDE -- | Generic integrator that records the whole time trace x(t) integ :: Stepper1 -> State -> Vector Double -> Int -> Double -> RHS -> Input -> (State, Vector Double) -- | Generic integrator for DDEs with a single delay integ' :: (State -> (Double, Double) -> (Double, Double) -> State) -> Int -> Int -> Int -> (State, Vector Double, Input) -> (State, Vector Double) -- | RK4 integrator shortcut for 2D DDEs with zero initial conditions integRk4_2D :: Int -> Double -> RHS -> Input -> (State, Vector Double) -- | Shortcut for Heun's 2nd order 2D DDEs with zero initial conditions integHeun2_2D :: Int -> Double -> RHS -> Input -> (State, Vector Double) -- | Vector of input data points newtype Input Input :: Vector Double -> Input [_input] :: Input -> Vector Double -- | Input u(t) is one-dimensional newtype InputSnapshot Inp :: Double -> InputSnapshot [_insnap] :: InputSnapshot -> Double -- | State of a dynamical system, it can be a vector of any length (x(t), -- y(t), ...). newtype State State :: Vector Double -> State [_state] :: State -> Vector Double -- | Contains only the required snapshot of history to make steppers (e.g. -- Heun) work. There could be several delay variables newtype HistorySnapshot Hist :: Vector Double -> HistorySnapshot [_histsnap] :: HistorySnapshot -> Vector Double -- | Stepper for DDEs with a single delay -- --
-- >>> _stepper stepSize rhs xyState xTau1_2 u1_2 --newtype Stepper1 Stepper1 :: (Double -> RHS -> State -> (Double, Double) -> (Double, Double) -> State) -> Stepper1 [_stepper] :: Stepper1 -> Double -> RHS -> State -> (Double, Double) -> (Double, Double) -> State -- | Fourth order Runge-Kutta stepper rk4 :: Stepper1 -- | Second order Heun's stepper heun2 :: Stepper1