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