library: libSpectrum
#include "TSpectrumTransform.h"

TSpectrumTransform


class description - header file - source file
viewCVS header - viewCVS source

class TSpectrumTransform: public TNamed

Inheritance Inherited Members Includes Libraries
Class Charts

Function Members (Methods)

Display options:
Show inherited
Show non-public
public:
TSpectrumTransform()
TSpectrumTransform(Int_t size)
TSpectrumTransform(const TSpectrumTransform&)
virtual~TSpectrumTransform()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
virtual TObject*TNamed::Clone(const char* newname = "") const
virtual Int_tTNamed::Compare(const TObject* obj) const
virtual voidTNamed::Copy(TObject& named) const
virtual voidTObject::Delete(Option_t* option = "")
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() const
virtual TObject*TObject::DrawClone(Option_t* option = "") const
virtual voidTObject::Dump() const
voidEnhance(const float* source, float* destVector)
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual voidTNamed::FillBuffer(char*& buffer)
voidFilterZonal(const float* source, float* destVector)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
virtual const char*TNamed::GetName() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
virtual const char*TNamed::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() const
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTNamed::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTNamed::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
virtual Bool_tTObject::Notify()
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
TSpectrumTransform&operator=(const TSpectrumTransform&)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidTNamed::Print(Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") const
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
voidSetDirection(Int_t direction)
virtual voidTObject::SetDrawOption(Option_t* option = "")
static voidTObject::SetDtorOnly(void* obj)
voidSetEnhanceCoeff(Float_t enhanceCoeff)
voidSetFilterCoeff(Float_t filterCoeff)
virtual voidTNamed::SetName(const char* name)
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
static voidTObject::SetObjectStat(Bool_t stat)
voidSetRegion(Int_t xmin, Int_t xmax)
virtual voidTNamed::SetTitle(const char* title = "")
voidSetTransformType(Int_t transType, Int_t degree)
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual Int_tTNamed::Sizeof() const
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
voidTransform(const float* source, float* destVector)
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = "0", Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = "0", Int_t option = 0, Int_t bufsize = 0) const
protected:
voidBitReverse(float* working_space, Int_t num)
voidBitReverseHaar(float* working_space, Int_t shift, Int_t num, Int_t start)
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
voidFourier(float* working_space, Int_t num, Int_t hartley, Int_t direction, Int_t zt_clear)
Int_tGeneralExe(float* working_space, Int_t zt_clear, Int_t num, Int_t degree, Int_t type)
Int_tGeneralInv(float* working_space, Int_t num, Int_t degree, Int_t type)
voidHaar(float* working_space, Int_t num, Int_t direction)
voidTObject::MakeZombie()
voidWalsh(float* working_space, Int_t num)

Data Members

public:
enum { kTransformHaar
kTransformWalsh
kTransformCos
kTransformSin
kTransformFourier
kTransformHartley
kTransformFourierWalsh
kTransformFourierHaar
kTransformWalshHaar
kTransformCosWalsh
kTransformCosHaar
kTransformSinWalsh
kTransformSinHaar
kTransformForward
kTransformInverse
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
protected:
Int_tfSizelength of transformed data
Int_tfTransformTypetype of transformation (Haar, Walsh, Cosine, Sine, Fourier, Hartley, Fourier-Walsh, Fourier-Haar, Walsh-Haar, Cosine-Walsh, Cosine-Haar, Sine-Walsh, Sine-Haar)
Int_tfDegreedegree of mixed transform, applies only for Fourier-Walsh, Fourier-Haar, Walsh-Haar, Cosine-Walsh, Cosine-Haar, Sine-Walsh, Sine-Haar transforms
Int_tfDirectionforward or inverse transform
Int_tfXminfirst channel of filtered or enhanced region
Int_tfXmaxlast channel of filtered or enhanced region
Float_tfFilterCoeffvalue set in the filtered region
Float_tfEnhanceCoeffmultiplication coefficient applied in enhanced region;
TStringTNamed::fNameobject identifier
TStringTNamed::fTitleobject title

Class Description

   THIS CLASS CONTAINS ORTHOGONAL TRANSFORM  FUNCTIONS.                  
                                                                         
