Auryn simulator  v0.8.1-206-gb56e451
Plastic Spiking Neural Network Simulator
Classes | Typedefs | Enumerations | Functions | Variables
auryn Namespace Reference

Classes

class  ABSConnection
 
class  AdExGroup
 Conductance based Adaptive Exponential neuron model - Brette and Gerstner (2005). More...
 
class  AIF2Group
 An adaptive integrate and fire group comparable to AIFGroup but with two independent adaptation timescales. More...
 
class  AIFGroup
 A simple extension of IFGroup with spike triggered adaptation. More...
 
class  AllToAllConnection
 Provides all to all connectivity. More...
 
class  AuditoryBeepGroup
 A special Poisson generator that mimicks thalamo-cortical input to auditory cortex layer 3/4. More...
 
class  AurynConnectionAllocationException
 
class  AurynDelayTooSmallException
 
class  AurynDelayVector
 AurynDelayVector is a AurynVectorFloat which keeps its own history in a ring buffer. More...
 
class  AurynGenericException
 
class  AurynMatrixBufferException
 
class  AurynMatrixComplexStateException
 
class  AurynMatrixDimensionalityException
 
class  AurynMatrixPushBackException
 
class  AurynMemoryAlignmentException
 
class  AurynMMFileException
 
class  AurynOpenFileException
 
class  AurynSpikeAttributeSizeException
 
class  AurynStateVectorException
 
class  AurynTimeOverFlowException
 
class  AurynVector
 Default Auryn vector template. More...
 
class  AurynVectorDimensionalityException
 
class  AurynVectorFloat
 Default AurynVectorFloat class for performance computation. More...
 
class  AurynVersion
 Container class providing Auryn version number. More...
 
class  BinarySpikeMonitor
 The standard Monitor object to record spikes from a SpikingGroup and write them to a binary file. More...
 
class  BinaryStateMonitor
 Records from an arbitray state vector of one unit from the source SpikingGroup to a binary file. More...
 
class  BurstRateMonitor
 Monitor class to record population firing rates. More...
 
class  Checker
 The abstract base class for all checkers. More...
 
class  ComplexMatrix
 Template for a sparse matrix with row major ordering and fast access of rows and capability to store float values per matrix entry. More...
 
class  Connection
 The abstract base class for all Connection objects in Auryn. More...
 
class  CorrelatedPoissonGroup
 A PoissonGroup with multiple subpopulations that co-modulate their firing rate according to an Ornstein Uhlenbeck process. More...
 
class  CubaIFGroup
 Current based neuron model with absolute refractoriness as used in Vogels and Abbott 2005. More...
 
class  CurrentInjector
 Stimulator class to add values in each timestep to arbitrary neuronal states. More...
 
class  DelayConnection
 DelayConnection implements a SparseConnection with adjustable delays. More...
 
class  DelayedSpikeMonitor
 SpikeMonitor that reads the delayed spikes as they are received by a postsynaptic neuron. More...
 
class  Device
 Abstract base class for all Device, Stimulator, etc objects. More...
 
class  DuplexConnection
 Duplex connection is the base class of most plastic connections. More...
 
class  EulerTrace
 Solves a set of identical linear differential equations with the Euler method. It is used to implement synaptic traces in most STDP models. More...
 
class  ExpCobaSynapse
 Implements an exponential conductance-based synapse model. More...
 
class  ExpCubaSynapse
 Implements an exponential current-based synapse model. More...
 
class  FanOutConnection
 Provides all to all connectivity. More...
 
class  FileInputGroup
 Reads spikes from a ras file and emits them as SpikingGroup in a simulation. More...
 
class  FileModulatedPoissonGroup
 A special Poisson generator that reads its instantaneous firing rate from a tiser file. Datapoints in the rate file are interpolated linearly. More...
 
class  IafPscDeltaGroup
 Conductance based neuron model with absolute refractoriness as used in Vogels and Abbott 2005. More...
 
class  IafPscExpGroup
 Simple LIF neuron model with absolute refractoriness and current based synapses. More...
 
class  IdentityConnection
 Provides a unity matrix like connectivity. More...
 
class  IFGroup
 Implements the standard integrate and file model used in Auryn. More...
 
class  InputChannelGroup
 A PoissonGroup with multiple subpopulations that co-modulate their firing rate according to an Ornstein Uhlenbeck process. More...
 
class  IzhikevichGroup
 This NeuronGroup implements the Izhikevich neuron model with conductance based AMPA and GABA synapses. More...
 
class  LinearComboSynapse
 Implements Auryn's default conductance based AMPA, Combo synapse without Combo voltage dependence. More...
 
class  LinearTrace
 Exponential synaptic trace which exactly solves in an event-based manner in non-follow scenarios. More...
 
class  Logger
 A generic logger class that logs to screen and a log-file. More...
 
class  LPTripletConnection
 Implements triplet STDP in which weight updates are low-pass filtered. More...
 
class  MinimalTripletConnection
 Implements minimal triplet STDP as described by Pfister and Gerstner 2006. More...
 
class  ModSynIFGroup
 Implements the standard integrate and file model used in Auryn. More...
 
class  Monitor
 Abstract base class for all Monitor objects. More...
 
class  MovingBumpGroup
 A special PoissonGroup that generates jumping Gaussian bumps in the firing rate profile. More...
 
