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 "AliMUONCDB.h"
37 #include "AliMUONTrackParam.h"
38 #include "AliMUONTrackExtrap.h"
39 #include "AliMUONESDInterface.h"
43 #include "AliRunLoader.h"
44 #include "AliHeader.h"
45 #include "AliLoader.h"
47 #include "AliESDEvent.h"
48 #include "AliESDVertex.h"
49 #include "AliCDBManager.h"
50 #include "AliESDMuonTrack.h"
55 #include "TLorentzVector.h"
59 #include "TParticle.h"
61 #include <Riostream.h>
62 #include <TGeoManager.h>
69 Bool_t MUONefficiency(char* filename = "generated/galice.root", char* esdFileName = "AliESDs.root",
70 char* geoFilename = "geometry.root", char* ocdbPath = "local://$ALICE_ROOT/OCDB",
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 Double_t MUON_MASS = 0.105658369;
86 Double_t UPSILON_MASS = 9.4603 ;
87 Double_t JPSI_MASS = 3.097;
89 // Upper and lower bound for counting entries in the mass peak
90 // +/- 300 MeV/c^2 in this case.
91 Float_t countingRange = 0.300 ;
93 Float_t massResonance = 5.;
94 Float_t invMassMinInPeak = 0. ;
95 Float_t invMassMaxInPeak = 0. ;
97 Float_t nBinsPerGev = 40 ;
98 Float_t invMassMin = 0; Float_t invMassMax = 20;
99 Float_t ptMinResonance = 0 ; Float_t ptMaxResonance = 20 ; Int_t ptBinsResonance = 100;
102 massResonance = JPSI_MASS ;
103 invMassMinInPeak = JPSI_MASS - countingRange ; invMassMaxInPeak = JPSI_MASS + countingRange ;
104 //limits for histograms
105 invMassMin = 0 ; invMassMax = 6.;
106 ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100;
109 massResonance = UPSILON_MASS;
110 invMassMinInPeak = UPSILON_MASS - countingRange ; invMassMaxInPeak = UPSILON_MASS + countingRange;
111 //limits for histograms
112 invMassMin = 0 ; invMassMax = 12.;
113 ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100;
116 // Single Tracks muon cuts
117 Float_t Chi2Cut = 100.;
118 Float_t PtCutMin = 0. ;
119 Float_t PtCutMax = 10000. ;
122 // Limits for histograms
123 Float_t ptMinMuon = 0. ; Float_t ptMaxMuon = 20.; Int_t ptBinsMuon = 100 ;
124 Float_t pMinMuon = 0. ; Float_t pMaxMuon = 200.; Int_t pBinsMuon = 100 ;
127 //Reset ROOT and connect tree file
131 Int_t PRINTLEVEL = 0 ;
133 //for kinematic, i.e. reference tracks
134 TNtuple *Ktuple = new TNtuple("Ktuple","Kinematics NTuple","ev:npart:id:idmo:idgdmo:p:pt:y:theta:pseudorap:vx:vy:vz");
137 TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
138 TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
139 TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
140 TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", pBinsMuon, pMinMuon, pMaxMuon);
144 TH2F *hInvMassAll_vs_Pt;
145 TH2F *hInvMassBgk_vs_Pt;
149 hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax - invMassMin)), invMassMin, invMassMax);
150 hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax);
151 hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
152 hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
154 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Resonance",(Int_t) (nBinsPerGev*3*countingRange*2),massResonance-3*countingRange,massResonance+3*countingRange);
156 TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
157 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
158 TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
159 TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
160 TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
161 TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
162 TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
163 TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
165 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");
166 TNtuple *ESDtupleBck = new TNtuple("ESDtupleBck","Reconstructed Mu+Mu- pairs for Background","ev:pt:y:theta:minv:pt1:y1:theta1:pt2:y2:theta2");
170 Int_t EventInMass = 0;
171 Int_t EventInMassMatch = 0;
181 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
182 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
183 Int_t fCharge1, fCharge2;
185 Int_t ntrackhits, nevents;
186 Int_t nprocessedevents = 0 ;
189 TLorentzVector fV1, fV2, fVtot;
191 // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
193 TGeoManager::Import(geoFilename);
195 Error("MUONefficiency", "getting geometry from file %s failed", geoFilename);
201 TFile* esdFile = TFile::Open(esdFileName);
202 if (!esdFile || !esdFile->IsOpen()) {
203 Error("MUONefficiency", "opening ESD file %s failed", esdFileName);
207 AliESDEvent* esd = new AliESDEvent();
208 TTree* tree = (TTree*) esdFile->Get("esdTree");
210 Error("MUONefficiency", "no ESD tree found");
213 esd->ReadFromTree(tree);
216 if (tree->GetEvent(0) <= 0) {
217 Error("MUONefficiency", "no ESD object found for event 0");
220 Int_t runNumber = esd->GetRunNumber();
222 // load necessary data from OCDB
223 AliCDBManager::Instance()->SetDefaultStorage(ocdbPath);
224 AliCDBManager::Instance()->SetRun(runNumber);
225 if (!AliMUONCDB::LoadField()) return kFALSE;
227 // set the magnetic field for track extrapolations
228 AliMUONTrackExtrap::SetField();
230 // open run loader and load gAlice, kinematics and header
231 AliRunLoader* runLoader = AliRunLoader::Open(filename);
233 Error("MUONefficiency", "getting run loader from file %s failed", filename);
237 runLoader->LoadgAlice();
238 gAlice = runLoader->GetAliRun();
240 Error("MUONefficiency", "no galice object found");
244 runLoader->LoadHeader();
245 if (runNumber != runLoader->GetHeader()->GetRun()) {
246 Error("MUONefficiency", "mismatch between run number from ESD and from runLoader");
250 nevents = runLoader->GetNumberOfEvents();
251 AliMUONTrackParam trackParam;
253 // to access the particle Stack
254 runLoader->LoadKinematics("READ");
256 Int_t numberOfGeneratedResonances = 0 ;
260 Int_t track1Trigger = 0 ;
261 Float_t track1TriggerChi2 = 0 ;
262 Int_t track2Trigger = 0 ;
263 Float_t track2TriggerChi2 = 0 ;
267 for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { // Start event loop
269 if (iEvent%1000 == 0 )
270 printf("\n Nb of events analysed: %d \n",iEvent);
273 runLoader->GetEvent(iEvent);
276 // get the stack and fill the kine tree
277 AliStack *theStack = runLoader->Stack();
278 if (PRINTLEVEL > 0) theStack->DumpPStack ();
280 Int_t nparticles = (Int_t)runLoader->TreeK()->GetEntries();
281 Int_t nprimarypart = theStack->GetNprimary();
282 Int_t ntracks = theStack->GetNtrack();
284 if (PRINTLEVEL || (iEvent%100==0)) printf("\n >>> Event %d \n",iEvent);
285 if (PRINTLEVEL) cout << nprimarypart << " Particles generated (total is " << ntracks << ")"<< endl ;
287 for(Int_t iparticle=0; iparticle<nparticles; iparticle++) { // Start loop over particles
288 particle = theStack->Particle(iparticle);
290 Int_t muId = particle->GetPdgCode();
291 Int_t muM = particle->GetFirstMother();
293 Float_t muP = particle->P();
294 Float_t muPt = TMath::Sqrt(particle->Px()*particle->Px()+particle->Py()*particle->Py());
295 Float_t muY = 0.5*TMath::Log((particle->Energy()+particle->Pz()+1.e-13)/(particle->Energy()-particle->Pz()+1.e-13));
297 TParticle *theMum = theStack->Particle(muM);
298 muM = theMum->GetPdgCode();
299 muGM = theMum->GetFirstMother() ;
301 TParticle *grandMa = theStack->Particle(muGM);
302 muGM = grandMa->GetPdgCode();
308 if (muId==ResType) numberOfGeneratedResonances++;
311 Float_t muT = particle->Theta()*180/TMath::Pi();
312 Float_t muE = particle->Eta();
314 Float_t muVx = particle->Vx();
315 Float_t muVy = particle->Vy();
316 Float_t muVz = particle->Vz();
318 // If a write error occurs, the number of bytes returned is -1.
319 // If no data are written, because e.g. the branch is disabled,
320 // the number of bytes returned is 0.
321 Int_t errCode = Ktuple->Fill(iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz);
322 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);
324 } // End loop over particles
328 // get the event summary data
329 if (tree->GetEvent(iEvent) <= 0) {
330 Error("CheckESD", "no ESD object found for event %d", iEvent);
334 // get the SPD reconstructed vertex (vertexer) and fill the histogram
335 AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
336 if (Vertex->GetNContributors()) {
337 fZVertex = Vertex->GetZv();
338 fYVertex = Vertex->GetYv();
339 fXVertex = Vertex->GetXv();
340 errXVtx = Vertex->GetXRes();
341 errYVtx = Vertex->GetYRes();
343 hPrimaryVertex->Fill(fZVertex);
345 Int_t triggerWord = esd->GetTriggerMask();
346 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
349 printf("\n Nb of events analysed: %d \n",iEvent);
350 cout << " number of tracks: " << nTracks <<endl;
353 // loop over all reconstructed tracks (also first track of combination)
354 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
357 if (!esd->GetMuonTrack(iTrack)->ContainTrackerData()) continue;
359 AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
361 // extrapolate to vertex if required and available
362 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
363 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam);
364 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
365 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack
366 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
367 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam);
368 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
369 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack
373 if (PRINTLEVEL > 5) cout << "MatchTrigger " << muonTrack->GetMatchTrigger() << " and Chi2 of matching tracks " << track1TriggerChi2 << endl ;
374 track1Trigger = muonTrack->GetMatchTrigger();
376 track1TriggerChi2 = muonTrack->GetChi2MatchTrigger();
378 track1TriggerChi2 = 0. ;
380 fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
382 muonTrack->LorentzP(fV1);
384 ntrackhits = muonTrack->GetNHit();
385 fitfmin = muonTrack->GetChi2();
387 // transverse momentum
388 Float_t pt1 = fV1.Pt();
391 Float_t p1 = fV1.P();
394 Float_t rapMuon1 = fV1.Rapidity();
398 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
399 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);
402 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { // condition for good track (Chi2Cut and PtCut)
403 if (PRINTLEVEL > 8) cout << "inside pt and chi2 cuts " << endl ;
405 // fill histos hPtMuon and hChi2PerDof
408 hChi2PerDof->Fill(ch1);
409 hRapMuon->Fill(rapMuon1);
412 hPtMuonPlus->Fill(pt1);
413 hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
415 hPtMuonMinus->Fill(pt1);
416 hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
419 // loop over second track of combination
420 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
423 if (!esd->GetMuonTrack(iTrack2)->ContainTrackerData()) continue;
425 AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
427 // extrapolate to vertex if required and available
428 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
429 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam);
430 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
431 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
432 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
433 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam);
434 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
435 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
438 track2Trigger = muonTrack2->GetMatchTrigger();
440 track2TriggerChi2 = muonTrack2->GetChi2MatchTrigger();
442 track2TriggerChi2 = 0. ;
444 fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
446 muonTrack2->LorentzP(fV2);
448 ntrackhits = muonTrack2->GetNHit();
449 fitfmin = muonTrack2->GetChi2();
451 // transverse momentum
452 Float_t pt2 = fV2.Pt();
455 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
458 // condition for good track (Chi2Cut and PtCut)
459 if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
461 // condition for opposite charges
462 if ((fCharge1 * fCharge2) == -1) {
464 if (PRINTLEVEL > 8) cout << "---------> Now filling the Ntuple " << endl ;
468 Float_t invMass = fVtot.M();
470 if (fCharge1 < 0){ //mu_minus is index 1 in the ntuple
471 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};
472 ESDtuple->Fill(ESDFill);
475 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};
476 ESDtuple->Fill(ESDFill);
479 // fill histos hInvMassAll and hInvMassRes
480 hInvMassAll->Fill(invMass);
481 hInvMassRes->Fill(invMass);
482 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
486 ptTrig = 0x08;// mask for Hpt unlike sign pair
487 else if (ResType == 443)
488 ptTrig = 0x04;// mask for Lpt unlike sign pair
491 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
493 if (invMass > invMassMinInPeak && invMass < invMassMaxInPeak) {
495 hRapResonance->Fill(fVtot.Rapidity());
496 hPtResonance->Fill(fVtot.Pt());
498 // match with trigger
499 if (muonTrack2->GetMatchTrigger()>=0 && (esd->GetTriggerMask() & ptTrig)) EventInMassMatch++;
503 } // if (fCharge1 * fCharge2) == -1)
504 } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
506 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
507 } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
509 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
511 hNumberOfTrack->Fill(nTracks);
512 // esdFile->Delete();
514 } // End of event loop
517 // Loop over events for bg event
519 Double_t thetaPlus, phiPlus;
520 Double_t thetaMinus, phiMinus;
521 Float_t PtMinus, PtPlus;
523 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { // Loop over events for bg event
524 // according to Christian a 3d phi-theta-pt random pick would take better care
525 // of all correlations
527 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
528 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
529 PtPlus = hPtMuonPlus->GetRandom();
530 PtMinus = hPtMuonMinus->GetRandom();
532 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
533 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
534 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
536 fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
537 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
539 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
540 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
541 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
543 fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
544 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
549 // fill histos hInvMassAll and hInvMassRes
550 hInvMassBg->Fill(fVtot.M());
551 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
553 // Ntuple for background... more convenient
554 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);
556 } // End loop over events for background
559 // File for histograms and histogram booking
560 TString outfilename = "MUONefficiency.root";
561 TFile *ntupleFile = new TFile(outfilename.Data(), "RECREATE");
565 ESDtupleBck->Write();
569 TFile *histoFile = new TFile("MUONhistos.root", "RECREATE");
570 hPrimaryVertex->Write();
572 hPtMuonPlus->Write();
573 hPtMuonMinus->Write();
575 hChi2PerDof->Write();
576 hInvMassAll->Write();
578 hInvMassAll_vs_Pt ->Write();
579 hInvMassBgk_vs_Pt->Write();
580 hInvMassRes->Write();
581 hNumberOfTrack->Write();
583 hRapResonance ->Write();
584 hPtResonance ->Write();
585 hThetaPhiPlus ->Write();
586 hThetaPhiMinus ->Write();
590 cout << "*************************************************" << endl;
592 cout << "MUONefficiency : " << nprocessedevents << " events processed" << endl;
594 cout << "Number of generated J/Psi (443) : " << numberOfGeneratedResonances << endl ;
596 cout << "Number of generated Upsilon (553) :" << numberOfGeneratedResonances << endl ;
597 cout << "Chi2Cut for muon tracks = " << Chi2Cut << endl;
598 cout << "PtCutMin for muon tracks = " << PtCutMin << endl;
599 cout << "PtCutMax for muon tracks = " << PtCutMax << endl;
601 cout << "Entries (unlike sign dimuons) : " << hInvMassAll->GetEntries();
603 if (hInvMassAll->GetEntries() > 0) {
604 hInvMassAll->Fit("gaus","q0");
605 TF1* f1 = hInvMassAll->GetFunction("gaus");
606 cout << Form(". Rough sigma = %7.2f MeV/c2",f1->GetParameter(2)*1000.0);
609 cout << endl << "Entries (unlike sign dimuons) in the mass range ["<<invMassMinInPeak<<";"<<invMassMaxInPeak<<"] : " << EventInMass <<endl;
611 if (ptTrig==0x800) cout << "Unlike Pair - All Pt" ;
612 if (ptTrig==0x400) cout << "Unlike Pair - High Pt" ;
613 if (ptTrig==0x200) cout << "Unlike Pair - Low Pt" ;
614 cout << " triggers : " << NbTrigger << endl;
616 cout << "Entries in the mass range with matching between reconstructed tracks and trigger tracks " << EventInMassMatch << endl;