]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/PartCorrDep/AliAnaNeutralMeson.h
fixing compilation problem
[u/mrichter/AliRoot.git] / PWG4 / PartCorrDep / AliAnaNeutralMeson.h
CommitLineData
9725fd2a 1#ifndef ALIANANEUTRALMESON_H\r
2#define ALIANANEUTRALMESON_H\r
3/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
4 * See cxx source for full Copyright notice */\r
5/* $Id: $ */\r
6\r
7//_________________________________________________________________________\r
8//a general class to fill two-photon and pi0+gamma invariant mass hisograms \r
9// to be used to extract pi0, eta and omega raw yield.\r
10// also for PHOS acceptance \r
11//--\r
12//-- by Renzhuo Wan (Iopp-wuhan) May 13,2009\r
13//_________________________________________________________________________\r
14\r
15//Root\r
16class TList;\r
17class TH3F ;\r
18class TH2F ;\r
19class TLorentzVector;\r
20//Analysis\r
21class AliAODEvent ;\r
22class AliESDEvent ;\r
23#include "AliAnaPartCorrBaseClass.h"\r
24class TParticle;\r
25\r
591cc579 26#ifdef __PHOSUTIL__\r
9725fd2a 27 class AliPHOSGeoUtils;\r
28#endif\r
591cc579 29
30#ifdef __EMCALUTIL__\r
31 class AliEMCALGeoUtils;\r
32#endif\r
33
9725fd2a 34\r
35class AliAnaNeutralMeson : public AliAnaPartCorrBaseClass {\r
36 \r
37 public: \r
38 \r
39 AliAnaNeutralMeson() ; // default ctor\r
40 AliAnaNeutralMeson(const char *name) ; // default ctor\r
41 AliAnaNeutralMeson(const AliAnaNeutralMeson & g) ; // cpy ctor\r
42 AliAnaNeutralMeson & operator = (const AliAnaNeutralMeson & api0) ;//cpy assignment\r
43 virtual ~AliAnaNeutralMeson() ;//virtual dtor\r
44 \r
45 //......\r
46 TList * GetCreateOutputObjects(); \r
47 void Print(const Option_t * opt) const;\r
48 \r
49 void Init();\r
50 void InitParameters();\r
51 void MakeAnalysisFillHistograms();\r
52 //void Terminate(TList * outList);\r
53 //......\r
54 void Get2GammaAsyDistPtM(AliAODPWG4Particle *p1, AliAODPWG4Particle *p2, Double_t &asy, Double_t &dist, \r
55 Double_t &pt, Double_t &mass); //get dist,asy, pair pt, pairmass btween two ptc\r
56 void GetAsyDistPtM(TLorentzVector photon1, TLorentzVector photon2, Double_t &asy, Double_t &dist, \r
57 Double_t &pt, Double_t &mass); //get dist,asy, pair pt, pairmass btween two ptc\r
58 Bool_t IsPhotonSelected(AliAODPWG4Particle *p1, Int_t ipid, Int_t nmod); //select the photon based on pid and nmod to be analyzed\r
59\r
60 void RealPi0Eta(TClonesArray * array) ; //perform the 2gamma IVM analysis for one event to extract pi0 and eta\r
61 void MixPi0Eta(TClonesArray * array1, TClonesArray *array2) ; //perform the 2gamma IVM analysis for mixed eventsto reconstruct the background\r
62 void RealOmega(TClonesArray *array) ; //perform the pi0+gamma->3gamma IVM analysis for one event to extract omega\r
63 void MixOmegaAB(TClonesArray * array1, TClonesArray * array2); ////perform the 2gamma IVM analysis for mixed events to reconstruct the bacground (r1_event1+r2_event1)+r3_event2 and (r1_event1+r2_event2)+r3_event2\r
64 void MixOmegaC(TClonesArray * array1, TClonesArray * array2, TClonesArray * array3); //perform the 2gamma IVM analysis for mixed events to reconstruct the bacground (r1_event1+r2_event2)+r3_event3\r
65 \r
66 void PhotonAcceptance(AliStack *stack); //to get MC photon pt and photon can be hit on dectector\r
67 void Pi0EtaAcceptance(AliStack *stack); //to get MC pi0 and eta pt and photon can be hit on dectector\r
68 void OmegaAcceptance(AliStack * stack); //to get MC omega pt and photon can be hit on dectector\r
69 \r
70 void GetMCDistAsy(TParticle * p1, TParticle * p2, Double_t & dist, Double_t &asy); //to get MC information between two particle\r
71 \r
72 void SetNAsyBinsMinMax(Int_t bins, Double_t min, Double_t max) {fNbinsAsy=bins; fMinAsy=min; fMaxAsy=max; } //set Asy bins, min and max\r
73 void SetNPtBinsMinMax(Int_t bins, Double_t min, Double_t max) {fNbinsPt=bins; fMinPt=min; fMaxPt=max; } //set pt bins, min and max\r
74 void SetNMassBinsMinMas(Int_t bins, Double_t min, Double_t max) {fNbinsM=bins; fMinM=min; fMaxM=max; } //set mass pt bins, min and max\r
75 void SetNEventsMixed(Int_t nevents) { fNmaxMixEv=nevents;} //events to be mixed \r
76 \r
77 TString GetCalorimeter() const { return fCalorimeter ; } //PHOS or EMCAL\r
78 void SetCalorimeter(TString det);// { fCalorimeter = det ; } \r
79 \r
80 void SetNPID(Int_t pid) {fNPID=pid;} //number of PID cut\r
81 void SetInvMassCut(Int_t mass) {fInvMassCut = mass;} // \r
82 void SetPi0MassPeakWidthCut(Double_t cut) {fPi0MassPeakWidthCut = cut;}\r
83 \r
84 Bool_t SetAnaPi0Eta(Bool_t ana) {fAnaPi0Eta = ana; return fAnaPi0Eta;} //analysis pi0 and eta ????\r
85 Bool_t SetAnaOmega(Bool_t ana) {fAnaOmega = ana; return fAnaOmega ;} //analysis omega ????\r
86 \r
87 Bool_t IsBadRun(Int_t /*iRun*/) const {return kFALSE;} //Tests if this run bad according to private list\r
88 \r
89 private:\r
90\r
91 Bool_t fAnaPi0Eta; //whether ana the 2gamma IVM to extract pi0 and eta spectrum\r
92 Bool_t fAnaOmega; //whether ana the 3gamma IVM to extract omega spectrum\r
93\r
94 Int_t fNPID; //PID bin\r
95 Int_t fNmaxMixEv ; //number events to be mixed\r
96 Int_t fNAsy; //number of Asy cut\r
97 Int_t fNMod; //number of PHOS or EMCAL modules to be analyzed\r
98 Int_t fNHistos; //number of histograms\r
99 Double_t fPerPtBin;//plot the IVM distribution per pt bin\r
100 \r
101 Double_t *fAsyCut; // asy cut(< 0.7, 0.8, 0.9, 1)\r
102 Int_t *fModCut; // mod cut(1,3,5) for PHOS, and 0 (default) for EMCAL\r
103\r
104 Double_t fDist; //distance between two particle in (eta, phi) plane\r
105 Double_t fAsy; //asysmetry |(e1-e2)|/(e1+e2)\r
106 Double_t fPt; //pair pt\r
107 Double_t fMass; //pair mass\r
108 Double_t fInvMassCut; // cut of invmass\r
109\r
110 Int_t fNbinsAsy; //Asy bin number, min and max\r
111 Double_t fMinAsy;\r
112 Double_t fMaxAsy;\r
113 Int_t fNbinsPt; //Pt bin number, min and max\r
114 Double_t fMinPt;\r
115 Double_t fMaxPt;\r
116 Int_t fNbinsM; //mass bin number, min and max\r
117 Double_t fMinM;\r
118 Double_t fMaxM;\r
119\r
120 //! containers for photons in stored events\r
121 TList * fEventsListPi0Eta ; //! events with cluster larger than 2 for pi0 and eta analysis\r
122 TList * fEventsListOmega ; //! events with cluster larger than 3 for omega analysis\r
123\r
124 Double_t fPi0Mass; //nominal pi0 mass\r
125 Double_t fPi0MassPeakWidthCut; //pi0 candidate\r
126 \r
127 TH1I * fhNClusters; //cluster multiplicity distribution per event\r
128 TH1D ** fhRecPhoton; //reconstructed photon pt distribution\r
129 TH2F ** fhRecPhotonEtaPhi; //photon (eta, phi) distruction\r
130\r
131 TH1F ** fReal2Gamma; //real 2gamma IVM in per pt bin (ipid, imod, iasy)\r
132 TH1F ** fMix2Gamma; //mix 2gamma IVM in per pt bin (ipid, imod, iasy, pt)\r
133 TH1F ** fRealOmega; //real omega IVM in per pt bin (ipid, imod, iasy, pt)\r
134 TH1F ** fMixOmegaA; //mixA omega IVM in per pt bin (ipid, imod, iasy, pt) (r1_event1+r2_event1)+r3_event2\r
135 TH1F ** fMixOmegaB; //mixB omega IVM in per pt bin (ipid, imod, iasy, pt) (r1_event1+r2_event2)+r3_event2\r
136 TH1F ** fMixOmegaC; //mixC omega IVM in per pt bin (ipid, imod, iasy, pt) (r1_event1+r2_event2)+r3_event3 \r
137\r
138 TH3F ** fRealTwoGammaAsyPtM; //real 2gamma IVM(asy, pt, m) without PID, mod cut\r
139 TH3F ** fMixTwoGammaAsyPtM; //mix 2gamma IVM(asy, pt, m) without PID, mod cut\r
140 TH3F ** fRealPi0GammaAsyPtM; //real omega IVM(asy, pt, m) without PID, mod cut\r
141 TH3F ** fMixAPi0GammaAsyPtM; //mixA omega IVM(asy, pt, m) without PID, mod cut\r
142 TH3F ** fMixBPi0GammaAsyPtM; //mixB omega IVM(asy, pt, m) without PID, mod cut\r
143 TH3F ** fMixCPi0GammaAsyPtM; //mixC omega IVM(asy, pt, m) without PID, mod cut\r
144\r
145 //Histograms filled only if MC data is requested \r
146 //pt distribution\r
147 //currently only PHOS included \r
148 TH1D * fhPrimPhotonPt ; //! Spectrum of Primary \r
149 TH1D * fhPrimPhotonAccPt ; //! Spectrum of primary with accepted daughters \r
150 TH1D * fhPrimPhotonY ; //! Rapidity distribution of primary particles\r
151 TH1D * fhPrimPhotonAccY ; //! Rapidity distribution of primary with accepted daughters\r
152 TH1D * fhPrimPhotonPhi ; //! Azimutal distribution of primary particles\r
153 TH1D * fhPrimPhotonAccPhi; //! Azimutal distribution of primary with accepted daughters \r
154 \r
155 TH1D * fhPrimPi0Pt ; //! Spectrum of Primary \r
156 TH1D * fhPrimPi0AccPt ; //! Spectrum of primary with accepted daughters \r
157 TH1D * fhPrimPi0Y ; //! Rapidity distribution of primary particles\r
158 TH1D * fhPrimPi0AccY ; //! Rapidity distribution of primary with accepted daughters\r
159 TH1D * fhPrimPi0Phi ; //! Azimutal distribution of primary particles\r
160 TH1D * fhPrimPi0AccPhi; //! Azimutal distribution of primary with accepted daughters \r
161 \r
162 TH1D * fhPrimEtaPt ; //! Spectrum of Primary \r
163 TH1D * fhPrimEtaAccPt ; //! Spectrum of primary with accepted daughters \r
164 TH1D * fhPrimEtaY ; //! Rapidity distribution of primary particles\r
165 TH1D * fhPrimEtaAccY ; //! Rapidity distribution of primary with accepted daughters\r
166 TH1D * fhPrimEtaPhi ; //! Azimutal distribution of primary particles\r
167 TH1D * fhPrimEtaAccPhi; //! Azimutal distribution of primary with accepted daughters\r
168 \r
169 TH1D * fhPrimOmegaPt ; //! Spectrum of Primary \r
170 TH1D * fhPrimOmegaAccPt ; //! Spectrum of primary with accepted daughters \r
171 TH1D * fhPrimOmegaY ; //! Rapidity distribution of primary particles\r
172 TH1D * fhPrimOmegaAccY ; //! Rapidity distribution of primary with accepted daughters\r
173 TH1D * fhPrimOmegaPhi ; //! Azimutal distribution of primary particles\r
174 TH1D * fhPrimOmegaAccPhi; //! Azimutal distribution of primary with accepted daughters\r
175\r
176 TString fCalorimeter ; //Select Calorimeter \r
177 \r
591cc579 178#ifdef __PHOSUTIL__\r
9725fd2a 179 AliPHOSGeoUtils * fPHOSGeo ; //! PHOS geometry pointer\r
180#endif \r
591cc579 181#ifdef __EMCALUTIL__\r
182 AliEMCALGeoUtils * fEMCALGeo ; //! EMCAL geometry pointer\r
183#endif \r
9725fd2a 184 \r
185 ClassDef(AliAnaNeutralMeson,1)\r
186} ;\r
187\r
188\r
189#endif //ALIANANEUTRALMESON\r
190\r
191\r
192\r