class  NaudGroup
 This file implements NaudGroup, Richard Naud's reduced two compartment model with active dendrites. More...
 
struct  neuron_pair
 Struct that defines a pair of neurons in SparseConnection. More...
 
class  NeuronGroup
 Abstract base class for all neuron groups. More...
 
class  NormalStimulator
 Stimulator class to inject timeseries of currents to patterns (subpopulations) of neurons. More...
 
class  PairInteractionConnection
 STDP Connection class to simulate arbitrary nearest-neighbor STDP windows. More...
 
class  ParrotGroup
 A SpikingGroup that copies the output of another source SpikingGroup. More...
 
struct  pattern_member
 Struct used to define neuronal assembly patterns in SparseConnection. More...
 
class  PatternMonitor
 Monitor class to record population firing rates. More...
 
class  PatternStimulator
 Stimulator class to inject timeseries of currents to patterns (subpopulations) of neurons. More...
 
class  PoissonGroup
 A SpikingGroup that creates poissonian spikes with a given rate. More...
 
class  PoissonSpikeInjector
 A PoissonGroup which can directly add its output spike to another SpikingGroup by piggy backing onto it. More...
 
class  PoissonStimulator
 Stimulator class to inject timeseries of currents NeuronGroups. More...
 
class  PopulationRateMonitor
 Monitor class to record population firing rates. More...
 
class  ProfilePoissonGroup
 A SpikingGroup that creates poissonian spikes with a given rate and spatial profile. More...
 
class  RateChecker
 A Checker class that tracks population firing rate as a moving average and breaks a run if it goes out of bound. More...
 
class  RateModulatedConnection
 Rate Modulated Connection implements a SparseConnection in which the weights depend. More...
 
class  RateMonitor
 Monitor class to record neural firing rates. More...
 
class  RealTimeMonitor
 Monitor class to record the system time in every timestep. More...
 
class  SimpleMatrix
 Template for a sparse matrix with row major ordering and fast access of rows. More...
 
class  SparseConnection
 The base class to create sparse random connections. More...
 
class  SpikeDelay
 Delay object for spikes which is synchronized between nodes using the SyncBuffer formalism implemented in System. More...
 
struct  SpikeEvent_type
 Auryn spike event for binary monitors. More...
 
class  SpikeMonitor
 The standard Monitor object to record spikes from a SpikingGroup and write them to a text file. More...
 
class  SpikeTimingStimGroup
 
class  SpikingGroup
 Abstract base class of all objects producing spikes. More...
 
class  SRM0Group
 Implements SRM0 neuron model with escape noise. More...
 
class  StateMonitor
 Records from an arbitray state vector of one unit from the source SpikingGroup to a file. More...
 
struct  StateValue_type
 Auryn spike event for binary monitors. More...
 
class  STDPConnection
 Double STDP All-to-All Connection. More...
 
class  STDPwdConnection
 Doublet STDP All-to-All as implemented in NEST as stdp_synapse_hom. More...
 
class  StimulusGroup
 Provides a poisson stimulus at random intervals in one or more predefined subsets of the group that are read from a file. More...
 
class  STPConnection
 This class implements short term plasticity according to the Tsodyks-Markram synapse. More...
 
class  StructuredPoissonGroup
 A special Poisson generator that can hide a fixed number of spatio-temporal patterns in the spike data. More...
 
class  SymmetricSTDPConnection
 Implements a symmetric STDP window with an optional presynaptic offset as used for inhibitory plasticity in Vogels et al. 2011. More...
 
class  SynapseModel
 Implements base class for modular synapse models. More...
 
class  SyncBuffer
 Buffer object to capsulate native MPI_Allgather for SpikingGroups. More...
 
class  System
 Class that implements system wide variables and methods to manage and run simulations. More...
 
class  TIFGroup
 Conductance based LIF neuron model with absolute refractoriness as used in Vogels and Abbott 2005. More...
 
class  Trace
 Abstract base class of synaptic traces. More...
 
class  TripletConnection
 Implements triplet STDP with metaplasticity as described by Pfister and Gerstner 2006. More...
 
class  TripletDecayConnection
 Implements triplet STDP with an exponential weight decay. More...
 
class  TripletScalingConnection
 
class  VoltageClampMonitor
 Implements a voltage clamp for one neuron and records the clamp current. More...
 
class  VoltageMonitor
 Records the membrane potential from one unit from the source neuron group to a file. More...
 
class  WeightChecker
 A Checker class that tracks the meain weight of a Connection and breaks a run if it goes out of bound. More...
 
class  WeightMatrixMonitor
 Saves the weight matrix of a given connection in regular time intervals. More...
 
class  WeightMonitor
 Monitors the evolution of a single or a set of weights. More...
 
class  WeightPatternMonitor
 Records mean weights from a connection specified by one or two pattern files. Can be used to easily monitor the mean synaptic weight in assemblies or feed-forward connections of populations of neurons. More...
 
class  WeightStatsMonitor
 Records mean and standard deviation of a weight matrix in predefined intervals. More...
 
class  WeightSumMonitor
 Records sum of all weights in synaptic weight matrix in predefined intervals. More...
 

Typedefs

typedef unsigned int NeuronID
 NeuronID is an unsigned integeger type used to index neurons in Auryn. More...
 
typedef NeuronID AurynInt
 
