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