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... | |
RFSignal * | getConvolution (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... | |
FFTWComplex * | doFFT (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) |
FFTWComplex * | makeMinimumPhase (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) |
FFTWComplex * | getUnevenDFT (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) |
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 d@ki cp.uc hica go.ed 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).
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.
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)
j | The index to evaluate the window at. |
n | The length of sample. |
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.
numGraphs | The number of TGraphs to sum. |
grPtr | A pointer to an array of numGraphs TGraph pointers. |
theWeights | An optional array of weights with which to sum the TGraphs. |
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.
numArrays | The number of arrays to sum. |
thePtrPtr | A pointer to a two-dimensional array of doubles [numArrays][eachLength]. |
eachLength | The length of each array. |
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);.
inputMag | TGraph containg the inputMagnitude |
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.
numGraphs | Number of graphs to average |
thePtrPtr | A pointer to a two-dimensional array of TGraphs *grPtrArray[numGraphs]. |
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
numGraphs | Number of graphs to average |
thePtrPtr | A pointer to a two-dimensional array of TGraphs *grPtrArray[numGraphs]. |
correlationBinPtr | A pointer to the pre-determined average correlation on an event-by-event basis |
binWiggleRoom | How 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 |
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.
grWave | The input graph. |
minTime | The lower edge of the time window. |
maxTime | The upper edge of the time window. |
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)
grA | A pointer to the first graph. |
grB | A pointer to the second graph. |
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.
g | What to look at |
freq | the frequency |
phase | pointer where phase will be stored or 0 to not store anything |
amp | pointer where amplitude will be stored or 0 to not store anything |
real | pointer where real will be stored or 0 to not store anything |
imag | pointer 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.
g | waveform to look at |
freq | base frequency |
nmultiples | number of multiples to consider |
phase | pointer to array where phase will be stored, or 0 to not store anything. |
amp | pointer to array where amp will be stored, or 0 to not store anything. |
real | pointer where real will be stored or 0 to not store anything |
imag | pointer 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).
grA | A pointer to the first graph. |
grB | A pointer to the second graph. |
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.
length | The length of the input array. |
theInput | The input array of length real 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.
length | The length of the output array |
theInput | The input array of complex numbers of (length/2 +1) |
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.
|
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.
theNum | 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.
grWave | A pointer to the input TGraph |
halfWidth | The 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. |
Definition at line 984 of file FFTtools.cxx.
TGraph * FFTtools::getConvolution | ( | const TGraph * | grA, |
const TGraph * | grB | ||
) |
Convolution.
grA | A pointer to the input TGraph A. |
grB | A pointer to the input TGraph B |
Definition at line 2309 of file FFTtools.cxx.
Convolution.
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.
gr1 | The first TGraph in the correlation. |
gr2 | The second TGraph in the correlation. |
firstIndex | The index of the first sample in the sub traces. |
lastIndex | The index of the last sample in the sub traces. |
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.
length | The length of the arrays |
oldY1 | The first array in the correlation. |
oldY2 | The second array in the correlation. |
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.
length | The length of the arrays |
oldY1 | The first array in the correlation. |
oldY2 | The second array in the correlation. |
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.
gr1 | The first input TGraph |
gr2 | The second input TGraph |
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.
grIn | The input graph. |
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.
grIn | The input graph. |
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.
grWave | A pointer to the input TGraph |
Definition at line 964 of file FFTtools.cxx.
TGraph * FFTtools::getHilbertTransform | ( | const TGraph * | grWave | ) |
The Hilbert transform of the input TGraph.
grWave | A pointer to the input TGraph |
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.
grIn1 | The first input TGraph |
grIn2 | The second input TGraph |
deltaT | The desired time step for the interpolated graphs |
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.
grIn | A pointer to the input TGraph. |
deltaT | The desired period (1/rate) of the interpolated waveform. |
Definition at line 32 of file FFTtools.cxx.
TGraph * FFTtools::getInterpolatedGraphFreqDom | ( | const TGraph * | grIn, |
Double_t | deltaT | ||
) |
Interpolation Routines that zeropads the FFT.
grIn | A pointer to the input TGraph. |
deltaT | The desired period (1/rate) of the interpolated waveform. |
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.
gr1 | The first input TGraph (must be zero meaned) |
gr2 | The second input TGraph (must be zero meaned) |
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.
gr1 | The first input TGraph (must be zero meaned) |
gr2 | The second input TGraph (must be zero meaned) |
zeroOffset | A pointer to an integer where the sample corresponding to zero offset will be stored |
useDtRange | A flag to enable the setting of a limited range of deltat's to try |
dtMin | The minimum delta-t to include in the correlation, the maximum delta-t to include in the correlation |
Definition at line 627 of file FFTtools.cxx.
Int_t FFTtools::getPeakBin | ( | const TGraph * | gr | ) |
Find the peak (maximum positive) bin in a TGraph.
gr | A pointer to the input TGraph |
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.
gr | A pointer to the input TGraph |
firstBin | The first bin to include. |
lastBin | The last bin to include. |
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.
gr | A pointer to the input TGraph. |
peak | A reference to a Double_t where the peak value will be stored. |
rms | A reference to a Double_t where the rms value will be stored. |
index | An 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.
gr | A pointer to the input TGraph. |
peak | A reference to a Double_t where the peak value will be stored. |
rms | A reference to a Double_t where the rms value will be stored. |
index | An 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.
gr | A pointer to the input TGraph. |
index | An optional pointer in which the peak bin can be stored. |
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.
gr | A pointer to the input TGraph. |
index | An optional pointer in which the peak bin can be stored. |
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.
gr | A pointer to the input TGraph. |
firstBin | The first bin to include. |
lastBin | The last bin to include. |
index | An optional pointer in which the peak bin can be stored. |
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.
gr | A pointer to the input TGraph |
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.
gr | A pointer to the input TGraph. |
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.
gr | The input graph. |
Definition at line 2176 of file FFTtools.cxx.
Double_t FFTtools::getWaveformSNR | ( | const TGraph * | gr | ) |
This returns the SNR ratio of the input waveform.
gr | The input graph. |
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.
gr | The input graph. |
peakToPeak | A reference to a double which will be set to half of the peak-to-peak |
rms | A reference to a double which will be set to the RMS of teh 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)
gr | A pointer to the input TGraph (normally a PSD) |
firstBin | The first bin to include in the integral. |
lastBin | The last bin to include in the integral. |
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.
gr | A pointer to the input TGraph |
firstBin | The first bin to include in the integral. |
lastBin | The last bin to include in 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.
deltaTInt | The time-step to interpolate to before doing the correlation and averaging |
numGraphs | Number of graphs to average |
thePtrPtr | A pointer to a two-dimensional array of TGraphs *grPtrArray[numGraphs]. |
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.
g | the time-domain signal, probably needs to be zero-padded |
phase | the phase (in radians) of the output |
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.
grWave | A pointer to the input TGraph |
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.
N,the | length of the amplitude response (typically the time domain length/ 2 + 1), |
F,the | amplitude response G(w) with w in normalized units (first component is dc, last is nyquist). This is in linear units. |
mindb,an | absolute 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 |
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.
g | the time-domain signal, probably needs to be zero-padded |
mindb,an | absolute 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 |
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.
grWave | The input time domain waveform |
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.
grWave | The input time domain waveform with units of millivolts-nanoseconds. |
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.
grWave | The input time domain waveform |
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.
grWave | The input time domain waveform with units of volts-seconds. |
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.
grWave | The input time domain waveform with units of volts-seconds. |
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.
grWave | The input time domain waveform with units of volts-seconds. |
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.
grWave | The input time domain waveform with units of volts-seconds. |
padFactor | The factor by which to zero pad the wave. |
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.
grWave | The input time domain waveform with units of volts-seconds. |
padFactor | The factor by which to zero pad the wave. |
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.
grWave | The input time domain waveform with units of volts-seconds. |
padFactor | The factor by which to zero pad the wave see padWave |
numFreqs | The number of frequency bins required in the output, which is related to the length of overlapping segments. |
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.
grWave | The input time domain waveform. |
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.
grWave | The input graph. |
numNotches | The number of notch regiosn |
minFreq | An array of lower frequency of the notches. |
maxFreq | An array upper frequency of the notch. |
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.
grA | A pointer to the input graph. |
padFactor | The factor by whcih to increas the length (eg. new length = factor * old length) |
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.
grA | A pointer to the input graph. |
newLength | The new length the graph should be (it gets zero padded on either side) |
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|.
N | the number of points in input / output arrays |
min_distance | the minimum distance between peaks. This should be close to the half-period of the carrier wave you expect. |
x | input x values (times) |
y | input y values (voltages) |
out | output y values (voltages) , or 0 to allocate new |
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).
grA | A pointer to the first graph. |
grB | A pointer to the second graph. |
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.
gr | A pointer to the input graph. |
makeNeg | An optional parameter which if true will return a negative waveform. |
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
N | length of input (and output) arrays |
W | window length, in units of the x values (not all windows may contain the same length) |
x | input x values |
y | input y values |
out | output y values, if 0 output will be allocated |
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.
x1 | x value of point 1 |
y1 | y value of point 1 |
x2 | x value of point 2 |
y2 | y value of point 2 |
x | is the point of interest |
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.
grWave | The input graph. |
minFreq | The lower frequency of the notch. |
maxFreq | The upper frequency of the notch. |
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.
grWave | The input graph. |
minFreq | The lowest frequency to pass. |
maxFreq | The highest frequency to pass. |
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.
gr | A pointer to the input TGraph. |
factor | The factor by which to smooth (eg. 2 returns a PSD with half as many frequency bins). |
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).
grA | A pointer to the first graph. |
grB | A pointer to the second graph. |
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)
gr | A pointer to the input TGraph (normally a PSD) |
firstBin | The first bin to include in the sum. |
lastBin | The last bin to include in the sum. |
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.
gr | A pointer to the input TGraph |
firstBin | The first bin to include in the sum. |
lastBin | The last bin to include in 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.
numPoints | The number of points in the input data (the output will be numPoints-1 in length). |
inputX | The input time values. |
inputY | The input voltage values. |
outputX | The output time values. |
outputY | The 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.
grWave | A pointer to the input graph. |
deltaT | The amount to move time axis by |
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
gin | evenly sampled graph to estimate power spectrum of |
segment_siz | ethe number of samples in each segment |
overlap_fraction | the fraction of overlap for each segment. 0 is no overlap. 1 is non-sensical full overlap which will probably cause some calamity. |
window | the window to use for each segment |
truncate_extra | If true, any partial segment is discarded. Otherwise, it will be zero-padded |
gout | if non-zero, this TGraph will be used for output |
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).
j | The index to evaluate the window at |
n | The length of sample. |
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.
|
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.
|
inline |
Wraps a value to be within period. Center assumed to be period/2
Definition at line 732 of file FFTtools.h.