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 **************************************************************************/
19 /// \file MUONefficiency.C
20 /// \brief add brief description
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
29 /// \author Christophe Suire, IPN Orsay
33 #if !defined(__CINT__) || defined(__MAKECINT__)
36 #include "AliMUONTrackParam.h"
37 #include "AliMUONTrackExtrap.h"
38 #include "AliESDMuonTrack.h"
42 #include "AliRunLoader.h"
43 #include "AliHeader.h"
44 #include "AliLoader.h"
46 #include "AliMagFMaps.h"
47 #include "AliESDEvent.h"
48 #include "AliESDVertex.h"
49 #include "AliTracker.h"
50 #include "AliCDBManager.h"
56 #include "TClonesArray.h"
57 #include "TLorentzVector.h"
61 #include "TParticle.h"
64 #include <Riostream.h>
65 #include <TGeoManager.h>
70 Bool_t MUONefficiency( char* filename = "galice.root", char* geoFilename = "geometry.root", char* esdFileName = "AliESDs.root",
71 Int_t ExtrapToVertex = -1, Int_t ResType = 553, Int_t FirstEvent = 0, Int_t LastEvent = 1000000 )
73 /// \param ExtrapToVertex (default -1)
74 /// - <0: no extrapolation;
75 /// - =0: extrapolation to (0,0,0);
76 /// - >0: extrapolation to ESDVertex if available, else to (0,0,0)
77 /// \param ResType 553 for Upsilon, 443 for J/Psi (default 553)
78 /// \param FirstEvent (default 0)
79 /// \param LastEvent (default 1.e6)
80 /// \param Chi2Cut to keep only tracks with chi2 per d.o.f. < Chi2Cut (default 100)
83 // MUONefficiency starts
85 // Set default CDB storage
86 AliCDBManager* man = AliCDBManager::Instance();
87 man->SetDefaultStorage("local://$ALICE_ROOT");
89 Double_t MUON_MASS = 0.105658369;
90 Double_t UPSILON_MASS = 9.4603 ;
91 Double_t JPSI_MASS = 3.097;
93 // Upper and lower bound for counting entries in the mass peak
94 // +/- 300 MeV/c^2 in this case.
95 Float_t countingRange = 0.300 ;
97 Float_t massResonance = 5.;
98 Float_t invMassMinInPeak = 0. ;
99 Float_t invMassMaxInPeak = 0. ;
101 Float_t nBinsPerGev = 40 ;
102 Float_t invMassMin = 0; Float_t invMassMax = 20;
103 Float_t ptMinResonance = 0 ; Float_t ptMaxResonance = 20 ; Int_t ptBinsResonance = 100;
106 massResonance = JPSI_MASS ;
107 invMassMinInPeak = JPSI_MASS - countingRange ; invMassMaxInPeak = JPSI_MASS + countingRange ;
108 //limits for histograms
109 invMassMin = 0 ; invMassMax = 6.;
110 ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100;
113 massResonance = UPSILON_MASS;
114 invMassMinInPeak = UPSILON_MASS - countingRange ; invMassMaxInPeak = UPSILON_MASS + countingRange;
115 //limits for histograms
116 invMassMin = 0 ; invMassMax = 12.;
117 ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100;
120 // Single Tracks muon cuts
121 Float_t Chi2Cut = 100.;
122 Float_t PtCutMin = 0. ;
123 Float_t PtCutMax = 10000. ;
126 // Limits for histograms
127 Float_t ptMinMuon = 0. ; Float_t ptMaxMuon = 20.; Int_t ptBinsMuon = 100 ;
128 Float_t pMinMuon = 0. ; Float_t pMaxMuon = 200.; Int_t pBinsMuon = 100 ;
131 //Reset ROOT and connect tree file
135 Int_t PRINTLEVEL = 0 ;
137 //for kinematic, i.e. reference tracks
138 TNtuple *Ktuple = new TNtuple("Ktuple","Kinematics NTuple","ev:npart:id:idmo:idgdmo:p:pt:y:theta:pseudorap:vx:vy:vz");
141 TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
142 TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
143 TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
144 TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", pBinsMuon, pMinMuon, pMaxMuon);
148 TH2F *hInvMassAll_vs_Pt;
149 TH2F *hInvMassBgk_vs_Pt;
153 hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax - invMassMin)), invMassMin, invMassMax);
154 hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax);
155 hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
156 hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
158 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 TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
161 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
162 TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
163 TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
164 TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
165 TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
166 TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
167 TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
169 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");
170 TNtuple *ESDtupleBck = new TNtuple("ESDtupleBck","Reconstructed Mu+Mu- pairs for Background","ev:pt:y:theta:minv:pt1:y1:theta1:pt2:y2:theta2");
174 Int_t EventInMass = 0;
175 Int_t EventInMassMatch = 0;
185 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
186 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
187 Int_t fCharge1, fCharge2;
189 Int_t ntrackhits, nevents;
190 Int_t nprocessedevents = 0 ;
193 TLorentzVector fV1, fV2, fVtot;
195 // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
197 TGeoManager::Import(geoFilename);
199 Error("MUONmass_ESD", "getting geometry from file %s failed", filename);
205 // waiting for mag field in CDB
206 printf("Loading field map...\n");
207 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
208 AliTracker::SetFieldMap(field, kFALSE);
210 // open run loader and load gAlice, kinematics and header
211 AliRunLoader* runLoader = AliRunLoader::Open(filename);
213 Error("MUONefficiency", "getting run loader from file %s failed", filename);
217 runLoader->LoadgAlice();
218 gAlice = runLoader->GetAliRun();
220 Error("MUONefficiency", "no galice object found");
225 TFile* esdFile = TFile::Open(esdFileName);
226 if (!esdFile || !esdFile->IsOpen()) {
227 Error("MUONefficiency", "opening ESD file %s failed", esdFileName);
231 AliESDEvent* esd = new AliESDEvent();
232 TTree* tree = (TTree*) esdFile->Get("esdTree");
234 Error("CheckESD", "no ESD tree found");
237 esd->ReadFromTree(tree);
239 runLoader->LoadHeader();
240 Int_t runNumber = runLoader->GetHeader()->GetRun();
241 AliCDBManager::Instance()->SetRun(runNumber);
243 nevents = runLoader->GetNumberOfEvents();
244 AliMUONTrackParam trackParam;
246 // to access the particle Stack
247 runLoader->LoadKinematics("READ");
249 Int_t numberOfGeneratedResonances = 0 ;
253 Int_t track1Trigger = 0 ;
254 Float_t track1TriggerChi2 = 0 ;
255 Int_t track2Trigger = 0 ;
256 Float_t track2TriggerChi2 = 0 ;
260 for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { // Start event loop
262 if (iEvent%1000 == 0 )
263 printf("\n Nb of events analysed: %d \n",iEvent);
266 runLoader->GetEvent(iEvent);
269 // get the stack and fill the kine tree
270 AliStack *theStack = runLoader->Stack();
271 if (PRINTLEVEL > 0) theStack->DumpPStack ();
273 Int_t nparticles = (Int_t)runLoader->TreeK()->GetEntries();
274 Int_t nprimarypart = theStack->GetNprimary();
275 Int_t ntracks = theStack->GetNtrack();
277 if (PRINTLEVEL || (iEvent%100==0)) printf("\n >>> Event %d \n",iEvent);
278 if (PRINTLEVEL) cout << nprimarypart << " Particles generated (total is " << ntracks << ")"<< endl ;
280 for(Int_t iparticle=0; iparticle<nparticles; iparticle++) { // Start loop over particles
281 particle = theStack->Particle(iparticle);
283 Int_t muId = particle->GetPdgCode();
284 Int_t muM = particle->GetFirstMother();
286 Float_t muP = particle->P();
287 Float_t muPt = TMath::Sqrt(particle->Px()*particle->Px()+particle->Py()*particle->Py());
288 Float_t muY = 0.5*TMath::Log((particle->Energy()+particle->Pz()+1.e-13)/(particle->Energy()-particle->Pz()+1.e-13));
290 TParticle *theMum = theStack->Particle(muM);
291 muM = theMum->GetPdgCode();
292 muGM = theMum->GetFirstMother() ;
294 TParticle *grandMa = theStack->Particle(muGM);
295 muGM = grandMa->GetPdgCode();
301 if (muId==ResType) numberOfGeneratedResonances++;
304 Float_t muT = particle->Theta()*180/TMath::Pi();
305 Float_t muE = particle->Eta();
307 Float_t muVx = particle->Vx();
308 Float_t muVy = particle->Vy();
309 Float_t muVz = particle->Vz();
311 // If a write error occurs, the number of bytes returned is -1.
312 // If no data are written, because e.g. the branch is disabled,
313 // the number of bytes returned is 0.
314 Int_t errCode = Ktuple->Fill(iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz);
315 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 } // End loop over particles
321 // get the event summary data
322 tree->GetEvent(iEvent);
324 Error("CheckESD", "no ESD object found for event %d", iEvent);
328 // get the SPD reconstructed vertex (vertexer) and fill the histogram
329 AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
330 if (Vertex->GetNContributors()) {
331 fZVertex = Vertex->GetZv();
332 fYVertex = Vertex->GetYv();
333 fXVertex = Vertex->GetXv();
334 errXVtx = Vertex->GetXRes();
335 errYVtx = Vertex->GetYRes();
337 hPrimaryVertex->Fill(fZVertex);
339 Int_t triggerWord = esd->GetTriggerMask();
340 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
343 printf("\n Nb of events analysed: %d \n",iEvent);
344 cout << " number of tracks: " << nTracks <<endl;
347 // set the magnetic field for track extrapolations
348 AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
349 // loop over all reconstructed tracks (also first track of combination)
350 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
352 AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
354 // extrapolate to vertex if required and available
355 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
356 trackParam.GetParamFromUncorrected(*muonTrack);
357 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
358 trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
359 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
360 trackParam.GetParamFromUncorrected(*muonTrack);
361 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
362 trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
366 if (PRINTLEVEL > 5) cout << "MatchTrigger " << muonTrack->GetMatchTrigger() << " and Chi2 of matching tracks " << track1TriggerChi2 << endl ;
367 track1Trigger = muonTrack->GetMatchTrigger();
369 track1TriggerChi2 = muonTrack->GetChi2MatchTrigger();
371 track1TriggerChi2 = 0. ;
373 fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
375 muonTrack->LorentzP(fV1);
377 ntrackhits = muonTrack->GetNHit();
378 fitfmin = muonTrack->GetChi2();
380 // transverse momentum
381 Float_t pt1 = fV1.Pt();
384 Float_t p1 = fV1.P();
387 Float_t rapMuon1 = fV1.Rapidity();
391 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
392 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);
395 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { // condition for good track (Chi2Cut and PtCut)
396 if (PRINTLEVEL > 8) cout << "inside pt and chi2 cuts " << endl ;
398 // fill histos hPtMuon and hChi2PerDof
401 hChi2PerDof->Fill(ch1);
402 hRapMuon->Fill(rapMuon1);
405 hPtMuonPlus->Fill(pt1);
406 hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
408 hPtMuonMinus->Fill(pt1);
409 hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
412 // loop over second track of combination
413 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
415 AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
417 // extrapolate to vertex if required and available
418 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
419 trackParam.GetParamFromUncorrected(*muonTrack2);
420 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
421 trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
422 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
423 trackParam.GetParamFromUncorrected(*muonTrack2);
424 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
425 trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
428 track2Trigger = muonTrack2->GetMatchTrigger();
430 track2TriggerChi2 = muonTrack2->GetChi2MatchTrigger();
432 track2TriggerChi2 = 0. ;
434 fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
436 muonTrack2->LorentzP(fV2);
438 ntrackhits = muonTrack2->GetNHit();
439 fitfmin = muonTrack2->GetChi2();
441 // transverse momentum
442 Float_t pt2 = fV2.Pt();
445 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
448 // condition for good track (Chi2Cut and PtCut)
449 if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
451 // condition for opposite charges
452 if ((fCharge1 * fCharge2) == -1) {
454 if (PRINTLEVEL > 8) cout << "---------> Now filling the Ntuple " << endl ;
458 Float_t invMass = fVtot.M();
460 if (fCharge1 < 0){ //mu_minus is index 1 in the ntuple
461 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};
462 ESDtuple->Fill(ESDFill);
465 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};
466 ESDtuple->Fill(ESDFill);
469 // fill histos hInvMassAll and hInvMassRes
470 hInvMassAll->Fill(invMass);
471 hInvMassRes->Fill(invMass);
472 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
476 ptTrig = 0x08;// mask for Hpt unlike sign pair
477 else if (ResType == 443)
478 ptTrig = 0x04;// mask for Lpt unlike sign pair
481 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
483 if (invMass > invMassMinInPeak && invMass < invMassMaxInPeak) {
485 hRapResonance->Fill(fVtot.Rapidity());
486 hPtResonance->Fill(fVtot.Pt());
488 // match with trigger
489 if (muonTrack2->GetMatchTrigger()>=0 && (esd->GetTriggerMask() & ptTrig)) EventInMassMatch++;
493 } // if (fCharge1 * fCharge2) == -1)
494 } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
496 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
497 } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
499 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
501 hNumberOfTrack->Fill(nTracks);
502 // esdFile->Delete();
504 } // End of event loop
507 // Loop over events for bg event
509 Double_t thetaPlus, phiPlus;
510 Double_t thetaMinus, phiMinus;
511 Float_t PtMinus, PtPlus;
513 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { // Loop over events for bg event
514 // according to Christian a 3d phi-theta-pt random pick would take better care
515 // of all correlations
517 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
518 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
519 PtPlus = hPtMuonPlus->GetRandom();
520 PtMinus = hPtMuonMinus->GetRandom();
522 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
523 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
524 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
526 fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
527 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
529 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
530 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
531 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
533 fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
534 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
539 // fill histos hInvMassAll and hInvMassRes
540 hInvMassBg->Fill(fVtot.M());
541 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
543 // Ntuple for background... more convenient
544 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 } // End loop over events for background
549 // File for histograms and histogram booking
550 TString outfilename = "MUONefficiency.root";
551 TFile *ntupleFile = new TFile(outfilename.Data(), "RECREATE");
555 ESDtupleBck->Write();
559 TFile *histoFile = new TFile("MUONhistos.root", "RECREATE");
560 hPrimaryVertex->Write();
562 hPtMuonPlus->Write();
563 hPtMuonMinus->Write();
565 hChi2PerDof->Write();
566 hInvMassAll->Write();
568 hInvMassAll_vs_Pt ->Write();
569 hInvMassBgk_vs_Pt->Write();
570 hInvMassRes->Write();
571 hNumberOfTrack->Write();
573 hRapResonance ->Write();
574 hPtResonance ->Write();
575 hThetaPhiPlus ->Write();
576 hThetaPhiMinus ->Write();
580 cout << "*************************************************" << endl;
582 cout << "MUONefficiency : " << nprocessedevents << " events processed" << endl;
584 cout << "Number of generated J/Psi (443) : " << numberOfGeneratedResonances << endl ;
586 cout << "Number of generated Upsilon (553) :" << numberOfGeneratedResonances << endl ;
587 cout << "Chi2Cut for muon tracks = " << Chi2Cut << endl;
588 cout << "PtCutMin for muon tracks = " << PtCutMin << endl;
589 cout << "PtCutMax for muon tracks = " << PtCutMax << endl;
590 cout << "Entries (unlike sign dimuons) in the mass range ["<<invMassMinInPeak<<";"<<invMassMaxInPeak<<"] : " << EventInMass <<endl;
591 if (ptTrig==0x800) cout << "Unlike Pair - All Pt" ;
592 if (ptTrig==0x400) cout << "Unlike Pair - High Pt" ;
593 if (ptTrig==0x200) cout << "Unlike Pair - Low Pt" ;
594 cout << " triggers : " << NbTrigger << endl;
596 cout << "Entries in the mass range with matching between reconstructed tracks and trigger tracks " << EventInMassMatch << endl;