]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONmassPlot_ESD.C
readers updated (mini header -> data header)
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot_ESD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 // ROOT includes
3 #include "TBranch.h"
4 #include "TClonesArray.h"
5 #include "TLorentzVector.h"
6 #include "TFile.h"
7 #include "TH1.h"
8 #include "TParticle.h"
9 #include "TTree.h"
10 #include <Riostream.h>
11
12 // STEER includes
13 #include "AliRun.h"
14 #include "AliRunLoader.h"
15 #include "AliHeader.h"
16 #include "AliLoader.h"
17 #include "AliStack.h"
18 #include "AliESD.h"
19
20 // MUON includes
21 #include "AliESDMuonTrack.h"
22 #endif
23 //
24 // Macro MUONmassPlot.C for ESD
25 // Ch. Finck, Subatech, April. 2004
26 //
27
28 // macro to make invariant mass plots
29 // for combinations of 2 muons with opposite charges,
30 // from root file "MUON.tracks.root" containing the result of track reconstruction.
31 // Histograms are stored on the "MUONmassPlot.root" file.
32 // introducing TLorentzVector for parameter calculations (Pt, P,rap,etc...)
33 // using Invariant Mass for rapidity.
34
35 // Arguments:
36 //   FirstEvent (default 0)
37 //   LastEvent (default 0)
38 //   ResType (default 553)
39 //      553 for Upsilon, anything else for J/Psi
40 //   Chi2Cut (default 100)
41 //      to keep only tracks with chi2 per d.o.f. < Chi2Cut
42 //   PtCutMin (default 1)
43 //      to keep only tracks with transverse momentum > PtCutMin
44 //   PtCutMax (default 10000)
45 //      to keep only tracks with transverse momentum < PtCutMax
46 //   massMin (default 9.17 for Upsilon) 
47 //      &  massMax (default 9.77 for Upsilon) 
48 //         to calculate the reconstruction efficiency for resonances with invariant mass
49 //         massMin < mass < massMax.
50
51 // Add parameters and histograms for analysis 
52
53 Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t LastEvent = 10000,
54                   char* esdFileName = "AliESDs.root", Int_t ResType = 553, 
55                   Float_t Chi2Cut = 100., Float_t PtCutMin = 1., Float_t PtCutMax = 10000.,
56                   Float_t massMin = 9.17,Float_t massMax = 9.77)
57 {
58   cout << "MUONmassPlot " << endl;
59   cout << "FirstEvent " << FirstEvent << endl;
60   cout << "LastEvent " << LastEvent << endl;
61   cout << "ResType " << ResType << endl;
62   cout << "Chi2Cut " << Chi2Cut << endl;
63   cout << "PtCutMin " << PtCutMin << endl;
64   cout << "PtCutMax " << PtCutMax << endl;
65   cout << "massMin " << massMin << endl;
66   cout << "massMax " << massMax << endl;
67
68  
69   //Reset ROOT and connect tree file
70   gROOT->Reset();
71
72
73   // File for histograms and histogram booking
74   TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
75   TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
76   TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
77   TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
78   TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
79   TH1F *hInvMassRes;
80
81   if (ResType == 553) {
82     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
83   } else {
84     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
85   }
86
87   TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
88   TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
89   TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
90   TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
91
92
93   // settings
94   Int_t EventInMass = 0;
95   Float_t muonMass = 0.105658389;
96 //   Float_t UpsilonMass = 9.46037;
97 //   Float_t JPsiMass = 3.097;
98
99   Double_t thetaX, thetaY, pYZ;
100   Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
101   Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
102   Int_t fCharge, fCharge2;
103
104   Int_t ntrackhits, nevents;
105   Double_t fitfmin;
106
107  
108   TLorentzVector fV1, fV2, fVtot;
109   
110   // open run loader and load gAlice, kinematics and header
111   AliRunLoader* runLoader = AliRunLoader::Open(filename);
112   if (!runLoader) {
113     Error("MUONmass_ESD", "getting run loader from file %s failed", 
114             filename);
115     return kFALSE;
116   }
117
118   runLoader->LoadgAlice();
119   gAlice = runLoader->GetAliRun();
120   if (!gAlice) {
121     Error("MUONmass_ESD", "no galice object found");
122     return kFALSE;
123   }
124   
125
126   // open the ESD file
127   TFile* esdFile = TFile::Open(esdFileName);
128   if (!esdFile || !esdFile->IsOpen()) {
129     Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
130     return kFALSE;
131   }
132   
133   runLoader->LoadHeader();
134   nevents = runLoader->GetNumberOfEvents();
135         
136   // Loop over events
137   for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) {
138
139     // get current event
140     runLoader->GetEvent(iEvent);
141     
142     // get the event summary data
143     char esdName[256]; 
144     sprintf(esdName, "ESD%d", iEvent);
145     AliESD* esd = (AliESD*) esdFile->Get(esdName);
146     if (!esd) {
147       Error("MUONmass_ESD", "no ESD object found for event %d", iEvent);
148       return kFALSE;
149     }
150
151     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; //
152
153     //    printf("\n Nb of events analysed: %d\r",iEvent);
154     //   cout << " number of tracks: " << nrectracks  <<endl;
155   
156     // loop over all reconstructed tracks (also first track of combination)
157     for (Int_t iTrack = 0; iTrack <  nTracks;  iTrack++) {
158
159       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
160
161       thetaX = muonTrack->GetThetaX();
162       thetaY = muonTrack->GetThetaY();
163
164       pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
165       fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
166       fPxRec1  = fPzRec1 * TMath::Tan(thetaX);
167       fPyRec1  = fPzRec1 * TMath::Tan(thetaY);
168       fCharge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
169
170       fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
171       fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
172
173       ntrackhits = muonTrack->GetNHit();
174       fitfmin    = muonTrack->GetChi2();
175
176       // transverse momentum
177       Float_t pt1 = fV1.Pt();
178
179       // total momentum
180       Float_t p1 = fV1.P();
181
182       // Rapidity
183       Float_t rapMuon1 = fV1.Rapidity();
184
185       // chi2 per d.o.f.
186       Float_t ch1 =  fitfmin / (2.0 * ntrackhits - 5);
187 //      printf(" px %f py %f pz %f NHits %d  Norm.chi2 %f charge %d\n", 
188 //           fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge);
189
190       // condition for good track (Chi2Cut and PtCut)
191
192       if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
193
194         // fill histos hPtMuon and hChi2PerDof
195         hPtMuon->Fill(pt1);
196         hPMuon->Fill(p1);
197         hChi2PerDof->Fill(ch1);
198         hRapMuon->Fill(rapMuon1);
199
200         // loop over second track of combination
201         for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
202           
203           AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack2);
204
205           thetaX = muonTrack->GetThetaX();
206           thetaY = muonTrack->GetThetaY();
207
208           pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
209           fPzRec2  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
210           fPxRec2  = fPzRec2 * TMath::Tan(thetaX);
211           fPyRec2  = fPzRec2 * TMath::Tan(thetaY);
212           fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
213
214           fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
215           fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
216
217           ntrackhits = muonTrack->GetNHit();
218           fitfmin    = muonTrack->GetChi2();
219
220           // transverse momentum
221           Float_t pt2 = fV2.Pt();
222
223           // chi2 per d.o.f.
224           Float_t ch2 = fitfmin  / (2.0 * ntrackhits - 5);
225
226           // condition for good track (Chi2Cut and PtCut)
227           if ((ch2 < Chi2Cut) && (pt2 > PtCutMin)  && (pt2 < PtCutMax)) {
228
229             // condition for opposite charges
230             if ((fCharge * fCharge2) == -1) {
231
232               // invariant mass
233               fVtot = fV1 + fV2;
234               Float_t invMass = fVtot.M();
235                     
236               // fill histos hInvMassAll and hInvMassRes
237               hInvMassAll->Fill(invMass);
238               hInvMassRes->Fill(invMass);
239
240               if (invMass > massMin && invMass < massMax) {
241                 EventInMass++;
242                 hRapResonance->Fill(fVtot.Rapidity());
243                 hPtResonance->Fill(fVtot.Pt());
244               }
245
246             } // if (fCharge * fCharge2) == -1)
247           } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
248         } //  for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
249       } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
250     } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
251
252     hNumberOfTrack->Fill(nTracks);
253   } // for (Int_t iEvent = FirstEvent;
254
255   histoFile->Write();
256   histoFile->Close();
257
258   cout << "MUONmassPlot " << endl;
259   cout << "FirstEvent " << FirstEvent << endl;
260   cout << "LastEvent " << LastEvent << endl;
261   cout << "ResType " << ResType << endl;
262   cout << "Chi2Cut " << Chi2Cut << endl;
263   cout << "PtCutMin " << PtCutMin << endl;
264   cout << "PtCutMax " << PtCutMax << endl;
265   cout << "massMin " << massMin << endl;
266   cout << "massMax " << massMax << endl;
267   cout << "EventInMass " << EventInMass << endl;
268
269   return kTRUE;
270 }
271