NCFS-Pack
A generic (astro)particle physics analysis framework
 
Loading...
Searching...
No Matches
NcDSP Class Reference

Various Digital Signal Processing (DSP) operations for (sequential) data samples. More...

#include "NcDSP.h"

Inheritance diagram for NcDSP:

Detailed Description

Various Digital Signal Processing (DSP) operations for (sequential) data samples.


Copyright(c) 2021 NCFS/IIHE, All Rights Reserved. *
*
Authors: The Netherlands Center for Fundamental Studies (NCFS). *
The Inter-university Institute for High Energies (IIHE). *
Website : http://www.iihe.ac.be *
Contact : Nick van Eijndhoven (nickve.nl@gmail.com) *
*
Contributors are mentioned in the code where appropriate. *
*
No part of this software may be used, copied, modified or distributed *
by any means nor transmitted or translated into machine language for *
commercial purposes without written permission by the IIHE representative. *
Permission to use the software strictly for non-commercial purposes *
is hereby granted without fee, provided that the above copyright notice *
appears in all copies and that both the copyright notice and this *
permission notice appear in the supporting documentation. *
This software is provided "as is" without express or implied warranty. *
The authors make no claims that this software is free of error, or is *
consistent with any particular standard of merchantability, or that it *
will meet your requirements for any particular application, other than *
indicated in the corresponding documentation. *
This software should not be relied on for solving a problem whose *
incorrect solution could result in injury to a person or loss of property. *
If you do use this software in such a manner, it is at your own risk. *
The authors disclaim all liability for direct or consequential damage *
resulting from your use of this software. *

// Class NcDSP
// Class to perform various Digital Signal Processing (DSP) operations on (sequential) data samples.
//
// For a description of most of the DSP techniques, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// The following discrete transformations (using the FFTW algorithms) are supported :
// Fourier (DFT), Hartley (DHT), Sine (DST) and Cosine (DCT).
//
// All transformation results have been normalized, such that
// the inverse transformation provides the original input spectrum.
//
// In addition to the above transformations, also convolution, correlation, periodogram,
// filter, Analog to Digital Converter (ADC), Digital to Analog Converter (DAC)
// and ADC-DAC chain transmission processors are provided.
//
// For details about the various operations and their options, please refer
// to the documentation in the corresponding member functions.
//
// Usage example for a discrete Fourier transform (DFT) :
// ------------------------------------------------------
//
// Int_t N=2048; // Number of samples
// Float_t fsample=3.2e9; // Sampling rate in Hz
// Float_t nu=3e8; // Signal frequency in Hz
//
// Float_t pi=acos(-1.);
// Float_t omega=2.*pi*nu; // Signal frequency in rad/sec
// Float_t step=1./fsample; // The time step of each sample
//
// // Fill the time domain sampled data
// Double_t tdata[N];
// Float_t t=0;
// Float_t y=0;
// for (Int_t i=0; i<N; i++)
// {
// y=cos(omega*t)+5.*sin(3.*omega*t);
// tdata[i]=y;
// t+=step;
// }
//
// NcDSP q;
// q.SetSamplingFrequency(fsample);
// q.Load(N,tdata);
//
// // Obtain the amplitudes in an array
// q.Fourier("R2C");
// TArrayD arr=q.GetData("AMP out");
//
// // Obtain the amplitudes in a histogram
// TH1F hist;
// q.Fourier("R2C",&hist,"AMP Hz");
//
// // Illustration of forward followed by inverse transformation
//
// // The original time domain input spectrum
// TH1F hist1;
// q.Fourier("R2C",&hist1,"t");
//
// // The forward transformation
// TH1F hist2;
// q.Fourier("R2C",&hist2,"AMP Hz");
//
// // Load the obtained frequency domain data as new input to obtain
// // the original time domain spectrum via the inverse transformation
// q.LoadResult();
//
// // The frequency domain input spectrum
// TH1F hist3;
// q.Fourier("C2R",&hist3,"AMP Hz");
//
// // Obtain the original time domain spectrum via the inverse transformation
// TH1F hist4;
// q.Fourier("C2R",&hist4,"t");
//
// TCanvas c1("c1","Forward");
// c1.Divide(1,2);
// c1.cd(1);
// hist1.Draw();
// c1.cd(2);
// hist2.Draw();
//
// TCanvas c2("c2","Inverse");
// c2.Divide(1,2);
// c2.cd(1);
// hist3.Draw();
// c2.cd(2);
// hist4.Draw();
//
// Usage example for Convolution :
// -------------------------------
//
// // The input signal
// const Int_t nx=17;
// Double_t x[nx]={0,0,1,1,1,1,1,1,1,1,1,1,1,1,1,0,0};
//
// // The system impulse response
// const Int_t nh=10;
// Double_t h[nh]={1,0.9,0.8,0.7,0.6,0.5,0.4,0.3,0.2,0.1};
//
// NcDSP q;
// q.SetSamplingFrequency(1e9); // Sampling frequency in Hz
// q.Load(nx,x); // Load the signal input data
// q.SetWaveform(nh,h); // Load the system response data
//
// printf(" Stored elements : Nsignal=%-i Nwaveform=%-i \n",q.GetN(),q.GetN(1));
//
// // Perform the Convolution
// TH1F hy;
// Int_t smin,smax;
// TArrayD y=q.Convolve(&hy,&smin,&smax);
//
// printf(" Convolution full coverage between sample numbers %-i and %-i \n",smin,smax);
//
// // Plot the various distributions
// TH1F hx("hx","Input signal",nx,0,nx);
// for (Int_t i=1; i<=nx; i++)
// {
// hx.SetBinContent(i,x[i-1]);
// }
//
// TH1F hh("hh","Impulse response",nh,0,nh);
// for (Int_t i=1; i<=nh; i++)
// {
// hh.SetBinContent(i,h[i-1]);
// }
//
// TCanvas c("c","Convolution");
// c.Divide(1,3);
// c.cd(1);
// hx.SetMarkerStyle(20);
// hx.Draw("PL");
// c.cd(2);
// hh.SetMarkerStyle(20);
// hh.Draw("PL");
// c.cd(3);
// hy.SetMarkerStyle(20);
// hy.Draw("PL");
//
// Usage example for cross-correlation :
// -------------------------------------
//
// // The input signal
// const Int_t nx=17;
// Double_t x[nx]={0,0,1,1,1,1,1,1,1,1,1,1,1,1,1,0,0};
//
// // The waveform pattern
// const Int_t nh=10;
// Double_t h[nh]={1,0.9,0.8,0.7,0.6,0.5,0.4,0.3,0.2,0.1};
//
// NcDSP q;
// q.SetSamplingFrequency(1e9); // Sampling frequency in Hz
// q.Load(nx,x); // Load the signal input data
// q.SetWaveform(nh,h); // Load the waveform data
//
// printf(" Stored elements : Nsignal=%-i Nwaveform=%-i \n",q.GetN(),q.GetN(1));
//
// // Perform the Cross-Correlation
// TH1F hy;
// Int_t smin,smax;
// TArrayD y=q.Correlate(&hy,&smin,&smax);
//
// printf(" Correlation full coverage between sample numbers %-i and %-i \n",smin,smax);
//
// // Plot the various distributions
// TH1F hx("hx","Input signal",nx,0,nx);
// for (Int_t i=1; i<=nx; i++)
// {
// hx.SetBinContent(i,x[i-1]);
// }
//
// TH1F hh("hh","Waveform pattern",nh,0,nh);
// for (Int_t i=1; i<=nh; i++)
// {
// hh.SetBinContent(i,h[i-1]);
// }
//
// TCanvas c("c","Correlation");
// c.Divide(1,3);
// c.cd(1);
// hx.SetMarkerStyle(20);
// hx.Draw("PL");
// c.cd(2);
// hh.SetMarkerStyle(20);
// hh.Draw("PL");
// c.cd(3);
// hy.SetMarkerStyle(20);
// hy.Draw("PL");
//
// Usage example for a Band Pass filter :
// --------------------------------------
//
// NcDSP q;
// Float_t fsample=1e9; // Sampling frequency in Hz
// q.SetSamplingFrequency(fsample);
// q.Load(nx,x); // Load some signal input data
// Double_t f1=200e6/fsample; // The lower bound (200 MHz) of the frequency band
// Double_t f2=300e6/fsample; // The upper bound (300 MHz) of the frequency band
// Int_t n=101; // The number of points describing the filter kernel
// TH1F hisf; // The frequency domain histogram
// TH1F hist; // The time domain histogram
// Int_t imin,imax; // The sample number boundaries of the signal for which the kernel is fully embedded
//
// // Perform the filtering
// q.FilterBandPass(f1,f2,n,&hisf,&hist,&imin,&imax);
//
// printf(" Filter kernel full coverage between sample numbers %-i and %-i \n",imin,imax);
//
// TCanvas cf("cf","cf");
// hisf.Draw();
// TCanvas ct("ct","ct");
// hist.Draw();
//
// Usage example for a Multi Band filter :
// ---------------------------------------
//
// NcDSP q;
// Float_t fsample=1e9; // Sampling frequency in Hz
// q.SetSamplingFrequency(fsample);
// q.Load(nx,x); // Load some signal input data
// Double_t f1=200e6/fsample; // The lower bound (200 MHz) of the 1st frequency band
// Double_t f2=300e6/fsample; // The upper bound (300 MHz) of the 1st frequency band
// Double_t f3=450e6/fsample; // The lower bound (450 MHz) of the 2nd frequency band
// Double_t f4=550e6/fsample; // The upper bound (550 MHz) of the 2nd frequency band
// TArrayF freqs(4); // Specification of the frequency bands
// freqs[0]=f1;
// freqs[1]=f2;
// freqs[2]=f3;
// freqs[3]=f4;
// Int_t n=101; // The number of points describing the filter kernel
// TH1F hisf; // The frequency domain histogram
// TH1F hist; // The time domain histogram
// Int_t imin,imax; // The sample number boundaries of the signal for which the kernel is fully embedded
//
// // Perform the filtering
// q.FilterMultiBand(freqs,n,&hisf,&hist,&imin,&imax);
//
// printf(" Filter kernel full coverage between sample numbers %-i and %-i \n",imin,imax);
//
// TCanvas cf("cf","cf");
// hisf.Draw();
// TCanvas ct("ct","ct");
// hist.Draw();
//
//
//--- Author: Nick van Eijndhoven, IIHE-VUB, Brussel, October 19, 2021 09:42Z
//- Modified: Nick van Eijndhoven, IIHE-VUB, Brussel, July 1, 2025 13:15Z

Definition at line 20 of file NcDSP.h.

Public Member Functions

 NcDSP (const char *name="", const char *title="")
 
 NcDSP (const NcDSP &q)
 
virtual ~NcDSP ()
 
TArrayL64 ADC (Int_t nbits, Double_t range, Double_t Vbias=0, TArray *Vsig=0, TH1 *hist=0, Int_t B=0, Int_t C=3) const
 
virtual TObject * Clone (const char *name="") const
 
TArrayD Convolve (TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Int_t shift=0)
 
TArrayD Correlate (TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Double_t *peak=0, TString norm="NONE")
 
void Cosine (Int_t type, TH1 *hist=0, TString sel="none")
 
TArrayD DAC (Int_t nbits, Double_t range, Double_t Vbias=0, TArray *adcs=0, TArray *peds=0, TH1 *hist=0, Int_t B=0, Int_t C=3) const
 
TArrayD Digitize (Int_t nbits, Double_t vcal, Int_t mode, TH1 *hist=0, Double_t *stp=0, Double_t *scale=0) const
 
TArrayD FilterBandPass (Double_t f1, Double_t f2, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Bool_t adaptn=kTRUE)
 
TArrayD FilterBandReject (Double_t f1, Double_t f2, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Bool_t adaptn=kTRUE)
 
TArrayD FilterHighPass (Double_t fcut, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Bool_t adaptn=kTRUE)
 
TArrayD FilterLowPass (Double_t fcut, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Bool_t adaptn=kTRUE)
 
TArrayD FilterMovingAverage (Int_t n, TString mode, TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, TH1 *hisf=0, Bool_t dB=kTRUE)
 
TArrayD FilterMultiBand (TArray &freqs, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Int_t *i1=0, Int_t *i2=0, Bool_t adaptn=kTRUE)
 
void Fourier (TString mode, TH1 *hist=0, TString sel="none")
 
TArrayD GetBandPassKernel (Double_t f1, Double_t f2, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Bool_t adaptn=kTRUE)
 
TArrayD GetBandRejectKernel (Double_t f1, Double_t f2, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Bool_t adaptn=kTRUE)
 
TArrayD GetData (TString mode) const
 
TArrayD GetHighPassKernel (Double_t fcut, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Bool_t adaptn=kTRUE)
 
TArrayD GetLowPassKernel (Double_t fcut, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Bool_t adaptn=kTRUE)
 
TArrayD GetMovingAverageKernel (Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0)
 
TArrayD GetMultiBandKernel (TArray &freqs, Int_t n, TH1 *hisf=0, Bool_t dB=kTRUE, TH1 *hist=0, Bool_t adaptn=kTRUE)
 
Int_t GetN (Int_t mode=0) const
 
Float_t GetSamplingFrequency () const
 
void Hartley (Int_t mode, TH1 *hist=0, TString sel="none")
 
void Hilbert (Int_t mode, TH1 *hist=0, TString sel="none")
 
void Load (Int_t n, Double_t *re, Double_t *im=0, Float_t f=-1)
 
void Load (NcSample *s, Int_t i, Float_t f=-1)
 
void Load (NcSample *s, TString name, Float_t f=-1)
 
void Load (TArray *re, TArray *im=0, Float_t f=-1)
 
void Load (TGraph *gr, Float_t f=-1)
 
void Load (TH1 *h, Float_t f=-1)
 
void LoadResult ()
 
TGraph Periodogram (TString tu, Double_t Tmin, Double_t Tmax, Int_t n, NcSample s, Int_t it, Int_t iy, Int_t idy=0, TF1 *Z=0, NcDevice *results=0) const
 
