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