1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // Macro (upgraded version of MUONmassPlot_ESD.C, better handling of Jpsi) to make :
19 // 1) Ntuple (Ktuple) containing Upsilon kinematics variables (from kinematics.root files)
20 // 2) Ntuple (ESDtuple) containing Upsilon kinematics variables from reconstruction and
21 // combinations of 2 muons with opposite charges (ESDtupleBck will be used later)
22 // 3) Some QA histograms
23 // Ntuple are stored in the file MUONefficiency.root and ESD tree and QA histograms in AliESDs.root
25 // Christophe Suire, IPN Orsay
30 // FirstEvent (default 0)
31 // LastEvent (default 1.e6)
32 // ResType (default 553)
33 // 553 for Upsilon, 443 for J/Psi
34 // Chi2Cut (default 100)
35 // to keep only tracks with chi2 per d.o.f. < Chi2Cut
39 #if !defined(__CINT__) || defined(__MAKECINT__)
44 #include "TClonesArray.h"
45 #include "TLorentzVector.h"
49 #include "TParticle.h"
52 #include <Riostream.h>
56 #include "AliRunLoader.h"
57 #include "AliHeader.h"
58 #include "AliLoader.h"
60 #include "AliMagFMaps.h"
62 #include "AliTracker.h"
65 #include "AliMUONTrackParam.h"
66 #include "AliMUONTrackExtrap.h"
67 #include "AliESDMuonTrack.h"
71 Bool_t MUONefficiency( Int_t ResType = 553, Int_t FirstEvent = 0, Int_t LastEvent = 1000000,
72 char* esdFileName = "AliESDs.root", char* filename = "galice.root")
73 { // MUONefficiency starts
75 Double_t MUON_MASS = 0.105658369;
76 Double_t UPSILON_MASS = 9.4603 ;
77 Double_t JPSI_MASS = 3.097;
79 // Upper and lower bound for counting entries in the mass peak
80 // +/- 300 MeV/c^2 in this case.
81 Float_t countingRange = 0.300 ;
83 Float_t massResonance = 5.;
84 Float_t invMassMinInPeak = 0. ;
85 Float_t invMassMaxInPeak = 0. ;
87 Float_t nBinsPerGev = 40 ;
88 Float_t invMassMin = 0; Float_t invMassMax = 20;
89 Float_t ptMinResonance = 0 ; Float_t ptMaxResonance = 20 ; Int_t ptBinsResonance = 100;
92 massResonance = JPSI_MASS ;
93 invMassMinInPeak = JPSI_MASS - countingRange ; invMassMaxInPeak = JPSI_MASS + countingRange ;
94 //limits for histograms
95 invMassMin = 0 ; invMassMax = 6.;
96 ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100;
99 massResonance = UPSILON_MASS;
100 invMassMinInPeak = UPSILON_MASS - countingRange ; invMassMaxInPeak = UPSILON_MASS + countingRange;
101 //limits for histograms
102 invMassMin = 0 ; invMassMax = 12.;
103 ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100;
106 // Single Tracks muon cuts
107 Float_t Chi2Cut = 100.;
108 Float_t PtCutMin = 0. ;
109 Float_t PtCutMax = 10000. ;
112 // Limits for histograms
113 Float_t ptMinMuon = 0. ; Float_t ptMaxMuon = 20.; Int_t ptBinsMuon = 100 ;
114 Float_t pMinMuon = 0. ; Float_t pMaxMuon = 200.; Int_t pBinsMuon = 100 ;
117 //Reset ROOT and connect tree file
121 Int_t PRINTLEVEL = 0 ;
123 //for kinematic, i.e. reference tracks
124 TNtuple *Ktuple = new TNtuple("Ktuple","Kinematics NTuple","ev:npart:id:idmo:idgdmo:p:pt:y:theta:pseudorap:vx:vy:vz");
127 TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
128 TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
129 TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
130 TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", pBinsMuon, pMinMuon, pMaxMuon);
134 TH2F *hInvMassAll_vs_Pt;
135 TH2F *hInvMassBgk_vs_Pt;
139 hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax - invMassMin)), invMassMin, invMassMax);
140 hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax);
141 hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
142 hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
144 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Resonance",(Int_t) (nBinsPerGev*3*countingRange*2),massResonance-3*countingRange,massResonance+3*countingRange);
146 TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
147 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
148 TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
149 TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
150 TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
151 TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
152 TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
153 TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
155 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");
156 TNtuple *ESDtupleBck = new TNtuple("ESDtupleBck","Reconstructed Mu+Mu- pairs for Background","ev:pt:y:theta:minv:pt1:y1:theta1:pt2:y2:theta2");
160 Int_t EventInMass = 0;
161 Int_t EventInMassMatch = 0;
169 Double_t thetaX, thetaY, pYZ;
170 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
171 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
172 Int_t fCharge1, fCharge2;
174 Int_t ntrackhits, nevents;
175 Int_t nprocessedevents = 0 ;
178 TLorentzVector fV1, fV2, fVtot;
181 // waiting for mag field in CDB
182 printf("Loading field map...\n");
183 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
184 AliTracker::SetFieldMap(field, kFALSE);
186 // open run loader and load gAlice, kinematics and header
187 AliRunLoader* runLoader = AliRunLoader::Open(filename);
189 Error("MUONefficiency", "getting run loader from file %s failed", filename);
193 runLoader->LoadgAlice();
194 gAlice = runLoader->GetAliRun();
196 Error("MUONefficiency", "no galice object found");
201 TFile* esdFile = TFile::Open(esdFileName);
202 if (!esdFile || !esdFile->IsOpen()) {
203 Error("MUONefficiency", "opening ESD file %s failed", esdFileName);
207 AliESD* esd = new AliESD();
208 TTree* tree = (TTree*) esdFile->Get("esdTree");
210 Error("CheckESD", "no ESD tree found");
213 tree->SetBranchAddress("ESD", &esd);
215 runLoader->LoadHeader();
216 nevents = runLoader->GetNumberOfEvents();
217 AliMUONTrackParam trackParam;
219 // to access the particle Stack
220 runLoader->LoadKinematics("READ");
222 Int_t numberOfGeneratedResonances = 0 ;
226 Int_t track1Trigger = 0 ;
227 Float_t track1TriggerChi2 = 0 ;
228 Int_t track2Trigger = 0 ;
229 Float_t track2TriggerChi2 = 0 ;
233 for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { // Start event loop
235 if (iEvent%1000 == 0 )
236 printf("\n Nb of events analysed: %d \n",iEvent);
239 runLoader->GetEvent(iEvent);
242 // get the stack and fill the kine tree
243 AliStack *theStack = runLoader->Stack();
244 if (PRINTLEVEL > 0) theStack->DumpPStack ();
246 Int_t nparticles = (Int_t)runLoader->TreeK()->GetEntries();
247 Int_t nprimarypart = theStack->GetNprimary();
248 Int_t ntracks = theStack->GetNtrack();
250 if (PRINTLEVEL || (iEvent%100==0)) printf("\n >>> Event %d \n",iEvent);
251 if (PRINTLEVEL) cout << nprimarypart << " Particles generated (total is " << ntracks << ")"<< endl ;
253 for(Int_t iparticle=0; iparticle<nparticles; iparticle++) { // Start loop over particles
254 particle = theStack->Particle(iparticle);
256 Int_t muId = particle->GetPdgCode();
257 Int_t muM = particle->GetFirstMother();
259 Float_t muP = particle->P();
260 Float_t muPt = TMath::Sqrt(particle->Px()*particle->Px()+particle->Py()*particle->Py());
261 Float_t muY = 0.5*TMath::Log((particle->Energy()+particle->Pz()+1.e-13)/(particle->Energy()-particle->Pz()+1.e-13));
263 TParticle *theMum = theStack->Particle(muM);
264 muM = theMum->GetPdgCode();
265 muGM = theMum->GetFirstMother() ;
267 TParticle *grandMa = theStack->Particle(muGM);
268 muGM = grandMa->GetPdgCode();
274 if (muId==ResType) numberOfGeneratedResonances++;
277 Float_t muT = particle->Theta()*180/TMath::Pi();
278 Float_t muE = particle->Eta();
280 Float_t muVx = particle->Vx();
281 Float_t muVy = particle->Vy();
282 Float_t muVz = particle->Vz();
284 // If a write error occurs, the number of bytes returned is -1.
285 // If no data are written, because e.g. the branch is disabled,
286 // the number of bytes returned is 0.
287 Int_t errCode = Ktuple->Fill(iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz);
288 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);
290 } // End loop over particles
294 // get the event summary data
295 tree->GetEvent(iEvent);
297 Error("CheckESD", "no ESD object found for event %d", iEvent);
301 // get the SPD reconstructed vertex (vertexer) and fill the histogram
302 AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
303 if (Vertex->GetNContributors()) {
304 fZVertex = Vertex->GetZv();
305 fYVertex = Vertex->GetYv();
306 fXVertex = Vertex->GetXv();
308 hPrimaryVertex->Fill(fZVertex);
310 Int_t triggerWord = esd->GetTriggerMask();
311 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
314 printf("\n Nb of events analysed: %d \n",iEvent);
315 cout << " number of tracks: " << nTracks <<endl;
318 // set the magnetic field for track extrapolations
319 AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
320 // loop over all reconstructed tracks (also first track of combination)
321 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
323 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
325 if (!Vertex->GetNContributors()) {
326 //re-extrapolate to vertex, if not kown before.
327 trackParam.GetParamFrom(*muonTrack);
328 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
329 trackParam.SetParamFor(*muonTrack);
333 if (PRINTLEVEL > 5) cout << "MatchTrigger " << muonTrack->GetMatchTrigger() << " and Chi2 of matching tracks " << track1TriggerChi2 << endl ;
334 track1Trigger = muonTrack->GetMatchTrigger();
336 track1TriggerChi2 = muonTrack->GetChi2MatchTrigger();
338 track1TriggerChi2 = 0. ;
340 thetaX = muonTrack->GetThetaX();
341 thetaY = muonTrack->GetThetaY();
343 pYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
344 fPzRec1 = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
345 fPxRec1 = fPzRec1 * TMath::Tan(thetaX);
346 fPyRec1 = fPzRec1 * TMath::Tan(thetaY);
347 fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
349 fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
350 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
352 ntrackhits = muonTrack->GetNHit();
353 fitfmin = muonTrack->GetChi2();
355 // transverse momentum
356 Float_t pt1 = fV1.Pt();
359 Float_t p1 = fV1.P();
362 Float_t rapMuon1 = fV1.Rapidity();
366 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
367 if (PRINTLEVEL > 5 ) printf(" px %f py %f pz %f pt %f NHits %d Norm.chi2 %f charge %d\n",fPxRec1, fPyRec1, fPzRec1, pt1, ntrackhits, ch1, fCharge1);
370 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { // condition for good track (Chi2Cut and PtCut)
371 if (PRINTLEVEL > 8) cout << "inside pt and chi2 cuts " << endl ;
373 // fill histos hPtMuon and hChi2PerDof
376 hChi2PerDof->Fill(ch1);
377 hRapMuon->Fill(rapMuon1);
380 hPtMuonPlus->Fill(pt1);
381 hThetaPhiPlus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
383 hPtMuonMinus->Fill(pt1);
384 hThetaPhiMinus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
387 // loop over second track of combination
388 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
390 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack2);
392 if (!Vertex->GetNContributors()) {
393 trackParam.GetParamFrom(*muonTrack);
394 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
395 trackParam.SetParamFor(*muonTrack);
398 track2Trigger = muonTrack->GetMatchTrigger();
400 track2TriggerChi2 = muonTrack->GetChi2MatchTrigger();
402 track2TriggerChi2 = 0. ;
404 thetaX = muonTrack->GetThetaX();
405 thetaY = muonTrack->GetThetaY();
407 pYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
408 fPzRec2 = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
409 fPxRec2 = fPzRec2 * TMath::Tan(thetaX);
410 fPyRec2 = fPzRec2 * TMath::Tan(thetaY);
411 fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
413 fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
414 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
416 ntrackhits = muonTrack->GetNHit();
417 fitfmin = muonTrack->GetChi2();
419 // transverse momentum
420 Float_t pt2 = fV2.Pt();
423 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
426 // condition for good track (Chi2Cut and PtCut)
427 if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
429 // condition for opposite charges
430 if ((fCharge1 * fCharge2) == -1) {
432 if (PRINTLEVEL > 8) cout << "---------> Now filling the Ntuple " << endl ;
436 Float_t invMass = fVtot.M();
438 if (fCharge1 < 0){ //mu_minus is index 1 in the ntuple
439 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};
440 ESDtuple->Fill(ESDFill);
443 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};
444 ESDtuple->Fill(ESDFill);
447 // fill histos hInvMassAll and hInvMassRes
448 hInvMassAll->Fill(invMass);
449 hInvMassRes->Fill(invMass);
450 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
454 ptTrig = 0x20;// mask for Hpt unlike sign pair
455 else if (ResType == 443)
456 ptTrig = 0x10;// mask for Lpt unlike sign pair
459 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
461 if (invMass > invMassMinInPeak && invMass < invMassMaxInPeak) {
463 hRapResonance->Fill(fVtot.Rapidity());
464 hPtResonance->Fill(fVtot.Pt());
466 // match with trigger
467 if (muonTrack->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig)) EventInMassMatch++;
471 } // if (fCharge1 * fCharge2) == -1)
472 } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
473 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
474 } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
475 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
477 hNumberOfTrack->Fill(nTracks);
478 // esdFile->Delete();
480 } // End of event loop
483 // Loop over events for bg event
485 Double_t thetaPlus, phiPlus;
486 Double_t thetaMinus, phiMinus;
487 Float_t PtMinus, PtPlus;
489 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { // Loop over events for bg event
490 // according to Christian a 3d phi-theta-pt random pick would take better care
491 // of all correlations
493 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
494 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
495 PtPlus = hPtMuonPlus->GetRandom();
496 PtMinus = hPtMuonMinus->GetRandom();
498 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
499 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
500 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
502 fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
503 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
505 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
506 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
507 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
509 fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
510 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
515 // fill histos hInvMassAll and hInvMassRes
516 hInvMassBg->Fill(fVtot.M());
517 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
519 // Ntuple for background... more convenient
520 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);
522 } // End loop over events for background
525 // File for histograms and histogram booking
526 TString outfilename = "MUONefficiency.root";
527 TFile *ntupleFile = new TFile(outfilename.Data(), "RECREATE");
531 ESDtupleBck->Write();
535 TFile *histoFile = new TFile("MUONhistos.root", "RECREATE");
536 hPrimaryVertex->Write();
538 hPtMuonPlus->Write();
539 hPtMuonMinus->Write();
541 hChi2PerDof->Write();
542 hInvMassAll->Write();
544 hInvMassAll_vs_Pt ->Write();
545 hInvMassBgk_vs_Pt->Write();
546 hInvMassRes->Write();
547 hNumberOfTrack->Write();
549 hRapResonance ->Write();
550 hPtResonance ->Write();
551 hThetaPhiPlus ->Write();
552 hThetaPhiMinus ->Write();
556 cout << "*************************************************" << endl;
558 cout << "MUONefficiency : " << nprocessedevents << " events processed" << endl;
560 cout << "Number of generated J/Psi (443) : " << numberOfGeneratedResonances << endl ;
562 cout << "Number of generated Upsilon (553) :" << numberOfGeneratedResonances << endl ;
563 cout << "Chi2Cut for muon tracks = " << Chi2Cut << endl;
564 cout << "PtCutMin for muon tracks = " << PtCutMin << endl;
565 cout << "PtCutMax for muon tracks = " << PtCutMax << endl;
566 cout << "Entries (unlike sign dimuons) in the mass range ["<<invMassMinInPeak<<";"<<invMassMaxInPeak<<"] : " << EventInMass <<endl;
567 if (ptTrig==0x800) cout << "Unlike Pair - All Pt" ;
568 if (ptTrig==0x400) cout << "Unlike Pair - High Pt" ;
569 if (ptTrig==0x200) cout << "Unlike Pair - Low Pt" ;
570 cout << " triggers : " << NbTrigger << endl;
572 cout << "Entries in the mass range with matching between reconstructed tracks and trigger tracks " << EventInMassMatch << endl;