TGraph Periodogram (TString tu, Double_t Tmin, Double_t Tmax, Int_t n, NcSample s, TString namet, TString namey, TString namedy="-", TF1 *Z=0, NcDevice *results=0) const
 
TGraph Periodogram (TString tu, Double_t Tmin, Double_t Tmax, Int_t n, TArray &t, TArray &y, TArray *dy=0, TF1 *Z=0, NcDevice *results=0) const
 
TGraph Periodogram (TString tu, Double_t Tmin, Double_t Tmax, Int_t n, TGraph &g, Bool_t err=kFALSE, TF1 *Z=0, NcDevice *results=0) const
 
TGraph Periodogram (TString tu, Double_t Tmin, Double_t Tmax, Int_t n, TH1 &h, Bool_t err=kFALSE, TF1 *Z=0, NcDevice *results=0) const
 
TArrayD SampleAndHold (Int_t ns, TH1 *hist=0, Int_t loc=-1, Int_t jmin=0, Int_t jmax=-1) const
 
TArrayD SampleAndHold (TF1 f, Double_t step, Double_t vmin, Double_t vmax, TH1 *hist=0, Int_t loc=-1) const
 
TArrayD SampleAndSum (Int_t ns, TH1 *hist, Int_t jmin=0, Int_t jmax=-1) const
 
TArrayD SampleAndSum (TF1 f, Double_t step, Double_t vmin, Double_t vmax, TH1 *hist=0) const
 
void SetSamplingFrequency (Float_t f)
 
void SetWaveform (Int_t n, Double_t *h, Float_t f=-1)
 
void SetWaveform (NcSample *s, Int_t i, Float_t f=-1)
 
void SetWaveform (NcSample *s, TString name, Float_t f=-1)
 
void SetWaveform (TArray *h, Float_t f=-1)
 
void SetWaveform (TGraph *gr, Float_t f=-1)
 
void SetWaveform (TH1 *h, Float_t f=-1)
 
void Sine (Int_t type, TH1 *hist=0, TString sel="none")
 
TArrayD Transmit (Int_t nbits, Double_t range, Double_t Vbias=0, TArray *Vsig=0, TArray *peds=0, TH1 *hist=0, Int_t B=0, Int_t C=3) const
 

Protected Member Functions

void HistogramFilterFFT (TArray *h, TH1 *hisf, Bool_t dB, Bool_t kernel, TH1 *hist=0)
 
void HistogramTrafoResult (TString name, Int_t mode, TH1 *hist, TString sel)
 
void Reset ()
 

Protected Attributes

TArrayD fHisto
 
TArrayD fImIn
 
TArrayD fImOut
 
Bool_t fKeepOutput
 
Int_t fN
 
TString fNorm
 
Int_t fNwf
 
TVirtualFFT * fProc
 
TArrayD fReIn
 
TArrayD fReOut
 
Float_t fSample
 
TArrayD fWaveform
 

Constructor & Destructor Documentation

◆ NcDSP() [1/2]

NcDSP::NcDSP ( const char * name = "",
const char * title = "" )
// Default constructor.

Definition at line 293 of file NcDSP.cxx.

◆ ~NcDSP()

NcDSP::~NcDSP ( )
virtual
// Default destructor.

Definition at line 307 of file NcDSP.cxx.

◆ NcDSP() [2/2]

NcDSP::NcDSP ( const NcDSP & q)
// Copy constructor.

Definition at line 322 of file NcDSP.cxx.

Member Function Documentation

◆ ADC()

TArrayL64 NcDSP::ADC ( Int_t nbits,
Double_t range,
Double_t Vbias = 0,
TArray * Vsig = 0,
TH1 * hist = 0,
Int_t B = 0,
Int_t C = 3 ) const
// Processing of an Analog to Digital Converter (ADC).
// Construct from analog input signals the discrete quantized data values of an "nbits" ADC,
// based on the "range" for the analog signal and a bias voltage "Vbias" (see below).
// The analog input signals may be provided by the (optional) TArray "Vsig".
// In case the array "Vsig" is not provided, the stored waveform is used to provide the analog input signals.
// The resulting (integer) quantized ADC values are returned in a TArrayL64 object,
// without modification of the original waveform data.
//
// Note : Make sure to use the same units for "range", "Vbias" and the analog input signals.
//
// The number of available quantization levels is given by N=2^|nbits|, of which the lowest
// level represents the value 0. This yields for the quantized values (adc) the range [0,N-1].
// The mapping of an analog input voltage (Vin) to one of the quantization levels depends on the
// specified "range" and "Vbias" and whether we have a linear or Log ADC (see below).
//
// As outlined below, the range of the analog input voltage can be specified as the Full Scale
// voltage range (Vfs) corresponding to adc=N-1 or as the Reference voltage (Vref) corresponding to
// the hypothetical adc=N.
//
// The Least Significant Bit (LSB) represents the smallest analog input voltage interval
// that can reliably be resolved. In other words LSB=Voltage(adc=1)-Voltage(adc=0).
// For a linear ADC we have : LSB=Vfs/(N-1) or equivalently LSB=Vref/N.
// For a Log_B ADC (see below) we have : LSB=Vref*pow(B,-C)*(pow(B,C/N)-1).
//
// The formulas for the Vin->adc mapping with Vin=(Vbias+Vsig) are :
// Linear ADC : adc=Vin/LSB.
// Log_B ADC (see below) : adc=(N/C)*Log_B(pow(B,C)*Vin/Vref).
//
// Which leads to the following relations between Vref and Vfs :
// Linear ADC : Vref=Vfs+LSB.
// Log_B ADC (see below) : Vfs=Vref*pow(B,-C)*pow(B,(N-1)*C/N).
//
// The Dynamic Range (DR) is defined as the ratio of the largest and smallest input voltages
// that can reliably be resolved.
// Expressed in decibel we have : DR=20*log_10(Vfs/LSB) dB.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// and the publication by Y. Sundarasaradula et al. :
// "A 6-bit, Two-step, Successive Approximation Logarithmic ADC for Biomedical Applications"
// which is online available at https://spiral.imperial.ac.uk/bitstream/10044/1/44156/2/2016_ICECS_LogADC_Camera.pdf.
//
// Input arguments :
// -----------------
// nbits : Digital quantization of the input values will be performed using nbits.
// range >0 : The full scale voltage range (Vfs) of the analog signal that corresponds to adc=N-1.
// <0 : |range| is the reference voltage (Vref) of the analog signal that corresponds to the hypothetical adc=N.
// Vbias : The bias voltage that will be added to the analog input signal before digitization.
// Vsig : (Optional) array to contain the analog input signals.
// B >1 : Base for a Log ADC (e.g. B=10 emulates a Log_10 ADC).
// =0 : The ADC will be linear
// =1 : The ADC will be a Log_e ADC.
// Note : When B>0 all (Vbias+Vsig) input values should be positive.
// C : Code efficiency factor for a Log ADC.
// Rule of thumb : pow(B,-C) is about the smallest signal/|range| ratio that can be resolved.
// So, for a log_10 ADC with C=3, the smallest signal that can be resolved is about |range|/1000.
// Note : It is required that C>0.
//
// Optional output arguments :
// ---------------------------
// hist : Histogram with the digitized result
//
// Notes :
// -------
// 1) In case no "Vsig" array is provided and no waveform is present, just the ADC specs will be printed
// and in the TArrayL64 only the adc value corresponding to Vbias is returned.
// This allows to quickly investigate various scenarios without any data treatment.
// 2) Providing a "Vsig" array with different small (random) amplitudes allows to mimic variations of the bias or
// characteristics of the various switched capacitor array elements for the individual samplings.
// This resembles a so called "pedestal run" in real life Data Acquisition (DAQ) applications.
// For a linear ADC, the contents of the returned TArrayL64 can be used afterwards for pedestal subtraction.
//
// The maximum number of bits that is supported is 60 to guarantee identical functioning
// on all machines.
//
// In case of inconsistent input parameters, no digitization is performed and an empty TArrayL64 is returned.
//
// The default values are Vbias=0, Vsig=0, hist=0, B=0 and C=3.

Definition at line 2798 of file NcDSP.cxx.

◆ Clone()

TObject * NcDSP::Clone ( const char * name = "") const
virtual
// Make a deep copy of the current object and provide the pointer to the copy.
// This memberfunction enables automatic creation of new objects of the
// correct type depending on the object type, a feature which may be very useful
// for containers like NcEvent when adding objects in case the
// container owns the objects.

Definition at line 6194 of file NcDSP.cxx.

◆ Convolve()

TArrayD NcDSP::Convolve ( TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Int_t shift = 0 )
// Convolve the loaded input data x[] with the data contained in the (system response)
// waveform h[] and return the resulting data y[] in a TArrayD object.
// The input data x[] have to be entered as real numbers by one of the Load() member functions,
// whereas the (system response) waveform h[] has to be specified by SetWaveform().
// The provided data of x[] and h[] are not modified.
//
// In formula : y[]=x[]*h[].
//
// The convolution of two (time) series expresses how the shape of one
// is modified by the other, which makes it a versatile tool to describe
// system responses, digital filtering processing, superposition of various influences
// in physical systems etc.
// The result y[] can be regarded as the weighted sum (or pdf) of X+H, where
// X and H are two independent random variables with as pdf x[] and h[], respectively.
// Note that here x[]*h[] is identical to h[]*x[].
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// hist : (optional) Histogram with the convolution result
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting convolved data array for which "h" was fully immersed in the
// input (signal) data x[].
// In other words : The indices range [i1,i2] in the resulting convolved data array y[]
// corresponds to all the y[] elements for which the convolution was completely performed
// by using all the elements of h[].
// So, values of y[j] with j<i1 or j>i2 contain incomplete convolutions, and as such should not
// be considered as reliable, especially when x[] and/or h[] contain large variations.
//
// The argument "shift" determines whether the Xaxis values of the (optional) result histogram
// are shifted or that the array y[] is truncated in order to provide (a view of) the resulting array y[]
// such that it can be directly compared to the original array x[] (see notes 2) and 3) below).
//
// The convention is :
// shift = 0 --> No shift and the full array y[] is returned and shown "as is" in the (optional) histogram.
// 1 --> The full array y[] is returned and shown in the (optional) histogram with y[i1] at the center.
// This mode is particulary useful for correlation studies between two distributions [see Correlate()],
// since it will show both negative and positive values at the Xaxis of the (optional) histogram,
// and the self-correlation peak coincides with the value 0 at the Xaxis.
// 2 --> The array y[] contains only the modified x[] elements and is shown with y[0] at the start of the (optional) histogram.
// This mode is particulary useful for digital filtering processes (see the various filtering member functions),
// since it allows to directly compare the filtered and the unfiltered x[] distribution.
//
// Notes :
// -------
// 1) The sampling (rate) of h[] has to be the same as for the input data x[].
// 2) Array sizes of x[nx] and h[nh] will result in a convolved data array of size y[(nx+nh-1)],
// except for shift=2 for which the convolved data array has size y[nx].
// 3) For an absolute comparison between the x[] and y[] values, one should realize that
// the array sizes of x[] and y[] are different and that x[k] should be compared with y[k+i1].
// 4) The values i1 and i2 (if requested) are indicated by vertical dashed blue lines
// in the time domain histogram.
//
// The default values are hist=0, i1=0, i2=0 and shift=0.

Definition at line 2126 of file NcDSP.cxx.

◆ Correlate()

TArrayD NcDSP::Correlate ( TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Double_t * peak = 0,
TString norm = "NONE" )
// (Cross) Correlate the data contained in the waveform h[] with the loaded input data x[]
// and return the resulting data y[] in a TArrayD object.
// The input data x[] have to be entered as real numbers by one of the Load() member functions,
// whereas the waveform h[] has to be specified by SetWaveform().
// The provided data of x[] and h[] are not modified.
//
// In formula : y[]=h[]*x[].
//
// The cross-correlation is a measure of similarity of two (time) series as a function
// of the displacement of one relative to the other.
// The result y[] can be regarded as the weighted sum (or pdf) of X-H, where
// X and H are two independent random variables with as pdf x[] and h[], respectively.
// Note that here h[]*x[] is different from x[]*h[].
//
// Mathematically it is seen that the cross-correlation h[]*x[] is equivalent
// to the convolution of h[] and x[], but with the ordering of the elements
// of the distribution h[] reversed.
//
// In formula : Cross-correlation h[m]*x[n] is equivalent to Convolution h[-m]*x[n].
//
// This feature is used here to centralize the computation in the member function Convolve(),
// so that also this Correlate() processor will automatically profit from possible CPU speed
// improvements in the Convolution processor.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// hist : (optional) Histogram with the correlation result
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting correlated data array for which "h" was fully immersed in the
// input (signal) data.
// These values i1 and i2 (if requested) are indicated by vertical
// dashed blue lines in the histogram.
//
// The (optional) argument "peak" provides the location of the correlation peak, indicating
// the shift (aka lag) of the h[] pattern w.r.t. the best matching pattern in x[], with the convention
// that the autocorrelation peak falls at 0.
//
// The (optional) argument "norm" allows to obtain the normalized correlation values.
// When normalized, the resulting correlation values will be in the interval [-1,1],
// which facilitates an interpretation of the amount of correlation.
// However, it should be realized that different normalizations may result in different cross-correlation profiles,
// as outlined below.
//
// norm = "NONE" --> No normalization
// "GNCC" --> Same as norm="NONE" but the resulting correlation values will be divided by (|h|*|x|) of the total vectors h[] and [x]
// "NCC" --> NCC normalization, i.e. Sum{h[i]*x[i]/(|h|*|x|)}, where |h| and |x| are the vector lengths of the overlapping parts
// "ZNCC" --> ZNCC normalization, i.e. (1/N)*Sum{(h[i]-hmean)*(x[i]-xmean)/(sx*sh)}, where sx and sh are the STDs of the overlapping parts
//
// The "NCC" and "ZNCC" normalizations are preferred when looking for pattern matches, irrespective of the amplitudes.
// This will provide consistent results for different amplitudes of x[] and/or h[], but may yield incorrect results when
// looking for a (time) lag between two pulses on top of some background noise.
// The "GNCC" normalization is preferred in case amplitude matches, next to pattern matches, are also relevant.
// This is the preferred method when looking for a (time) lag between two pulses on top of some background noise,
// but in case of a bad SNR, the results may become unreliable.
//
// Note : The sampling (rate) of h has to be the same as for the loaded input data x[].
//
// The default values are hist=0, i1=0, i2=0, peak=0 and norm="NONE".

