Basins.basin_entropyMethod
basin_entropy(basin; eps_x=20, eps_y=20)

This algorithm computes the basin entropy of a computed basin of attraction on a regular grid. The function return the basin entropy and the boundary basin entropy.

[A. Daza, A. Wagemakers, B. Georgeot, D. Guéry-Odelin and M. A. F. Sanjuán, Basin entropy: a new tool to analyze uncertainty in dynamical systems, Sci. Rep., 6, 31416, (2016).]

Arguments

  • basin : the matrix containing the information of the basin.

Keyword arguments

  • eps_x, eps_y : define the size of the box that will be used to sample the basin
Basins.basin_stabilityMethod
basin_stability(basin)

This algorithm computes the basin stability of a computed basin of attraction on a regular grid. This function returns a vector with the relative volume of the basins.

[P. Menck, J. Heitzig, N. Marwan et al. How basin stability complements the linear-stability paradigm. Nature Phys 9, 89–92 (2013).]

Arguments

  • basin : the matrix containing the information of the basin.
Basins.basins_generalMethod
basins_general(xg, yg, integ; T=1., idxs=1:2)

Compute an estimate of the basin of attraction on a two-dimensional plane using a stroboscopic map.

Arguments

  • xg, yg : 1-dim range vector that defines the grid of the initial conditions to test.
  • integ : integrator handle of the dynamical system.

Keyword Arguments

  • dt : Time step of the integrator. It recommended to use values above 1.
  • idxs = 1:D : Optionally you can choose which variables to save. Defaults to the entire state.
  • Ncheck : A parameter that sets the number of consecutives hits of an attractor before deciding the basin

of the initial condition.

Example:

