Classes | Enumerations | Functions
FFTtools Namespace Reference

Classes

class  AnalyticSignal
 
class  Averager
 
class  BandlimitedSampledSignal
 
class  BlackmanHarrisWindow
 
class  BlackmanWindow
 
class  BoxFilter
 
class  ButterworthFilter
 
class  ChebyshevIFilter
 
class  CompositeSignal
 
class  CWT
 
class  DifferenceFilter
 
class  DigitalFilter
 
class  DigitalFilterSeries
 
class  FFTWindow
 
class  FFTWindowType
 
class  FIRFilter
 
class  GaussianFilter
 
class  GaussianWindow
 
class  HammingWindow
 
class  HannWindow
 
class  IIRFilter
 
class  InverseFilter
 
class  KaiserWindow
 
class  LanczosFilter
 
class  MinimumPhaseFilter
 
class  RCFilter
 
class  RectangularWindow
 
class  SavitzkyGolayFilter
 
class  SincFilter
 
class  SineFitter
 
class  SineFitterLimits
 
class  SineSubtract
 
class  SineSubtractResult
 
class  ThermalNoise
 
class  ThiranFilter
 
class  TransformedZPKFilter
 
class  TriangularWindow
 
class  TukeyWindow
 

Enumerations

enum  FilterTopology { LOWPASS, HIGHPASS, BANDPASS, NOTCH }
 
enum  DirectConvolveEdgeBehavior { ZEROES_OUTSIDE, REPEAT_OUTSIDE }
 

Functions

TGraph * normalizeWaveform (TGraph *inGraph)
 
double * getCorrelationFromFFT (int length, const FFTWComplex *theFFT1, const FFTWComplex *theFFT2)
 
TGraph * getInterpolatedGraph (const TGraph *grIn, Double_t deltaT)
 Interpolation Routines that use ROOT::Math::Interpolator. More...
 
TGraph * getInterpolatedGraphFreqDom (const TGraph *grIn, Double_t deltaT)
 Interpolation Routines that zeropads the FFT. More...
 
TGraph * getConvolution (const TGraph *grA, const TGraph *grB)
 Convolution. More...
 
RFSignalgetConvolution (RFSignal *grA, RFSignal *grB)
 Convolution. More...
 
double getAbs (FFTWComplex &theNum)
 Returns the magnitude of a complex number. More...
 
double * doInvFFT (int length, const FFTWComplex *theInput)
 Computes an inverse FFT. More...
 
FFTWComplexdoFFT (int length, const double *theInput)
 Computes an FFT of an array of real numbers. More...
 
void doFFT (int length, const double *properly_aligned_input, FFTWComplex *properly_aligned_output)
 
void doInvFFTClobber (int length, FFTWComplex *properly_aligned_input_that_will_likely_be_clobbered, double *properly_aligned_output)
 
void doInvFFTNoClobber (int length, const FFTWComplex *properly_aligned_input, double *properly_aligned_output)
 
TGraph * convertMagnitudeToTimeDomain (const TGraph *inputMag)
 Converts inputMag (linear units) to time domain by means of hilbert transform assuming R_signal = 1/sqrt(2) (R_mag - i R^_mag);. More...
 
double * getCorrelation (const TGraph *gr1, const TGraph *gr2, int firstIndex, int lastIndex)
 Computes the correlation of two subsets of TGraphs. More...
 
double * getCorrelation (int length, const float *oldY1, const float *oldY2)
 Computes the correlation of two arrays. More...
 
double * getCorrelation (int length, const double *oldY1, const double *oldY2)
 Computes the correlation of two arrays. More...
 
TGraph * correlateAndAverage (Int_t numGraphs, TGraph **grPtrPtr)
 This is designed for when you want to average a number of graphs of the same thing together. It uses a correlation to find the deltaT between graphs and then shifts the graphs and coherently sums them. The return is the average of the input graphs. More...
 
TGraph * correlateAndAveragePredeterminedZone (Int_t numGraphs, TGraph **grPtrPtr, Int_t *correlationBinPtr, Int_t binWiggleRoom)
 
TGraph * interpolateCorrelateAndAverage (Double_t deltaTInt, Int_t numGraphs, TGraph **grPtrPtr)
 This is designed for when you want to average a number of graphs of the same thing together. It uses a correlation to find the deltaT between graphs and then shifts the graphs and coherently sums them. The return is the average of the input graphs. More...
 
Double_t * combineValuesUsingFFTs (Int_t numArrays, Double_t **thePtrPtr, Int_t eachLength)
 Returns the time domain result of a frequency domain sum of a number of arrays. As of writing this documentation, I'm not sure why this would be interesting. More...
 
TGraph * makePowerSpectrum (const TGraph *grWave)
 Returns the power spectral density. Note the PSD is unormalised (or if you prefer is normalised to the sum squared amplitude of the time domain). See this short note for my terminology. More...
 
TGraph * makePowerSpectrumPeriodogram (const TGraph *grWave)
 Returns the power spectral density. Note the PSD returned is the periodogram (or if you prefer is normalised to the mean squared amplitude of the time domain). See this short note for my terminology. More...
 
TGraph * makePowerSpectrumVoltsSeconds (const TGraph *grWave)
 Returns the power spectral density. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one. More...
 
TGraph * makePowerSpectrumVoltsSecondsBartlett (const TGraph *grWave)
 Returns the power spectral density of the input waveform convolved with a Bartlett Window. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one. More...
 
TGraph * makePSVSBartlettPaddedOverlap (const TGraph *grWave, Int_t padFactor=4, Int_t numFreqs=64)
 Returns the power spectral density of the input waveform. In this one we first zero pad the waveform and then split it up into overlapping segments and convolve each segment with the Bartlett window before summing the resulting PSD's. As the name suggests this function expects the input waveform to be a volts-seconds one. No idea if this one actually works, or where I read oabout this crazy method which is supposed to reduce the variance of the PSD estimator. More...
 
TGraph * makePowerSpectrumVoltsSecondsdB (const TGraph *grWave)
 Returns the power spectral density in dB units. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one. More...
 
TGraph * makePowerSpectrumMilliVoltsNanoSeconds (const TGraph *grWave)
 Returns the power spectral density in dB units. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a millivolts-nanoseconds one. More...
 
TGraph * makePowerSpectrumMilliVoltsNanoSecondsdB (const TGraph *grWave)
 
TGraph * makePowerSpectrumVoltsSecondsPadded (const TGraph *grWave, Int_t padFactor=4)
 Returns the power spectral density of the input waveform zero-padded by some factor. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one. More...
 
TGraph * makePowerSpectrumVoltsSecondsPaddeddB (const TGraph *grWave, Int_t padFactor=4)
 Returns the power spectral density in dB units of the input waveform zero-padded by some factor. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one. More...
 
