]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONefficiency.C
- Adding possibility to handle Jpsi (default is Upsilon) in efficiency
[u/mrichter/AliRoot.git] / MUON / MUONefficiency.C
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 // 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
24
25 // Christophe Suire, IPN Orsay
26
27
28
29 // Arguments:
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
36
37
38
39 #if !defined(__CINT__) || defined(__MAKECINT__)
40 // ROOT includes
41 #include "TTree.h"
42 #include "TNtuple.h"
43 #include "TBranch.h"
44 #include "TClonesArray.h"
45 #include "TLorentzVector.h"
46 #include "TFile.h"
47 #include "TH1.h"
48 #include "TH2.h"
49 #include "TParticle.h"
50 #include "TTree.h"
51 #include "TString.h"
52 #include <Riostream.h>
53
54 // STEER includes
55 #include "AliRun.h"
56 #include "AliRunLoader.h"
57 #include "AliHeader.h"
58 #include "AliLoader.h"
59 #include "AliStack.h"
60 #include "AliMagF.h"
61 #include "AliESD.h"
62
63 // MUON includes
64 #include "AliESDMuonTrack.h"
65 #endif
66
67
68 Bool_t MUONefficiency( Int_t ResType = 553, Int_t FirstEvent = 0, Int_t LastEvent = 1000000,
69                        char* esdFileName = "AliESDs.root", char* filename = "galice.root")
70 { // MUONefficiency starts
71
72   Double_t MUON_MASS = 0.105658369;
73   Double_t UPSILON_MASS = 9.4603 ;
74   Double_t JPSI_MASS = 3.097;
75
76   // Upper and lower bound for counting entries in the mass peak
77   // +/- 300 MeV/c^2 in this case.
78   Float_t countingRange = 0.300 ;  
79   
80   Float_t massResonance = 5.;
81   Float_t invMassMinInPeak = 0. ; 
82   Float_t invMassMaxInPeak = 0. ; 
83   
84   Float_t nBinsPerGev = 40 ; 
85   Float_t invMassMin = 0;   Float_t invMassMax = 20; 
86   Float_t ptMinResonance = 0 ; Float_t ptMaxResonance = 20 ; Int_t ptBinsResonance = 100; 
87
88   if (ResType==443) {
89     massResonance = JPSI_MASS ;
90     invMassMinInPeak =  JPSI_MASS - countingRange  ; invMassMaxInPeak = JPSI_MASS + countingRange ; 
91     //limits for histograms
92     invMassMin = 0 ; invMassMax = 6.;
93     ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100; 
94   }
95   if (ResType==553) {
96     massResonance = UPSILON_MASS;
97     invMassMinInPeak = UPSILON_MASS - countingRange ; invMassMaxInPeak = UPSILON_MASS + countingRange; 
98     //limits for histograms 
99     invMassMin = 0 ; invMassMax = 12.;
100     ptMinResonance = 0 ; ptMaxResonance = 20 ; ptBinsResonance = 100; 
101   }
102   
103   // Single Tracks muon cuts
104   Float_t Chi2Cut = 100.;
105   Float_t PtCutMin = 0. ;
106   Float_t PtCutMax = 10000. ; 
107
108
109   // Limits for histograms 
110   Float_t ptMinMuon = 0. ; Float_t ptMaxMuon = 20.; Int_t ptBinsMuon = 100 ;
111   Float_t pMinMuon = 0.  ; Float_t pMaxMuon = 200.; Int_t pBinsMuon = 100 ;
112  
113
114   //Reset ROOT and connect tree file
115   gROOT->Reset();
116   
117   // Printing Level 
118   Int_t PRINTLEVEL = 0 ;
119   Int_t SELECT =  0 ; // not used
120   
121   //for kinematic, i.e. reference tracks
122   TNtuple *Ktuple = new TNtuple("Ktuple","Kinematics NTuple","ev:npart:id:idmo:idgdmo:p:pt:y:theta:pseudorap:vx:vy:vz");
123
124   //for reconstruction  
125   TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
126   TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
127   TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", ptBinsMuon, ptMinMuon, ptMaxMuon);
128   TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", pBinsMuon, pMinMuon, pMaxMuon);
129   
130   TH1F *hInvMassAll;
131   TH1F *hInvMassBg;
132   TH2F *hInvMassAll_vs_Pt;
133   TH2F *hInvMassBgk_vs_Pt;
134   TH1F *hInvMassRes;
135
136
137   hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax - invMassMin)), invMassMin, invMassMax);
138   hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", (Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax);
139   hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
140   hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",(Int_t) (nBinsPerGev*(invMassMax- invMassMin)), invMassMin, invMassMax,ptBinsResonance,ptMinResonance,ptMaxResonance);
141   
142   hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Resonance",(Int_t) (nBinsPerGev*3*countingRange*2),massResonance-3*countingRange,massResonance+3*countingRange);
143
144   
145   TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
146   TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
147   TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
148   TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
149   TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
150   TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
151   TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
152   
153   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");
154   TNtuple *ESDtupleBck = new TNtuple("ESDtupleBck","Reconstructed Mu+Mu- pairs for Background","ev:pt:y:theta:minv:pt1:y1:theta1:pt2:y2:theta2");
155
156
157   // Variables
158   Int_t EventInMass = 0;
159   Int_t EventInMassMatch = 0;
160   Int_t NbTrigger = 0;
161   Int_t ptTrig = 0;
162
163   Double_t thetaX, thetaY, pYZ;
164   Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
165   Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
166   Int_t fCharge1, fCharge2;
167
168   Int_t ntrackhits, nevents; 
169   Int_t nprocessedevents = 0 ;
170   Double_t fitfmin;
171
172   TLorentzVector fV1, fV2, fVtot;
173
174   // set off mag field 
175   AliMagF::SetReadField(kFALSE);
176
177   // open run loader and load gAlice, kinematics and header
178   AliRunLoader* runLoader = AliRunLoader::Open(filename);
179   if (!runLoader) {
180     Error("MUONefficiency", "getting run loader from file %s failed", filename);
181     return kFALSE;
182   }
183
184   runLoader->LoadgAlice();
185   gAlice = runLoader->GetAliRun();
186   if (!gAlice) {
187     Error("MUONefficiency", "no galice object found");
188     return kFALSE;
189   }
190   
191   // open the ESD file
192   TFile* esdFile = TFile::Open(esdFileName);
193   if (!esdFile || !esdFile->IsOpen()) {
194     Error("MUONefficiency", "opening ESD file %s failed", esdFileName);
195     return kFALSE;
196   }
197   
198   AliESD* esd = new AliESD();
199   TTree* tree = (TTree*) esdFile->Get("esdTree");
200   if (!tree) {
201     Error("CheckESD", "no ESD tree found");
202     return kFALSE;
203   }
204   tree->SetBranchAddress("ESD", &esd);
205
206   runLoader->LoadHeader();
207   nevents = runLoader->GetNumberOfEvents();
208   
209   // to access the particle  Stack
210   runLoader->LoadKinematics("READ");
211
212   Int_t numberOfGeneratedResonances = 0 ;
213
214   TParticle *particle; 
215   
216   Int_t track1Id = 0 ;
217   Int_t track1PDGId = 0 ;
218   Int_t track1MotherId = 0 ;
219   Int_t track1MotherPDGId = 0 ;
220   Int_t track1Trigger = 0 ;
221   Float_t track1TriggerChi2 = 0 ;
222   Int_t track2Id = 0 ;
223   Int_t track2PDGId = 0 ;
224   Int_t track2MotherId = 0 ;
225   Int_t track2MotherPDGId = 0 ;
226   Int_t track2Trigger = 0 ;
227   Float_t track2TriggerChi2 = 0 ;
228
229
230   // Loop over events
231   for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) { // Start event loop
232
233     if (iEvent%1000 == 0 )
234       printf("\n Nb of events analysed: %d \n",iEvent);
235
236     // get current event
237     runLoader->GetEvent(iEvent);
238     nprocessedevents++;
239
240     // get the stack and fill the kine tree
241     AliStack *theStack = runLoader->Stack();
242     if (PRINTLEVEL > 0) theStack->DumpPStack ();    
243     
244     Int_t nparticles = (Int_t)runLoader->TreeK()->GetEntries();
245     Int_t nprimarypart = theStack->GetNprimary();
246     Int_t ntracks = theStack->GetNtrack();
247   
248     if (PRINTLEVEL || (iEvent%100==0)) printf("\n  >>> Event %d \n",iEvent);
249     if (PRINTLEVEL) cout << nprimarypart << " Particles generated (total is " << ntracks << ")"<< endl ;
250     
251     
252     for(Int_t iparticle=0; iparticle<nparticles; iparticle++) { // Start loop over particles
253       particle = theStack->Particle(iparticle);
254       
255       Int_t   muId = particle->GetPdgCode(); 
256       Int_t   muM  = particle->GetFirstMother();
257       Int_t   muGM = 0;
258       Float_t muP  = particle->P();
259       Float_t muPt = TMath::Sqrt(particle->Px()*particle->Px()+particle->Py()*particle->Py());
260       Float_t muY  = 0.5*TMath::Log((particle->Energy()+particle->Pz()+1.e-13)/(particle->Energy()-particle->Pz()+1.e-13));
261       if (muM >= 0) {
262         TParticle *theMum = theStack->Particle(muM);
263         muM  =  theMum->GetPdgCode();
264         muGM  = theMum->GetFirstMother() ;
265         if (muGM >= 0){
266           TParticle *grandMa = theStack->Particle(muGM);
267           muGM = grandMa->GetPdgCode();
268         }
269         else muGM=0;
270       }
271       else muM=0;
272     
273       if (muId==ResType) numberOfGeneratedResonances++;
274
275   
276       Float_t muT  = particle->Theta()*180/TMath::Pi();
277       Float_t muE  = particle->Eta();
278       
279       Float_t muVx = particle->Vx();
280       Float_t muVy = particle->Vy();
281       Float_t muVz = particle->Vz();
282       
283       // If a write error occurs, the number of bytes returned is -1.
284       // If no data are written, because e.g. the branch is disabled,
285       // the number of bytes returned is 0.
286       Int_t errCode = Ktuple->Fill(iEvent,nparticles,muId,muM,muGM,muP,muPt,muY,muT,muE,muVx,muVy,muVz);
287       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);
288
289     } // End loop over particles
290     
291
292     
293     // get the event summary data
294     tree->GetEvent(iEvent);
295     if (!esd) {
296       Error("CheckESD", "no ESD object found for event %d", iEvent);
297       return kFALSE;
298     }
299
300     Int_t triggerWord = esd->GetTrigger();
301     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
302
303     if (PRINTLEVEL > 0){
304       printf("\n Nb of events analysed: %d \n",iEvent);
305       cout << " number of tracks: " << nTracks  <<endl;
306     }
307
308     // loop over all reconstructed tracks (also first track of combination)
309     for (Int_t iTrack = 0; iTrack <  nTracks;  iTrack++) {
310
311       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
312
313       //if (PRINTLEVEL > 5) cout << "1st muonTrack->GetTrackID() : " << track1Id << endl; 
314
315       if(SELECT && track1Id) {
316         particle = theStack->Particle(track1Id);
317         track1PDGId = particle->GetPdgCode() ;
318         track1MotherId = particle->GetFirstMother();
319         if (track1MotherId >=0 )
320           track1MotherPDGId = ((TParticle*) theStack->Particle(track1MotherId))->GetPdgCode();
321         if (PRINTLEVEL > 0) cout << "track1MotherPDGId = " << track1MotherPDGId << endl ;
322       }
323
324       // Trigger
325       if (PRINTLEVEL > 5) cout << "MatchTrigger " << muonTrack->GetMatchTrigger() << " and Chi2 of matching tracks " << track1TriggerChi2 <<   endl ;
326       track1Trigger = muonTrack->GetMatchTrigger();
327       if (track1Trigger)
328         track1TriggerChi2 = muonTrack->GetChi2MatchTrigger();
329       else 
330         track1TriggerChi2 = 0. ;
331
332       thetaX = muonTrack->GetThetaX();
333       thetaY = muonTrack->GetThetaY();
334
335       pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
336       fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
337       fPxRec1  = fPzRec1 * TMath::Tan(thetaX);
338       fPyRec1  = fPzRec1 * TMath::Tan(thetaY);
339       fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
340       
341       fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
342       fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
343
344       ntrackhits = muonTrack->GetNHit();
345       fitfmin    = muonTrack->GetChi2();
346
347       // transverse momentum
348       Float_t pt1 = fV1.Pt();
349       
350       // total momentum
351       Float_t p1 = fV1.P();
352       
353       // Rapidity
354       Float_t rapMuon1 = fV1.Rapidity();
355       
356       // chi2 per d.o.f.
357       
358       Float_t ch1 =  fitfmin / (2.0 * ntrackhits - 5);
359       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);
360       
361       
362       if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) { // condition for good track (Chi2Cut and PtCut)
363         if (PRINTLEVEL > 8) cout << "inside pt and chi2 cuts " << endl ; 
364         
365         // fill histos hPtMuon and hChi2PerDof
366         hPtMuon->Fill(pt1);
367         hPMuon->Fill(p1);
368         hChi2PerDof->Fill(ch1);
369         hRapMuon->Fill(rapMuon1);
370
371         if (fCharge1 > 0) {
372           hPtMuonPlus->Fill(pt1);
373           hThetaPhiPlus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
374         } else {
375           hPtMuonMinus->Fill(pt1);
376           hThetaPhiMinus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
377         }
378
379         // loop over second track of combination
380         for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
381           
382           AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack2);
383           //Int_t track2Id = muonTrack->GetTrackID();
384           //if (PRINTLEVEL > 5) cout << "2nd muonTrack->GetTrackID() : " << track2Id << endl;
385
386           if(SELECT && track2Id) {
387             particle = theStack->Particle(track2Id);
388             track2PDGId = particle->GetPdgCode();
389             track2MotherId = particle->GetFirstMother();
390             if (track2MotherId >=0 ) 
391               track2MotherPDGId = ((TParticle*) theStack->Particle(track2MotherId))->GetPdgCode();
392           }
393
394           track2Trigger = muonTrack->GetMatchTrigger();
395           if (track2Trigger) 
396             track2TriggerChi2 = muonTrack->GetChi2MatchTrigger();
397           else 
398             track2TriggerChi2 = 0. ;
399
400           thetaX = muonTrack->GetThetaX();
401           thetaY = muonTrack->GetThetaY();
402
403           pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
404           fPzRec2  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
405           fPxRec2  = fPzRec2 * TMath::Tan(thetaX);
406           fPyRec2  = fPzRec2 * TMath::Tan(thetaY);
407           fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
408
409           fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
410           fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
411
412           ntrackhits = muonTrack->GetNHit();
413           fitfmin    = muonTrack->GetChi2();
414
415           // transverse momentum
416           Float_t pt2 = fV2.Pt();
417
418           // chi2 per d.o.f.
419           Float_t ch2 = fitfmin  / (2.0 * ntrackhits - 5);
420
421           if (PRINTLEVEL > 5)  cout << "track1MotherId : "<< track1MotherId << "    track2MotherId : " << track2MotherId << endl ;
422           if (PRINTLEVEL > 5) cout << "track1MotherPDGId : " << track1MotherPDGId << "    track2MotherPDGId : "  << track2MotherPDGId << endl ;
423
424           // Select Condition
425           if (!SELECT || (track2MotherId == track1MotherId && track2MotherPDGId == ResType && TMath::Abs(track1PDGId)==13 && TMath::Abs(track2PDGId)==13 )) {
426
427             // condition for good track (Chi2Cut and PtCut)
428             if ((ch2 < Chi2Cut) && (pt2 > PtCutMin)  && (pt2 < PtCutMax)) {
429             
430               // condition for opposite charges
431               if ((fCharge1 * fCharge2) == -1) {
432
433                 if (PRINTLEVEL > 8) cout << "---------> Now filling the Ntuple "  <<  endl ;
434                 
435                 // invariant mass
436                 fVtot = fV1 + fV2;
437                 Float_t invMass = fVtot.M();
438                 
439                 if (fCharge1 < 0){ //mu_minus is index 1 in the ntuple
440                   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};
441                   ESDtuple->Fill(ESDFill);
442                 }
443                 else{
444                   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};
445                   ESDtuple->Fill(ESDFill);
446                 }
447                 
448                 // fill histos hInvMassAll and hInvMassRes
449                 hInvMassAll->Fill(invMass);
450                 hInvMassRes->Fill(invMass);
451                 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
452
453                 //trigger info 
454                 if (ResType == 553)
455                   ptTrig = 0x400;// mask for Hpt unlike sign pair
456                 else if (ResType == 443)
457                   ptTrig = 0x800;// mask for Apt unlike sign pair
458                 else 
459                   ptTrig = 0x200;// mask for Lpt unlike sign pair
460                 
461
462                 if (esd->GetTrigger() &  ptTrig) NbTrigger++;
463                 
464                 if (invMass > invMassMinInPeak && invMass < invMassMaxInPeak) {
465                   EventInMass++;
466                   hRapResonance->Fill(fVtot.Rapidity());
467                   hPtResonance->Fill(fVtot.Pt());
468
469                   // match with trigger
470                   if (muonTrack->GetMatchTrigger() && (esd->GetTrigger() & ptTrig))  EventInMassMatch++;
471
472                 }
473                 
474               } // if (fCharge1 * fCharge2) == -1)
475             } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
476           } // if (track2MotherId == track1MotherId && track2MotherPDGId == ResType)
477         } //  for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
478       } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
479     } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
480     
481     hNumberOfTrack->Fill(nTracks);
482     //    esdFile->Delete();
483   
484   } // End of event loop
485
486
487   // Loop over events for bg event
488
489   Double_t thetaPlus,  phiPlus;
490   Double_t thetaMinus, phiMinus;
491   Float_t PtMinus, PtPlus;
492   
493   for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {  // Loop over events for bg event
494     // according to Christian a 3d phi-theta-pt random pick  would take better care 
495     // of all correlations 
496
497     hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
498     hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
499     PtPlus = hPtMuonPlus->GetRandom();
500     PtMinus = hPtMuonMinus->GetRandom();
501
502     fPxRec1  = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
503     fPyRec1  = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
504     fPzRec1  = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
505
506     fE1 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
507     fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
508
509     fPxRec2  = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
510     fPyRec2  = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
511     fPzRec2  = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
512
513     fE2 = TMath::Sqrt(MUON_MASS * MUON_MASS + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
514     fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
515
516     // invariant mass
517     fVtot = fV1 + fV2;
518       
519     // fill histos hInvMassAll and hInvMassRes
520     hInvMassBg->Fill(fVtot.M());
521     hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
522
523     // Ntuple for background... more convenient
524     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);
525
526   } // End loop over events for background
527
528
529   // File for histograms and histogram booking
530   TString outfilename = "MUONefficiency.root";
531   TFile *ntupleFile = new TFile(outfilename.Data(), "RECREATE");  
532   
533   Ktuple->Write();
534   ESDtuple->Write();
535   ESDtupleBck->Write();
536
537   ntupleFile->Close();
538   
539   TFile *histoFile = new TFile("MUONhistos.root", "RECREATE");  
540   hPtMuon->Write();
541   hPtMuonPlus->Write();
542   hPtMuonMinus->Write();
543   hPMuon->Write();
544   hChi2PerDof->Write();
545   hInvMassAll->Write();
546   hInvMassBg->Write();
547   hInvMassAll_vs_Pt ->Write();
548   hInvMassBgk_vs_Pt->Write();
549   hInvMassRes->Write();
550   hNumberOfTrack->Write();
551   hRapMuon ->Write();
552   hRapResonance ->Write();
553   hPtResonance ->Write();
554   hThetaPhiPlus ->Write();
555   hThetaPhiMinus ->Write();
556   histoFile->Close();
557
558   cout << "" << endl ;
559   cout << "*************************************************" << endl;
560  
561   cout << "MUONefficiency : " << nprocessedevents  << " events processed" << endl;
562   if (ResType==443)
563     cout << "Number of generated J/Psi (443)  : " <<  numberOfGeneratedResonances  << endl ;
564   if (ResType==553)
565     cout << "Number of generated Upsilon (553)  :" <<  numberOfGeneratedResonances  << endl ;
566   cout << "Chi2Cut for muon tracks = " << Chi2Cut << endl;
567   cout << "PtCutMin for muon tracks = " << PtCutMin << endl;
568   cout << "PtCutMax for muon tracks = " << PtCutMax << endl;
569   cout << "Entries (unlike sign dimuons) in the mass range  ["<<invMassMinInPeak<<";"<<invMassMaxInPeak<<"] : " << EventInMass <<endl;
570   if (ptTrig==0x800) cout << "Unlike Pair - All Pt" ;   
571   if (ptTrig==0x400) cout << "Unlike Pair - High Pt" ;   
572   if (ptTrig==0x200) cout << "Unlike Pair - Low Pt" ; 
573   cout << " triggers : " << NbTrigger << endl;
574   
575   cout << "Entries in the mass range with matching between reconstructed tracks and trigger tracks " << EventInMassMatch << endl;
576
577
578   return kTRUE;
579 }
580