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>
53 #include <TGeoManager.h>
58 #include "AliRunLoader.h"
59 #include "AliHeader.h"
60 #include "AliLoader.h"
62 #include "AliMagFMaps.h"
63 #include "AliESDEvent.h"
64 #include "AliESDVertex.h"
65 #include "AliTracker.h"
66 #include "AliCDBManager.h"
69 #include "AliMUONTrackParam.h"
70 #include "AliMUONTrackExtrap.h"
71 #include "AliESDMuonTrack.h"
75 // ExtrapToVertex (default -1)
76 // <0: no extrapolation;
77 // =0: extrapolation to (0,0,0);
78 // >0: extrapolation to ESDVertex if available, else to (0,0,0)
79 // ResType (default 553)
80 // 553 for Upsilon, anything else for J/Psi
82 Bool_t MUONefficiency( char* filename = "galice.root", char* geoFilename = "geometry.root", char* esdFileName = "AliESDs.root",
83 Int_t ExtrapToVertex = -1, Int_t ResType = 553, Int_t FirstEvent = 0, Int_t LastEvent = 1000000 )
84 { // MUONefficiency starts
86 // Set default CDB storage
87 AliCDBManager* man = AliCDBManager::Instance();
88 man->SetDefaultStorage("local://$ALICE_ROOT");
90 Double_t MUON_MASS = 0.105658369;
91 Double_t UPSILON_MASS = 9.4603 ;
92 Double_t JPSI_MASS = 3.097;
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 ;
98 Float_t massResonance = 5.;
99 Float_t invMassMinInPeak = 0. ;
100 Float_t invMassMaxInPeak = 0. ;
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;
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;
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;
121 // Single Tracks muon cuts
122 Float_t Chi2Cut = 100.;
123 Float_t PtCutMin = 0. ;
124 Float_t PtCutMax = 10000. ;
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 ;
132 //Reset ROOT and connect tree file
136 Int_t PRINTLEVEL = 0 ;
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");
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);
149 TH2F *hInvMassAll_vs_Pt;
150 TH2F *hInvMassBgk_vs_Pt;
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);
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);
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.);
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");
175 Int_t EventInMass = 0;
176 Int_t EventInMassMatch = 0;
186 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
187 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
188 Int_t fCharge1, fCharge2;
190 Int_t ntrackhits, nevents;
191 Int_t nprocessedevents = 0 ;
194 TLorentzVector fV1, fV2, fVtot;
196 // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
198 TGeoManager::Import(geoFilename);
200 Error("MUONmass_ESD", "getting geometry from file %s failed", filename);
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);
211 // open run loader and load gAlice, kinematics and header
212 AliRunLoader* runLoader = AliRunLoader::Open(filename);
214 Error("MUONefficiency", "getting run loader from file %s failed", filename);
218 runLoader->LoadgAlice();
219 gAlice = runLoader->GetAliRun();
221 Error("MUONefficiency", "no galice object found");
226 TFile* esdFile = TFile::Open(esdFileName);
227 if (!esdFile || !esdFile->IsOpen()) {
228 Error("MUONefficiency", "opening ESD file %s failed", esdFileName);
232 AliESDEvent* esd = new AliESDEvent();
233 TTree* tree = (TTree*) esdFile->Get("esdTree");
235 Error("CheckESD", "no ESD tree found");
238 esd->ReadFromTree(tree);
240 runLoader->LoadHeader();
241 nevents = runLoader->GetNumberOfEvents();
242 AliMUONTrackParam trackParam;
244 // to access the particle Stack
245 runLoader->LoadKinematics("READ");
247 Int_t numberOfGeneratedResonances = 0 ;
251 Int_t track1Trigger = 0 ;
252 Float_t track1TriggerChi2 = 0 ;
253 Int_t track2Trigger = 0 ;
254 Float_t track2TriggerChi2 = 0 ;
258 for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { // Start event loop
260 if (iEvent%1000 == 0 )
261 printf("\n Nb of events analysed: %d \n",iEvent);
264 runLoader->GetEvent(iEvent);
267 // get the stack and fill the kine tree
268 AliStack *theStack = runLoader->Stack();
269 if (PRINTLEVEL > 0) theStack->DumpPStack ();
271 Int_t nparticles = (Int_t)runLoader->TreeK()->GetEntries();
272 Int_t nprimarypart = theStack->GetNprimary();
273 Int_t ntracks = theStack->GetNtrack();
275 if (PRINTLEVEL || (iEvent%100==0)) printf("\n >>> Event %d \n",iEvent);
276 if (PRINTLEVEL) cout << nprimarypart << " Particles generated (total is " << ntracks << ")"<< endl ;
278 for(Int_t iparticle=0; iparticle<nparticles; iparticle++) { // Start loop over particles
279 particle = theStack->Particle(iparticle);
281 Int_t muId = particle->GetPdgCode();
282 Int_t muM = particle->GetFirstMother();
284 Float_t muP = particle->P();
285 Float_t muPt = TMath::Sqrt(particle->Px()*particle->Px()+particle->Py()*particle->Py());
286 Float_t muY = 0.5*TMath::Log((particle->Energy()+particle->Pz()+1.e-13)/(particle->Energy()-particle->Pz()+1.e-13));
288 TParticle *theMum = theStack->Particle(muM);
289 muM = theMum->GetPdgCode();
290 muGM = theMum->GetFirstMother() ;
292 TParticle *grandMa = theStack->Particle(muGM);
293 muGM = grandMa->GetPdgCode();
299 if (muId==ResType) numberOfGeneratedResonances++;
302 Float_t muT = particle->Theta()*180/TMath::Pi();
303 Float_t muE = particle->Eta();
305 Float_t muVx = particle->Vx();
306 Float_t muVy = particle->Vy();
307 Float_t muVz = particle->Vz();
309 // If a write error occurs, the number of bytes returned is -1.
310 // If no data are written, because e.g. the branch is disabled,
311 // the number of bytes returned is 0.
312 Int_t errCode = Ktuple->Fill(iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz);
313 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);
315 } // End loop over particles
319 // get the event summary data
320 tree->GetEvent(iEvent);
322 Error("CheckESD", "no ESD object found for event %d", iEvent);
326 // get the SPD reconstructed vertex (vertexer) and fill the histogram
327 AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
328 if (Vertex->GetNContributors()) {
329 fZVertex = Vertex->GetZv();
330 fYVertex = Vertex->GetYv();
331 fXVertex = Vertex->GetXv();
332 errXVtx = Vertex->GetXRes();
333 errYVtx = Vertex->GetYRes();
335 hPrimaryVertex->Fill(fZVertex);
337 Int_t triggerWord = esd->GetTriggerMask();
338 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
341 printf("\n Nb of events analysed: %d \n",iEvent);
342 cout << " number of tracks: " << nTracks <<endl;
345 // set the magnetic field for track extrapolations
346 AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
347 // loop over all reconstructed tracks (also first track of combination)
348 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
350 AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
352 // extrapolate to vertex if required and available
353 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
354 trackParam.GetParamFromUncorrected(*muonTrack);
355 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
356 trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
357 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
358 trackParam.GetParamFromUncorrected(*muonTrack);
359 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
360 trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
364 if (PRINTLEVEL > 5) cout << "MatchTrigger " << muonTrack->GetMatchTrigger() << " and Chi2 of matching tracks " << track1TriggerChi2 << endl ;
365 track1Trigger = muonTrack->GetMatchTrigger();
367 track1TriggerChi2 = muonTrack->GetChi2MatchTrigger();
369 track1TriggerChi2 = 0. ;
371 fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
373 muonTrack->LorentzP(fV1);
375 ntrackhits = muonTrack->GetNHit();
376 fitfmin = muonTrack->GetChi2();
378 // transverse momentum
379 Float_t pt1 = fV1.Pt();
382 Float_t p1 = fV1.P();
385 Float_t rapMuon1 = fV1.Rapidity();
389 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
390 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);
393 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { // condition for good track (Chi2Cut and PtCut)
394 if (PRINTLEVEL > 8) cout << "inside pt and chi2 cuts " << endl ;
396 // fill histos hPtMuon and hChi2PerDof
399 hChi2PerDof->Fill(ch1);
400 hRapMuon->Fill(rapMuon1);
403 hPtMuonPlus->Fill(pt1);
404 hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
406 hPtMuonMinus->Fill(pt1);
407 hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
410 // loop over second track of combination
411 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
413 AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
415 // extrapolate to vertex if required and available
416 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
417 trackParam.GetParamFromUncorrected(*muonTrack2);
418 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
419 trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
420 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
421 trackParam.GetParamFromUncorrected(*muonTrack2);
422 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
423 trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
426 track2Trigger = muonTrack2->GetMatchTrigger();
428 track2TriggerChi2 = muonTrack2->GetChi2MatchTrigger();
430 track2TriggerChi2 = 0. ;
432 fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
434 muonTrack2->LorentzP(fV2);
436 ntrackhits = muonTrack2->GetNHit();
437 fitfmin = muonTrack2->GetChi2();
439 // transverse momentum
440 Float_t pt2 = fV2.Pt();
443 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
446 // condition for good track (Chi2Cut and PtCut)
447 if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
449 // condition for opposite charges
450 if ((fCharge1 * fCharge2) == -1) {
452 if (PRINTLEVEL > 8) cout << "---------> Now filling the Ntuple " << endl ;
456 Float_t invMass = fVtot.M();
458 if (fCharge1 < 0){ //mu_minus is index 1 in the ntuple
459 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};
460 ESDtuple->Fill(ESDFill);
463 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};
464 ESDtuple->Fill(ESDFill);
467 // fill histos hInvMassAll and hInvMassRes
468 hInvMassAll->Fill(invMass);
469 hInvMassRes->Fill(invMass);
470 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
474 ptTrig = 0x08;// mask for Hpt unlike sign pair
475 else if (ResType == 443)
476 ptTrig = 0x04;// mask for Lpt unlike sign pair
479 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
481 if (invMass > invMassMinInPeak && invMass < invMassMaxInPeak) {
483 hRapResonance->Fill(fVtot.Rapidity());
484 hPtResonance->Fill(fVtot.Pt());
486 // match with trigger
487 if (muonTrack2->GetMatchTrigger()>=0 && (esd->GetTriggerMask() & ptTrig)) EventInMassMatch++;
491 } // if (fCharge1 * fCharge2) == -1)
492 } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
494 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
495 } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
497 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
499 hNumberOfTrack->Fill(nTracks);
500 // esdFile->Delete();
502 } // End of event loop
505 // Loop over events for bg event
507 Double_t thetaPlus, phiPlus;
508 Double_t thetaMinus, phiMinus;
509 Float_t PtMinus, PtPlus;
511 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { // Loop over events for bg event
512 // according to Christian a 3d phi-theta-pt random pick would take better care
513 // of all correlations
515 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
516 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
517 PtPlus = hPtMuonPlus->GetRandom();
518 PtMinus = hPtMuonMinus->GetRandom();
520 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
521 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
522 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
524 fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
525 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
527 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
528 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
529 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
531 fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
532 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
537 // fill histos hInvMassAll and hInvMassRes
538 hInvMassBg->Fill(fVtot.M());
539 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
541 // Ntuple for background... more convenient
542 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);
544 } // End loop over events for background
547 // File for histograms and histogram booking
548 TString outfilename = "MUONefficiency.root";
549 TFile *ntupleFile = new TFile(outfilename.Data(), "RECREATE");
553 ESDtupleBck->Write();
557 TFile *histoFile = new TFile("MUONhistos.root", "RECREATE");
558 hPrimaryVertex->Write();
560 hPtMuonPlus->Write();
561 hPtMuonMinus->Write();
563 hChi2PerDof->Write();
564 hInvMassAll->Write();
566 hInvMassAll_vs_Pt ->Write();
567 hInvMassBgk_vs_Pt->Write();
568 hInvMassRes->Write();
569 hNumberOfTrack->Write();
571 hRapResonance ->Write();
572 hPtResonance ->Write();
573 hThetaPhiPlus ->Write();
574 hThetaPhiMinus ->Write();
578 cout << "*************************************************" << endl;
580 cout << "MUONefficiency : " << nprocessedevents << " events processed" << endl;
582 cout << "Number of generated J/Psi (443) : " << numberOfGeneratedResonances << endl ;
584 cout << "Number of generated Upsilon (553) :" << numberOfGeneratedResonances << endl ;
585 cout << "Chi2Cut for muon tracks = " << Chi2Cut << endl;
586 cout << "PtCutMin for muon tracks = " << PtCutMin << endl;
587 cout << "PtCutMax for muon tracks = " << PtCutMax << endl;
588 cout << "Entries (unlike sign dimuons) in the mass range ["<<invMassMinInPeak<<";"<<invMassMaxInPeak<<"] : " << EventInMass <<endl;
589 if (ptTrig==0x800) cout << "Unlike Pair - All Pt" ;
590 if (ptTrig==0x400) cout << "Unlike Pair - High Pt" ;
591 if (ptTrig==0x200) cout << "Unlike Pair - Low Pt" ;
592 cout << " triggers : " << NbTrigger << endl;
594 cout << "Entries in the mass range with matching between reconstructed tracks and trigger tracks " << EventInMassMatch << endl;