Definition at line 2402 of file NcDSP.cxx.

◆ Cosine()

void NcDSP::Cosine ( Int_t type,
TH1 * hist = 0,
TString sel = "none" )
// Perform a normalized 1-dimensional Discrete Cosine Transformation (DCT).
// Actually, this is just a regular Discrete Fourier Transformation (DFT)
// with only real input values which contain an even symmetry.
// Consequently, the resulting transformed array is also only real with even symmetry.
//
// Conventions :
// -------------
// N = The number of data elements
// Time domain array .... : X[]=X[0],...,X[N-1]
// Frequency domain array : Q[]=Q[0],...,Q[N-1]
//
// Cosine transform type 1 : Q[k]=(1/sqrt(2*(N-1)))*[X[0]+pow(-1,k)*X[N-1]+2*sum(n=1,n=N-2){X[n]*cos(pi*n*k/(N-1))}]
//
// Cosine transform type 2 : Q[k]=(1/sqrt(2N))*2*sum(n=0,n=N-1){X[n]*cos(pi*(n+0.5)*k/N)}]
//
// Cosine transform type 3 : Q[k]=(1/sqrt(2N))[X[0]+2*sum(n=0,n=N-1){X[n]*cos(pi*n*(k+0.5)/N)}]
//
// Cosine transform type 4 : Q[k]=(1/sqrt(2N))*2*sum(n=0,n=N-1){X[n]*cos(pi*(n+0.5)*(k+0.5)/N)}
//
// Notes :
// -------
// 1) The type 1 transformation is only defined for N>1 and is its own inverse.
// 2) The type 4 transformation is its own inverse.
// 3) The type 3 transformation is the inverse of type 2 (and vice versa).
// 4) The type 2 transformation is often called "the" DCT.
// 5) The type 3 transformation id often called "the" inverse DCT (IDCT).
//
// Input arguments :
// -----------------
// type : The type of transformation (i.e. 1,2,3 or 4) to be performed
// The inverse transformations are specified with the corresponding negative type value
// hist : (optional) Histogram with selected results
// sel : String to specify the representation of the result histogram
// "k" --> X-axis represents the index k in the frequency domain
// "f" --> X-axis represents the fraction f of the sampling rate in the frequency domain
// "Hz" --> X-axis represents the actual frequency in Hz in the frequency domain
// "n" --> X-axis represents the index n in the time domain
// "t" --> X-axis represents the actual time in seconds in the time domain
// "2" --> X-axis spans the full number of data points, instead of the usual (N/2)+1
//
// Note : The options "Hz" and "t" can only be used if the actual data acquisition sampling rate
// has been provided via the Load() memberfunction.
//
// Examples :
// ----------
// sel="f" will show the (N/2)+1 amplitudes as a function of the fractional sampling rate.
// sel="k 2" will show all N amplitudes as a function of the index k in the frequency domain.
//
// The default values are hist=0 and sel="none"

Definition at line 1334 of file NcDSP.cxx.

◆ DAC()

TArrayD NcDSP::DAC ( Int_t nbits,
Double_t range,
Double_t Vbias = 0,
TArray * adcs = 0,
TArray * peds = 0,
TH1 * hist = 0,
Int_t B = 0,
Int_t C = 3 ) const
// Processing of a Digital to Analog Converter (DAC).
// Reconstruct the analog signals based on the discrete quantized digital data from an "nbits" ADC,
// based on the "range" for the analog signal and a bias voltage "Vbias" or array "peds" of pedestal values (see below).
// The digital input signals may be provided by the (optional) TArray "adcs".
// In case the array "adcs" is not provided, the stored waveform is used to provide the digital input signals.
// The resulting analog values are returned in a TArrayD object, without modification of the original waveform data.
//
// Note : Make sure to use the same units for "range" and "Vbias".
//
// The number of digital quantization levels is given by N=2^|nbits|, of which the lowest
// level represents the value 0. This implies a range [0,N-1] for the various digital adc values.
// The correspondance of an analog input voltage (Vin) to one of the quantization levels depends on the
// specified "range" and "Vbias" and whether we have a linear or Log DAC (see below).
//
// As outlined below, the range of the analog input voltage can be specified as the Full Scale
// voltage range (Vfs) corresponding to adc=N-1 or as the Reference voltage (Vref) corresponding to
// the hypothetical adc=N.
//
// The Least Significant Bit (LSB) represents the smallest analog input voltage interval
// that can reliably be resolved. In other words LSB=Voltage(adc=1)-Voltage(adc=0).
// For a linear DAC we have : LSB=Vfs/(N-1) or equivalently LSB=Vref/N.
// For a Log_B DAC (see below) we have : LSB=Vref*pow(B,-C)*(pow(B,C/N)-1).
//
// The formulas for the adc->Vin mapping with Vin=(Vbias+Vsig) are :
// Linear DAC : Vin=adc*LSB.
// Log_B DAC (see below) : Vin=Vref*pow(B,-C)*pow(B,C*adc/N)
//
// Which implies the following relations between Vref and Vfs :
// Linear DAC : Vref=Vfs+LSB.
// Log_B DAC (see below) : Vfs=Vref*pow(B,-C)*pow(B,C*(N-1)/N).
//
// The Dynamic Range (DR) is defined as the ratio of the largest and smallest analog voltages
// that can reliably be resolved.
// Expressed in decibel we have : DR=20*log_10(Vfs/LSB) dB.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// and the publication by Y. Sundarasaradula et al. :
// "A 6-bit, Two-step, Successive Approximation Logarithmic ADC for Biomedical Applications"
// which is online available at https://spiral.imperial.ac.uk/bitstream/10044/1/44156/2/2016_ICECS_LogADC_Camera.pdf.
//
// Input arguments :
// -----------------
// nbits : Digitial quantization was performed using an nbits ADC.
// range >0 : The full scale voltage range (Vfs) of the analog signal that corresponds to adc=N-1.
// <0 : |range| is the reference voltage (Vref) of the analog signal that corresponds to the hypothetical adc=N.
// Vbias : The bias voltage that was added to the analog input signal before digitization.
// If specified, the resulting analog signals will be corrected for the bias voltage.
// For a linear DAC the correction via Vbias will only be performed if no pedestal array "peds"
// is specified (see below). If the array "peds" is specified, Vbias will be ignored for a linear DAC.
// For a Log DAC the pedestal values will never be used and bias correction may only be obtained via Vbias.
// adcs : (Optional) array to contain the digital input signals.
// peds : (Optional) array to contain the pedestal values for the individual digital input signals.
// The array "peds" must contain (at least) the same number of values as the number of digital input signals.
// If the array "peds" is specified, the pedestals will be subtracted from the corresponding digital input signals
// before the conversion to analog signals is performed in case of a linear DAC.
// For a Log DAC the pedestal values will never be used and bias correction may only be obtained via Vbias.
// B >1 : Base for a Log DAC (e.g. B=10 emulates a Log_10 DAC).
// =0 : The DAC will be linear
// =1 : The DAC will be a Log_e DAC.
// C : Code efficiency factor that was used for a Log ADC.
// Rule of thumb : pow(B,-C) is about the smallest signal/|range| ratio that can be resolved.
// So, for a Log_10 ADC with C=3, the smallest signal that can be resolved is about |range|/1000.
// Note : It is required that C>0.
//
// Optional output arguments :
// ---------------------------
// hist : Histogram with the analog result
//
// Notes :
// -------
// 1) In case no "adcs" array is provided and no waveform is present, just the DAC specs will be printed
// and in the TArrayD only the adc value corresponding to Vbias is returned.
// This allows to quickly investigate various scenarios without any data treatment.
//
// The maximum number of bits that is supported is 60 to guarantee identical functioning
// on all machines.
//
// In case of inconsistent input parameters, no processing is performed and an empty TArrayD is returned.
//
// The default values are Vbias=0, adcs=0, peds=0, hist=0, B=0 and C=3.

Definition at line 3025 of file NcDSP.cxx.

◆ Digitize()

TArrayD NcDSP::Digitize ( Int_t nbits,
Double_t vcal,
Int_t mode,
TH1 * hist = 0,
Double_t * stp = 0,
Double_t * scale = 0 ) const
// **************************************************************************************
// *** This function has become obsolete and is only kept for backward compatibility. ***
// *** Please refer to the new, more flexible memberfunction Transmit(). ***
// *** The user can also invoke the ADC and DAC processors individually by means of ***
// *** the corresponding memberfunctions ADC() and DAC(). ***
// **************************************************************************************
//
// Digitize the values of the stored waveform according to an "nbits" ADC.
// The resulting digitized values are returned in a TArrayD object,
// without modification of the original waveform data.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// nbits >0 : Digitization of the values will be performed using nbits.
// <0 : Digitization of the Log10 of the values will be performed using |nbits|.
// After digitization of the Log10 value, the digitized result (digval) is
// used to store the corresponding linear value via value=pow(10,digval).
// So, nbits<0 emulates a Log10 ADC to enhance the dynamic range.
// Note : When nbits<0 all values to be digitized should be positive.
// vcal : Range calibration value of the ADC according to "mode" as indicated below.
// mode =0 : Range for the digitized result (digval) will be set to [0,vcal] (or [vcal,0] if vcal<0).
// 1 : Full scale range for the digitized result (digval) will be set to [-|vcal|,|vcal|].
// 2 : A step size of |vcal| is used providing a digval interval of [0,scale] (or [-scale,0] if vcal<0).
// 3 : A step size of |vcal| is used providing a digval interval of [-scale,scale].
//
// Optional output arguments :
// ---------------------------
// hist : Histogram with the digitized result
// stp : The value of "step size".
// scale : The value of "scale".
//
// Notes :
// -------
// 1) The step size corresponds closely to the Least Significant Bit (LSB) precision for the
// digitized result (digval).
// For an n-bit ADC we have stepsize=range/(-1+2^n), whereas LSB=1/(2^n).
// 2) In case of a Log10 ADC, the value of "vcal" relates to the Log10 values.
// So, for a Log10 ADC, the "vcal" interval [-2,2] represents linear values [0.01,100].
// A step size of 0.0315 for an 8-bit Log10 ADC provides a Log10 interval of [-4,4]
// for mode=3, which corresponds to a linear value interval of [0.0001,10000].
// For mode=2 this would result in a Log10 interval of [0,8.0325], which corresponds
// to a linear value interval of [1,1.0777e8].
// A step size of 0.01 between the Log10 values corresponds to a multiplication factor
// of 10^0.01 (i.e. about 1.023) for the linear values at each step.
// 3) In case no waveform is present, just the specs of the specified ADC performance will be printed,
// but no digitization is performed.
// This allows to quickly investigate various scenarios without any data treatment.
//
// The maximum number of bits that is supported is 60 to guarantee identical functioning
// on all machines.
//
// In case of inconsistent input parameters, no digitization is performed and an empty TArrayD is returned.
//
// The default values are hist=0, stp=0 and scale=0.

Definition at line 2577 of file NcDSP.cxx.

◆ FilterBandPass()

TArrayD NcDSP::FilterBandPass ( Double_t f1,
Double_t f2,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Bool_t adaptn = kTRUE )
// Perform a Band Pass filter on the loaded input data x[] in the frequency band
// specified by "f1" and "f2" (see below) and a filter kernel consisting of "n" points (see below).
// The time domain result is returned in a TArrayD object and (optionally) in the histogram "hist",
// whereas the frequency domain result is returned in the (optional) histogram "hisf".
// The "hisf" amplitudes may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The original input data x[] are not modified.
// The input data x[] have to be entered as real numbers by one of the Load() member functions.
//
// Note :
// ------
// When the input data x[] contains en even (odd) number of samples, then an odd (even) value of "n"
// will result again in an even (odd) number of samples for the filtered data in the time domain.
//
// The implementation here is based on a combination of Blackman windowed-sinc Low Pass and High Pass filters,
// which is spectrally inverted.
//
// A Blackman Band Pass filter is an excellent frequency domain filter for passing
// only a certain band of frequencies.
// Large values of "n" will result in sharp transitions at the edges of the pass band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in less sharp transitions at the edges
// of the pass band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) at the edges of the pass band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) The returned TArrayD contains the filtered time domain data,
// which can be directly used as input for subsequent DSP processing.
// 2) Repeated application of this filter will further decrease the amplitudes outside the selected band.
// 3) Actually this filter represents a convolution with a time domain filter kernel
// that consists of a sinc pulse window of "n" points and an area of 1.
// 4) The performance in the time domain is very poor, since the (inverse) Fourier transform
// of a frequency domain rectangular pulse filter kernel results in a sinc function.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// f1 : The lower bound of the frequency band expressed as a fraction of the sampling frequency
// f2 : The upper bound of the frequency band expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<f1<0.5 and 0<f2<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below).
// hisf : (optional) Histogram with the filtered result (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filtered result in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting filtered data array for which the filter kernel was fully immersed
// in the input (signal) data.
// These values i1 and i2 (if requested) are indicated by vertical
// dashed blue lines in the time domain histogram.
//
// The default values are hisf=0, dB=kTRUE, hist=0, i1=0, i2=0 and adaptn=TRUE.

Definition at line 4236 of file NcDSP.cxx.

◆ FilterBandReject()

