|
pyclustering
0.10.1
pyclustring is a Python, C++ data mining library.
|
3 @brief Oscillatory Neural Network based on Kuramoto model in frequency domain.
4 @details Implementation based on paper @cite book::chemical_oscillatorions_waves.
6 @authors Andrei Novikov (pyclustering@yandex.ru)
8 @copyright BSD-3-Clause
17 from scipy.integrate
import odeint
24 @brief Represents output dynamic of Sync in frequency domain.
31 @brief Constructor of Sync dynamic in frequency domain.
33 @param[in] amplitude (list): Dynamic of oscillators on each step of simulation.
34 @param[in] time (list): Simulation time where each time-point corresponds to amplitude-point.
45 @brief (list) Returns output dynamic of the Sync network (amplitudes of each oscillator in the network) during simulation.
55 @brief (list) Returns time-points corresponds to dynamic-points points.
64 @brief (uint) Returns number of simulation steps that are stored in dynamic.
73 @brief Indexing of the dynamic.
83 raise NameError(
'Out of range ' + index +
': only indexes 0 and 1 are supported.')
88 @brief Allocate clusters in line with ensembles of synchronous oscillators where each synchronous ensemble corresponds to only one cluster.
90 @param[in] tolerance (double): Maximum error for allocation of synchronous ensemble oscillators.
92 @return (list) Grours of indexes of synchronous oscillators, for example, [ [index_osc1, index_osc3], [index_osc2], [index_osc4, index_osc5] ].
101 @brief Extracts number of oscillations of specified oscillator.
103 @param[in] index (uint): Index of oscillator whose dynamic is considered.
104 @param[in] amplitude_threshold (double): Amplitude threshold when oscillation is taken into account, for example,
105 when oscillator amplitude is greater than threshold then oscillation is incremented.
107 @return (uint) Number of oscillations of specified oscillator.
117 @brief Visualizer of output dynamic of sync network in frequency domain.
124 @brief Shows output dynamic (output of each oscillator) during simulation.
126 @param[in] fsync_output_dynamic (fsync_dynamic): Output dynamic of the fSync network.
128 @see show_output_dynamics
138 @brief Shows several output dynamics (output of each oscillator) during simulation.
139 @details Each dynamic is presented on separate plot.
141 @param[in] fsync_output_dynamics (list): list of output dynamics 'fsync_dynamic' of the fSync network.
143 @see show_output_dynamic
153 @brief Model of oscillatory network that uses Landau-Stuart oscillator and Kuramoto model as a synchronization mechanism.
154 @details Dynamic of each oscillator in the network is described by following differential Landau-Stuart equation with feedback:
157 \dot{z}_{i} = (i\omega_{i} + \rho^{2}_{i} - |z_{i}|^{2} )z_{i} + \frac{1}{N}\sum_{j=0}^{N}k_{ij}(z_{j} - z_{i});
160 Where left part of the equation is Landau-Stuart equation and the right is a Kuramoto model for synchronization.
161 For solving this equation Runge-Kutta 4 method is used by default.
165 # Prepare oscillatory network parameters.
166 amount_oscillators = 3;
168 radiuses = [1.0, 2.0, 3.0];
169 coupling_strength = 1.0;
171 # Create oscillatory network
172 oscillatory_network = fsync_network(amount_oscillators, frequency, radiuses, coupling_strength);
174 # Simulate network during 200 steps on 10 time-units of time-axis.
175 output_dynamic = oscillatory_network.simulate(200, 10, True); # True is to collect whole output dynamic.
177 # Visualize output result
178 fsync_visualizer.show_output_dynamic(output_dynamic);
181 Example of output dynamic of the network:
182 @image html fsync_sync_examples.png
186 __DEFAULT_FREQUENCY_VALUE = 1.0;
187 __DEFAULT_RADIUS_VALUE = 1.0;
188 __DEFAULT_COUPLING_STRENGTH = 1.0;
191 def __init__(self, num_osc, factor_frequency = 1.0, factor_radius = 1.0, factor_coupling = 1.0, type_conn = conn_type.ALL_TO_ALL, representation = conn_represent.MATRIX):
193 @brief Constructor of oscillatory network based on synchronization Kuramoto model and Landau-Stuart oscillator.
195 @param[in] num_osc (uint): Amount oscillators in the network.
196 @param[in] factor_frequency (double|list): Frequency of oscillators, it can be specified as common value for all oscillators by
197 single double value and for each separately by list.
198 @param[in] factor_radius (double|list): Radius of oscillators that affects amplitude, it can be specified as common value for all oscillators by
199 single double value and for each separately by list.
200 @param[in] factor_coupling (double): Coupling strength between oscillators.
201 @param[in] type_conn (conn_type): Type of connection between oscillators in the network (all-to-all, grid, bidirectional list, etc.).
202 @param[in] representation (conn_represent): Internal representation of connection in the network: matrix or list.
206 super().
__init__(num_osc, type_conn, representation);
208 self.
__frequency = factor_frequency
if isinstance(factor_frequency, list)
else [ fsync_network.__DEFAULT_FREQUENCY_VALUE * factor_frequency
for _
in range(num_osc) ];
209 self.
__radius = factor_radius
if isinstance(factor_radius, list)
else [ fsync_network.__DEFAULT_RADIUS_VALUE * factor_radius
for _
in range(num_osc) ];
214 self.
__amplitude = [ random.random()
for _
in range(num_osc) ];
217 def simulate(self, steps, time, collect_dynamic = False):
219 @brief Performs static simulation of oscillatory network.
221 @param[in] steps (uint): Number simulation steps.
222 @param[in] time (double): Time of simulation.
223 @param[in] collect_dynamic (bool): If True - returns whole dynamic of oscillatory network, otherwise returns only last values of dynamics.
225 @return (list) Dynamic of oscillatory network. If argument 'collect_dynamic' is True, than return dynamic for the whole simulation time,
226 otherwise returns only last values (last step of simulation) of output dynamic.
229 @see simulate_dynamic()
233 dynamic_amplitude, dynamic_time = ([], [])
if collect_dynamic
is False else ([self.
__amplitude], [0]);
236 int_step = step / 10.0;
238 for t
in numpy.arange(step, time + step, step):
241 if collect_dynamic
is True:
242 dynamic_amplitude.append([ numpy.real(amplitude)[0]
for amplitude
in self.
__amplitude ]);
243 dynamic_time.append(t);
245 if collect_dynamic
is False:
246 dynamic_amplitude.append([ numpy.real(amplitude)[0]
for amplitude
in self.
__amplitude ]);
247 dynamic_time.append(time);
249 output_sync_dynamic =
fsync_dynamic(dynamic_amplitude, dynamic_time);
250 return output_sync_dynamic;
253 def __calculate(self, t, step, int_step):
255 @brief Calculates new amplitudes for oscillators in the network in line with current step.
257 @param[in] t (double): Time of simulation.
258 @param[in] step (double): Step of solution at the end of which states of oscillators should be calculated.
259 @param[in] int_step (double): Step differentiation that is used for solving differential equation.
261 @return (list) New states (phases) for oscillators.
265 next_amplitudes = [0.0] * self.
_num_osc;
267 for index
in range (0, self.
_num_osc, 1):
268 z = numpy.array(self.
__amplitude[index], dtype = numpy.complex128, ndmin = 1);
269 result = odeint(self.
__calculate_amplitude, z.view(numpy.float64), numpy.arange(t - step, t, int_step), (index , ));
270 next_amplitudes[index] = (result[len(result) - 1]).view(numpy.complex128);
272 return next_amplitudes;
275 def __oscillator_property(self, index):
277 @brief Calculate Landau-Stuart oscillator constant property that is based on frequency and radius.
279 @param[in] index (uint): Oscillator index whose property is calculated.
281 @return (double) Oscillator property.
285 return numpy.array(1j * self.
__frequency[index] + self.
__radius[index]**2, dtype = numpy.complex128, ndmin = 1);
288 def __landau_stuart(self, amplitude, index):
290 @brief Calculate Landau-Stuart state.
292 @param[in] amplitude (double): Current amplitude of oscillator.
293 @param[in] index (uint): Oscillator index whose state is calculated.
295 @return (double) Landau-Stuart state.
299 return (self.
__properties[index] - numpy.absolute(amplitude) ** 2) * amplitude;
302 def __synchronization_mechanism(self, amplitude, index):
304 @brief Calculate synchronization part using Kuramoto synchronization mechanism.
306 @param[in] amplitude (double): Current amplitude of oscillator.
307 @param[in] index (uint): Oscillator index whose synchronization influence is calculated.
309 @return (double) Synchronization influence for the specified oscillator.
313 sync_influence = 0.0;
317 amplitude_neighbor = numpy.array(self.
__amplitude[k], dtype = numpy.complex128, ndmin = 1);
318 sync_influence += amplitude_neighbor - amplitude;
323 def __calculate_amplitude(self, amplitude, t, argv):
325 @brief Returns new amplitude value for particular oscillator that is defined by index that is in 'argv' argument.
326 @details The method is used for differential calculation.
328 @param[in] amplitude (double): Current amplitude of oscillator.
329 @param[in] t (double): Current time of simulation.
330 @param[in] argv (uint): Index of the current oscillator.
332 @return (double) New amplitude of the oscillator.
336 z = amplitude.view(numpy.complex);
339 return dzdt.view(numpy.float64);
def __calculate_amplitude(self, amplitude, t, argv)
Returns new amplitude value for particular oscillator that is defined by index that is in 'argv' argu...
def __init__(self, num_osc, factor_frequency=1.0, factor_radius=1.0, factor_coupling=1.0, type_conn=conn_type.ALL_TO_ALL, representation=conn_represent.MATRIX)
Constructor of oscillatory network based on synchronization Kuramoto model and Landau-Stuart oscillat...
def extract_number_oscillations(self, index, amplitude_threshold)
Extracts number of oscillations of specified oscillator.
def __landau_stuart(self, amplitude, index)
Calculate Landau-Stuart state.
def draw_dynamics_set(dynamics, xtitle=None, ytitle=None, xlim=None, ylim=None, xlabels=False, ylabels=False)
Draw lists of dynamics of neurons (oscillators) in the network.
def __synchronization_mechanism(self, amplitude, index)
Calculate synchronization part using Kuramoto synchronization mechanism.
Visualizer of output dynamic of sync network in frequency domain.
def extract_number_oscillations(osc_dyn, index=0, amplitude_threshold=1.0)
Extracts number of oscillations of specified oscillator.
def allocate_sync_ensembles(dynamic, tolerance=0.1, threshold=1.0, ignore=None)
Allocate clusters in line with ensembles of synchronous oscillators where each synchronous ensemble c...
def show_output_dynamic(fsync_output_dynamic)
Shows output dynamic (output of each oscillator) during simulation.
def draw_dynamics(t, dyn, x_title=None, y_title=None, x_lim=None, y_lim=None, x_labels=True, y_labels=True, separate=False, axes=None)
Draw dynamics of neurons (oscillators) in the network.
Represents output dynamic of Sync in frequency domain.
def __len__(self)
(uint) Returns number of simulation steps that are stored in dynamic.
def __oscillator_property(self, index)
Calculate Landau-Stuart oscillator constant property that is based on frequency and radius.
Common network description that consists of information about oscillators and connection between them...
def has_connection(self, i, j)
Returns True if there is connection between i and j oscillators and False - if connection doesn't exi...
def time(self)
(list) Returns time-points corresponds to dynamic-points points.
def __getitem__(self, index)
Indexing of the dynamic.
Model of oscillatory network that uses Landau-Stuart oscillator and Kuramoto model as a synchronizati...
Neural and oscillatory network module. Consists of models of bio-inspired networks.
def allocate_sync_ensembles(self, tolerance=0.1)
Allocate clusters in line with ensembles of synchronous oscillators where each synchronous ensemble c...
def __calculate(self, t, step, int_step)
Calculates new amplitudes for oscillators in the network in line with current step.
def show_output_dynamics(fsync_output_dynamics)
Shows several output dynamics (output of each oscillator) during simulation.
def simulate(self, steps, time, collect_dynamic=False)
Performs static simulation of oscillatory network.
def __init__(self, amplitude, time)
Constructor of Sync dynamic in frequency domain.
Utils that are used by modules of pyclustering.
def output(self)
(list) Returns output dynamic of the Sync network (amplitudes of each oscillator in the network) duri...