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