typedef unsigned int StateID
 StateID is an unsigned integeger type used to index synaptic states in Auryn. More...
 
typedef unsigned long AurynLong
 An unsigned long type used to count synapses or similar. More...
 
typedef NeuronID AurynTime
 Defines Auryns discrete time unit of the System clock. Change to AurynLong if 120h of simtime are not sufficient. More...
 
typedef std::string string
 Standard library string type which is imported into Auryn namespace. More...
 
typedef float AurynFloat
 Low precision floating point datatype. More...
 
typedef double AurynDouble
 Higher precision floating point datatype. More...
 
typedef AurynFloat AurynWeight
 Unit of synaptic weights. More...
 
typedef AurynFloat AurynState
 Type for Auryn state variables (default single precision since it needs to be compatible with auryn_vector_float). More...
 
typedef std::vector< NeuronIDSpikeContainer
 Spike container type. Used for storing spikes. More...
 
typedef std::vector< float > AttributeContainer
 Attribute container type. Used for storing spike attributes that are needed for efficient STP implementations. More...
 
typedef std::vector< pattern_membertype_pattern
 
typedef AurynVectorFloat AurynStateVector
 Defines AurynStateVector type as synonymous to AurynVectorFloat. More...
 
typedef AurynVector< AurynWeight, AurynLongAurynSynStateVector
 Defines AurynSynStateVector for synaptic states. More...
 
typedef AurynStateVector auryn_vector_float
 Legacy definition of AurynStateVector. More...
 
typedef AurynVector< unsigned short, NeuronIDauryn_vector_ushort
 Legacy definition of AurynVector<unsigned short> More...
 
typedef SimpleMatrix< AurynWeight * > BackwardMatrix
 
typedef ComplexMatrix< AurynWeightForwardMatrix
 

Enumerations

enum  TransmitterType {
  GLUT, GABA, AMPA, NMDA,
  MEM, CURSYN
}
 Specifies the different transmitter types that Auryn knows. More...
 
enum  StimulusGroupModeType {
  MANUAL, RANDOM, SEQUENTIAL, SEQUENTIAL_REV,
  STIMFILE
}
 Specifies stimulus order used in StimulusGroup. More...
 
enum  LogMessageType {
  EVERYTHING, VERBOSE, INFO, NOTIFICATION,
  SETTINGS, PROGRESS, WARNING, ERROR,
  NONE
}
 Enum type for significance level of a given message send to the Logger. More...
 
enum  NodeDistributionMode { AUTO, ROUNDROBIN, BLOCKLOCK, RANKLOCK }
 Specifies howto distribute different neurons across ranks when simulation is run in parallel. More...
 
enum  RecordingMode { SINGLE, DATARANGE, ELEMENTLIST, GROUPS }
 
enum  PatternMode { ALLTOALL, ASSEMBLIES_ONLY }
 

Functions

int auryn_AlignOffset (const int N, const void *vp, const int inc, const int align)
 Determines memory alignment (adapted from ATLAS library) More...
 
NeuronID calculate_vector_size (NeuronID i)
 Rounds vector size to multiple of four to allow using the SSE optimizations. More...
 
void auryn_vector_float_mul (auryn_vector_float *a, auryn_vector_float *b)
 
void auryn_vector_float_add_constant (auryn_vector_float *a, float b)
 Computes a := a + b. More...
 
void auryn_vector_float_scale (const float a, auryn_vector_float *b)
 
void auryn_vector_float_saxpy (const float a, auryn_vector_float *x, auryn_vector_float *y)
 
void auryn_vector_float_add (auryn_vector_float *a, auryn_vector_float *b)
 Internal version of to add GSL vectors. More...
 
void auryn_vector_float_sub (auryn_vector_float *a, auryn_vector_float *b)
 Computes a := a-b. More...
 
void auryn_vector_float_sub (auryn_vector_float *a, auryn_vector_float *b, auryn_vector_float *r)
 Computes r := a-b. More...
 
void auryn_vector_float_clip (auryn_vector_float *v, const float a, const float b)
 
void auryn_vector_float_clip (auryn_vector_float *v, const float a)
 
auryn_vector_floatauryn_vector_float_alloc (const NeuronID n)
 
void auryn_vector_float_free (auryn_vector_float *v)
 
void auryn_vector_float_set_all (auryn_vector_float *v, AurynFloat x)
 
void auryn_vector_float_set_zero (auryn_vector_float *v)
 
AurynFloat auryn_vector_float_get (auryn_vector_float *v, const NeuronID i)
 
AurynFloatauryn_vector_float_ptr (auryn_vector_float *v, const NeuronID i)
 
void auryn_vector_float_set (auryn_vector_float *v, const NeuronID i, AurynFloat x)
 
void auryn_vector_float_copy (auryn_vector_float *src, auryn_vector_float *dst)
 Copies vector src to dst assuming they have the same size. More...
 
auryn_vector_ushortauryn_vector_ushort_alloc (const NeuronID n)
 
void auryn_vector_ushort_free (auryn_vector_ushort *v)
 
void auryn_vector_ushort_set_all (auryn_vector_ushort *v, unsigned short x)
 
void auryn_vector_ushort_set_zero (auryn_vector_ushort *v)
 
unsigned short auryn_vector_ushort_get (auryn_vector_ushort *v, const NeuronID i)
 
