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