snudda.detect.detect

class snudda.detect.detect.SnuddaDetect(config_file=None, network_path=None, snudda_data=None, position_file=None, voxel_size=3e-06, hyper_voxel_size=100, verbose=False, logfile_name=None, logfile=None, save_file=None, work_history_file=None, slurm_id=0, volume_id=None, role=None, rc=None, axon_stump_id_flag=False, simulation_origo=None, h5libver=None, random_seed=None, debug_flag=False)[source]

SnuddaDetect places synapses in the network based on touch detection.

Constructor.

Parameters
  • network_path (str) – Network directory

  • snudda_data (str, optional) – Path to SNUDDA_DATA (if not specified will try to read from config file)

  • config_file (str, optional) – Network config file (default network-config.json in network_path)

  • position_file (str, optional) – Network position file (default network-neuron-positions in network_path)

  • voxel_size (float, optional) – Width of voxel (default 3e-6m)

  • hyper_voxel_size (int, optional) – Number of voxels per side (default: 100, ie 100x100x100 voxels total)

  • verbose (bool, optional) – Verbose mode (default False)

  • logfile_name (str, optional) – Name of log file

  • logfile (_io.TextIOWrapper, optional) – Pointer to already open log file

  • save_file (str, optional) – Name of output file (default voxels/network-putative-synapses.hdf5 in network_path)

  • work_history_file (str, optional) – Work log file (default network-detect-worklog.hdf5 in network_path)

  • slurm_id (int, optional) – SlurmID of job

  • volume_id (str, optional) – Volume ID to do touch detection on

  • role (str, optional) – Parallel role, i.e. “master” or “worker”

  • rc (ipyparallel.Client, optional) – iPyParallel client, if given program will run in parallel

  • axon_stump_id_flag (bool, optional) – Recalculate segment IDs to account for axon stump? (default False)

  • simulation_origo (np.array, optional) – Origo for touch detection hypervoxels and voxels, voxel coordinates must always positive.

  • h5libver (string, optional) – h5py library version (default “latest”)

  • random_seed (int, optional) – Random seed

  • debug_flag (bool, optional) – Save additional information for debugging (Default: False)

static create_lookup_table(data, n_rows, data_type, num_neurons, max_synapse_type)[source]

This creates a lookup table where all synapses in the hyper voxel between the same pair of neurons are grouped together in the synapse matrix. Returns a matrix where first column is a UID = srcID*nNeurons + destID and the following two columns are start row and end row (-1) in matrix

Parameters
  • data – either synapse matrix, or gap junction matrix

  • n_rows – number of rows in matrix that are used (matrix itself can be larger)

  • data_type – “synapses” or “gap_junctions”

  • num_neurons – number of neurons

  • max_synapse_type – the synapse types are numbered, this number must not be too small.

Returns a matrix where first column is a UID = src_ID*num_neurons + dest_ID and the following two columns are start row and end row (-1) in matrix

delete_old_merge()[source]

Cleans up data files from previous detection run.

detect(restart_detection_flag=True, rc=None)[source]

Synapse placement based on touch detection. Space is divided into hyper voxels, containing 100x100x100 voxels. Each hyper voxel is processed separately.

Parameters
  • restart_detection_flag (bool, optional) – Restart detection or resume previous partial run.

  • rc (ipyparallel.Client, optional) – Remote Client, used for parallel execution

detect_gap_junctions()[source]

Helper function, triggers detection of gap junctions. Called by process_hyper_voxel.

detect_synapses()[source]

Helper function, triggers detection of synapses. Called by process_hyper_voxel.

distribute_neurons(neuron_idx=None, distribution_seeds=None, min_coord=None, max_coord=None)[source]

This creates a list for each hyper voxel of the neurons that has any neurites within its border (here defined as vertices inside region)