unsigned short * auryn_vector_ushort_ptr (auryn_vector_ushort *v, const NeuronID i)
 
void auryn_vector_ushort_set (auryn_vector_ushort *v, const NeuronID i, unsigned short x)
 
void auryn_vector_ushort_copy (auryn_vector_ushort *src, auryn_vector_ushort *dst)
 
void auryn_env_init (int ac, char *av[], string dir=".", string logfile_prefix="", LogMessageType filelog_level=NOTIFICATION, LogMessageType consolelog_level=PROGRESS)
 Initalizes Auryn base environment (used internally) More...
 
void auryn_kernel_init (string dir=".", string simulation_name="default")
 Initalizes the Auryn kernel (used internally) More...
 
void auryn_init (int ac, char *av[], string dir=".", string simulation_name="default", string logfile_prefix="", LogMessageType filelog_level=NOTIFICATION, LogMessageType consolelog_level=PROGRESS)
 Initalizes MPI and the Auryn simulation environment. More...
 
void auryn_kernel_free ()
 Frees the current auryn kernel (used interally) More...
 
void auryn_env_free ()
 Frees logger and MPI. More...
 
void auryn_free ()
 Cleanly shuts down Auryn simulation environment. More...
 
void auryn_abort (int errcode=0)
 Terminates Auryn simulation abnormally. More...
 

Variables

double auryn_timestep = 1e-4
 Simulation timestep in seconds. More...
 
mpi::communicator * mpicommunicator
 Global pointer to instance of mpi::mpicommunicator which needs to be initialized in every simulation main program. More...
 
mpi::environment * mpienv
 Global pointer to instance of mpi::environment which needs to be initialized in every simulation main program. More...
 
Loggerlogger
 Global pointer to instance of Logger which needs to be initialized in every simulation main program. More...
 
Systemsys
 Global pointer to instance of System which needs to be initialized in every simulation main program. More...
 

Typedef Documentation

◆ AttributeContainer

typedef std::vector<float> auryn::AttributeContainer

Attribute container type. Used for storing spike attributes that are needed for efficient STP implementations.

◆ auryn_vector_float

Legacy definition of AurynStateVector.

Default legacy Auryn state vector type

◆ auryn_vector_ushort

Legacy definition of AurynVector<unsigned short>

Default legacy Auryn ushort vector type

◆ AurynDouble

typedef double auryn::AurynDouble

Higher precision floating point datatype.

◆ AurynFloat

typedef float auryn::AurynFloat

Low precision floating point datatype.

◆ AurynInt

◆ AurynLong

typedef unsigned long auryn::AurynLong

An unsigned long type used to count synapses or similar.

◆ AurynState

Type for Auryn state variables (default single precision since it needs to be compatible with auryn_vector_float).

◆ AurynStateVector

Defines AurynStateVector type as synonymous to AurynVectorFloat.

Auryn state vectors are used to implement vectorized code for SpikingGroup and NeuronGroup. An AurynStateVector in a SpikingGroup typically has the local rank size of that group and each neuronal state variable corresponds to a state vector that houses this state for all neurons on that rank. AurynStateVectors are defined as AurynVectorFloat. This typically needs to change when AurynState or AurynFloat types are changed.

◆ AurynSynStateVector

Defines AurynSynStateVector for synaptic states.

◆ AurynTime

Defines Auryns discrete time unit of the System clock. Change to AurynLong if 120h of simtime are not sufficient.

◆ AurynWeight

Unit of synaptic weights.

◆ BackwardMatrix

Definition of BackwardMatrix - a sparsematrix of pointers to weight values.

◆ ForwardMatrix

◆ NeuronID

typedef unsigned int auryn::NeuronID

NeuronID is an unsigned integeger type used to index neurons in Auryn.

◆ SpikeContainer

typedef std::vector<NeuronID> auryn::SpikeContainer

Spike container type. Used for storing spikes.

◆ StateID

typedef unsigned int auryn::StateID

StateID is an unsigned integeger type used to index synaptic states in Auryn.

◆ string

typedef std::string auryn::string

Standard library string type which is imported into Auryn namespace.

◆ type_pattern

typedef std::vector<pattern_member> auryn::type_pattern

Enumeration Type Documentation

◆ LogMessageType

Enum type for significance level of a given message send to the Logger.

Enumerator
EVERYTHING 
VERBOSE 
INFO 
NOTIFICATION 
SETTINGS 
PROGRESS 
WARNING 
ERROR 
NONE 
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41
Definition: Logger.h:41

◆ NodeDistributionMode

Specifies howto distribute different neurons across ranks when simulation is run in parallel.

Enumerator
AUTO 

Tries to make a smart choice.

ROUNDROBIN 

Default mode of distribution.

BLOCKLOCK 

Tries to implement block lock.

RANKLOCK 