TGraph * makeRawPowerSpectrum (const TGraph *grWave)
 Returns the power spectral density in completely unormalised unit (as in Parseval's theorem is not obeyed and there is an extra factor of N not removed form the PSD). See this short note for my terminology. More...
 
TGraph * getCorrelationGraph (const TGraph *gr1, const TGraph *gr2, Int_t *zeroOffset=0)
 Returns the correlation of two TGraphs. More...
 
TGraph * getNormalisedCorrelationGraph (const TGraph *gr1, const TGraph *gr2, Int_t *zeroOffset=0)
 Returns the normalised correlation of two TGraphs. More...
 
TGraph * getNormalisedCorrelationGraphTimeDomain (const TGraph *gr1, const TGraph *gr2, Int_t *zeroOffset=0, Int_t useDtRange=0, Double_t dtMin=-1000, Double_t dtMax=1000)
 Returns the normalised correlation of two TGraphs. More...
 
TGraph * getInterpolatedCorrelationGraph (const TGraph *grIn1, const TGraph *grIn2, Double_t deltaT)
 Returns the correlation of two interpolated TGraphs. More...
 
TGraph * makeInverseInverseSpectrum (const TGraph *grWave)
 Returns the inverse FFT of the FFT of the input TGraph. Seems pointless. More...
 
TGraph * combineGraphsUsingFFTs (Int_t numGraphs, TGraph **grPtr, const double *theWeights=0)
 Returns the time domain result of a frequency domain sum of a number of TGraphs. In the sum each TGraph is weighted by a value. As of writing this documentation, I'm not sure why this would be interesting. More...
 
TGraph * getBoxCar (const TGraph *grWave, Int_t halfWidth)
 Smooth graph using box car smoothing. More...
 
TGraph * getHilbertTransform (const TGraph *grWave)
 The Hilbert transform of the input TGraph. More...
 
double * getHilbertTransform (int N, const double *y)
 
TGraph * getHilbertEnvelope (const TGraph *grWave)
 The Hilbert envelope of the input TGraph. This is defined as e_i=sqrt(v_i^2 + h_i^2), where e_i, v_i and h_i are the i-th sample of the envelope, input graph and hilbert transform of the input graph repsectively. More...
 
Double_t sumPower (const TGraph *gr, Int_t firstBin=-1, Int_t lastBin=-1)
 The linear sum of the power in a TGraph (normally a PSD) More...
 
Double_t integratePower (const TGraph *gr, Int_t firstBin=-1, Int_t lastBin=-1)
 The integral of the power in a TGraph (normally a PSD) (i.e the sum of bin content*bin width) More...
 
Double_t sumVoltageSquared (const TGraph *gr, Int_t firstBin, Int_t lastBin)
 The sum of the voltage squared in a waveform. More...
 
Double_t integrateVoltageSquared (const TGraph *gr, Int_t firstBin=-1, Int_t lastBin=-1)
 The integral of the v^2*dt in a waveform. Now works for unevenly sampled waeforms. More...
 
Int_t getPeakBin (const TGraph *gr)
 Find the peak (maximum positive) bin in a TGraph. More...
 
Int_t getPeakBin (const TGraph *gr, Int_t firstBin, Int_t lastBin)
 Find the peak (maximum positive) bin in a TGraph. More...
 
Double_t getPeakXvalue (const TGraph *gr)
 Find the x value associated with the peak (maximum positive) in a TGraph. More...
 
Double_t getPeakVal (const TGraph *gr, int *index=0)
 Find the peak (maximum positive) in a TGraph. More...
 
Double_t getPeakVal (const TGraph *gr, Int_t firstBin, Int_t lastBin, int *index=0)
 Find the peak (maximum positive) in a TGraph. More...
 
Double_t getPeakSqVal (const TGraph *gr, int *index=0)
 Find the peak (v^2) in a TGraph. More...
 
void getPeakRmsSqVal (const TGraph *gr, Double_t &peak, Double_t &rms, Int_t *index=0)
 Find the peak (v^2) and RMS (of v^2) in a TGraph. More...
 
void getPeakRmsRectified (const TGraph *gr, Double_t &peak, Double_t &rms, Int_t *index=0)
 Find the peak (v) and RMS (of v) of a rectified TGraph. More...
 
TGraph * getSimplePowerEnvelopeGraph (const TGraph *gr)
 Returns the simple power envelope of a waveform. The waveform is first squared and then only local peaks are taken to form the power envelope. More...
 
TGraph * smoothFFT (const TGraph *gr, Int_t factor)
 Returns a smoothed FFT where N-bins are averaged to reduce variance. More...
 
TGraph * subtractGraphs (const TGraph *grA, const TGraph *grB)
 Returns the difference between two graphs (A-B). More...
 
TGraph * divideGraphs (const TGraph *grA, const TGraph *grB)
 Returns the ratio between two graphs (A/B). More...
 
TGraph * translateGraph (const TGraph *grWave, const Double_t deltaT)
 Returns a graph translated by deltaT. Such that t'=t+dt. More...
 
TGraph * ratioSubtractOneGraphs (const TGraph *grA, const TGraph *grB)
 Returns the one minus the ratio between two graphs (A/B). More...
 
TGraph * dbGraphs (const TGraph *grA, const TGraph *grB)
 Returns the ratio of two graphs as 10 *log(A/B) More...
 
TGraph * padWave (const TGraph *grA, Int_t padFactor)
 Zero pad a wave making it a factor of N longer. More...
 
TGraph * padWaveToLength (const TGraph *grA, Int_t newLength)
 Zero pad a wave making it up to newLength points. More...
 
TGraph * rectifyWave (const TGraph *gr, Int_t makeNeg=0)
 Rectify a waveform, optionally returning an all negative waveform. More...
 
Double_t bartlettWindow (Int_t j, Int_t n)
 The Bartlett window function (it's basically a triangle that peaks in the middle at 1 and goes to zero at the end points) More...
 
Double_t welchWindow (Int_t j, Int_t n)
 The Welch window function similar to the Bartlett window except that it falls off from the middle as dx^2 (ie. less steeply). More...
 
void takeDerivative (Int_t numPoints, const Double_t *inputX, const Double_t *inputY, Double_t *outputX, Double_t *outputY)
 This fucntions just calculates the simple bin by bin dv/dt derivative of the input data. More...
 
TGraph * getDerviative (const TGraph *grIn)
 This returns a TGraph which is the derivative of the input graph. More...
 
TGraph * getDerivative (const TGraph *grIn)
 This returns a TGraph which is the derivative of the input graph. More...
 
TGraph * simplePassBandFilter (const TGraph *grWave, Double_t minFreq, Double_t maxFreq)
 This returns a TGraph which has had a simple pass band filter applied. More...
 
TGraph * simpleNotchFilter (const TGraph *grWave, Double_t minFreq, Double_t maxFreq)
 This returns a TGraph which has had a simple notch band filter applied. More...
 
TGraph * multipleSimpleNotchFilters (const TGraph *grWave, Int_t numNotches, const Double_t minFreq[], const Double_t maxFreq[])
 This returns a TGraph which has had N simple notch band filters applied. More...
 
TGraph * cropWave (const TGraph *grWave, Double_t minTime, Double_t maxTime)
 This returns a TGraph which has been cropped in. More...
 
Double_t getWaveformSNR (const TGraph *gr)
 This returns the SNR ratio of the input waveform. More...
 
Double_t getWaveformSNR (const TGraph *gr, Double_t &peakToPeak, Double_t &rms)
 This returns the SNR ratio of the input waveform. More...
 
Double_t getWaveformPeak (const TGraph *gr)
 This returns the largest (i.e most positive, or least negative) value. More...
 
Double_t getEnvelopeSNR (const TGraph *gr)
 
Double_t getEnvelopeSNR (const TGraph *gr, Double_t &peakToPeak, Double_t &rms, Double_t &timeOfPeak)
 
Double_t simpleInterploate (Double_t x1, Double_t y1, Double_t x2, Double_t y2, Double_t x)
 Linear interpolate to find the value at some point near two known points. More...
 
double * FFTCorrelation (int waveformlength, const FFTWComplex *A, const FFTWComplex *B, FFTWComplex *work=0, int min_i=0, int max_i=0, int order=1)
 
void inPlaceShift (int N, double *x)
 
void rotate (TGraph *g, int rot)
 
void polySubtract (TGraph *g, int order=1)
 
double * directConvolve (int N, const double *x, int M, const double *h, double *y=0, int delay=0, DirectConvolveEdgeBehavior edge_behavior=ZEROES_OUTSIDE)
 
void wrap (size_t N, double *vals, double period=360)
 
void wrap (size_t N, double *vals, double period, double center)
 
void wrap (size_t N, float *vals, float period, float center)
 
void wrap (size_t N, float *vals, float period=360)
 
int fast_floor (double val)
 
double wrap (double val, double period, double center)
 
double wrap (double val, double period=360)
 
void unwrap (size_t N, double *vals, double period=360)
 
void unwrap (size_t N, float *vals, float period=360)
 
double getDt (const TGraph *g, int realN=0)
 
double evalEvenGraph (const TGraph *g, double x)
 
void applyWindow (TGraph *g, const FFTWindowType *w)
 
double randomRayleigh (double sigma=1, TRandom *rng=0)
 
double sinc (double x, double eps=0)
 
int loadWisdom (const char *file)
 
int saveWisdom (const char *file)
 
TGraph * welchPeriodogram (const TGraph *gin, int segment_size, double overlap_fraction=0.5, const FFTWindowType *window=&GAUSSIAN_WINDOW, bool truncate_extra=true, TGraph *gout=0)
 
double * lombScarglePeriodogramSlow (int N, const double *x, const double *y, int nfreqs, const double *freqs, double *answer=0)
 
TGraph * lombScarglePeriodogram (const TGraph *g, double dt=0, double oversample_factor=4, double high_factor=2, TGraph *replaceme=0, int extirpolation_factor=4)
 
TGraph * lombScarglePeriodogram (int N, double dt, const double *__restrict x, const double *__restrict y, double oversample_factor=4, double high_factor=2, TGraph *replaceme=0, int extirpolation_factor=4)
 
void stokesParameters (int N, const double *hpol, const double *hpol_hat, const double *vpol, const double *vpol_hat, double *Iavg=0, double *Qavg=0, double *Uavg=0, double *Vavg=0, double *Iins=0, double *Qins=0, double *Uins=0, double *Vins=0, bool only_max_instantaneous=true)
 
void dftAtFreq (const TGraph *g, double freq, double *phase, double *amp=0, double *real=0, double *imag=0)
 
void dftAtFreqAndMultiples (const TGraph *g, double freq, int nmultiples, double *phase, double *amp=0, double *real=0, double *imag=0)
 
FFTWComplexmakeMinimumPhase (int N, const double *G, double mindb=-100)
 
TGraph * makeMinimumPhase (const TGraph *g, double mindb=-100)
 
TGraph * makeFixedDelay (const TGraph *g, double delay=1)
 
double checkCausality (int N, const FFTWComplex *signal)
 
double * rmsEnvelope (int N, double W, const double *x, const double *y, double *out=0)
 
double * peakEnvelope (int N, double min_distance, const double *x, const double *y, double *out=0)
 
double * getCrossCov (int length, const double *oldY1, const double *oldY2)
 
double * getCrossCorr (int length, const double *oldY1, const double *oldY2)
 
TGraph * getCovGraph (const TGraph *gr1, const TGraph *gr2, int *zeroOffset=0)
 
TGraph * getCorrGraph (const TGraph *gr1, const TGraph *gr2, int *zeroOffset=0)
 
TGraph * getInterpCovGraph (const TGraph *gr1, const TGraph *gr2, double deltaT)
 
TGraph * getInterpCorrGraph (const TGraph *gr1, const TGraph *gr2, double deltaT)
 
TGraph * getInterpolatedGraphInvert (const TGraph *g, double dt=0, int nout=0)
 
TGraph * getInterpolatedGraphWeightedInvert (const TGraph *g, double dt=0, int nout=0)
 
TGraphErrors * getInterpolatedGraphSparseInvert (const TGraph *g, double dt=0, int nout=0, double max_dist=32, double eps=0, double weight_exp=0, double mu=1e-3, int regularization_order=0, double error_scale=1, TH2 *A=0)
 
void getInterpolatedGraphSparseInvert (const TGraph *in, TGraph *out, double max_dist=32, double eps=0, double weight_exp=0, double mu=1e-3, int regularization_order=0, double error_scale=1, TH2 *A=0)
 
TGraph * getInterpolatedGraphInvertLapped (const TGraph *g, double dt=0, int lapsize=64, int nout=0)
 
TGraph * getInterpolatedGraphInvertSplit (const TGraph *g, double dt=0, int splitsize=64, int overlap=8, int out=0)
 
TGraph * getInterpolatedGraphDFT (const TGraph *g, double dt=0, int nout=0, double maxF=0)
 
FFTWComplexgetUnevenDFT (const TGraph *g, double df, int npts)
 
TGraph * supersample (const TGraph *g, int supersample_factor, int sinc_radius=16)
 
TGraph * getInterpolatedGraphSincKernel (const TGraph *g, double dt=0, int nout=0)
 
TGraph * getInterpolatedGraphLagrange (const TGraph *g, double dt=0, double supersample=1)
 
double shannonWhitakerInterpolateValue (double t, const TGraph *regular_graph, int max_lobe=0, const FFTWindowType *win=0)
 
double shannonWhitakerInterpolateValueAndError (double t, const TGraphErrors *regular_graph, double *err, int max_lobe=0, const FFTWindowType *win=0)
 
double linearInterpolateValueAndError (double t, const TGraph *regular_graph, double *err=0)
 
double shannonWhitakerInterpolateDerivative (double t, const TGraph *regular_graph)
 

Detailed Description

Stuff I use that should maybe be in FFTtools?

A set of basic filter operations that serve as examples and also should be quite useful SimplePassBandFilter just cuts everything outside the pass band in fourier space like a brick wall, with all the ensuing acausal goodness.

Digital Filter Implementation.

Cosmin Deaconu cozzy.nosp@m.d@ki.nosp@m.cp.uc.nosp@m.hica.nosp@m.go.ed.nosp@m.u

This is the namespace that holds all the various useful FFT related functions. All of the functions canbe called using something like FFTtools::doYourFunkyThing(arg1,arg2).

Enumeration Type Documentation

direct convolution for kernel can either assume zero's outside boundaries or repetition of first and last values

Definition at line 693 of file FFTtools.h.

Function Documentation

void FFTtools::applyWindow ( TGraph *  g,
const FFTWindowType w 
)

applies window to graph

Definition at line 2713 of file FFTtools.cxx.

Double_t FFTtools::bartlettWindow ( Int_t  j,
Int_t  n 
)

The Bartlett window function (it's basically a triangle that peaks in the middle at 1 and goes to zero at the end points)

Parameters
jThe index to evaluate the window at.
nThe length of sample.
Returns
The value of the window function at index j.

Definition at line 21 of file FFTtools.cxx.

double FFTtools::checkCausality ( int  N,
const FFTWComplex signal 
)

Checks if a signal is causal or not by comparing the imaginary part to the hilbert transform of the real part and returning the square difference (scaled by N). A perfectly causal system will return 0, but numerical inaccuracy may make that not true.

Definition at line 3372 of file FFTtools.cxx.

TGraph * FFTtools::combineGraphsUsingFFTs ( Int_t  numGraphs,
TGraph **  grPtr,
const double *  theWeights = 0 
)

Returns the time domain result of a frequency domain sum of a number of TGraphs. In the sum each TGraph is weighted by a value. As of writing this documentation, I'm not sure why this would be interesting.

Parameters
numGraphsThe number of TGraphs to sum.
grPtrA pointer to an array of numGraphs TGraph pointers.
theWeightsAn optional array of weights with which to sum the TGraphs.
Returns
A pointer to a TGraph containing the inverse FFT of the weighted summed FFT of the input graphs.

Definition at line 570 of file FFTtools.cxx.

Double_t * FFTtools::combineValuesUsingFFTs ( Int_t  numArrays,
Double_t **  thePtrPtr,
Int_t  eachLength 
)

Returns the time domain result of a frequency domain sum of a number of arrays. As of writing this documentation, I'm not sure why this would be interesting.

Parameters
numArraysThe number of arrays to sum.
thePtrPtrA pointer to a two-dimensional array of doubles [numArrays][eachLength].
eachLengthThe length of each array.
Returns
The time domain result of the frequency domain summation.

Definition at line 538 of file FFTtools.cxx.

TGraph * FFTtools::convertMagnitudeToTimeDomain ( const TGraph *  inputMag)

Converts inputMag (linear units) to time domain by means of hilbert transform assuming R_signal = 1/sqrt(2) (R_mag - i R^_mag);.

Parameters
inputMagTGraph containg the inputMagnitude
Returns
TGraph of the time domain signal

Definition at line 2256 of file FFTtools.cxx.

TGraph * FFTtools::correlateAndAverage ( Int_t  numGraphs,
TGraph **  grPtrPtr 
)

This is designed for when you want to average a number of graphs of the same thing together. It uses a correlation to find the deltaT between graphs and then shifts the graphs and coherently sums them. The return is the average of the input graphs.

Parameters
numGraphsNumber of graphs to average
thePtrPtrA pointer to a two-dimensional array of TGraphs *grPtrArray[numGraphs].
Returns
A pointer to a graph containing the averaged waveform.

Definition at line 2461 of file FFTtools.cxx.

TGraph * FFTtools::correlateAndAveragePredeterminedZone ( Int_t  numGraphs,
TGraph **  grPtrPtr,
Int_t *  correlationBinPtr,
Int_t  binWiggleRoom = 30 
)

This is designed for when you want to average a number of graphs of the same thing together. It uses a correlation to find the deltaT between graphs and then shifts the graphs and coherently sums them. The return is the average of the input graphs This is intended to be used for a set of multiple waveforms where we have found some event-wide peak bin from previous correlation Using this method assures that the bin associated with the maximum correlation value is found within a limited range of the predetermined event wide peak bin

Parameters
numGraphsNumber of graphs to average
thePtrPtrA pointer to a two-dimensional array of TGraphs *grPtrArray[numGraphs].
correlationBinPtrA pointer to the pre-determined average correlation on an event-by-event basis
binWiggleRoomHow many bins either side of the bin associated with the peak value in the reduced range used to find the max correlation, with the default set to 20 bins
Returns
A pointer to a graph containing the averaged waveform.

Definition at line 2535 of file FFTtools.cxx.

TGraph * FFTtools::cropWave ( const TGraph *  grWave,
Double_t  minTime,
Double_t  maxTime 
)

This returns a TGraph which has been cropped in.

Parameters
grWaveThe input graph.
minTimeThe lower edge of the time window.
maxTimeThe upper edge of the time window.
Returns
The derviative of grWave.

Definition at line 2012 of file FFTtools.cxx.

TGraph * FFTtools::dbGraphs ( const TGraph *  grA,
const TGraph *  grB 
)

Returns the ratio of two graphs as 10 *log(A/B)

Parameters
grAA pointer to the first graph.
grBA pointer to the second graph.
Returns
A pointer to a TGraph containing 10 * log(A/B)

Definition at line 1682 of file FFTtools.cxx.

void FFTtools::dftAtFreq ( const TGraph *  g,
double  freq,
double *  phase,
double *  amp = 0,
double *  real = 0,
double *  imag = 0 
)

Compute the DFT term at a given frequency. Probably not want you want to do for a series of frequencies where one can take advantage of trig relations, etc.

Parameters
gWhat to look at
freqthe frequency
phasepointer where phase will be stored or 0 to not store anything
amppointer where amplitude will be stored or 0 to not store anything
realpointer where real will be stored or 0 to not store anything
imagpointer where imag will be stored or 0 to not store anything

Definition at line 3206 of file FFTtools.cxx.

void FFTtools::dftAtFreqAndMultiples ( const TGraph *  g,
double  freq,
int  nmultiples,
double *  phase,
double *  amp = 0,
double *  real = 0,
double *  imag = 0 
)

Compute the DFT term at a given frequency and n multiples of it. Timebase does not have to be even. Phase and amp should be big enough to store nmultiples or 0 if you don't want to fill one.

Parameters
gwaveform to look at
freqbase frequency
nmultiplesnumber of multiples to consider
phasepointer to array where phase will be stored, or 0 to not store anything.
amppointer to array where amp will be stored, or 0 to not store anything.
realpointer where real will be stored or 0 to not store anything
imagpointer where imag will be stored or 0 to not store anything

Definition at line 3034 of file FFTtools.cxx.

double * FFTtools::directConvolve ( int  N,
const double *  x,
int  M,
const double *  h,
double *  y = 0,
int  delay = 0,
DirectConvolveEdgeBehavior  edge_behavior = ZEROES_OUTSIDE 
)

convolution without FFT of x with kernel h. y should have same size as x (it is cropped symmetrically) if delay = 0, h will be centered around its middle sample, if - M/2 , will be purely causal, etc.

Definition at line 2757 of file FFTtools.cxx.

TGraph * FFTtools::divideGraphs ( const TGraph *  grA,
const TGraph *  grB 
)

Returns the ratio between two graphs (A/B).

Parameters
grAA pointer to the first graph.
grBA pointer to the second graph.
Returns
A pointer to a TGraph containing the ratio (A/B) of the input graphs.

Definition at line 1632 of file FFTtools.cxx.

FFTWComplex * FFTtools::doFFT ( int  length,
const double *  theInput 
)

Computes an FFT of an array of real numbers.

Parameters
lengthThe length of the input array.
theInputThe input array of length real numbers.
Returns
An array of length/2 + 1 complex numbers

Definition at line 436 of file FFTtools.cxx.

void FFTtools::doFFT ( int  length,
const double *  properly_aligned_input,
FFTWComplex properly_aligned_output 
)

Version of doFFT don't require copying of memory. If these are not aligned properly (i.e. allocated with fftw_malloc, memalign or equivalent), bad things might happen.

Definition at line 351 of file FFTtools.cxx.

double * FFTtools::doInvFFT ( int  length,
const FFTWComplex theInput 
)

Computes an inverse FFT.

Parameters
lengthThe length of the output array
theInputThe input array of complex numbers of (length/2 +1)
Returns
An array of length real numbers

Definition at line 482 of file FFTtools.cxx.

void FFTtools::doInvFFTClobber ( int  length,
FFTWComplex properly_aligned_input_that_will_likely_be_clobbered,
double *  properly_aligned_output 
)

Version of doInvFFt don't require copying of memory. If these are not aligned properly (i.e. allocated with fftw_malloc, memalign or equivalent), bad things might happen. Note that the input may be clobbered in this case.

Definition at line 408 of file FFTtools.cxx.

void FFTtools::doInvFFTNoClobber ( int  length,
const FFTWComplex properly_aligned_input,
double *  properly_aligned_output 
)

Version of doInvFFt that only requires copying of input. The input is copied to a properly-aligned temporary on the stack. If the output is not aligned properly (i.e. allocated with fftw_malloc, memalign or equivalent), bad things might happen. Things will probably be a bit faster if the input is aligned properly, as long as memcpy dispatches correctly.

Definition at line 372 of file FFTtools.cxx.

double FFTtools::evalEvenGraph ( const TGraph *  g,
double  x 
)

linearly interpolates value x in g . Similar to doing g->Eval() but much faster (since doesn't need to sort or binary search)

Definition at line 2881 of file FFTtools.cxx.

int FFTtools::fast_floor ( double  val)
inline

Faster than the libm floor, which cares about errno and the floating point environment and all that crap

Definition at line 721 of file FFTtools.h.

double * FFTtools::FFTCorrelation ( int  waveformlength,
const FFTWComplex A,
const FFTWComplex B,
FFTWComplex work = 0,
int  min_i = 0,
int  max_i = 0,
int  order = 1 
)

Correlation between two FFTs bandpasses between min_i and max_i using butterworth filter of order order work can be used for temporary to avoid allocation of new memory

Definition at line 2641 of file FFTtools.cxx.

double FFTtools::getAbs ( FFTWComplex theNum)

Returns the magnitude of a complex number.

Parameters
theNumThe complex number.
Returns
The magnitude of the complex number.

Definition at line 1388 of file FFTtools.cxx.

TGraph * FFTtools::getBoxCar ( const TGraph *  grWave,
Int_t  halfWidth 
)

Smooth graph using box car smoothing.

Parameters
grWaveA pointer to the input TGraph
halfWidthThe halfWidth in samples of the size of the box over which to smooth the input graph (eg. halfWidth=1 means that sample i is the average if i-1, i and i+1.
Returns
A pointer to a TGraph containing the smoothed graph.

Definition at line 984 of file FFTtools.cxx.

TGraph * FFTtools::getConvolution ( const TGraph *  grA,
const TGraph *  grB 
)

Convolution.

Parameters
grAA pointer to the input TGraph A.
grBA pointer to the input TGraph B
Returns
A pointer to the convolution of A and B stored in a TGraph, it is the users responsibility to delete this after use.

Definition at line 2309 of file FFTtools.cxx.

RFSignal * FFTtools::getConvolution ( RFSignal grA,
RFSignal grB 
)

Convolution.

Parameters
grAA pointer to the input RFSignal A.
grBA pointer to the input RFSignal B
Returns
A pointer to the convolution of A and B stored in a TGraph, it is the users responsibility to delete this after use.

Definition at line 2399 of file FFTtools.cxx.

double * FFTtools::getCorrelation ( const TGraph *  gr1,
const TGraph *  gr2,
int  firstIndex,
int  lastIndex 
)

Computes the correlation of two subsets of TGraphs.

Parameters
gr1The first TGraph in the correlation.
gr2The second TGraph in the correlation.
firstIndexThe index of the first sample in the sub traces.
lastIndexThe index of the last sample in the sub traces.
Returns
The correlation as an array of lastIndex-firstIndex real numbers.

Definition at line 902 of file FFTtools.cxx.

double * FFTtools::getCorrelation ( int  length,
const float *  oldY1,
const float *  oldY2 
)

Computes the correlation of two arrays.

Parameters
lengthThe length of the arrays
oldY1The first array in the correlation.
oldY2The second array in the correlation.
Returns
The correlation as an array of length real numbers.

Definition at line 852 of file FFTtools.cxx.

double * FFTtools::getCorrelation ( int  length,
const double *  oldY1,
const double *  oldY2 
)

Computes the correlation of two arrays.

Parameters
lengthThe length of the arrays
oldY1The first array in the correlation.
oldY2The second array in the correlation.
Returns
The correlation as an array of length real numbers.

Definition at line 868 of file FFTtools.cxx.

TGraph * FFTtools::getCorrelationGraph ( const TGraph *  gr1,
const TGraph *  gr2,
Int_t *  zeroOffset = 0 
)

Returns the correlation of two TGraphs.

Parameters
gr1The first input TGraph
gr2The second input TGraph
Returns
A pointer to a TGraph containing the correlation of gr1 and gr2.

Definition at line 754 of file FFTtools.cxx.

TGraph * FFTtools::getDerivative ( const TGraph *  grIn)

This returns a TGraph which is the derivative of the input graph.

Parameters
grInThe input graph.
Returns
The derviative of grIn.

Definition at line 1922 of file FFTtools.cxx.

TGraph * FFTtools::getDerviative ( const TGraph *  grIn)

This returns a TGraph which is the derivative of the input graph.

Parameters
grInThe input graph.
Returns
The derviative of grIn.

Definition at line 1917 of file FFTtools.cxx.

double FFTtools::getDt ( const TGraph *  g,
int  realN = 0 
)

computes reasonable dt from unevenly sampled graph

Definition at line 2734 of file FFTtools.cxx.

TGraph * FFTtools::getHilbertEnvelope ( const TGraph *  grWave)

The Hilbert envelope of the input TGraph. This is defined as e_i=sqrt(v_i^2 + h_i^2), where e_i, v_i and h_i are the i-th sample of the envelope, input graph and hilbert transform of the input graph repsectively.

Parameters
grWaveA pointer to the input TGraph
Returns
A pointer to a TGraph containing the Hilbert envelope function.

Definition at line 964 of file FFTtools.cxx.

TGraph * FFTtools::getHilbertTransform ( const TGraph *  grWave)

The Hilbert transform of the input TGraph.

Parameters
grWaveA pointer to the input TGraph
Returns
A pointer to a TGraph containing the Hilbert transform

Definition at line 951 of file FFTtools.cxx.

double * FFTtools::getHilbertTransform ( int  N,
const double *  y 
)

Hilbert transform on an array y is input,N is length

Definition at line 936 of file FFTtools.cxx.

TGraph * FFTtools::getInterpolatedCorrelationGraph ( const TGraph *  grIn1,
const TGraph *  grIn2,
Double_t  deltaT 
)

Returns the correlation of two interpolated TGraphs.

Parameters
grIn1The first input TGraph
grIn2The second input TGraph
deltaTThe desired time step for the interpolated graphs
Returns
A pointer to a TGraph containing the correlation of gr1 and gr2, after each is interpolated to have a timestep of deltaT.

std::cout << grCor->GetN() << "\n";

Definition at line 838 of file FFTtools.cxx.

TGraph * FFTtools::getInterpolatedGraph ( const TGraph *  grIn,
Double_t  deltaT 
)

Interpolation Routines that use ROOT::Math::Interpolator.

Parameters
grInA pointer to the input TGraph.
deltaTThe desired period (1/rate) of the interpolated waveform.
Returns
A pointer to ther interpolated TGraph, it is the users responsibility to delete this after use.

Definition at line 32 of file FFTtools.cxx.

TGraph * FFTtools::getInterpolatedGraphFreqDom ( const TGraph *  grIn,
Double_t  deltaT 
)

Interpolation Routines that zeropads the FFT.

Parameters
grInA pointer to the input TGraph.
deltaTThe desired period (1/rate) of the interpolated waveform.
Returns
A pointer to ther interpolated TGraph, it is the users responsibility to delete this after use.

Definition at line 84 of file FFTtools.cxx.

TGraph * FFTtools::getNormalisedCorrelationGraph ( const TGraph *  gr1,
const TGraph *  gr2,
Int_t *  zeroOffset = 0 
)

Returns the normalised correlation of two TGraphs.

Parameters
gr1The first input TGraph (must be zero meaned)
gr2The second input TGraph (must be zero meaned)
Returns
A pointer to a TGraph containing the correlation of gr1 and gr2 where each point is normalised by the number of valid samples in the correlation and by the product of the RMS of the input graphs.

Definition at line 709 of file FFTtools.cxx.

TGraph * FFTtools::getNormalisedCorrelationGraphTimeDomain ( const TGraph *  gr1,
const TGraph *  gr2,
Int_t *  zeroOffset = 0,
Int_t  useDtRange = 0,
Double_t  dtMin = -1000,
Double_t  dtMax = 1000 
)

Returns the normalised correlation of two TGraphs.

Parameters
gr1The first input TGraph (must be zero meaned)
gr2The second input TGraph (must be zero meaned)
zeroOffsetA pointer to an integer where the sample corresponding to zero offset will be stored
useDtRangeA flag to enable the setting of a limited range of deltat's to try
dtMinThe minimum delta-t to include in the correlation, the maximum delta-t to include in the correlation
Returns
A pointer to a TGraph containing the correlation of gr1 and gr2 where each point is normalised by the number of valid samples in the correlation and by the product of the RMS of the input graphs.

Definition at line 627 of file FFTtools.cxx.

Int_t FFTtools::getPeakBin ( const TGraph *  gr)

Find the peak (maximum positive) bin in a TGraph.

Parameters
grA pointer to the input TGraph
Returns
The index of the bin with peak value.

Definition at line 1456 of file FFTtools.cxx.

Int_t FFTtools::getPeakBin ( const TGraph *  gr,
Int_t  firstBin,
Int_t  lastBin 
)

Find the peak (maximum positive) bin in a TGraph.

Parameters
grA pointer to the input TGraph
firstBinThe first bin to include.
lastBinThe last bin to include.
Returns
The index of the bin with peak value.

Definition at line 1490 of file FFTtools.cxx.

void FFTtools::getPeakRmsRectified ( const TGraph *  gr,
Double_t &  peak,
Double_t &  rms,
Int_t *  index = 0 
)

Find the peak (v) and RMS (of v) of a rectified TGraph.

Parameters
grA pointer to the input TGraph.
peakA reference to a Double_t where the peak value will be stored.
rmsA reference to a Double_t where the rms value will be stored.
indexAn optional pointer in which the peak bin can be stored.

Definition at line 1586 of file FFTtools.cxx.

void FFTtools::getPeakRmsSqVal ( const TGraph *  gr,
Double_t &  peak,
Double_t &  rms,
Int_t *  index = 0 
)

Find the peak (v^2) and RMS (of v^2) in a TGraph.

Parameters
grA pointer to the input TGraph.
peakA reference to a Double_t where the peak value will be stored.
rmsA reference to a Double_t where the rms value will be stored.
indexAn optional pointer in which the peak bin can be stored.

Definition at line 1570 of file FFTtools.cxx.

Double_t FFTtools::getPeakSqVal ( const TGraph *  gr,
int *  index = 0 
)

Find the peak (v^2) in a TGraph.

Parameters
grA pointer to the input TGraph.
indexAn optional pointer in which the peak bin can be stored.
Returns
The peak (v^2)value.

Definition at line 1550 of file FFTtools.cxx.

Double_t FFTtools::getPeakVal ( const TGraph *  gr,
int *  index = 0 
)

Find the peak (maximum positive) in a TGraph.

Parameters
grA pointer to the input TGraph.
indexAn optional pointer in which the peak bin can be stored.
Returns
The peak value.

Definition at line 1510 of file FFTtools.cxx.

Double_t FFTtools::getPeakVal ( const TGraph *  gr,
Int_t  firstBin,
Int_t  lastBin,
int *  index = 0 
)

Find the peak (maximum positive) in a TGraph.

Parameters
grA pointer to the input TGraph.
firstBinThe first bin to include.
lastBinThe last bin to include.
indexAn optional pointer in which the peak bin can be stored.
Returns
The peak value.

Definition at line 1529 of file FFTtools.cxx.

Double_t FFTtools::getPeakXvalue ( const TGraph *  gr)

Find the x value associated with the peak (maximum positive) in a TGraph.

Parameters
grA pointer to the input TGraph
Returns
The x value of the bin with peak value.

Definition at line 1473 of file FFTtools.cxx.

TGraph * FFTtools::getSimplePowerEnvelopeGraph ( const TGraph *  gr)

Returns the simple power envelope of a waveform. The waveform is first squared and then only local peaks are taken to form the power envelope.

Parameters
grA pointer to the input TGraph.
Returns
A pointer to a TGraph containing the power envelope

Definition at line 1817 of file FFTtools.cxx.

Double_t FFTtools::getWaveformPeak ( const TGraph *  gr)

This returns the largest (i.e most positive, or least negative) value.

Parameters
grThe input graph.
Returns
The peak

Definition at line 2176 of file FFTtools.cxx.

Double_t FFTtools::getWaveformSNR ( const TGraph *  gr)

This returns the SNR ratio of the input waveform.

Parameters
grThe input graph.
Returns
The SNR of the input waveform, where S is half of the peak-to-peak and N is the RMS of the first 25 ns.

Definition at line 2081 of file FFTtools.cxx.

Double_t FFTtools::getWaveformSNR ( const TGraph *  gr,
Double_t &  peakToPeak,
Double_t &  rms 
)

This returns the SNR ratio of the input waveform.

Parameters
grThe input graph.
peakToPeakA reference to a double which will be set to half of the peak-to-peak
rmsA reference to a double which will be set to the RMS of teh first 25 samples
Returns
The SNR of the input waveform, where S is half of the peak-to-peak and N is the RMS of the first 25 samples.

Definition at line 2092 of file FFTtools.cxx.

void FFTtools::inPlaceShift ( int  N,
double *  x 
)

in place array rotation

Definition at line 2721 of file FFTtools.cxx.

Double_t FFTtools::integratePower ( const TGraph *  gr,
Int_t  firstBin = -1,
Int_t  lastBin = -1 
)

The integral of the power in a TGraph (normally a PSD) (i.e the sum of bin content*bin width)

Parameters
grA pointer to the input TGraph (normally a PSD)
firstBinThe first bin to include in the integral.
lastBinThe last bin to include in the integral.
Returns
The integral value.

Definition at line 1406 of file FFTtools.cxx.

Double_t FFTtools::integrateVoltageSquared ( const TGraph *  gr,
Int_t  firstBin = -1,
Int_t  lastBin = -1 
)

The integral of the v^2*dt in a waveform. Now works for unevenly sampled waeforms.

Parameters
grA pointer to the input TGraph
firstBinThe first bin to include in the integral.
lastBinThe last bin to include in the integral.
Returns
The value of the integral.

Definition at line 1436 of file FFTtools.cxx.

TGraph * FFTtools::interpolateCorrelateAndAverage ( Double_t  deltaTInt,
Int_t  numGraphs,
TGraph **  grPtrPtr 
)

This is designed for when you want to average a number of graphs of the same thing together. It uses a correlation to find the deltaT between graphs and then shifts the graphs and coherently sums them. The return is the average of the input graphs.

Parameters
deltaTIntThe time-step to interpolate to before doing the correlation and averaging
numGraphsNumber of graphs to average
thePtrPtrA pointer to a two-dimensional array of TGraphs *grPtrArray[numGraphs].
Returns
A pointer to a graph containing the averaged waveform.

Definition at line 2447 of file FFTtools.cxx.

TGraph * FFTtools::lombScarglePeriodogram ( const TGraph *  g,
double  dt = 0,
double  oversample_factor = 4,
double  high_factor = 2,
TGraph *  replaceme = 0,
int  extirpolation_factor = 4 
)

fast periodogoram (as in Press & Rybicki) . Implementation in Periodogram.cxx

Definition at line 250 of file Periodogram.cxx.

double * FFTtools::lombScarglePeriodogramSlow ( int  N,
const double *  x,
const double *  y,
int  nfreqs,
const double *  freqs,
double *  answer = 0 
)

"normal" lomb-scargle periodogram, not using the N log N algorithm in Press & Rybicki)

Definition at line 327 of file Periodogram.cxx.

TGraph * FFTtools::makeFixedDelay ( const TGraph *  g,
double  delay = 1 
)

Compute a fixed group delay version of a time-domain signal.

Parameters
gthe time-domain signal, probably needs to be zero-padded
phasethe phase (in radians) of the output
Returns
fixed phase version of g

Definition at line 3341 of file FFTtools.cxx.

TGraph * FFTtools::makeInverseInverseSpectrum ( const TGraph *  grWave)

Returns the inverse FFT of the FFT of the input TGraph. Seems pointless.

Parameters
grWaveA pointer to the input TGraph
Returns
A pointer to a TGraph containing the inverse FFT of the FFT of grWave.

Definition at line 919 of file FFTtools.cxx.

FFTWComplex * FFTtools::makeMinimumPhase ( int  N,
const double *  G,
double  mindb = -100 
)

Given a desired amplitude response G, makes a minimum phase version of the waveform. This takes advantage of the fact that for a minimum phase signal, the phase response equals the hilbert transform of the log magnitude.

Parameters
N,thelength of the amplitude response (typically the time domain length/ 2 + 1),
F,theamplitude response G(w) with w in normalized units (first component is dc, last is nyquist). This is in linear units.
mindb,anabsolute zero in the frequency response cannot be represented with a minimum phase filter. Anything with a (power) magnitude smaller than mindb will be made equal to mindb
Returns
A minimum phase version of the waveform.

Definition at line 3278 of file FFTtools.cxx.

TGraph * FFTtools::makeMinimumPhase ( const TGraph *  g,
double  mindb = -100 
)

Compute a minimum phase version of a time-domain signal.

Parameters
gthe time-domain signal, probably needs to be zero-padded
mindb,anabsolute zero in the frequency response cannot be represented with a minimum phase filter. Anything with a (power) magnitude smaller than mindb will be made equal to mindb
Returns
minimum phase version of g

Definition at line 3317 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrum ( const TGraph *  grWave)

Returns the power spectral density. Note the PSD is unormalised (or if you prefer is normalised to the sum squared amplitude of the time domain). See this short note for my terminology.

Parameters
grWaveThe input time domain waveform
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1031 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumMilliVoltsNanoSeconds ( const TGraph *  grWave)

Returns the power spectral density in dB units. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a millivolts-nanoseconds one.

Parameters
grWaveThe input time domain waveform with units of millivolts-nanoseconds.
Returns
A pointer to a TGraph containing the power spectrum in dB units, the frequency units are MHz.

Definition at line 1147 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumPeriodogram ( const TGraph *  grWave)

Returns the power spectral density. Note the PSD returned is the periodogram (or if you prefer is normalised to the mean squared amplitude of the time domain). See this short note for my terminology.

Parameters
grWaveThe input time domain waveform
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1070 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumVoltsSeconds ( const TGraph *  grWave)

Returns the power spectral density. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one.

Parameters
grWaveThe input time domain waveform with units of volts-seconds.
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1106 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumVoltsSecondsBartlett ( const TGraph *  grWave)

Returns the power spectral density of the input waveform convolved with a Bartlett Window. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one.

Parameters
grWaveThe input time domain waveform with units of volts-seconds.
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1017 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumVoltsSecondsdB ( const TGraph *  grWave)

Returns the power spectral density in dB units. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one.

Parameters
grWaveThe input time domain waveform with units of volts-seconds.
Returns
A pointer to a TGraph containing the power spectrum in dB units, the frequency units are MHz.

Definition at line 1286 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumVoltsSecondsPadded ( const TGraph *  grWave,
Int_t  padFactor = 4 
)

Returns the power spectral density of the input waveform zero-padded by some factor. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one.

Parameters
grWaveThe input time domain waveform with units of volts-seconds.
padFactorThe factor by which to zero pad the wave.
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1336 of file FFTtools.cxx.

TGraph * FFTtools::makePowerSpectrumVoltsSecondsPaddeddB ( const TGraph *  grWave,
Int_t  padFactor = 4 
)

Returns the power spectral density in dB units of the input waveform zero-padded by some factor. Note the PSD returned is normalised and divided by frequency bin width (or if you prefer it is normalised to the time-integral squared amplitude of the time domain and then divided by frequency bin width). See this short note for my terminology. As the name suggests this function expects the input waveform to be a volts-seconds one.

Parameters
grWaveThe input time domain waveform with units of volts-seconds.
padFactorThe factor by which to zero pad the wave.
Returns
A pointer to a TGraph containing the power spectrum in dB units with MHz as the frequency unit.

Definition at line 1346 of file FFTtools.cxx.

TGraph * FFTtools::makePSVSBartlettPaddedOverlap ( const TGraph *  grWave,
Int_t  padFactor = 4,
Int_t  numFreqs = 64 
)

Returns the power spectral density of the input waveform. In this one we first zero pad the waveform and then split it up into overlapping segments and convolve each segment with the Bartlett window before summing the resulting PSD's. As the name suggests this function expects the input waveform to be a volts-seconds one. No idea if this one actually works, or where I read oabout this crazy method which is supposed to reduce the variance of the PSD estimator.

Parameters
grWaveThe input time domain waveform with units of volts-seconds.
padFactorThe factor by which to zero pad the wave see padWave
numFreqsThe number of frequency bins required in the output, which is related to the length of overlapping segments.
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1856 of file FFTtools.cxx.

TGraph * FFTtools::makeRawPowerSpectrum ( const TGraph *  grWave)

Returns the power spectral density in completely unormalised unit (as in Parseval's theorem is not obeyed and there is an extra factor of N not removed form the PSD). See this short note for my terminology.

Parameters
grWaveThe input time domain waveform.
Returns
A pointer to a TGraph containing the power spectrum.

Definition at line 1354 of file FFTtools.cxx.

TGraph * FFTtools::multipleSimpleNotchFilters ( const TGraph *  grWave,
Int_t  numNotches,
const Double_t  minFreq[],
const Double_t  maxFreq[] 
)

This returns a TGraph which has had N simple notch band filters applied.

Parameters
grWaveThe input graph.
numNotchesThe number of notch regiosn
minFreqAn array of lower frequency of the notches.
maxFreqAn array upper frequency of the notch.
Returns
The derviative of grWave.

Definition at line 2038 of file FFTtools.cxx.

TGraph * FFTtools::normalizeWaveform ( TGraph *  inGraph)

stuff I want to add to FFTtools

Definition at line 752 of file AnitaTemplates.cc.

TGraph * FFTtools::padWave ( const TGraph *  grA,
Int_t  padFactor 
)

Zero pad a wave making it a factor of N longer.

Parameters
grAA pointer to the input graph.
padFactorThe factor by whcih to increas the length (eg. new length = factor * old length)
Returns
A pointer to the zero padded TGraph.

Definition at line 1745 of file FFTtools.cxx.

TGraph * FFTtools::padWaveToLength ( const TGraph *  grA,
Int_t  newLength 
)

Zero pad a wave making it up to newLength points.

Parameters
grAA pointer to the input graph.
newLengthThe new length the graph should be (it gets zero padded on either side)
Returns
A pointer to the zero padded TGraph.

Definition at line 1769 of file FFTtools.cxx.

double* FFTtools::peakEnvelope ( int  N,
double  min_distance,
const double *  x,
const double *  y,
double *  out = 0 
)

Peak envelope Approximates envelope by connecting peaks of |y|.

Parameters
Nthe number of points in input / output arrays
min_distancethe minimum distance between peaks. This should be close to the half-period of the carrier wave you expect.
xinput x values (times)
yinput y values (voltages)
outoutput y values (voltages) , or 0 to allocate new
Returns
output values
void FFTtools::polySubtract ( TGraph *  g,
int  order = 1 
)

does fit to polynomial of order order then subtracts it

Definition at line 2745 of file FFTtools.cxx.

double FFTtools::randomRayleigh ( double  sigma = 1,
TRandom *  rng = 0 
)

random variable from rayleigh distribution

Definition at line 2792 of file FFTtools.cxx.

TGraph * FFTtools::ratioSubtractOneGraphs ( const TGraph *  grA,
const TGraph *  grB 
)

Returns the one minus the ratio between two graphs (A/B).

Parameters
grAA pointer to the first graph.
grBA pointer to the second graph.
Returns
A pointer to a TGraph containing the one minus the ratio (1 - A/B) of the input graphs.

Definition at line 1652 of file FFTtools.cxx.

TGraph * FFTtools::rectifyWave ( const TGraph *  gr,
Int_t  makeNeg = 0 
)

Rectify a waveform, optionally returning an all negative waveform.

Parameters
grA pointer to the input graph.
makeNegAn optional parameter which if true will return a negative waveform.
Returns
A pointer to the rectified TGraph

Definition at line 1799 of file FFTtools.cxx.

double* FFTtools::rmsEnvelope ( int  N,
double  W,
const double *  x,
const double *  y,
double *  out = 0 
)

Avg RMS envelope Approximates envelope by computing rms in each window

Parameters
Nlength of input (and output) arrays
Wwindow length, in units of the x values (not all windows may contain the same length)
xinput x values
yinput y values
outoutput y values, if 0 output will be allocated
Returns
output
void FFTtools::rotate ( TGraph *  g,
int  rot 
)

inplace graph rotation

Definition at line 2627 of file FFTtools.cxx.

double FFTtools::shannonWhitakerInterpolateDerivative ( double  t,
const TGraph *  regular_graph 
)

Derivative of Shannon-Whitaker interpolation of a regular graph

Definition at line 733 of file RFInterpolate.cxx.

double FFTtools::shannonWhitakerInterpolateValue ( double  t,
const TGraph *  regular_graph,
int  max_lobe = 0,
const FFTWindowType win = 0 
)

Shannon-Whitaker interpolation of a regular graph

Definition at line 575 of file RFInterpolate.cxx.

Double_t FFTtools::simpleInterploate ( Double_t  x1,
Double_t  y1,
Double_t  x2,
Double_t  y2,
Double_t  x 
)

Linear interpolate to find the value at some point near two known points.

Parameters
x1x value of point 1
y1y value of point 1
x2x value of point 2
y2y value of point 2
xis the point of interest
Returns
y

Definition at line 2304 of file FFTtools.cxx.

TGraph * FFTtools::simpleNotchFilter ( const TGraph *  grWave,
Double_t  minFreq,
Double_t  maxFreq 
)

This returns a TGraph which has had a simple notch band filter applied.

Parameters
grWaveThe input graph.
minFreqThe lower frequency of the notch.
maxFreqThe upper frequency of the notch.
Returns
The derviative of grWave.

Definition at line 1975 of file FFTtools.cxx.

TGraph * FFTtools::simplePassBandFilter ( const TGraph *  grWave,
Double_t  minFreq,
Double_t  maxFreq 
)

This returns a TGraph which has had a simple pass band filter applied.

Parameters
grWaveThe input graph.
minFreqThe lowest frequency to pass.
maxFreqThe highest frequency to pass.
Returns
The derviative of grWave.

Definition at line 1938 of file FFTtools.cxx.

double FFTtools::sinc ( double  x,
double  eps = 0 
)

sinc function, if |x| < eps, sinc(x) = 1

Definition at line 2609 of file FFTtools.cxx.

TGraph * FFTtools::smoothFFT ( const TGraph *  gr,
Int_t  factor 
)

Returns a smoothed FFT where N-bins are averaged to reduce variance.

Parameters
grA pointer to the input TGraph.
factorThe factor by which to smooth (eg. 2 returns a PSD with half as many frequency bins).
Returns
A pointer to a TGraph containing the smoothed PSD.

Definition at line 1715 of file FFTtools.cxx.

TGraph * FFTtools::subtractGraphs ( const TGraph *  grA,
const TGraph *  grB 
)

Returns the difference between two graphs (A-B).

Parameters
grAA pointer to the first graph.
grBA pointer to the second graph.
Returns
A pointer to a TGraph containing the difference (A-B) between the input graphs.

Definition at line 1599 of file FFTtools.cxx.

Double_t FFTtools::sumPower ( const TGraph *  gr,
Int_t  firstBin = -1,
Int_t  lastBin = -1 
)

The linear sum of the power in a TGraph (normally a PSD)

Parameters
grA pointer to the input TGraph (normally a PSD)
firstBinThe first bin to include in the sum.
lastBinThe last bin to include in the sum.
Returns
The summed power

Definition at line 1393 of file FFTtools.cxx.

Double_t FFTtools::sumVoltageSquared ( const TGraph *  gr,
Int_t  firstBin,
Int_t  lastBin 
)

The sum of the voltage squared in a waveform.

Parameters
grA pointer to the input TGraph
firstBinThe first bin to include in the sum.
lastBinThe last bin to include in the sum.
Returns
The value of the sum.

Definition at line 1423 of file FFTtools.cxx.

void FFTtools::takeDerivative ( Int_t  numPoints,
const Double_t *  inputX,
const Double_t *  inputY,
Double_t *  outputX,
Double_t *  outputY 
)

This fucntions just calculates the simple bin by bin dv/dt derivative of the input data.

Parameters
numPointsThe number of points in the input data (the output will be numPoints-1 in length).
inputXThe input time values.
inputYThe input voltage values.
outputXThe output time values.
outputYThe output voltage values.

Definition at line 1905 of file FFTtools.cxx.

TGraph * FFTtools::translateGraph ( const TGraph *  grWave,
const Double_t  deltaT 
)

Returns a graph translated by deltaT. Such that t'=t+dt.

Parameters
grWaveA pointer to the input graph.
deltaTThe amount to move time axis by
Returns
A pointer to a TGraph containing the translated graph

Definition at line 1619 of file FFTtools.cxx.

void FFTtools::unwrap ( size_t  N,
double *  vals,
double  period = 360 
)

unwraps periodic array of doubles.

Definition at line 2857 of file FFTtools.cxx.

void FFTtools::unwrap ( size_t  N,
float *  vals,
float  period = 360 
)

unwraps periodic array of floats.

Definition at line 2845 of file FFTtools.cxx.

TGraph * FFTtools::welchPeriodogram ( const TGraph *  gin,
int  segment_size,
double  overlap_fraction = 0.5,
const FFTWindowType window = &GAUSSIAN_WINDOW,
bool  truncate_extra = true,
TGraph *  gout = 0 
)

Make a welch periodogram of evenly sampled input. A Bartlett periodogram can be emulated using overlap_fraction = 0 and window = &RECTANGULAR_WINDOW

Parameters
ginevenly sampled graph to estimate power spectrum of
segment_sizethe number of samples in each segment
overlap_fractionthe fraction of overlap for each segment. 0 is no overlap. 1 is non-sensical full overlap which will probably cause some calamity.
windowthe window to use for each segment
truncate_extraIf true, any partial segment is discarded. Otherwise, it will be zero-padded
goutif non-zero, this TGraph will be used for output
Returns
the Welch Periodogram

Definition at line 258 of file Periodogram.cxx.

Double_t FFTtools::welchWindow ( Int_t  j,
Int_t  n 
)

The Welch window function similar to the Bartlett window except that it falls off from the middle as dx^2 (ie. less steeply).

Parameters
jThe index to evaluate the window at
nThe length of sample.
Returns
The value of the window function at index j.

Definition at line 27 of file FFTtools.cxx.

void FFTtools::wrap ( size_t  N,
double *  vals,
double  period = 360 
)

wraps periodic array of doubles. center assumed to be period/2

Definition at line 2851 of file FFTtools.cxx.

void FFTtools::wrap ( size_t  N,
double *  vals,
double  period,
double  center 
)

wraps periodic array of doubles.

Definition at line 2872 of file FFTtools.cxx.

void FFTtools::wrap ( size_t  N,
float *  vals,
float  period,
float  center 
)

wraps periodic array of floats

Definition at line 2864 of file FFTtools.cxx.

void FFTtools::wrap ( size_t  N,
float *  vals,
float  period = 360 
)

wraps periodic array of floats. center assumed to be period/2

Definition at line 2839 of file FFTtools.cxx.

double FFTtools::wrap ( double  val,
double  period,
double  center 
)
inline

Wraps a value to be within period centered at center. This could be implemented in terms of the other wraps, but in practice it's better to inline since the period and center are usually known at compile them.

Definition at line 728 of file FFTtools.h.

double FFTtools::wrap ( double  val,
double  period = 360 
)
inline

Wraps a value to be within period. Center assumed to be period/2

Definition at line 732 of file FFTtools.h.