Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches

Advanced Spectra Processing.

Author
Miroslav Morhac

Legacy Code

TSpectrum is a legacy interface: there will be no bug fixes nor new developments. Therefore it is not recommended to use it in new long-term production code. But, depending on the context, using TSpectrum might still be a valid solution. For modeling a spectrum fitting and estimating the background one can use RooFit while for deconvolution and unfolding one can use TUnfold.

This class contains advanced spectra processing functions for:

  • One-dimensional background estimation
  • One-dimensional smoothing
  • One-dimensional deconvolution
  • One-dimensional peak search

The algorithms in this class have been published in the following references:

  1. M.Morhac et al.: Background elimination methods for multidimensional coincidence gamma-ray spectra. Nuclear Instruments and Methods in Physics Research A 401 (1997) 113-132.
  2. M.Morhac et al.: Efficient one- and two-dimensional Gold deconvolution and its application to gamma-ray spectra decomposition. Nuclear Instruments and Methods in Physics Research A 401 (1997) 385-408.
  3. M.Morhac et al.: Identification of peaks in multidimensional coincidence gamma-ray spectra. Nuclear Instruments and Methods in Research Physics A 443(2000), 108-125.

Definition at line 18 of file TSpectrum.h.

Public Types

enum  {
  kIsOnHeap = 0x01000000 , kNotDeleted = 0x02000000 , kZombie = 0x04000000 , kInconsistent = 0x08000000 ,
  kBitMask = 0x00ffffff
}
 
enum  { kSingleKey = (1ULL << ( 0 )) , kOverwrite = (1ULL << ( 1 )) , kWriteDelete = (1ULL << ( 2 )) }
 
enum  {
  kBackOrder2 =0 , kBackOrder4 =1 , kBackOrder6 =2 , kBackOrder8 =3 ,
  kBackIncreasingWindow =0 , kBackDecreasingWindow =1 , kBackSmoothing3 =3 , kBackSmoothing5 =5 ,
  kBackSmoothing7 =7 , kBackSmoothing9 =9 , kBackSmoothing11 =11 , kBackSmoothing13 =13 ,
  kBackSmoothing15 =15
}
 
enum  EDeprecatedStatusBits { kObjInCanvas = (1ULL << ( 3 )) }
 
enum  EStatusBits {
  kCanDelete = (1ULL << ( 0 )) , kMustCleanup = (1ULL << ( 3 )) , kIsReferenced = (1ULL << ( 4 )) , kHasUUID = (1ULL << ( 5 )) ,
  kCannotPick = (1ULL << ( 6 )) , kNoContextMenu = (1ULL << ( 8 )) , kInvalidObject = (1ULL << ( 13 ))
}
 

Public Member Functions

 TSpectrum ()
 Constructor.
 
 TSpectrum (Int_t maxpositions, Double_t resolution=1)
 
 ~TSpectrum () override
 Destructor.
 
void AbstractMethod (const char *method) const
 Call this function within a function that you don't want to define as purely virtual, in order not to force all users deriving from that class to implement that maybe (on their side) unused function; but at the same time, emit a run-time warning if they try to call it, telling that it is not implemented in the derived class: action must thus be taken on the user side to override it.
 
virtual void AppendPad (Option_t *option="")
 Append graphics object to current pad.
 
virtual TH1Background (const TH1 *hist, Int_t nIter=20, Option_t *option="")
 One-dimensional background estimation function.
 
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)
 This function calculates background spectrum from source spectrum.
 
virtual void Browse (TBrowser *b)
 Browse object. May be overridden for another default action.
 
ULong_t CheckedHash ()
 Check and record whether this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object.
 
virtual const char * ClassName () const
 Returns name of class to which the object belongs.
 
void Clear (Option_t *option="") override
 Set name and title to empty strings ("").
 
TObjectClone (const char *newname="") const override
 Make a clone of an object using the Streamer facility.
 
Int_t Compare (const TObject *obj) const override
 Compare two TNamed objects.
 
void Copy (TObject &named) const override
 Copy this to obj.
 
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.
 
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.
 
virtual void Delete (Option_t *option="")
 Delete this object.
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 Computes distance from point (px,py) to the object.
 
virtual void Draw (Option_t *option="")
 Default Draw method for all objects.
 
virtual void DrawClass () const
 Draw class inheritance tree of the class to which this object belongs.
 
virtual TObjectDrawClone (Option_t *option="") const
 Draw a clone of this object in the current selected pad with: gROOT->SetSelectedPad(c1).
 
virtual void Dump () const
 Dump contents of object on stdout.
 
virtual void Error (const char *method, const char *msgfmt,...) const
 Issue error message.
 
virtual void Execute (const char *method, const char *params, Int_t *error=nullptr)
 Execute method on this object with the given parameter string, e.g.
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=nullptr)
 Execute method on this object with parameters stored in the TObjArray.
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 Execute action corresponding to an event at (px,py).
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 Issue fatal error message.
 
virtual void FillBuffer (char *&buffer)
 Encode TNamed into output buffer.
 
virtual TObjectFindObject (const char *name) const
 Must be redefined in derived classes.
 
virtual TObjectFindObject (const TObject *obj) const
 Must be redefined in derived classes.
 
virtual Option_tGetDrawOption () const
 Get option used by the graphics system to draw this object.
 
TH1GetHistogram () const
 
virtual const char * GetIconName () const
 Returns mime type name of object.
 
const char * GetName () const override
 Returns name of object.
 
Int_t GetNPeaks () const
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const
 Returns string containing info about the object at position (px,py).
 
virtual Option_tGetOption () const
 
Double_tGetPositionX () const
 
Double_tGetPositionY () const
 
const char * GetTitle () const override
 Returns title of object.
 
virtual UInt_t GetUniqueID () const
 Return the unique object id.
 
virtual Bool_t HandleTimer (TTimer *timer)
 Execute action in response of a timer timing out.
 
ULong_t Hash () const override
 Return hash value for this object.
 
Bool_t HasInconsistentHash () const
 Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e.
 
virtual void Info (const char *method, const char *msgfmt,...) const
 Issue info message.
 
virtual Bool_t InheritsFrom (const char *classname) const
 Returns kTRUE if object inherits from class "classname".
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 Returns kTRUE if object inherits from TClass cl.
 
virtual void Inspect () const
 Dump contents of this object in a graphics canvas.
 
void InvertBit (UInt_t f)
 
TClassIsA () const override
 
Bool_t IsDestructed () const
 IsDestructed.
 
virtual Bool_t IsEqual (const TObject *obj) const
 Default equal comparison (objects are equal if they have the same address in memory).
 
virtual Bool_t IsFolder () const
 Returns kTRUE in case object contains browsable objects (like containers or lists of other objects).
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
Bool_t IsSortable () const override
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
void ls (Option_t *option="") const override
 List TNamed name and title.
 
void MayNotUse (const char *method) const
 Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary).
 
virtual Bool_t Notify ()
 This method must be overridden to handle object notification (the base implementation is no-op).
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 Use this method to declare a method obsolete.
 
void operator delete (void *, size_t)
 Operator delete for sized deallocation.
 
void operator delete (void *ptr)
 Operator delete.
 
void operator delete (void *ptr, void *vp)
 Only called by placement new when throwing an exception.
 
void operator delete[] (void *, size_t)
 Operator delete [] for sized deallocation.
 
void operator delete[] (void *ptr)
 Operator delete [].
 
void operator delete[] (void *ptr, void *vp)
 Only called by placement new[] when throwing an exception.
 
void * operator new (size_t sz)
 
void * operator new (size_t sz, void *vp)
 
void * operator new[] (size_t sz)
 
void * operator new[] (size_t sz, void *vp)
 
virtual void Paint (Option_t *option="")
 This method must be overridden if a class wants to paint itself.
 
virtual void Pop ()
 Pop on object drawn in a pad to the top of the display list.
 
void Print (Option_t *option="") const override
 Print the array of positions.
 
virtual Int_t Read (const char *name)
 Read contents of object with specified name from the current directory.
 
virtual void RecursiveRemove (TObject *obj)
 Recursively remove this object from a list.
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 Save this object in the file specified by filename.
 
virtual void SavePrimitive (std::ostream &out, Option_t *option="")
 Save a primitive as a C++ statement(s) on output stream "out".
 
virtual Int_t Search (const TH1 *hist, Double_t sigma=2, Option_t *option="", Double_t threshold=0.05)
 One-dimensional peak search function.
 
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.
 
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.
 
void SetBit (UInt_t f)
 
void SetBit (UInt_t f, Bool_t set)
 Set or unset the user status bits as specified in f.
 
virtual void SetDrawOption (Option_t *option="")
 Set drawing option for object.
 
virtual void SetName (const char *name)
 Set the name of the TNamed.
 
virtual void SetNameTitle (const char *name, const char *title)
 Set all the TNamed parameters (name and title).
 
void SetResolution (Double_t resolution=1)
 NOT USED resolution: determines resolution of the neighbouring peaks default value is 1 correspond to 3 sigma distance between peaks.
 
virtual void SetTitle (const char *title="")
 Set the title of the TNamed.
 
virtual void SetUniqueID (UInt_t uid)
 Set the unique object id.
 
virtual Int_t Sizeof () const
 Return size of the TNamed part of the TObject.
 
const char * SmoothMarkov (Double_t *source, Int_t ssize, Int_t averWindow)
 One-dimensional markov spectrum smoothing function.
 
void Streamer (TBuffer &) override
 Stream an object of class TObject.
 
void StreamerNVirtual (TBuffer &ClassDef_StreamerNVirtual_b)
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 Issue system error message.
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
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.
 
virtual void UseCurrentStyle ()
 Set current style settings in this object This function is called when either TCanvas::UseCurrentStyle or TROOT::ForceStyle have been invoked.
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 Issue warning message.
 
virtual Int_t Write (const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
 Write this object to the current directory.
 
virtual Int_t Write (const char *name=nullptr, Int_t option=0, Int_t bufsize=0) const
 Write this object to the current directory.
 

Static Public Member Functions

static TClassClass ()
 
static const char * Class_Name ()
 
static constexpr Version_t Class_Version ()
 
static const char * DeclFileName ()
 
static Longptr_t GetDtorOnly ()
 Return destructor only flag.
 
static Bool_t GetObjectStat ()
 Get status of object stat flag.
 
static void SetAverageWindow (Int_t w=3)
 Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes).
 
static void SetDeconIterations (Int_t n=3)
 Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::SearchHighRes).
 
static void SetDtorOnly (void *obj)
 Set destructor only flag.
 
static void SetObjectStat (Bool_t stat)
 Turn on/off tracking of objects in the TObjectTable.
 
static TH1StaticBackground (const TH1 *hist, Int_t niter=20, Option_t *option="")
 Static function, interface to TSpectrum::Background.
 
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.
 

Protected Types

enum  { kOnlyPrepStep = (1ULL << ( 3 )) }
 

Protected Member Functions

virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 Interface to ErrorHandler (protected).
 
void MakeZombie ()
 
void SavePrimitiveNameTitle (std::ostream &out, const char *variable_name)
 Save object name and title into the output stream "out".
 

Static Protected Member Functions

static void SavePrimitiveConstructor (std::ostream &out, TClass *cl, const char *variable_name, const char *constructor_agrs="", Bool_t empty_line=kTRUE)
 Save object constructor in the output stream "out".
 
static void SavePrimitiveDraw (std::ostream &out, const char *variable_name, Option_t *option=nullptr)
 Save invocation of primitive Draw() method Skipped if option contains "nodraw" string.
 
static TString SavePrimitiveVector (std::ostream &out, const char *prefix, Int_t len, Double_t *arr, Int_t flag=0)
 Save array in the output stream "out" as vector.
 

Protected Attributes

TH1fHistogram
 resulting histogram
 
Int_t fMaxPeaks
 Maximum number of peaks to be found.
 
TString fName
 
Int_t fNPeaks
 number of peaks found
 
Double_tfPosition
 [fNPeaks] array of current peak positions
 
Double_tfPositionX
 [fNPeaks] X position of peaks
 
Double_tfPositionY
 [fNPeaks] Y position of peaks
 
Double_t fResolution
 NOT USED resolution of the neighboring peaks
 
TString fTitle
 

Static Protected Attributes

static Int_t fgAverageWindow = 3
 Average window of searched peaks.
 
static Int_t fgIterations = 3
 Maximum number of decon iterations (default=3)
 

Private Member Functions

 TSpectrum (const TSpectrum &)
 
TSpectrumoperator= (const TSpectrum &)
 

Static Private Member Functions

static void AddToTObjectTable (TObject *)
 Private helper function which will dispatch to TObjectTable::AddObj.
 

Private Attributes

UInt_t fBits
 bit field status word
 
