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