Stimfit 0.12.7
Classes | Typedefs | Enumerations | Functions | Variables
Generic stimfit classes and functions

Classes

class  Channel
 A Channel contains several data Sections representing observations of the same physical quantity. More...
class  Recording
 Represents the data within a file. More...
class  Section
 Represents a continuously sampled sweep of data points. More...
class  stf::wxProgressInfo
 Progress Info interface adapter; maps to wxProgressDialog. More...
class  stf::Table
 A table used for printing information. More...
struct  stf::UserInput
 Represents user input from dialogs that can be used in plugins. More...
struct  stf::Plugin
 User-defined plugin. More...
struct  stf::Extension
 User-defined Python extension. More...
struct  stf::ifstreamMan
 Resource manager for ifstream objects. More...
struct  stf::ofstreamMan
 Resource manager for ofstream objects. More...
class  stf::Event
 Describes the attributes of an event. More...
struct  stf::PyMarker
 A marker that can be set from Python. More...
struct  stf::SectionAttributes
struct  stf::SectionPointer
struct  stf::parInfo
 Information about parameters used in storedFunc. More...
struct  stf::storedFunc
 Function used for least-squares fitting. More...

Typedefs

typedef boost::function
< Recording(const Recording
&, const Vector_double
&, std::map< std::string,
double > &)> 
stf::PluginFunc
 Get a Recording, do something with it, return the new Recording.
typedef boost::function
< double(double, const
Vector_double &)> 
stf::Func
 A function taking a double and a vector and returning a double.
typedef boost::function
< Vector_double(double, const
Vector_double &)> 
stf::Jac
 The jacobian of a stf::Func.
typedef boost::function
< double(double, double,
double, double, double)> 
stf::Scale
 Scaling function for fit parameters.
typedef boost::function< Table(const
Vector_double &, const
std::vector< stf::parInfo >
, double)> 
stf::Output
 Print the output of a fit into a stf::Table.
typedef boost::function< void(const
Vector_double &, double,
double, double, double, double,
Vector_double &)> 
stf::Init
 Initialising function for the parameters in stf::Func to start a fit.

Enumerations

enum  stf::cursor_type {
  stf::measure_cursor, stf::peak_cursor, stf::base_cursor, stf::decay_cursor,
  stf::latency_cursor, stf::zoom_cursor, stf::event_cursor, stf::undefined_cursor
}
 Mouse cursor types. More...
enum  stf::direction { stf::up, stf::down, stf::both, stf::undefined_direction }
 The direction of peak calculations. More...
enum  stf::zoom_channels { stf::zoomch1, stf::zoomch2, stf::zoomboth }
 Determines which channels to scale. More...
enum  stf::latency_mode {
  stf::manualMode = 0, stf::peakMode = 1, stf::riseMode = 2, stf::halfMode = 3,
  stf::footMode = 4, stf::undefinedMode
}
 Latency cursor settings. More...
enum  stf::latency_window_mode { stf::defaultMode = 0, stf::windowMode = 1 }
 Latency window settings. More...
enum  stf::extraction_mode { stf::criterion, stf::correlation, stf::deconvolution }
 Deconvolution. More...

Functions

std::string stf::wx2std (const wxString &wxs)
wxString stf::std2wx (const std::string &sst)
int stf::round (double toRound)
 Does what it says.
template<typename T >
stf::linFit (const std::vector< T > &x, const std::vector< T > &y, T &m, T &c)
 Performs a linear fit.
double StfDll stf::lmFit (const Vector_double &data, double dt, const stf::storedFunc &fitFunc, const Vector_double &opts, bool use_scaling, Vector_double &p, std::string &info, int &warning)
 Uses the Levenberg-Marquardt algorithm to perform a non-linear least-squares fit.
double stf::flin (double x, const Vector_double &p)
 Linear function.
