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