]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrDep/AliAnaPi0.h
AliAnaPi0: Add possibility to cut on opening angle of pairs, to reject background...
[u/mrichter/AliRoot.git] / PWG4 / PartCorrDep / AliAnaPi0.h
1 #ifndef ALIANAPI0_H
2 #define ALIANAPI0_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5 /* $Id: $ */
6
7 //_________________________________________________________________________
8 // Class to fill two-photon invariant mass hisograms 
9 // to be used to extract pi0 raw yield.
10 //
11 //-- Author: Dmitri Peressounko (RRC "KI")
12 //-- Adapted to PartCorr frame by Lamia Benhabib (SUBATECH)
13 //-- and Gustavo Conesa (INFN-Frascati)
14
15 //Root
16 class TList;
17 class TH3D ;
18 class TH2D ;
19
20 //Analysis
21 class AliAODEvent ;
22 class AliESDEvent ;
23 #include "AliAnaPartCorrBaseClass.h"
24
25 class AliAnaPi0 : public AliAnaPartCorrBaseClass {
26   
27   public: 
28   
29   AliAnaPi0() ; // default ctor
30   AliAnaPi0(const char *name) ; // default ctor
31   AliAnaPi0(const AliAnaPi0 & g) ; // cpy ctor
32   AliAnaPi0 & operator = (const AliAnaPi0 & api0) ;//cpy assignment
33   virtual ~AliAnaPi0() ;//virtual dtor
34   
35   TList * GetCreateOutputObjects(); 
36   
37   void Print(const Option_t * opt) const;
38   
39   //void Init();
40   void InitParameters();
41   
42   //void MakeAnalysisFillAOD() {;} //Not needed
43   void MakeAnalysisFillHistograms();
44   
45   //    void SetBadRunsList(){;} ;     //Set list of runs which can be used for this analysis
46   //To be defined in future.
47   
48   //void SetEtalonHisto(TH3D * h);//Provide etalon of binning for histograms
49   
50   //Setters for parameters of event buffers
51   void SetNCentrBin(Int_t n=5){fNCentrBin=n ;} //number of bins in centrality 
52   void SetNZvertBin(Int_t n=5){fNZvertBin=n ;} //number of bins for vertex position
53   void SetNRPBin(Int_t n=6)   {fNrpBin=n ;}    //number of bins in reaction plain
54   void SetNMaxEvMix(Int_t n=20){fNmaxMixEv=n ;}//Maximal number of events for mixing
55   
56   //Setters for event selection
57   void SetZvertexCut(Float_t zcut=40.){fZvtxCut=zcut ;} //cut on vertex position
58   
59   TString GetCalorimeter()   const {return fCalorimeter ; }
60   void SetCalorimeter(TString det)    {fCalorimeter = det ; }
61         
62   void Terminate(TList* outputList);
63   void ReadHistograms(TList * outputList); //Fill histograms with histograms in ouput list, needed in Terminate.
64         
65   Int_t GetModuleNumber(AliAODPWG4Particle * particle);
66   void SetNumberOfModules(Int_t nmod) {fNModules = nmod;}
67         
68   Int_t GetNPID()   const {return fNPID ; }
69   void SetNPID(Int_t n)    {fNPID = n ; }
70         
71   void SwitchOnAngleSelection()    {fUseAngleCut = kTRUE ; }
72   void SwitchOffAngleSelection()   {fUseAngleCut = kFALSE ; }
73
74   private:
75   Bool_t IsBadRun(Int_t /*iRun*/) const {return kFALSE;} //Tests if this run bad according to private list
76   
77   private:
78   Int_t    fNCentrBin ;   // Number of bins in event container for centrality
79   Int_t    fNZvertBin ;   // Number of bins in event container for vertex position
80   Int_t    fNrpBin ;      // Number of bins in event container for reaction plain
81   Int_t    fNPID ;                // Number of possible PID combinations
82   Int_t    fNmaxMixEv ;   // Maximal number of events stored in buffer for mixing
83   Float_t  fZvtxCut ;     // Cut on vertex position
84   TString  fCalorimeter ; // Select Calorimeter for IM
85   Int_t    fNModules ;    // Number of EMCAL/PHOS modules, set as many histogras as modules 
86   Bool_t   fUseAngleCut ; // Select pairs depending on their opening angle
87   TList ** fEventsList ;  //! Containers for photons in stored events
88   
89   //Histograms
90   
91   //TH3D * fhEtalon ; //Etalon histo, all distributions will have same binning as this one
92   
93   TH3D ** fhReMod ;  //!REAL two-photon invariant mass distribution for different calorimeter modules.
94         
95   TH3D ** fhRe1 ;  //!REAL two-photon invariant mass distribution for different centralities and PID 
96   TH3D ** fhMi1 ;  //!MIXED two-photon invariant mass distribution for different centralities and PID
97   TH3D ** fhRe2 ;  //!REAL two-photon invariant mass distribution for different centralities and PID 
98   TH3D ** fhMi2 ;  //!MIXED two-photon invariant mass distribution for different centralities and PID
99   TH3D ** fhRe3 ;  //!REAL two-photon invariant mass distribution for different centralities and PID 
100   TH3D ** fhMi3 ;  //!MIXED two-photon invariant mass distribution for different centralities and PID
101   TH3D * fhEvents; //!Number of events per centrality, RP, zbin
102
103   TH2D * fhRealOpeningAngle ;    //! Opening angle of pair versus pair energy
104   TH2D * fhRealCosOpeningAngle ; //! Cosinus of opening angle of pair version pair energy
105
106   //Acceptance
107   TH1D * fhPrimPt ;    //! Spectrum of Primary 
108   TH1D * fhPrimAccPt ; //! Spectrum of primary with accepted daughters 
109   TH1D * fhPrimY ;     //! Rapidity distribution of primary particles
110   TH1D * fhPrimAccY ;  //! Rapidity distribution of primary with accepted daughters
111   TH1D * fhPrimPhi ;   //! Azimutal distribution of primary particles
112   TH1D * fhPrimAccPhi; //! Azimutal distribution of primary with accepted daughters     
113   TH2D * fhPrimOpeningAngle ;    //! Opening angle of pair versus pair energy, primaries
114   TH2D * fhPrimCosOpeningAngle ; //! Cosinus of opening angle of pair version pair energy, primaries
115         
116   ClassDef(AliAnaPi0,8)
117 } ;
118
119
120 #endif //ALIANAPI0_H
121
122
123