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