   These functions were written by:                                      
   Miroslav Morhac                                                       
   Institute of Physics                                                  
   Slovak Academy of Sciences                                            
   Dubravska cesta 9, 842 28 BRATISLAVA                                  
   SLOVAKIA                                                              
                                                                         
   email:fyzimiro@savba.sk,    fax:+421 7 54772479                       
                                                                         
  The original code in C has been repackaged as a C++ class by R.Brun    
                                                                         
  The algorithms in this class have been published in the following      
  references:                                                            
                                                                         
  [1] C.V. Hampton, B. Lian, Wm. C. McHarris: Fast-Fourier-transform     
      spectral enhancement techniques for gamma-ray spectroscopy.NIM A353
      (1994) 280-284.                                                    
  [2] Morhac M., Matousek V., New adaptive Cosine-Walsh  transform and   
      its application to nuclear data compression, IEEE Transactions on  
      Signal Processing 48 (2000) 2693.                                  //  
  [3] Morhac M., Matousek V., Data compression using new fast adaptive   
      Cosine-Haar transforms, Digital Signal Processing 8 (1998) 63.     
  [4] Morhac M., Matousek V.: Multidimensional nuclear data compression  
      using fast adaptive Walsh-Haar transform. Acta Physica Slovaca 51  
     (2001) 307.                                                         
____________________________________________________________________________
TSpectrumTransform()
default constructor
TSpectrumTransform(Int_t size)
the constructor creates TSpectrumTransform object. Its size must be > than zero and must be power of 2.
It sets default transform type to be Cosine transform. Transform parameters can be changed using setter functions.
~TSpectrumTransform()
destructor
void Haar(float *working_space, int num, int direction)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion calculates Haar transform of a part of data                   
      Function parameters:                                                    
              -working_space-pointer to vector of transformed data            
              -num-length of processed data                                   
              -direction-forward or inverse transform                         
                                                                              

void Walsh(float *working_space, int num)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion calculates Walsh transform of a part of data                  
      Function parameters:                                                    
              -working_space-pointer to vector of transformed data            
              -num-length of processed data                                   
                                                                              

void BitReverse(float *working_space, int num)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion carries out bir-reverse reordering of data                    
      Function parameters:                                                    
              -working_space-pointer to vector of processed data              
              -num-length of processed data                                   
                                                                              

void Fourier(float *working_space, int num, int hartley, int direction, int zt_clear)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion calculates Fourier based transform of a part of data          
      Function parameters:                                                    
              -working_space-pointer to vector of transformed data            
              -num-length of processed data                                   
              -hartley-1 if it is Hartley transform, 0 othewise               
              -direction-forward or inverse transform                         
                                                                              

void BitReverseHaar(float *working_space, int shift, int num, int start)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion carries out bir-reverse reordering for Haar transform         
      Function parameters:                                                    
              -working_space-pointer to vector of processed data              
              -shift-shift of position of processing                          
              -start-initial position of processed data                       
              -num-length of processed data                                   
                                                                              

int GeneralExe(float *working_space, int zt_clear, int num, int degree, int type)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion calculates generalized (mixed) transforms of different degrees
      Function parameters:                                                    
              -working_space-pointer to vector of transformed data            
              -zt_clear-flag to clear imaginary data before staring           
              -num-length of processed data                                   
              -degree-degree of transform (see manual)                        
              -type-type of mixed transform (see manual)                      
                                                                              

int GeneralInv(float *working_space, int num, int degree, int type)
   AUXILIARY FUNCION                                                          
                                                                              
   This funcion calculates inverse generalized (mixed) transforms             
      Function parameters:                                                    
              -working_space-pointer to vector of transformed data            
              -num-length of processed data                                   
              -degree-degree of transform (see manual)                        
              -type-type of mixed transform (see manual)                      
                                                                              

void Transform(const float *source, float *destVector)
        ONE-DIMENSIONAL TRANSFORM FUNCTION                                 
        This function transforms the source spectrum. The calling program 
        should fill in input parameters.                                    
        Transformed data are written into dest spectrum.                  
                                                                           