void stf::flin_init (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Dummy function to be passed to stf::storedFunc for linear functions.
stf::storedFunc stf::initLinFunc ()
 initializes a linear function
Vector_double stf::get_scale (Vector_double &data, double oldx)
 Compute and perform normalisation.
Vector_double stf::LM_default_opts ()
 Return default LM options.
double stf::base (double &var, const std::vector< double > &data, std::size_t llb, std::size_t ulb)
 Calculate the average of all sampling points between and including llb and ulb.
double stf::peak (const std::vector< double > &data, double base, std::size_t llp, std::size_t ulp, int pM, stf::direction, double &maxT)
 Find the peak value of data between llp and ulp.
double stf::threshold (const std::vector< double > &data, std::size_t llp, std::size_t ulp, double slope, double &thrT, long windowLength)
 Find the value within data between llp and ulp at which slope is exceeded.
double stf::risetime (const std::vector< double > &data, double base, double ampl, double left, double right, double frac, std::size_t &tLoId, std::size_t &tHiId, double &tLoReal)
 Find 20 to 80% rise time of an event in data.
double stf::t_half (const std::vector< double > &data, double base, double ampl, double left, double right, double center, std::size_t &t50LeftId, std::size_t &t50RightId, double &t50LeftReal)
 Find the full width at half-maximal amplitude of an event within data.
double stf::maxRise (const std::vector< double > &data, double left, double right, double &maxRiseT, double &maxRiseY, long windowLength)
 Find the maximal slope during the rising phase of an event within data.
double stf::maxDecay (const std::vector< double > &data, double left, double right, double &maxDecayT, double &maxDecayY, long windowLength)
 Find the maximal slope during the decaying phase of an event within data.
double stf::fexp (double x, const Vector_double &p)
 Sum of n exponential functions.
Vector_double stf::fexp_jac (double x, const Vector_double &p)
 Computes the Jacobian of stf::fexp().
void stf::fexp_init (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::fexp() to data.
void stf::fexp_init2 (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Yet another initialiser for fitting stf::fexp() to data.
double stf::fexpde (double x, const Vector_double &p)
 Monoexponential function with delay.
void stf::fexpde_init (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::fexpde() to data.
double stf::fexpbde (double x, const Vector_double &p)
 Biexponential function with delay.
void stf::fexpbde_init (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::fexpde() to data.
double stf::falpha (double x, const Vector_double &p)
 Alpha function.
Vector_double stf::falpha_jac (double x, const Vector_double &p)
 Computes the Jacobian of stf::falpha().
double stf::fHH (double x, const Vector_double &p)
 Hodgkin-Huxley sodium conductance function.
double stf::fgauss (double x, const Vector_double &p)
 Computes the sum of an arbitrary number of Gaussians.
Vector_double stf::fgauss_jac (double x, const Vector_double &p)
 Computes the Jacobian of a sum of Gaussians.
double stf::fgnabiexp (double x, const Vector_double &p)
 power of 1 sodium conductance function.
void stf::falpha_init (const Vector_double &data, double base, double peak, double RTLoHI, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::falpha() to data.
void stf::fgauss_init (const Vector_double &data, double base, double peak, double RTLoHI, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::fgauss() to data.
void stf::fHH_init (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::falpha() to data.
void stf::fgnabiexp_init (const Vector_double &data, double base, double peak, double RTLoHi, double HalfWidth, double dt, Vector_double &pInit)
 Initialises parameters for fitting stf::falpha() to data.
double stf::xscale (double param, double xscale, double xoff, double yscale, double yoff)
 Scales a parameter that linearly depends on x.
double stf::xunscale (double param, double xscale, double xoff, double yscale, double yoff)
 Unscales a parameter that linearly depends on x.
double stf::yscale (double param, double xscale, double xoff, double yscale, double yoff)
 Scales a parameter that linearly depends on y.
double stf::yscaleoffset (double param, double xscale, double xoff, double yscale, double yoff)
 Scales a parameter that linearly depends on y and adds an offset.
double stf::yunscale (double param, double xscale, double xoff, double yscale, double yoff)
 Unscales a parameter that linearly depends on y.
double stf::yunscaleoffset (double param, double xscale, double xoff, double yscale, double yoff)
 Unscales a parameter that linearly depends on y and removes the offset.
std::vector< stf::parInfostf::getParInfoExp (int n_exp)
 Creates stf::parInfo structs for n-exponential functions.
stf::Table stf::outputWTau (const Vector_double &p, const std::vector< stf::parInfo > &parsInfo, double chisqr)
 Calculates a weighted time constant.
std::size_t stf::whereis (const Vector_double &data, double value)
 Finds the index of data where value is encountered for the first time.
std::vector< stf::storedFuncstf::GetFuncLib ()
 Returns the library of functions for non-linear regression.
Vector_double stf::nojac (double x, const Vector_double &p)
 Dummy function, serves as a placeholder to initialize functions without a Jacobian.
double stf::noscale (double param, double xscale, double xoff, double yscale, double yoff)
 Dummy function, serves as a placeholder to initialize parameters without a scaling function.
Table stf::defaultOutput (const Vector_double &pars, const std::vector< parInfo > &parsInfo, double chisqr)
 Default fit output function, constructing a stf::Table from the parameters, their description and chisqr.
template<typename T >
stf::SQR (T a)
 Calculates the square of a number.
Vector_double stf::filter (const Vector_double &toFilter, std::size_t filter_start, std::size_t filter_end, const Vector_double &a, int SR, stf::Func func, bool inverse=false)
 Convolves a data set with a filter function.
std::map< double, int > stf::histogram (const Vector_double &data, int nbins=-1)
 Computes a histogram.
StfDll Vector_double stf::deconvolve (const Vector_double &data, const Vector_double &templ, int SR, double hipass, double lopass, stfio::ProgressInfo &progDlg)
 Deconvolves a template from a signal.
template<class T >
std::vector< T > stf::cubicSpline (const std::vector< T > &y, T oldF, T newF)
 Interpolates a dataset using cubic splines.
wxString stf::sectionToString (const Section &section)
 Converts a Section to a wxString.
template<class T >
std::vector< T > stf::diff (const std::vector< T > &input, T x_scale)
 Differentiate data.
double stf::integrate_simpson (const Vector_double &input, std::size_t a, std::size_t b, double x_scale)
 Integration using Simpson's rule.
double stf::integrate_trapezium (const Vector_double &input, std::size_t a, std::size_t b, double x_scale)
 Integration using the trapezium rule.
int stf::linsolv (int m, int n, int nrhs, Vector_double &A, Vector_double &B)
 Solves a linear equation system using LAPACK.
Vector_double stf::quad (const Vector_double &data, std::size_t begin, std::size_t end)
 Solve quadratic equations for 3 adjacent sampling points.
StfDll Vector_double stf::detectionCriterion (const Vector_double &data, const Vector_double &templ, stfio::ProgressInfo &progDlg)
 Computes the event detection criterion according to Clements & Bekkers (1997).
StfDll std::vector< int > stf::peakIndices (const Vector_double &data, double threshold, int minDistance)
 Searches for positive-going peaks.
StfDll Vector_double stf::linCorr (const Vector_double &va1, const Vector_double &va2, stfio::ProgressInfo &progDlg)
 Computes the linear correlation between two arrays.
double stf::fgaussColqu (double x, const Vector_double &p)
 Computes a Gaussian that can be used as a filter kernel.
double stf::fboltz (double x, const Vector_double &p)
 Computes a Boltzmann function.
double stf::fbessel (double x, int n)
 Computes a Bessel polynomial.
double stf::fbessel4 (double x, const Vector_double &p)
 Computes a 4th-order Bessel polynomial that can be used as a filter kernel.
wxString stf::CreatePreview (const wxString &fName)
 Creates a preview of a text file.
int stf::fac (int arg)
 Computes the faculty of an integer.
int stf::pow2 (int arg)
 Computes \( 2^{arg} \). Uses the bitwise-shift operator (<<).
wxString stf::noPath (const wxString &fName)
 Strips the directory off a full path name, returns only the filename.

Variables

const double stf::PI = 3.14159265358979323846
 Add decimals if you are not satisfied.

Typedef Documentation

typedef boost::function<double(double, const Vector_double&)> stf::Func

A function taking a double and a vector and returning a double.

Type definition for a function (or, to be precise, any 'callable entity') that takes a double (the x-value) and a vector of parameters and returns the function's result (the y-value).


Enumeration Type Documentation

Mouse cursor types.

Enumerator:
measure_cursor 

Measurement cursor (crosshair).

peak_cursor 

Peak calculation limits cursor.

base_cursor 

Baseline calculation limits cursor.

decay_cursor 

Fit limits cursor.

latency_cursor 

Latency cursor.

zoom_cursor 

Zoom rectangle cursor.

event_cursor 

Event mode cursor.

undefined_cursor 

Undefined cursor.

The direction of peak calculations.

Enumerator:
up 

Find positive-going peaks.

down 

Find negative-going peaks.

both 

Find negative- or positive-going peaks, whichever is larger.

undefined_direction 

Undefined direction.

Deconvolution.

Enumerator:
criterion 

Clements & Bekkers criterion.

correlation 

Jonas et al. correlation coefficient.

deconvolution 

Pernia-Andrade et al. deconvolution.

Latency cursor settings.

Enumerator:
manualMode 

Set the corresponding latency cursor manually (by clicking on the graph).

peakMode 

Set the corresponding latency cursor to the peak.

riseMode 

Set the corresponding latency cursor to the maximal slope of rise.

halfMode 

Set the corresponding latency cursor to the half-maximal amplitude.

footMode 

Set the corresponding latency cursor to the beginning of an event.

undefinedMode 

undefined mode.

Latency window settings.

Enumerator:
defaultMode 

Use the current peak cursor window for the active channel.

windowMode 

Use a window of 100 sampling points around the peak.

Determines which channels to scale.

Enumerator:
zoomch1 

Scaling applies to channel 1 only.

zoomch2 

Scaling applies to channel 2 only.

zoomboth 

Scaling applies to both channels.


Function Documentation

double stf::base ( double &  var,
const std::vector< double > &  data,
std::size_t  llb,
std::size_t  ulb 
)

Calculate the average of all sampling points between and including llb and ulb.

Parameters:
varWill contain the variance on exit.
dataThe data waveform to be analysed.
llbAveraging will be started at this index.
ulbIndex of the last data point included in the average (legacy of the PASCAL version).
llpLower limit of the peak window (see stf::peak()).
ulpUpper limit of the peak window (see stf::peak()).
Returns:
The baseline value.

Referenced by wxStfDoc::GetBase().

wxString stf::CreatePreview ( const wxString fName)

Creates a preview of a text file.

Parameters:
fNameFull path name of the file.
Returns:
A string showing at most the initial 100 lines of the text file.
template<class T >
std::vector< T > stf::cubicSpline ( const std::vector< T > &  y,
oldF,
newF 
)

Interpolates a dataset using cubic splines.

Parameters:
yThe valarray to be interpolated.
oldFThe original sampling frequency.
newFThe new frequency of the interpolated array.
Returns:
The interpolated data set.
StfDll Vector_double stf::deconvolve ( const Vector_double &  data,
const Vector_double &  templ,
int  SR,
double  hipass,
double  lopass,
stfio::ProgressInfo progDlg 
)

Deconvolves a template from a signal.

Parameters:
dataThe input signal
templThe template
SRThe sampling rate in kHz.
hipassHighpass filter cutoff frequency in kHz
lopassLowpass filter cutoff frequency in kHz
Returns:
The result of the deconvolution
StfDll Vector_double stf::detectionCriterion ( const Vector_double &  data,
const Vector_double &  templ,
stfio::ProgressInfo progDlg 
)

Computes the event detection criterion according to Clements & Bekkers (1997).

Parameters:
dataThe valarray from which to extract events.
templA template waveform that is used for event detection.
Returns:
The detection criterion for every value of data.
int stf::fac ( int  arg)

Computes the faculty of an integer.

Parameters:
argArgument of the function.
Returns:
The faculty of arg.
double stf::falpha ( double  x,
const Vector_double &  p 
)

Alpha function.

\[f(x)=p_0 p_1^2 x \mathrm{e}^{-p_1 x} + p_2\]

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the amplitude,
p[1] is the rate and
p[2] is the offset.
Returns:
The evaluated function.
void stf::falpha_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHI,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::falpha() to data.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 3. On exit, will contain initial parameter estimates.
Vector_double stf::falpha_jac ( double  x,
const Vector_double &  p 
)

Computes the Jacobian of stf::falpha().

\begin{eqnarray*} j_0(x) &=& \frac{df(x)}{dp_0} = p_1^2 x \mathrm{e}^{-p_1 x} \\ j_1(x) &=& \frac{df(x)}{dp_1} = p_0 p_1 x \left( 2 \mathrm{e}^{-p_1 x} - p_1 x \mathrm{e}^{-p_1 x} \right) \\ j_2(x) &=& \frac{df(x)}{dp_2} = 1.0 \end{eqnarray*}

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the amplitude,
p[1] is the rate and
p[2] is the offset.
Returns:
A valarray j with the evaluated Jacobian, where
j[0] contains the derivative with respect to p[0],
j[1] contains the derivative with respect to p[1] and
j[2] contains the derivative with respect to p[2].
double stf::fbessel ( double  x,
int  n 
)

Computes a Bessel polynomial.

\[ f(x, n) = \sum_{k=0}^n \frac{ \left( 2n - k \right) ! }{ \left( n - k \right) ! k! } \frac{x^k}{ 2^{n-k} } \]

Parameters:
xArgument of the function.
nOrder of the polynomial.
Returns:
The evaluated function.
double stf::fbessel4 ( double  x,
const Vector_double &  p 
)

Computes a 4th-order Bessel polynomial that can be used as a filter kernel.

\[ f(x) = \frac{b(0,4)}{b(\frac{0.355589x}{p_0},4)} \]

where \( b(a,n) \) is the bessel polynomial stf::fbessel().

Parameters:
xArgument of the function.
pFunction parameters, where
p[0] is the corner frequency (-3 dB attenuation)
Returns:
The evaluated function.
double stf::fboltz ( double  x,
const Vector_double &  p 
)

Computes a Boltzmann function.

\[f(x)=\frac{1}{1+\mathrm{e}^{\frac{p_0-x}{p_1}}}\]

Parameters:
xArgument of the function.
pFunction parameters, where
p[0] is the midpoint and
p[1] is the slope of the function.
Returns:
The evaluated function.
double stf::fexp ( double  x,
const Vector_double &  p 
)

Sum of n exponential functions.

\[f(x)=p_{2n} + \sum_{i=0}^{2 n - 1}p_{2i}\mathrm{e}^{\frac{x}{p_{2i + 1}}}\]

Parameters:
xFunction argument.
pA valarray of parameters of size 2n+1, where
n is the number of exponential terms,
p[2i] is the amplitude term,
p[2i+1] is the time constant,
p[2n], the last element, contains the offset and
i denotes the i -th exponential term (running from 0 to n-1).
Returns:
The evaluated function.
void stf::fexp_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHi,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::fexp() to data.

This needs to be made more robust.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 2n+1, where n is the number of exponential functions. On exit, will contain initial parameter estimates.
void stf::fexp_init2 ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHi,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Yet another initialiser for fitting stf::fexp() to data.

In this case, one of the amplitude terms will have another sign than the others, making it more suitable for fitting PSCs or PSPs. However, this often fails to work in practice.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 2n+1, where n is the number of exponential functions. On exit, will contain initial parameter estimates.
Vector_double stf::fexp_jac ( double  x,
const Vector_double &  p 
)

Computes the Jacobian of stf::fexp().

\begin{eqnarray*} j_{2i}(x) &=& \frac{df(x)}{dp_{2i}} = \mathrm{e}^{\frac{-x}{p_{2i+1}}} \\ j_{2i+1}(x) &=& \frac{df(x)}{dp_{2i+1}} = \frac{p_{2i}}{p_{2i+1}^2} x \mathrm{e}^{\frac{-x}{p_{2i+1}}} \\ j_n(x) &=& \frac{df(x)}{dp_{n}} = 1 \end{eqnarray*}

Parameters:
xFunction argument.
pA valarray of parameters of size 2n+1, where
n is the number of exponential terms,
p[2i] is the amplitude term,
p[2i+1] is the time constant,
p[2n], the last element, contains the offset and
i denotes the i -th exponential term (running from 0 to n-1).
Returns:
A valarray j of size 2n+1 with the evaluated Jacobian, where
j[2i] contains the derivative with respect to p[2i],
j[2i+1] contains the derivative with respect to p[2i+1] and
j[2n], the last element, contains the derivative with respect to p[2n].
double stf::fexpbde ( double  x,
const Vector_double &  p 
)

Biexponential function with delay.

\begin{eqnarray*} f(x)= \begin{cases} p_0 & \mbox{if }x < p_1 \\ p_3 \left( \mathrm{e}^{\frac{p_1 - x}{p_2}} - \mathrm{e}^{\frac{p_1 - x}{p_4}} \right) + p_0 & \mbox{if }x \geq p_1 \end{cases} \end{eqnarray*}

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the baseline,
p[1] is the delay,
p[2] is the later (slower) time constant,
p[3] is the amplitude and
p[4] is the earlier (faster) time constant,
Returns:
The evaluated function.
void stf::fexpbde_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHi,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::fexpde() to data.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 4. On exit, will contain initial parameter estimates.
double stf::fexpde ( double  x,
const Vector_double &  p 
)

Monoexponential function with delay.

\begin{eqnarray*} f(x)= \begin{cases} p_0, & \mbox{if }x < p_3 \\ \left( p_0 - p_2 \right) \mathrm{e}^{\frac{p_3 - x}{p_1}} + p_2, & \mbox{if }x \geq p_3 \end{cases} \end{eqnarray*}

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the baseline,
p[1] is the time constant,
p[2] is the amplitude and
p[3] is the delay.
Returns:
The evaluated function.
void stf::fexpde_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHi,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::fexpde() to data.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 4. On exit, will contain initial parameter estimates.
double stf::fgauss ( double  x,
const Vector_double &  p 
)

Computes the sum of an arbitrary number of Gaussians.

\[ f(x) = \sum_{i=0}^{n-1}p_{3i}\mathrm{e}^{- \left( \frac{x-p_{3i+1}}{p_{3i+2}} \right) ^2} \]

Parameters:
xArgument of the function.
pA valarray of function parameters of size 3n, where
p[3i] is the amplitude of the Gaussian
p[3i+1] is the position of the center of the peak,
p[3i+2] is the width of the Gaussian,
n is the number of Gaussian functions and
i is the 0-based index of the i-th Gaussian.
Returns:
The evaluated function.
void stf::fgauss_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHI,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::fgauss() to data.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 3. On exit, will contain initial parameter estimates.
double stf::fgaussColqu ( double  x,
const Vector_double &  p 
)

Computes a Gaussian that can be used as a filter kernel.

\[ f(x) = \mathrm{e}^{-0.3466 \left( \frac{x}{p_{0}} \right) ^2} \]

Parameters:
xArgument of the function.
pFunction parameters, where
p[0] is the corner frequency (-3 dB according to Colquhoun)
Returns:
The evaluated function.
double stf::fgnabiexp ( double  x,
const Vector_double &  p 
)

power of 1 sodium conductance function.

\[f(x)=p_0\left(1-\mathrm{e}^{\frac{-x}{p_1}}\right)\mathrm{e}^{\frac{-x}{p_2}} + p_3\]

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the amplitude \(g'_{Na}\),
p[1] is \(\tau_m\),
p[2] is \(\tau_h\) and
p[3] is the offset.
Returns:
The evaluated function.
void stf::fgnabiexp_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHi,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::falpha() to data.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 4. On exit, will contain initial parameter estimates.
double stf::fHH ( double  x,
const Vector_double &  p 
)

Hodgkin-Huxley sodium conductance function.

\[f(x)=p_0\left(1-\mathrm{e}^{\frac{-x}{p_1}}\right)^3\mathrm{e}^{\frac{-x}{p_2}} + p_3\]

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the amplitude \(g'_{Na}\),
p[1] is \(\tau_m\),
p[2] is \(\tau_h\) and
p[3] is the offset.
Returns:
The evaluated function.
void stf::fHH_init ( const Vector_double &  data,
double  base,
double  peak,
double  RTLoHi,
double  HalfWidth,
double  dt,
Vector_double &  pInit 
)

Initialises parameters for fitting stf::falpha() to data.

Parameters:
dataThe waveform of the data for the fit.
baseBaseline of data.
peakPeak value of data.
dtThe sampling interval.
pInitOn entry, pass a valarray of size 4. On exit, will contain initial parameter estimates.
Vector_double stf::filter ( const Vector_double &  toFilter,
std::size_t  filter_start,
std::size_t  filter_end,
const Vector_double &  a,
int  SR,
stf::Func  func,
bool  inverse = false 
)

Convolves a data set with a filter function.

Parameters:
toFilterThe valarray to be filtered.
filter_startThe index from which to start filtering.
filter_endThe index at which to stop filtering.
aA valarray of parameters for the filter function.
SRThe sampling rate.
funcThe filter function in the frequency domain.
inversetrue if (1- func) should be used as the filter function, false otherwise
Returns:
The convolved data set.
double stf::flin ( double  x,
const Vector_double &  p 
)

Linear function.

\[f(x)=p_0 x + p_1\]

Parameters:
xFunction argument.
pA valarray of parameters, where
p[0] is the slope and
p[1] is the y intersection.
Returns:
The evaluated function.
Vector_double stf::get_scale ( Vector_double &  data,
double  oldx 
)

Compute and perform normalisation.

Parameters:
dataData vector; will be scaled upon return
oldxoriginal x interval
Returns:
A vector with
[0] x scale [1] x offset [2] y scale [3] y offset
std::vector<stf::storedFunc> stf::GetFuncLib ( )

Returns the library of functions for non-linear regression.

Returns:
A vector of non-linear regression functions.
std::vector<stf::parInfo> stf::getParInfoExp ( int  n_exp)

Creates stf::parInfo structs for n-exponential functions.

Parameters:
n_expNumber of exponential terms.
Returns:
A vector of parameter information structs.
std::map<double, int> stf::histogram ( const Vector_double &  data,
int  nbins = -1 
)

Computes a histogram.

Parameters:
dataThe signal
nbinsNumber of bins in the histogram.
Returns:
A map with lower bin limits as keys, number of observations as values.
stf::storedFunc stf::initLinFunc ( )

initializes a linear function

Returns:
An stf::storedFunc that can be used to store a linear function after a fit
double stf::integrate_simpson ( const Vector_double &  input,
std::size_t  a,
std::size_t  b,
double  x_scale 
)

Integration using Simpson's rule.

Parameters:
inputThe valarray to be integrated.
aStart of the integration interval.
bEnd of the integration interval.
x_scaleSampling interval.
Returns:
The integral of input between a and b.
double stf::integrate_trapezium ( const Vector_double &  input,
std::size_t  a,
std::size_t  b,
double  x_scale 
)

Integration using the trapezium rule.

Parameters:
inputThe valarray to be integrated.
aStart of the integration interval.
bEnd of the integration interval.
x_scaleSampling interval.
Returns:
The integral of input between a and b.
StfDll Vector_double stf::linCorr ( const Vector_double &  va1,
const Vector_double &  va2,
stfio::ProgressInfo progDlg 
)

Computes the linear correlation between two arrays.

Parameters:
va1First array.
va2Second array.
Returns:
The linear correlation between the two arrays for each data point of va1.
template<typename T >
T stf::linFit ( const std::vector< T > &  x,
const std::vector< T > &  y,
T &  m,
T &  c 
)

Performs a linear fit.

Parameters:
xThe x- values of the data that are to be fitted.
yThe y- values of the data that are to be fitted.
mOn exit, the slope of the regression line.
cOn exit, the y-intercept of the regression line.
Returns:
A valarray containing the waveform of the fitted function.
int stf::linsolv ( int  m,
int  n,
int  nrhs,
Vector_double &  A,
Vector_double &  B 
)

Solves a linear equation system using LAPACK.

Uses column-major order for matrices. For an example, see Section::SetIsIntegrated()

Parameters:
mNumber of rows of the matrix A.
nNumber of columns of the matrix A.
nrhsNumber of columns of the matrix B.
AOn entry, the left-hand-side matrix. On exit, the factors L and U from the factorization A = P*L*U; the unit diagonal elements of L are not stored.
BOn entry, the right-hand-side matrix. On exit, the solution to the linear equation system.
Returns:
At present, always returns 0.
Vector_double stf::LM_default_opts ( )

Return default LM options.

Returns:
Default LM options
double StfDll stf::lmFit ( const Vector_double &  data,
double  dt,
const stf::storedFunc fitFunc,
const Vector_double &  opts,
bool  use_scaling,
Vector_double &  p,
std::string &  info,
int &  warning 
)

Uses the Levenberg-Marquardt algorithm to perform a non-linear least-squares fit.

Parameters:
dataA valarray containing the data.
dtThe sampling interval of data.
fitFuncAn stf::storedFunc to be fitted to data.
optsOptions controlling Lourakis' implementation of the algorithm.
use_scalingWhether to scale x and y-amplitudes to 1.0
pfunc's parameters. Should be set to an initial guess on entry. Will contain the best-fit values on exit.
infoInformation about why the fit stopped iterating
warningA warning code on return.
Returns:
The sum of squred errors between data and the best-fit function.
double stf::maxDecay ( const std::vector< double > &  data,
double  left,
double  right,
double &  maxDecayT,
double &  maxDecayY,
long  windowLength 
)

Find the maximal slope during the decaying phase of an event within data.

Parameters:
dataThe data waveform to be analysed.
leftDelimits the search to the left.
rightDelimits the search to the right.
maxDecayTThe interpolated time point of the maximal slope of decay in units of sampling points.
maxDecayYThe interpolated value of data at maxDecayT.
windowLengthis the distance (in number of samples) used to compute the slope, the default value is 1.
Returns:
The maximal slope during the decaying phase.

Referenced by wxStfDoc::GetMaxDecay().

double stf::maxRise ( const std::vector< double > &  data,
double  left,
double  right,
double &  maxRiseT,
double &  maxRiseY,
long  windowLength 
)

Find the maximal slope during the rising phase of an event within data.

Parameters:
dataThe data waveform to be analysed.
leftDelimits the search to the left.
rightDelimits the search to the right.
maxRiseTThe interpolated time point of the maximal slope of rise in units of sampling points.
maxRiseYThe interpolated value of data at maxRiseT.
windowLengthis the distance (in number of samples) used to compute the slope, the default value is 1.
Returns:
The maximal slope during the rising phase.

Referenced by wxStfDoc::GetMaxRise().

wxString stf::noPath ( const wxString fName)

Strips the directory off a full path name, returns only the filename.

Parameters:
fNameThe full path of a file.
Returns:
The file name without the directory.
stf::Table stf::outputWTau ( const Vector_double &  p,
const std::vector< stf::parInfo > &  parsInfo,
double  chisqr 
)

Calculates a weighted time constant.

Parameters:
pParameters of an exponential function (see stf::fexp()).
parsInfoInformation about the parameters p.
chisqrThe sum of squared errors, as returned from a least-squares fit.
Returns:
A formatted table of results.
double stf::peak ( const std::vector< double > &  data,
double  base,
std::size_t  llp,
std::size_t  ulp,
int  pM,
stf::direction  ,
double &  maxT 
)

Find the peak value of data between llp and ulp.

Note that peaks will be detected by measuring from base, but the return value is given from 0. Data points at both llp and ulp will be included in the search (legacy of Stimfit for PASCAL).

Parameters:
dataThe data waveform to be analysed.
baseThe baseline value.
llpLower limit of the peak window.
ulpUpper limit of the peak window.
pMIf pM > 1, a sliding (boxcar) average of width pM will be used to measure the peak.
dirCan be
stf::up for positive-going peaks,
stf::down for negative-going peaks or
stf::both for negative- or positive-going peaks, whichever is larger.
maxTOn exit, the index of the peak value. May be interpolated if pM > 1.
Returns:
The peak value, measured from 0.

Referenced by wxStfDoc::GetPeak().

StfDll std::vector<int> stf::peakIndices ( const Vector_double &  data,
double  threshold,
int  minDistance 
)

Searches for positive-going peaks.

Parameters:
dataThe valarray to be searched for peaks.
thresholdMinimal amplitude of a peak.
minDistanceMinimal distance between subsequent peaks.
Returns:
A vector of indices where peaks have occurred in data.
int stf::pow2 ( int  arg) [inline]

Computes \( 2^{arg} \). Uses the bitwise-shift operator (<<).

Parameters:
argArgument of the function.
Returns:
\( 2^{arg} \).
Vector_double stf::quad ( const Vector_double &  data,
std::size_t  begin,
std::size_t  end 
)

Solve quadratic equations for 3 adjacent sampling points.

Parameters:
dataThe data vector
beginStart of interval to be used
endEnd of interval to be used
Returns:
Parameters of quadratic equation
double stf::risetime ( const std::vector< double > &  data,
double  base,
double  ampl,
double  left,
double  right,
double  frac,
std::size_t &  tLoId,
std::size_t &  tHiId,
double &  tLoReal 
)

Find 20 to 80% rise time of an event in data.

Although t80real is not explicitly returned, it can be calculated from t20Real+risetime.

Parameters:
dataThe data waveform to be analysed.
baseThe baseline value.
amplThe amplitude of the event (typically, peak-base).
leftDelimits the search to the left.
rightDelimits the search to the right.
t20IdOn exit, the index wich is closest to the 20%-point.
t80IdOn exit, the index wich is closest to the 80%-point.
t20Realthe linearly interpolated 20%-timepoint in units of sampling points.
Returns:
The rise time.
int stf::round ( double  toRound) [inline]

Does what it says.

Parameters:
toRoundThe double to be rounded.
Returns:
The rounded integer.
wxString stf::sectionToString ( const Section section)

Converts a Section to a wxString.

Parameters:
sectionThe Section to be written to a string.
Returns:
A string containing the x- and y-values of the section in two columns.
template<typename T >
T stf::SQR ( a) [inline]

Calculates the square of a number.

Parameters:
aArgument of the function.
Returns:
a ^2
double stf::t_half ( const std::vector< double > &  data,
double  base,
double  ampl,
double  left,
double  right,
double  center,
std::size_t &  t50LeftId,
std::size_t &  t50RightId,
double &  t50LeftReal 
)

Find the full width at half-maximal amplitude of an event within data.

Although t50RightReal is not explicitly returned, it can be calculated from t50LeftReal+t_half.

Parameters:
dataThe data waveform to be analysed.
baseThe baseline value.
amplThe amplitude of the event (typically, peak-base).
leftDelimits the search to the left.
rightDelimits the search to the right.
centerThe estimated center of an event from which to start searching to the left and to the right (typically, the index of the peak).
t50LeftIdOn exit, the index wich is closest to the left 50%-point.
t50RightIdOn exit, the index wich is closest to the right 50%-point.
t50LeftRealthe linearly interpolated left 50%-timepoint in units of sampling points.
Returns:
The full width at half-maximal amplitude.
double stf::threshold ( const std::vector< double > &  data,
std::size_t  llp,
std::size_t  ulp,
double  slope,
double &  thrT,
long  windowLength 
)

Find the value within data between llp and ulp at which slope is exceeded.

Parameters:
dataThe data waveform to be analysed.
llpLower limit of the peak window.
ulpUpper limit of the peak window.
thrTOn exit, The interpolated time point of the threshold crossing in units of sampling points, or a negative value if the threshold wasn't found.
windowLengthis the distance (in number of samples) used to compute the difference, the default value is 1.
Returns:
The interpolated threshold value.

Referenced by wxStfDoc::GetThreshold().

std::size_t stf::whereis ( const Vector_double &  data,
double  value 
)

Finds the index of data where value is encountered for the first time.

Parameters:
dataThe waveform to be searched.
valueThe value to be found.
Returns:
The index of data right after value has been crossed.
double stf::xscale ( double  param,
double  xscale,
double  xoff,
double  yscale,
double  yoff 
)

Scales a parameter that linearly depends on x.

Parameters:
Theparameter to scale
xscalex scaling factor
xoffx offset
yscaley scaling factor
yoffy offset
Returns:
Scaled parameter
double stf::xunscale ( double  param,
double  xscale,
double  xoff,
double  yscale,
double  yoff 
)

Unscales a parameter that linearly depends on x.

Parameters:
Theparameter to scale
xscalex scaling factor
xoffx offset
yscaley scaling factor
yoffy offset
Returns:
Unscaled parameter
double stf::yscale ( double  param,
double  xscale,
double  xoff,
double  yscale,
double  yoff 
)

Scales a parameter that linearly depends on y.

Parameters:
Theparameter to scale
xscalex scaling factor
xoffx offset
yscaley scaling factor
yoffy offset
double stf::yscaleoffset ( double  param,
double  xscale,
double  xoff,
double  yscale,
double  yoff 
)

Scales a parameter that linearly depends on y and adds an offset.

Parameters:
Theparameter to scale
xscalex scaling factor
xoffx offset
yscaley scaling factor
yoffy offset
double stf::yunscale ( double  param,
double  xscale,
double  xoff,
double  yscale,
double  yoff 
)

Unscales a parameter that linearly depends on y.

Parameters:
Theparameter to scale
xscalex scaling factor
xoffx offset
yscaley scaling factor
yoffy offset
Returns:
Unscaled parameter
double stf::yunscaleoffset ( double  param,
double  xscale,
double  xoff,
double  yscale,
double  yoff 
)

Unscales a parameter that linearly depends on y and removes the offset.

Parameters:
Theparameter to scale
xscalex scaling factor
xoffx offset
yscaley scaling factor
yoffy offset
Returns:
Unscaled parameter