]>
Commit | Line | Data |
---|---|---|
1 | /************************************************************************** | |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | /* $Id$ */ | |
17 | ||
18 | /// \ingroup macros | |
19 | /// \file MUONefficiency.C | |
20 | /// \brief add brief description | |
21 | /// | |
22 | /// Macro (upgraded version of MUONmassPlot_ESD.C, better handling of Jpsi) to make : | |
23 | /// - Ntuple (Ktuple) containing Upsilon kinematics variables (from kinematics.root files) | |
24 | /// - Ntuple (ESDtuple) containing Upsilon kinematics variables from reconstruction and | |
25 | /// combinations of 2 muons with opposite charges (ESDtupleBck will be used later) | |
26 | /// - Some QA histograms | |
27 | /// Ntuple are stored in the file MUONefficiency.root and ESD tree and QA histograms in AliESDs.root | |
28 | /// | |
29 | /// \author Christophe Suire, IPN Orsay | |
30 | ||
31 | ||
32 | ||
33 | #if !defined(__CINT__) || defined(__MAKECINT__) | |
34 | ||
35 | // MUON includes | |
36 | #include "AliMUONTrackParam.h" | |
37 | #include "AliMUONTrackExtrap.h" | |
38 | #include "AliMUONESDInterface.h" | |
39 | ||
40 | // STEER includes | |
41 | #include "AliRun.h" | |
42 | #include "AliRunLoader.h" | |
43 | #include "AliHeader.h" | |
44 | #include "AliLoader.h" | |
45 | #include "AliStack.h" | |
46 | #include "AliMagFMaps.h" | |
47 | #include "AliESDEvent.h" | |
48 | #include "AliESDVertex.h" | |
49 | #include "AliTracker.h" | |
50 | #include "AliCDBManager.h" | |
51 | #include "AliESDMuonTrack.h" | |
52 | ||
53 | // ROOT includes | |
54 | #include "TTree.h" | |
55 | #include "TNtuple.h" | |
56 | #include "TBranch.h" | |
57 | #include "TClonesArray.h" | |
58 | #include "TLorentzVector.h" | |
59 | #include "TFile.h" | |
60 | #include "TH1.h" | |
61 | #include "TH2.h" | |
62 | #include "TParticle.h" | |
63 | #include "TTree.h" | |
64 | #include "TString.h" | |
65 | #include <Riostream.h> | |
66 | #include <TGeoManager.h> | |
67 | #include <TROOT.h> | |
68 | ||
69 | #endif | |
70 | ||
71 | Bool_t MUONefficiency( char* filename = "galice.root", char* geoFilename = "geometry.root", char* esdFileName = "AliESDs.root", | |
72 | Int_t ExtrapToVertex = -1, Int_t ResType = 553, Int_t FirstEvent = 0, Int_t LastEvent = 1000000 ) | |
73 | { | |
74 | /// \param ExtrapToVertex (default -1) | |
75 | /// - <0: no extrapolation; | |
76 | /// - =0: extrapolation to (0,0,0); | |
77 | /// - >0: extrapolation to ESDVertex if available, else to (0,0,0) | |
78 | /// \param ResType 553 for Upsilon, 443 for J/Psi (default 553) | |
79 | /// \param FirstEvent (default 0) | |
80 | /// \param LastEvent (default 1.e6) | |
81 | /// \param Chi2Cut to keep only tracks with chi2 per d.o.f. < Chi2Cut (default 100) | |
82 | ||
83 | ||
84 | // MUONefficiency starts | |
85 | ||
86 | // Set default CDB storage | |
87 | AliCDBManager* man = AliCDBManager::Instance(); | |
88 | man->SetDefaultStorage("local://$ALICE_ROOT"); | |
89 | ||
90 | Double_t MUON_MASS = 0.105658369; | |
91 | Double_t UPSILON_MASS = 9.4603 ; | |
92 | Double_t JPSI_MASS = 3.097; | |
93 | ||
94 | // Upper and lower bound for counting entries in the mass peak | |
95 | // +/- 300 MeV/c^2 in this case. | |
96 | Float_t countingRange = 0.300 ; | |
97 | ||
98 | Float_t massResonance = 5.; | |
99 | Float_t invMassMinInPeak = 0. ; | |
100 | Float_t invMassMaxInPeak = 0. ; | |
101 | ||
102 | Float_t nBinsPerGev = 40 ; | |
103 | Float_t invMassMin = 0; Float_t invMassMax = 20; | |
104 | Float_t ptMinResonance = 0 ; Float_t ptMaxResonance = 20 ; Int_t ptBinsResonance = 100; | |
105 | ||
106 | if (ResType==443) { | |
107 | massResonance = JPSI_MASS ; | |
108 | invMassMinInPeak = JPSI_MASS - countingRange ; invMassMaxInPeak = JPSI_MASS + countingRange ; | |
109 | //limits for histograms | |
110 | invMassMin = 0 ; invMassMax = 6.; | |
111 | ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100; | |
112 | } | |
113 | if (ResType==553) { | |
114 | massResonance = UPSILON_MASS; | |
115 | invMassMinInPeak = UPSILON_MASS - countingRange ; invMassMaxInPeak = UPSILON_MASS + countingRange; | |
116 | //limits for histograms | |
117 | invMassMin = 0 ; invMassMax = 12.; | |
118 | ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100; | |
119 | } | |
120 | ||
121 | // Single Tracks muon cuts | |
122 | Float_t Chi2Cut = 100.; | |
123 | Float_t PtCutMin = 0. ; | |
124 | Float_t PtCutMax = 10000. ; | |
125 | ||
126 | ||
127 | // Limits for histograms | |
128 | Float_t ptMinMuon = 0. ; Float_t ptMaxMuon = 20.; Int_t ptBinsMuon = 100 ; | |
129 | Float_t pMinMuon = 0. ; Float_t pMaxMuon = 200.; Int_t pBinsMuon = 100 ; | |
130 | ||
131 | ||
132 | //Reset ROOT and connect tree file | |
133 | gROOT->Reset(); | |
134 | ||
135 | // Printing Level | |
136 | Int_t PRINTLEVEL = 0 ; | |
137 | ||
138 | //for kinematic, i.e. reference tracks | |
139 | TNtuple *Ktuple = new TNtuple("Ktuple","Kinematics NTuple","ev:npart:id:idmo:idgdmo:p:pt:y:theta:pseudorap:vx:vy:vz"); | |
140 | ||
141 | //for reconstruction | |
142 | TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon); | |
143 | TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon); | |
144 | TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon); | |
145 | TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", pBinsMuon, pMinMuon, pMaxMuon); | |
146 | ||
147 | TH1F *hInvMassAll; | |
148 | TH1F *hInvMassBg; | |
149 | TH2F *hInvMassAll_vs_Pt; | |
150 | TH2F *hInvMassBgk_vs_Pt; | |
151 | TH1F *hInvMassRes; | |
152 | ||
153 | ||
154 | hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax - invMassMin)), invMassMin, invMassMax); | |
155 | hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax); | |
156 | hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance); | |
157 | hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance); | |
158 | ||
159 | hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Resonance",(Int_t) (nBinsPerGev*3*countingRange*2),massResonance-3*countingRange,massResonance+3*countingRange); | |
160 | ||
161 | TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15); | |
162 | TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.); | |
163 | TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5); | |
164 | TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2); | |
165 | TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2); | |
166 | TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.); | |
167 | TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.); | |
168 | TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.); | |
169 | ||
170 | TNtuple *ESDtuple = new TNtuple("ESDtuple","Reconstructed Mu+Mu- pairs and Upsilon","ev:tw:pt:y:theta:minv:pt1:y1:theta1:q1:trig1:pt2:y2:theta2:q2:trig2"); | |
171 | TNtuple *ESDtupleBck = new TNtuple("ESDtupleBck","Reconstructed Mu+Mu- pairs for Background","ev:pt:y:theta:minv:pt1:y1:theta1:pt2:y2:theta2"); | |
172 | ||
173 | ||
174 | // Variables | |
175 | Int_t EventInMass = 0; | |
176 | Int_t EventInMassMatch = 0; | |
177 | Int_t NbTrigger = 0; | |
178 | Int_t ptTrig = 0; | |
179 | ||
180 | Double_t fXVertex=0; | |
181 | Double_t fYVertex=0; | |
182 | Double_t fZVertex=0; | |
183 | Double_t errXVtx=0; | |
184 | Double_t errYVtx=0; | |
185 | ||
186 | Double_t fPxRec1, fPyRec1, fPzRec1, fE1; | |
187 | Double_t fPxRec2, fPyRec2, fPzRec2, fE2; | |
188 | Int_t fCharge1, fCharge2; | |
189 | ||
190 | Int_t ntrackhits, nevents; | |
191 | Int_t nprocessedevents = 0 ; | |
192 | Double_t fitfmin; | |
193 | ||
194 | TLorentzVector fV1, fV2, fVtot; | |
195 | ||
196 | // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex) | |
197 | if (!gGeoManager) { | |
198 | TGeoManager::Import(geoFilename); | |
199 | if (!gGeoManager) { | |
200 | Error("MUONmass_ESD", "getting geometry from file %s failed", filename); | |
201 | return kFALSE; | |
202 | } | |
203 | } | |
204 | ||
205 | // set mag field | |
206 | // waiting for mag field in CDB | |
207 | printf("Loading field map...\n"); | |
208 | AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG); | |
209 | AliTracker::SetFieldMap(field, kFALSE); | |
210 | ||
211 | // open run loader and load gAlice, kinematics and header | |
212 | AliRunLoader* runLoader = AliRunLoader::Open(filename); | |
213 | if (!runLoader) { | |
214 | Error("MUONefficiency", "getting run loader from file %s failed", filename); | |
215 | return kFALSE; | |
216 | } | |
217 | ||
218 | runLoader->LoadgAlice(); | |
219 | gAlice = runLoader->GetAliRun(); | |
220 | if (!gAlice) { | |
221 | Error("MUONefficiency", "no galice object found"); | |
222 | return kFALSE; | |
223 | } | |
224 | ||
225 | // open the ESD file | |
226 | TFile* esdFile = TFile::Open(esdFileName); | |
227 | if (!esdFile || !esdFile->IsOpen()) { | |
228 | Error("MUONefficiency", "opening ESD file %s failed", esdFileName); | |
229 | return kFALSE; | |
230 | } | |
231 | ||
232 | AliESDEvent* esd = new AliESDEvent(); | |
233 | TTree* tree = (TTree*) esdFile->Get("esdTree"); | |
234 | if (!tree) { | |
235 | Error("CheckESD", "no ESD tree found"); | |
236 | return kFALSE; | |
237 | } | |
238 | esd->ReadFromTree(tree); | |
239 | ||
240 | runLoader->LoadHeader(); | |
241 | Int_t runNumber = runLoader->GetHeader()->GetRun(); | |
242 | AliCDBManager::Instance()->SetRun(runNumber); | |
243 | ||
244 | nevents = runLoader->GetNumberOfEvents(); | |
245 | AliMUONTrackParam trackParam; | |
246 | ||
247 | // to access the particle Stack | |
248 | runLoader->LoadKinematics("READ"); | |
249 | ||
250 | Int_t numberOfGeneratedResonances = 0 ; | |
251 | ||
252 | TParticle *particle; | |
253 | ||
254 | Int_t track1Trigger = 0 ; | |
255 | Float_t track1TriggerChi2 = 0 ; | |
256 | Int_t track2Trigger = 0 ; | |
257 | Float_t track2TriggerChi2 = 0 ; | |
258 | ||
259 | ||
260 | // Loop over events | |
261 | for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { // Start event loop | |
262 | ||
263 | if (iEvent%1000 == 0 ) | |
264 | printf("\n Nb of events analysed: %d \n",iEvent); | |
265 | ||
266 | // get current event | |
267 | runLoader->GetEvent(iEvent); | |
268 | nprocessedevents++; | |
269 | ||
270 | // get the stack and fill the kine tree | |
271 | AliStack *theStack = runLoader->Stack(); | |
272 | if (PRINTLEVEL > 0) theStack->DumpPStack (); | |
273 | ||
274 | Int_t nparticles = (Int_t)runLoader->TreeK()->GetEntries(); | |
275 | Int_t nprimarypart = theStack->GetNprimary(); | |
276 | Int_t ntracks = theStack->GetNtrack(); | |
277 | ||
278 | if (PRINTLEVEL || (iEvent%100==0)) printf("\n >>> Event %d \n",iEvent); | |
279 | if (PRINTLEVEL) cout << nprimarypart << " Particles generated (total is " << ntracks << ")"<< endl ; | |
280 | ||
281 | for(Int_t iparticle=0; iparticle<nparticles; iparticle++) { // Start loop over particles | |
282 | particle = theStack->Particle(iparticle); | |
283 | ||
284 | Int_t muId = particle->GetPdgCode(); | |
285 | Int_t muM = particle->GetFirstMother(); | |
286 | Int_t muGM = 0; | |
287 | Float_t muP = particle->P(); | |
288 | Float_t muPt = TMath::Sqrt(particle->Px()*particle->Px()+particle->Py()*particle->Py()); | |
289 | Float_t muY = 0.5*TMath::Log((particle->Energy()+particle->Pz()+1.e-13)/(particle->Energy()-particle->Pz()+1.e-13)); | |
290 | if (muM >= 0) { | |
291 | TParticle *theMum = theStack->Particle(muM); | |
292 | muM = theMum->GetPdgCode(); | |
293 | muGM = theMum->GetFirstMother() ; | |
294 | if (muGM >= 0){ | |
295 | TParticle *grandMa = theStack->Particle(muGM); | |
296 | muGM = grandMa->GetPdgCode(); | |
297 | } | |
298 | else muGM=0; | |
299 | } | |
300 | else muM=0; | |
301 | ||
302 | if (muId==ResType) numberOfGeneratedResonances++; | |
303 | ||
304 | ||
305 | Float_t muT = particle->Theta()*180/TMath::Pi(); | |
306 | Float_t muE = particle->Eta(); | |
307 | ||
308 | Float_t muVx = particle->Vx(); | |
309 | Float_t muVy = particle->Vy(); | |
310 | Float_t muVz = particle->Vz(); | |
311 | ||
312 | // If a write error occurs, the number of bytes returned is -1. | |
313 | // If no data are written, because e.g. the branch is disabled, | |
314 | // the number of bytes returned is 0. | |
315 | Int_t errCode = Ktuple->Fill(iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz); | |
316 | if (PRINTLEVEL || errCode < 1) printf("iEvent %d,nparticles %d,muId %d,muM %d,muGM %d,muP %.2f,muPt %.2f,muY %.2f,muT %.2f,muE %.2f,muVx %.2f,muVy %.2f,muVz %.2f \n", iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz); | |
317 | ||
318 | } // End loop over particles | |
319 | ||
320 | ||
321 | ||
322 | // get the event summary data | |
323 | tree->GetEvent(iEvent); | |
324 | if (!esd) { | |
325 | Error("CheckESD", "no ESD object found for event %d", iEvent); | |
326 | return kFALSE; | |
327 | } | |
328 | ||
329 | // get the SPD reconstructed vertex (vertexer) and fill the histogram | |
330 | AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex(); | |
331 | if (Vertex->GetNContributors()) { | |
332 | fZVertex = Vertex->GetZv(); | |
333 | fYVertex = Vertex->GetYv(); | |
334 | fXVertex = Vertex->GetXv(); | |
335 | errXVtx = Vertex->GetXRes(); | |
336 | errYVtx = Vertex->GetYRes(); | |
337 | } | |
338 | hPrimaryVertex->Fill(fZVertex); | |
339 | ||
340 | Int_t triggerWord = esd->GetTriggerMask(); | |
341 | Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; | |
342 | ||
343 | if (PRINTLEVEL > 0){ | |
344 | printf("\n Nb of events analysed: %d \n",iEvent); | |
345 | cout << " number of tracks: " << nTracks <<endl; | |
346 | } | |
347 | ||
348 | // set the magnetic field for track extrapolations | |
349 | AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap()); | |
350 | // loop over all reconstructed tracks (also first track of combination) | |
351 | for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) { | |
352 | ||
353 | AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack))); | |
354 | ||
355 | // extrapolate to vertex if required and available | |
356 | if (ExtrapToVertex > 0 && Vertex->GetNContributors()) { | |
357 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam); | |
358 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx); | |
359 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack | |
360 | } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){ | |
361 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam); | |
362 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.); | |
363 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack | |
364 | } | |
365 | ||
366 | // Trigger | |
367 | if (PRINTLEVEL > 5) cout << "MatchTrigger " << muonTrack->GetMatchTrigger() << " and Chi2 of matching tracks " << track1TriggerChi2 << endl ; | |
368 | track1Trigger = muonTrack->GetMatchTrigger(); | |
369 | if (track1Trigger) | |
370 | track1TriggerChi2 = muonTrack->GetChi2MatchTrigger(); | |
371 | else | |
372 | track1TriggerChi2 = 0. ; | |
373 | ||
374 | fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum())); | |
375 | ||
376 | muonTrack->LorentzP(fV1); | |
377 | ||
378 | ntrackhits = muonTrack->GetNHit(); | |
379 | fitfmin = muonTrack->GetChi2(); | |
380 | ||
381 | // transverse momentum | |
382 | Float_t pt1 = fV1.Pt(); | |
383 | ||
384 | // total momentum | |
385 | Float_t p1 = fV1.P(); | |
386 | ||
387 | // Rapidity | |
388 | Float_t rapMuon1 = fV1.Rapidity(); | |
389 | ||
390 | // chi2 per d.o.f. | |
391 | ||
392 | Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5); | |
393 | if (PRINTLEVEL > 5 ) printf(" px %f py %f pz %f pt %f NHits %d Norm.chi2 %f charge %d\n",fV1.Px(), fV1.Py(), fV1.Pz(), pt1, ntrackhits, ch1, fCharge1); | |
394 | ||
395 | ||
396 | if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { // condition for good track (Chi2Cut and PtCut) | |
397 | if (PRINTLEVEL > 8) cout << "inside pt and chi2 cuts " << endl ; | |
398 | ||
399 | // fill histos hPtMuon and hChi2PerDof | |
400 | hPtMuon->Fill(pt1); | |
401 | hPMuon->Fill(p1); | |
402 | hChi2PerDof->Fill(ch1); | |
403 | hRapMuon->Fill(rapMuon1); | |
404 | ||
405 | if (fCharge1 > 0) { | |
406 | hPtMuonPlus->Fill(pt1); | |
407 | hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi()); | |
408 | } else { | |
409 | hPtMuonMinus->Fill(pt1); | |
410 | hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi()); | |
411 | } | |
412 | ||
413 | // loop over second track of combination | |
414 | for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) { | |
415 | ||
416 | AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2))); | |
417 | ||
418 | // extrapolate to vertex if required and available | |
419 | if (ExtrapToVertex > 0 && Vertex->GetNContributors()) { | |
420 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam); | |
421 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx); | |
422 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack | |
423 | } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){ | |
424 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam); | |
425 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.); | |
426 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack | |
427 | } | |
428 | ||
429 | track2Trigger = muonTrack2->GetMatchTrigger(); | |
430 | if (track2Trigger) | |
431 | track2TriggerChi2 = muonTrack2->GetChi2MatchTrigger(); | |
432 | else | |
433 | track2TriggerChi2 = 0. ; | |
434 | ||
435 | fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum())); | |
436 | ||
437 | muonTrack2->LorentzP(fV2); | |
438 | ||
439 | ntrackhits = muonTrack2->GetNHit(); | |
440 | fitfmin = muonTrack2->GetChi2(); | |
441 | ||
442 | // transverse momentum | |
443 | Float_t pt2 = fV2.Pt(); | |
444 | ||
445 | // chi2 per d.o.f. | |
446 | Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5); | |
447 | ||
448 | ||
449 | // condition for good track (Chi2Cut and PtCut) | |
450 | if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) { | |
451 | ||
452 | // condition for opposite charges | |
453 | if ((fCharge1 * fCharge2) == -1) { | |
454 | ||
455 | if (PRINTLEVEL > 8) cout << "---------> Now filling the Ntuple " << endl ; | |
456 | ||
457 | // invariant mass | |
458 | fVtot = fV1 + fV2; | |
459 | Float_t invMass = fVtot.M(); | |
460 | ||
461 | if (fCharge1 < 0){ //mu_minus is index 1 in the ntuple | |
462 | Float_t ESDFill[16] = {iEvent,triggerWord,fVtot.Pt(),fVtot.Rapidity(),fVtot.Theta()/TMath::Pi()*180,invMass,fV1.Pt(),fV1.Rapidity(),fV1.Theta()/TMath::Pi()*180,fCharge1,track1TriggerChi2,fV2.Pt(),fV2.Rapidity(),fV2.Theta()/TMath::Pi()*180,fCharge2,track2TriggerChi2}; | |
463 | ESDtuple->Fill(ESDFill); | |
464 | } | |
465 | else{ | |
466 | Float_t ESDFill[16] = {iEvent,triggerWord,fVtot.Pt(),fVtot.Rapidity(),fVtot.Theta()/TMath::Pi()*180,invMass,fV2.Pt(),fV2.Rapidity(),fV2.Theta()/TMath::Pi()*180,fCharge2,track2TriggerChi2,fV1.Pt(),fV1.Rapidity(),fV1.Theta()/TMath::Pi()*180,fCharge1,track1TriggerChi2}; | |
467 | ESDtuple->Fill(ESDFill); | |
468 | } | |
469 | ||
470 | // fill histos hInvMassAll and hInvMassRes | |
471 | hInvMassAll->Fill(invMass); | |
472 | hInvMassRes->Fill(invMass); | |
473 | hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt()); | |
474 | ||
475 | //trigger info | |
476 | if (ResType == 553) | |
477 | ptTrig = 0x08;// mask for Hpt unlike sign pair | |
478 | else if (ResType == 443) | |
479 | ptTrig = 0x04;// mask for Lpt unlike sign pair | |
480 | ||
481 | ||
482 | if (esd->GetTriggerMask() & ptTrig) NbTrigger++; | |
483 | ||
484 | if (invMass > invMassMinInPeak && invMass < invMassMaxInPeak) { | |
485 | EventInMass++; | |
486 | hRapResonance->Fill(fVtot.Rapidity()); | |
487 | hPtResonance->Fill(fVtot.Pt()); | |
488 | ||
489 | // match with trigger | |
490 | if (muonTrack2->GetMatchTrigger()>=0 && (esd->GetTriggerMask() & ptTrig)) EventInMassMatch++; | |
491 | ||
492 | } | |
493 | ||
494 | } // if (fCharge1 * fCharge2) == -1) | |
495 | } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) | |
496 | delete muonTrack2; | |
497 | } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++) | |
498 | } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) ) | |
499 | delete muonTrack; | |
500 | } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++) | |
501 | ||
502 | hNumberOfTrack->Fill(nTracks); | |
503 | // esdFile->Delete(); | |
504 | ||
505 | } // End of event loop | |
506 | ||
507 | ||
508 | // Loop over events for bg event | |
509 | ||
510 | Double_t thetaPlus, phiPlus; | |
511 | Double_t thetaMinus, phiMinus; | |
512 | Float_t PtMinus, PtPlus; | |
513 | ||
514 | for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { // Loop over events for bg event | |
515 | // according to Christian a 3d phi-theta-pt random pick would take better care | |
516 | // of all correlations | |
517 | ||
518 | hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus); | |
519 | hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus); | |
520 | PtPlus = hPtMuonPlus->GetRandom(); | |
521 | PtMinus = hPtMuonMinus->GetRandom(); | |
522 | ||
523 | fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus); | |
524 | fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus); | |
525 | fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus); | |
526 | ||
527 | fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1); | |
528 | fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1); | |
529 | ||
530 | fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus); | |
531 | fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus); | |
532 | fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus); | |
533 | ||
534 | fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2); | |
535 | fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2); | |
536 | ||
537 | // invariant mass | |
538 | fVtot = fV1 + fV2; | |
539 | ||
540 | // fill histos hInvMassAll and hInvMassRes | |
541 | hInvMassBg->Fill(fVtot.M()); | |
542 | hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() ); | |
543 | ||
544 | // Ntuple for background... more convenient | |
545 | ESDtupleBck->Fill(iEvent,fVtot.Pt(),fVtot.Rapidity(),fVtot.Theta()/TMath::Pi()*180,fVtot.M(),fV2.Pt(),fV2.Rapidity(),fV2.Theta()/TMath::Pi()*180,fV1.Pt(),fV1.Rapidity(),fV1.Theta()/TMath::Pi()*180); | |
546 | ||
547 | } // End loop over events for background | |
548 | ||
549 | ||
550 | // File for histograms and histogram booking | |
551 | TString outfilename = "MUONefficiency.root"; | |
552 | TFile *ntupleFile = new TFile(outfilename.Data(), "RECREATE"); | |
553 | ||
554 | Ktuple->Write(); | |
555 | ESDtuple->Write(); | |
556 | ESDtupleBck->Write(); | |
557 | ||
558 | ntupleFile->Close(); | |
559 | ||
560 | TFile *histoFile = new TFile("MUONhistos.root", "RECREATE"); | |
561 | hPrimaryVertex->Write(); | |
562 | hPtMuon->Write(); | |
563 | hPtMuonPlus->Write(); | |
564 | hPtMuonMinus->Write(); | |
565 | hPMuon->Write(); | |
566 | hChi2PerDof->Write(); | |
567 | hInvMassAll->Write(); | |
568 | hInvMassBg->Write(); | |
569 | hInvMassAll_vs_Pt ->Write(); | |
570 | hInvMassBgk_vs_Pt->Write(); | |
571 | hInvMassRes->Write(); | |
572 | hNumberOfTrack->Write(); | |
573 | hRapMuon ->Write(); | |
574 | hRapResonance ->Write(); | |
575 | hPtResonance ->Write(); | |
576 | hThetaPhiPlus ->Write(); | |
577 | hThetaPhiMinus ->Write(); | |
578 | histoFile->Close(); | |
579 | ||
580 | cout << "" << endl ; | |
581 | cout << "*************************************************" << endl; | |
582 | ||
583 | cout << "MUONefficiency : " << nprocessedevents << " events processed" << endl; | |
584 | if (ResType==443) | |
585 | cout << "Number of generated J/Psi (443) : " << numberOfGeneratedResonances << endl ; | |
586 | if (ResType==553) | |
587 | cout << "Number of generated Upsilon (553) :" << numberOfGeneratedResonances << endl ; | |
588 | cout << "Chi2Cut for muon tracks = " << Chi2Cut << endl; | |
589 | cout << "PtCutMin for muon tracks = " << PtCutMin << endl; | |
590 | cout << "PtCutMax for muon tracks = " << PtCutMax << endl; | |
591 | cout << "Entries (unlike sign dimuons) in the mass range ["<<invMassMinInPeak<<";"<<invMassMaxInPeak<<"] : " << EventInMass <<endl; | |
592 | if (ptTrig==0x800) cout << "Unlike Pair - All Pt" ; | |
593 | if (ptTrig==0x400) cout << "Unlike Pair - High Pt" ; | |
594 | if (ptTrig==0x200) cout << "Unlike Pair - Low Pt" ; | |
595 | cout << " triggers : " << NbTrigger << endl; | |
596 | ||
597 | cout << "Entries in the mass range with matching between reconstructed tracks and trigger tracks " << EventInMassMatch << endl; | |
598 | ||
599 | ||
600 | return kTRUE; | |
601 | } |