1 #if !defined(__CINT__) || defined(__MAKECINT__)
5 #include "TClonesArray.h"
6 #include "TLorentzVector.h"
10 #include "TParticle.h"
12 #include <Riostream.h>
17 #include "AliRunLoader.h"
18 #include "AliHeader.h"
19 #include "AliLoader.h"
21 #include "AliMagFMaps.h"
23 #include "AliTracker.h"
26 #include "AliMUONTrackParam.h"
27 #include "AliMUONTrackExtrap.h"
28 #include "AliESDMuonTrack.h"
31 // Macro MUONmassPlot.C for ESD
32 // Ch. Finck, Subatech, April. 2004
35 // macro to make invariant mass plots
36 // for combinations of 2 muons with opposite charges,
37 // from root file "MUON.tracks.root" containing the result of track reconstruction.
38 // Histograms are stored on the "MUONmassPlot.root" file.
39 // introducing TLorentzVector for parameter calculations (Pt, P,rap,etc...)
40 // using Invariant Mass for rapidity.
43 // FirstEvent (default 0)
44 // LastEvent (default 0)
45 // ResType (default 553)
46 // 553 for Upsilon, anything else for J/Psi
47 // Chi2Cut (default 100)
48 // to keep only tracks with chi2 per d.o.f. < Chi2Cut
49 // PtCutMin (default 1)
50 // to keep only tracks with transverse momentum > PtCutMin
51 // PtCutMax (default 10000)
52 // to keep only tracks with transverse momentum < PtCutMax
53 // massMin (default 9.17 for Upsilon)
54 // & massMax (default 9.77 for Upsilon)
55 // to calculate the reconstruction efficiency for resonances with invariant mass
56 // massMin < mass < massMax.
58 // Add parameters and histograms for analysis
60 Bool_t MUONmassPlot(char* filename = "galice.root", Int_t FirstEvent = 0, Int_t LastEvent = 10000,
61 char* esdFileName = "AliESDs.root", Int_t ResType = 553,
62 Float_t Chi2Cut = 100., Float_t PtCutMin = 1., Float_t PtCutMax = 10000.,
63 Float_t massMin = 9.17,Float_t massMax = 9.77)
65 cout << "MUONmassPlot " << endl;
66 cout << "FirstEvent " << FirstEvent << endl;
67 cout << "LastEvent " << LastEvent << endl;
68 cout << "ResType " << ResType << endl;
69 cout << "Chi2Cut " << Chi2Cut << endl;
70 cout << "PtCutMin " << PtCutMin << endl;
71 cout << "PtCutMax " << PtCutMax << endl;
72 cout << "massMin " << massMin << endl;
73 cout << "massMax " << massMax << endl;
76 //Reset ROOT and connect tree file
79 // File for histograms and histogram booking
80 TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
81 TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
82 TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", 100, 0., 20.);
83 TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", 100, 0., 20.);
84 TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
85 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
86 TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
87 TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
88 TH2F *hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
89 TH2F *hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
91 TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
94 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
96 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
99 TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
100 TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
101 TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
102 TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
103 TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
104 TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
108 Int_t EventInMass = 0;
109 Int_t EventInMassMatch = 0;
112 Float_t muonMass = 0.105658389;
113 // Float_t UpsilonMass = 9.46037;
114 // Float_t JPsiMass = 3.097;
116 Double_t thetaX, thetaY, pYZ;
117 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
118 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
119 Int_t fCharge, fCharge2;
121 Int_t ntrackhits, nevents;
127 TLorentzVector fV1, fV2, fVtot;
130 // waiting for mag field in CDB
131 printf("Loading field map...\n");
132 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
133 AliTracker::SetFieldMap(field, kFALSE);
135 // open run loader and load gAlice, kinematics and header
136 AliRunLoader* runLoader = AliRunLoader::Open(filename);
138 Error("MUONmass_ESD", "getting run loader from file %s failed",
144 Error("MUONmass_ESD", "no galice object found");
150 TFile* esdFile = TFile::Open(esdFileName);
151 if (!esdFile || !esdFile->IsOpen()) {
152 Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
156 AliESD* esd = new AliESD();
157 TTree* tree = (TTree*) esdFile->Get("esdTree");
159 Error("CheckESD", "no ESD tree found");
162 tree->SetBranchAddress("ESD", &esd);
166 runLoader->LoadHeader();
167 nevents = runLoader->GetNumberOfEvents();
169 AliMUONTrackParam trackParam;
172 for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) {
175 runLoader->GetEvent(iEvent);
177 // get the event summary data
178 tree->GetEvent(iEvent);
180 Error("CheckESD", "no ESD object found for event %d", iEvent);
184 // get the SPD reconstructed vertex (vertexer) and fill the histogram
185 AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex();
187 if (Vertex->GetNContributors()) {
188 fZVertex = Vertex->GetZv();
189 fYVertex = Vertex->GetYv();
190 fXVertex = Vertex->GetXv();
193 hPrimaryVertex->Fill(fZVertex);
195 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
197 // printf("\n Nb of events analysed: %d\r",iEvent);
198 // cout << " number of tracks: " << nTracks <<endl;
200 // set the magnetic field for track extrapolations
201 AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
202 // loop over all reconstructed tracks (also first track of combination)
203 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
205 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
207 if (!Vertex->GetNContributors()) {
208 //re-extrapolate to vertex, if not kown before.
209 trackParam.GetParamFrom(*muonTrack);
210 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
211 trackParam.SetParamFor(*muonTrack);
213 thetaX = muonTrack->GetThetaX();
214 thetaY = muonTrack->GetThetaY();
216 pYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
217 fPzRec1 = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
218 fPxRec1 = fPzRec1 * TMath::Tan(thetaX);
219 fPyRec1 = fPzRec1 * TMath::Tan(thetaY);
220 fCharge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
222 fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
223 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
225 ntrackhits = muonTrack->GetNHit();
226 fitfmin = muonTrack->GetChi2();
228 // transverse momentum
229 Float_t pt1 = fV1.Pt();
232 Float_t p1 = fV1.P();
235 Float_t rapMuon1 = fV1.Rapidity();
238 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
239 // printf(" px %f py %f pz %f NHits %d Norm.chi2 %f charge %d\n",
240 // fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge);
242 // condition for good track (Chi2Cut and PtCut)
244 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
246 // fill histos hPtMuon and hChi2PerDof
249 hChi2PerDof->Fill(ch1);
250 hRapMuon->Fill(rapMuon1);
252 hPtMuonPlus->Fill(pt1);
253 hThetaPhiPlus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
255 hPtMuonMinus->Fill(pt1);
256 hThetaPhiMinus->Fill(TMath::ATan2(fPyRec1,fPxRec1)*180./TMath::Pi(),TMath::ATan2(pt1,fPzRec1)*180./3.1415);
258 // loop over second track of combination
259 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
261 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack2);
263 if (!Vertex->GetNContributors()) {
264 trackParam.GetParamFrom(*muonTrack);
265 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
266 trackParam.SetParamFor(*muonTrack);
269 thetaX = muonTrack->GetThetaX();
270 thetaY = muonTrack->GetThetaY();
272 pYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
273 fPzRec2 = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
274 fPxRec2 = fPzRec2 * TMath::Tan(thetaX);
275 fPyRec2 = fPzRec2 * TMath::Tan(thetaY);
276 fCharge2 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
278 fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
279 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
281 ntrackhits = muonTrack->GetNHit();
282 fitfmin = muonTrack->GetChi2();
284 // transverse momentum
285 Float_t pt2 = fV2.Pt();
288 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
290 // condition for good track (Chi2Cut and PtCut)
291 if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
293 // condition for opposite charges
294 if ((fCharge * fCharge2) == -1) {
298 Float_t invMass = fVtot.M();
300 // fill histos hInvMassAll and hInvMassRes
301 hInvMassAll->Fill(invMass);
302 hInvMassRes->Fill(invMass);
303 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
306 ptTrig = 0x20;// mask for Hpt unlike sign pair
308 ptTrig = 0x10;// mask for Lpt unlike sign pair
310 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
311 if (invMass > massMin && invMass < massMax) {
313 if (muonTrack->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger
316 hRapResonance->Fill(fVtot.Rapidity());
317 hPtResonance->Fill(fVtot.Pt());
320 } // if (fCharge * fCharge2) == -1)
321 } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
322 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
323 } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
324 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
326 hNumberOfTrack->Fill(nTracks);
327 // esdFile->Delete();
328 } // for (Int_t iEvent = FirstEvent;
330 // Loop over events for bg event
332 Double_t thetaPlus, phiPlus;
333 Double_t thetaMinus, phiMinus;
334 Float_t PtMinus, PtPlus;
336 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
338 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
339 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
340 PtPlus = hPtMuonPlus->GetRandom();
341 PtMinus = hPtMuonMinus->GetRandom();
343 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
344 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
345 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
347 fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
348 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
350 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
351 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
352 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
354 fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
355 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
360 // fill histos hInvMassAll and hInvMassRes
361 hInvMassBg->Fill(fVtot.M());
362 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
369 cout << "EventInMass " << EventInMass << endl;
370 cout << "NbTrigger " << NbTrigger << endl;
371 cout << "EventInMass match with trigger " << EventInMassMatch << endl;