-- 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.2.0 module Numeric.DDE.Types -- | 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), ...). newtype RHS state RHS :: ((state, HistorySnapshot state, InputSnapshot) -> state) -> RHS state [_state] :: RHS state -> (state, HistorySnapshot state, InputSnapshot) -> state -- | Contains only the required snapshot of history to make steppers (e.g. -- Heun) work. There could be several delay variables newtype HistorySnapshot state Hist :: state -> HistorySnapshot state [_histsnap] :: HistorySnapshot state -> state -- | 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 -- | DDE stepper (all delays are equal). -- -- Stepper is a function of the following arguments: -- --
-- 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 --module Numeric.DDE -- | Generic integrator that records the whole time trace x(t) -- (single delay time). integ :: (Functor state, Storable (state Double), VectorSpace (state Double), Num (Scalar (state Double))) => Stepper -> state Double -> Vector (state Double) -> Int -> Scalar (state Double) -> RHS (state Double) -> Input -> (state Double, Vector (state Double)) -- | Generic integrator for DDEs (single delay time). Records all dynamical -- variables. integ' :: Storable state => (state -> (HistorySnapshot state, HistorySnapshot state) -> (Double, Double) -> state) -> Int -> Int -> Int -> (state, Vector state, Input) -> (state, Vector state) -- | RK4 integrator shortcut for 1D DDEs with zero initial conditions integRk4 :: Int -> Double -> RHS (V1 Double) -> Input -> (V1 Double, Vector (V1 Double)) -- | Shortcut for Heun's 2nd order 1D DDEs with zero initial conditions integHeun2 :: Int -> Double -> RHS (V1 Double) -> Input -> (V1 Double, Vector (V1 Double)) -- | RK4 integrator shortcut for 2D DDEs with zero initial conditions integRk4_2D :: Int -> Double -> RHS (V2 Double) -> Input -> (V2 Double, Vector (V2 Double)) -- | Shortcut for Heun's 2nd order 2D DDEs with zero initial conditions integHeun2_2D :: Int -> Double -> RHS (V2 Double) -> Input -> (V2 Double, Vector (V2 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 -- | Contains only the required snapshot of history to make steppers (e.g. -- Heun) work. There could be several delay variables newtype HistorySnapshot state Hist :: state -> HistorySnapshot state [_histsnap] :: HistorySnapshot state -> state -- | 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), ...). newtype RHS state RHS :: ((state, HistorySnapshot state, InputSnapshot) -> state) -> RHS state [_state] :: RHS state -> (state, HistorySnapshot state, InputSnapshot) -> state -- | DDE stepper (all delays are equal). -- -- Stepper is a function of the following arguments: -- --