TArrayD NcDSP::FilterBandReject ( Double_t f1,
Double_t f2,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Bool_t adaptn = kTRUE )
// Perform a Band Reject filter on the loaded input data x[] in the frequency band
// specified by "f1" and "f2" (see below) and a filter kernel consisting of "n" points (see below).
// The time domain result is returned in a TArrayD object and (optionally) in the histogram "hist",
// whereas the frequency domain result is returned in the (optional) histogram "hisf".
// The "hisf" amplitudes may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The original input data x[] are not modified.
// The input data x[] have to be entered as real numbers by one of the Load() member functions.
//
// Note :
// ------
// When the input data x[] contains en even (odd) number of samples, then an odd (even) value of "n"
// will result again in an even (odd) number of samples for the filtered data in the time domain.
//
// The implementation here is based on a combination of Blackman windowed-sinc Low Pass and High Pass filters.
//
// A Blackman Band Reject filter is an excellent frequency domain filter for suppressing
// a certain band of frequencies.
// Large values of "n" will result in sharp transitions at the edges of the rejected band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in less sharp transitions at the edges
// of the rejected band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) at the edges of the rejected band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) The returned TArrayD contains the filtered time domain data,
// which can be directly used as input for subsequent DSP processing.
// 2) Repeated application of this filter will further decrease the amplitudes inside the selected band.
// 3) Actually this filter represents a convolution with a time domain filter kernel
// that consists of a sinc pulse window of "n" points and an area of 1.
// 4) The performance in the time domain is very poor, since the (inverse) Fourier transform
// of a frequency domain rectangular pulse filter kernel results in a sinc function.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// f1 : The lower bound of the frequency band expressed as a fraction of the sampling frequency
// f2 : The upper bound of the frequency band expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<f1<0.5 and 0<f2<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below).
// hisf : (optional) Histogram with the filtered result (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filtered result in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting filtered data array for which the filter kernel was fully immersed
// in the input (signal) data.
// These values i1 and i2 (if requested) are indicated by vertical
// dashed blue lines in the time domain histogram.
//
// The default values are hisf=0, dB=kTRUE, hist=0, i1=0, i2=0 and adaptn=kTRUE.

Definition at line 4385 of file NcDSP.cxx.

◆ FilterHighPass()

TArrayD NcDSP::FilterHighPass ( Double_t fcut,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Bool_t adaptn = kTRUE )
// Perform a High Pass filter on the loaded input data x[] with a frequency cut-off
// specified by "fcut" (see below) and a filter kernel consisting of "n" points (see below).
// The time domain result is returned in a TArrayD object and (optionally) in the histogram "hist",
// whereas the frequency domain result is returned in the (optional) histogram "hisf".
// The "hisf" amplitudes may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The original input data x[] are not modified.
// The input data x[] have to be entered as real numbers by one of the Load() member functions.
//
// Note :
// ------
// When the input data x[] contains en even (odd) number of samples, then an odd (even) value of "n"
// will result again in an even (odd) number of samples for the filtered data in the time domain.
//
// The implementation here is based on a spectrally inverted Blackman windowed-sinc Low Pass filter.
//
// A Blackman High Pass filter is an excellent frequency domain filter for seperating
// one band of frequencies from another.
// Large values of "n" will result in a sharp transition between the pass band and the stop band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in a less sharp transition between the
// pass band and the stop band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) between the pass band and the stop band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) The returned TArrayD contains the filtered time domain data,
// which can be directly used as input for subsequent DSP processing.
// 2) Repeated application of this filter will further decrease the amplitudes in the stop band.
// 3) Actually this filter represents a convolution with a time domain filter kernel
// that consists of a sinc pulse window of "n" points and an area of 1.
// 4) The performance in the time domain is very poor, since the (inverse) Fourier transform
// of the frequency domain step pulse filter kernel results in a sinc function.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// fcut : The cut-off frequency expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<fcut<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below).
// hisf : (optional) Histogram with the filtered result (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filtered result in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting filtered data array for which the filter kernel was fully immersed
// in the input (signal) data.
// These values i1 and i2 (if requested) are indicated by vertical
// dashed blue lines in the time domain histogram.
//
// The default values are hisf=0, dB=kTRUE, hist=0, i1=0, i2=0 and adaptn=kTRUE.

Definition at line 4091 of file NcDSP.cxx.

◆ FilterLowPass()

TArrayD NcDSP::FilterLowPass ( Double_t fcut,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Bool_t adaptn = kTRUE )
// Perform a Low Pass filter on the loaded input data x[] with a frequency cut-off
// specified by "fcut" (see below) and a filter kernel consisting of "n" points (see below).
// The time domain result is returned in a TArrayD object and (optionally) in the histogram "hist",
// whereas the frequency domain result is returned in the (optional) histogram "hisf".
// The "hisf" amplitudes may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The original input data x[] are not modified.
// The input data x[] have to be entered as real numbers by one of the Load() member functions.
//
// Note :
// ------
// When the input data x[] contains en even (odd) number of samples, then an odd (even) value of "n"
// will result again in an even (odd) number of samples for the filtered data in the time domain.
//
// The implementation here is based on the Blackman windowed-sinc filtering procedure.
//
// A Blackman Low Pass filter is an excellent frequency domain filter for seperating
// one band of frequencies from another.
// Large values of "n" will result in a sharp transition between the pass band and the stop band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in a less sharp transition between the
// pass band and the stop band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) between the pass band and the stop band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) The returned TArrayD contains the filtered time domain data,
// which can be directly used as input for subsequent DSP processing.
// 2) Repeated application of this filter will further decrease the amplitudes in the stop band.
// 3) Actually this filter represents a convolution with a time domain filter kernel
// that consists of a sinc pulse window of "n" points and an area of 1.
// 4) The performance in the time domain is very poor, since the (inverse) Fourier transform
// of the frequency domain step pulse filter kernel results in a sinc function.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// fcut : The cut-off frequency expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<fcut<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below)
// hisf : (optional) Histogram with the filtered result (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filtered result in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting filtered data array for which the filter kernel was fully immersed
// in the input (signal) data.
// These values i1 and i2 (if requested) are indicated by vertical
// dashed blue lines in the time domain histogram.
//
// The default values are hisf=0, dB-kTRUE, hist=0, i1=0, i2=0 and adaptn=kTRUE.

Definition at line 3946 of file NcDSP.cxx.

◆ FilterMovingAverage()

TArrayD NcDSP::FilterMovingAverage ( Int_t n,
TString mode,
TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
TH1 * hisf = 0,
Bool_t dB = kTRUE )
// Perform a Moving Average filter on the loaded input data x[] with averaging over "n" samples.
// The time domain result is returned in a TArrayD object and (optionally) in the histogram "hist".
// The frequency domain result is returned in the (optional) histogram "hisf", for which
// the amplitudes may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The original input data x[] are not modified.
// The input data x[] have to be entered as real numbers by one of the Load() member functions.
//
// A Moving Average filter is the optimal time domain filter for reducing random (noise) fluctuations,
// while retaining sharp step responses.
// Large values of "n" will result in a large noise reduction, but the edges of the steps
// become less sharp.
// On the contrary, small values of "n" will result in sharp step edges, but less noise reduction.
//
// Rule of thumb :
// ---------------
// The noise is reduced by a factor of sqrt(n) and the rise or fall of an edge
// is smeared over "n" samples.
//
// This filter may be invoked in two different modes, namely Recursion or Convolution.
// The resulting time domain data are equivalent, but are contained in TArrayD objects
// of different length, as outlined below..
//
// The Recursion mode :
// --------------------
// The straightforward procedure is expressed as :
//
// y[i]=(1/n)*Sum(j=0,n-1){x[i+j]} (1)
//
// To reduce rounding errors and computational speed, the following recursion relation is used :
//
// y[k]=y[k-1]-(x[k-1]/n)+(x[k+n-1]/n) (2)
//
// whereby the first term y[0] is calulated via equation (1).
//
// In case the number of elements of x[] is m, the number of elements of y[] will "only" be (m-n+1),
// such that the last element of y[] still has been completely computed from all corresponding "n" values of x[].
//
// The Convolution mode :
// ----------------------
// Based on the above, this filter actually represents a convolution with a filter kernel
// that consists of a rectangular pulse of "n" points of height 1/n (i.e. an area of 1).
// In this respect we can regard this filtering as a convolution y[]=x[]*h[] where h[]
// represents the "n" point rectangular pulse filter kernel.
// In this mode, an array x[] with m elements, will result in an array y[] of (m+n-1) elements.
// For further details see the documentation of Convolve().
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// n : The number of samples that will be averaged over
// mode : To invoke Recursion (mode="rec") or Convolution (mode="conv") processing.
// hist : (optional) Histogram with the filtered result in the time domain
// i1 : Optional argument (see below)
// i2 : Optional argument (see below)
// hisf : (optional) Histogram with the filtered result (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting filtered data array for which the filter kernel was fully immersed
// in the input (signal) data x[].
// In other words : The indices range [i1,i2] in the resulting filtered data array y[]
// corresponds to all the y[] elements for which the averaging was completely performed
// over all "n" samples.
// So, values of y[j] with j<i1 or j>i2 contain incomplete summations, and as such should not
// be considered as reliable, especially when x[] contains large variations.
// For an absolute comparison between the x[] and y[] values, one should realize that
// the array sizes of x[] and y[] are different and that x[k] should be compared with y[k+i1].
// The values i1 and i2 (if requested) are indicated by vertical dashed blue lines
// in the time domain histogram.
//
// Notes :
// -------
// 1) The returned TArrayD contains the filtered time domain data,
// which can be directly used as input for subsequent DSP processing.
// 2) The performance in the frequency domain is very poor, since the Fourier transform
// of the time domain rectangular pulse filter kernel results in a sinc function.
// 3) A Bayesian Block analysis (see the class NcBlocks) in general provides
// even better results by allowing "n" to vary dynamically based on the data pattern.
// However, the computation time involved for a Bayesian Block analysis
// is vastly larger than for a Moving Average filter.
//
// The default values are hist=0, i1=0 i2=0, hisf=0 and dB=kTRUE.

Definition at line 3698 of file NcDSP.cxx.

◆ FilterMultiBand()

TArrayD NcDSP::FilterMultiBand ( TArray & freqs,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Int_t * i1 = 0,
Int_t * i2 = 0,
Bool_t adaptn = kTRUE )
// Perform a MultiBand filter on the loaded input data x[] in various frequency bands as
// specified by the array "freqs" (see below) and filter kernels consisting of "n" points (see below).
// The time domain result is returned in a TArrayD object and (optionally) in the histogram "hist",
// whereas the frequency domain result is returned in the (optional) histogram "hisf".
// The "hisf" amplitudes may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The original input data x[] are not modified.
// The input data x[] have to be entered as real numbers by one of the Load() member functions.
//
// Note :
// ------
// When the input data x[] contains en even (odd) number of samples, then an odd (even) value of "n"
// will result again in an even (odd) number of samples for the filtered data in the time domain.
//
// The procedure is based on a convolution of the various provided Blackman
// single Low Pass and/or High Pass and/or Band Pass and/or Band Reject filters.
//
// Large values of "n" will result in sharp transitions at the edges of the specified bands,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in less sharp transitions at the edges
// of the specified bands, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) at the edges of the specified bands (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) The size of the TArrayD "freqs" should match exactly twice the number of
// frequency bands to be specified (see below).
// 2) The returned TArrayD contains the filtered time domain data,
// which can be directly used as input for subsequent DSP processing.
// 3) Specifying the same frequency band more than once will further suppress the unwanted frequencies.
// 4) The performance in the time domain is very poor, since the (inverse) Fourier transform
// of a frequency domain rectangular pulse filter kernel results in a sinc function.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// freqs : Array containing the lower and upper bounds of the frequency bands expressed as
// fractions of the sampling frequency.
// The array elements represent the various pairs [flow,fup] to define the frequency bands,
// ordered as (f1low,f1up,f2low,f2up,...).
// The following conventions are used :
// flow>0 and fup>0 --> Apply a Band Pass filter over [flow,fup]
// flow<0 and fup<0 --> Apply a Band Reject filter over [|flow|,|fup|]
// flow<0 and fup>0 --> Apply a Low Pass filter with fcut=fup
// flow>0 and fup<0 --> Apply a High Pass filter with fcut=flow
// n : The number of values in the corresponding filter kernels
// For best functionality this must be an odd integer (see adaptn below).
// hisf : (optional) Histogram with the filtered result (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filtered result in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The (optional) arguments i1 and i2 provide the range [i1,i2] in the
// resulting filtered data array for which the filter kernel was fully immersed
// in the input (signal) data.
// These values i1 and i2 (if requested) are indicated by vertical
// dashed blue lines in the time domain histogram.
//
// The default values are hisf=0, dB=kTRUE, hist=0, i1=0, i2=0 and adaptn=kTRUE.

Definition at line 4533 of file NcDSP.cxx.

◆ Fourier()

void NcDSP::Fourier ( TString mode,
TH1 * hist = 0,
TString sel = "none" )
// Perform a normalized 1-dimensional Discrete Fourier Transformation (DFT).
//
// Conventions :
// -------------
// N = The number of data elements
// Time domain array : X[]=X[0],...,X[N-1]
// Frequency domain array : Q[]=Q[0],...,Q[N-1]
//
// Fourier transform : Q[k]=(1/sqrt(N))*sum(n=0,n=N-1){X[n]*exp(-i*2pi*(k/N)*n)}
//
// Inverse Fourier transform : X[n]=(1/sqrt(N))*sum(k=0,k=N-1){Q[k]*exp(i*2pi*(n/N)*k)}
//
// Input arguments :
// -----------------
// mode : "R2C" --> Perform a real-input to complex-output discrete Fourier transformation
// "C2R" --> Perform the inverse transformation of "R2C"
// "C2C" --> Perform a complex-input to complex-output discrete Fourier transformation
// "C2CI" --> Perform the inverse of "C2C"
// hist : (optional) Histogram with selected results
// sel : String to specify the contents and representation of the result histogram
// "RE" --> Y-axis shows the values of the real (re) components
// "IM" --> Y-axis shows the values of the imaginary (im) components
// "AMP" --> Y-axis shows the values of the amplitudes, i.e. sqrt(re*re+im*im)
// "dB" --> Y-axis shows the values of the amplitudes, i.e. sqrt(re*re+im*im), in decibel
// "PHIR" --> Y-axis shows the values of the phases, i.e. arctan(im/re), in radians
// "PHID" --> Y-axis shows the values of the phases, i.e. arctan(im/re), in degrees
// "CPHI" --> Y-axis shows the values of cos(phase), i.e. cos(arctan(im/re))
// "SPHI" --> Y-axis shows the values of sin(phase), i.e. sin(arctan(im/re))
// "k" --> X-axis represents the index k in the frequency domain
// "f" --> X-axis represents the fraction f of the sampling rate in the frequency domain
// "Hz" --> X-axis represents the actual frequency in Hz in the frequency domain
// "n" --> X-axis represents the index n in the time domain
// "t" --> X-axis represents the actual time in seconds in the time domain
// "2" --> X-axis spans the full number of data points, instead of the usual (N/2)+1
//
// Note : The options "Hz" and "t" can only be used if the actual data acquisition sampling rate
// has been provided via the Load() memberfunction.
//
// Examples :
// ----------
// sel="AMP f" will show the (N/2)+1 amplitudes as a function of the fractional sampling rate.
// sel="dB f" will show the (N/2)+1 amplitudes in dB as a function of the fractional sampling rate.
// sel="RE k 2" will show all N real components as a function of the index k in the frequency domain.
//
// Note : If only the X-axis variable is specified, the Y-axis will show "RE" in the corresponding domain.
// So, sel="t" is equivalent to sel="RE t" etc...
//
// The default values are hist=0 and sel="none".

