3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * See cxx source for full Copyright notice */
7 //_________________________________________________________________________
8 // Class to fill two-photon invariant mass histograms
9 // to be used to extract pi0 raw yield.
10 // Input is produced by AliAnaPhoton (or any other analysis producing output AliAODPWG4Particles),
11 // it will do nothing if executed alone
13 //-- Author: Dmitri Peressounko (RRC "KI")
14 //-- Adapted to PartCorr frame by Lamia Benhabib (SUBATECH)
15 //-- and Gustavo Conesa (INFN-Frascati)
24 #include "AliAnaPartCorrBaseClass.h"
27 class AliAODPWG4Particle ;
29 class AliAnaPi0 : public AliAnaPartCorrBaseClass {
32 AliAnaPi0() ; // default ctor
33 virtual ~AliAnaPi0() ;//virtual dtor
35 AliAnaPi0(const AliAnaPi0 & g) ; // cpy ctor
36 AliAnaPi0 & operator = (const AliAnaPi0 & api0) ;//cpy assignment
40 //-------------------------------
41 // General analysis frame methods
42 //-------------------------------
44 TObjString * GetAnalysisCuts();
45 TList * GetCreateOutputObjects();
46 void Terminate(TList* outputList);
47 void ReadHistograms(TList * outputList); //Fill histograms with histograms in ouput list, needed in Terminate.
48 void Print(const Option_t * opt) const;
49 //void MakeAnalysisFillAOD() {;} //Not needed
50 void MakeAnalysisFillHistograms();
52 void InitParameters();
55 TString GetCalorimeter() const { return fCalorimeter; }
56 void SetCalorimeter(TString & det) { fCalorimeter = det ; }
57 void SetNumberOfModules(Int_t nmod) { fNModules = nmod; }
59 //-------------------------------
61 //-------------------------------
63 virtual Int_t GetEventIndex(AliAODPWG4Particle * part, Double_t * vert) ;
65 void CountAndGetAverages(Int_t &nClus,Int_t &nCell, Float_t &eClusTot,Float_t &eCellTot, Float_t &eDenClus,Float_t &eDenCell) ;
67 //Switchs for event multiplicity bin option, by default, centrality
68 void SwitchOnTrackMultBins() {fUseTrackMultBins = kTRUE ; }
69 void SwitchOffTrackMultBins() {fUseTrackMultBins = kFALSE ; }
71 void SwitchOnPhotonMultBins() {fUsePhotonMultBins = kTRUE ; }
72 void SwitchOffPhotonMultBins() {fUsePhotonMultBins = kFALSE ; }
74 void SwitchOnClusterEBins() {fUseAverClusterEBins = kTRUE ; }
75 void SwitchOffClusterEBins() {fUseAverClusterEBins = kFALSE ; }
77 void SwitchOnCellEBins() {fUseAverCellEBins = kTRUE ; }
78 void SwitchOffCellEBins() {fUseAverCellEBins = kFALSE ; }
80 void SwitchOnClusterEDenBins() {fUseAverClusterEDenBins = kTRUE ; }
81 void SwitchOffClusterEDenBins() {fUseAverClusterEDenBins = kFALSE ; }
83 // void SwitchOnClusterPairRBins() {fUseAverClusterPairRBins = kTRUE ; }
84 // void SwitchOffClusterPairRBins() {fUseAverClusterPairRBins = kFALSE ; }
86 // void SwitchOnClusterPairRWeightBins() {fUseAverClusterPairRWeightBins = kTRUE ; }
87 // void SwitchOffClusterPairRWeightBins(){fUseAverClusterPairRWeightBins = kFALSE ; }
89 // void SwitchOnEMaxBins() {fUseEMaxBins = kTRUE ; }
90 // void SwitchOffEMaxBins() {fUseEMaxBins = kFALSE ; }
92 //-------------------------------
93 //Opening angle pair selection
94 //-------------------------------
95 void SwitchOnAngleSelection() {fUseAngleCut = kTRUE ; }
96 void SwitchOffAngleSelection() {fUseAngleCut = kFALSE ; }
97 void SwitchOnAngleEDepSelection() {fUseAngleEDepCut = kTRUE ; }
98 void SwitchOffAngleEDepSelection() {fUseAngleEDepCut = kFALSE ; }
99 void SetAngleCut(Float_t a) {fAngleCut = a ; }
100 void SetAngleMaxCut(Float_t a) {fAngleMaxCut = a ; }
102 //-------------------------------
103 // Use mixing code of this class
104 //-------------------------------
105 void SwitchOnOwnMix() {fDoOwnMix = kTRUE ; }
106 void SwitchOffOwnMix() {fDoOwnMix = kFALSE ; }
108 //------------------------------------------
109 //Do analysis only with clusters in same SM or different combinations of SM
110 //------------------------------------------
111 void SwitchOnSameSM() {fSameSM = kTRUE ; }
112 void SwitchOffSameSM() {fSameSM = kFALSE ; }
114 void SwitchOnSMCombinations() {fFillSMCombinations = kTRUE ; }
115 void SwitchOffSMCombinations() {fFillSMCombinations = kFALSE ; }
117 //-------------------------------
118 //Histogram filling options off by default
119 //-------------------------------
120 void SwitchOnInvPtWeight() {fMakeInvPtPlots = kTRUE ; }
121 void SwitchOffInvPtWeight() {fMakeInvPtPlots = kFALSE ; }
123 void SwitchOnFillBadDistHisto() {fFillBadDistHisto = kTRUE;}
124 void SwitchOffFillBadDistHisto() {fFillBadDistHisto = kFALSE;}
126 //-------------------------------------------
127 //Cuts for multiple analysis, off by default
128 //-------------------------------------------
129 void SwitchOnMultipleCutAnalysis() {fMultiCutAna = kTRUE ;}
130 void SwitchOffMultipleCutAnalysis() {fMultiCutAna = kFALSE;}
132 void SetNPtCuts (Int_t size) {if(size <= 10)fNPtCuts = size; }
133 void SetNAsymCuts (Int_t size) {if(size <= 10)fNAsymCuts = size; }
134 void SetNNCellCuts(Int_t size) {if(size <= 10)fNCellNCuts = size; }
135 void SetNPIDBits (Int_t size) {if(size <= 10)fNPIDBits = size; }
137 void SetPtCutsAt (Int_t pos,Float_t val) {if(pos < 10)fPtCuts[pos] = val;}
138 void SetAsymCutsAt (Int_t pos,Float_t val) {if(pos < 10)fAsymCuts[pos] = val;}
139 void SetNCellCutsAt(Int_t pos,Int_t val) {if(pos < 10)fCellNCuts[pos] = val;}
140 void SetPIDBitsAt (Int_t pos,Int_t val) {if(pos < 10)fPIDBits[pos] = val;}
142 //MC analysis related methods
143 void FillAcceptanceHistograms();
144 void FillMCVersusRecDataHistograms(const Int_t index1, const Int_t index2,
145 const Float_t pt1, const Float_t pt2,
146 const Int_t ncells1, const Int_t ncells2,
147 const Double_t mass, const Double_t pt, const Double_t asym,
148 const Double_t deta, const Double_t dphi);
150 void SwitchOnMultipleCutAnalysisInSimulation() {fMultiCutAnaSim = kTRUE;}
151 void SwitchOffMultipleCutAnalysisInSimulation() {fMultiCutAnaSim = kFALSE;}
153 void SwitchOnConversionChecker() { fCheckConversion = kTRUE ; }
154 void SwitchOffConversionChecker() { fCheckConversion = kFALSE ; }
158 Bool_t fDoOwnMix; // Do combinatorial background not the one provided by the frame
159 TString fCalorimeter ; // Select Calorimeter for IM
160 Int_t fNModules ; // Number of EMCAL/PHOS modules, set as many histogras as modules
161 Bool_t fUseAngleCut ; // Select pairs depending on their opening angle
162 Bool_t fUseAngleEDepCut ; // Select pairs depending on their opening angle
163 Float_t fAngleCut ; // Select pairs with opening angle larger than a threshold
164 Float_t fAngleMaxCut ; // Select pairs with opening angle smaller than a threshold
165 TList ** fEventsList ; //![GetNCentrBin()*GetNZvertBin()*GetNRPBin()] Containers for photons in stored events
167 //Multiple cuts analysis
168 Bool_t fMultiCutAna; // Do analysis with several or fixed cut
169 Bool_t fMultiCutAnaSim; // Do analysis with several or fixed cut, in the simulation related part
170 Int_t fNPtCuts; // Number of pt cuts
171 Float_t fPtCuts[10]; // Array with different pt cuts
172 Int_t fNAsymCuts; // Number of assymmetry cuts
173 Float_t fAsymCuts[10]; // Array with different assymetry cuts
174 Int_t fNCellNCuts; // Number of cuts with number of cells in cluster
175 Int_t fCellNCuts[10]; // Array with different cell number cluster cuts
176 Int_t fNPIDBits ; // Number of possible PID bit combinations
177 Int_t fPIDBits[10]; // Array with different PID bits
179 //Switchs of different analysis options
180 Bool_t fMakeInvPtPlots; // D plots with inverse pt weight
181 Bool_t fSameSM; // Select only pairs in same SM;
182 Bool_t fFillSMCombinations; // Fill histograms with different cluster pairs in SM combinations
183 Bool_t fCheckConversion; // Fill histograms with tagged photons as conversion
184 Bool_t fUseTrackMultBins; // Use track multiplicity and not centrality bins
185 Bool_t fUsePhotonMultBins; // Use photon multiplicity and not centrality bins
186 Bool_t fUseAverCellEBins; // Use cell average energy and not centrality bins
187 Bool_t fUseAverClusterEBins; // Use cluster average energy and not centrality bins
188 Bool_t fUseAverClusterEDenBins; // Use cluster average energy density and not centrality bins
189 // Bool_t fUseAverClusterPairRBins; // Use cluster average energy and not centrality bins
190 // Bool_t fUseAverClusterPairRWeightBins; // Use cluster average energy and not centrality bins
191 // Bool_t fUseEMaxBins; // Use Emax bins
192 Bool_t fFillBadDistHisto; // Do plots for different distances to bad channels
196 //Event characterization
197 TH1F* fhAverTotECluster; //! Average number of clusters in SM
198 TH1F* fhAverTotECell; //! Average number of cells in SM
199 TH2F* fhAverTotECellvsCluster; //! Average number of cells in SM
200 TH1F* fhEDensityCluster; //! Deposited energy in event per cluster
201 TH1F* fhEDensityCell; //! Deposited energy in event per cell vs cluster
202 TH2F* fhEDensityCellvsCluster; //! Deposited energy in event per cell vs cluster
203 // TH1F* fhClusterPairDist; //! Distance between clusters
204 // TH1F* fhClusterPairDistWeight; //! Distance between clusters weighted by pair energy
205 // TH1F* fhAverClusterPairDist; //! Average distance between cluster pairs
206 // TH1F* fhAverClusterPairDistWeight;//! Average distance between cluster pairs weighted with pair energy
207 // TH2F* fhAverClusterPairDistvsAverE; //! Average distance between cluster pairs vs average cluster energy
208 // TH2F* fhAverClusterPairDistWeightvsAverE; //! Average distance between cluster pairs weighted with pair energy vs average cluster energy
209 // TH2F* fhAverClusterPairDistvsN; //! Average distance between cluster pairs vs number of clusters
210 // TH2F* fhAverClusterPairDistWeightvsN; //! Average distance between cluster pairs weighted with pair energy vs number of clusters
211 // TH2F* fhMaxEvsClustMult; //!
212 // TH2F* fhMaxEvsClustEDen; //!
215 TH2D ** fhReMod ; //![fNModules] REAL two-photon invariant mass distribution for different calorimeter modules.
216 TH2D ** fhReSameSideEMCALMod ; //![fNModules-2] REAL two-photon invariant mass distribution for different clusters in different calorimeter modules.
217 TH2D ** fhReSameSectorEMCALMod ; //![fNModules/2] REAL two-photon invariant mass distribution for different clusters in different calorimeter modules.
218 TH2D ** fhReDiffPHOSMod ; //![fNModules] REAL two-photon invariant mass distribution for different clusters in different calorimeter modules.
219 TH2D ** fhMiMod ; //![fNModules] MIXED two-photon invariant mass distribution for different calorimeter modules.
220 TH2D ** fhMiSameSideEMCALMod ; //![fNModules-2] REAL two-photon invariant mass distribution for different clusters in different calorimeter modules.
221 TH2D ** fhMiSameSectorEMCALMod ; //![fNModules/2] REAL two-photon invariant mass distribution for different clusters in different calorimeter modules.
222 TH2D ** fhMiDiffPHOSMod ; //![fNModules-1] REAL two-photon invariant mass distribution for different clusters in different calorimeter modules.
224 // Pairs with at least one cluster tagged as conversion
225 TH2D * fhReConv ; //! REAL two-photon invariant mass distribution one of the pair was 2 clusters with small mass
226 TH2D * fhMiConv ; //! MIXED two-photon invariant mass distribution one of the pair was 2 clusters with small mass
227 TH2D * fhReConv2 ; //! REAL two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass
228 TH2D * fhMiConv2 ; //! MIXED two-photon invariant mass distribution both pair photons recombined from 2 clusters with small mass
230 TH2D ** fhRe1 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL two-photon invariant mass distribution for different centralities and Asymmetry
231 TH2D ** fhMi1 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
232 TH2D ** fhRe2 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL two-photon invariant mass distribution for different centralities and Asymmetry
233 TH2D ** fhMi2 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
234 TH2D ** fhRe3 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL two-photon invariant mass distribution for different centralities and Asymmetry
235 TH2D ** fhMi3 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry
237 //Histograms weighted by inverse pT
238 TH2D ** fhReInvPt1 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
239 TH2D ** fhMiInvPt1 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
240 TH2D ** fhReInvPt2 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
241 TH2D ** fhMiInvPt2 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
242 TH2D ** fhReInvPt3 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] REAL two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
243 TH2D ** fhMiInvPt3 ; //![GetNCentrBin()*fNPIDBits*fNAsymCuts] MIXED two-photon invariant mass distribution for different centralities and Asymmetry, inverse pT
245 //Multiple cuts: Assymmetry, pt, n cells, PID
246 TH2D ** fhRePtNCellAsymCuts ; //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
247 TH2D ** fhRePtNCellAsymCutsSM0 ; //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
248 TH2D ** fhRePtNCellAsymCutsSM1 ; //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
249 TH2D ** fhRePtNCellAsymCutsSM2 ; //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
250 TH2D ** fhRePtNCellAsymCutsSM3 ; //![fNPtCuts*fNAsymCuts*fNCellNCuts] REAL two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
251 TH2D ** fhMiPtNCellAsymCuts ; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Mixed two-photon invariant mass distribution for different pt cut, n cell cuts and assymetry
252 TH2D ** fhRePIDBits ; //![fNPIDBits] REAL two-photon invariant mass distribution for different PID bits
253 TH3D ** fhRePtMult ; //![fNAsymCuts] REAL two-photon invariant mass distribution for different track multiplicity and assymetry cuts
255 // Asymmetry vs pt, in pi0/eta regions
256 TH2D * fhRePtAsym ; //! REAL two-photon pt vs asymmetry
257 TH2D * fhRePtAsymPi0 ; //! REAL two-photon pt vs asymmetry, close to pi0 mass
258 TH2D * fhRePtAsymEta ; //! REAL two-photon pt vs asymmetry, close to eta mass
260 //Centrality, Event plane bins
261 TH3D * fhEvents; //! Number of events per centrality, RP, zbin
262 TH1D * fhCentrality; //! Histogram with centrality bins with at least one pare
263 TH1D * fhCentralityNoPair; //! Histogram with centrality bins with no pair
265 TH1D * fhEventPlaneAngle; //! Histogram with Event plane angle
266 TH2D * fhEventPlaneResolution; //! Histogram with Event plane resolution vs centrality
268 // Pair opening angle
269 TH2D * fhRealOpeningAngle ; //! Opening angle of pair versus pair energy
270 TH2D * fhRealCosOpeningAngle ; //! Cosinus of opening angle of pair version pair energy
271 TH2D * fhMixedOpeningAngle ; //! Opening angle of pair versus pair energy
272 TH2D * fhMixedCosOpeningAngle ; //! Cosinus of opening angle of pair version pair energy
274 //MC analysis histograms
276 TH1D * fhPrimPi0Pt ; //! Spectrum of Primary
277 TH1D * fhPrimPi0AccPt ; //! Spectrum of primary with accepted daughters
278 TH2D * fhPrimPi0Y ; //! Rapidity distribution of primary particles vs pT
279 TH2D * fhPrimPi0AccY ; //! Rapidity distribution of primary with accepted daughters vs pT
280 TH2D * fhPrimPi0Phi ; //! Azimutal distribution of primary particles vs pT
281 TH2D * fhPrimPi0AccPhi; //! Azimutal distribution of primary with accepted daughters vs pT
282 TH2D * fhPrimPi0OpeningAngle ; //! Opening angle of pair versus pair energy, primaries
283 TH2D * fhPrimPi0CosOpeningAngle ; //! Cosinus of opening angle of pair version pair energy, primaries
285 TH1D * fhPrimEtaPt ; //! Spectrum of Primary
286 TH1D * fhPrimEtaAccPt ; //! Spectrum of primary with accepted daughters
287 TH2D * fhPrimEtaY ; //! Rapidity distribution of primary particles vs pT
288 TH2D * fhPrimEtaAccY ; //! Rapidity distribution of primary with accepted daughters vs pT
289 TH2D * fhPrimEtaPhi ; //! Azimutal distribution of primary particles vs pT
290 TH2D * fhPrimEtaAccPhi; //! Azimutal distribution of primary with accepted daughters vs pT
293 TH2D * fhPrimPi0PtOrigin ; //! Spectrum of generated pi0 vs mother
294 TH2D * fhPrimEtaPtOrigin ; //! Spectrum of generated eta vs mother
297 //Array of histograms ordered as follows: 0-Photon, 1-electron, 2-pi0, 3-eta, 4-a-proton, 5-a-neutron, 6-stable particles,
298 // 7-other decays, 8-string, 9-final parton, 10-initial parton, intermediate, 11-colliding proton, 12-unrelated
299 TH2D * fhMCOrgMass[13]; //! Mass vs pt of real pairs, check common origin of pair
300 TH2D * fhMCOrgAsym[13]; //! Asymmetry vs pt of real pairs, check common origin of pair
301 TH2D * fhMCOrgDeltaEta[13]; //! Delta Eta vs pt of real pairs, check common origin of pair
302 TH2D * fhMCOrgDeltaPhi[13]; //! Delta Phi vs pt of real pairs, check common origin of pair
304 //Multiple cuts in simulation, origin pi0 or eta
305 TH2D ** fhMCPi0MassPtRec; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs reconstructed pt of original pair
306 TH2D ** fhMCPi0MassPtTrue; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed mass vs generated pt of original pair
307 TH2D ** fhMCPi0PtTruePtRec; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real pi0 pairs, reconstructed pt vs generated pt of pair
308 TH2D ** fhMCEtaMassPtRec; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs reconstructed pt of original pair
309 TH2D ** fhMCEtaMassPtTrue; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed mass vs generated pt of original pair
310 TH2D ** fhMCEtaPtTruePtRec; //![fNPtCuts*fNAsymCuts*fNCellNCuts] Real eta pairs, reconstructed pt vs generated pt of pair
312 TH2D * fhMCPi0PtOrigin ; //! Mass of reoconstructed pi0 pairs in calorimeter vs mother
313 TH2D * fhMCEtaPtOrigin ; //! Mass of reoconstructed pi0 pairs in calorimeter vs mother
315 ClassDef(AliAnaPi0,19)