]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/CaloTrackCorrelations/AliAnaPhotonConvInCalo.h
add option to fill or not histograms in eta gaps, and momentum imbalance in associate...
[u/mrichter/AliRoot.git] / PWGGA / CaloTrackCorrelations / AliAnaPhotonConvInCalo.h
1 #ifndef ALIANAPHOTONINCALO_H
2 #define ALIANAPHOTONINCALO_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5
6 //_________________________________________________________________________
7 //
8 // Conversions pairs analysis
9 // Check if cluster comes from a conversion in the material in front of the calorimeter
10 // Do invariant mass of all pairs, if mass is close to 0, then it is conversion.
11 // Input are selected clusters with AliAnaPhoton
12 //
13 //-- Author: Gustavo Conesa (LPSC-IN2P3-CNRS)
14
15 // --- ROOT system ---
16 class TH2F;
17 class TH1F;
18 class TString ;
19 class TObjString;
20 class TList ;
21
22 // --- ANALYSIS system ---
23 #include "AliAnaCaloTrackCorrBaseClass.h"
24
25 class AliAnaPhotonConvInCalo : public AliAnaCaloTrackCorrBaseClass {
26
27  public: 
28   AliAnaPhotonConvInCalo() ;              // default ctor
29   virtual ~AliAnaPhotonConvInCalo() { ; } // virtual dtor
30         
31   //---------------------------------------
32   // General analysis frame methods
33   //---------------------------------------
34   
35   TObjString * GetAnalysisCuts();
36   
37   TList      * GetCreateOutputObjects();
38   
39   void         InitParameters();
40
41   void         MakeAnalysisFillAOD()  ;
42
43   void         MakeAnalysisFillHistograms() ; 
44   
45   void         Print(const Option_t * opt)const ;
46     
47   //---------------------------------------
48   // Analysis parameters setters getters
49   //---------------------------------------
50     
51   Float_t      GetMassCut()                     const { return fMassCut            ; }
52   void         SetMassCut(Float_t m)                  { fMassCut    = m            ; }
53         
54   Bool_t       AreConvertedPairsInAOD()         const { return fAddConvertedPairsToAOD   ; }
55   void         SwitchOnAdditionConvertedPairsToAOD()  { fAddConvertedPairsToAOD = kTRUE  ; }
56   void         SwitchOffAdditionConvertedPairsToAOD() { fAddConvertedPairsToAOD = kFALSE ; }  
57         
58   Bool_t       AreConvertedPairsRemoved()       const { return fRemoveConvertedPair      ; }
59   void         SwitchOnConvertedPairsRemoval()        { fRemoveConvertedPair  = kTRUE    ; }
60   void         SwitchOffConvertedPairsRemoval()       { fRemoveConvertedPair  = kFALSE   ; }    
61   
62   void         SetConvAsymCut(Float_t c)              { fConvAsymCut = c           ; }
63   Float_t      GetConvAsymCut()                 const { return fConvAsymCut        ; }
64   
65   void         SetConvDEtaCut(Float_t c)              { fConvDEtaCut = c           ; }
66   Float_t      GetConvDEtaCut()                 const { return fConvDEtaCut        ; }
67   
68   void         SetConvDPhiCut(Float_t min, Float_t max)  { fConvDPhiMinCut = min   ;  
69                                                            fConvDPhiMaxCut = max   ; }
70   Float_t      GetConvDPhiMinCut()              const { return fConvDPhiMinCut     ; }
71   Float_t      GetConvDPhiMaxCut()              const { return fConvDPhiMaxCut     ; }
72   
73   private:
74  
75   Bool_t   fRemoveConvertedPair;          // Remove conversion pairs
76   Bool_t   fAddConvertedPairsToAOD;       // Put Converted pairs in AOD
77   Float_t  fMassCut;                      // Mass cut for the conversion pairs selection  
78   Float_t  fConvAsymCut;                  // Select conversion pairs when asymmetry is smaller than cut
79   Float_t  fConvDEtaCut;                  // Select conversion pairs when deta of pair smaller than cut
80   Float_t  fConvDPhiMinCut;               // Select conversion pairs when dphi of pair lager than cut
81   Float_t  fConvDPhiMaxCut;               // Select conversion pairs when dphi of pair smaller than cut
82
83   // Histograms
84   TH1F * fhPtPhotonConv   ;               //! Number of identified photon vs transerse momentum 
85   TH2F * fhEtaPhiPhotonConv  ;            //! Pseudorapidity vs Phi of identified  photon for transerse momentum > 0.5, for converted
86   TH2F * fhEtaPhi05PhotonConv  ;          //! Pseudorapidity vs Phi of identified  photon for transerse momentum < 0.5, for converted
87   TH2F * fhConvDeltaEta;                  //! Small mass photons, correlation in eta
88   TH2F * fhConvDeltaPhi;                  //! Small mass photons, correlation in phi
89   TH2F * fhConvDeltaEtaPhi;               //! Small mass photons, correlation in phi and eta
90   TH2F * fhConvAsym;                      //! Small mass photons, correlation in energy asymmetry
91   TH2F * fhConvPt;                        //! Small mass photons, pT of pair
92   
93   //Vertex distance
94   TH2F * fhConvDistEta;                   //! Approx distance to vertex vs cluster Eta 
95   TH2F * fhConvDistEn;                    //! Approx distance to vertex vs Energy
96   TH2F * fhConvDistMass;                  //! Approx distance to vertex vs Mass
97   TH2F * fhConvDistEtaCutEta;             //! Approx distance to vertex vs cluster Eta, dEta < 0.05 
98   TH2F * fhConvDistEnCutEta;              //! Approx distance to vertex vs Energy, dEta < 0.05
99   TH2F * fhConvDistMassCutEta;            //! Approx distance to vertex vs Mass, dEta < 0.05
100   TH2F * fhConvDistEtaCutMass;            //! Approx distance to vertex vs cluster Eta, dEta < 0.05, m < 10 MeV 
101   TH2F * fhConvDistEnCutMass;             //! Approx distance to vertex vs Energy, dEta < 0.05, m < 10 MeV
102   TH2F * fhConvDistEtaCutAsy;             //! Approx distance to vertex vs cluster Eta, dEta < 0.05, m < 10 MeV, A < 0.1
103   TH2F * fhConvDistEnCutAsy;              //! Approx distance to vertex vs energy, dEta < 0.05, m < 10 MeV, A < 0.1
104
105   //Conversion pairs analysis histograms
106   TH1F * fhPtConversionTagged;            //! Number of identified gamma from Conversion , tagged as conversion 
107   TH1F * fhPtAntiNeutronTagged;           //! Number of identified gamma from AntiNeutrons gamma, tagged as conversion 
108   TH1F * fhPtAntiProtonTagged;            //! Number of identified gamma from AntiProtons gamma, tagged as conversion 
109   TH1F * fhPtUnknownTagged;               //! Number of identified gamma from unknown, tagged as conversion 
110   
111   TH2F * fhConvDeltaEtaMCConversion;      //! Small mass cluster pairs, correlation in eta, origin of both clusters is conversion
112   TH2F * fhConvDeltaPhiMCConversion;      //! Small mass cluster pairs, correlation in phi, origin of both clusters is conversion
113   TH2F * fhConvDeltaEtaPhiMCConversion;   //! Small mass cluster pairs, correlation in eta-phi, origin of both clusters is conversion
114   TH2F * fhConvAsymMCConversion;          //! Small mass cluster pairs, correlation in energy asymmetry, origin of both clusters is conversion
115   TH2F * fhConvPtMCConversion;            //! Small mass cluster pairs, pt of pair, origin of both clusters is conversion
116   TH2F * fhConvDispersionMCConversion;    //! Small mass cluster pairs, dispersion of cluster 1 vs cluster 2
117   TH2F * fhConvM02MCConversion;           //! Small mass cluster pairs, m02 of cluster 1 vs cluster 2 
118
119   TH2F * fhConvDeltaEtaMCAntiNeutron;     //! Small mass cluster pairs, correlation in eta, origin of both clusters is anti neutron
120   TH2F * fhConvDeltaPhiMCAntiNeutron;     //! Small mass cluster pairs, correlation in phi, origin of both clusters is anti neutron
121   TH2F * fhConvDeltaEtaPhiMCAntiNeutron;  //! Small mass cluster pairs, correlation in eta-phi, origin of both clusters is anti neutron
122   TH2F * fhConvAsymMCAntiNeutron;         //! Small mass cluster pairs, correlation in energy asymmetry, origin of both clusters is anti neutron
123   TH2F * fhConvPtMCAntiNeutron;           //! Small mass cluster pairs, pt of pair, origin of both clusters is anti neutron
124   TH2F * fhConvDispersionMCAntiNeutron;   //! Small mass cluster pairs, dispersion of cluster 1 vs cluster 2, origin of both clusters is anti neutron
125   TH2F * fhConvM02MCAntiNeutron;          //! Small mass cluster pairs, m02 of cluster 1 vs cluster 2, origin of both clusters is anti neutron
126
127   TH2F * fhConvDeltaEtaMCAntiProton;      //! Small mass cluster pairs, correlation in eta, origin of both clusters is anti proton
128   TH2F * fhConvDeltaPhiMCAntiProton;      //! Small mass cluster pairs, correlation in phi, origin of both clusters is anti proton
129   TH2F * fhConvDeltaEtaPhiMCAntiProton;   //! Small mass cluster pairs, correlation in eta-phi, origin of both clusters is anti proton
130   TH2F * fhConvAsymMCAntiProton;          //! Small mass cluster pairs, correlation in energy asymmetry, origin of both clusters is anti proton
131   TH2F * fhConvPtMCAntiProton;            //! Small mass cluster pairs, pt of pairs, origin of both clusters is anti proton
132   TH2F * fhConvDispersionMCAntiProton;    //! Small mass cluster pairs, dispersion of cluster 1 vs cluster 2, origin of both clusters is anti proton
133   TH2F * fhConvM02MCAntiProton;           //! Small mass cluster pairs, m02 of cluster 1 vs cluster 2, origin of both clusters is anti proton
134
135   TH2F * fhConvDeltaEtaMCString;          //! Small mass cluster pairs, correlation in eta, origin of both clusters is string
136   TH2F * fhConvDeltaPhiMCString;          //! Small mass cluster pairs, correlation in phi, origin of both clusters is string
137   TH2F * fhConvDeltaEtaPhiMCString;       //! Small mass cluster pairs, correlation in eta-phi, origin of both clusters is string
138   TH2F * fhConvAsymMCString;              //! Small mass cluster pairs, correlation in energy asymmetry, origin of both clusters is string
139   TH2F * fhConvPtMCString;                //! Small mass cluster pairs, pt of pairs, origin of both clusters is string
140   TH2F * fhConvDispersionMCString;        //! Small mass cluster pairs, dispersion of cluster 1 vs cluster 2, origin of both clusters is string
141   TH2F * fhConvM02MCString;               //! Small mass cluster pairs, m02 of cluster 1 vs cluster 2, origin of both clusters is string
142   TH2F * fhConvDistMCConversion;          //! Calculated conversion distance vs real distance to vertex       
143   TH2F * fhConvDistMCConversionCuts;      //! Calculated conversion distance vs real distance to vertex        
144
145   AliAnaPhotonConvInCalo(              const AliAnaPhotonConvInCalo & g) ; // cpy ctor
146   AliAnaPhotonConvInCalo & operator = (const AliAnaPhotonConvInCalo & g) ; // cpy assignment
147   
148   ClassDef(AliAnaPhotonConvInCalo,1)
149
150 } ;
151  
152 #endif//ALIANAPHOTONINCALO_H
153
154
155