UInt_t fUniqueID
 object unique identifier
 

Static Private Attributes

static Longptr_t fgDtorOnly = 0
 object for which to call dtor only (i.e. no delete)
 
static Bool_t fgObjectStat = kTRUE
 if true keep track of objects in TObjectTable
 

#include <TSpectrum.h>

Inheritance diagram for TSpectrum:
TNamed TObject

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
inherited
Enumerator
kIsOnHeap 

object is on heap

kNotDeleted 

object has not been deleted

kZombie 

object ctor failed

kInconsistent 

class overload Hash but does call RecursiveRemove in destructor

kBitMask 

Definition at line 89 of file TObject.h.

◆ anonymous enum

anonymous enum
inherited
Enumerator
kSingleKey 

write collection with single key

kOverwrite 

overwrite existing object with same name

kWriteDelete 

write object, then delete previous key with same name

Definition at line 99 of file TObject.h.

◆ anonymous enum

anonymous enum
protectedinherited
Enumerator
kOnlyPrepStep 

Used to request that the class specific implementation of TObject::Write just prepare the objects to be ready to be written but do not actually write them into the TBuffer.

This is just for example by TBufferMerger to request that the TTree inside the file calls TTree::FlushBaskets (outside of the merging lock) and TBufferMerger will later ask for the write (inside the merging lock). To take advantage of this feature the class needs to overload TObject::Write and use this enum value accordingly. (See TTree::Write and TObject::Write) Do not use, this feature will be migrate to the Merge function (See TClass and TTree::Merge)

Definition at line 106 of file TObject.h.

◆ anonymous enum

Enumerator
kBackOrder2 
kBackOrder4 
kBackOrder6 
kBackOrder8 
kBackIncreasingWindow 
kBackDecreasingWindow 
kBackSmoothing3 
kBackSmoothing5 
kBackSmoothing7 
kBackSmoothing9 
kBackSmoothing11 
kBackSmoothing13 
kBackSmoothing15 

Definition at line 36 of file TSpectrum.h.

◆ EDeprecatedStatusBits

Enumerator
kObjInCanvas 

for backward compatibility only, use kMustCleanup

Definition at line 84 of file TObject.h.

◆ EStatusBits

Enumerator
kCanDelete 

if object in a list can be deleted

kMustCleanup 

if object destructor must call RecursiveRemove()

kIsReferenced 

if object is referenced by a TRef or TRefArray

kHasUUID 

if object has a TUUID (its fUniqueID=UUIDNumber)

kCannotPick 

if object in a pad cannot be picked

kNoContextMenu 

if object does not want context menu

kInvalidObject 

if object ctor succeeded but object should not be used

Definition at line 70 of file TObject.h.

Constructor & Destructor Documentation

◆ TSpectrum() [1/3]

TSpectrum::TSpectrum ( const TSpectrum & )
private

◆ TSpectrum() [2/3]

TSpectrum::TSpectrum ( )

Constructor.

Definition at line 42 of file TSpectrum.cxx.

◆ TSpectrum() [3/3]