Parameters
  • neuron_idx – NeuronIDs to process

  • distribution_seeds – Random seed (used for neurons without axon)

  • min_coord (float, float, float) – Minimum x,y,z coordinates

  • max_coord (float, float, float) – Maximum x,y,z coordinates

Updates self.hyper_voxels. Also returns min_coord, max_coord

distribute_neurons_parallel(d_view=None)[source]

Locates which hyper voxel each neuron is present in.

export_voxel_visualisation_csv(neuron_id)[source]

Export CSV file with voxel data, used for visualisation.

fill_voxels_axon(voxel_space, voxel_space_ctr, voxel_axon_dist, coords, links, neuron_id)[source]

Mark all voxels containing axons.

voxel_space : n x n x n x k matrix holding the voxel content, normally self.axon_voxels (neuron IDs) voxel_space_ctr : n x n x n matrix holding count of how many items each voxel holds coords : neuron vertices, n x 3 matrix links : how do vertices link up to for axon segments n x 2 matrix neuron_id : ID of the neurons

fill_voxels_dend(voxel_space, voxel_space_ctr, voxel_sec_id, voxel_sec_x, voxel_soma_dist, coords, links, seg_id, seg_x, neuron_id)[source]

Mark all voxels containing dendrites.

voxel_space : n x n x n x k matrix holding the voxel content, normally self.dend_voxels (neuron IDs) voxel_space_ctr : n x n x n matrix holding count of how many items each voxel holds voxel_sec_id : n x n x n x k matrix, holding section ID of each item voxel_sec_x : n x n x n x k matrix, holding section X of each item voxel_soma_dist : n x n x n x k matrix, holding distance to soma along dendrite coords : neuron vertices, n x 3 matrix links : how do vertices link up to for dendrite segments n x 2 matrix seg_id : segment ID of the links end points seg_x : segment X of the links end points neuron_id : ID of the neurons

fill_voxels_soma(voxel_space, voxel_space_ctr, voxel_sec_id, voxel_sec_x, soma_coord, neuron_id, verbose=False)[source]

Marks all the dendrite voxels that all the somas in the hyper voxel occupy.

voxel_space : n x n x n x k matrix holding the voxel content, normally self.dend_voxels (neuron IDs) voxel_space_ctr : n x n x n matrix holding count of how many items each voxel holds voxel_sec_id : n x n x n x k matrix, holding section ID of each item voxel_sec_x : n x n x n x k matrix, holding section X of each item soma_coord : (x,y,z,r) location of soma to place, and radius neuron_id : ID of the neurons verbose (bool) : how much to print out

find_min_max_coord(volume_id=None, neuron_idx=None)[source]

Finds the minimum and maximum coordinates in entire model for all neuron components.

Parameters

volume_id – Volume ID to check, None means all

Returns

min_coord, max_coord

find_min_max_coord_parallel(volume_id=None, d_view=None)[source]

Finds the minimum and maximum coordinates in entire model for all neuron components.

Parameters
  • volume_id – Volume ID to check, None means all

  • d_view – Direct view object

Returns

min_coord, max_coord

generate_hyper_voxel_random_seeds()[source]

Generates a seed sequence for each hyper voxel based on the master seed for touch detection.

generate_neuron_distribution_random_seeds()[source]

Generate seed sequence for neuron distribution from master seed.

get_hyper_voxel_axon_points(neuron_position, rotation, axon_density_bounds_xyz, num_points=1000)[source]

Helper function to give points inside axon bounding box, that are inside hyper voxel

Parameters
  • neuron_position (float,float,float) – coordinates of neuron

  • rotation – rotation matrix

  • axon_density_bounds_xyz – boundary box for axon

  • num_points – number of points to place

get_neuron_distribution_history()[source]

Returns info about what neurons each hyper voxel contains etc.

Returns

containing

hyper_voxels (dictionary): dictionary with keys ‘neurons’, ‘neuronCtr’, ‘origo’, ‘randomSeed’ hyper_voxel_id_lookup (3D matrix with int): hypervoxel ID, spatially arranged n_hyper_voxels (int): number of hypervoxels simulation_origo (float, float, float): origo of entire simulation

