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