        Function parameters:                                               
        source-pointer to the vector of source spectrum, its length should 
             be size except for inverse FOURIER, FOUR-WALSH, FOUR-HAAR     
             transform. These need 2*size length to supply real and        
             imaginary coefficients.                                       
        destVector-pointer to the vector of dest data, its length should be
             size except for direct FOURIER, FOUR-WALSH, FOUR-HAAR. These    
             need 2*size length to store real and imaginary coefficients    
                                                                         

Transform methods

 

Goal: to analyze experimental data using orthogonal transforms

         orthogonal transforms can be successfully used for the processing of nuclear spectra (not only)

         they can be used to remove high frequency noise, to increase signal-to-background ratio as well as to enhance low intensity components [1], to carry out e.g. Fourier analysis etc.

         we have implemented the function for the calculation of the commonly used orthogonal transforms as well as functions for the filtration and enhancement of experimental data

 

Function:

void TSpectrumTransform::Transform(const float *fSource, float *fDest)

 

This function transforms the source spectrum according to the given input parameters. Transformed data are written into dest spectrum. Before the Transform function is called the class must be created by constructor and the type of the transform as well as some other parameters must be set using a set of setter functions.

 

Member variables of TSpectrumTransform class:

 fSource-pointer to the vector of source spectrum. Its length should be equal to the “fSize” parameter except for inverse FOURIER, FOUR-WALSH, FOUR-HAAR transforms. These need “2*fSize” length to supply real and imaginary coefficients.                   

fDest-pointer to the vector of destination spectrum. Its length should be equal to the “fSize” parameter except for inverse FOURIER, FOUR-WALSH, FOUR-HAAR transforms. These need “2*fSize” length to store real and imaginary coefficients.

        fSize-basic length of the source and dest spectrum. It should be power of 2.

fType-type of transform

            Classic transforms:

                        kTransformHaar

                        kTransformWalsh

                        kTransformCos

                        kTransformSin

                        kTransformFourier

                        kTransformHartey

            Mixed transforms:

                        kTransformFourierWalsh

                        kTransformFourierHaar

                        kTransformWalshHaar

                        kTransformCosWalsh

                        kTransformCosHaar

                        kTransformSinWalsh

                        kTransformSinHaar

fDirection-direction-transform direction (forward, inverse)

                        kTransformForward

                        kTransformInverse

fDegree-applies only for mixed transforms [2], [3], [4].

                  Allowed range  .

References:

[1] C.V. Hampton, B. Lian, Wm. C. McHarris: Fast-Fourier-transform spectral enhancement techniques for gamma-ray spectroscopy. NIM A353 (1994) 280-284.

[2] Morháč M., Matoušek V., New adaptive Cosine-Walsh  transform and its application to nuclear data compression, IEEE Transactions on Signal Processing 48 (2000) 2693. 

[3] Morháč M., Matoušek V., Data compression using new fast adaptive Cosine-Haar transforms, Digital Signal Processing 8 (1998) 63.

[4] Morháč M., Matoušek V.: Multidimensional nuclear data compression using fast adaptive Walsh-Haar transform. Acta Physica Slovaca 51 (2001) 307.

 

Example  – script Transform.c:

Fig. 1 Original gamma-ray spectrum

 

Fig. 2 Transformed spectrum from Fig. 1 using Cosine transform

 

Script:

// Example to illustrate Transform function (class TSpectrumTransform).

// To execute this example, do

// root > .x Transform.C

  

#include <TSpectrum>

#include <TSpectrumTransform>

 

void Transform() {

   Int_t i;

   Double_t nbins = 4096;

   Double_t xmin  = 0;

   Double_t xmax  = (Double_t)nbins;

   Float_t * source = new float[nbins];

   Float_t * dest = new float[nbins];  

   TH1F *h = new TH1F("h","Transformed spectrum using Cosine transform",nbins,xmin,xmax);

   TFile *f = new TFile("spectra\\TSpectrum.root");

   h=(TH1F*) f->Get("transform1;1");  

   for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);        

   TCanvas *Transform1 = gROOT->GetListOfCanvases()->FindObject("Transform1");

   if (!Transform1) Transform1 = new TCanvas("Transform","Transform1",10,10,1000,700);

   TSpectrum *s = new TSpectrum();

   TSpectrumTransform *t = new TSpectrumTransform(4096);