Return type

(tuple)

includes_gap_junctions()[source]

Checks if any gap junctions are defined in self.connectivity_distribution. Returns True or False.

load_neuron(neuron_info)[source]

Load neuron.

Parameters

neuron_info – dictionary with neuron information, i.e. ‘name’, ‘parameterID’, ‘morphologyID’, ‘modulationID’, ‘rotation’, ‘position’

no_axon_points_sphere(soma_centre, r_cum_distribution, num_points)[source]

Helper function placing axon segments with spherical probability distribution. This picks points around soma centre. num_points are randomized, points outside the hyper sphere are rejected, so fewer than nPoints might be returned.

Parameters
  • soma_centre (float,float,float) – x,y,z coordinates of soma centre

  • r_cum_distribution – cumulative distribution

  • num_points – number of points to place

Returns

voxel_coordinates synapse_distance_to_soma

no_axon_points_xyz(neuron_position, rotation, axon_density_func, axon_density_bounds_xyz)[source]

Placing fake axon segments based on probability distribution speicified with x,y,z.

Parameters
  • neuron_position (float,float,float) – location of neuron soma

  • rotation (3x3 rotation matrix) – rotation of neuron

  • axon_density_func – axon density function in x,y,z coordinates (SWC coords), must handle x,y,z as vectors e.g. axon_density_func = eval(“lambda x,y,z: ” + axonPstr)

  • axon_density_bounds_xyz – [xmin,xmax,ymin,ymax,zmin,zmax] using the coordinates in the SWC file

parallel_process_hyper_voxels(rc=None, d_view=None)[source]

Distributes touch detection in hyper voxels to workers.

Parameters

rc (ipyparallel.Client, optional) – Remote client, for parallel execution

place_synapses_no_axon(hyper_id, voxel_space, voxel_space_ctr, voxel_axon_dist)[source]

Places fake axon segments for neurons without axons.

Parameters
  • hyper_id (int) – Hypervoxel ID

  • voxel_space – Axon voxel space (n_bins x n_bins x n_bins x max_syn, voxel space matrix)

  • voxel_space_ctr – Synapse counter (int) for voxels (n_bins x n_bins x n_bins)

  • voxel_axon_dist – Axonal distance from soma to synapses (n_bins x n_bins x n_bins)

plot_hyper_voxel(plot_neurons=False, draw_axons=True, draw_dendrites=True, draw_axon_voxels=True, draw_dendrite_voxels=True, detect_done=True, elev_azim=None, show_axis=True, title=None, fig_file_name=None, dpi=300)[source]

Plot hyper voxel.

Parameters
  • plot_neurons – Should neurons be plotted

  • draw_axons – Draw axons

  • draw_dendrites – Draw dendrites

  • draw_axon_voxels – Draw axon voxels marked

  • draw_dendrite_voxels – Draw dendrite voxels marked

  • detect_done

  • elev_azim – View angle

  • show_axis – Show x,y,z axis?

  • title – Title of plot

  • fig_file_name – Fig name to save figure to

  • dpi – Resolution

plot_neurons_in_hyper_voxel(neuron_id, neuron_colour, axon_alpha=None, dend_alpha=None, show_plot=True, dpi=300)[source]

Plot neurons in hyper voltage

Parameters
  • neuron_id – ID of neurons to plot

  • neuron_colour – Colur of neurons to plot

  • axon_alpha – Alpha value of neuron axons

  • dend_alpha – Alpha value of neuron dendrites

  • show_plot – Should we dispaly the plot or keep it hidden

  • dpi – Resolution of output file

process_hyper_voxel(hyper_id)[source]

Process hyper voxel, ie do touch detection, and save results.

Parameters

hyper_id – ID of hyper voxel to process

read_neuron_positions(position_file)[source]