Definition at line 1092 of file NcDSP.cxx.

◆ GetBandPassKernel()

TArrayD NcDSP::GetBandPassKernel ( Double_t f1,
Double_t f2,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Bool_t adaptn = kTRUE )
// Provide via the returned TArrayD an n-point time domain Band Pass Filter kernel
// for the frequency band [f1,f2] expressed as fractions of the sampling frequency.
// The optional argument "hisf" may be used to obtain a histogram of the kernel
// in the frequency domain with the amplitude in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
//
// The implementation here is based on a combination of Blackman windowed-sinc Low Pass and High Pass filters,
// which is spectrally inverted.
//
// A Blackman Band Pass filter is an excellent frequency domain filter for passing
// only a certain band of frequencies.
// Large values of "n" will result in sharp transitions at the edges of the pass band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in less sharp transitions at the edges
// of the pass band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) at the edges of the pass band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) Repeated application of this filter will further decrease the amplitudes outside the selected band.
// 2) Actually this filter kernel consists of a time domain sinc pulse window of "n" points
// and an area of 1.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// f1 : The lower bound of the frequency band expressed as a fraction of the sampling frequency
// f2 : The upper bound of the frequency band expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<f1<0.5 and 0<f2<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below).
// hisf : (optional) Histogram with the filter kernel (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filter kernel in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The default values are hisf=0, dB=kTRUE, hist=0 and adaptn=kTRUE.

Definition at line 5030 of file NcDSP.cxx.

◆ GetBandRejectKernel()

TArrayD NcDSP::GetBandRejectKernel ( Double_t f1,
Double_t f2,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Bool_t adaptn = kTRUE )
// Provide via the returned TArrayD an n-point time domain Band Reject Filter kernel
// for the frequency band [f1,f2] expressed as fractions of the sampling frequency.
// The optional argument "hisf" may be used to obtain a histogram of the kernel
// in the frequency domain with the amplitude in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
//
// The implementation here is based on a combination of Blackman windowed-sinc Low Pass and High Pass filters.
//
// A Blackman Band Reject filter is an excellent frequency domain filter for suppressing
// a certain band of frequencies.
// Large values of "n" will result in sharp transitions at the edges of the rejected band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in less sharp transitions at the edges
// of the rejected band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) at the edges of the rejected band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) Repeated application of this filter will further decrease the amplitudes inside the selected band.
// 2) Actually this filter kernel consists of a time domain sinc pulse window of "n" points
// and an area of 1.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// f1 : The lower bound of the frequency band expressed as a fraction of the sampling frequency
// f2 : The upper bound of the frequency band expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<f1<0.5 and 0<f2<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer.
// hisf : (optional) Histogram with the filter kernel (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filter kernel in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The default values are hisf=0, dB=kTRUE, hist=0 and adaptn=kTRUE.

Definition at line 5153 of file NcDSP.cxx.

◆ GetData()

TArrayD NcDSP::GetData ( TString mode) const
// Provide a selected set of data.
//
// Input argument :
// ----------------
// sel : String to specify the contents of the provided data array
// "RE" --> The values of the real (re) components
// "IM" --> The values of the imaginary (im) components
// "AMP" --> The amplitudes, i.e. sqrt(re*re+im*im)
// "PHIR" --> The phases, i.e. arctan(im/re), in radians
// "PHID" --> The phases, i.e. arctan(im/re), in degrees
// "in" --> The values of the input data are provided
// "out" --> The values of the output data are provided
// "Hist" --> The contents of the selected transformation histogram (if any) are provided.
// "Wave" --> The values of the stored (system response) waveform data are provided
//
// Examples :
// ----------
// sel="AMP out" will provide all the N amplitudes of the resulting data after transformation.
// sel="RE in" will provide all the N real components of the input data.
// sel="Wave" will provide all the amplitudes of the stored (system response) waveform

Definition at line 1028 of file NcDSP.cxx.

◆ GetHighPassKernel()

TArrayD NcDSP::GetHighPassKernel ( Double_t fcut,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Bool_t adaptn = kTRUE )
// Provide via the returned TArrayD an n-point time domain High Pass Filter kernel
// with cut-off frequency "fcut" expressed as a fraction of the sampling frequency.
// The optional argument "hisf" may be used to obtain a histogram of the kernel
// in the frequency domain with the amplitude in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
//
// The implementation here is based on a spectrally inverted Blackman windowed-sinc Low Pass filter.
//
// A Blackman High Pass filter is an excellent frequency domain filter for seperating
// one band of frequencies from another.
// Large values of "n" will result in a sharp transition between the pass band and the stop band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in a less sharp transition between the
// pass band and the stop band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) between the pass band and the stop band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) Repeated application of this filter will further decrease the amplitudes in the stop band.
// 2) Actually this filter kernel consists of a time domain sinc pulse window of "n" points
// and an area of 1.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// fcut : The cut-off frequency expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<fcut<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below).
// hisf : (optional) Histogram with the filter kernel (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filter kernel in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The default values are hisf=0, dB=kTRUE, hist=0 and adaptn=kTRUE.

Definition at line 4911 of file NcDSP.cxx.

◆ GetLowPassKernel()

TArrayD NcDSP::GetLowPassKernel ( Double_t fcut,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Bool_t adaptn = kTRUE )
// Provide via the returned TArrayD an n-point time domain Low Pass Filter kernel
// with cut-off frequency "fcut" expressed as a fraction of the sampling frequency.
// The optional argument "hisf" may be used to obtain a histogram of the kernel
// in the frequency domain with the amplitude in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
//
// The implementation here is based on the Blackman windowed-sinc filtering procedure.
//
// A Blackman Low Pass filter is an excellent frequency domain filter for seperating
// one band of frequencies from another.
// Large values of "n" will result in a sharp transition between the pass band and the stop band,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in a less sharp transition between the
// pass band and the stop band, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) between the pass band and the stop band (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) Repeated application of this filter will further decrease the amplitudes in the stop band.
// 2) Actually this filter kernel consists of a time domain sinc pulse window of "n" points
// and an area of 1.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// fcut : The cut-off frequency expressed as a fraction of the sampling frequency
// For proper functionality one should choose 0<fcut<0.5, because of the
// underlying Fourier transform.
// n : The number of values in the filter kernel
// For best functionality this must be an odd integer (see adaptn below)
// hisf : (optional) Histogram with the filter kernel (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filter kernel in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The default values are hisf=0, dB=kTRUE, hist=0 and adaptn=kTRUE.

Definition at line 4776 of file NcDSP.cxx.

◆ GetMovingAverageKernel()

TArrayD NcDSP::GetMovingAverageKernel ( Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0 )
// Provide via the returned TArrayD an n-point time domain Moving Average Filter kernel.
// The optional argument "hisf" may be used to obtain a histogram of the frequency domain kernel
// with the amplitude in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
//
// A Moving Average filter is the optimal time domain filter for reducing random (noise) fluctuations,
// while retaining sharp step responses.
// Large values of "n" will result in a large noise reduction, but the edges of the steps
// become less sharp.
// On the contrary, small values of "n" will result in sharp step edges, but less noise reduction.
//
// Rule of thumb :
// ---------------
// The noise is reduced by a factor of sqrt(n) and the rise or fall of an edge
// is smeared over "n" samples.
//
// Notes :
// -------
// 1) Actually this filter kernel consists of a time domain rectangular pulse of "n" points
// and an area of 1.
// 2) The performance in the frequency domain is very poor, since the Fourier transform
// of the time domain rectangular pulse filter kernel results in a sinc function.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// n : The number of samples that will be averaged over
// hisf : (optional) Histogram with the filter kernel (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filter kernel in the time domain
//
// The default values are hisf=0, dB=kTRUE and hist=0.

Definition at line 4688 of file NcDSP.cxx.

◆ GetMultiBandKernel()

TArrayD NcDSP::GetMultiBandKernel ( TArray & freqs,
Int_t n,
TH1 * hisf = 0,
Bool_t dB = kTRUE,
TH1 * hist = 0,
Bool_t adaptn = kTRUE )
// Provide via the returned TArrayD a time domain Multi Band Filter kernel in various frequency bands
// specified by the array "freqs" (see below) and filter kernels consisting of "n" points each (see below).
// The optional argument "hisf" may be used to obtain a histogram of the kernel
// in the frequency domain with the amplitude in decibel (dB=kTRUE) or linear (dB=kFALSE).
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
//
// The procedure is based on a convolution of the various provided Blackman
// single Low Pass and/or High Pass and/or Band Pass and/or Band Reject filters.
//
// Large values of "n" will result in sharp transitions at the edges of the specified bands,
// but may result in long(er) computation times.
// On the contrary, small values of "n" will result in less sharp transitions at the edges
// of the specified bands, but result in short(er) computation times.
//
// Rule of thumb :
// ---------------
// The transition bandwidth (BW) at the edges of the specified bands (aka roll-off)
// may be approximated as BW=4/n, where BW is expressed as a fraction of the sampling frequency.
//
// Notes :
// -------
// 1) The size of the TArrayD "freqs" should match exactly twice the number of
// frequency bands to be specified (see below).
// 2) Specifying the same frequency band more than once will further suppress the unwanted frequencies.
//
// For details about the procedure, please refer to the excellent textbook :
// "The Scientist and Engineer's Guide to Digital Signal Processing" by Steven W. Smith,
// which is online available at http://www.dspguide.com/pdfbook.htm
//
// Input arguments :
// -----------------
// freqs : Array containing the lower and upper bounds of the frequency bands expressed as
// fractions of the sampling frequency.
// The array elements represent the various pairs [flow,fup] to define the frequency bands,
// ordered as (f1low,f1up,f2low,f2up,...).
// The following conventions are used :
// flow>0 and fup>0 --> Apply a Band Pass filter over [flow,fup]
// flow<0 and fup<0 --> Apply a Band Reject filter over [|flow|,|fup|]
// flow<0 and fup>0 --> Apply a Low Pass filter with fcut=fup
// flow>0 and fup<0 --> Apply a High Pass filter with fcut=flow
// In case flow=0 or fup=0 the pair [flow,fup] is neglected.
// n : The number of points in the corresponding filter kernels
// For best functionality this must be an odd integer (see adaptn below)
// hisf : (optional) Histogram with the filter kernel (amplitudes) in the frequency domain
// dB : Amplitudes of "hisf" are represented in decibel (kTRUE) or linear (kFALSE)
// hist : (optional) Histogram with the filter kernel in the time domain
// adaptn : If adaptn=kTRUE an even value of "n" will be increased by 1 to obtain an odd value
//
// The default values are hisf=0, dB=kTRUE, hist=0 and adaptn=kTRUE.

Definition at line 5278 of file NcDSP.cxx.

◆ GetN()

Int_t NcDSP::GetN ( Int_t mode = 0) const
// Provide the number of data elements (to be) processed.
//
// mode : 0 --> Provide the number of input data elements entered via Load()
// 1 --> Provide the number of input data elements entered via SetWaveform()
//
// The default value is mode=0 for backward compatibility.

Definition at line 1008 of file NcDSP.cxx.

◆ GetSamplingFrequency()

Float_t NcDSP::GetSamplingFrequency ( ) const
// Provide the current value of the DAQ sampling frequency in Hz.

Definition at line 384 of file NcDSP.cxx.

◆ Hartley()

void NcDSP::Hartley ( Int_t mode,
TH1 * hist = 0,
TString sel = "none" )
// Perform a normalized 1-dimensional Discrete Hartley Transformation (DHT).
// Actually, a DHT is closely related to a Discrete Fourier Transformation (DFT)
// with only real input values.
// Consequently, the resulting transformed array is also only real.
//
// Indicating in the frequency domain the DHT data elements by H[k] and the
// DFT data elements by F[k], we have the following relations :
//
// 1) Re(F[k])=(H[k]+H[N-k])/2 and Im(F[k])=(H[N-k]-H[k])/2
// 2) H[k]=Re((1+i)*F[k])
//
// Conventions :
// -------------
// N = The number of data elements
// Time domain array : X[]=X[0],...,X[N-1]
// Frequency domain array : Q[]=Q[0],...,Q[N-1]
//
// Hartley transform : Q[k]=(1/sqrt(N))*sum(n=0,n=N-1){X[n]*[cos(2pi*(k/N)*n)+sin(2pi*(k/N)*n)]}
//
// Inverse Hartley transform : X[n]=(1/sqrt(N))*sum(n=0,k=N-1){Q[k]*[cos(2pi*(n/N)*k)+sin(2pi*(n/N)*k)]}
//
// So, it is seen that the Hartley transform is its own inverse.
//
// Input arguments :
// -----------------
// mode : >0 --> Perform the forward X[n]->Q[k] Hartley transformation
// <0 --> Perform the backward Q[k]->X[n] Hartley transformation
// hist : (optional) Histogram with selected results
// sel : String to specify the representation of the result histogram
// "k" --> X-axis represents the index k in the frequency domain
// "f" --> X-axis represents the fraction f of the sampling rate in the frequency domain
// "Hz" --> X-axis represents the actual frequency in Hz in the frequency domain
// "n" --> X-axis represents the index n in the time domain
// "t" --> X-axis represents the actual time in seconds in the time domain
// "2" --> X-axis spans the full number of data points, instead of the usual (N/2)+1
//
// Note : The options "Hz" and "t" can only be used if the actual data acquisition sampling rate
// has been provided via the Load() memberfunction.
//
// Examples :
// ----------
// sel="f" will show the (N/2)+1 amplitudes as a function of the fractional sampling rate in the frequency domain.
// sel="k 2" will show all N amplitudes as a function of the index k in the frequency domain.
//
// The default values are hist=0 and sel="none"

