]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/MUONmassPlot_ESD.C
Can no longer use $ALICE_ROOT/OCDB so switch to raw://
[u/mrichter/AliRoot.git] / MUON / MUONmassPlot_ESD.C
CommitLineData
5473f16a 1#if !defined(__CINT__) || defined(__MAKECINT__)
2// ROOT includes
48d6b312 3#include "TTree.h"
5473f16a 4#include "TBranch.h"
5#include "TClonesArray.h"
6#include "TLorentzVector.h"
7#include "TFile.h"
8#include "TH1.h"
ac3c5325 9#include "TH2.h"
5473f16a 10#include "TParticle.h"
11#include "TTree.h"
12#include <Riostream.h>
8cde4af5 13#include <TGeoManager.h>
21939432 14#include <TROOT.h>
5473f16a 15
16// STEER includes
0ff94351 17#include "AliLog.h"
a99c3449 18#include "AliCDBManager.h"
1c2bdf00 19#include "AliESDEvent.h"
995a61aa 20#include "AliESDVertex.h"
103e6575 21#include "AliESDMuonTrack.h"
5473f16a 22
23// MUON includes
a99c3449 24#include "AliMUONCDB.h"
211c52eb 25#include "AliMUONTrackParam.h"
37827b29 26#include "AliMUONTrackExtrap.h"
103e6575 27#include "AliMUONESDInterface.h"
05b4bc73 28#include "AliMUONConstants.h"
5473f16a 29#endif
5473f16a 30
e54bf126 31/// \ingroup macros
32/// \file MUONmassPlot_ESD.C
a99c3449 33/// \brief Macro MUONmassPlot_ESD.C for ESD
e54bf126 34///
35/// \author Ch. Finck, Subatech, April. 2004
36///
37///
38/// Macro to make invariant mass plots
39/// for combinations of 2 muons with opposite charges,
40/// from root file "MUON.tracks.root" containing the result of track reconstruction.
41/// Histograms are stored on the "MUONmassPlot.root" file.
42/// introducing TLorentzVector for parameter calculations (Pt, P,rap,etc...)
43/// using Invariant Mass for rapidity.
44///
45/// Add parameters and histograms for analysis
5473f16a 46
a99c3449 47Bool_t MUONmassPlot(const char* esdFileName = "AliESDs.root", const char* geoFilename = "geometry.root",
48 const char* ocdbPath = "local://$ALICE_ROOT/OCDB",
05b4bc73 49 Int_t FirstEvent = 0, Int_t LastEvent = -1, Int_t ExtrapToVertex = -1,
a99c3449 50 Int_t ResType = 553, Float_t Chi2Cut = 100., Float_t PtCutMin = 1.,
51 Float_t PtCutMax = 10000., Float_t massMin = 9.17,Float_t massMax = 9.77)
5473f16a 52{
a99c3449 53 /// \param FirstEvent (default 0)
54 /// \param LastEvent (default 10000)
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 ResType 553 for Upsilon, anything else for J/Psi (default 553)
60 /// \param Chi2Cut to keep only tracks with chi2 per d.o.f. < Chi2Cut (default 100)
61 /// \param PtCutMin to keep only tracks with transverse momentum > PtCutMin (default 1)
62 /// \param PtCutMax to keep only tracks with transverse momentum < PtCutMax (default 10000)
63 /// \param massMin (default 9.17 for Upsilon)
64 /// \param massMax (default 9.77 for Upsilon);
65 /// to calculate the reconstruction efficiency for resonances with invariant mass
66 /// massMin < mass < massMax.
e54bf126 67
5473f16a 68 cout << "MUONmassPlot " << endl;
69 cout << "FirstEvent " << FirstEvent << endl;
05b4bc73 70 cout << "LastEvent " << ((LastEvent>=0) ? Form("%d",LastEvent) : "all") << endl;
5473f16a 71 cout << "ResType " << ResType << endl;
72 cout << "Chi2Cut " << Chi2Cut << endl;
73 cout << "PtCutMin " << PtCutMin << endl;
74 cout << "PtCutMax " << PtCutMax << endl;
75 cout << "massMin " << massMin << endl;
76 cout << "massMax " << massMax << endl;
77
78
79 //Reset ROOT and connect tree file
80 gROOT->Reset();
81
5473f16a 82 // File for histograms and histogram booking
83 TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
84 TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
ac3c5325 85 TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", 100, 0., 20.);
86 TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", 100, 0., 20.);
5473f16a 87 TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
88 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
89 TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
ac3c5325 90 TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
f57d136a 91 TH2F *hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
92 TH2F *hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
93 TH1F *hInvMassRes;
6ab68e5f 94 TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
5473f16a 95
96 if (ResType == 553) {
97 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
98 } else {
99 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
100 }
101
102 TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
103 TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
104 TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
105 TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
ac3c5325 106 TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
107 TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
5473f16a 108
109
110 // settings
111 Int_t EventInMass = 0;
6678cd54 112 Int_t EventInMassMatch = 0;
113 Int_t NbTrigger = 0;
114
5473f16a 115 Float_t muonMass = 0.105658389;
116// Float_t UpsilonMass = 9.46037;
117// Float_t JPsiMass = 3.097;
118
22ccc301 119 Int_t fCharge1, fCharge2;
5473f16a 120 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
121 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
5473f16a 122
a99c3449 123 Int_t ntrackhits;
5473f16a 124 Double_t fitfmin;
e3c6ae4c 125 Double_t fZVertex=0;
211c52eb 126 Double_t fYVertex=0;
127 Double_t fXVertex=0;
690d2205 128 Double_t errXVtx=0;
129 Double_t errYVtx=0;
5473f16a 130
131 TLorentzVector fV1, fV2, fVtot;
a99c3449 132
8cde4af5 133 // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
134 if (!gGeoManager) {
135 TGeoManager::Import(geoFilename);
136 if (!gGeoManager) {
cec5541f 137 Error("MUONmass_ESD", "getting geometry from file %s failed", geoFilename);
8cde4af5 138 return kFALSE;
139 }
140 }
141
5473f16a 142 // open the ESD file
143 TFile* esdFile = TFile::Open(esdFileName);
144 if (!esdFile || !esdFile->IsOpen()) {
145 Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
146 return kFALSE;
147 }
1c2bdf00 148 AliESDEvent* esd = new AliESDEvent();
48d6b312 149 TTree* tree = (TTree*) esdFile->Get("esdTree");
150 if (!tree) {
a99c3449 151 Error("MUONmass_ESD", "no ESD tree found");
48d6b312 152 return kFALSE;
153 }
8b4a9b89 154 esd->ReadFromTree(tree);
f9ebb3bd 155
a99c3449 156 // get run number
157 if (tree->GetEvent(0) <= 0) {
158 Error("MUONmass_ESD", "no ESD object found for event 0");
159 return kFALSE;
160 }
161 Int_t runNumber = esd->GetRunNumber();
162
163 // load necessary data from OCDB
164 AliCDBManager::Instance()->SetDefaultStorage(ocdbPath);
05b4bc73 165 AliCDBManager::Instance()->SetSpecificStorage("GRP/GRP/Data","local://.");
a99c3449 166 AliCDBManager::Instance()->SetRun(runNumber);
167 if (!AliMUONCDB::LoadField()) return kFALSE;
168
169 // set the magnetic field for track extrapolations
170 AliMUONTrackExtrap::SetField();
211c52eb 171
172 AliMUONTrackParam trackParam;
05b4bc73 173 AliMUONTrackParam trackParamAtAbsEnd;
174
5473f16a 175 // Loop over events
05b4bc73 176 Int_t nevents = (LastEvent >= 0) ? TMath::Min(LastEvent, (Int_t)tree->GetEntries()-1) : (Int_t)tree->GetEntries()-1;
177 for (Int_t iEvent = FirstEvent; iEvent <= nevents; iEvent++) {
5473f16a 178
5473f16a 179 // get the event summary data
a99c3449 180 if (tree->GetEvent(iEvent) <= 0) {
181 Error("MUONmass_ESD", "no ESD object found for event %d", iEvent);
5473f16a 182 return kFALSE;
183 }
a99c3449 184
f9ebb3bd 185 // get the SPD reconstructed vertex (vertexer) and fill the histogram
43939bd8 186 AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
43939bd8 187 if (Vertex->GetNContributors()) {
e690d4d0 188 fZVertex = Vertex->GetZ();
189 fYVertex = Vertex->GetY();
190 fXVertex = Vertex->GetX();
690d2205 191 errXVtx = Vertex->GetXRes();
192 errYVtx = Vertex->GetYRes();
211c52eb 193 }
f9ebb3bd 194 hPrimaryVertex->Fill(fZVertex);
195
48d6b312 196 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
5473f16a 197
198 // printf("\n Nb of events analysed: %d\r",iEvent);
48d6b312 199 // cout << " number of tracks: " << nTracks <<endl;
5473f16a 200
201 // loop over all reconstructed tracks (also first track of combination)
202 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
203
b1fea02e 204 // skip ghosts
205 if (!esd->GetMuonTrack(iTrack)->ContainTrackerData()) continue;
206
8cde4af5 207 AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
5473f16a 208
8cde4af5 209 // extrapolate to vertex if required and available
210 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
103e6575 211 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam);
690d2205 212 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
103e6575 213 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack
8cde4af5 214 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
103e6575 215 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam);
690d2205 216 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
103e6575 217 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack
211c52eb 218 }
8cde4af5 219
05b4bc73 220 // compute track position at the end of the absorber
221 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParamAtAbsEnd);
222 AliMUONTrackExtrap::ExtrapToZ(&trackParamAtAbsEnd, AliMUONConstants::AbsZEnd());
223 Double_t xAbs = trackParamAtAbsEnd.GetNonBendingCoor();
224 Double_t yAbs = trackParamAtAbsEnd.GetBendingCoor();
225 Double_t dAbs1 = TMath::Sqrt(xAbs*xAbs + yAbs*yAbs);
226 Double_t aAbs1 = TMath::ATan(-dAbs1/AliMUONConstants::AbsZEnd()) * TMath::RadToDeg();
227
22ccc301 228 fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
229
230 muonTrack->LorentzP(fV1);
8cde4af5 231
5473f16a 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",
22ccc301 247// fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge1);
5473f16a 248
249 // condition for good track (Chi2Cut and PtCut)
250
a99c3449 251// if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
5473f16a 252
253 // fill histos hPtMuon and hChi2PerDof
254 hPtMuon->Fill(pt1);
255 hPMuon->Fill(p1);
256 hChi2PerDof->Fill(ch1);
257 hRapMuon->Fill(rapMuon1);
22ccc301 258 if (fCharge1 > 0) {
ac3c5325 259 hPtMuonPlus->Fill(pt1);
22ccc301 260 hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
ac3c5325 261 } else {
262 hPtMuonMinus->Fill(pt1);
22ccc301 263 hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
ac3c5325 264 }
5473f16a 265 // loop over second track of combination
266 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
267
b1fea02e 268 // skip ghosts
269 if (!esd->GetMuonTrack(iTrack2)->ContainTrackerData()) continue;
270
8cde4af5 271 AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
272
273 // extrapolate to vertex if required and available
274 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
103e6575 275 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam);
690d2205 276 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx);
103e6575 277 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
8cde4af5 278 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
103e6575 279 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam);
690d2205 280 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.);
103e6575 281 AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
211c52eb 282 }
8cde4af5 283
05b4bc73 284 // compute track position at the end of the absorber
285 AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParamAtAbsEnd);
286 AliMUONTrackExtrap::ExtrapToZ(&trackParamAtAbsEnd, AliMUONConstants::AbsZEnd());
287 xAbs = trackParamAtAbsEnd.GetNonBendingCoor();
288 yAbs = trackParamAtAbsEnd.GetBendingCoor();
289 Double_t dAbs2 = TMath::Sqrt(xAbs*xAbs + yAbs*yAbs);
290 Double_t aAbs2 = TMath::ATan(-dAbs2/AliMUONConstants::AbsZEnd()) * TMath::RadToDeg();
291
8cde4af5 292 fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
5473f16a 293
22ccc301 294 muonTrack2->LorentzP(fV2);
5473f16a 295
8cde4af5 296 ntrackhits = muonTrack2->GetNHit();
297 fitfmin = muonTrack2->GetChi2();
5473f16a 298
299 // transverse momentum
300 Float_t pt2 = fV2.Pt();
301
302 // chi2 per d.o.f.
303 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
304
305 // condition for good track (Chi2Cut and PtCut)
a99c3449 306// if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
05b4bc73 307// if (aAbs1 < 2. || aAbs2 < 2.) {
308// if (aAbs1 > 2. && aAbs2 > 2. && (aAbs1 < 2.1 || aAbs2 < 2.1)) {
309// if (aAbs1 > 2. && aAbs2 > 2. && (dAbs1 < 17.8 || dAbs2 < 17.8)) {
310// if (dAbs1 > 17.8 && dAbs2 > 17.8 && (dAbs1 < 18. || dAbs2 < 18.)) {
311// if (dAbs1 > 18. && dAbs2 > 18. && (dAbs1 < 18.2 || dAbs2 < 18.2)) {
312// if (dAbs1 > 18.2 && dAbs2 > 18.2 && (aAbs1 < 2.1 || aAbs2 < 2.1)) {
313// if (dAbs1 > 18. && dAbs2 > 18.) {
314// if (aAbs1 > 2.1 && aAbs2 > 2.1) {
315// if (muonTrack->GetMatchTrigger() && muonTrack2->GetMatchTrigger()) {
5473f16a 316
317 // condition for opposite charges
22ccc301 318 if ((fCharge1 * fCharge2) == -1) {
5473f16a 319
320 // invariant mass
321 fVtot = fV1 + fV2;
322 Float_t invMass = fVtot.M();
323
324 // fill histos hInvMassAll and hInvMassRes
325 hInvMassAll->Fill(invMass);
326 hInvMassRes->Fill(invMass);
ac3c5325 327 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
6678cd54 328 Int_t ptTrig;
329 if (ResType == 553)
8d4fefab 330 ptTrig = 0x20;// mask for Hpt unlike sign pair
6678cd54 331 else
8d4fefab 332 ptTrig = 0x10;// mask for Lpt unlike sign pair
6678cd54 333
f57d136a 334 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
5473f16a 335 if (invMass > massMin && invMass < massMax) {
336 EventInMass++;
8cde4af5 337 if (muonTrack2->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger
6678cd54 338 EventInMassMatch++;
339
5473f16a 340 hRapResonance->Fill(fVtot.Rapidity());
341 hPtResonance->Fill(fVtot.Pt());
342 }
343
22ccc301 344 } // if (fCharge1 * fCharge2) == -1)
a99c3449 345// } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
8cde4af5 346 delete muonTrack2;
5473f16a 347 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
a99c3449 348// } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
8cde4af5 349 delete muonTrack;
5473f16a 350 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
351
352 hNumberOfTrack->Fill(nTracks);
48d6b312 353 // esdFile->Delete();
5473f16a 354 } // for (Int_t iEvent = FirstEvent;
355
ac3c5325 356// Loop over events for bg event
357
358 Double_t thetaPlus, phiPlus;
359 Double_t thetaMinus, phiMinus;
360 Float_t PtMinus, PtPlus;
361
362 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
363
364 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
365 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
366 PtPlus = hPtMuonPlus->GetRandom();
367 PtMinus = hPtMuonMinus->GetRandom();
368
369 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
370 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
371 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
372
373 fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
374 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
375
376 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
377 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
378 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
379
380 fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
381 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
382
383 // invariant mass
384 fVtot = fV1 + fV2;
385
386 // fill histos hInvMassAll and hInvMassRes
387 hInvMassBg->Fill(fVtot.M());
388 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
389 }
390
5473f16a 391 histoFile->Write();
392 histoFile->Close();
393
6678cd54 394 cout << endl;
5473f16a 395 cout << "EventInMass " << EventInMass << endl;
6678cd54 396 cout << "NbTrigger " << NbTrigger << endl;
397 cout << "EventInMass match with trigger " << EventInMassMatch << endl;
5473f16a 398
399 return kTRUE;
400}
401