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