Definition at line 1240 of file NcDSP.cxx.

◆ Hilbert()

void NcDSP::Hilbert ( Int_t mode,
TH1 * hist = 0,
TString sel = "none" )
// Perform a 1-dimensional Discrete Hilbert Transformation (DHIT).
// Actually, a DHIT is closely related to a Discrete Fourier Transformation (DFT)
// with only real input values, followed by a phase shift and inverse DFT.
// Consequently, the transformation results in a time domain array of only real values.
//
// Practical application :
// -----------------------
// Since the time domain input consists only of real values, the input signal can
// be fully described by cos() terms only.
// Consider in the time domain the following time series (x) that results from
// a carrier wave that is modulated both in amplitude (A) and phase (P) :
//
// x(t)=A(t)cos(P(t))
//
// The Hilbert transformation allows to extract A(t) and P(t) from the time series x(t)
// via the so called complex analytic signal Z[] as outlined below.
// The values of A(t) provide the so-called the Hilbert envelope, and will always by positive or 0.
// The instantaneous angular frequency omega(t) is obtained via omega(t)=dP(t)/dt, whereas
// the values of cos(P(t)) provide the so-called Temporal Fine Structure (TFS).
//
// Conventions :
// -------------
// N = The number of data elements
// Time domain input array .............................. : X[]=X[0],...,X[N-1]
// Fourier transform array DFT(X) in the frequency domain : Q[]=Q[0],...,Q[N-1]
// Hilbert transform array DHIT(X) in the time domain ... : HT[]=HT[0],...,HT[N-1]
// Analytic signal array in the time domain ............. : Z[]=Z[0],...,Z[N-1]
//
// In this memberfunction, the Hilbert transformation action is executed in the frequency domain,
// after which the time domain result HT[] is obtained via an inverse Fourier transformation.
//
// Indicating in the frequency domain the (intermediate) DHIT data elements by H[k],
// we have the following relations :
//
// H[k] = -i*Q[k] for positive frequencies (which in our convention corresponds to 1<=k<=N/2)
// 0 for frequency 0 (which in our convention corresponds to k=0)
// i*Q[k] for negative frequencies (which in our convention corresponds to k>N/2)
//
// In other words : If Q[k] is represented as Q[k]=a+b*i then we have :
//
// H[k] = b-a*i for positive frequencies
// 0 for frequency 0
// -b+a*i for negative frequencies
//
// From this we obtain via an inverse DFT on H[k] the Hilbert transform HT[n] in the time domain.
//
// The instantaneous amplitude A(t) and phase P(t) can be obtained via the so-called analytic signal Z(t) :
//
// Z[n]=X[n]+i*HT[n] which can be represented as Z[n]=c+d*i
//
// where A[n]=|Z[n]|=sqrt(c*c+d*d) and P[n]=arctan(d/c).
//
// Inverse transformation :
// ------------------------
// From the above procedure it is seen that DHIT(HT[])=-X[], which is used to perform the inverse transformation HT[]->X[].
// For the result of the inverse transformation HT[]->X[] we keep the definition Z[n]=X[n]+i*HT[n], which allows
// to also extract the instantaneous amplitude A(t) and phase P(t) from the inverse transformation result.
// The inverse Hilbert transformation has been implemented for completeness, but in practice this will not be used often.
//
// Input arguments :
// -----------------
// mode : >0 --> Perform the forward X[n]->Z[n] Hilbert transformation
// <0 --> Perform the backward HT[n]->X[n] Hilbert transformation
// hist : (optional) Histogram with selected results
// sel : String to specify the representation of the result histogram
// "X" --> Y-axis shows the values of the real X[] components
// "HT" --> Y-axis shows the values of the real HT[] components
// "AMP" --> Y-axis shows the values of the amplitudes, i.e. sqrt(re*re+im*im)
// "dB" --> Y-axis shows the values of the amplitudes, i.e. sqrt(re*re+im*im), in decibel
// "PHIR" --> Y-axis shows the values of the phases, i.e. arctan(im/re), in radians
// "PHID" --> Y-axis shows the values of the phases, i.e. arctan(im/re), in degrees
// "CPHI" --> Y axis shows the values of cos(phase), i.e. cos(arctan(im/re))
// This is also called the Temporal Fine Structure (TFS).
// "OMEGA" --> Y axis shows the values of angular frequency in rad/sec, i.e. d(arctan(im/re))/dt
// "n" --> X-axis represents the index n in the time domain
// "t" --> X-axis represents the actual time in seconds in the time domain
//
// Note : The option "t" can only be used if the actual data acquisition sampling rate
// has been provided via the Load() memberfunction.
//
// Examples :
// ----------
// sel="AMP t" will show the amplitudes as a function of the actual time in seconds.
// sel="PHID n" will show the phases in degrees as a function of the index n in the time domain.
// sel="HT n" will show the values of HT[] as a function of the index n in the time domain.
//
// Note :
// ------
// sel="t" for the forward transformation is equivalent to sel="HT t"
// sel="t" for the inverse transformation is equivalent to sel="X t"
// sel="n" for the forward transformation is equivalent to sel="HT n"
// sel="n" for the inverse transformation is equivalent to sel="X n"
//
// The default values are hist=0 and sel="none"

Definition at line 1568 of file NcDSP.cxx.

◆ HistogramFilterFFT()

void NcDSP::HistogramFilterFFT ( TArray * h,
TH1 * hisf,
Bool_t dB,
Bool_t kernel,
TH1 * hist = 0 )
protected
// Internal member function to provide a filter kernel or filter result histogram in the frequency domain,
// based on the time domain filter kernel or time domain filter result data contained in "h".
// The amplitude may be represented in decibel (dB=kTRUE) or linear (dB=kFALSE).
// In case of filter kernel data (kernel=kTRUE), the histogram will be rescaled
// such that the maximum value is at 1 for linear amplitudes or 0 for amplitudes in dB.
// The optional argument "hist" may be used to obtain a (zero padded) histogram of the time domain kernel.
// However, in case kernel=kFALSE, the "hist" histogram will be left empty.

Definition at line 5420 of file NcDSP.cxx.

◆ HistogramTrafoResult()

void NcDSP::HistogramTrafoResult ( TString name,
Int_t mode,
TH1 * hist,
TString sel )
protected
// Internal memberfunction to provide a histogram with the requested transformation result.
// The histogram contents are stored in the fHisto array, which can be accessed via GetData().
//
// Input arguments :
// -----------------
// mode : >0 --> Result of a forward transformation
// <0 --> Result of a backward transformation
// hist : Histogram with selected results
// sel : String to specify the representation of the result histogram
// "RE" --> Y-axis shows the values of the real (re) components
// "IM" --> Y-axis shows the values of the imaginary (im) components
// "AMP" --> Y-axis shows the values of the amplitudes, i.e. sqrt(re*re+im*im)
// "dB" --> Y-axis shows the values of the amplitudes, i.e. sqrt(re*re+im*im), in decibel
// "PHIR" --> Y-axis shows the values of the phases, i.e. arctan(im/re), in radians
// "PHID" --> Y-axis shows the values of the phases, i.e. arctan(im/re), in degrees
// "CPHI" --> Y axis shows the values of cos(phase), i.e. cos(arctan(im/re))
// "SPHI" --> Y axis shows the values of cos(phase), i.e. cos(arctan(im/re))
// "OMEGA" --> Y axis shows the values of the angular frequency in rad/sec, i.e. d(arctan(im/re))/dt
// "k" --> X-axis represents the index k in the frequency domain
// "f" --> X-axis represents the fraction f of the sampling rate in the frequency domain
// "Hz" --> X-axis represents the actual frequency in Hz in the frequency domain
// "n" --> X-axis represents the index n in the time domain
// "t" --> X-axis represents the actual time in seconds in the time domain
//
// Note : The options "Hz", "t" and "OMEGA" can only be used if the actual data acquisition sampling rate
// has been provided via the Load() memberfunction.
//
// Examples :
// ----------
// sel="AMP t" will show the amplitudes as a function of the actual time in seconds.
// sel="PHID n" will show the phases in degrees as a function of the index n in the time domain.

Definition at line 1768 of file NcDSP.cxx.

◆ Load() [1/6]

void NcDSP::Load ( Int_t n,
Double_t * re,
Double_t * im = 0,
Float_t f = -1 )
// Provide new input data to be processed and reset the FFTW processor.
//
// Note : The (optional) waveform data stored via SetWaveform() will not be modified.
//
// Input arguments :
// -----------------
// n : The number of data elements
// re : Array with real data elements
// im : Array with imaginary data elements
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// If provided, the array size must be at least of dimension "n".
//
// The defaults are im=0 and f=-1.

Definition at line 395 of file NcDSP.cxx.

◆ Load() [2/6]

void NcDSP::Load ( NcSample * s,
Int_t i,
Float_t f = -1 )
// Provide new input data from the NcSample "s" to be processed and reset the FFTW processor.
//
// Note : The (optional) waveform data stored via SetWaveform() will not be modified.
//
// Input arguments :
// -----------------
// i : The data of the i-th variable (1=first) of the NcSample are used
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// Note : The Store Mode of the NcSample must be activated.
//
// The default value is f=-1.

Definition at line 493 of file NcDSP.cxx.

◆ Load() [3/6]

void NcDSP::Load ( NcSample * s,
TString name,
Float_t f = -1 )
// Provide new input data from the NcSample "s" to be processed and reset the FFTW processor.
//
// Note : The (optional) waveform data stored via SetWaveform() will not be modified.
//
// Input arguments :
// -----------------
// name : Name of the NcSample variable of which the data are used
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// Note : The Store Mode of the NcSample must be activated.
//
// The default value is f=-1.

Definition at line 545 of file NcDSP.cxx.

◆ Load() [4/6]

void NcDSP::Load ( TArray * re,
TArray * im = 0,
Float_t f = -1 )
// Provide new input data to be processed and reset the FFTW processor.
//
// Note : The (optional) waveform data stored via SetWaveform() will not be modified.
//
// Input arguments :
// -----------------
// re : Array with real data elements
// im : Array with imaginary data elements
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// If both arrays are provided, the size of the smallest non-empty array
// will be used to compose the input data.
//
// The defaults are im=0 and f=-1.

Definition at line 434 of file NcDSP.cxx.

◆ Load() [5/6]

void NcDSP::Load ( TGraph * gr,
Float_t f = -1 )
// Provide new input data from a TGraph object to be processed and reset the FFTW processor.
//
// Note : The (optional) waveform data stored via SetWaveform() will not be modified.
//
// Input arguments :
// -----------------
// gr : Pointer to the TGraph object of which the bin contents are used
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// The default value is f=-1.

Definition at line 640 of file NcDSP.cxx.

◆ Load() [6/6]

void NcDSP::Load ( TH1 * h,
Float_t f = -1 )
// Provide new input data from a 1-Dimensional histogram to be processed and reset the FFTW processor.
//
// Note : The (optional) waveform data stored via SetWaveform() will not be modified.
//
// Input arguments :
// -----------------
// h : Pointer to the 1-D histogram of which the bin contents are used
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// The default value is f=-1.

Definition at line 591 of file NcDSP.cxx.

◆ LoadResult()

void NcDSP::LoadResult ( )
// Load the current transformation result as new input data in order to enable
// inverse transformations acting on previous transformation results.
//
// Notes :
// -------
// 1) Invokation of this memberfunction will reset all output arrays, so that
// previously obtained results will internally be lost.
// Use the GetData() memberfunction to retrieve data that might be needed later again.
// 2) Invokation of one of the other Load() memberfunctions will also reset the input arrays.
// Use the GetData() memberfunction to retrieve data that might be needed later again.
// 3) The (optional) waveform data stored via SetWaveform() will not be modified.

Definition at line 691 of file NcDSP.cxx.

◆ Periodogram() [1/5]

TGraph NcDSP::Periodogram ( TString tu,
Double_t Tmin,
Double_t Tmax,
Int_t n,
NcSample s,
Int_t it,
Int_t iy,
Int_t idy = 0,
TF1 * Z = 0,
NcDevice * results = 0 ) const
// Provide via the returned TGraph a periodogram of the specified data as stored in the NcSample object "s".
// Note that the storage mode for the NcSample object has to be activated.
// If the index "idy" is provided, the corresponding values will be regarded as uncertainties
// and will be used as weights in the fitting procedure (see below).
//
// This memberfunction provides a tool to search for periodicity in time series,
// even in case of unequally spaced data.
// The implementation consists of an analytic algorithm known as the
// Generalised Lomb-Scargle (GLS) periodogram and is equivalent to fitting the function
//
// y(t)=a*Z(t)*cos(2*pi*f*t)+b*Z(t)*sin(2*pi*f*t)+c
//
// where :
// a,b = Amplitudes
// Z(t) = Amplitude modulation
// c = Offset
// f = Frequency
//
// By performing the fit at various frequencies, the corresponding amplitudes can be used
// to provide for each frequency the corresponding power P.
// It is the distribution P(f) that is reflected in the periodogram, where for the implementation
// here the power is normalised to the interval [0,1].
//
// Details of the used algorithm can be found at :
//
// M. Zechmeister and M. Kuerster, Astronomy & Astrophysics 496 (2009) 577.
//
// Input arguments :
// -----------------
// tu : Units of the time recording.
// Supported are : "ps", "ns", "sec", "hour", "day".
// Tmin : Lower bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency upper bound of fmax=1/Tmin.
// Tmax : Upper bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency lower bound of fmin=1/Tmax.
// Specification of Tmax<=Tmin will result in fmin=0.
// n : Number of steps to be used for the frequency scan.
// This will result in n+1 frequency samplings to cover the interval [fmin,fmax].
// s : The NcSample object containing the recorded data.
// it : Index (1=first) of the NcSample variable with the times of the data recordings.
// iy : Index (1=first) of the NcSample variable with the observed values y(t).
// idy : Optional index (1=first) of the NcSample varaible with the uncertainties of the y-values.
// idy=0 will result in unweighted treatment of the data.
// Z : Optional function to describe the amplitude modulation Z(t).
// Z=0 is equivalent to Z(t)=1, i.e. unmodulated amplitudes.
// results : Optional NcDevice to contain the resulting (f,a,b,c,P) values.
// Each set (f,a,b,c,P) is stored as a "hit" in the NcDevice, ordered w.r.t. increasing f.
// Ordering w.r.t. any value can be obtained via the SortHits() facility of NcDevice.
//
// Notes :
// -------
// 1) Each entry of the provided NcSample object represents a single data recording of the time series.
// 2) In case the index idy>0 is provided but a stored "y uncertainty" is 0, the corresponding NcSample entry will be neglected.
// 3) The settings idxdy=0 and Z=0 basically result in the original Lomb-Scargle periodogram.
//
// The default values are : idxdy=0, Z=0 and results=0.

