Logo ROOT  
Reference Guide
TSpectrum.h
Go to the documentation of this file.
1// @(#)root/spectrum:$Id$
2// Author: Miroslav Morhac 27/05/99
3
4/*************************************************************************
5 * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11#ifndef ROOT_TSpectrum
12#define ROOT_TSpectrum
13
14#include "TNamed.h"
15
16class TH1;
17
18class TSpectrum : public TNamed {
19private:
20
21 TSpectrum(const TSpectrum&); // Not implemented
22 TSpectrum& operator=(const TSpectrum&); // Not implemented
23
24protected:
25 Int_t fMaxPeaks; ///< Maximum number of peaks to be found
26 Int_t fNPeaks; ///< number of peaks found
27 Double_t *fPosition; ///< [fNPeaks] array of current peak positions
28 Double_t *fPositionX; ///< [fNPeaks] X position of peaks
29 Double_t *fPositionY; ///< [fNPeaks] Y position of peaks
30 Double_t fResolution; ///< *NOT USED* resolution of the neighboring peaks
31 TH1 *fHistogram; ///< resulting histogram
32static Int_t fgAverageWindow; ///< Average window of searched peaks
33static Int_t fgIterations; ///< Maximum number of decon iterations (default=3)
34
35public:
36 enum {
50 };
51
52 TSpectrum();
53 TSpectrum(Int_t maxpositions, Double_t resolution=1); // resolution is *NOT USED*
54 virtual ~TSpectrum();
55 virtual TH1 *Background(const TH1 *hist,Int_t niter=20, Option_t *option="");
56 TH1 *GetHistogram() const {return fHistogram;}
57 Int_t GetNPeaks() const {return fNPeaks;}
58 Double_t *GetPositionX() const {return fPositionX;}
59 Double_t *GetPositionY() const {return fPositionY;}
60 virtual void Print(Option_t *option="") const;
61 virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05);
62 static void SetAverageWindow(Int_t w=3); //set average window
63 static void SetDeconIterations(Int_t n=3); //set max number of decon iterations
64 void SetResolution(Double_t resolution=1); // *NOT USED*
65
66 //new functions January 2006
67 const char *Background(Double_t *spectrum, Int_t ssize,Int_t numberIterations,Int_t direction, Int_t filterOrder,bool smoothing,Int_t smoothWindow,bool compton);
68 const char *SmoothMarkov(Double_t *source, Int_t ssize, Int_t averWindow);
69 const char *Deconvolution(Double_t *source, const Double_t *response,Int_t ssize, Int_t numberIterations,Int_t numberRepetitions, Double_t boost );
70 const char *DeconvolutionRL(Double_t *source, const Double_t *response,Int_t ssize, Int_t numberIterations,Int_t numberRepetitions, Double_t boost );
71 const char *Unfolding(Double_t *source,const Double_t **respMatrix,Int_t ssizex, Int_t ssizey,Int_t numberIterations,Int_t numberRepetitions, Double_t boost);
72 Int_t SearchHighRes(Double_t *source,Double_t *destVector, Int_t ssize,Double_t sigma, Double_t threshold,bool backgroundRemove,Int_t deconIterations,bool markov, Int_t averWindow);
73 Int_t Search1HighRes(Double_t *source,Double_t *destVector, Int_t ssize,Double_t sigma, Double_t threshold,bool backgroundRemove,Int_t deconIterations,bool markov, Int_t averWindow);
74
75 static Int_t StaticSearch(const TH1 *hist, Double_t sigma=2, Option_t *option="goff", Double_t threshold=0.05);
76 static TH1 *StaticBackground(const TH1 *hist,Int_t niter=20, Option_t *option="");
77
78 ClassDef(TSpectrum,3) //Peak Finder, background estimator, Deconvolution
79};
80
81#endif
82
int Int_t
Definition: RtypesCore.h:41
double Double_t
Definition: RtypesCore.h:55
const char Option_t
Definition: RtypesCore.h:62
#define ClassDef(name, id)
Definition: Rtypes.h:326
The TH1 histogram class.
Definition: TH1.h:56
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
Advanced Spectra Processing.
Definition: TSpectrum.h:18
Double_t fResolution
NOT USED resolution of the neighboring peaks
Definition: TSpectrum.h:30
static TH1 * StaticBackground(const TH1 *hist, Int_t niter=20, Option_t *option="")
Static function, interface to TSpectrum::Background.
Definition: TSpectrum.cxx:2587
Int_t fMaxPeaks
Maximum number of peaks to be found.
Definition: TSpectrum.h:25
Double_t * fPositionX
[fNPeaks] X position of peaks
Definition: TSpectrum.h:28
const char * SmoothMarkov(Double_t *source, Int_t ssize, Int_t averWindow)
One-dimensional markov spectrum smoothing function.
Definition: TSpectrum.cxx:1195
const char * Unfolding(Double_t *source, const Double_t **respMatrix, Int_t ssizex, Int_t ssizey, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional unfolding function.
Definition: TSpectrum.cxx:1889
Double_t * GetPositionY() const
Definition: TSpectrum.h:59
TH1 * GetHistogram() const
Definition: TSpectrum.h:56
void SetResolution(Double_t resolution=1)
NOT USED resolution: determines resolution of the neighbouring peaks default value is 1 correspond to...
Definition: TSpectrum.cxx:351
virtual ~TSpectrum()
Destructor.
Definition: TSpectrum.cxx:87
@ kBackOrder8
Definition: TSpectrum.h:40
@ kBackSmoothing11
Definition: TSpectrum.h:47
@ kBackSmoothing15
Definition: TSpectrum.h:49
@ kBackOrder2
Definition: TSpectrum.h:37
@ kBackOrder6
Definition: TSpectrum.h:39
@ kBackSmoothing7
Definition: TSpectrum.h:45
@ kBackIncreasingWindow
Definition: TSpectrum.h:41
@ kBackDecreasingWindow
Definition: TSpectrum.h:42
@ kBackSmoothing13
Definition: TSpectrum.h:48
@ kBackSmoothing9
Definition: TSpectrum.h:46
@ kBackOrder4
Definition: TSpectrum.h:38
@ kBackSmoothing5
Definition: TSpectrum.h:44
@ kBackSmoothing3
Definition: TSpectrum.h:43
TSpectrum()
Constructor.
Definition: TSpectrum.cxx:50
virtual TH1 * Background(const TH1 *hist, Int_t niter=20, Option_t *option="")
One-dimensional background estimation function.
Definition: TSpectrum.cxx:152
Int_t SearchHighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
One-dimensional high-resolution peak search function.
Definition: TSpectrum.cxx:2126
virtual Int_t Search(const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
One-dimensional peak search function.
Definition: TSpectrum.cxx:266
Int_t fNPeaks
number of peaks found
Definition: TSpectrum.h:26
Int_t Search1HighRes(Double_t *source, Double_t *destVector, Int_t ssize, Double_t sigma, Double_t threshold, bool backgroundRemove, Int_t deconIterations, bool markov, Int_t averWindow)
Old name of SearcHighRes introduced for back compatibility.
Definition: TSpectrum.cxx:2563
TSpectrum & operator=(const TSpectrum &)
const char * Deconvolution(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
Definition: TSpectrum.cxx:1459
static Int_t fgAverageWindow
Average window of searched peaks.
Definition: TSpectrum.h:32
Int_t GetNPeaks() const
Definition: TSpectrum.h:57
static void SetDeconIterations(Int_t n=3)
Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::Search...
Definition: TSpectrum.cxx:108
Double_t * fPosition
[fNPeaks] array of current peak positions
Definition: TSpectrum.h:27
static void SetAverageWindow(Int_t w=3)
Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes).
Definition: TSpectrum.cxx:99
Double_t * GetPositionX() const
Definition: TSpectrum.h:58
const char * DeconvolutionRL(Double_t *source, const Double_t *response, Int_t ssize, Int_t numberIterations, Int_t numberRepetitions, Double_t boost)
One-dimensional deconvolution function.
Definition: TSpectrum.cxx:1665
static Int_t fgIterations
Maximum number of decon iterations (default=3)
Definition: TSpectrum.h:33
TSpectrum(const TSpectrum &)
TH1 * fHistogram
resulting histogram
Definition: TSpectrum.h:31
virtual void Print(Option_t *option="") const
Print the array of positions.
Definition: TSpectrum.cxx:219
static Int_t StaticSearch(const TH1 *hist, Double_t sigma=2, Option_t *option="goff", Double_t threshold=0.05)
Static function, interface to TSpectrum::Search.
Definition: TSpectrum.cxx:2577
Double_t * fPositionY
[fNPeaks] Y position of peaks
Definition: TSpectrum.h:29
const Double_t sigma
const Int_t n
Definition: legend1.C:16
LVector boost(const LVector &v, const BoostVector &b)
Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost The only re...
Definition: VectorUtil.h:327