   t->SetTransformType(t->kTransformCos,0);

   t->SetDirection(t->kTransformForward);

   t->Transform(source,dest);

   for (i = 0; i < nbins; i++) h->SetBinContent(i + 1,dest[i]);  

   h->SetLineColor(kRed);     

   h->Draw("L");

}


void FilterZonal(const float *source, float *destVector)
        ONE-DIMENSIONAL FILTER ZONAL FUNCTION                               
        This function transforms the source spectrum. The calling program  
        should fill in input parameters. Then it sets transformed          
        coefficients in the given region (fXmin, fXmax) to the given         
        fFilterCoeff and transforms it back.
        Filtered data are written into dest spectrum.                     
                                                                           
        Function parameters:                                               
        source-pointer to the vector of source spectrum, its length should 
             be size except for inverse FOURIER, FOUR-WALSH, FOUR-HAAR    
             transform. These need 2*size length to supply real and       
             imaginary coefficients.                                      
        destVector-pointer to the vector of dest data, its length should be
           size except for direct FOURIER, FOUR-WALSH, FOUR-HAAR. These  
           need 2*size length to store real and imaginary coefficients   
                                                                          

       

Example of zonal filtering

 

Function:

void TSpectrumTransform::FilterZonal(const float *fSource, float *fDest)

 

This function transforms the source spectrum (for details see Transform function). Before the FilterZonal function is called the class must be created by constructor and the type of the transform as well as some other parameters must be set using a set of setter funcions. The FilterZonal function sets transformed coefficients in the given region (fXmin, fXmax) to the given fFilterCoeff and transforms it back. Filtered data are written into dest spectrum.

 

Example – script Filter.c:

Fig. 1 Original spectrum (black line) and filtered spectrum (red line) using Cosine transform and zonal filtration (channels 2048-4095 were set to 0)

 

Script:

// Example to illustrate FilterZonal function (class TSpectrumTransform).

// To execute this example, do

// root > .x Filter.C

  

#include <TSpectrum>

#include <TSpectrumTransform>

 

void Filter() {

   Int_t i;

   Double_t nbins = 4096;

   Double_t xmin  = 0;

   Double_t xmax  = (Double_t)nbins;

   Float_t * source = new float[nbins];

   Float_t * dest = new float[nbins];  

   TH1F *h = new TH1F("h","Zonal filtering using Cosine transform",nbins,xmin,xmax);

   TH1F *d = new TH1F("d","",nbins,xmin,xmax);        

   TFile *f = new TFile("spectra\\TSpectrum.root");

   h=(TH1F*) f->Get("transform1;1");  

   for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);    

   TCanvas *Transform1 = gROOT->GetListOfCanvases()->FindObject("Transform1");

   if (!Transform1) Transform1 = new TCanvas("Transform","Transform1",10,10,1000,700);

   h->SetAxisRange(700,1024);  

   h->Draw("L");  

   TSpectrum *s = new TSpectrum();

   TSpectrumTransform *t = new TSpectrumTransform(4096);

   t->SetTransformType(t->kTransformCos,0);

   t->SetRegion(2048, 4095);

   t->FilterZonal(source,dest);    

   for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,dest[i]);

   d->SetLineColor(kRed);  

   d->Draw("SAME L");

}


void Enhance(const float *source, float *destVector)
        ONE-DIMENSIONAL ENHANCE ZONAL FUNCTION                             
        This function transforms the source spectrum. The calling program  
      should fill in input parameters. Then it multiplies transformed      
      coefficients in the given region (fXmin, fXmax) by the given          
      fEnhanceCoeff and transforms it back                                   
        Processed data are written into dest spectrum.                      
                                                                            
        Function parameters:                                                
        source-pointer to the vector of source spectrum, its length should  
             be size except for inverse FOURIER, FOUR-WALSh, FOUR-HAAR      
             transform. These need 2*size length to supply real and         
             imaginary coefficients.                                        
        destVector-pointer to the vector of dest data, its length should be 
           size except for direct FOURIER, FOUR-WALSh, FOUR-HAAR. These     
           need 2*size length to store real and imaginary coefficients      
                                                                           

Example of enhancement

 

Function:

void TSpectrumTransform::Enhance(const float *fSource, float *fDest)

 

This function transforms the source spectrum (for details see Transform function). Before the Enhance function is called the class must be created by constructor and the type of the transform as well as some other parameters must be set using a set of setter funcions. The Enhance function multiplies transformed coefficients in the given region (fXmin, fXmax) by the given fEnhancCoeff and transforms it back. Enhanced data are written into dest spectrum.

 

Example  – script Enhance.c:

 

Fig. 1 Original spectrum (black line) and enhanced spectrum (red line) using Cosine transform (channels 0-1024 were multiplied by 2)

 

Script:

 

// Example to illustrate Enhance function (class TSpectrumTransform).

// To execute this example, do

// root > .x Enhance.C

 

void Enhance() {

   Int_t i;

   Double_t nbins = 4096;

   Double_t xmin  = 0;

   Double_t xmax  = (Double_t)nbins;

   Float_t * source = new float[nbins];

   Float_t * dest = new float[nbins];  

   TH1F *h = new TH1F("h","Enhancement using Cosine transform",nbins,xmin,xmax);

   TH1F *d = new TH1F("d","",nbins,xmin,xmax);        

   TFile *f = new TFile("spectra\\TSpectrum.root");

   h=(TH1F*) f->Get("transform1;1");  

   for (i = 0; i < nbins; i++) source[i]=h->GetBinContent(i + 1);    

   TCanvas *Transform1 = gROOT->GetListOfCanvases()->FindObject("Transform1");

   if (!Transform1) Transform1 = new TCanvas("Transform","Transform1",10,10,1000,700);

   h->SetAxisRange(700,1024);  

   h->Draw("L");  

   TSpectrum *s = new TSpectrum();

   TSpectrumTransform *t = new TSpectrumTransform(4096);

   t->SetTransformType(t->kTransformCos,0);

   t->SetRegion(0, 1024);

   t->SetEnhanceCoeff(2);

   t->Enhance(source,dest);       

   for (i = 0; i < nbins; i++) d->SetBinContent(i + 1,dest[i]);

   d->SetLineColor(kRed);  

   d->Draw("SAME L");

}


void SetTransformType(Int_t transType, Int_t degree)
   SETTER FUNCION                                                      
                                                     
   This funcion sets the following parameters for transform:
         -transType - type of transform (Haar, Walsh, Cosine, Sine, Fourier, Hartley, Fourier-Walsh, Fourier-Haar, Walsh-Haar, Cosine-Walsh, Cosine-Haar, Sine-Walsh, Sine-Haar)
         -degree - degree of mixed transform, applies only for Fourier-Walsh, Fourier-Haar, Walsh-Haar, Cosine-Walsh, Cosine-Haar, Sine-Walsh, Sine-Haar transforms
////////////////////////////////////////////////////////////////////////////      
void SetRegion(Int_t xmin, Int_t xmax)
   SETTER FUNCION                                                      
                                                     
   This funcion sets the filtering or enhancement region:
         -xmin, xmax
////////////////////////////////////////////////////////////////////////////         
void SetDirection(Int_t direction)
   SETTER FUNCION                                                      
                                                     
   This funcion sets the direction of the transform:
         -direction (forward or inverse)
////////////////////////////////////////////////////////////////////////////   
void SetFilterCoeff(Float_t filterCoeff)
   SETTER FUNCION                                                      
                                                     
   This funcion sets the filter coefficient:
         -filterCoeff - after the transform the filtered region (xmin, xmax) is replaced by this coefficient. Applies only for filtereng operation.
////////////////////////////////////////////////////////////////////////////   
void SetEnhanceCoeff(Float_t enhanceCoeff)
   SETTER FUNCION                                                      
                                                     
   This funcion sets the enhancement coefficient:
         -enhanceCoeff - after the transform the enhanced region (xmin, xmax) is multiplied by this coefficient. Applies only for enhancement operation.
////////////////////////////////////////////////////////////////////////////   
TSpectrumTransform()

Author: Miroslav Morhac 25/09/06
Last update: root/spectrum:$Name: $:$Id: TSpectrumTransform.cxx,v 1.4 2006/10/09 09:29:59 brun Exp $


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.