Definition at line 5810 of file NcDSP.cxx.

◆ Periodogram() [2/5]

TGraph NcDSP::Periodogram ( TString tu,
Double_t Tmin,
Double_t Tmax,
Int_t n,
NcSample s,
TString namet,
TString namey,
TString namedy = "-",
TF1 * Z = 0,
NcDevice * results = 0 ) const
// Provide via the returned TGraph a periodogram of the specified data as stored in the NcSample object "s".
// Note that the storage mode for the NcSample object has to be activated.
// If a valid "namedy" is provided, the corresponding values will be regarded as uncertainties
// and will be used as weights in the fitting procedure (see below).
//
// This memberfunction provides a tool to search for periodicity in time series,
// even in case of unequally spaced data.
// The implementation consists of an analytic algorithm known as the
// Generalised Lomb-Scargle (GLS) periodogram and is equivalent to fitting the function
//
// y(t)=a*Z(t)*cos(2*pi*f*t)+b*Z(t)*sin(2*pi*f*t)+c
//
// where :
// a,b = Amplitudes
// Z(t) = Amplitude modulation
// c = Offset
// f = Frequency
//
// By performing the fit at various frequencies, the corresponding amplitudes can be used
// to provide for each frequency the corresponding power P.
// It is the distribution P(f) that is reflected in the periodogram, where for the implementation
// here the power is normalised to the interval [0,1].
//
// Details of the used algorithm can be found at :
//
// M. Zechmeister and M. Kuerster, Astronomy & Astrophysics 496 (2009) 577.
//
// Input arguments :
// -----------------
// tu : Units of the time recording.
// Supported are : "ps", "ns", "sec", "hour", "day".
// Tmin : Lower bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency upper bound of fmax=1/Tmin.
// Tmax : Upper bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency lower bound of fmin=1/Tmax.
// Specification of Tmax<=Tmin will result in fmin=0.
// n : Number of steps to be used for the frequency scan.
// This will result in n+1 frequency samplings to cover the interval [fmin,fmax].
// s : The NcSample object containing the recorded data.
// namet : Name of the NcSample variable with the times of the data recordings.
// namey : Name of the NcSample variable with the observed values y(t).
// namedy : Name of the optional NcSample variable with the uncertainties of the y-values.
// namedy="-" will result in unweighted treatment of the data.
// Z : Optional function to describe the amplitude modulation Z(t).
// Z=0 is equivalent to Z(t)=1, i.e. unmodulated amplitudes.
// results : Optional NcDevice to contain the resulting (f,a,b,c,P) values.
// Each set (f,a,b,c,P) is stored as a "hit" in the NcDevice, ordered w.r.t. increasing f.
// Ordering w.r.t. any value can be obtained via the SortHits() facility of NcDevice.
//
// Notes :
// -------
// 1) Each entry of the provided NcSample object represents a single data recording of the time series.
// 2) In case a valid namedy is provided but a stored "y uncertainty" is 0, the corresponding NcSample entry will be neglected.
// 3) The settings namedy="-" and Z=0 basically result in the original Lomb-Scargle periodogram.
//
// The default values are : namedy="-", Z=0 and results=0.

Definition at line 5913 of file NcDSP.cxx.

◆ Periodogram() [3/5]

TGraph NcDSP::Periodogram ( TString tu,
Double_t Tmin,
Double_t Tmax,
Int_t n,
TArray & t,
TArray & y,
TArray * dy = 0,
TF1 * Z = 0,
NcDevice * results = 0 ) const
// Provide via the returned TGraph a periodogram of the specified data.
// If the array "dy" is provided, these uncertainties will be used as weights
// in the fitting procedure (see below).
//
// This memberfunction provides a tool to search for periodicity in time series,
// even in case of unequally spaced data.
// The implementation consists of an analytic algorithm known as the
// Generalised Lomb-Scargle (GLS) periodogram and is equivalent to fitting the function
//
// y(t)=a*Z(t)*cos(2*pi*f*t)+b*Z(t)*sin(2*pi*f*t)+c
//
// where :
// a,b = Amplitudes
// Z(t) = Amplitude modulation
// c = Offset
// f = Frequency
//
// By performing the fit at various frequencies, the corresponding amplitudes can be used
// to provide for each frequency the corresponding power P.
// It is the distribution P(f) that is reflected in the periodogram, where for the implementation
// here the power is normalised to the interval [0,1].
//
// Details of the used algorithm can be found at :
//
// M. Zechmeister and M. Kuerster, Astronomy & Astrophysics 496 (2009) 577.
//
// Input arguments :
// -----------------
// tu : Units of the time recording.
// Supported are : "ps", "ns", "sec", "hour", "day".
// Tmin : Lower bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency upper bound of fmax=1/Tmin.
// Tmax : Upper bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency lower bound of fmin=1/Tmax.
// Specification of Tmax<=Tmin will result in fmin=0.
// n : Number of steps to be used for the frequency scan.
// This will result in n+1 frequency samplings to cover the interval [fmin,fmax].
// t : Array with the times of the data recordings.
// y : Array with the observed values y(t).
// dy : Optional array with the uncertainties of the y-values.
// dy=0 will result in unweighted treatment of the data.
// Z : Optional function to describe the amplitude modulation Z(t).
// Z=0 is equivalent to Z(t)=1, i.e. unmodulated amplitudes.
// results : Optional NcDevice to contain the resulting (f,a,b,c,P) values.
// Each set (f,a,b,c,P) is stored as a "hit" in the NcDevice, ordered w.r.t. increasing f.
// Ordering w.r.t. any value can be obtained via the SortHits() facility of NcDevice.
//
// Notes :
// -------
// 1) All provided input arrays should have the same size.
// 2) The values stored in the arrays t, y and (optionally) dy should match
// such that (t[i],y[i],dy[i]) represents a single data recording of the time series.
// 3) In case the array dy is provided and some dy[i]=0, the corresponding data recording will be neglected.
// 4) The settings dy=0 and Z=0 basically result in the original Lomb-Scargle periodogram.
//
// The default values are : dy=0, Z=0 and results=0.

Definition at line 5552 of file NcDSP.cxx.

◆ Periodogram() [4/5]

TGraph NcDSP::Periodogram ( TString tu,
Double_t Tmin,
Double_t Tmax,
Int_t n,
TGraph & g,
Bool_t err = kFALSE,
TF1 * Z = 0,
NcDevice * results = 0 ) const
// Provide via the returned TGraph a periodogram of the data contained in the graph "g".
// If err=kTRUE, the uncertainties of the graph values will be used as weights in the fitting procedure (see below).
// Note: Treatment of the uncertainties can only be performed if "g" is a TGraphErrors object. Not for a TGraph.
//
// This memberfunction provides a tool to search for periodicity in time series,
// even in case of unequally spaced data.
// The implementation consists of an analytic algorithm known as the
// Generalised Lomb-Scargle (GLS) periodogram and is equivalent to fitting the function
//
// y(t)=a*Z(t)*cos(2*pi*f*t)+b*Z(t)*sin(2*pi*f*t)+c
//
// where :
// a,b = Amplitudes
// Z(t) = Amplitude modulation
// c = Offset
// f = Frequency
//
// By performing the fit at various frequencies, the corresponding amplitudes can be used
// to provide for each frequency the corresponding power P.
// It is the distribution P(f) that is reflected in the periodogram, where for the implementation
// here the power is normalised to the interval [0,1].
//
// Details of the used algorithm can be found at :
//
// M. Zechmeister and M. Kuerster, Astronomy & Astrophysics 496 (2009) 577.
//
// Input arguments :
// -----------------
// tu : Units of the time recording.
// Supported are : "ps", "ns", "sec", "hour", "day".
// Tmin : Lower bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency upper bound of fmax=1/Tmin.
// Tmax : Upper bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency lower bound of fmin=1/Tmax.
// Specification of Tmax<=Tmin will result in fmin=0.
// n : Number of steps to be used for the frequency scan.
// This will result in n+1 frequency samplings to cover the interval [fmin,fmax].
// g : The graph containing the recorded data (t,y).
// err : Flag to trigger the use of the uncertainties of the y-values (kTRUE) or not (kFALSE).
// err=kFALSE will result in unweighted treatment of the data.
// Z : Optional function to describe the amplitude modulation Z(t).
// Z=0 is equivalent to Z(t)=1, i.e. unmodulated amplitudes.
// results : Optional NcDevice to contain the resulting (f,a,b,c,P) values.
// Each set (f,a,b,c,P) is stored as a "hit" in the NcDevice, ordered w.r.t. increasing f.
// Ordering w.r.t. any value can be obtained via the SortHits() facility of NcDevice.
//
// Notes :
// -------
// 1) Each point of the provided graph represents a single data recording of the time series.
// 2) In case err=kTRUE but a stored "y uncertainty" is 0, the corresponding point will be neglected.
// 3) The settings err=kFALSE and Z=0 basically result in the original Lomb-Scargle periodogram.
//
// The default values are : err=kFALSE, Z=0 and results=0.

Definition at line 6092 of file NcDSP.cxx.

◆ Periodogram() [5/5]

TGraph NcDSP::Periodogram ( TString tu,
Double_t Tmin,
Double_t Tmax,
Int_t n,
TH1 & h,
Bool_t err = kFALSE,
TF1 * Z = 0,
NcDevice * results = 0 ) const
// Provide via the returned TGraph a periodogram of the data contained in the histogram "h".
// If err=kTRUE, the uncertainties of the histogram values will be used as weights in the fitting procedure (see below).
//
// This memberfunction provides a tool to search for periodicity in time series,
// even in case of unequally spaced data.
// The implementation consists of an analytic algorithm known as the
// Generalised Lomb-Scargle (GLS) periodogram and is equivalent to fitting the function
//
// y(t)=a*Z(t)*cos(2*pi*f*t)+b*Z(t)*sin(2*pi*f*t)+c
//
// where :
// a,b = Amplitudes
// Z(t) = Amplitude modulation
// c = Offset
// f = Frequency
//
// By performing the fit at various frequencies, the corresponding amplitudes can be used
// to provide for each frequency the corresponding power P.
// It is the distribution P(f) that is reflected in the periodogram, where for the implementation
// here the power is normalised to the interval [0,1].
//
// Details of the used algorithm can be found at :
//
// M. Zechmeister and M. Kuerster, Astronomy & Astrophysics 496 (2009) 577.
//
// Input arguments :
// -----------------
// tu : Units of the time recording.
// Supported are : "ps", "ns", "sec", "hour", "day".
// Tmin : Lower bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency upper bound of fmax=1/Tmin.
// Tmax : Upper bound of the period (in units of "tu") to be investigated.
// This will lead to a frequency lower bound of fmin=1/Tmax.
// Specification of Tmax<=Tmin will result in fmin=0.
// n : Number of steps to be used for the frequency scan.
// This will result in n+1 frequency samplings to cover the interval [fmin,fmax].
// h : The histogram containing the recorded data (t,y).
// err : Flag to trigger the use of the uncertainties of the y-values (kTRUE) or not (kFALSE).
// err=kFALSE will result in unweighted treatment of the data.
// Z : Optional function to describe the amplitude modulation Z(t).
// Z=0 is equivalent to Z(t)=1, i.e. unmodulated amplitudes.
// results : Optional NcDevice to contain the resulting (f,a,b,c,P) values.
// Each set (f,a,b,c,P) is stored as a "hit" in the NcDevice, ordered w.r.t. increasing f.
// Ordering w.r.t. any value can be obtained via the SortHits() facility of NcDevice.
//
// Notes :
// -------
// 1) Each bin of the provided histogram represents a single data recording of the time series.
// 2) In case err=kTRUE but a stored "y uncertainty" is 0, the corresponding bin will be neglected.
// 3) The settings err=kFALSE and Z=0 basically result in the original Lomb-Scargle periodogram.
//
// The default values are : err=kFALSE, Z=0 and results=0.

Definition at line 5999 of file NcDSP.cxx.

◆ Reset()

void NcDSP::Reset ( )
protected
// Internal member function to reset all data and the FFTW processor.

Definition at line 344 of file NcDSP.cxx.

◆ SampleAndHold() [1/2]

TArrayD NcDSP::SampleAndHold ( Int_t ns,
TH1 * hist = 0,
Int_t loc = -1,
Int_t jmin = 0,
Int_t jmax = -1 ) const
// Perform a Sample-And-Hold operation on the data contained in the stored waveform
// over the sampled interval [jmin,jmax], using "n" original samples as the new sampling step size.
// By convention, the first sample is at j=0.
// The result is returned in a TArrayD object and (optionally) in the histogram "hist",
// without modification of the original waveform data.
//
// If the waveform can be regarded as a pulse generator in time, this mimics a
// sample and hold device with a lock time of "n" time units,
// or equivalently a sampling frequency which is reduced by a factor 1/n.
//
// The input argument "loc" determines whether the resulting data
// will be recorded at the start (loc<0), center (loc=0) or end (loc>0)
// of the new sampling step size.
// However, in case the recording location of the last sampling step would
// exceed "jmax", the data will be recorded at the value of "jmax".
//
// If jmax<=jmin the full data array of the stored waveform will be used.
//
// The default values are hist=0, loc=-1, jmin=0 and jmax=-1.

