CARLsim  5.0.0
CARLsim: a GPU-accelerated SNN simulator
poisson_rate.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2014 Regents of the University of California. All rights reserved.
3  *
4  * Redistribution and use in source and binary forms, with or without
5  * modification, are permitted provided that the following conditions
6  * are met:
7  *
8  * 1. Redistributions of source code must retain the above copyright
9  * notice, this list of conditions and the following disclaimer.
10  *
11  * 2. Redistributions in binary form must reproduce the above copyright
12  * notice, this list of conditions and the following disclaimer in the
13  * documentation and/or other materials provided with the distribution.
14  *
15  * 3. The names of its contributors may not be used to endorse or promote
16  * products derived from this software without specific prior written
17  * permission.
18  *
19  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
20  * "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
21  * LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
22  * A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR
23  * CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
24  * EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
25  * PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
26  * PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
27  * LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
28  * NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
29  * SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
30  *
31  * *********************************************************************************************** *
32  * CARLsim
33  * created by: (MDR) Micah Richert, (JN) Jayram M. Nageswaran
34  * maintained by: (MA) Mike Avery <averym@uci.edu>, (MB) Michael Beyeler <mbeyeler@uci.edu>,
35  * (KDC) Kristofor Carlson <kdcarlso@uci.edu>
36  * (TSC) Ting-Shuo Chou <tingshuc@uci.edu>
37  *
38  * CARLsim available from http://socsci.uci.edu/~jkrichma/CARLsim/
39  * Ver 11/26/2014
40  */
41 
42 #ifndef _POISSON_RATE_H_
43 #define _POISSON_RATE_H_
44 
45 #include <vector>
46 
84 class PoissonRate {
85 public:
94  PoissonRate(int nNeur, bool onGPU=false);
95 
102  ~PoissonRate();
103 
112  int getNumNeurons();
113 
123  float getRate(int neurId);
124 
132  std::vector<float> getRates();
133 
143  float* getRatePtrCPU();
144 
154  float* getRatePtrGPU();
155 
163  bool isOnGPU();
164 
174  void setRate(int neurId, float rate);
175 
184  void setRates(float rate);
185 
194  void setRates(const std::vector<float>& rates);
195 
196 
197 private:
198  // This class provides a pImpl for the CARLsim User API.
199  // \see https://marcmutz.wordpress.com/translated-articles/pimp-my-pimpl/
200  class Impl;
201  Impl* _impl;
202 };
203 
204 #endif
PoissonRate::setRates
void setRates(float rate)
Assigns the same mean firing rate to all neurons.
Definition: poisson_rate.cpp:229
PoissonRate::getRatePtrCPU
float * getRatePtrCPU()
Returns pointer to CPU-allocated firing rate array (deprecated)
Definition: poisson_rate.cpp:225
PoissonRate
Class for generating Poisson spike trains.
Definition: poisson_rate.h:84
PoissonRate::getNumNeurons
int getNumNeurons()
Returns the number of neurons for which to generate Poisson spike trains.
Definition: poisson_rate.cpp:222
PoissonRate::getRate
float getRate(int neurId)
Returns the mean firing rate of a specific neuron ID.
Definition: poisson_rate.cpp:223
PoissonRate::~PoissonRate
~PoissonRate()
PoissonRate destructor.
Definition: poisson_rate.cpp:220
PoissonRate::setRate
void setRate(int neurId, float rate)
Sets the mean firing rate of a particular neuron ID.
Definition: poisson_rate.cpp:228
PoissonRate::PoissonRate
PoissonRate(int nNeur, bool onGPU=false)
PoissonRate constructor.
Definition: poisson_rate.cpp:219
PoissonRate::Impl
Definition: poisson_rate.cpp:54
PoissonRate::getRates
std::vector< float > getRates()
Returns a vector of firing rates, one element per neuron.
Definition: poisson_rate.cpp:224
PoissonRate::getRatePtrGPU
float * getRatePtrGPU()
Returns pointer to GPU-allocated firing rate array (deprecated)
Definition: poisson_rate.cpp:226
PoissonRate::isOnGPU
bool isOnGPU()
Checks whether the firing rates are allocated on CPU or GPU.
Definition: poisson_rate.cpp:227