Loads neuron positions from network’s position_file.

Parameters

position_file – path to network position file (network-neuron-positions.hdf5)

read_prototypes(config_file=None, axon_stump_id_flag=False)[source]

Read in neuron prototypes. A neuron prototype can have multiple parameters, and morphology variations.

Parameters
  • config_file (str) – path to network config file

  • axon_stump_id_flag (bool) – Should segments be renumbered as if axon is replaced by axon stump

remove_empty(hyper_id)[source]

Removes empty hypervoxels from the list

Parameters

hyper_id (list) – List of hyper voxel IDs

Returns

List of remaining, and removed, hyper voxel ID

Return type

hyper_id_kept, hyper_id_removed (list, list)

resize_hyper_voxel_synapses_matrix(new_size=None)[source]

Increase the maximal size of the synapse matrix used for the hypervoxel.

Parameters

new_size (int) – Number of rows in synapse matrix

save_neuron_distribution_history(hyper_voxels, min_coord, max_coord)[source]

Save neuron distribution history to file.

Parameters
  • hyper_voxels (3D matrix with int) – Hypervoxel IDs, spatially arranged

  • min_coord (float,float,float) – minimal coordinate for all neurons/neurites in simulation

  • max_coord (float,float,float) – maximal coordinate for all neurons/neurites in simulation

setup_hyper_voxel(hyper_voxel_origo, hyper_voxel_id)[source]

Initialise all variables for a hypervoxel (containing NxNxN voxels) before touch detection.

Parameters
  • hyper_voxel_origo (float,float,float) – Origo of new hyper voxel

  • hyper_voxel_id (int) – ID of hypervoxel to process

setup_log(logfile_name=None)[source]

Initiates log file.

Parameters

logfile_name (str) – Path to log file

setup_parallel(d_view=None)[source]

Prepares workers for parallel execution if d_view is not None.

setup_process_hyper_voxel_state_history()[source]

Initialises the variables for tracking progress of touch detection.

Returns

The items returned are (all_hyper_id_list, num_completed,

remaining, voxel_overflow_counter). The remaining hypervoxel id list is sorted descending by size.

Return type

progress_data (list, int, list, int)

setup_work_history(work_history_file=None)[source]

Sets up work history. By logging progress we are able to restart an earlier partial touch detection run.

Parameters

work_history_file (str, optional) – Path to work history file

sort_gap_junctions()[source]

Sort gap junctions in self.hyper_voxel_gap_junctions. New sort order is columns 1 (dest), 0 (src).

sort_remaining_by_size(remaining)[source]

Sorts the remaining hypervoxel ID list by descending size (number of neurons in hypervoxel)

Parameters

remaining (list) – List of hypervoxels

Returns

Sorted list of hypervoxels

Return type

sorted_remaining (list)

sort_synapses()[source]

Sort synapses stored in self.hyper_voxel_synapses. New sort order is columns 1 (dest), 0 (src), 6 (synapse type).

update_process_hyper_voxel_state(hyper_id, num_syn, num_gj, exec_time, voxel_overflow_counter)[source]

Updates the process log with new hypervoxel state

Parameters
  • hyper_id (int) – Hypervoxel id completed

  • num_syn (int) – Number of synapses detected in hyper voxel

  • num_gj (int) – Number of gap junctions detected in hyper voxel

  • exec_time – Execution time (currently not used!)

  • voxel_overflow_counter – How many synapses/gap junctions did we miss due to memory overflow? (Should be 0)

write_hyper_voxel_to_hdf5()[source]

Saves hyper voxel synapses to data file.

write_log(text, flush=True, is_error=False, force_print=False)[source]

Writes to log file. Use setup_log first. Text is only written to screen if self.verbose=True, or is_error = True, or force_print = True.

test (str) : Text to write flush (bool) : Should all writes be flushed to disk directly? is_error (bool) : Is this an error, always written. force_print (bool) : Force printing, even if self.verbose=False.