Definition at line 3430 of file NcDSP.cxx.

◆ SampleAndHold() [2/2]

TArrayD NcDSP::SampleAndHold ( TF1 f,
Double_t step,
Double_t vmin,
Double_t vmax,
TH1 * hist = 0,
Int_t loc = -1 ) const
// Perform a Sample-And-Hold operation on the specified function "f"
// in the interval [vmin,vmax], using "step" as the sampling step size.
// The result is returned in a TArrayD object and (optionally) in the histogram "hist".
//
// If "f" can be regarded as a pulse generator in time, this mimics a
// sample and hold device with a lock time of "step" time units,
// or equivalently a sampling frequency of 1/step.
//
// The input argument "loc" determines whether the resulting data
// will be recorded at the start (loc<0), center (loc=0) or end (loc>0)
// of the sampling step size.
// However, in case the recording location of the last sampling step would
// exceed "vmax", the data will be recorded at the value of "vmax".
//
// The default values are hist=0 and loc=-1.

Definition at line 3357 of file NcDSP.cxx.

◆ SampleAndSum() [1/2]

TArrayD NcDSP::SampleAndSum ( Int_t ns,
TH1 * hist,
Int_t jmin = 0,
Int_t jmax = -1 ) const
// Perform a Sample-And-Sum operation on the data contained in the stored waveform
// over the sampled interval [jmin,jmax], using "n" original samples as the new sampling step size.
// By convention, the first sample is at j=0.
// The result is returned in a TArrayD object and (optionally) in the histogram "hist",
// without modification of the original waveform data.
//
// This procedure resembles a Sample-And-Hold operation, but instead of locking
// the data recording during the stepsize of "n" samplings, the data that appear
// within "n" samplings are summed.
//
// If the waveform data can be regarded as sampling in time, this mimics a
// Switched Capacitor Array (SCA) with a time gate of "n" time units,
// or equivalently a sampling frequency which is reduced by a factor 1/n.
//
// If jmax<=jmin the full data array of the stored waveform will be used.
//
// The default values are hist=0, jmin=0 and jmax=-1.

Definition at line 3604 of file NcDSP.cxx.

◆ SampleAndSum() [2/2]

TArrayD NcDSP::SampleAndSum ( TF1 f,
Double_t step,
Double_t vmin,
Double_t vmax,
TH1 * hist = 0 ) const
// Perform a Sample-And-Sum operation on the specified function "f"
// in the interval [vmin,vmax], using "step" as the sampling step size.
// The result is returned in a TArrayD object and (optionally) in the histogram "hist".
//
// This procedure resembles a Sample-And-Hold operation, but instead of locking
// the data recording during the stepsize, the data that appear within "step"
// are summed c.q. integrated.
//
// If "f" can be regarded as sampling in time, this mimics a
// Switched Capacitor Array (SCA) with a time gate of "step" time units,
// or equivalently a sampling frequency which is reduced by a factor 1/step.
//
// The default value is hist=0.

Definition at line 3534 of file NcDSP.cxx.

◆ SetSamplingFrequency()

void NcDSP::SetSamplingFrequency ( Float_t f)
// Set the actual DAQ sampling frequency in Hz.
// This sampling frequency may be overwritten by invokation of
// one of the various Load() member functions.
// The sampling frequency is set to 0 in the constructor of this class.

Definition at line 370 of file NcDSP.cxx.

◆ SetWaveform() [1/6]

void NcDSP::SetWaveform ( Int_t n,
Double_t * h,
Float_t f = -1 )
// Set the (system response) waveform for Convolution, Correlation etc.
//
// Note : The input data stored via Load() will not be modified.
//
// Input arguments :
// -----------------
// n : The number of data elements
// h : Array with the waveform data
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// The array size of "h" must be at least of dimension "n".
//
// The default value is f=-1.

Definition at line 716 of file NcDSP.cxx.

◆ SetWaveform() [2/6]

void NcDSP::SetWaveform ( NcSample * s,
Int_t i,
Float_t f = -1 )
// Set the (system response) waveform for Convolution, Correlation etc.
//
// Note : The input data stored via Load() will not be modified.
//
// Input arguments :
// -----------------
// s : NcSample with the waveform data
// i : The data of the i-th variable (1=first) of the NcSample are used
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// Note : The Store Mode of the NcSample must be activated.
//
// The default value is f=-1.

Definition at line 804 of file NcDSP.cxx.

◆ SetWaveform() [3/6]

void NcDSP::SetWaveform ( NcSample * s,
TString name,
Float_t f = -1 )
// Set the (system response) waveform for Convolution, Correlation etc.
//
// Note : The input data stored via Load() will not be modified.
//
// Input arguments :
// -----------------
// s : NcSample with the waveform data
// name : Name of the NcSample variable of which the data are used
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// Note : The Store Mode of the NcSample must be activated.
//
// The default value is f=-1.

Definition at line 858 of file NcDSP.cxx.

◆ SetWaveform() [4/6]

void NcDSP::SetWaveform ( TArray * h,
Float_t f = -1 )
// Set the (system response) waveform for Convolution, Correlation etc.
//
// Note : The input data stored via Load() will not be modified.
//
// Input arguments :
// -----------------
// h : Array with the waveform data
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// The default value is f=-1.

Definition at line 759 of file NcDSP.cxx.

◆ SetWaveform() [5/6]

void NcDSP::SetWaveform ( TGraph * gr,
Float_t f = -1 )
// Set the (system response) waveform for Convolution, Correlation etc.
//
// Note : The input data stored via Load() will not be modified.
//
// Input arguments :
// -----------------
// gr : Pointer to the TGraph object of which the bin contents are used for the waveform data
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// The default value is f=-1.

Definition at line 956 of file NcDSP.cxx.

◆ SetWaveform() [6/6]

void NcDSP::SetWaveform ( TH1 * h,
Float_t f = -1 )
// Set the (system response) waveform for Convolution, Correlation etc.
//
// Note : The input data stored via Load() will not be modified.
//
// Input arguments :
// -----------------
// h : Pointer to the 1-D histogram of which the bin contents are used for the waveform data
// f : (optional) Actual data acquisition sampling frequency in Hz
// In case f<0 the current sampling frequency is not modified.
//
// The default value is f=-1.

Definition at line 906 of file NcDSP.cxx.

◆ Sine()

void NcDSP::Sine ( Int_t type,
TH1 * hist = 0,
TString sel = "none" )
// Perform a normalized 1-dimensional Discrete Sine Transformation (DST).
// Actually, this is just a regular Discrete Fourier Transformation (DFT)
// with only real input values which contain an odd symmetry.
// Consequently, the resulting transformed array is also only real with odd symmetry.
//
// Conventions :
// -------------
// N = The number of data elements
// Time domain array .... : X[]=X[0],...,X[N-1]
// Frequency domain array : Q[]=Q[0],...,Q[N-1]
//
// Sine transform type 1 : Q[k]=(1/sqrt(2*(N+1)))*[2*sum(n=0,n=N-1){X[n]*sin(pi*(n+1)*(k+1)/(N+1))}]
//
// Sine transform type 2 : Q[k]=(1/sqrt(2N))*2*sum(n=0,n=N-1){X[n]*sin(pi*(n+0.5)*(k+1)/N)}
//
// Sine transform type 3 : Q[k]=(1/sqrt(2N))*[pow(-1,k)*X[N-1]+2*sum(n=0,n=N-2){X[n]*sin(pi*(n+1)*(k+0.5)/N)}]
//
// Sine transform type 4 : Q[k]=(1/sqrt(2N))*2*sum(n=0,n=N-1){X[n]*sin(pi*(n+0.5)*(k+0.5)/N)}
//
// Notes :
// -------
// 1) The type 1 transformation is its own inverse.
// 2) The type 4 transformation is its own inverse.
// 3) The type 3 transformation is the inverse of type 2 (and vice versa).
//
// Input arguments :
// -----------------
// type : The type of transformation (i.e. 1,2,3 or 4) to be performed
// The inverse transformations are specified with the corresponding negative type value
// hist : (optional) Histogram with selected results
// sel : String to specify the representation of the result histogram
// "k" --> X-axis represents the index k in the frequency domain
// "f" --> X-axis represents the fraction f of the sampling rate in the frequency domain
// "Hz" --> X-axis represents the actual frequency in Hz in the frequency domain
// "n" --> X-axis represents the index n in the time domain
// "t" --> X-axis represents the actual time in seconds in the time domain
// "2" --> X-axis spans the full number of data points, instead of the usual (N/2)+1
//
// Note : The options "Hz" and "t" can only be used if the actual data acquisition sampling rate
// has been provided via the Load() memberfunction.
//
// Examples :
// ----------
// sel="f" will show the (N/2)+1 amplitudes as a function of the fractional sampling rate.
// sel="k 2" will show all N amplitudes as a function of the index k in the frequency domain.
//
// The default values are hist=0 and sel="none"

Definition at line 1452 of file NcDSP.cxx.

◆ Transmit()

TArrayD NcDSP::Transmit ( Int_t nbits,
Double_t range,
Double_t Vbias = 0,
TArray * Vsig = 0,
TArray * peds = 0,
TH1 * hist = 0,
Int_t B = 0,
Int_t C = 3 ) const
// Mimic signal transmission according to an "nbits" ADC-DAC chain.
// Analog input signals are digitized via the discrete quantization levels of an "nbits" ADC,
// based on the "range" for the analog signal and a bias voltage "Vbias" or array "peds" of pedestal values (see below).
// The analog input signals may be provided by the (optional) TArray "Vsig".
// In case the array "Vsig" is not provided, the stored waveform is used to provide the analog input signals.
// After digitization, the digital signals are converted into analog signals via the corresponding "nbits" DAC.
// In this way the effect of digitization on the original input signals can be investigated, and can provide
// a guideline in selecting the most suitable ADC-DAC system for data transmission.
// The resulting analog values are returned in a TArrayD object, without modification of the original waveform data.
//
// Note : Make sure to use the same units for "range", "Vbias" and the analog input signals.
//
// For further details, please refer to the documentation of the memberfunctions ADC() and DAC().
//
// Input arguments :
// -----------------
// nbits : Digital quantization was performed using an nbits ADC.
// range >0 : The full scale voltage range (Vfs) of the analog signal that corresponds to adc=N-1.
// <0 : |range| is the reference voltage (Vref) of the analog signal that corresponds to the hypothetical adc=N.
// Vbias : The bias voltage that was added to the analog input signal before digitization.
// If specified, the resulting analog signals will be corrected for the bias voltage.
// For a linear DAC the correction via Vbias will only be performed if no pedestal array "peds"
// is specified (see below). If the array "peds" is specified, Vbias will be ignored for a linear DAC.
// For a Log DAC the pedestal values will never be used and bias correction may only be obtained via Vbias.
// Vsig : (Optional) array to contain the digital input signals.
// peds : (Optional) array to contain the pedestal values for the individual digital signals.
// The array "peds" must contain (at least) the same number of values as the number of analog input signals.
// If the array "peds" is specified, the pedestals will be subtracted from the corresponding digital signals
// before the conversion to analog signals is performed in case of a linear DAC.
// For a Log DAC the pedestal values will never be used and bias correction may only be obtained via Vbias.
// B >1 : Base for a log ADC (e.g. B=10 emulates a Log_10 ADC).
// =0 : The ADC will be linear
// =1 : The ADC will be a Log_e ADC.
// Note : When B>0 all (Vbias+Vsig) input values should be positive.
// C : Code efficiency factor that was used for a Log ADC.
// Rule of thumb : pow(B,-C) is about the smallest signal/|range| ratio that can be resolved.
// So, for a Log_10 ADC with C=3, the smallest signal that can be resolved is about |range|/1000.
// Note : It is required that C>0.
//
// Optional output arguments :
// ---------------------------
// hist : Histogram with the analog result
//
// Note :
// ------
// In case no "Vsig" array is provided and no waveform is present, just the ADC and DAC specs will be printed
// and in the TArrayD only the adc value corresponding to Vbias is returned.
// This allows to quickly investigate various scenarios without any data treatment.
//
// The maximum number of bits that is supported is 60 to guarantee identical functioning
// on all machines.
//
// In case of inconsistent input parameters, no processing is performed and an empty TArrayD is returned.
//
// The default values are Vbias=0, Vsig=0, peds=0, hist=0, B=0 and C=3.

Definition at line 3265 of file NcDSP.cxx.

Member Data Documentation

◆ fHisto

TArrayD NcDSP::fHisto
protected

Definition at line 87 of file NcDSP.h.

◆ fImIn

TArrayD NcDSP::fImIn
protected

Definition at line 84 of file NcDSP.h.

◆ fImOut

TArrayD NcDSP::fImOut
protected

Definition at line 86 of file NcDSP.h.

◆ fKeepOutput

Bool_t NcDSP::fKeepOutput
protected

Definition at line 91 of file NcDSP.h.

◆ fN

Int_t NcDSP::fN
protected

Definition at line 81 of file NcDSP.h.

◆ fNorm

TString NcDSP::fNorm
protected

Definition at line 89 of file NcDSP.h.

◆ fNwf

Int_t NcDSP::fNwf
protected

Definition at line 82 of file NcDSP.h.

◆ fProc

TVirtualFFT* NcDSP::fProc
protected

Definition at line 80 of file NcDSP.h.

◆ fReIn

TArrayD NcDSP::fReIn
protected

Definition at line 83 of file NcDSP.h.

◆ fReOut

TArrayD NcDSP::fReOut
protected

Definition at line 85 of file NcDSP.h.

◆ fSample

Float_t NcDSP::fSample
protected

Definition at line 90 of file NcDSP.h.

◆ fWaveform

TArrayD NcDSP::fWaveform
protected

Definition at line 88 of file NcDSP.h.


The documentation for this class was generated from the following files: