1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: Boris Polishchuk (Boris.Polishchuk@cern.ch) *
6 * Permission to use, copy, modify and distribute this software and its *
7 * documentation strictly for non-commercial purposes is hereby granted *
8 * without fee, provided that the above copyright notice appears in all *
9 * copies and that both the copyright notice and this permission notice *
10 * appear in the supporting documentation. The authors make no claims *
11 * about the suitability of this software for any purpose. It is *
12 * provided "as is" without express or implied warranty. *
13 **************************************************************************/
15 ////////////////////////////////////////////////////////////////////////////
16 //------------------------------------------------------------------------//
17 // Class used to do omega(782) -> pi0 pi+ pi- analysis. //
18 //------------------------------------------------------------------------//
19 ////////////////////////////////////////////////////////////////////////////
21 #include "AliAnalysisTaskOmegaPi0PiPi.h"
25 #include "AliESDEvent.h"
26 #include "TRefArray.h"
27 #include "TLorentzVector.h"
28 #include "AliESDCaloCluster.h"
30 #include "AliESDtrack.h"
33 ClassImp(AliAnalysisTaskOmegaPi0PiPi)
35 AliAnalysisTaskOmegaPi0PiPi::AliAnalysisTaskOmegaPi0PiPi() :
36 AliAnalysisTaskSE(),fOutputContainer(0x0),fhM2piSel(0x0),fhDxy(0x0),fhMggSel(0x0),
37 fhImpXY(0x0),fhM3pi0to2(0x0),fhM3pi2to4(0x0),fhM3pi4to6(0x0),fhM3pi6to8(0x0),
38 fhM3pi8to10(0x0),fhM3pi10to12(0x0),fhM3pi12(0x0),fhM3pi0to8(0x0)
40 //Default constructor.
44 AliAnalysisTaskOmegaPi0PiPi::AliAnalysisTaskOmegaPi0PiPi(const char* name) :
45 AliAnalysisTaskSE(name),fOutputContainer(0x0),fhM2piSel(0x0),fhDxy(0x0),fhMggSel(0x0),
46 fhImpXY(0x0),fhM3pi0to2(0x0),fhM3pi2to4(0x0),fhM3pi4to6(0x0),fhM3pi6to8(0x0),
47 fhM3pi8to10(0x0),fhM3pi10to12(0x0),fhM3pi12(0x0),fhM3pi0to8(0x0)
49 //Constructor used to create a named task.
51 DefineOutput(1,TList::Class());
55 AliAnalysisTaskOmegaPi0PiPi::~AliAnalysisTaskOmegaPi0PiPi()
60 fOutputContainer->Delete() ;
61 delete fOutputContainer ;
66 void AliAnalysisTaskOmegaPi0PiPi::UserCreateOutputObjects()
68 //Allocates histograms and puts it to the output container.
70 fOutputContainer = new TList();
72 fhM2piSel = new TH1F("hM2pi_sel","V0 inv. mass, Dxy<2cm",100,0.,1.);
73 fOutputContainer->Add(fhM2piSel);
75 fhDxy = new TH1F("hDxy","Distance from V0 to primary vertex in XY-plane", 500, 0.,500.);
76 fOutputContainer->Add(fhDxy);
78 fhMggSel = new TH1F("hMgg_sel","Inv. mass of two clusters, pT>1GeV",100,0.,0.3);
79 fOutputContainer->Add(fhMggSel);
81 fhImpXY = new TH1F("hImpXY","Impact parameters in XY-plane",1000,0.,1.);
82 fOutputContainer->Add(fhImpXY);
85 //M(3pi) vs pT(3pi), ptrack < 2GeV
86 fhM3pi0to2 = new TH2F("hM3pi_0_2","pTrack: 0-2GeV",100,0.,10., 200,0.4,1.);
87 fOutputContainer->Add(fhM3pi0to2);
89 //M(3pi) vs pT(3pi), 2GeV < ptrack < 4GeV
90 fhM3pi2to4 = new TH2F("hM3pi_2_4","pTrack: 2-4GeV",100,0.,10., 200,0.4,1.);
91 fOutputContainer->Add(fhM3pi2to4);
93 //M(3pi) vs pT(3pi), 4GeV < ptrack < 6GeV
94 fhM3pi4to6 = new TH2F("hM3pi_4_6","pTrack: 4-6GeV",100,0.,10., 200,0.4,1.);
95 fOutputContainer->Add(fhM3pi4to6);
97 //M(3pi) vs pT(3pi), 6GeV < ptrack < 8GeV
98 fhM3pi6to8 = new TH2F("hM3pi_6_8","pTrack: 6-8GeV",100,0.,10., 200,0.4,1.);
99 fOutputContainer->Add(fhM3pi6to8);
101 //M(3pi) vs pT(3pi), 8GeV < ptrack < 10GeV
102 fhM3pi8to10 = new TH2F("hM3pi_8_10","pTrack: 8-10GeV",100,0.,10., 200,0.4,1.);
103 fOutputContainer->Add(fhM3pi8to10);
105 //M(3pi) vs pT(3pi), 10GeV < ptrack < 12GeV
106 fhM3pi10to12 = new TH2F("hM3pi_10_12","pTrack: 10-12GeV",100,0.,10., 200,0.4,1.);
107 fOutputContainer->Add(fhM3pi10to12);
109 //M(3pi) vs pT(3pi), ptrack > 12GeV
110 fhM3pi12 = new TH2F("hM3pi_12","pTrack>12GeV",100,0.,10., 200,0.4,1.);
111 fOutputContainer->Add(fhM3pi12);
113 //M(3pi) vs pT(3pi), ptrack < 8GeV
114 fhM3pi0to8 = new TH2F("hM3pi_0_8","pTrack: 0-8GeV",100,0.,10., 200,0.4,1.);
115 fOutputContainer->Add(fhM3pi0to8);
120 void AliAnalysisTaskOmegaPi0PiPi::UserExec(Option_t* /* option */)
122 //Main function that does all the job.
124 AliVEvent* event = InputEvent();
125 AliESDEvent* esd = (AliESDEvent*)event ;
127 Double_t v[3] ; //vertex ;
128 esd->GetVertex()->GetXYZ(v) ;
131 AliDebug(2,Form("Vertex: (%.3f,%.3f,%.3f)\n",vtx.X(),vtx.Y(),vtx.Z()));
133 TRefArray * caloClustersArr = new TRefArray();
134 esd->GetPHOSClusters(caloClustersArr);
136 const Int_t kNumberOfTracks = esd->GetNumberOfTracks();
137 const Int_t kNumberOfV0s = esd->GetNumberOfV0s();
139 AliDebug(2,Form("Tracks: %d. V0s: %d\n",kNumberOfTracks,kNumberOfV0s));
141 Float_t xyImp,zImp,imp1,imp2;
143 const Int_t kNumberOfPhosClusters = caloClustersArr->GetEntries() ;
144 AliDebug(2,Form("CaloClusters: %d\n", kNumberOfPhosClusters));
145 if(kNumberOfPhosClusters<2){
146 delete caloClustersArr;
149 TLorentzVector pc1; //4-momentum of PHOS cluster 1
150 TLorentzVector pc2; //4-momentum of PHOS cluster 2
153 TLorentzVector ppos; //4-momentum of positive track
154 TLorentzVector pneg; //4-momentum of negative track
155 TLorentzVector ptr12;
157 TLorentzVector pomega;
159 Double_t p1,p2; // 3-momentum of tracks
161 const Double_t kMpi = 0.1396;
163 for(Int_t iClu=0; iClu<kNumberOfPhosClusters; iClu++) {
164 AliESDCaloCluster *c1 = (AliESDCaloCluster *) caloClustersArr->At(iClu);
165 c1->GetMomentum(pc1,v);
167 for (Int_t jClu=iClu; jClu<kNumberOfPhosClusters; jClu++) {
168 AliESDCaloCluster *c2 = (AliESDCaloCluster *) caloClustersArr->At(jClu);
169 if(c2->IsEqual(c1)) continue;
170 c2->GetMomentum(pc2,v);
173 AliDebug(2,Form("pc12.M(): %.3f\n",pc12.M()));
175 if(pc12.M()<0.115 || pc12.M()>0.155) continue; //not a pi0 candidate!
176 if(pc12.Pt()<1.) continue; //pT(pi0) > 1GeV
178 fhMggSel->Fill(pc12.M());
180 for(Int_t iTrack=0; iTrack<kNumberOfTracks; iTrack++) {
181 AliESDtrack* tr1 = esd->GetTrack(iTrack);
183 tr1->GetImpactParameters(xyImp,zImp);
184 imp1 = TMath::Abs(xyImp);
186 Short_t charge = tr1->Charge();
187 if(imp1>0.004) continue; // not from the primary vertex!
189 etrack = TMath::Sqrt(p1*p1 + kMpi*kMpi);
190 ppos.SetPxPyPzE(tr1->Px(),tr1->Py(),tr1->Pz(),etrack);
192 for(Int_t jTrack=iTrack; jTrack<kNumberOfTracks; jTrack++) {
193 AliESDtrack* tr2 = esd->GetTrack(jTrack);
195 if(tr2->Charge()==charge) continue; //same sign track (WRONG!!)
196 tr2->GetImpactParameters(xyImp,zImp);
197 imp2=TMath::Abs(xyImp);
198 if(imp2>0.004) continue; // not from the primary vertex!
200 etrack = TMath::Sqrt(p2*p2 + kMpi*kMpi);
201 pneg.SetPxPyPzE(tr2->Px(),tr2->Py(),tr2->Pz(),etrack);
205 // printf("ptr12.M()=%.3f, xyImp1=%f, xyImp2=%f, ch1=%d, ch2=%d,
206 // p1=%.3f, p2=%.3f, m1=%.3f, m2=%.3f\n",
207 // ptr12.M(),imp1,imp2,tr1->Charge(),tr2->Charge(),
208 // tr1->P(),tr2->P(),tr1->M(),tr2->M());
210 pomega = pc12 + ptr12;
211 AliDebug(2,Form("pomega.M(): %f\n",pomega.M()));
214 fhM3pi0to2->Fill(pomega.Pt(),pomega.M());
216 if( (p1>=2. && p1<4.) && (p2>=2. && p2<4.) )
217 fhM3pi2to4->Fill(pomega.Pt(),pomega.M());
219 if( (p1>=4. && p1<6.) && (p2>=4. && p2<6.) )
220 fhM3pi4to6->Fill(pomega.Pt(),pomega.M());
222 if( (p1>=6. && p1<8.) && (p2>=6. && p2<8.) )
223 fhM3pi6to8->Fill(pomega.Pt(),pomega.M());
225 if( (p1>=8. && p1<10.) && (p2>=8. && p2<10.) )
226 fhM3pi8to10->Fill(pomega.Pt(),pomega.M());
228 if( (p1>=10. && p1<12.) && (p2>=10. && p2<12.) )
229 fhM3pi10to12->Fill(pomega.Pt(),pomega.M());
231 if( (p1>=12.) && (p2>=12.) )
232 fhM3pi12->Fill(pomega.Pt(),pomega.M());
235 fhM3pi0to8->Fill(pomega.Pt(),pomega.M());
240 for(Int_t iVert=0; iVert<kNumberOfV0s; iVert++) {
241 AliESDv0* v0 = esd->GetV0(iVert);
242 AliESDtrack* ptrack = esd->GetTrack(v0->GetPindex()); //positive track
243 AliESDtrack* ntrack = esd->GetTrack(v0->GetNindex()); //negative track
245 etrack = TMath::Sqrt(ptrack->P()*ptrack->P() + kMpi*kMpi);
246 ppos.SetPxPyPzE(ptrack->Px(),ptrack->Py(),ptrack->Pz(),etrack);
248 etrack = TMath::Sqrt(ntrack->P()*ntrack->P() + kMpi*kMpi);
249 pneg.SetPxPyPzE(ntrack->Px(),ntrack->Py(),ntrack->Pz(),etrack);
253 Double_t dx = vtx.X() - v0->Xv();
254 Double_t dy = vtx.Y() - v0->Yv();
255 Double_t dxy = TMath::Sqrt(dx*dx + dy*dy);
258 AliDebug(2,Form("V0: dxy=%.3f\n",dxy));
260 if(dxy<2.) fhM2piSel->Fill(ptr12.M());
266 delete caloClustersArr;
267 PostData(1,fOutputContainer);