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