TSpectrum::TSpectrum ( Int_t maxpositions,
Double_t resolution = 1 )
  • maxpositions: maximum number of peaks
  • resolution: NOT USED determines resolution of the neighbouring peaks default value is 1 correspond to 3 sigma distance between peaks. Higher values allow higher resolution (smaller distance between peaks. May be set later through SetResolution.

Definition at line 62 of file TSpectrum.cxx.

◆ ~TSpectrum()

TSpectrum::~TSpectrum ( )
override

Destructor.

Definition at line 79 of file TSpectrum.cxx.

Member Function Documentation

◆ AbstractMethod()

void TObject::AbstractMethod ( const char * method) const
inherited

Call this function within a function that you don't want to define as purely virtual, in order not to force all users deriving from that class to implement that maybe (on their side) unused function; but at the same time, emit a run-time warning if they try to call it, telling that it is not implemented in the derived class: action must thus be taken on the user side to override it.

In other word, this method acts as a "runtime purely virtual" warning instead of a "compiler purely virtual" error.

Warning
This interface is a legacy function that is no longer recommended to be used by new development code.
Note
The name "AbstractMethod" does not imply that it's an abstract method in the strict C++ sense.

Definition at line 1149 of file TObject.cxx.

◆ AddToTObjectTable()

void TObject::AddToTObjectTable ( TObject * op)
staticprivateinherited

Private helper function which will dispatch to TObjectTable::AddObj.

Included here to avoid circular dependency between header files.

Definition at line 195 of file TObject.cxx.

◆ AppendPad()

void TObject::AppendPad ( Option_t * option = "")
virtualinherited

Append graphics object to current pad.

In case no current pad is set yet, create a default canvas with the name "c1".

Definition at line 204 of file TObject.cxx.

◆ Background() [1/2]

TH1 * TSpectrum::Background ( const TH1 * h,
Int_t nIter = 20,
Option_t * option = "" )
virtual

One-dimensional background estimation function.

This function calculates the background spectrum in the input histogram h. The background is returned as a histogram.

Parameters:

  • h: input 1-d histogram
  • nIter, (default value = 20) Increasing number of iterations makes the result smoother and lower.
  • option: may contain one of the following options:
    • to set the direction parameter "BackIncreasingWindow". By default the direction is BackDecreasingWindow
    • filterOrder-order of clipping filter, (default "BackOrder2") -possible values= "BackOrder4" "BackOrder6" "BackOrder8"
    • "nosmoothing"- if selected, the background is not smoothed By default the background is smoothed.
    • smoothWindow-width of smoothing window, (default is "BackSmoothing3") -possible values= "BackSmoothing5" "BackSmoothing7" "BackSmoothing9" "BackSmoothing11" "BackSmoothing13" "BackSmoothing15"
    • "Compton" if selected the estimation of Compton edge will be included.
    • "same" : if this option is specified, the resulting background histogram is superimposed on the picture in the current pad.

NOTE that the background is only evaluated in the current range of h. ie, if h has a bin range (set via h->GetXaxis()->SetRange(binmin,binmax), the returned histogram will be created with the same number of bins as the input histogram h, but only bins from binmin to binmax will be filled with the estimated background.

Definition at line 144 of file TSpectrum.cxx.

◆ Background() [2/2]

const char * TSpectrum::Background ( Double_t * spectrum,
Int_t ssize,
Int_t numberIterations,
Int_t direction,
Int_t filterOrder,
bool smoothing,
Int_t smoothWindow,
bool compton )

This function calculates background spectrum from source spectrum.

The result is placed in the vector pointed by spectrum pointer. The goal is to separate the useful information (peaks) from useless information (background).

  • method is based on Sensitive Nonlinear Iterative Peak (SNIP) clipping algorithm.
  • new value in the channel "i" is calculated

\[ v_p(i) = min \left\{ v_{p-1}(i)^{\frac{\left[v_{p-1}(i+p)+v_{p-1}(i-p)\right]}{2}} \right\} \]

where p = 1, 2, ..., numberIterations. In fact it represents second order difference filter (-1,2,-1).

One can also change the direction of the change of the clipping window, the order of the clipping filter, to include smoothing, to set width of smoothing window and to include the estimation of Compton edges. On successful completion it returns 0. On error it returns pointer to the string describing error.

Parameters:

  • spectrum: pointer to the vector of source spectrum
  • ssize: length of the spectrum vector
  • numberIterations: maximal width of clipping window,
  • direction: direction of change of clipping window. Possible values: kBackIncreasingWindow, kBackDecreasingWindow
  • filterOrder: order of clipping filter. Possible values: kBackOrder2, kBackOrder4, kBackOrder6, kBackOrder8
  • smoothing: logical variable whether the smoothing operation in the estimation of background will be included. Possible values: kFALSE, kTRUE
  • smoothWindow: width of smoothing window. Possible values: kBackSmoothing3, kBackSmoothing5, kBackSmoothing7, kBackSmoothing9, kBackSmoothing11, kBackSmoothing13, kBackSmoothing15.
  • compton: logical variable whether the estimation of Compton edge will be included. Possible values: kFALSE, kTRUE.

References:

  1. C. G Ryan et al.: SNIP, a statistics-sensitive background treatment for the quantitative analysis of PIXE spectra in geoscience applications. NIM, B34 (1988), 396-402.
  2. M. Morhac;, J. Kliman, V. Matouoek, M. Veselsky, I. Turzo: Background elimination methods for multidimensional gamma-ray spectra. NIM, A401 (1997) 113-132.
  3. D. D. Burgess, R. J. Tervo: Background estimation for gamma-ray spectroscopy. NIM 214 (1983), 431-434.

Example 1 script Background_incr.C:

Example of the estimation of background for number of iterations=6. Original spectrum is shown in black color, estimated background in red color.

{
Int_t i;
const Int_t nbins = 1024;
Double_t xmax = nbins;
Double_t source[nbins];
gROOT->ForceStyle();
TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back1");
back->GetXaxis()->SetRange(1, nbins);
back->SetTitle("Estimation of background with increasing window");
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
// Estimate the background
// Draw the estimated background
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}
#define d(i)
Definition RSha256.hxx:102
#define f(i)
Definition RSha256.hxx:104
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
@ kRed
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
float xmin
float xmax
#define gROOT
Definition TROOT.h:417
virtual void SetRange(Int_t first=0, Int_t last=0)
Set the viewing range for the axis using bin numbers.
Definition TAxis.cxx:1061
A file, usually with extension .root, that stores data and code in the form of serialized objects in ...
Definition TFile.h:130
1-D histogram with a float per channel (see TH1 documentation)
Definition TH1.h:878
void SetTitle(const char *title) override
Change/set the title.
Definition TH1.cxx:6836
TAxis * GetXaxis()
Definition TH1.h:571
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3097
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5143
Advanced Spectra Processing.
Definition TSpectrum.h:18
@ kBackOrder2
Definition TSpectrum.h:37
@ kBackIncreasingWindow
Definition TSpectrum.h:41
@ kBackSmoothing3
Definition TSpectrum.h:43
virtual TH1 * Background(const TH1 *hist, Int_t nIter=20, Option_t *option="")
One-dimensional background estimation function.
Basic string class.
Definition TString.h:138
const char * Data() const
Definition TString.h:384

Example 2 script Background_decr.C:

In Example 1. one can notice that at the edges of the peaks the estimated background goes under the peaks. An alternative approach is to decrease the clipping window from a given value numberIterations to the value of one, which is presented in this example.

Example of the estimation of background for numberIterations=6 using decreasing clipping window algorithm. Original spectrum is shown in black color, estimated background in red color.

{
Int_t i;
const Int_t nbins = 1024;
Double_t xmax = nbins;
Double_t source[nbins];
gROOT->ForceStyle();
TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back1");
back->SetTitle("Estimation of background with decreasing window");
back->GetXaxis()->SetRange(1, nbins);
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
// Estimate the background
// Draw the estimated background
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}
@ kBackDecreasingWindow
Definition TSpectrum.h:42

Example 3 script Background_width.C:

The question is how to choose the width of the clipping window, i.e., numberIterations parameter. The influence of this parameter on the estimated background is illustrated in Example 3.

Example of the influence of clipping window width on the estimated background for numberIterations=4 (red line), 6 (orange line) 8 (green line) using decreasing clipping window algorithm.

in general one should set this parameter so that the value 2*numberIterations+1 was greater than the widths of preserved objects (peaks).

{
Int_t i;
const Int_t nbins = 1024;
Double_t xmax = nbins;
Double_t source[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back1");
TH1F *d1 = new TH1F("d1", "", nbins, xmin, xmax);
TH1F *d2 = new TH1F("d2", "", nbins, xmin, xmax);
TH1F *d3 = new TH1F("d3", "", nbins, xmin, xmax);
back->GetXaxis()->SetRange(1, nbins);
back->SetTitle("Influence of clipping window width on the estimated background");
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d1->SetBinContent(i + 1, source[i]);
d1->SetLineColor(kRed);
d1->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d2->SetBinContent(i + 1, source[i]);
d2->SetLineColor(kOrange);
d2->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d3->SetBinContent(i + 1, source[i]);
d3->SetLineColor(kGreen);
d3->Draw("SAME L");
}
@ kOrange
Definition Rtypes.h:68
@ kGreen
Definition Rtypes.h:67

Example 4 script Background_width2.C:

another example for very complex spectrum is given here.

Example of the influence of clipping window width on the estimated background for numberIterations=10 (red line), 20 (blue line), 30 (green line) and 40 (magenta line) using decreasing clipping window algorithm.

{
Int_t i;
const Int_t nbins = 4096;
Double_t xmax = 4096;
Double_t source[nbins];
gROOT->ForceStyle();
TH1F *d1 = new TH1F("d1", "", nbins, xmin, xmax);
TH1F *d2 = new TH1F("d2", "", nbins, xmin, xmax);
TH1F *d3 = new TH1F("d3", "", nbins, xmin, xmax);
TH1F *d4 = new TH1F("d4", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back2");
back->SetTitle("Influence of clipping window width on the estimated background");
back->SetAxisRange(0, 1000);
back->SetMaximum(7000);
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d1->SetBinContent(i + 1, source[i]);
d1->SetLineColor(kRed);
d1->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d2->SetBinContent(i + 1, source[i]);
d2->SetLineColor(kBlue);
d2->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d3->SetBinContent(i + 1, source[i]);
d3->SetLineColor(kGreen);
d3->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d4->SetBinContent(i + 1, source[i]);
d4->SetLineColor(kMagenta);
d4->Draw("SAME L");
}
@ kMagenta
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
virtual void SetMaximum(Double_t maximum=-1111)
Definition TH1.h:652
virtual void SetAxisRange(Double_t xmin, Double_t xmax, Option_t *axis="X")
Set the "axis" range.
Definition Haxis.cxx:201

Example 5 script Background_order.C:

Second order difference filter removes linear (quasi-linear) background and preserves symmetrical peaks. However if the shape of the background is more complex one can employ higher-order clipping filters.

Example of the influence of clipping filter difference order on the estimated background for fNnumberIterations=40, 2-nd order red line, 4-th order blue line, 6-th order green line and 8-th order magenta line, and using decreasing clipping window algorithm.

{
Int_t i;
const Int_t nbins = 4096;
Double_t xmax = 4096;
Double_t source[nbins];
gROOT->ForceStyle();
TH1F *d1 = new TH1F("d1", "", nbins, xmin, xmax);
TH1F *d2 = new TH1F("d2", "", nbins, xmin, xmax);
TH1F *d3 = new TH1F("d3", "", nbins, xmin, xmax);
TH1F *d4 = new TH1F("d4", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back2");
back->SetTitle("Influence of clipping filter difference order on the estimated background");
back->SetAxisRange(1220, 1460);
back->SetMaximum(3000);
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d1->SetBinContent(i + 1, source[i]);
d1->SetLineColor(kRed);
d1->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d2->SetBinContent(i + 1, source[i]);
d2->SetLineColor(kBlue);
d2->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d3->SetBinContent(i + 1, source[i]);
d3->SetLineColor(kGreen);
d3->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d4->SetBinContent(i + 1, source[i]);
d4->SetLineColor(kMagenta);
d4->Draw("SAME L");
}
@ kBackOrder8
Definition TSpectrum.h:40
@ kBackOrder6
Definition TSpectrum.h:39
@ kBackOrder4
Definition TSpectrum.h:38

Example 6 script Background_smooth.C:

The estimate of the background can be influenced by noise present in the spectrum. We proposed the algorithm of the background estimate with simultaneous smoothing. In the original algorithm without smoothing, the estimated background snatches the lower spikes in the noise. Consequently, the areas of peaks are biased by this error.

Principle of background estimation algorithm with simultaneous smoothing.
{
Int_t i;
const Int_t nbins = 4096;
Double_t xmax = nbins;
Double_t source[nbins];
gROOT->ForceStyle();
TH1F *d1 = new TH1F("d1", "", nbins, xmin, xmax);
TH1F *d2 = new TH1F("d2", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back1");
back->SetTitle("Estimation of background with noise");
back->SetAxisRange(3460, 3830);
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d1->SetBinContent(i + 1, source[i]);
d1->SetLineColor(kRed);
d1->Draw("SAME L");
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d2->SetBinContent(i + 1, source[i]);
d2->SetLineColor(kBlue);
d2->Draw("SAME L");
}
constexpr Bool_t kTRUE
Definition RtypesCore.h:107

Example 8 script Background_compton.C:

Sometimes it is necessary to include also the Compton edges into the estimate of the background. This example presents the synthetic spectrum with Compton edges. The background was estimated using the 8-th order filter with the estimation of the Compton edges using decreasing clipping window algorithm (numberIterations=10) with smoothing (smoothingWindow=5).

Example of the estimate of the background with Compton edges (red line) for numberIterations=10, 8-th order difference filter, using decreasing clipping window algorithm and smoothing (smoothingWindow=5).

{
Int_t i;
const Int_t nbins = 256;
Double_t xmax = nbins;
Double_t source[nbins];
gROOT->ForceStyle();
TH1F *d1 = new TH1F("d1", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *back = (TH1F *)f->Get("back3");
back->SetTitle("Estimation of background with Compton edges under peaks");
back->GetXaxis()->SetRange(1, nbins);
back->Draw("L");
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = back->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
d1->SetBinContent(i + 1, source[i]);
d1->SetLineColor(kRed);
d1->Draw("SAME L");
}
@ kBackSmoothing5
Definition TSpectrum.h:44

Definition at line 504 of file TSpectrum.cxx.

◆ Browse()

◆ CheckedHash()

ULong_t TObject::CheckedHash ( )
inlineinherited

Check and record whether this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object.

The intent is for this routine to be called instead of directly calling the function Hash during "insert" operations. See TObject::HasInconsistenTObjectHash();

(*) The setup is consistent when all classes in the class hierarchy that overload TObject::Hash do call ROOT::CallRecursiveRemoveIfNeeded in their destructor. i.e. it is safe to call the Hash virtual function during the RecursiveRemove operation.

Definition at line 332 of file TObject.h.

◆ Class()

static TClass * TSpectrum::Class ( )
static
Returns
TClass describing this class

◆ Class_Name()

static const char * TSpectrum::Class_Name ( )
static
Returns
Name of this class

◆ Class_Version()

static constexpr Version_t TSpectrum::Class_Version ( )
inlinestaticconstexpr
Returns
Version of this class

Definition at line 78 of file TSpectrum.h.

◆ ClassName()

const char * TObject::ClassName ( ) const
virtualinherited

Returns name of class to which the object belongs.

Definition at line 227 of file TObject.cxx.

◆ Clear()

void TNamed::Clear ( Option_t * option = "")
overridevirtualinherited

Set name and title to empty strings ("").

Reimplemented from TObject.

Reimplemented in TStreamerInfo, TVirtualStreamerInfo, TProcessID, TTask, TPrincipal, and TVirtualFitter.

Definition at line 63 of file TNamed.cxx.

◆ Clone()

TObject * TNamed::Clone ( const char * newname = "") const
overridevirtualinherited

Make a clone of an object using the Streamer facility.

If newname is specified, this will be the name of the new object.

Reimplemented from TObject.

Reimplemented in TStreamerInfo, and TTreeIndex.

Definition at line 73 of file TNamed.cxx.

◆ Compare()

Int_t TNamed::Compare ( const TObject * obj) const
overridevirtualinherited

Compare two TNamed objects.

Returns 0 when equal, -1 when this is smaller and +1 when bigger (like strcmp).

Reimplemented from TObject.

Reimplemented in TStructNodeProperty.

Definition at line 84 of file TNamed.cxx.

◆ Copy()

void TNamed::Copy ( TObject & named) const
overridevirtualinherited

Copy this to obj.

Reimplemented from TObject.

Reimplemented in TSystemDirectory, TSystemFile, TProfile, TProfile2D, TProfile3D, TPieSlice, TStyle, TText, and TXTRU.

Definition at line 93 of file TNamed.cxx.

◆ DeclFileName()

static const char * TSpectrum::DeclFileName ( )
inlinestatic
Returns
Name of the file containing the class declaration

Definition at line 78 of file TSpectrum.h.

◆ Deconvolution()

const char * TSpectrum::Deconvolution ( Double_t * source,
const Double_t * response,
Int_t ssize,
Int_t numberIterations,
Int_t numberRepetitions,
Double_t boost )

One-dimensional deconvolution function.

This function calculates deconvolution from source spectrum according to response spectrum using Gold deconvolution algorithm. The result is placed in the vector pointed by source pointer. On successful completion it returns 0. On error it returns pointer to the string describing error. If desired after every numberIterations one can apply boosting operation (exponential function with exponent given by boost coefficient) and repeat it numberRepetitions times.

Parameters:

  • source: pointer to the vector of source spectrum
  • response: pointer to the vector of response spectrum
  • ssize: length of source and response spectra
  • numberIterations, for details we refer to the reference given below
  • numberRepetitions, for repeated boosted deconvolution
  • boost, boosting coefficient

The goal of this function is the improvement of the resolution in spectra, decomposition of multiplets. The mathematical formulation of the convolution system is:

\[ y(i) = \sum_{k=0}^{N-1} h(i-k)x(k), i=0,1,2,...,N-1 \]

where h(i) is the impulse response function, x, y are input and output vectors, respectively, N is the length of x and h vectors. In matrix form we have:

\[ \begin{bmatrix} y(0) \\ y(1) \\ \dots \\ y(2N-2) \\ y(2N-1) \end{bmatrix} = \begin{bmatrix} h(0) & 0 & 0 & \dots & 0 \\ h(1) & h(0) & 0 & \dots & \dots \\ \dots & h(1) & h(0) & \dots & \dots \\ \dots & \dots & h(1) & \dots & \dots \\ \dots & \dots & \dots & \dots & \dots \\ h(N-1) & \dots & \dots &\dots & 0 \\ 0 & h(N-1) & \dots & \dots & h(0) \\ 0 & 0 & h(N-1) & \dots & h(1) \\ \dots & \dots & \dots & \dots & \dots \\ 0 & 0 & 0 & \dots & h(N-1) \end{bmatrix} \begin{bmatrix} x(0) \\ x(1) \\ \dots \\ x(N-1) \end{bmatrix} \]

Let us assume that we know the response and the output vector (spectrum) of the above given system. The deconvolution represents solution of the overdetermined system of linear equations, i.e., the calculation of the vector x. From numerical stability point of view the operation of deconvolution is extremely critical (ill-posed problem) as well as time consuming operation. The Gold deconvolution algorithm proves to work very well, other methods (Fourier, VanCittert etc) oscillate. It is suitable to process positive definite data (e.g. histograms).

Gold deconvolution algorithm:

\[ y = Hx \\ H^{T}y = H^{T}Hx \\ y^{'} = H^{'}x \\ x_{i}^{(k+1)} = \frac{y_{i}^{'}}{\sum_{m=0}^{N-1}H_{im}^{'}x_{m}{(k)}}x_{i}{(k)}, i=0,1,2,...,N-1 \\ k = 1,2,3,...,L x^{0} = [1,1, ..., 1]^T \]

Where L is given number of iterations (numberIterations parameter).

Boosted deconvolution:

  1. Set the initial solution: \( x^{(0)} = [1,1,...,1]^{T} \)
  2. Set required number of repetitions R and iterations L.
  3. Set r = 1.
  4. Using Gold deconvolution algorithm for k=1,2,...,L find \( x^{(L)} \)
  5. If r = R stop calculation, else
    1. Apply boosting operation, i.e., set \( x^{(0)}(i) = [x^{(L)}(i)]^{p} \) i=0,1,...N-1 and p is boosting coefficient >0.
    2. r = r + 1
    3. continue in 4.

References:

  1. Gold R., ANL-6984, Argonne National Laboratories, Argonne Ill, 1964.
  2. Coote G.E., Iterative smoothing and deconvolution of one- and two-dimensional elemental distribution data, NIM B 130 (1997) 118.
  3. M. Morhac;, J. Kliman, V. Matouoek, M. Veselsky, I. Turzo: Efficient one- and two-dimensional Gold deconvolution and its application to gamma-ray spectra decomposition. NIM, A401 (1997) 385-408.
  4. Morhac; M., Matouoek V., Kliman J., Efficient algorithm of multidimensional deconvolution and its application to nuclear data processing, Digital Signal Processing 13 (2003) 144.

Example 8 - script Deconvolution.C :

response function (usually peak) should be shifted left to the first non-zero channel (bin).

Principle how the response matrix is composed inside of the Deconvolution function.
void Deconvolution()
{
Int_t i;
const Int_t nbins = 256;
Double_t xmax = nbins;
Double_t source[nbins];
Double_t response[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("decon1");
h->SetTitle("Deconvolution");
TH1F *d = (TH1F *)f->Get("decon_response");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
response[i] = d->GetBinContent(i + 1);
h->SetMaximum(30000);
h->Draw("L");
TSpectrum *s = new TSpectrum();
s->Deconvolution(source, response, 256, 1000, 1, 1);
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}
#define h(i)
Definition RSha256.hxx:106
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.

Examples of Gold deconvolution method:

First let us study the influence of the number of iterations on the deconvolved spectrum (Fig. 12).

Fig. 12 Study of Gold deconvolution algorithm.The original source spectrum is drawn with black color, spectrum after 100 iterations with red color, spectrum after 1000 iterations with blue color, spectrum after 10000 iterations with green color and spectrum after 100000 iterations with magenta color.

For relatively narrow peaks in the above given example the Gold deconvolution method is able to decompose overlapping peaks practically to delta - functions. In the next example we have chosen a synthetic data (spectrum, 256 channels) consisting of 5 very closely positioned, relatively wide peaks (sigma =5), with added noise (Fig. 13). Thin lines represent pure Gaussians (see Table 1); thick line is a resulting spectrum with additive noise (10% of the amplitude of small peaks).

Fig. 13 Testing example of synthetic spectrum composed of 5 Gaussians with added noise.
Peak # Position Height Area
1 50 500 10159
2 70 3000 60957
3 80 1000 20319
4 100 5000 101596
5 110 500 10159

Table 1 Positions, heights and areas of peaks in the spectrum shown in Fig. 13.

In ideal case, we should obtain the result given in Fig. 14. The areas of the Gaussian components of the spectrum are concentrated completely to delta-functions. When solving the overdetermined system of linear equations with data from Fig. 13 in the sense of minimum least squares criterion without any regularisation we obtain the result with large oscillations (Fig. 15). From mathematical point of view, it is the optimal solution in the unconstrained space of independent variables. From physical point of view we are interested only in a meaningful solution. Therefore, we have to employ regularisation techniques (e.g. Gold deconvolution) and/or to confine the space of allowed solutions to subspace of positive solutions.

Fig. 14 The same spectrum like in Fig. 13, outlined bars show the contents of present components (peaks).
Fig. 15 Least squares solution of the system of linear equations without regularisation.

Example 9 - script Deconvolution_wide.C

When we employ Gold deconvolution algorithm we obtain the result given in Fig. 16. One can observe that the resulting spectrum is smooth. On the other hand the method is not able to decompose completely the peaks in the spectrum.

Example of Gold deconvolution for closely positioned wide peaks. The original source spectrum is drawn with black color, the spectrum after the deconvolution (10000 iterations) with red color.

{
Int_t i;
const Int_t nbins = 256;
Double_t xmax = nbins;
Double_t source[nbins];
Double_t response[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("decon3");
h->SetTitle("Deconvolution of closely positioned overlapping peaks using Gold deconvolution method");
TH1F *d = (TH1F *)f->Get("decon_response_wide");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
response[i] = d->GetBinContent(i + 1);
h->SetMaximum(50000);
h->Draw("L");
TSpectrum *s = new TSpectrum();
s->Deconvolution(source, response, 256, 10000, 1, 1);
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}

Example 10 - script Deconvolution_wide_boost.C :

Further let us employ boosting operation into deconvolution (Fig. 17).

The original source spectrum is drawn with black color, the spectrum after the deconvolution with red color. Number of iterations = 200, number of repetitions = 50 and boosting coefficient = 1.2.

One can observe that peaks are decomposed practically to delta functions. Number of peaks is correct, positions of big peaks as well as their areas are relatively well estimated. However there is a considerable error in the estimation of the position of small right hand peak.

{
Int_t i;
const Int_t nbins = 256;
Double_t xmax = nbins;
Double_t source[nbins];
Double_t response[nbins];
gROOT->ForceStyle();
TH1F *h = new TH1F("h", "Deconvolution", nbins, xmin, xmax);
TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
h = (TH1F *)f->Get("decon3");
h->SetTitle("Deconvolution of closely positioned overlapping peaks using boosted Gold deconvolution method");
d = (TH1F *)f->Get("decon_response_wide");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
response[i] = d->GetBinContent(i + 1);
h->SetMaximum(200000);
h->Draw("L");
TSpectrum *s = new TSpectrum();
s->Deconvolution(source, response, 256, 200, 50, 1.2);
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}

Definition at line 1451 of file TSpectrum.cxx.

◆ DeconvolutionRL()

const char * TSpectrum::DeconvolutionRL ( Double_t * source,
const Double_t * response,
Int_t ssize,
Int_t numberIterations,
Int_t numberRepetitions,
Double_t boost )

One-dimensional deconvolution function.

This function calculates deconvolution from source spectrum according to response spectrum using Richardson-Lucy deconvolution algorithm. The result is placed in the vector pointed by source pointer. On successful completion it returns 0. On error it returns pointer to the string describing error. If desired after every numberIterations one can apply boosting operation (exponential function with exponent given by boost coefficient) and repeat it numberRepetitions times (see Gold deconvolution).

Parameters:

  • source: pointer to the vector of source spectrum
  • response: pointer to the vector of response spectrum
  • ssize: length of source and response spectra
  • numberIterations, for details we refer to the reference given above
  • numberRepetitions, for repeated boosted deconvolution
  • boost, boosting coefficient

Richardson-Lucy deconvolution algorithm:

For discrete systems it has the form:

\[ x^{(n)}(i) = x^{(n-1)}(i) \sum_{j=0}^{N-1}h(i,j)\frac{y(j)}{\sum_{k=0}^{M-1}h(j,k)x^{(n-1)}(k)} \\ i \in \left<0,M-1\right> \]

for positive input data and response matrix this iterative method forces the deconvoluted spectra to be non-negative. The Richardson-Lucy iteration converges to the maximum likelihood solution for Poisson statistics in the data.

References:

  1. Abreu M.C. et al., A four-dimensional deconvolution method to correct NA38 experimental data, NIM A 405 (1998) 139.
  2. Lucy L.B., A.J. 79 (1974) 745.
  3. Richardson W.H., J. Opt. Soc. Am. 62 (1972) 55.

Examples of Richardson-Lucy deconvolution method:

Example 11 - script DeconvolutionRL_wide.C :

When we employ Richardson-Lucy deconvolution algorithm to our data from Fig. 13 we obtain the following result. One can observe improvements as compared to the result achieved by Gold deconvolution. Nevertheless it is unable to decompose the multiplet.

Example of Richardson-Lucy deconvolution for closely positioned wide peaks. The original source spectrum is drawn with black color, the spectrum after the deconvolution (10000 iterations) with red color.

{
Int_t i;
const Int_t nbins = 256;
Double_t xmax = nbins;
Double_t source[nbins];
Double_t response[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("decon3");
h->SetTitle("Deconvolution of closely positioned overlapping peaks using Richardson-Lucy deconvolution method");
TH1F *d = (TH1F *)f->Get("decon_response_wide");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
response[i] = d->GetBinContent(i + 1);
h->SetMaximum(30000);
h->Draw("L");
TSpectrum *s = new TSpectrum();
s->DeconvolutionRL(source, response, 256, 10000, 1, 1);
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}
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.

Example 12 - script DeconvolutionRL_wide_boost.C :

Further let us employ boosting operation into deconvolution.

The original source spectrum is drawn with black color, the spectrum after the deconvolution with red color. Number of iterations = 200, number of repetitions = 50 and boosting coefficient = 1.2.

One can observe improvements in the estimation of peak positions as compared to the results achieved by Gold deconvolution.

{
Int_t i;
const Int_t nbins = 256;
Double_t xmax = nbins;
Double_t source[nbins];
Double_t response[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("decon3");
"Deconvolution of closely positioned overlapping peaks using boosted Richardson-Lucy deconvolution method");
TH1F *d = (TH1F *)f->Get("decon_response_wide");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
for (i = 0; i < nbins; i++)
response[i] = d->GetBinContent(i + 1);
h->Draw("L");
TSpectrum *s = new TSpectrum();
s->DeconvolutionRL(source, response, 256, 200, 50, 1.2);
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, source[i]);
d->SetLineColor(kRed);
d->Draw("SAME L");
}

Definition at line 1657 of file TSpectrum.cxx.

◆ Delete()

void TObject::Delete ( Option_t * option = "")
virtualinherited

◆ DistancetoPrimitive()

◆ DoError()

void TObject::DoError ( int level,
const char * location,
const char * fmt,
va_list va ) const
protectedvirtualinherited

Interface to ErrorHandler (protected).

Reimplemented in TTreeViewer, and TThread.

Definition at line 1059 of file TObject.cxx.

◆ Draw()

◆ DrawClass()

void TObject::DrawClass ( ) const
virtualinherited

Draw class inheritance tree of the class to which this object belongs.

If a class B inherits from a class A, description of B is drawn on the right side of description of A. Member functions overridden by B are shown in class A with a blue line crossing-out the corresponding member function. The following picture is the class inheritance tree of class TPaveLabel:

Reimplemented in TSystemDirectory, TSystemFile, and TGFrame.

Definition at line 308 of file TObject.cxx.

◆ DrawClone()

TObject * TObject::DrawClone ( Option_t * option = "") const
virtualinherited

Draw a clone of this object in the current selected pad with: gROOT->SetSelectedPad(c1).

If pad was not selected - gPad will be used.

Note
For histograms, use the more specialised TH1::DrawCopy().

Reimplemented in TSystemDirectory, TSystemFile, TGFrame, TAxis, and TCanvas.

Definition at line 319 of file TObject.cxx.

◆ Dump()

void TObject::Dump ( ) const
virtualinherited

Dump contents of object on stdout.

Using the information in the object dictionary (class TClass) each data member is interpreted. If a data member is a pointer, the pointer value is printed

The following output is the Dump of a TArrow object:

fAngle 0 Arrow opening angle (degrees)
fArrowSize 0.2 Arrow Size
fOption.*fData
fX1 0.1 X of 1st point
fY1 0.15 Y of 1st point
fX2 0.67 X of 2nd point
fY2 0.83 Y of 2nd point
fBits 50331648 bit field status word
fLineColor 1 line color
fLineStyle 1 line style
fLineWidth 1 line width
fFillColor 19 fill area color
fFillStyle 1001 fill area style
#define X(type, name)
Option_t Option_t TPoint TPoint angle
Option_t Option_t width
Option_t Option_t style
UInt_t fUniqueID
object unique identifier
Definition TObject.h:46
UInt_t fBits
bit field status word
Definition TObject.h:47
TLine * line

Reimplemented in TSystemFile, TCollection, TClass, TGFrame, and TGPack.

Definition at line 367 of file TObject.cxx.

◆ Error()

void TObject::Error ( const char * location,
const char * fmt,
... ) const
virtualinherited

Issue error message.

Use "location" to specify the method where the error occurred. Accepts standard printf formatting arguments.

Reimplemented in TFitResult.

Definition at line 1098 of file TObject.cxx.

◆ Execute() [1/2]

void TObject::Execute ( const char * method,
const char * params,
Int_t * error = nullptr )
virtualinherited

Execute method on this object with the given parameter string, e.g.

"3.14,1,\"text\"".

Reimplemented in TMethodCall, TCling, TInterpreter, ROOT::R::TRInterface, and TContextMenu.

Definition at line 378 of file TObject.cxx.

◆ Execute() [2/2]

void TObject::Execute ( TMethod * method,
TObjArray * params,
Int_t * error = nullptr )
virtualinherited

Execute method on this object with parameters stored in the TObjArray.

The TObjArray should contain an argv vector like:

argv[0] ... argv[n] = the list of TObjString parameters
Collectable string class.
Definition TObjString.h:28
const Int_t n
Definition legend1.C:16

Reimplemented in TCling, TMethodCall, TInterpreter, ROOT::R::TRInterface, and TContextMenu.

Definition at line 398 of file TObject.cxx.

◆ ExecuteEvent()

◆ Fatal()

void TObject::Fatal ( const char * location,
const char * fmt,
... ) const
virtualinherited

Issue fatal error message.

Use "location" to specify the method where the fatal error occurred. Accepts standard printf formatting arguments.

Definition at line 1126 of file TObject.cxx.

◆ FillBuffer()

void TNamed::FillBuffer ( char *& buffer)
virtualinherited

Encode TNamed into output buffer.

Reimplemented in TKeySQL, TSQLFile, TKeyXML, TXMLFile, TDirectoryFile, TFile, and TKey.

Definition at line 103 of file TNamed.cxx.

◆ FindObject() [1/2]

TObject * TObject::FindObject ( const char * name) const
virtualinherited

Must be redefined in derived classes.

This function is typically used with TCollections, but can also be used to find an object by name inside this object.

Reimplemented in TListOfEnums, TMap, TDirectory, TFolder, TROOT, TListOfTypes, TListOfTypes, TBtree, TCollection, THashList, THashTable, TList, TObjArray, TListOfDataMembers, TListOfDataMembers, TListOfEnums, TListOfEnumsWithLock, TListOfFunctions, TListOfFunctionTemplates, TListOfFunctionTemplates, TViewPubDataMembers, TViewPubFunctions, TPad, TGeometry, THbookFile, TGraph, TGraph2D, TH1, RooAbsCollection, and RooLinkedList.

Definition at line 425 of file TObject.cxx.

◆ FindObject() [2/2]

TObject * TObject::FindObject ( const TObject * obj) const
virtualinherited

Must be redefined in derived classes.

This function is typically used with TCollections, but can also be used to find an object inside this object.

Reimplemented in TMap, TDirectory, TFolder, TROOT, TListOfTypes, TBtree, TCollection, THashList, THashTable, TList, TObjArray, TListOfDataMembers, TListOfEnums, TListOfEnumsWithLock, TListOfFunctions, TListOfFunctionTemplates, TViewPubDataMembers, TViewPubFunctions, TPad, TGeometry, THbookFile, TGraph, TGraph2D, TH1, RooAbsCollection, and RooLinkedList.

Definition at line 435 of file TObject.cxx.

◆ GetDrawOption()

Option_t * TObject::GetDrawOption ( ) const
virtualinherited

Get option used by the graphics system to draw this object.

Note that before calling object.GetDrawOption(), you must have called object.Draw(..) before in the current pad.

Reimplemented in TBrowser, TFitEditor, TGedFrame, TGFileBrowser, TRootBrowser, and TRootBrowserLite.

Definition at line 445 of file TObject.cxx.

◆ GetDtorOnly()

Longptr_t TObject::GetDtorOnly ( )
staticinherited

Return destructor only flag.

Definition at line 1196 of file TObject.cxx.

◆ GetHistogram()

TH1 * TSpectrum::GetHistogram ( ) const
inline

Definition at line 56 of file TSpectrum.h.

◆ GetIconName()

const char * TObject::GetIconName ( ) const
virtualinherited

Returns mime type name of object.

Used by the TBrowser (via TGMimeTypes class). Override for class of which you would like to have different icons for objects of the same class.

Reimplemented in TSystemFile, TGeoVolume, TASImage, TGMainFrame, TKey, ROOT::Experimental::XRooFit::xRooNode, TBranch, TVirtualBranchBrowsable, TMethodBrowsable, and TBranchElement.

Definition at line 472 of file TObject.cxx.

◆ GetName()

const char * TNamed::GetName ( ) const
inlineoverridevirtualinherited

Returns name of object.

This default method returns the class name. Classes that give objects a name should override this method.

Reimplemented from TObject.

Definition at line 49 of file TNamed.h.

◆ GetNPeaks()

Int_t TSpectrum::GetNPeaks ( ) const
inline

Definition at line 57 of file TSpectrum.h.

◆ GetObjectInfo()

char * TObject::GetObjectInfo ( Int_t px,
Int_t py ) const
virtualinherited

Returns string containing info about the object at position (px,py).

This method is typically overridden by classes of which the objects can report peculiarities for different positions. Returned string will be re-used (lock in MT environment).

Reimplemented in TGeoNode, TGeoVolume, TGeoTrack, TASImage, TColorWheel, TAxis3D, TNode, TGL5DDataSet, TGLHistPainter, TGLParametricEquation, TGLTH3Composition, TF1, TF2, TGraph, TH1, THistPainter, TPaletteAxis, TFileDrawMap, TParallelCoordVar, and TVirtualHistPainter.

Definition at line 491 of file TObject.cxx.

◆ GetObjectStat()

Bool_t TObject::GetObjectStat ( )
staticinherited

Get status of object stat flag.

Definition at line 1181 of file TObject.cxx.

◆ GetOption()

virtual Option_t * TObject::GetOption ( ) const
inlinevirtualinherited

◆ GetPositionX()

Double_t * TSpectrum::GetPositionX ( ) const
inline

Definition at line 58 of file TSpectrum.h.

◆ GetPositionY()

Double_t * TSpectrum::GetPositionY ( ) const
inline

Definition at line 59 of file TSpectrum.h.

◆ GetTitle()

const char * TNamed::GetTitle ( ) const
inlineoverridevirtualinherited

Returns title of object.

This default method returns the class title (i.e. description). Classes that give objects a title should override this method.

Reimplemented from TObject.

Definition at line 50 of file TNamed.h.

◆ GetUniqueID()

UInt_t TObject::GetUniqueID ( ) const
virtualinherited

Return the unique object id.

Definition at line 480 of file TObject.cxx.

◆ HandleTimer()

Bool_t TObject::HandleTimer ( TTimer * timer)
virtualinherited

Execute action in response of a timer timing out.

This method must be overridden if an object has to react to timers.

Reimplemented in TGWindow, TGuiBldDragManager, TGraphTime, TGLEventHandler, TGCommandPlugin, TGDNDManager, TGFileContainer, TGPopupMenu, TGScrollBar, TGShutter, TGTextEdit, TGTextEditor, TGTextEntry, TGTextView, TGToolTip, TGHtml, and TTreeViewer.

Definition at line 516 of file TObject.cxx.

◆ Hash()

ULong_t TNamed::Hash ( ) const
inlineoverridevirtualinherited

Return hash value for this object.

Note: If this routine is overloaded in a derived class, this derived class should also add

void CallRecursiveRemoveIfNeeded(TObject &obj)
call RecursiveRemove for obj if gROOT is valid and obj.TestBit(kMustCleanup) is true.
Definition TROOT.h:406

Otherwise, when RecursiveRemove is called (by ~TObject or example) for this type of object, the transversal of THashList and THashTable containers will will have to be done without call Hash (and hence be linear rather than logarithmic complexity). You will also see warnings like

ULong_t Hash() const override
Return hash value for this object.
Definition TNamed.h:51
Mother of all ROOT objects.
Definition TObject.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1098
void RecursiveRemove(TObject *obj) override
Recursively remove this object from the list of Cleanups.
Definition TROOT.cxx:2651

Reimplemented from TObject.

Definition at line 51 of file TNamed.h.

◆ HasInconsistentHash()

Bool_t TObject::HasInconsistentHash ( ) const
inlineinherited

Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e.

missing call to RecursiveRemove in destructor).

Note: Since the consistency is only tested for during inserts, this routine will return true for object that have never been inserted whether or not they have a consistent setup. This has no negative side-effect as searching for the object with the right or wrong Hash will always yield a not-found answer (Since anyway no hash can be guaranteed unique, there is always a check)

Definition at line 366 of file TObject.h.

◆ Info()

void TObject::Info ( const char * location,
const char * fmt,
... ) const
virtualinherited

Issue info message.

Use "location" to specify the method where the warning occurred. Accepts standard printf formatting arguments.

Definition at line 1072 of file TObject.cxx.

◆ InheritsFrom() [1/2]

Bool_t TObject::InheritsFrom ( const char * classname) const
virtualinherited

Returns kTRUE if object inherits from class "classname".

Reimplemented in TClass.

Definition at line 549 of file TObject.cxx.

◆ InheritsFrom() [2/2]

Bool_t TObject::InheritsFrom ( const TClass * cl) const
virtualinherited

Returns kTRUE if object inherits from TClass cl.

Reimplemented in TClass.

Definition at line 557 of file TObject.cxx.

◆ Inspect()

void TObject::Inspect ( ) const
virtualinherited

Dump contents of this object in a graphics canvas.

Same action as Dump but in a graphical form. In addition pointers to other objects can be followed.

The following picture is the Inspect of a histogram object:

Reimplemented in TSystemFile, TInspectorObject, TGFrame, and ROOT::Experimental::XRooFit::xRooNode.

Definition at line 570 of file TObject.cxx.

◆ InvertBit()

void TObject::InvertBit ( UInt_t f)
inlineinherited

Definition at line 206 of file TObject.h.

◆ IsA()

TClass * TSpectrum::IsA ( ) const
inlineoverridevirtual
Returns
TClass describing current object

Reimplemented from TNamed.

Definition at line 78 of file TSpectrum.h.

◆ IsDestructed()

Bool_t TObject::IsDestructed ( ) const
inlineinherited

IsDestructed.

Note
This function must be non-virtual as it can be used on destructed (but not yet modified) memory. This is used for example in TClonesArray to record the element that have been destructed but not deleted and thus are ready for re-use (by operator new with placement).
Returns
true if this object's destructor has been run.

Definition at line 186 of file TObject.h.

◆ IsEqual()

Bool_t TObject::IsEqual ( const TObject * obj) const
virtualinherited

Default equal comparison (objects are equal if they have the same address in memory).

More complicated classes might want to override this function.

Reimplemented in TObjString, TQCommand, TPair, and TGObject.

Definition at line 589 of file TObject.cxx.

◆ IsFolder()

◆ IsOnHeap()

R__ALWAYS_INLINE Bool_t TObject::IsOnHeap ( ) const
inlineinherited

Definition at line 160 of file TObject.h.

◆ IsSortable()

Bool_t TNamed::IsSortable ( ) const
inlineoverridevirtualinherited

Reimplemented from TObject.

Reimplemented in TStructNodeProperty.

Definition at line 52 of file TNamed.h.

◆ IsZombie()

R__ALWAYS_INLINE Bool_t TObject::IsZombie ( ) const
inlineinherited

Definition at line 161 of file TObject.h.

◆ ls()

void TNamed::ls ( Option_t * option = "") const
overridevirtualinherited

List TNamed name and title.

Reimplemented from TObject.

Reimplemented in ROOT::Experimental::XRooFit::xRooBrowser, TVirtualStreamerInfo, TROOT, TStreamerElement, TStreamerBase, TStreamerSTL, TText, TStreamerInfo, TTask, and TNode.

Definition at line 112 of file TNamed.cxx.

◆ MakeZombie()

void TObject::MakeZombie ( )
inlineprotectedinherited

Definition at line 55 of file TObject.h.

◆ MayNotUse()

void TObject::MayNotUse ( const char * method) const
inherited

Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary).

Definition at line 1160 of file TObject.cxx.

◆ Notify()

Bool_t TObject::Notify ( )
virtualinherited

This method must be overridden to handle object notification (the base implementation is no-op).

Different objects in ROOT use the Notify method for different purposes, in coordination with other objects that call this method at the appropriate time.

For example, TLeaf uses it to load class information; TBranchRef to load contents of referenced branches TBranchRef; most notably, based on Notify, TChain implements a callback mechanism to inform interested parties when it switches to a new sub-tree.

Reimplemented in TMessageHandler, TNotifyLink< Type >, TNotifyLink< RNoCleanupNotifierHelper >, TNotifyLink< ROOT::Detail::TBranchProxy >, TNotifyLink< TTreeReader >, TFileHandler, TSignalHandler, TStdExceptionHandler, TProcessEventTimer, TTimer, TIdleTimer, TSingleShotCleaner, TCollection, TRefTable, TBrowserTimer, TInterruptHandler, TTermInputHandler, TThreadTimer, TGLRedrawTimer, TViewTimer, TGContainerKeyboardTimer, TGContainerScrollTimer, TGInputHandler, TViewUpdateTimer, TPopupDelayTimer, TRepeatTimer, TSBRepeatTimer, TGTextEditHist, TInsCharCom, TDelCharCom, TBreakLineCom, TInsTextCom, TDelTextCom, TBlinkTimer, TTipDelayTimer, TGuiBldDragManagerRepeatTimer, TARInterruptHandler, TASLogHandler, TASInterruptHandler, TASSigPipeHandler, TASInputHandler, TSocketHandler, TTimeOutTimer, TBranchElement, TBranchRef, TLeafObject, TSelector, TTree, TSelectorDraw, TSelectorEntries, TTreeFormula, TTreeFormulaManager, TTreeReader, h1analysis, h1analysisTreeReader, and TSysEvtHandler.

Definition at line 618 of file TObject.cxx.

◆ Obsolete()

void TObject::Obsolete ( const char * method,
const char * asOfVers,
const char * removedFromVers ) const
inherited

Use this method to declare a method obsolete.

Specify as of which version the method is obsolete and as from which version it will be removed.

Definition at line 1169 of file TObject.cxx.

◆ operator delete() [1/3]

void TObject::operator delete ( void * ptr,
size_t size )
inherited

Operator delete for sized deallocation.

Definition at line 1234 of file TObject.cxx.

◆ operator delete() [2/3]

void TObject::operator delete ( void * ptr)
inherited

Operator delete.

Definition at line 1212 of file TObject.cxx.

◆ operator delete() [3/3]

void TObject::operator delete ( void * ptr,
void * vp )
inherited

Only called by placement new when throwing an exception.

Definition at line 1266 of file TObject.cxx.

◆ operator delete[]() [1/3]

void TObject::operator delete[] ( void * ptr,
size_t size )
inherited

Operator delete [] for sized deallocation.

Definition at line 1245 of file TObject.cxx.

◆ operator delete[]() [2/3]

void TObject::operator delete[] ( void * ptr)
inherited

Operator delete [].

Definition at line 1223 of file TObject.cxx.

◆ operator delete[]() [3/3]

void TObject::operator delete[] ( void * ptr,
void * vp )
inherited

Only called by placement new[] when throwing an exception.

Definition at line 1274 of file TObject.cxx.

◆ operator new() [1/2]

void * TObject::operator new ( size_t sz)
inlineinherited

Definition at line 189 of file TObject.h.

◆ operator new() [2/2]

void * TObject::operator new ( size_t sz,
void * vp )
inlineinherited

Definition at line 191 of file TObject.h.

◆ operator new[]() [1/2]

void * TObject::operator new[] ( size_t sz)
inlineinherited

Definition at line 190 of file TObject.h.

◆ operator new[]() [2/2]

void * TObject::operator new[] ( size_t sz,
void * vp )
inlineinherited

Definition at line 192 of file TObject.h.

◆ operator=()

TSpectrum & TSpectrum::operator= ( const TSpectrum & )
private

◆ Paint()

void TObject::Paint ( Option_t * option = "")
virtualinherited

This method must be overridden if a class wants to paint itself.

The difference between Paint() and Draw() is that when a object draws itself it is added to the display list of the pad in which it is drawn (and automatically redrawn whenever the pad is redrawn). While paint just draws the object without adding it to the pad display list.

Reimplemented in TPie, TSQLFile, TXMLFile, ROOT::RGeoPainter, TGaxis, TGraph, TGraphTime, THStack, TMultiGraph, TScatter, TScatter2D, TTreePerfStats, ROOT::Experimental::RTreeMapPainter, TEfficiency, TRatioPlot, TGeoBoolNode, TGeoUnion, TGeoIntersection, TGeoSubtraction, TMarker3DBox, TGL5DDataSet, TGLHistPainter, TGLParametricEquation, TGLTH3Composition, TGraph2DPainter, TSpectrum2Painter, TFileDrawMap, TDirectory, TExec, TMacro, TStyle, TBits, TCollection, TGeoNode, TGeoPhysicalNode, TGeoShape, TGeoVolume, TGeoOverlap, TGeoPainter, TGeoTrack, TGeoVGShape, TASImage, TASPaletteEditor::PaintPalette, TASPaletteEditor::LimitLine, TAnnotation, TButton, TCanvas, TClassTree, TColorWheel, TPad, TArrow, TBox, TCrown, TDiamond, TEllipse, TFrame, TLatex, TLegend, TLine, TMarker, TMathText, TPave, TPaveLabel, TPaveStats, TPavesText, TPaveText, TPolyLine, TText, TWbox, TGraphEdge, TGraphNode, TEveArrow, TEveCaloViz, TEveDigitSet, TEveGeoTopNode, TEveGeoShape, TEvePlot3D, TEvePointSet, TEveProjectionAxes, TEveScene, TEveShape, TEveStraightLineSet, TEveText, TEveTriangleSet, TAxis3D, TNode, TNodeDiv, TPolyLine3D, TPolyMarker3D, TShape, TF1, TF2, TF3, TGraph2D, TH1, TPolyMarker, TSpline, THistPainter, TPaletteAxis, TFile, TGenerator, TParticle, TPrimary, TParallelCoordVar, TVirtualPad, TVirtualGeoPainter, TVirtualGeoTrack, TVirtualHistPainter, TParallelCoordRange, TSpider, TGraphPolargram, and TParallelCoord.

Definition at line 631 of file TObject.cxx.

◆ Pop()

void TObject::Pop ( )
virtualinherited

Pop on object drawn in a pad to the top of the display list.

I.e. it will be drawn last and on top of all other primitives.

Reimplemented in TPad, TFrame, and TVirtualPad.

Definition at line 640 of file TObject.cxx.

◆ Print()

void TSpectrum::Print ( Option_t * option = "") const
overridevirtual

Print the array of positions.

Reimplemented from TNamed.

Definition at line 211 of file TSpectrum.cxx.

◆ Read()

Int_t TObject::Read ( const char * name)
virtualinherited

Read contents of object with specified name from the current directory.

First the key with the given name is searched in the current directory, next the key buffer is deserialized into the object. The object must have been created before via the default constructor. See TObject::Write().

Reimplemented in TKeyXML, TBuffer, TKey, and TKeySQL.

Definition at line 673 of file TObject.cxx.

◆ RecursiveRemove()

◆ ResetBit()

void TObject::ResetBit ( UInt_t f)
inlineinherited

Definition at line 203 of file TObject.h.

◆ SaveAs()

void TObject::SaveAs ( const char * filename = "",
Option_t * option = "" ) const
virtualinherited

Save this object in the file specified by filename.

  • if "filename" contains ".root" the object is saved in filename as root binary file.
  • if "filename" contains ".xml" the object is saved in filename as a xml ascii file.
  • if "filename" contains ".cc" the object is saved in filename as C code independent from ROOT. The code is generated via SavePrimitive(). Specific code should be implemented in each object to handle this option. Like in TF1::SavePrimitive().
  • otherwise the object is written to filename as a CINT/C++ script. The C++ code to rebuild this object is generated via SavePrimitive(). The "option" parameter is passed to SavePrimitive. By default it is an empty string. It can be used to specify the Draw option in the code generated by SavePrimitive.

    The function is available via the object context menu.

Reimplemented in TSpline, TFolder, TGeoVolume, TClassTree, TPad, TPaveClass, TGObject, TSpline3, TSpline5, ROOT::Experimental::XRooFit::xRooNode, TTreePerfStats, TVirtualPad, TGraph, and TH1.

Definition at line 708 of file TObject.cxx.

◆ SavePrimitive()

void TObject::SavePrimitive ( std::ostream & out,
Option_t * option = "" )
virtualinherited

Save a primitive as a C++ statement(s) on output stream "out".

Reimplemented in TGeoTessellated, TGraphEdge, TGraphNode, TGeoIdentity, TStyle, TCurlyArc, TCurlyLine, TGedMarkerSelect, TGedPatternSelect, TGColorSelect, TGFont, TGVerticalLayout, TGHorizontalLayout, TGRowLayout, TGColumnLayout, TGMatrixLayout, TGTileLayout, TGListLayout, TGListDetailsLayout, TGTextLBEntry, TGNumberEntryField, TGNumberEntry, TGTableLayoutHints, TGTableLayout, TGTextEdit, TGTextView, TGXYLayoutHints, TGXYLayout, TRootContainer, TGHtml, TEfficiency, TExec, TMacro, TGeoArb8, TGeoTrap, TGeoGtra, TGeoBBox, TGeoBoolNode, TGeoUnion, TGeoIntersection, TGeoSubtraction, TGeoCompositeShape, TGeoCone, TGeoConeSeg, TGeoElementRN, TGeoDecayChannel, TGeoEltu, TGeoHalfSpace, TGeoHype, TGeoMaterial, TGeoMixture, TGeoTranslation, TGeoRotation, TGeoCombiTrans, TGeoHMatrix, TGeoMedium, TGeoPara, TGeoParaboloid, TGeoPatternX, TGeoPatternY, TGeoPatternZ, TGeoPatternParaX, TGeoPatternParaY, TGeoPatternParaZ, TGeoPatternTrapZ, TGeoPatternCylR, TGeoPatternCylPhi, TGeoPatternSphR, TGeoPatternSphTheta, TGeoPatternSphPhi, TGeoPcon, TGeoPgon, TGeoScaledShape, TGeoShapeAssembly, TGeoSphere, TGeoTorus, TGeoTrd1, TGeoTrd2, TGeoTube, TGeoTubeSeg, TGeoCtub, TGeoVolume, TGeoXtru, TASImage, TAnnotation, TButton, TCanvas, TGroupButton, TPad, TPaveClass, TSlider, TSliderBox, TArc, TArrow, TBox, TCrown, TCutG, TDiamond, TEllipse, TFrame, TGaxis, TGraphPolar, TGraphPolargram, TLatex, TLegend, TLine, TMarker, TMathText, TPave, TPaveLabel, TPaveStats, TPavesText, TPaveText, TPolyLine, TText, TWbox, TGraphStruct, TAxis3D, THelix, TMarker3DBox, TPolyLine3D, TPolyMarker3D, TGHorizontal3DLine, TGVertical3DLine, TGButton, TGTextButton, TGPictureButton, TGCheckButton, TGRadioButton, TGButtonGroup, TGVButtonGroup, TGHButtonGroup, TGContainer, TGCanvas, TGComboBox, TGLineStyleComboBox, TGLineWidthComboBox, TGDockableFrame, TGDoubleVSlider, TGDoubleHSlider, TGFrame, TGCompositeFrame, TGVerticalFrame, TGHorizontalFrame, TGMainFrame, TGTransientFrame, TGGroupFrame, TGFSComboBox, TGFileContainer, TGGC, TGIcon, TGLabel, TGLayoutHints, TGListBox, TGListTree, TGListView, TGLVContainer, TGMdiFrame, TGMdiMainFrame, TGMdiMenuBar, TGPopupMenu, TGMenuTitle, TGMenuBar, TGProgressBar, TGHProgressBar, TGVProgressBar, TGHScrollBar, TGVScrollBar, TGShapedFrame, TGShutterItem, TGShutter, TGVSlider, TGHSlider, TGSplitFrame, TGVSplitter, TGHSplitter, TGVFileSplitter, TGStatusBar, TGTabLayout, TGTab, TGTextEntry, TGToolBar, TGTripleVSlider, TGTripleHSlider, TRootEmbeddedCanvas, TF1, TF12, TF2, TF3, TGraph, TGraph2D, TGraph2DAsymmErrors, TGraph2DErrors, TGraphAsymmErrors, TGraphBentErrors, TGraphErrors, TGraphMultiErrors, TH1, TH2Poly, THStack, TMultiGraph, TPolyMarker, TProfile, TProfile2D, TProfile3D, TScatter, TScatter2D, TSpline3, TSpline5, TPaletteAxis, TChain, TTreePerfStats, TParallelCoord, TParallelCoordVar, TPie, and TPieSlice.

Definition at line 858 of file TObject.cxx.

◆ SavePrimitiveConstructor()

void TObject::SavePrimitiveConstructor ( std::ostream & out,
TClass * cl,
const char * variable_name,
const char * constructor_agrs = "",
Bool_t empty_line = kTRUE )
staticprotectedinherited

Save object constructor in the output stream "out".

Can be used as first statement when implementing SavePrimitive() method for the object

Definition at line 777 of file TObject.cxx.

◆ SavePrimitiveDraw()

void TObject::SavePrimitiveDraw ( std::ostream & out,
const char * variable_name,
Option_t * option = nullptr )
staticprotectedinherited

Save invocation of primitive Draw() method Skipped if option contains "nodraw" string.

Definition at line 845 of file TObject.cxx.

◆ SavePrimitiveNameTitle()

void TNamed::SavePrimitiveNameTitle ( std::ostream & out,
const char * variable_name )
protectedinherited

Save object name and title into the output stream "out".

Definition at line 135 of file TNamed.cxx.

◆ SavePrimitiveVector()

TString TObject::SavePrimitiveVector ( std::ostream & out,
const char * prefix,
Int_t len,
Double_t * arr,
Int_t flag = 0 )
staticprotectedinherited

Save array in the output stream "out" as vector.

Create unique variable name based on prefix value Returns name of vector which can be used in constructor or in other places of C++ code If flag === kTRUE, just add empty line If flag === 111, check if array is empty and return nullptr or <vectorname>.data()

Definition at line 796 of file TObject.cxx.

◆ Search()

Int_t TSpectrum::Search ( const TH1 * hin,
Double_t sigma = 2,
Option_t * option = "",
Double_t threshold = 0.05 )
virtual

One-dimensional peak search function.

This function searches for peaks in source spectrum in hin The number of found peaks and their positions are written into the members fNpeaks and fPositionX. The search is performed in the current histogram range.

Parameters:

  • hin: pointer to the histogram of source spectrum
  • sigma: sigma of searched peaks, for details we refer to manual
  • threshold: (default=0.05) peaks with amplitude less than threshold*highest_peak are discarded. 0<threshold<1

By default, the background is removed before deconvolution. Specify the option "nobackground" to not remove the background.

By default the "Markov" chain algorithm is used. Specify the option "noMarkov" to disable this algorithm Note that by default the source spectrum is replaced by a new spectrum

By default a polymarker object is created and added to the list of functions of the histogram. The histogram is drawn with the specified option and the polymarker object drawn on top of the histogram. The polymarker coordinates correspond to the npeaks peaks found in the histogram.

A pointer to the polymarker object can be retrieved later via:

TList *functions = hin->GetListOfFunctions();
TPolyMarker *pm = (TPolyMarker*)functions->FindObject("TPolyMarker");
A doubly linked list.
Definition TList.h:38
A PolyMarker is defined by an array on N points in a 2-D space.
Definition TPolyMarker.h:31

Specify the option "goff" to disable the storage and drawing of the polymarker.

To disable the final drawing of the histogram with the search results (in case you want to draw it yourself) specify "nodraw" in the options parameter.

Definition at line 258 of file TSpectrum.cxx.

◆ Search1HighRes()

Int_t TSpectrum::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.

This function will be removed after the June 2006 release

Definition at line 2555 of file TSpectrum.cxx.

◆ SearchHighRes()

Int_t TSpectrum::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.

This function searches for peaks in source spectrum. It is based on deconvolution method. First the background is removed (if desired), then Markov smoothed spectrum is calculated (if desired), then the response function is generated according to given sigma and deconvolution is carried out. The order of peaks is arranged according to their heights in the spectrum after background elimination. The highest peak is the first in the list. On success it returns number of found peaks.

Parameters:

  • source: pointer to the vector of source spectrum.
  • destVector: pointer to the vector of resulting deconvolved spectrum.
  • ssize: length of source spectrum.
  • sigma: sigma of searched peaks, for details we refer to manual.
  • threshold: threshold value in % for selected peaks, peaks with amplitude less than threshold*highest_peak/100 are ignored, see manual.
  • backgroundRemove: logical variable, set if the removal of background before deconvolution is desired.
  • deconIterations-number of iterations in deconvolution operation.
  • markov: logical variable, if it is true, first the source spectrum is replaced by new spectrum calculated using Markov chains method.
  • averWindow: averaging window of searched peaks, for details we refer to manual (applies only for Markov method).

Peaks searching:

The goal of this function is to identify automatically the peaks in spectrum with the presence of the continuous background and statistical fluctuations - noise.

The common problems connected with correct peak identification are:

  • non-sensitivity to noise, i.e., only statistically relevant peaks should be identified.
  • non-sensitivity of the algorithm to continuous background.
  • ability to identify peaks close to the edges of the spectrum region. Usually peak finders fail to detect them.
  • resolution, decomposition of Double_tts and multiplets. The algorithm should be able to recognise close positioned peaks.
  • ability to identify peaks with different sigma.
Fig. 27 An example of one-dimensional synthetic spectrum with found peaks denoted by markers.

References:

  1. M.A. Mariscotti: A method for identification of peaks in the presence of background and its application to spectrum analysis. NIM 50 (1967), 309-320.
  2. M. Morhac;, J. Kliman, V. Matouoek, M. Veselsky, I. Turzo.:Identification of peaks in multidimensional coincidence gamma-ray spectra. NIM, A443 (2000) 108-125.
  3. Z.K. Silagadze, A new algorithm for automatic photopeak searches. NIM A 376 (1996), 451.

Examples of peak searching method:

The SearchHighRes function provides users with the possibility to vary the input parameters and with the access to the output deconvolved data in the destination spectrum. Based on the output data one can tune the parameters.

Example 15 - script SearchHR1.C:

One-dimensional spectrum with found peaks denoted by markers, 3 iterations steps in the deconvolution.

Script:

void SearchHR1()
{
Double_t fPositionX[100];
Double_t fPositionY[100];
Int_t fNPeaks = 0;
const Int_t nbins = 1024;
Double_t xmax = nbins;
Double_t source[nbins], dest[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("back2");
h->SetTitle("High resolution peak searching, number of iterations = 3");
h->GetXaxis()->SetRange(1, nbins);
TH1F *d = new TH1F("d", "", nbins, xmin, xmax);
h->Draw("L");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
h->Draw("L");
TSpectrum *s = new TSpectrum();
nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
for (i = 0; i < nfound; i++) {
a = xpeaks[i];
bin = 1 + Int_t(a + 0.5);
fPositionX[i] = h->GetBinCenter(bin);
fPositionY[i] = h->GetBinContent(bin);
}
TPolyMarker *pm = (TPolyMarker *)h->GetListOfFunctions()->FindObject("TPolyMarker");
if (pm) {
h->GetListOfFunctions()->Remove(pm);
delete pm;
}
pm = new TPolyMarker(nfound, fPositionX, fPositionY);
h->GetListOfFunctions()->Add(pm);
pm->SetMarkerStyle(23);
pm->SetMarkerColor(kRed);
pm->SetMarkerSize(1.3);
for (i = 0; i < nbins; i++)
d->SetBinContent(i + 1, dest[i]);
d->SetLineColor(kRed);
d->Draw("SAME");
printf("Found %d candidate peaks\n", nfound);
for (i = 0; i < nfound; i++)
printf("posx= %f, posy= %f\n", fPositionX[i], fPositionY[i]);
}
#define a(i)
Definition RSha256.hxx:99
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t dest
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.
Double_t * GetPositionX() const
Definition TSpectrum.h:58

Example 16 - script SearchHR3.C:

Influence of number of iterations (3-red, 10-blue, 100- green, 1000-magenta), sigma=8, smoothing width=3.

void SearchHR3()
{
Double_t fPositionX[100];
Double_t fPositionY[100];
Int_t fNPeaks = 0;
const Int_t nbins = 1024;
Double_t xmax = nbins;
Double_t source[nbins], dest[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("back2");
h->SetTitle("Influence of # of iterations in deconvolution in peak searching");
h->GetXaxis()->SetRange(1, nbins);
TH1F *d1 = new TH1F("d1", "", nbins, xmin, xmax);
TH1F *d2 = new TH1F("d2", "", nbins, xmin, xmax);
TH1F *d3 = new TH1F("d3", "", nbins, xmin, xmax);
TH1F *d4 = new TH1F("d4", "", nbins, xmin, xmax);
TSpectrum *s = new TSpectrum();
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
nfound = s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 3, kTRUE, 3);
for (i = 0; i < nfound; i++) {
a = xpeaks[i];
bin = 1 + Int_t(a + 0.5);
fPositionX[i] = h->GetBinCenter(bin);
fPositionY[i] = h->GetBinContent(bin);
}
TPolyMarker *pm = (TPolyMarker *)h->GetListOfFunctions()->FindObject("TPolyMarker");
if (pm) {
h->GetListOfFunctions()->Remove(pm);
delete pm;
}
pm = new TPolyMarker(nfound, fPositionX, fPositionY);
h->GetListOfFunctions()->Add(pm);
pm->SetMarkerStyle(23);
pm->SetMarkerColor(kRed);
pm->SetMarkerSize(1.3);
h->Draw("L");
for (i = 0; i < nbins; i++)
d1->SetBinContent(i + 1, dest[i]);
d1->SetLineColor(kRed);
d1->Draw("SAME");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 10, kTRUE, 3);
for (i = 0; i < nbins; i++)
d2->SetBinContent(i + 1, dest[i]);
d2->SetLineColor(kBlue);
d2->Draw("SAME");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 100, kTRUE, 3);
for (i = 0; i < nbins; i++)
d3->SetBinContent(i + 1, dest[i]);
d3->SetLineColor(kGreen);
d3->Draw("SAME");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
s->SearchHighRes(source, dest, nbins, 8, 2, kTRUE, 1000, kTRUE, 3);
for (i = 0; i < nbins; i++)
d4->SetBinContent(i + 1, dest[i]);
d4->SetLineColor(kMagenta);
d4->Draw("SAME");
printf("Found %d candidate peaks\n", nfound);
}

Definition at line 2118 of file TSpectrum.cxx.

◆ SetAverageWindow()

void TSpectrum::SetAverageWindow ( Int_t w = 3)
static

Static function: Set average window of searched peaks (see TSpectrum::SearchHighRes).

Definition at line 91 of file TSpectrum.cxx.

◆ SetBit() [1/2]

void TObject::SetBit ( UInt_t f)
inlineinherited

Definition at line 202 of file TObject.h.

◆ SetBit() [2/2]

void TObject::SetBit ( UInt_t f,
Bool_t set )
inherited

Set or unset the user status bits as specified in f.

Definition at line 888 of file TObject.cxx.

◆ SetDeconIterations()

void TSpectrum::SetDeconIterations ( Int_t n = 3)
static

Static function: Set max number of decon iterations in deconvolution operation (see TSpectrum::SearchHighRes).

Definition at line 100 of file TSpectrum.cxx.

◆ SetDrawOption()

void TObject::SetDrawOption ( Option_t * option = "")
virtualinherited

Set drawing option for object.

This option only affects the drawing style and is stored in the option field of the TObjOptLink supporting a TPad's primitive list (TList). Note that it does not make sense to call object.SetDrawOption(option) before having called object.Draw().

Reimplemented in TSystemDirectory, TSystemFile, TPad, TGFrame, TAxis, TBrowser, TPaveStats, TGedFrame, TRootBrowserLite, and RooPlot.

Definition at line 871 of file TObject.cxx.

◆ SetDtorOnly()

void TObject::SetDtorOnly ( void * obj)
staticinherited

Set destructor only flag.

Definition at line 1204 of file TObject.cxx.

◆ SetName()

void TNamed::SetName ( const char * name)
virtualinherited

Set the name of the TNamed.

WARNING: if the object is a member of a THashTable or THashList container the container must be Rehash()'ed after SetName(). For example the list of objects in the current directory is a THashList.

Reimplemented in TEveScene, TColor, TSystemDirectory, TSystemFile, TNode, TRotMatrix, TShape, TEfficiency, TFormula, TGraph2D, TH1, RooAbsArg, RooAbsData, RooDataHist, RooDataSet, RooFitResult, RooPlot, ROOT::Experimental::XRooFit::xRooNode, TChain, TEventList, TTree, TGraph, and TDirectory.

Definition at line 149 of file TNamed.cxx.

◆ SetNameTitle()

void TNamed::SetNameTitle ( const char * name,
const char * title )
virtualinherited

Set all the TNamed parameters (name and title).

WARNING: if the name is changed and the object is a member of a THashTable or THashList container the container must be Rehash()'ed after SetName(). For example the list of objects in the current directory is a THashList.

Reimplemented in TContextMenu, TNode, TGraph2D, TH1, RooAbsArg, RooAbsData, RooDataHist, RooDataSet, RooFitResult, RooPlot, and TGraph.

Definition at line 163 of file TNamed.cxx.

◆ SetObjectStat()

void TObject::SetObjectStat ( Bool_t stat)
staticinherited

Turn on/off tracking of objects in the TObjectTable.

Definition at line 1188 of file TObject.cxx.

◆ SetResolution()

void TSpectrum::SetResolution ( Double_t resolution = 1)

NOT USED resolution: determines resolution of the neighbouring peaks default value is 1 correspond to 3 sigma distance between peaks.

Higher values allow higher resolution (smaller distance between peaks. May be set later through SetResolution.

Definition at line 343 of file TSpectrum.cxx.

◆ SetTitle()

void TNamed::SetTitle ( const char * title = "")
virtualinherited

◆ SetUniqueID()

void TObject::SetUniqueID ( UInt_t uid)
virtualinherited

Set the unique object id.

Definition at line 899 of file TObject.cxx.

◆ Sizeof()

Int_t TNamed::Sizeof ( ) const
virtualinherited

Return size of the TNamed part of the TObject.

Reimplemented in TSQLFile, TXMLFile, TDirectory, TDirectoryFile, TFile, and TKey.

Definition at line 182 of file TNamed.cxx.

◆ SmoothMarkov()

const char * TSpectrum::SmoothMarkov ( Double_t * source,
Int_t ssize,
Int_t averWindow )

One-dimensional markov spectrum smoothing function.

This function calculates smoothed spectrum from source spectrum based on Markov chain method. The result is placed in the array pointed by source pointer. On successful completion it returns 0. On error it returns pointer to the string describing error.

Parameters:

  • source: pointer to the array of source spectrum
  • ssize: length of source array
  • averWindow: width of averaging smoothing window

The goal of this function is the suppression of the statistical fluctuations. The algorithm is based on discrete Markov chain, which has very simple invariant distribution:

\[ U_2 = \frac{p_{1,2}}{p_{2,1}}U_1, U_3 = \frac{p_{2,3}}{p_{3,2}}U_2U_1, ... , U_n = \frac{p_{n-1,n}}{p_{n,n-1}}U_{n-1}...U_2U_1 \]

\( U_1\) being defined from the normalization condition \( \sum_{i=1}^{n} U_i=1\). \( n \) is the length of the smoothed spectrum and

\[ p_{i,i\pm 1} = A_i\sum_{k=1}^{m} exp\left[ \frac{y(i\pm k)-y(i)}{y(i\pm k)+y(i)}\right] \]

Reference:

  1. Z.K. Silagadze, A new algorithm for automatic photopeak searches. NIM A 376 (1996), 451.

Example 14 - script Smoothing.C

void Smoothing()
{
Int_t i;
const Int_t nbins = 1024;
Double_t xmax = nbins;
Double_t source[nbins];
gROOT->ForceStyle();
TString dir = gROOT->GetTutorialDir();
TString file = dir + "/legacy/spectrum/TSpectrum.root";
TFile *f = new TFile(file.Data());
TH1F *h = (TH1F *)f->Get("back1");
h->SetTitle("Smoothed spectrum for m=3");
for (i = 0; i < nbins; i++)
source[i] = h->GetBinContent(i + 1);
h->SetAxisRange(1, 1024);
h->Draw("L");
TSpectrum *s = new TSpectrum();
TH1F *smooth = new TH1F("smooth", "smooth", nbins, 0., nbins);
smooth->SetLineColor(kRed);
s->SmoothMarkov(source, 1024, 3); // 3, 7, 10
for (i = 0; i < nbins; i++)
smooth->SetBinContent(i + 1, source[i]);
smooth->Draw("L SAME");
}
const char * SmoothMarkov(Double_t *source, Int_t ssize, Int_t averWindow)
One-dimensional markov spectrum smoothing function.

Definition at line 1187 of file TSpectrum.cxx.

◆ StaticBackground()

TH1 * TSpectrum::StaticBackground ( const TH1 * hist,
Int_t niter = 20,
Option_t * option = "" )
static

Static function, interface to TSpectrum::Background.

Definition at line 2579 of file TSpectrum.cxx.

◆ StaticSearch()

Int_t TSpectrum::StaticSearch ( const TH1 * hist,
Double_t sigma = 2,
Option_t * option = "goff",
Double_t threshold = 0.05 )
static

Static function, interface to TSpectrum::Search.

Definition at line 2569 of file TSpectrum.cxx.

◆ Streamer()

void TSpectrum::Streamer ( TBuffer & R__b)
overridevirtual

Stream an object of class TObject.

Reimplemented from TNamed.

◆ StreamerNVirtual()

void TSpectrum::StreamerNVirtual ( TBuffer & ClassDef_StreamerNVirtual_b)
inline

Definition at line 78 of file TSpectrum.h.

◆ SysError()

void TObject::SysError ( const char * location,
const char * fmt,
... ) const
virtualinherited

Issue system error message.

Use "location" to specify the method where the system error occurred. Accepts standard printf formatting arguments.

Definition at line 1112 of file TObject.cxx.

◆ TestBit()

R__ALWAYS_INLINE Bool_t TObject::TestBit ( UInt_t f) const
inlineinherited

Definition at line 204 of file TObject.h.

◆ TestBits()

Int_t TObject::TestBits ( UInt_t f) const
inlineinherited

Definition at line 205 of file TObject.h.

◆ Unfolding()

const char * TSpectrum::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.

This function unfolds source spectrum according to response matrix columns. The result is placed in the vector pointed by source pointer. The coefficients of the resulting vector represent contents of the columns (weights) in the input vector. On successful completion it returns 0. On error it returns pointer to the string describing error. If desired after every numberIterations one can apply boosting operation (exponential function with exponent given by boost coefficient) and repeat it numberRepetitions times. For details we refer to [1].

Parameters:

  • source: pointer to the vector of source spectrum
  • respMatrix: pointer to the matrix of response spectra
  • ssizex: length of source spectrum and # of rows of the response matrix. ssizex must be >= ssizey.
  • ssizey: length of destination coefficients and # of columns of the response matrix.
  • numberIterations: number of iterations
  • numberRepetitions: number of repetitions for boosted deconvolution. It must be greater or equal to one.
  • boost: boosting coefficient, applies only if numberRepetitions is greater than one.

Unfolding:

The goal is the decomposition of spectrum to a given set of component spectra.

The mathematical formulation of the discrete linear system is:

\[ y(i) = \sum_{k=0}^{N_y-1} h(i,k)x(k), i = 0,1,2,...,N_x-1 \]

\[ \begin{bmatrix} y(0) \\ y(1) \\ \dots \\ y(N_x-1) \end{bmatrix} = \begin{bmatrix} h(0,0) & h(0,1) & \dots & h(0,N_y-1) \\ h(1,0) & h(1,1) & \dots & h(1,N_y-1) \\ \dots \\ h(N_x-1,0) & h(N_x-1,1) & \dots & h(N_x-1,N_y-1) \end{bmatrix} \begin{bmatrix} x(0) \\ x(1) \\ \dots \\ x(N_y-1) \end{bmatrix} \]

References:

  1. Jandel M., Morhac; M., Kliman J., Krupa L., Matouoek V., Hamilton J. H., Ramaya A. V.: Decomposition of continuum gamma-ray spectra using synthesised response matrix. NIM A 516 (2004), 172-183.

Example of unfolding:

Example 13 - script Unfolding.C:

Fig. 20 Response matrix composed of neutron spectra of pure chemical elements.
Fig. 21 Source neutron spectrum to be decomposed
Fig. 22 Spectrum after decomposition, contains 10 coefficients, which correspond to contents of chemical components (dominant 8-th and 10-th components, i.e. O, Si)

Script:

// Example to illustrate unfolding function (class TSpectrum).
// To execute this example, do
// root > .x Unfolding.C
void Unfolding() {
Int_t i, j;
Int_t nbinsx = 2048;
Int_t nbinsy = 10;
double xmin = 0;
double xmax = nbinsx;
double ymin = 0;
double ymax = nbinsy;
double *source = new double[nbinsx];
double **response = new double *[nbinsy];
for (i=0;i<nbinsy;i++) response[i] = new double[nbinsx];
TH1F *h = new TH1F("h","",nbinsx,xmin,xmax);
TH1F *d = new TH1F("d","Decomposition - unfolding",nbinsx,xmin,xmax);
TH2F *decon_unf_resp = new TH2F("decon_unf_resp","Root File",nbinsy,ymin,ymax,nbinsx,xmin,xmax);
TFile *f = new TFile("TSpectrum.root");
h = (TH1F*) f->Get("decon_unf_in;1");
TFile *fr = new TFile("TSpectrum.root");
decon_unf_resp = (TH2F*) fr->Get("decon_unf_resp;1");
for (i = 0; i < nbinsx; i++) source[i] = h->GetBinContent(i + 1);
for (i = 0; i < nbinsy; i++){
for (j = 0; j< nbinsx; j++){
response[i][j] = decon_unf_resp->GetBinContent(i + 1, j + 1);
}
}
TCanvas *Decon1 = (TCanvas *)gROOT->GetListOfCanvases()->FindObject("Decon1");
if (!Decon1) Decon1 = new TCanvas("Decon1","Decon1",10,10,1000,700);
h->Draw("L");
TSpectrum *s = new TSpectrum();
s->Unfolding(source,(const double**)response,nbinsx,nbinsy,1000,1,1);
for (i = 0; i < nbinsy; i++) d->SetBinContent(i + 1,source[i]);
d->SetLineColor(kRed);
d->SetAxisRange(0,nbinsy);
d->Draw("");
}
float ymin
float ymax
The Canvas class.
Definition TCanvas.h:23
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:345
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.
TSpectrum()
Constructor.
Definition TSpectrum.cxx:42

Definition at line 1881 of file TSpectrum.cxx.

◆ UseCurrentStyle()

void TObject::UseCurrentStyle ( )
virtualinherited

Set current style settings in this object This function is called when either TCanvas::UseCurrentStyle or TROOT::ForceStyle have been invoked.

Reimplemented in TCanvas, TPad, TFrame, TPaveStats, TPaveText, TAxis3D, TGraph, TH1, and TTree.

Definition at line 909 of file TObject.cxx.

◆ Warning()

void TObject::Warning ( const char * location,
const char * fmt,
... ) const
virtualinherited

Issue warning message.

Use "location" to specify the method where the warning occurred. Accepts standard printf formatting arguments.

Definition at line 1084 of file TObject.cxx.

◆ Write() [1/2]

Int_t TObject::Write ( const char * name = nullptr,
Int_t option = 0,
Int_t bufsize = 0 )
virtualinherited

Write this object to the current directory.

For more see the const version of this method.

Reimplemented in TSQLFile, TXMLFile, TDirectory, TBuffer, ROOT::TBufferMergerFile, TDirectoryFile, TFile, TParallelMergingFile, TCollection, TMap, and TTree.

Definition at line 989 of file TObject.cxx.

◆ Write() [2/2]

Int_t TObject::Write ( const char * name = nullptr,
Int_t option = 0,
Int_t bufsize = 0 ) const
virtualinherited

Write this object to the current directory.

The data structure corresponding to this object is serialized. The corresponding buffer is written to the current directory with an associated key with name "name".

Writing an object to a file involves the following steps:

  • Creation of a support TKey object in the current directory. The TKey object creates a TBuffer object.
  • The TBuffer object is filled via the class::Streamer function.
  • If the file is compressed (default) a second buffer is created to hold the compressed buffer.
  • Reservation of the corresponding space in the file by looking in the TFree list of free blocks of the file.
  • The buffer is written to the file.

Bufsize can be given to force a given buffer size to write this object. By default, the buffersize will be taken from the average buffer size of all objects written to the current file so far.

If a name is specified, it will be the name of the key. If name is not given, the name of the key will be the name as returned by GetName().

The option can be a combination of: kSingleKey, kOverwrite or kWriteDelete Using the kOverwrite option a previous key with the same name is overwritten. The previous key is deleted before writing the new object. Using the kWriteDelete option a previous key with the same name is deleted only after the new object has been written. This option is safer than kOverwrite but it is slower. NOTE: Neither kOverwrite nor kWriteDelete reduces the size of a TFile– the space is simply freed up to be overwritten; in the case of a TTree, it is more complicated. If one opens a TTree, appends some entries, then writes it out, the behaviour is effectively the same. If, however, one creates a new TTree and writes it out in this way, only the metadata is replaced, effectively making the old data invisible without deleting it. TTree::Delete() can be used to mark all disk space occupied by a TTree as free before overwriting its metadata this way. The kSingleKey option is only used by TCollection::Write() to write a container with a single key instead of each object in the container with its own key.

An object is read from the file into memory via TKey::Read() or via TObject::Read().

The function returns the total number of bytes written to the file. It returns 0 if the object cannot be written.

Reimplemented in TSQLFile, TXMLFile, TDirectory, TBuffer, TDirectoryFile, TFile, TParallelMergingFile, TCollection, TMap, and TTree.

Definition at line 964 of file TObject.cxx.

Member Data Documentation

◆ fBits

UInt_t TObject::fBits
privateinherited

bit field status word

Definition at line 47 of file TObject.h.

◆ fgAverageWindow

Int_t TSpectrum::fgAverageWindow = 3
staticprotected

Average window of searched peaks.

Definition at line 32 of file TSpectrum.h.

◆ fgDtorOnly

Longptr_t TObject::fgDtorOnly = 0
staticprivateinherited

object for which to call dtor only (i.e. no delete)

Definition at line 49 of file TObject.h.

◆ fgIterations

Int_t TSpectrum::fgIterations = 3
staticprotected

Maximum number of decon iterations (default=3)

Definition at line 33 of file TSpectrum.h.

◆ fgObjectStat

Bool_t TObject::fgObjectStat = kTRUE
staticprivateinherited

if true keep track of objects in TObjectTable

Definition at line 50 of file TObject.h.

◆ fHistogram

TH1* TSpectrum::fHistogram
protected

resulting histogram

Definition at line 31 of file TSpectrum.h.

◆ fMaxPeaks

Int_t TSpectrum::fMaxPeaks
protected

Maximum number of peaks to be found.

Definition at line 25 of file TSpectrum.h.

◆ fName

TString TNamed::fName
protectedinherited

Definition at line 32 of file TNamed.h.

◆ fNPeaks

Int_t TSpectrum::fNPeaks
protected

number of peaks found

Definition at line 26 of file TSpectrum.h.

◆ fPosition

Double_t* TSpectrum::fPosition
protected

[fNPeaks] array of current peak positions

Definition at line 27 of file TSpectrum.h.

◆ fPositionX

Double_t* TSpectrum::fPositionX
protected

[fNPeaks] X position of peaks

Definition at line 28 of file TSpectrum.h.

◆ fPositionY

Double_t* TSpectrum::fPositionY
protected

[fNPeaks] Y position of peaks

Definition at line 29 of file TSpectrum.h.

◆ fResolution

Double_t TSpectrum::fResolution
protected

NOT USED resolution of the neighboring peaks

Definition at line 30 of file TSpectrum.h.

◆ fTitle

TString TNamed::fTitle
protectedinherited

Definition at line 33 of file TNamed.h.

◆ fUniqueID

UInt_t TObject::fUniqueID
privateinherited

object unique identifier

Definition at line 46 of file TObject.h.


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