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