Locks to single rank (this is a special case of BLOCKLOCK.

50  {
51  AUTO,
52  ROUNDROBIN,
53  BLOCKLOCK,
54  RANKLOCK
55 };
Locks to single rank (this is a special case of BLOCKLOCK.
Definition: SpikingGroup.h:54
Tries to make a smart choice.
Definition: SpikingGroup.h:51
Default mode of distribution.
Definition: SpikingGroup.h:52
Tries to implement block lock.
Definition: SpikingGroup.h:53

◆ PatternMode

Determines how pattern files are interpreted for loading. ALLTOALL will add connections between all possible pattern combinations, whereas ASSEMBLIES_ONLY restricts recording to inside each pattern.

Enumerator
ALLTOALL 
ASSEMBLIES_ONLY 
Definition: WeightMonitor.h:55
Definition: WeightMonitor.h:55

◆ RecordingMode

RecordingMode determines the default recording behavior of the monitor. The modes SINGLE, DATARANGE and ELEMENTLIST (default) record weight values from single synapses while GROUPS records the statistics over groups of synapses

Enumerator
SINGLE 

The entire Monitor will record from a single synapse specified at initialization.

DATARANGE 

The Monitor will record from a range of synapses specified at initialization.

ELEMENTLIST 

The Monitor records from selected synapses stored in a list. This is the default behavior.

GROUPS 

This mode is added in versions >0.4.1 and allows to record summary statistics of synapses between neural groups/patterns.

44  {
45  SINGLE,
46  DATARANGE,
47  ELEMENTLIST,
48  GROUPS
50 };
The entire Monitor will record from a single synapse specified at initialization. ...
Definition: WeightMonitor.h:45
The Monitor will record from a range of synapses specified at initialization.
Definition: WeightMonitor.h:46
The Monitor records from selected synapses stored in a list. This is the default behavior.
Definition: WeightMonitor.h:47
Definition: WeightMonitor.h:48

◆ StimulusGroupModeType

Specifies stimulus order used in StimulusGroup.

Enumerator
MANUAL 
RANDOM 
SEQUENTIAL 
SEQUENTIAL_REV 
STIMFILE 
Definition: auryn_definitions.h:148
Definition: auryn_definitions.h:148
Definition: auryn_definitions.h:148
Definition: auryn_definitions.h:148
Definition: auryn_definitions.h:148

◆ TransmitterType

Specifies the different transmitter types that Auryn knows.

Enumerator
GLUT 

Standard Glutamatergic (excitatory) transmission.

GABA 

Standard Gabaergic (inhibitory) transmission.

AMPA 

Only targets AMPA channels.

NMDA 

Only targets NMDA.

MEM 

Current based synapse. Adds the transmitted quantity directly to membrane voltage.

CURSYN 

Current based synapse with dynamics.

138  {
139  GLUT,
140  GABA,
141  AMPA,
142  NMDA,
143  MEM,
144  CURSYN
145  };
Current based synapse. Adds the transmitted quantity directly to membrane voltage.
Definition: auryn_definitions.h:143
Standard Glutamatergic (excitatory) transmission.
Definition: auryn_definitions.h:139
Only targets NMDA.
Definition: auryn_definitions.h:142
Current based synapse with dynamics.
Definition: auryn_definitions.h:144
Standard Gabaergic (inhibitory) transmission.
Definition: auryn_definitions.h:140
Only targets AMPA channels.
Definition: auryn_definitions.h:141

Function Documentation

◆ auryn_abort()

void auryn::auryn_abort ( int  errcode = 0)

Terminates Auryn simulation abnormally.

This issues a term signal to all MPI processes in case of error, but first closes the Auryn kernel and terminates the logger to ensure all information of the issuing rank are written to disk.

Parameters
errcodeThe errorcode to be returned by the MPI processes
113  {
114  delete sys;
115  delete logger;
116 #ifdef AURYN_CODE_USE_MPI
117  mpienv->abort(errcode);
118 #endif // AURYN_CODE_USE_MPI
119  // In the MPI case the above line should have killed this process already ...
120  std::exit(errcode);
121  }
mpi::environment * mpienv
Global pointer to instance of mpi::environment which needs to be initialized in every simulation main...
Definition: auryn_global.cpp:33
Logger * logger
Global pointer to instance of Logger which needs to be initialized in every simulation main program...
Definition: auryn_global.cpp:36
System * sys
Global pointer to instance of System which needs to be initialized in every simulation main program...
Definition: auryn_global.cpp:37

◆ auryn_AlignOffset()

int auryn::auryn_AlignOffset ( const int  N,
const void *  vp,
const int  inc,
const int  align 
)

Determines memory alignment (adapted from ATLAS library)

Parameters
Nmax return value
*vpPointer to be aligned
incsize of element, in bytes
alignrequired alignment, in bytes
38 {
39  const int p = align/inc;
40  const size_t k=(size_t)vp, j=k/inc;
41  int iret;
42  if (k == (j)*inc && p*inc == align)
43  {
44  iret = ((j+p-1) / p)*p - j;
45  if (iret <= N) return(iret);
46  }
47  return(N);
48 }
#define N
Definition: sim_delay_connection.cpp:23

◆ auryn_env_free()

void auryn::auryn_env_free ( )

Frees logger and MPI.

This function frees the MPI modules and logger. It is usually called by auryn_free()

See also
auryn_kernel_init
99  {
100  delete logger;
101 #ifdef AURYN_CODE_USE_MPI
102  delete mpicommunicator;
103  delete mpienv;
104 #endif // AURYN_CODE_USE_MPI
105  }
mpi::environment * mpienv
Global pointer to instance of mpi::environment which needs to be initialized in every simulation main...
Definition: auryn_global.cpp:33
mpi::communicator * mpicommunicator
Global pointer to instance of mpi::mpicommunicator which needs to be initialized in every simulation ...
Definition: auryn_global.cpp:32
Logger * logger
Global pointer to instance of Logger which needs to be initialized in every simulation main program...
Definition: auryn_global.cpp:36

◆ auryn_env_init()

void auryn::auryn_env_init ( int  ac,
char *  av[],
string  dir = ".",
string  logfile_prefix = "",
LogMessageType  filelog_level = NOTIFICATION,
LogMessageType  consolelog_level = PROGRESS 
)

Initalizes Auryn base environment (used internally)

This function is called as port of auryn_init and is typically called internally. It instantiates the global Logger object and initializes the MPI environment.

41  {
42 #ifdef AURYN_CODE_USE_MPI
43  // init MPI environment
44  mpienv = new mpi::environment(ac, av);
45  mpicommunicator = new mpi::communicator();
46  const unsigned int local_rank = mpicommunicator->rank();
47 #else
48  const unsigned int local_rank = 0;
49 #endif // AURYN_CODE_USE_MPI
50 
51  // Init logger environment
52  try
53  {
54  string log_prefix_ = boost::filesystem::basename(av[0]);
55  log_prefix_.erase(std::remove(log_prefix_.begin(),log_prefix_.end(),' '),log_prefix_.end()); // remove spaces
56  std::transform(log_prefix_.begin(), log_prefix_.end(), log_prefix_.begin(), ::tolower); // convert to lower case
57 
58  if ( !logfile_prefix.empty() ) log_prefix_ = logfile_prefix;
59 
60  char strbuf_tmp [255];
61  sprintf(strbuf_tmp, "%s/%s.%d.log", dir.c_str(), log_prefix_.c_str(), local_rank);
62  string auryn_simulation_logfile = strbuf_tmp;
63  logger = new Logger(auryn_simulation_logfile, local_rank, consolelog_level, filelog_level);
64  }
65  catch ( AurynOpenFileException excpt )
66  {
67  std::cerr << "Cannot proceed without log file. Exiting all ranks ..." << std::endl;
68  auryn_abort(10);
69  }
70 
71  }
mpi::environment * mpienv
Global pointer to instance of mpi::environment which needs to be initialized in every simulation main...
Definition: auryn_global.cpp:33
mpi::communicator * mpicommunicator
Global pointer to instance of mpi::mpicommunicator which needs to be initialized in every simulation ...
Definition: auryn_global.cpp:32
Logger * logger
Global pointer to instance of Logger which needs to be initialized in every simulation main program...
Definition: auryn_global.cpp:36
void auryn_abort(int errcode)
Terminates Auryn simulation abnormally.
Definition: auryn_global.cpp:113
Here is the call graph for this function:

◆ auryn_free()

void auryn::auryn_free ( )

Cleanly shuts down Auryn simulation environment.

Deletes global variables mpienv, mpicommunicator, sys and logger and ensures that all data is written to disk.

108  {
110  auryn_env_free();
111  }
void auryn_kernel_free()
Frees the current auryn kernel (used interally)
Definition: auryn_global.cpp:93
void auryn_env_free()
Frees logger and MPI.
Definition: auryn_global.cpp:98
Here is the call graph for this function:

◆ auryn_init()

void auryn::auryn_init ( int  ac,
char *  av[],
string  dir = ".",
string  simulation_name = "default",
string  logfile_prefix = "",
LogMessageType  filelog_level = NOTIFICATION,
LogMessageType  consolelog_level = PROGRESS 
)

Initalizes MPI and the Auryn simulation environment.

This function has to be called only once and before any Auryn class is used. It is thus typically invoked in the main function of your simulation program. The function initalizes up the global objects mpienv, mpicommunicator, sys and logger. It takes the command line parameters from your main method (needed to initializes the MPI environment) as well as some additional arguments. You can for instance pass a default output directory which can be used by Device instances or any class which writes to disk. Moreover, you can pass a simulation_name to auryn_init which will appear in logging output. Finally, you can explicitly pass a logfile_prefix string which Auryn will decorate with the mpi rank number and a "log" extension. If this option is omitted, Auryn will derive a log file name from the name of your simulation binary. At the end of your simulation, make sure to call auryn_free() to cleanly terminate Auryn and avoid data loss.

Parameters
acNumber of command line parameters argc passed to main
avCommand line parameters passed as argv to main
dirThe default output directory for files generated by your simulation
simulation_nameA name for your simulation which will appear in logfiles
logfile_prefixA file prefix (without path) which Auryn will use to generate a log file name.
85  {
86  // Init MPI and Logger
87  auryn_env_init(ac, av, dir, logfile_prefix, filelog_level, consolelog_level);
88  // Init Auryn Kernel
89  auryn_kernel_init(dir, simulation_name);
90  }
void auryn_kernel_init(string dir, string simulation_name)
Initalizes the Auryn kernel (used internally)
Definition: auryn_global.cpp:73
void auryn_env_init(int ac, char *av[], string dir, string logfile_prefix, LogMessageType filelog_level, LogMessageType consolelog_level)
Initalizes Auryn base environment (used internally)
Definition: auryn_global.cpp:40
Here is the call graph for this function:

◆ auryn_kernel_free()

void auryn::auryn_kernel_free ( )

Frees the current auryn kernel (used interally)

This function frees the current auryn kernel. It is usually called by auryn_free(), but in some cases it may make sense to use this independently.

See also
auryn_kernel_init
94  {
95  delete sys;
96  }
System * sys
Global pointer to instance of System which needs to be initialized in every simulation main program...
Definition: auryn_global.cpp:37

◆ auryn_kernel_init()

void auryn::auryn_kernel_init ( string  dir = ".",
string  simulation_name = "default" 
)

Initalizes the Auryn kernel (used internally)

This function is called as port of auryn_init and is typically called internally. However, it can be used if several simulations are to be run from a single simulation binary. Then manual calles of auryn_kernel_init and auryn_kernel_free can be used to reinitalize the kernel without closing down the MPI environment.

74  {
75 #ifdef AURYN_CODE_USE_MPI
76  auryn::sys = new System(mpicommunicator);
77 #else
78  auryn::sys = new System();
79 #endif // AURYN_CODE_USE_MPI
80  sys->set_output_dir(dir);
81  sys->set_simulation_name(simulation_name);
82  }
void set_output_dir(std::string path)
Set output dir for fn function.
Definition: System.cpp:664
mpi::communicator * mpicommunicator
Global pointer to instance of mpi::mpicommunicator which needs to be initialized in every simulation ...
Definition: auryn_global.cpp:32
void set_simulation_name(std::string name)
Sets the simulation name.
Definition: System.cpp:654
System * sys
Global pointer to instance of System which needs to be initialized in every simulation main program...
Definition: auryn_global.cpp:37
Here is the call graph for this function:

◆ auryn_vector_float_add()

void auryn::auryn_vector_float_add ( auryn_vector_float a,
auryn_vector_float b 
)

Internal version of to add GSL vectors.

Add vectors a and b and store the result in a.

85 {
86  a->add(b);
87 }
Here is the call graph for this function:

◆ auryn_vector_float_add_constant()

void auryn::auryn_vector_float_add_constant ( auryn_vector_float a,
float  b 
)

Computes a := a + b.

Internal version of auryn_vector_float_add between a constant and a vector

70 {
71  a->add(b);
72 }
Here is the call graph for this function:

◆ auryn_vector_float_alloc()

auryn_vector_float * auryn::auryn_vector_float_alloc ( const NeuronID  n)

Allocates an auryn_vector_float

107  {
108  return new auryn_vector_float(n);
109 }
int n
Definition: mkpat.py:5
AurynStateVector auryn_vector_float
Legacy definition of AurynStateVector.
Definition: auryn_definitions.h:345

◆ auryn_vector_float_clip() [1/2]

void auryn::auryn_vector_float_clip ( auryn_vector_float v,
const float  a,
const float  b 
)

Internal version to clip all the elements of a vector between [a:b]

99  {
100  v->clip(a, b);
101 }
Here is the call graph for this function:

◆ auryn_vector_float_clip() [2/2]

void auryn::auryn_vector_float_clip ( auryn_vector_float v,
const float  a 
)

Internal version to clip all the elements of a vector between [a:0]

103  {
104  v->clip(a,0.0);
105 }
Here is the call graph for this function:

◆ auryn_vector_float_copy()

void auryn::auryn_vector_float_copy ( auryn_vector_float src,
auryn_vector_float dst 
)

Copies vector src to dst assuming they have the same size.

Otherwise this will lead to undefined results. No checking of size is performed for performance reasons.

135  {
136  dst->copy(src);
137 }
Here is the call graph for this function:

◆ auryn_vector_float_free()

void auryn::auryn_vector_float_free ( auryn_vector_float v)

Frees an auryn_vector_float

111  {
112  delete v;
113 }

◆ auryn_vector_float_get()

AurynFloat auryn::auryn_vector_float_get ( auryn_vector_float v,
const NeuronID  i 
)

Auryn vector getter

123  {
124  return v->data[i];
125 }

◆ auryn_vector_float_mul()

void auryn::auryn_vector_float_mul ( auryn_vector_float a,
auryn_vector_float b 
)

Internal version of auryn_vector_float_mul of gsl operations

65 {
66  a->mul(b);
67 }
Here is the call graph for this function:

◆ auryn_vector_float_ptr()

AurynFloat * auryn::auryn_vector_float_ptr ( auryn_vector_float v,
const NeuronID  i 
)

Auryn vector gets pointer to designed element.

127  {
128  return v->data+i;
129 }

◆ auryn_vector_float_saxpy()

void auryn::auryn_vector_float_saxpy ( const float  a,
auryn_vector_float x,
auryn_vector_float y 
)

Computes y := a*x+y

Internal SAXPY version

80 {
81  y->saxpy(a,x);
82 }
Here is the call graph for this function:

◆ auryn_vector_float_scale()

void auryn::auryn_vector_float_scale ( const float  a,
auryn_vector_float b 
)

Internal version to scale a vector with a constant b

75 {
76  b->scale(a);
77 }
Here is the call graph for this function:

◆ auryn_vector_float_set()

void auryn::auryn_vector_float_set ( auryn_vector_float v,
const NeuronID  i,
AurynFloat  x 
)

Auryn vector setter

131  {
132  v->data[i] = x;
133 }

◆ auryn_vector_float_set_all()

void auryn::auryn_vector_float_set_all ( auryn_vector_float v,
AurynFloat  x 
)

Sets all elements in an auryn_vector_float to value x

115  {
116  v->set_all(x);
117 }
Here is the call graph for this function:

◆ auryn_vector_float_set_zero()

void auryn::auryn_vector_float_set_zero ( auryn_vector_float v)

Initializes an auryn_vector_float with zeros

119  {
120  v->set_zero();
121 }
Here is the call graph for this function:

◆ auryn_vector_float_sub() [1/2]

void auryn::auryn_vector_float_sub ( auryn_vector_float a,
auryn_vector_float b 
)

Computes a := a-b.

Internal version of to subtract GSL vectors.

90 {
91  a->sub(b);
92 }
Here is the call graph for this function:

◆ auryn_vector_float_sub() [2/2]

void auryn::auryn_vector_float_sub ( auryn_vector_float a,
auryn_vector_float b,
auryn_vector_float r 
)

Computes r := a-b.

95 {
96  r->diff(a,b);
97 }
r
Definition: mkpat.py:9
Here is the call graph for this function:

◆ auryn_vector_ushort_alloc()

auryn_vector_ushort * auryn::auryn_vector_ushort_alloc ( const NeuronID  n)

Allocates an auryn_vector_ushort

140  {
141  return new auryn_vector_ushort(n);
142 }
int n
Definition: mkpat.py:5
AurynVector< unsigned short, NeuronID > auryn_vector_ushort
Legacy definition of AurynVector<unsigned short>
Definition: auryn_definitions.h:348

◆ auryn_vector_ushort_copy()

void auryn::auryn_vector_ushort_copy ( auryn_vector_ushort src,
auryn_vector_ushort dst 
)

Copies vector src to dst assuming they have the same size. Otherwise this will lead to undefined results. No checking of size is performed for performance reasons.

168  {
169  for ( NeuronID i = 0 ; i < dst->size ; ++i )
170  dst->data[i] = src->data[i];
171 }
unsigned int NeuronID
NeuronID is an unsigned integeger type used to index neurons in Auryn.
Definition: auryn_definitions.h:151

◆ auryn_vector_ushort_free()

void auryn::auryn_vector_ushort_free ( auryn_vector_ushort v)

Frees an auryn_vector_ushort

144  {
145  delete v;
146 }

◆ auryn_vector_ushort_get()

unsigned short auryn::auryn_vector_ushort_get ( auryn_vector_ushort v,
const NeuronID  i 
)

Auryn vector getter

156  {
157  return v->data[i];
158 }

◆ auryn_vector_ushort_ptr()

unsigned short * auryn::auryn_vector_ushort_ptr ( auryn_vector_ushort v,
const NeuronID  i 
)

Auryn vector gets pointer to designed element.

160  {
161  return v->data+i;
162 }

◆ auryn_vector_ushort_set()

void auryn::auryn_vector_ushort_set ( auryn_vector_ushort v,
const NeuronID  i,
unsigned short  x 
)

Auryn vector setter

164  {
165  v->data[i] = x;
166 }

◆ auryn_vector_ushort_set_all()

void auryn::auryn_vector_ushort_set_all ( auryn_vector_ushort v,
unsigned short  x 
)

Sets all elements in an auryn_vector_ushort to value x

148  {
149  v->set_all(x);
150 }
Here is the call graph for this function:

◆ auryn_vector_ushort_set_zero()

void auryn::auryn_vector_ushort_set_zero ( auryn_vector_ushort v)

Initializes an auryn_vector_ushort with zeros

152  {
153  v->set_zero();
154 }
Here is the call graph for this function:

◆ calculate_vector_size()

NeuronID auryn::calculate_vector_size ( NeuronID  i)

Rounds vector size to multiple of four to allow using the SSE optimizations.

51 {
52  #if defined(CODE_USE_SIMD_INSTRUCTIONS_EXPLICITLY)
54  return i;
55  NeuronID div = i/SIMD_NUM_OF_PARALLEL_FLOAT_OPERATIONS; // rounds down
56  NeuronID new_size = (div+1)*SIMD_NUM_OF_PARALLEL_FLOAT_OPERATIONS; // is multiple of SIMD...
57  return new_size;
58  #else
59  return i;
60  #endif
61 }
#define SIMD_NUM_OF_PARALLEL_FLOAT_OPERATIONS
Use Intel Cilk Plus – only has an effect when CODE_USE_SIMD_INSTRUCTIONS_EXPLICITLY is enabled...
Definition: auryn_definitions.h:64
unsigned int NeuronID
NeuronID is an unsigned integeger type used to index neurons in Auryn.
Definition: auryn_definitions.h:151

Variable Documentation

◆ auryn_timestep

double auryn::auryn_timestep = 1e-4

Simulation timestep in seconds.

Todo:
For consistency the global variable for auryn_timestep should be moved to either auryn_definitions or into the System class

◆ logger

Logger * auryn::logger

Global pointer to instance of Logger which needs to be initialized in every simulation main program.

◆ mpicommunicator

mpi::communicator * auryn::mpicommunicator

Global pointer to instance of mpi::mpicommunicator which needs to be initialized in every simulation main program.

◆ mpienv

mpi::environment * auryn::mpienv

Global pointer to instance of mpi::environment which needs to be initialized in every simulation main program.

◆ sys

System * auryn::sys

Global pointer to instance of System which needs to be initialized in every simulation main program.