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 GNU Public License 10 @cond GNU_PUBLIC_LICENSE 11 PyClustering is free software: you can redistribute it and/or modify 12 it under the terms of the GNU General Public License as published by 13 the Free Software Foundation, either version 3 of the License, or 14 (at your option) any later version. 16 PyClustering is distributed in the hope that it will be useful, 17 but WITHOUT ANY WARRANTY; without even the implied warranty of 18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 19 GNU General Public License for more details. 21 You should have received a copy of the GNU General Public License 22 along with this program. If not, see <http://www.gnu.org/licenses/>. 32 from scipy.integrate
import odeint
39 @brief Represents output dynamic of Sync in frequency domain. 46 @brief Constructor of Sync dynamic in frequency domain. 48 @param[in] amplitude (list): Dynamic of oscillators on each step of simulation. 49 @param[in] time (list): Simulation time where each time-point corresponds to amplitude-point. 60 @brief (list) Returns output dynamic of the Sync network (amplitudes of each oscillator in the network) during simulation. 70 @brief (list) Returns time-points corresponds to dynamic-points points. 79 @brief (uint) Returns number of simulation steps that are stored in dynamic. 88 @brief Indexing of the dynamic. 98 raise NameError(
'Out of range ' + index +
': only indexes 0 and 1 are supported.')
103 @brief Allocate clusters in line with ensembles of synchronous oscillators where each synchronous ensemble corresponds to only one cluster. 105 @param[in] tolerance (double): Maximum error for allocation of synchronous ensemble oscillators. 107 @return (list) Grours of indexes of synchronous oscillators, for example, [ [index_osc1, index_osc3], [index_osc2], [index_osc4, index_osc5] ]. 116 @brief Extracts number of oscillations of specified oscillator. 118 @param[in] index (uint): Index of oscillator whose dynamic is considered. 119 @param[in] amplitude_threshold (double): Amplitude threshold when oscillation is taken into account, for example, 120 when oscillator amplitude is greater than threshold then oscillation is incremented. 122 @return (uint) Number of oscillations of specified oscillator. 132 @brief Visualizer of output dynamic of sync network in frequency domain. 139 @brief Shows output dynamic (output of each oscillator) during simulation. 141 @param[in] fsync_output_dynamic (fsync_dynamic): Output dynamic of the fSync network. 143 @see show_output_dynamics 153 @brief Shows several output dynamics (output of each oscillator) during simulation. 154 @details Each dynamic is presented on separate plot. 156 @param[in] fsync_output_dynamics (list): list of output dynamics 'fsync_dynamic' of the fSync network. 158 @see show_output_dynamic 168 @brief Model of oscillatory network that uses Landau-Stuart oscillator and Kuramoto model as a synchronization mechanism. 169 @details Dynamic of each oscillator in the network is described by following differential Landau-Stuart equation with feedback: 172 \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}); 175 Where left part of the equation is Landau-Stuart equation and the right is a Kuramoto model for synchronization. 176 For solving this equation Runge-Kutta 4 method is used by default. 180 # Prepare oscillatory network parameters. 181 amount_oscillators = 3; 183 radiuses = [1.0, 2.0, 3.0]; 184 coupling_strength = 1.0; 186 # Create oscillatory network 187 oscillatory_network = fsync_network(amount_oscillators, frequency, radiuses, coupling_strength); 189 # Simulate network during 200 steps on 10 time-units of time-axis. 190 output_dynamic = oscillatory_network.simulate(200, 10, True); # True is to collect whole output dynamic. 192 # Visualize output result 193 fsync_visualizer.show_output_dynamic(output_dynamic); 196 Example of output dynamic of the network: 197 @image html fsync_sync_examples.png 201 __DEFAULT_FREQUENCY_VALUE = 1.0;
202 __DEFAULT_RADIUS_VALUE = 1.0;
203 __DEFAULT_COUPLING_STRENGTH = 1.0;
206 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):
208 @brief Constructor of oscillatory network based on synchronization Kuramoto model and Landau-Stuart oscillator. 210 @param[in] num_osc (uint): Amount oscillators in the network. 211 @param[in] factor_frequency (double|list): Frequency of oscillators, it can be specified as common value for all oscillators by 212 single double value and for each separately by list. 213 @param[in] factor_radius (double|list): Radius of oscillators that affects amplitude, it can be specified as common value for all oscillators by 214 single double value and for each separately by list. 215 @param[in] factor_coupling (double): Coupling strength between oscillators. 216 @param[in] type_conn (conn_type): Type of connection between oscillators in the network (all-to-all, grid, bidirectional list, etc.). 217 @param[in] representation (conn_represent): Internal representation of connection in the network: matrix or list. 221 super().
__init__(num_osc, type_conn, representation);
223 self.
__frequency = factor_frequency
if isinstance(factor_frequency, list)
else [ fsync_network.__DEFAULT_FREQUENCY_VALUE * factor_frequency
for _
in range(num_osc) ];
224 self.
__radius = factor_radius
if isinstance(factor_radius, list)
else [ fsync_network.__DEFAULT_RADIUS_VALUE * factor_radius
for _
in range(num_osc) ];
229 self.
__amplitude = [ random.random()
for _
in range(num_osc) ];
232 def simulate(self, steps, time, collect_dynamic = False):
234 @brief Performs static simulation of oscillatory network. 236 @param[in] steps (uint): Number simulation steps. 237 @param[in] time (double): Time of simulation. 238 @param[in] collect_dynamic (bool): If True - returns whole dynamic of oscillatory network, otherwise returns only last values of dynamics. 240 @return (list) Dynamic of oscillatory network. If argument 'collect_dynamic' is True, than return dynamic for the whole simulation time, 241 otherwise returns only last values (last step of simulation) of output dynamic. 244 @see simulate_dynamic() 248 dynamic_amplitude, dynamic_time = ([], [])
if collect_dynamic
is False else ([self.
__amplitude], [0]);
251 int_step = step / 10.0;
253 for t
in numpy.arange(step, time + step, step):
256 if collect_dynamic
is True:
257 dynamic_amplitude.append([ numpy.real(amplitude)[0]
for amplitude
in self.
__amplitude ]);
258 dynamic_time.append(t);
260 if collect_dynamic
is False:
261 dynamic_amplitude.append([ numpy.real(amplitude)[0]
for amplitude
in self.
__amplitude ]);
262 dynamic_time.append(time);
264 output_sync_dynamic =
fsync_dynamic(dynamic_amplitude, dynamic_time);
265 return output_sync_dynamic;
268 def __calculate(self, t, step, int_step):
270 @brief Calculates new amplitudes for oscillators in the network in line with current step. 272 @param[in] t (double): Time of simulation. 273 @param[in] step (double): Step of solution at the end of which states of oscillators should be calculated. 274 @param[in] int_step (double): Step differentiation that is used for solving differential equation. 276 @return (list) New states (phases) for oscillators. 280 next_amplitudes = [0.0] * self.
_num_osc;
282 for index
in range (0, self.
_num_osc, 1):
283 z = numpy.array(self.
__amplitude[index], dtype = numpy.complex128, ndmin = 1);
284 result = odeint(self.
__calculate_amplitude, z.view(numpy.float64), numpy.arange(t - step, t, int_step), (index , ));
285 next_amplitudes[index] = (result[len(result) - 1]).view(numpy.complex128);
287 return next_amplitudes;
290 def __oscillator_property(self, index):
292 @brief Calculate Landau-Stuart oscillator constant property that is based on frequency and radius. 294 @param[in] index (uint): Oscillator index whose property is calculated. 296 @return (double) Oscillator property. 300 return numpy.array(1j * self.
__frequency[index] + self.
__radius[index]**2, dtype = numpy.complex128, ndmin = 1);
303 def __landau_stuart(self, amplitude, index):
305 @brief Calculate Landau-Stuart state. 307 @param[in] amplitude (double): Current amplitude of oscillator. 308 @param[in] index (uint): Oscillator index whose state is calculated. 310 @return (double) Landau-Stuart state. 314 return (self.
__properties[index] - numpy.absolute(amplitude) ** 2) * amplitude;
317 def __synchronization_mechanism(self, amplitude, index):
319 @brief Calculate synchronization part using Kuramoto synchronization mechanism. 321 @param[in] amplitude (double): Current amplitude of oscillator. 322 @param[in] index (uint): Oscillator index whose synchronization influence is calculated. 324 @return (double) Synchronization influence for the specified oscillator. 328 sync_influence = 0.0;
332 amplitude_neighbor = numpy.array(self.
__amplitude[k], dtype = numpy.complex128, ndmin = 1);
333 sync_influence += amplitude_neighbor - amplitude;
338 def __calculate_amplitude(self, amplitude, t, argv):
340 @brief Returns new amplitude value for particular oscillator that is defined by index that is in 'argv' argument. 341 @details The method is used for differential calculation. 343 @param[in] amplitude (double): Current amplitude of oscillator. 344 @param[in] t (double): Current time of simulation. 345 @param[in] argv (uint): Index of the current oscillator. 347 @return (double) New amplitude of the oscillator. 351 z = amplitude.view(numpy.complex);
354 return dzdt.view(numpy.float64);
def output(self)
(list) Returns output dynamic of the Sync network (amplitudes of each oscillator in the network) duri...
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...
Model of oscillatory network that uses Landau-Stuart oscillator and Kuramoto model as a synchronizati...
def __synchronization_mechanism(self, amplitude, index)
Calculate synchronization part using Kuramoto synchronization mechanism.
Utils that are used by modules of pyclustering.
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 show_output_dynamics(fsync_output_dynamics)
Shows several output dynamics (output of each oscillator) during simulation.
def show_output_dynamic(fsync_output_dynamic)
Shows output dynamic (output of each oscillator) during simulation.
def extract_number_oscillations(osc_dyn, index=0, amplitude_threshold=1.0)
Extracts number of oscillations of specified oscillator.
def __len__(self)
(uint) Returns number of simulation steps that are stored in dynamic.
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 has_connection(self, i, j)
Returns True if there is connection between i and j oscillators and False - if connection doesn't exi...
def __oscillator_property(self, index)
Calculate Landau-Stuart oscillator constant property that is based on frequency and radius...
Represents output dynamic of Sync in frequency domain.
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 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.
Common network description that consists of information about oscillators and connection between them...
def time(self)
(list) Returns time-points corresponds to dynamic-points points.
def __getitem__(self, index)
Indexing of the dynamic.
def allocate_sync_ensembles(self, tolerance=0.1)
Allocate clusters in line with ensembles of synchronous oscillators where each synchronous ensemble c...
Visualizer of output dynamic of sync network in frequency domain.
def __calculate(self, t, step, int_step)
Calculates new amplitudes for oscillators in the network in line with current step.
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.
Neural and oscillatory network module.