```jl using DynamicalSystems, ChaosTools ds = Systems.magneticpendulum(γ=1, d=0.2, α=0.2, ω=0.8, N=3) integ = integrator(ds, u0=[0,0,0,0], reltol=1e-9) xg=range(-4,4,length=150) yg=range(-4,4,length=150) @time bsn = basinsgeneral(xg, yg, integ; dt=1., idxs=1:2)

Basins.basins_map2DMethod
basins_map2D(xg, yg, integ; kwargs...)

Compute an estimate of the basin of attraction on a two-dimensional plane using a map of the plane onto itself. The dynamical system should be a discrete two dimensional system such as: * Discrete 2D map. * 2D poincaré map. * A 2D stroboscopic map.

[H. E. Nusse and J. A. Yorke, Dynamics: numerical explorations, Springer, New York, 1997]

Arguments

  • xg, yg : 1-dim range vector that defines the grid of the initial conditions to test.
  • integ : A integrator handle of the dynamical system. For a Poincaré map, the handle is a pmap

as defined in poincaremap

Keyword Arguments

  • T : Period of the stroboscopic map.
  • Ncheck : A parameter that sets the number of consecutives hits of an attractor before deciding the basin

of the initial condition.

Example:

using DynamicalSystems, ChaosTools
ds = Systems.rikitake(μ = 0.47, α = 1.0)
xg=range(-6.,6.,length=150); yg=range(-6.,6.,length=150)
pmap = poincaremap(ds, (3, 0.), Tmax=1e6; idxs = 1:2, rootkw = (xrtol = 1e-8, atol = 1e-8), reltol=1e-9)
bsn = basins_map2D(xg, yg, pmap)
Basins.box_counting_dimMethod
box_counting_dim(xg, yg, basin; kwargs...)

This function compute the box-counting dimension of the boundary. It is related to the uncertainty dimension.

[C. Grebogi, S. W. McDonald, E. Ott, J. A. Yorke, Final state sensitivity: An obstruction to predictability, Physics Letters A, 99, 9, 1983]

Arguments

  • basin : the matrix containing the information of the basin.
  • xg, yg : 1-dim range vector that defines the grid of the initial conditions to test.

Keyword arguments

  • kwargs... these arguments are passed to generalized_dim see ChaosTools.jl for the docs.
Basins.compute_saddleMethod
compute_saddle(bsn_nfo::BasinInfo, integ; max_iter=10)

The algorithm the saddle that lies in a boundary of the basin of attraction of a dynamical system. When the boundary is fractal, this set is known as the chaotic saddle and is responsible for the transient dynamics of the system. The saddle is computed with the Saddle-Straddle algorithm. It is necessary to define two generalized basin, that is we must separate the basin into two sets (see also keyword arguments).

[H. E. Nusse and J. A. Yorke, Dynamics: numerical explorations, Springer, New York, 2012]

Arguments

  • bsn_nfo : structure that holds the information of the basin as well as the map function. This structure is set when the basin is first computed with basin_stroboscopic_map or basin_poincare_map.
  • integ : the matrix containing the information of the basin.
  • bas_A : vector with the indices of the attractors that will represent the generalized basin A
  • bas_B : vector with the indices of the attractors that will represent the generalized basin B. Notice that bas_A ∪ bas_B = [1:N] and bas_A ∩ bas_B = ∅

Keyword arguments

  • N : number of points of the saddle to compute
  • init_tol: length of the straddle segment

Example

# compute 1000 points between generalized A = [1] and basin B = [2,3]
sa,sb = compute_saddle(bsn, integ_df, [1], [2,3],1000)
# sa is the "left" set and sb the right "set"
Basins.detect_wada_grid_methodMethod
detect_wada_grid_method(integ, bsn_nfo::BasinInfo; max_iter=10)

The algorithm test for Wada basin in a dynamical system. It uses the dynamical system to look if all the atractors are represented in the boundary.

[A. Daza, A. Wagemakers, M. A. F. Sanjuán and J. A. Yorke, Testing for Basins of Wada, Sci. Rep., 5, 16579 (2015)]

Arguments

  • integ : the matrix containing the information of the basin.
  • bsn_nfo : structure that holds the information of the basin as well as the map function. This structure is set when the basin is first computed with basin_stroboscopic_map or basin_poincare_map.

Keyword arguments

  • max_iter : set the maximum depth of subdivisions to look for an atractor. The number of points doubles at each step.
Basins.detect_wada_merge_methodMethod
detect_wada_merge_method(xg,yg, basin)

The algorithm gives the maximum and minimum Haussdorff distances between merged basins. These two distances can help to decide if the basin has the Wada property.

[A. Daza, A. Wagemakers and M. A. F. Sanjuán, Ascertaining when a basin is Wada: the merging method, Sci. Rep., 8, 9954 (2018)]

Arguments

  • basin : the matrix containing the information of the basin.
  • xg, yg : 1-dim range vector that defines the grid of the initial conditions to test.

Example

max_dist,min_dist = detect_wada_merge_method(basin,xg,yg)
# grid resolution
epsilon = xg[2]-xg[1]
# if dmax is large then this is not Wada
@show dmax = max_dist/epsilon
@show dmin = min_dist/epsilon
Basins.draw_basin!Method
draw_basin!(xg, yg, integ, iter_f!::Function, reinit_f!::Function)

Compute an estimate of the basin of attraction on a two-dimensional plane. This is a low level function, for higher level functions see: basin_poincare_map, basin_discrete_map, basin_stroboscopic_map

Arguments

  • xg, yg : 1-dim range vector that defines the grid of the initial conditions to test.
  • integ : integrator handle of the dynamical system.
  • iter_f! : function that iterates the map or the system, see step! from DifferentialEquations.jl and

examples for a Poincaré map of a continuous system.

  • reinit_f! : function that sets the initial condition to test on a two dimensional projection of the phase space.
Basins.uncertainty_exponentMethod
uncertainty_exponent(bsn::BasinInfo, integ; precision=1e-3) -> α,e,ε,f_ε

This function estimates the uncertainty exponent of the boundary. It is related to the uncertainty dimension.

[C. Grebogi, S. W. McDonald, E. Ott, J. A. Yorke, Final state sensitivity: An obstruction to predictability, Physics Letters A, 99, 9, 1983]

Arguments

  • bsn : structure that holds the information of the basin.
  • integ : handle to the iterator of the dynamical system.

Keyword arguments

  • precision variance of the estimator of the uncertainty function.