]>
Commit | Line | Data |
---|---|---|
5473f16a | 1 | #if !defined(__CINT__) || defined(__MAKECINT__) |
2 | // ROOT includes | |
48d6b312 | 3 | #include "TTree.h" |
5473f16a | 4 | #include "TBranch.h" |
5 | #include "TClonesArray.h" | |
6 | #include "TLorentzVector.h" | |
7 | #include "TFile.h" | |
8 | #include "TH1.h" | |
ac3c5325 | 9 | #include "TH2.h" |
5473f16a | 10 | #include "TParticle.h" |
11 | #include "TTree.h" | |
12 | #include <Riostream.h> | |
8cde4af5 | 13 | #include <TGeoManager.h> |
21939432 | 14 | #include <TROOT.h> |
5473f16a | 15 | |
16 | // STEER includes | |
0ff94351 | 17 | #include "AliLog.h" |
0ff94351 | 18 | #include "AliMagFMaps.h" |
1c2bdf00 | 19 | #include "AliESDEvent.h" |
995a61aa | 20 | #include "AliESDVertex.h" |
0ff94351 | 21 | #include "AliTracker.h" |
103e6575 | 22 | #include "AliESDMuonTrack.h" |
5473f16a | 23 | |
24 | // MUON includes | |
211c52eb | 25 | #include "AliMUONTrackParam.h" |
37827b29 | 26 | #include "AliMUONTrackExtrap.h" |
103e6575 | 27 | #include "AliMUONESDInterface.h" |
5473f16a | 28 | #endif |
5473f16a | 29 | |
e54bf126 | 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 | |
5473f16a | 45 | |
cec5541f | 46 | Bool_t MUONmassPlot(Int_t ExtrapToVertex = -1, char* geoFilename = "geometry.root", |
8cde4af5 | 47 | Int_t FirstEvent = 0, Int_t LastEvent = 10000, char* esdFileName = "AliESDs.root", Int_t ResType = 553, |
5473f16a | 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 | { | |
e54bf126 | 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 | ||
5473f16a | 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 | ||
5473f16a | 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.); | |
ac3c5325 | 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.); | |
5473f16a | 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.); | |
ac3c5325 | 88 | TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.); |
f57d136a | 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; | |
6ab68e5f | 92 | TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15); |
5473f16a | 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.); | |
ac3c5325 | 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.); | |
5473f16a | 106 | |
107 | ||
108 | // settings | |
109 | Int_t EventInMass = 0; | |
6678cd54 | 110 | Int_t EventInMassMatch = 0; |
111 | Int_t NbTrigger = 0; | |
112 | ||
5473f16a | 113 | Float_t muonMass = 0.105658389; |
114 | // Float_t UpsilonMass = 9.46037; | |
115 | // Float_t JPsiMass = 3.097; | |
116 | ||
22ccc301 | 117 | Int_t fCharge1, fCharge2; |
5473f16a | 118 | Double_t fPxRec1, fPyRec1, fPzRec1, fE1; |
119 | Double_t fPxRec2, fPyRec2, fPzRec2, fE2; | |
5473f16a | 120 | |
121 | Int_t ntrackhits, nevents; | |
122 | Double_t fitfmin; | |
e3c6ae4c | 123 | Double_t fZVertex=0; |
211c52eb | 124 | Double_t fYVertex=0; |
125 | Double_t fXVertex=0; | |
690d2205 | 126 | Double_t errXVtx=0; |
127 | Double_t errYVtx=0; | |
5473f16a | 128 | |
129 | TLorentzVector fV1, fV2, fVtot; | |
cb75342e | 130 | |
8cde4af5 | 131 | // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex) |
132 | if (!gGeoManager) { | |
133 | TGeoManager::Import(geoFilename); | |
134 | if (!gGeoManager) { | |
cec5541f | 135 | Error("MUONmass_ESD", "getting geometry from file %s failed", geoFilename); |
8cde4af5 | 136 | return kFALSE; |
137 | } | |
138 | } | |
139 | ||
0ff94351 | 140 | // set mag field |
141 | // waiting for mag field in CDB | |
142 | printf("Loading field map...\n"); | |
b97b210c | 143 | AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG); |
0ff94351 | 144 | AliTracker::SetFieldMap(field, kFALSE); |
cb75342e | 145 | |
5473f16a | 146 | // open the ESD file |
147 | TFile* esdFile = TFile::Open(esdFileName); | |
148 | if (!esdFile || !esdFile->IsOpen()) { | |
149 | Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName); | |
150 | return kFALSE; | |
151 | } | |
152 | ||
1c2bdf00 | 153 | AliESDEvent* esd = new AliESDEvent(); |
48d6b312 | 154 | TTree* tree = (TTree*) esdFile->Get("esdTree"); |
155 | if (!tree) { | |
156 | Error("CheckESD", "no ESD tree found"); | |
157 | return kFALSE; | |
158 | } | |
8b4a9b89 | 159 | // tree->SetBranchAddress("ESD", &esd); |
160 | esd->ReadFromTree(tree); | |
f9ebb3bd | 161 | |
cec5541f | 162 | nevents = (Int_t)tree->GetEntries(); |
211c52eb | 163 | |
164 | AliMUONTrackParam trackParam; | |
0ff94351 | 165 | |
5473f16a | 166 | // Loop over events |
167 | for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { | |
168 | ||
5473f16a | 169 | // get the event summary data |
48d6b312 | 170 | tree->GetEvent(iEvent); |
5473f16a | 171 | if (!esd) { |
48d6b312 | 172 | Error("CheckESD", "no ESD object found for event %d", iEvent); |
5473f16a | 173 | return kFALSE; |
174 | } | |
175 | ||
f9ebb3bd | 176 | // get the SPD reconstructed vertex (vertexer) and fill the histogram |
43939bd8 | 177 | AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex(); |
43939bd8 | 178 | if (Vertex->GetNContributors()) { |
211c52eb | 179 | fZVertex = Vertex->GetZv(); |
180 | fYVertex = Vertex->GetYv(); | |
181 | fXVertex = Vertex->GetXv(); | |
690d2205 | 182 | errXVtx = Vertex->GetXRes(); |
183 | errYVtx = Vertex->GetYRes(); | |
211c52eb | 184 | } |
f9ebb3bd | 185 | hPrimaryVertex->Fill(fZVertex); |
186 | ||
48d6b312 | 187 | Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; |
5473f16a | 188 | |
189 | // printf("\n Nb of events analysed: %d\r",iEvent); | |
48d6b312 | 190 | // cout << " number of tracks: " << nTracks <<endl; |
5473f16a | 191 | |
37827b29 | 192 | // set the magnetic field for track extrapolations |
193 | AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap()); | |
5473f16a | 194 | // loop over all reconstructed tracks (also first track of combination) |
195 | for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) { | |
196 | ||
b1fea02e | 197 | // skip ghosts |
198 | if (!esd->GetMuonTrack(iTrack)->ContainTrackerData()) continue; | |
199 | ||
8cde4af5 | 200 | AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack))); |
5473f16a | 201 | |
8cde4af5 | 202 | // extrapolate to vertex if required and available |
203 | if (ExtrapToVertex > 0 && Vertex->GetNContributors()) { | |
103e6575 | 204 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam); |
690d2205 | 205 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx); |
103e6575 | 206 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack |
8cde4af5 | 207 | } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){ |
103e6575 | 208 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam); |
690d2205 | 209 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.); |
103e6575 | 210 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack |
211c52eb | 211 | } |
8cde4af5 | 212 | |
22ccc301 | 213 | fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum())); |
214 | ||
215 | muonTrack->LorentzP(fV1); | |
8cde4af5 | 216 | |
5473f16a | 217 | ntrackhits = muonTrack->GetNHit(); |
218 | fitfmin = muonTrack->GetChi2(); | |
219 | ||
220 | // transverse momentum | |
221 | Float_t pt1 = fV1.Pt(); | |
222 | ||
223 | // total momentum | |
224 | Float_t p1 = fV1.P(); | |
225 | ||
226 | // Rapidity | |
227 | Float_t rapMuon1 = fV1.Rapidity(); | |
228 | ||
229 | // chi2 per d.o.f. | |
230 | Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5); | |
231 | // printf(" px %f py %f pz %f NHits %d Norm.chi2 %f charge %d\n", | |
22ccc301 | 232 | // fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge1); |
5473f16a | 233 | |
234 | // condition for good track (Chi2Cut and PtCut) | |
235 | ||
236 | if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { | |
237 | ||
238 | // fill histos hPtMuon and hChi2PerDof | |
239 | hPtMuon->Fill(pt1); | |
240 | hPMuon->Fill(p1); | |
241 | hChi2PerDof->Fill(ch1); | |
242 | hRapMuon->Fill(rapMuon1); | |
22ccc301 | 243 | if (fCharge1 > 0) { |
ac3c5325 | 244 | hPtMuonPlus->Fill(pt1); |
22ccc301 | 245 | hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi()); |
ac3c5325 | 246 | } else { |
247 | hPtMuonMinus->Fill(pt1); | |
22ccc301 | 248 | hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi()); |
ac3c5325 | 249 | } |
5473f16a | 250 | // loop over second track of combination |
251 | for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) { | |
252 | ||
b1fea02e | 253 | // skip ghosts |
254 | if (!esd->GetMuonTrack(iTrack2)->ContainTrackerData()) continue; | |
255 | ||
8cde4af5 | 256 | AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2))); |
257 | ||
258 | // extrapolate to vertex if required and available | |
259 | if (ExtrapToVertex > 0 && Vertex->GetNContributors()) { | |
103e6575 | 260 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam); |
690d2205 | 261 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx); |
103e6575 | 262 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack |
8cde4af5 | 263 | } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){ |
103e6575 | 264 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam); |
690d2205 | 265 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.); |
103e6575 | 266 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack |
211c52eb | 267 | } |
8cde4af5 | 268 | |
8cde4af5 | 269 | fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum())); |
5473f16a | 270 | |
22ccc301 | 271 | muonTrack2->LorentzP(fV2); |
5473f16a | 272 | |
8cde4af5 | 273 | ntrackhits = muonTrack2->GetNHit(); |
274 | fitfmin = muonTrack2->GetChi2(); | |
5473f16a | 275 | |
276 | // transverse momentum | |
277 | Float_t pt2 = fV2.Pt(); | |
278 | ||
279 | // chi2 per d.o.f. | |
280 | Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5); | |
281 | ||
282 | // condition for good track (Chi2Cut and PtCut) | |
283 | if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) { | |
284 | ||
285 | // condition for opposite charges | |
22ccc301 | 286 | if ((fCharge1 * fCharge2) == -1) { |
5473f16a | 287 | |
288 | // invariant mass | |
289 | fVtot = fV1 + fV2; | |
290 | Float_t invMass = fVtot.M(); | |
291 | ||
292 | // fill histos hInvMassAll and hInvMassRes | |
293 | hInvMassAll->Fill(invMass); | |
294 | hInvMassRes->Fill(invMass); | |
ac3c5325 | 295 | hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt()); |
6678cd54 | 296 | Int_t ptTrig; |
297 | if (ResType == 553) | |
8d4fefab | 298 | ptTrig = 0x20;// mask for Hpt unlike sign pair |
6678cd54 | 299 | else |
8d4fefab | 300 | ptTrig = 0x10;// mask for Lpt unlike sign pair |
6678cd54 | 301 | |
f57d136a | 302 | if (esd->GetTriggerMask() & ptTrig) NbTrigger++; |
5473f16a | 303 | if (invMass > massMin && invMass < massMax) { |
304 | EventInMass++; | |
8cde4af5 | 305 | if (muonTrack2->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger |
6678cd54 | 306 | EventInMassMatch++; |
307 | ||
5473f16a | 308 | hRapResonance->Fill(fVtot.Rapidity()); |
309 | hPtResonance->Fill(fVtot.Pt()); | |
310 | } | |
311 | ||
22ccc301 | 312 | } // if (fCharge1 * fCharge2) == -1) |
5473f16a | 313 | } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) |
8cde4af5 | 314 | delete muonTrack2; |
5473f16a | 315 | } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++) |
316 | } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) ) | |
8cde4af5 | 317 | delete muonTrack; |
5473f16a | 318 | } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++) |
319 | ||
320 | hNumberOfTrack->Fill(nTracks); | |
48d6b312 | 321 | // esdFile->Delete(); |
5473f16a | 322 | } // for (Int_t iEvent = FirstEvent; |
323 | ||
ac3c5325 | 324 | // Loop over events for bg event |
325 | ||
326 | Double_t thetaPlus, phiPlus; | |
327 | Double_t thetaMinus, phiMinus; | |
328 | Float_t PtMinus, PtPlus; | |
329 | ||
330 | for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { | |
331 | ||
332 | hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus); | |
333 | hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus); | |
334 | PtPlus = hPtMuonPlus->GetRandom(); | |
335 | PtMinus = hPtMuonMinus->GetRandom(); | |
336 | ||
337 | fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus); | |
338 | fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus); | |
339 | fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus); | |
340 | ||
341 | fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1); | |
342 | fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1); | |
343 | ||
344 | fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus); | |
345 | fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus); | |
346 | fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus); | |
347 | ||
348 | fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2); | |
349 | fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2); | |
350 | ||
351 | // invariant mass | |
352 | fVtot = fV1 + fV2; | |
353 | ||
354 | // fill histos hInvMassAll and hInvMassRes | |
355 | hInvMassBg->Fill(fVtot.M()); | |
356 | hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() ); | |
357 | } | |
358 | ||
5473f16a | 359 | histoFile->Write(); |
360 | histoFile->Close(); | |
361 | ||
6678cd54 | 362 | cout << endl; |
5473f16a | 363 | cout << "EventInMass " << EventInMass << endl; |
6678cd54 | 364 | cout << "NbTrigger " << NbTrigger << endl; |
365 | cout << "EventInMass match with trigger " << EventInMassMatch << endl; | |
5473f16a | 366 | |
367 | return kTRUE; | |
368 | } | |
369 |