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