]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONmassPlot_ESD.C
Adding histo with ITS Z vertex information (Sebastian)
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot_ESD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 // ROOT includes
3 #include "TTree.h"
4 #include "TBranch.h"
5 #include "TClonesArray.h"
6 #include "TLorentzVector.h"
7 #include "TFile.h"
8 #include "TH1.h"
9 #include "TH2.h"
10 #include "TParticle.h"
11 #include "TTree.h"
12 #include <Riostream.h>
13
14 // STEER includes
15 #include "AliRun.h"
16 #include "AliRunLoader.h"
17 #include "AliHeader.h"
18 #include "AliLoader.h"
19 #include "AliStack.h"
20 #include "AliMagF.h"
21 #include "AliESD.h"
22
23 // MUON includes
24 #include "AliESDMuonTrack.h"
25 #endif
26 //
27 // Macro MUONmassPlot.C for ESD
28 // Ch. Finck, Subatech, April. 2004
29 //
30
31 // macro to make invariant mass plots
32 // for combinations of 2 muons with opposite charges,
33 // from root file "MUON.tracks.root" containing the result of track reconstruction.
34 // Histograms are stored on the "MUONmassPlot.root" file.
35 // introducing TLorentzVector for parameter calculations (Pt, P,rap,etc...)
36 // using Invariant Mass for rapidity.
37
38 // Arguments:
39 //   FirstEvent (default 0)
40 //   LastEvent (default 0)
41 //   ResType (default 553)
42 //      553 for Upsilon, anything else for J/Psi
43 //   Chi2Cut (default 100)
44 //      to keep only tracks with chi2 per d.o.f. < Chi2Cut
45 //   PtCutMin (default 1)
46 //      to keep only tracks with transverse momentum > PtCutMin
47 //   PtCutMax (default 10000)
48 //      to keep only tracks with transverse momentum < PtCutMax
49 //   massMin (default 9.17 for Upsilon) 
50 //      &  massMax (default 9.77 for Upsilon) 
51 //         to calculate the reconstruction efficiency for resonances with invariant mass
52 //         massMin < mass < massMax.
53
54 // Add parameters and histograms for analysis 
55
56 Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t LastEvent = 10000,
57                   char* esdFileName = "AliESDs.root", Int_t ResType = 553, 
58                   Float_t Chi2Cut = 100., Float_t PtCutMin = 1., Float_t PtCutMax = 10000.,
59                   Float_t massMin = 9.17,Float_t massMax = 9.77)
60 {
61   cout << "MUONmassPlot " << endl;
62   cout << "FirstEvent " << FirstEvent << endl;
63   cout << "LastEvent " << LastEvent << endl;
64   cout << "ResType " << ResType << endl;
65   cout << "Chi2Cut " << Chi2Cut << endl;
66   cout << "PtCutMin " << PtCutMin << endl;
67   cout << "PtCutMax " << PtCutMax << endl;
68   cout << "massMin " << massMin << endl;
69   cout << "massMax " << massMax << endl;
70
71  
72   //Reset ROOT and connect tree file
73   gROOT->Reset();
74
75   // File for histograms and histogram booking
76   TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
77   TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
78   TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", 100, 0., 20.);
79   TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", 100, 0., 20.);
80   TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
81   TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
82   TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
83   TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
84 TH2F *hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
85 TH2F *hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
86 TH1F *hInvMassRes;
87   TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",120,-12,12);
88
89   if (ResType == 553) {
90     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
91   } else {
92     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
93   }
94
95   TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
96   TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
97   TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
98   TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
99   TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
100   TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
101
102
103   // settings
104   Int_t EventInMass = 0;
105   Int_t EventInMassMatch = 0;
106   Int_t NbTrigger = 0;
107
108   Float_t muonMass = 0.105658389;
109 //   Float_t UpsilonMass = 9.46037;
110 //   Float_t JPsiMass = 3.097;
111
112   Double_t thetaX, thetaY, pYZ;
113   Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
114   Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
115   Int_t fCharge, fCharge2;
116
117   Int_t ntrackhits, nevents;
118   Double_t fitfmin;
119   Double_t fZVertex;
120
121  
122   TLorentzVector fV1, fV2, fVtot;
123
124   // set off mag field 
125   AliMagF::SetReadField(kFALSE);
126
127   // open run loader and load gAlice, kinematics and header
128   AliRunLoader* runLoader = AliRunLoader::Open(filename);
129   if (!runLoader) {
130     Error("MUONmass_ESD", "getting run loader from file %s failed", 
131             filename);
132     return kFALSE;
133   }
134
135   if (!gAlice) {
136     Error("MUONmass_ESD", "no galice object found");
137     return kFALSE;
138   }
139   
140
141   // open the ESD file
142   TFile* esdFile = TFile::Open(esdFileName);
143   if (!esdFile || !esdFile->IsOpen()) {
144     Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
145     return kFALSE;
146   }
147   
148   AliESD* esd = new AliESD();
149   TTree* tree = (TTree*) esdFile->Get("esdTree");
150   if (!tree) {
151     Error("CheckESD", "no ESD tree found");
152     return kFALSE;
153   }
154   tree->SetBranchAddress("ESD", &esd);
155   
156   
157   AliESDVertex* Vertex = (AliESDVertex*) esd->AliESD::GetVertex();
158
159   runLoader->LoadHeader();
160   nevents = runLoader->GetNumberOfEvents();
161         
162   // Loop over events
163   for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) {
164
165     // get current event
166     runLoader->GetEvent(iEvent);
167    
168     // get the event summary data
169     tree->GetEvent(iEvent);
170     if (!esd) {
171       Error("CheckESD", "no ESD object found for event %d", iEvent);
172       return kFALSE;
173     }
174
175     // get the SPD reconstructed vertex (vertexer) and fill the histogram
176     fZVertex = Vertex->GetZv();
177     hPrimaryVertex->Fill(fZVertex);
178
179     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
180
181     //    printf("\n Nb of events analysed: %d\r",iEvent);
182     //      cout << " number of tracks: " << nTracks  <<endl;
183   
184     // loop over all reconstructed tracks (also first track of combination)
185     for (Int_t iTrack = 0; iTrack <  nTracks;  iTrack++) {
186
187       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
188
189       thetaX = muonTrack->GetThetaX();
190       thetaY = muonTrack->GetThetaY();
191
192       pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
193       fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
194       fPxRec1  = fPzRec1 * TMath::Tan(thetaX);
195       fPyRec1  = fPzRec1 * TMath::Tan(thetaY);
196       fCharge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
197
198       fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
199       fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
200
201       ntrackhits = muonTrack->GetNHit();
202       fitfmin    = muonTrack->GetChi2();
203
204       // transverse momentum
205       Float_t pt1 = fV1.Pt();
206
207       // total momentum
208       Float_t p1 = fV1.P();
209
210       // Rapidity
211       Float_t rapMuon1 = fV1.Rapidity();
212
213       // chi2 per d.o.f.
214       Float_t ch1 =  fitfmin / (2.0 * ntrackhits - 5);
215 //      printf(" px %f py %f pz %f NHits %d  Norm.chi2 %f charge %d\n", 
216 //           fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge);
217
218       // condition for good track (Chi2Cut and PtCut)
219
220       if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
221
222         // fill histos hPtMuon and hChi2PerDof
223         hPtMuon->Fill(pt1);
224         hPMuon->Fill(p1);
225         hChi2PerDof->Fill(ch1);
226         hRapMuon->Fill(rapMuon1);
227         if (fCharge > 0) {
228           hPtMuonPlus->Fill(pt1);
229           hThetaPhiPlus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
230         } else {
231           hPtMuonMinus->Fill(pt1);
232           hThetaPhiMinus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
233         }
234         // loop over second track of combination
235         for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
236           
237           AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack2);
238
239           thetaX = muonTrack->GetThetaX();
240           thetaY = muonTrack->GetThetaY();
241
242           pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
243           fPzRec2  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
244           fPxRec2  = fPzRec2 * TMath::Tan(thetaX);
245           fPyRec2  = fPzRec2 * TMath::Tan(thetaY);
246           fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
247
248           fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
249           fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
250
251           ntrackhits = muonTrack->GetNHit();
252           fitfmin    = muonTrack->GetChi2();
253
254           // transverse momentum
255           Float_t pt2 = fV2.Pt();
256
257           // chi2 per d.o.f.
258           Float_t ch2 = fitfmin  / (2.0 * ntrackhits - 5);
259
260           // condition for good track (Chi2Cut and PtCut)
261           if ((ch2 < Chi2Cut) && (pt2 > PtCutMin)  && (pt2 < PtCutMax)) {
262
263             // condition for opposite charges
264             if ((fCharge * fCharge2) == -1) {
265
266               // invariant mass
267               fVtot = fV1 + fV2;
268               Float_t invMass = fVtot.M();
269                     
270               // fill histos hInvMassAll and hInvMassRes
271               hInvMassAll->Fill(invMass);
272               hInvMassRes->Fill(invMass);
273               hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
274               Int_t ptTrig;
275               if (ResType == 553) 
276                 ptTrig =  0x400;// mask for Hpt unlike sign pair
277               else 
278                 ptTrig =  0x200;// mask for Lpt unlike sign pair
279
280               if (esd->GetTrigger() &  ptTrig) NbTrigger++; 
281               if (invMass > massMin && invMass < massMax) {
282                 EventInMass++;
283                 if (muonTrack->GetMatchTrigger() && (esd->GetTrigger() & ptTrig))// match with trigger
284                   EventInMassMatch++;
285
286                 hRapResonance->Fill(fVtot.Rapidity());
287                 hPtResonance->Fill(fVtot.Pt());
288               }
289
290             } // if (fCharge * fCharge2) == -1)
291           } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
292         } //  for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
293       } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
294     } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
295
296     hNumberOfTrack->Fill(nTracks);
297     //    esdFile->Delete();
298   } // for (Int_t iEvent = FirstEvent;
299
300 // Loop over events for bg event
301
302   Double_t thetaPlus,  phiPlus;
303   Double_t thetaMinus, phiMinus;
304   Float_t PtMinus, PtPlus;
305   
306   for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
307
308     hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
309     hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
310     PtPlus = hPtMuonPlus->GetRandom();
311     PtMinus = hPtMuonMinus->GetRandom();
312
313     fPxRec1  = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
314     fPyRec1  = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
315     fPzRec1  = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
316
317     fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
318     fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
319
320     fPxRec2  = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
321     fPyRec2  = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
322     fPzRec2  = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
323
324     fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
325     fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
326
327     // invariant mass
328     fVtot = fV1 + fV2;
329       
330     // fill histos hInvMassAll and hInvMassRes
331     hInvMassBg->Fill(fVtot.M());
332     hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
333   }
334
335   histoFile->Write();
336   histoFile->Close();
337
338   cout << endl;
339   cout << "EventInMass " << EventInMass << endl;
340   cout << "NbTrigger " << NbTrigger << endl;
341   cout << "EventInMass match with trigger " << EventInMassMatch << endl;
342
343   return kTRUE;
344 }
345