graph_tool.dynamics.KuramotoState#

class graph_tool.dynamics.KuramotoState(g, omega=1, w=1, sigma=0, t0=0, s=None)[source]#

Bases: ContinuousStateBase

The Kuramoto model.

Parameters:
gGraph

Graph to be used for the dynamics

omegaVertexPropertyMap or float (optional, default: 1)

Intrinsic frequencies for each node. If a scalar is given, it will be used for all nodes.

wEdgePropertyMap or float (optional, default: 1)

Coupling strength of each edge. If a scalar is given, it will be used for all edges.

sigmafloat (optional, default: .0)

Stochastic noise magnitude.

sVertexPropertyMap (optional, default: None)

Initial global state. If not provided, a random state will be chosen.

Notes

This implements Kuramoto’s model for synchronization [kuramoto_self-entrainment_1975] [rodrigues_kuramoto_2016].

Each node has an angle \(\theta_i\), which evolves in time obeying the differential equation:

\[\frac{\mathrm{d}\theta_i}{\mathrm{d}t} = \omega_i + \sum_{j}A_{ij}w_{ij}\sin(\theta_j-\theta_i) + \sigma\xi_i(t),\]

where \(\xi_i(t)\) is a Gaussian noise with zero mean and unit variance.

References

[kuramoto_self-entrainment_1975]

Y. Kuramoto, “Self-entrainment of a population of coupled non-linear oscillators”, International Symposium on Mathematical Problems in Theoretical Physics. Lecture Notes in Physics, vol 39. Springer, Berlin, Heidelberg (1975), DOI: 10.1007/BFb0013365 [sci-hub, @tor]

[rodrigues_kuramoto_2016]

Francisco A. Rodrigues, Thomas K. DM.Peron, Peng Ji, Jürgen Kurth, “The Kuramoto model in complex networks”, Physics Reports 610 1-98 (2016) DOI: 10.1016/j.physrep.2015.10.008 [sci-hub, @tor], arXiv: 1511.07139

Examples

>>> g = gt.collection.data["karate"]
>>> omega = g.new_vp("double", np.random.normal(0, 1, g.num_vertices())) 
>>> state = gt.KuramotoState(g, omega=omega, w=1.5)
>>> thetas = []
>>> ts = linspace(0, 40, 1000)
>>> for t in ts:
...     ret = state.solve(t, first_step=0.0001)
...     thetas.append(state.get_state().fa % (2 * pi))

>>> figure(figsize=(6, 4))
<...>
>>> for v in g.vertices():
...    plot(ts, [t[int(v)] - t.mean() for t in thetas])
[...]
>>> xlabel(r"Time")
Text(...)
>>> ylabel(r"$\theta_i - \left<\theta\right>$")
Text(...)
>>> tight_layout()
>>> savefig("karate-kuramoto.svg")
../_images/karate-kuramoto.svg

Kuramoto oscillator dynamics on the Karate Club network.#

Methods

copy()

Copy state.

get_diff(dt)

Returns the current time derivative for all the nodes.

get_r_phi()

Returns the phase coherence \(r\) and average phase \(\phi\), defined as

get_state()

Returns the internal VertexPropertyMap with the current state.

solve(t, *args, **kwargs)

Integrate the system up to time t.

solve_euler(t[, dt])

Integrate the system up o time t using a simple Euler's method with step size dt.

copy()#

Copy state.

get_diff(dt)#

Returns the current time derivative for all the nodes. The parameter dt is the time interval in consideration, which is used only if the ODE has a stochastic component.

If enabled during compilation, this algorithm runs in parallel.

get_r_phi()[source]#

Returns the phase coherence \(r\) and average phase \(\phi\), defined as

\[re^{i\phi} = \frac{1}{N}\sum_j e^{i\theta_j}.\]
get_state()#

Returns the internal VertexPropertyMap with the current state.

solve(t, *args, **kwargs)#

Integrate the system up to time t. The remaining parameters are passed to scipy.integrate.solve_ivp(). This solver is not suitable for stochastic ODEs.

solve_euler(t, dt=0.001)#

Integrate the system up o time t using a simple Euler’s method with step size dt. This solver is suitable for stochastic ODEs.