]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONmassPlot_ESD.C
Protection against special particle types.
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot_ESD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 // ROOT includes
3 #include "TTree.h"
4 #include "TBranch.h"
5 #include "TClonesArray.h"
6 #include "TLorentzVector.h"
7 #include "TFile.h"
8 #include "TH1.h"
9 #include "TH2.h"
10 #include "TParticle.h"
11 #include "TTree.h"
12 #include <Riostream.h>
13 #include <TGeoManager.h>
14
15 // STEER includes
16 #include "AliRun.h"
17 #include "AliLog.h"
18 #include "AliRunLoader.h"
19 #include "AliHeader.h"
20 #include "AliLoader.h"
21 #include "AliStack.h"
22 #include "AliMagFMaps.h"
23 #include "AliESD.h"
24 #include "AliTracker.h"
25
26 // MUON includes
27 #include "AliMUONTrackParam.h"
28 #include "AliMUONTrackExtrap.h"
29 #include "AliESDMuonTrack.h"
30 #endif
31 //
32 // Macro MUONmassPlot.C for ESD
33 // Ch. Finck, Subatech, April. 2004
34 //
35
36 // macro to make invariant mass plots
37 // for combinations of 2 muons with opposite charges,
38 // from root file "MUON.tracks.root" containing the result of track reconstruction.
39 // Histograms are stored on the "MUONmassPlot.root" file.
40 // introducing TLorentzVector for parameter calculations (Pt, P,rap,etc...)
41 // using Invariant Mass for rapidity.
42
43 // Arguments:
44 //   ExtrapToVertex (default -1)
45 //      <0: no extrapolation;
46 //      =0: extrapolation to (0,0,0);
47 //      >0: extrapolation to ESDVertex if available, else to (0,0,0)
48 //   FirstEvent (default 0)
49 //   LastEvent (default 0)
50 //   ResType (default 553)
51 //      553 for Upsilon, anything else for J/Psi
52 //   Chi2Cut (default 100)
53 //      to keep only tracks with chi2 per d.o.f. < Chi2Cut
54 //   PtCutMin (default 1)
55 //      to keep only tracks with transverse momentum > PtCutMin
56 //   PtCutMax (default 10000)
57 //      to keep only tracks with transverse momentum < PtCutMax
58 //   massMin (default 9.17 for Upsilon) 
59 //      &  massMax (default 9.77 for Upsilon) 
60 //         to calculate the reconstruction efficiency for resonances with invariant mass
61 //         massMin < mass < massMax.
62
63 // Add parameters and histograms for analysis 
64
65 Bool_t MUONmassPlot(Int_t ExtrapToVertex = -1, char* geoFilename = "geometry.root", char* filename = "galice.root",
66                   Int_t FirstEvent = 0, Int_t LastEvent = 10000, char* esdFileName = "AliESDs.root", Int_t ResType = 553, 
67                   Float_t Chi2Cut = 100., Float_t PtCutMin = 1., Float_t PtCutMax = 10000.,
68                   Float_t massMin = 9.17,Float_t massMax = 9.77)
69 {
70   cout << "MUONmassPlot " << endl;
71   cout << "FirstEvent " << FirstEvent << endl;
72   cout << "LastEvent " << LastEvent << endl;
73   cout << "ResType " << ResType << endl;
74   cout << "Chi2Cut " << Chi2Cut << endl;
75   cout << "PtCutMin " << PtCutMin << endl;
76   cout << "PtCutMax " << PtCutMax << endl;
77   cout << "massMin " << massMin << endl;
78   cout << "massMax " << massMax << endl;
79
80  
81   //Reset ROOT and connect tree file
82   gROOT->Reset();
83
84   // File for histograms and histogram booking
85   TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
86   TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
87   TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", 100, 0., 20.);
88   TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", 100, 0., 20.);
89   TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
90   TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
91   TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
92   TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
93   TH2F *hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
94   TH2F *hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
95   TH1F *hInvMassRes;
96   TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
97
98   if (ResType == 553) {
99     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
100   } else {
101     hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
102   }
103
104   TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
105   TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
106   TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
107   TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
108   TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
109   TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
110
111
112   // settings
113   Int_t EventInMass = 0;
114   Int_t EventInMassMatch = 0;
115   Int_t NbTrigger = 0;
116
117   Float_t muonMass = 0.105658389;
118 //   Float_t UpsilonMass = 9.46037;
119 //   Float_t JPsiMass = 3.097;
120
121   Int_t fCharge1, fCharge2;
122   Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
123   Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
124
125   Int_t ntrackhits, nevents;
126   Double_t fitfmin;
127   Double_t fZVertex=0;
128   Double_t fYVertex=0;
129   Double_t fXVertex=0;
130  
131   TLorentzVector fV1, fV2, fVtot;
132
133   // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
134   if (!gGeoManager) {
135     TGeoManager::Import(geoFilename);
136     if (!gGeoManager) {
137       Error("MUONmass_ESD", "getting geometry from file %s failed", filename);
138       return kFALSE;
139     }
140   }
141   
142   // set mag field
143   // waiting for mag field in CDB 
144   printf("Loading field map...\n");
145   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
146   AliTracker::SetFieldMap(field, kFALSE);
147
148   // open run loader and load gAlice, kinematics and header
149   AliRunLoader* runLoader = AliRunLoader::Open(filename);
150   if (!runLoader) {
151     Error("MUONmass_ESD", "getting run loader from file %s failed", filename);
152     return kFALSE;
153   }
154
155   if (!gAlice) {
156     Error("MUONmass_ESD", "no galice object found");
157     return kFALSE;
158   }
159   
160
161   // open the ESD file
162   TFile* esdFile = TFile::Open(esdFileName);
163   if (!esdFile || !esdFile->IsOpen()) {
164     Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
165     return kFALSE;
166   }
167   
168   AliESD* esd = new AliESD();
169   TTree* tree = (TTree*) esdFile->Get("esdTree");
170   if (!tree) {
171     Error("CheckESD", "no ESD tree found");
172     return kFALSE;
173   }
174   tree->SetBranchAddress("ESD", &esd);
175   
176   
177
178   runLoader->LoadHeader();
179   nevents = runLoader->GetNumberOfEvents();
180   
181   AliMUONTrackParam trackParam;
182
183   // Loop over events
184   for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) {
185
186     // get current event
187     runLoader->GetEvent(iEvent);
188    
189     // get the event summary data
190     tree->GetEvent(iEvent);
191     if (!esd) {
192       Error("CheckESD", "no ESD object found for event %d", iEvent);
193       return kFALSE;
194     }
195
196     // get the SPD reconstructed vertex (vertexer) and fill the histogram
197     AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
198     if (Vertex->GetNContributors()) {
199       fZVertex = Vertex->GetZv();
200       fYVertex = Vertex->GetYv();
201       fXVertex = Vertex->GetXv();
202     }
203     hPrimaryVertex->Fill(fZVertex);
204
205     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
206
207     //    printf("\n Nb of events analysed: %d\r",iEvent);
208     //      cout << " number of tracks: " << nTracks  <<endl;
209   
210     // set the magnetic field for track extrapolations
211     AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
212     // loop over all reconstructed tracks (also first track of combination)
213     for (Int_t iTrack = 0; iTrack <  nTracks;  iTrack++) {
214
215       AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
216
217       // extrapolate to vertex if required and available
218       if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
219         trackParam.GetParamFromUncorrected(*muonTrack);
220         AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
221         trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
222       } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
223         trackParam.GetParamFromUncorrected(*muonTrack);
224         AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0.);
225         trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
226       }
227
228       fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
229       
230       muonTrack->LorentzP(fV1);
231       
232       ntrackhits = muonTrack->GetNHit();
233       fitfmin    = muonTrack->GetChi2();
234
235       // transverse momentum
236       Float_t pt1 = fV1.Pt();
237
238       // total momentum
239       Float_t p1 = fV1.P();
240
241       // Rapidity
242       Float_t rapMuon1 = fV1.Rapidity();
243
244       // chi2 per d.o.f.
245       Float_t ch1 =  fitfmin / (2.0 * ntrackhits - 5);
246 //      printf(" px %f py %f pz %f NHits %d  Norm.chi2 %f charge %d\n", 
247 //           fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge1);
248
249       // condition for good track (Chi2Cut and PtCut)
250
251       if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
252
253         // fill histos hPtMuon and hChi2PerDof
254         hPtMuon->Fill(pt1);
255         hPMuon->Fill(p1);
256         hChi2PerDof->Fill(ch1);
257         hRapMuon->Fill(rapMuon1);
258         if (fCharge1 > 0) {
259           hPtMuonPlus->Fill(pt1);
260           hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
261         } else {
262           hPtMuonMinus->Fill(pt1);
263           hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
264         }
265         // loop over second track of combination
266         for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
267           
268           AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
269           
270           // extrapolate to vertex if required and available
271           if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
272             trackParam.GetParamFromUncorrected(*muonTrack2);
273             AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
274             trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
275           } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
276             trackParam.GetParamFromUncorrected(*muonTrack2);
277             AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0.);
278             trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
279           }
280           
281           fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
282
283           muonTrack2->LorentzP(fV2);
284
285           ntrackhits = muonTrack2->GetNHit();
286           fitfmin    = muonTrack2->GetChi2();
287
288           // transverse momentum
289           Float_t pt2 = fV2.Pt();
290
291           // chi2 per d.o.f.
292           Float_t ch2 = fitfmin  / (2.0 * ntrackhits - 5);
293
294           // condition for good track (Chi2Cut and PtCut)
295           if ((ch2 < Chi2Cut) && (pt2 > PtCutMin)  && (pt2 < PtCutMax)) {
296
297             // condition for opposite charges
298             if ((fCharge1 * fCharge2) == -1) {
299
300               // invariant mass
301               fVtot = fV1 + fV2;
302               Float_t invMass = fVtot.M();
303                     
304               // fill histos hInvMassAll and hInvMassRes
305               hInvMassAll->Fill(invMass);
306               hInvMassRes->Fill(invMass);
307               hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
308               Int_t ptTrig;
309               if (ResType == 553) 
310                 ptTrig =  0x20;// mask for Hpt unlike sign pair
311               else 
312                 ptTrig =  0x10;// mask for Lpt unlike sign pair
313
314               if (esd->GetTriggerMask() &  ptTrig) NbTrigger++; 
315               if (invMass > massMin && invMass < massMax) {
316                 EventInMass++;
317                 if (muonTrack2->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger
318                   EventInMassMatch++;
319
320                 hRapResonance->Fill(fVtot.Rapidity());
321                 hPtResonance->Fill(fVtot.Pt());
322               }
323
324             } // if (fCharge1 * fCharge2) == -1)
325           } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
326           delete muonTrack2;
327         } //  for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
328       } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
329       delete muonTrack;
330     } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
331
332     hNumberOfTrack->Fill(nTracks);
333     //    esdFile->Delete();
334   } // for (Int_t iEvent = FirstEvent;
335
336 // Loop over events for bg event
337
338   Double_t thetaPlus,  phiPlus;
339   Double_t thetaMinus, phiMinus;
340   Float_t PtMinus, PtPlus;
341   
342   for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
343
344     hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
345     hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
346     PtPlus = hPtMuonPlus->GetRandom();
347     PtMinus = hPtMuonMinus->GetRandom();
348
349     fPxRec1  = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
350     fPyRec1  = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
351     fPzRec1  = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
352
353     fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
354     fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
355
356     fPxRec2  = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
357     fPyRec2  = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
358     fPzRec2  = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
359
360     fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
361     fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
362
363     // invariant mass
364     fVtot = fV1 + fV2;
365       
366     // fill histos hInvMassAll and hInvMassRes
367     hInvMassBg->Fill(fVtot.M());
368     hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
369   }
370
371   histoFile->Write();
372   histoFile->Close();
373
374   cout << endl;
375   cout << "EventInMass " << EventInMass << endl;
376   cout << "NbTrigger " << NbTrigger << endl;
377   cout << "EventInMass match with trigger " << EventInMassMatch << endl;
378
379   return kTRUE;
380 }
381