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