1 #if !defined(__CINT__) || defined(__MAKECINT__)
5 #include "TClonesArray.h"
6 #include "TLorentzVector.h"
10 #include "TParticle.h"
12 #include <Riostream.h>
13 #include <TGeoManager.h>
19 #include "AliRunLoader.h"
20 #include "AliHeader.h"
21 #include "AliLoader.h"
23 #include "AliMagFMaps.h"
24 #include "AliESDEvent.h"
25 #include "AliESDVertex.h"
26 #include "AliTracker.h"
29 #include "AliMUONTrackParam.h"
30 #include "AliMUONTrackExtrap.h"
31 #include "AliESDMuonTrack.h"
34 // Macro MUONmassPlot.C for ESD
35 // Ch. Finck, Subatech, April. 2004
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.
46 // ExtrapToVertex (default -1)
47 // <0: no extrapolation;
48 // =0: extrapolation to (0,0,0);
49 // >0: extrapolation to ESDVertex if available, else to (0,0,0)
50 // FirstEvent (default 0)
51 // LastEvent (default 0)
52 // ResType (default 553)
53 // 553 for Upsilon, anything else for J/Psi
54 // Chi2Cut (default 100)
55 // to keep only tracks with chi2 per d.o.f. < Chi2Cut
56 // PtCutMin (default 1)
57 // to keep only tracks with transverse momentum > PtCutMin
58 // PtCutMax (default 10000)
59 // to keep only tracks with transverse momentum < PtCutMax
60 // massMin (default 9.17 for Upsilon)
61 // & massMax (default 9.77 for Upsilon)
62 // to calculate the reconstruction efficiency for resonances with invariant mass
63 // massMin < mass < massMax.
65 // Add parameters and histograms for analysis
67 Bool_t MUONmassPlot(char* filename = "galice_sim.root", Int_t ExtrapToVertex = -1, char* geoFilename = "geometry.root",
68 Int_t FirstEvent = 0, Int_t LastEvent = 10000, char* esdFileName = "AliESDs.root", Int_t ResType = 553,
69 Float_t Chi2Cut = 100., Float_t PtCutMin = 1., Float_t PtCutMax = 10000.,
70 Float_t massMin = 9.17,Float_t massMax = 9.77)
72 cout << "MUONmassPlot " << endl;
73 cout << "FirstEvent " << FirstEvent << endl;
74 cout << "LastEvent " << LastEvent << endl;
75 cout << "ResType " << ResType << endl;
76 cout << "Chi2Cut " << Chi2Cut << endl;
77 cout << "PtCutMin " << PtCutMin << endl;
78 cout << "PtCutMax " << PtCutMax << endl;
79 cout << "massMin " << massMin << endl;
80 cout << "massMax " << massMax << endl;
83 //Reset ROOT and connect tree file
86 // File for histograms and histogram booking
87 TFile *histoFile = new TFile("MUONmassPlot.root", "RECREATE");
88 TH1F *hPtMuon = new TH1F("hPtMuon", "Muon Pt (GeV/c)", 100, 0., 20.);
89 TH1F *hPtMuonPlus = new TH1F("hPtMuonPlus", "Muon+ Pt (GeV/c)", 100, 0., 20.);
90 TH1F *hPtMuonMinus = new TH1F("hPtMuonMinus", "Muon- Pt (GeV/c)", 100, 0., 20.);
91 TH1F *hPMuon = new TH1F("hPMuon", "Muon P (GeV/c)", 100, 0., 200.);
92 TH1F *hChi2PerDof = new TH1F("hChi2PerDof", "Muon track chi2/d.o.f.", 100, 0., 20.);
93 TH1F *hInvMassAll = new TH1F("hInvMassAll", "Mu+Mu- invariant mass (GeV/c2)", 480, 0., 12.);
94 TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.);
95 TH2F *hInvMassAll_vs_Pt = new TH2F("hInvMassAll_vs_Pt","hInvMassAll_vs_Pt",480,0.,12.,80,0.,20.);
96 TH2F *hInvMassBgk_vs_Pt = new TH2F("hInvMassBgk_vs_Pt","hInvMassBgk_vs_Pt",480,0.,12.,80,0.,20.);
98 TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15);
100 if (ResType == 553) {
101 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around Upsilon", 60, 8., 11.);
103 hInvMassRes = new TH1F("hInvMassRes", "Mu+Mu- invariant mass (GeV/c2) around J/Psi", 80, 0., 5.);
106 TH1F *hNumberOfTrack = new TH1F("hNumberOfTrack","nb of track /evt ",20,-0.5,19.5);
107 TH1F *hRapMuon = new TH1F("hRapMuon"," Muon Rapidity",50,-4.5,-2);
108 TH1F *hRapResonance = new TH1F("hRapResonance"," Resonance Rapidity",50,-4.5,-2);
109 TH1F *hPtResonance = new TH1F("hPtResonance", "Resonance Pt (GeV/c)", 100, 0., 20.);
110 TH2F *hThetaPhiPlus = new TH2F("hThetaPhiPlus", "Theta vs Phi +", 760, -190., 190., 400, 160., 180.);
111 TH2F *hThetaPhiMinus = new TH2F("hThetaPhiMinus", "Theta vs Phi -", 760, -190., 190., 400, 160., 180.);
115 Int_t EventInMass = 0;
116 Int_t EventInMassMatch = 0;
119 Float_t muonMass = 0.105658389;
120 // Float_t UpsilonMass = 9.46037;
121 // Float_t JPsiMass = 3.097;
123 Int_t fCharge1, fCharge2;
124 Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
125 Double_t fPxRec2, fPyRec2, fPzRec2, fE2;
127 Int_t ntrackhits, nevents;
133 TLorentzVector fV1, fV2, fVtot;
135 // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
137 TGeoManager::Import(geoFilename);
139 Error("MUONmass_ESD", "getting geometry from file %s failed", filename);
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);
150 // open run loader and load gAlice, kinematics and header
151 AliRunLoader* runLoader = AliRunLoader::Open(filename);
153 Error("MUONmass_ESD", "getting run loader from file %s failed", filename);
157 runLoader->LoadgAlice();
159 Error("MUONmass_ESD", "no galice object found");
165 TFile* esdFile = TFile::Open(esdFileName);
166 if (!esdFile || !esdFile->IsOpen()) {
167 Error("MUONmass_ESD", "opening ESD file %s failed", esdFileName);
171 AliESDEvent* esd = new AliESDEvent();
172 TTree* tree = (TTree*) esdFile->Get("esdTree");
174 Error("CheckESD", "no ESD tree found");
177 // tree->SetBranchAddress("ESD", &esd);
178 esd->ReadFromTree(tree);
181 runLoader->LoadHeader();
182 nevents = runLoader->GetNumberOfEvents();
184 AliMUONTrackParam trackParam;
187 for (Int_t iEvent = FirstEvent; iEvent <= TMath::Min(LastEvent, nevents - 1); iEvent++) {
190 runLoader->GetEvent(iEvent);
192 // get the event summary data
193 tree->GetEvent(iEvent);
195 Error("CheckESD", "no ESD object found for event %d", iEvent);
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();
206 hPrimaryVertex->Fill(fZVertex);
208 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
210 // printf("\n Nb of events analysed: %d\r",iEvent);
211 // cout << " number of tracks: " << nTracks <<endl;
213 // set the magnetic field for track extrapolations
214 AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
215 // loop over all reconstructed tracks (also first track of combination)
216 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
218 AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack)));
220 // extrapolate to vertex if required and available
221 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
222 trackParam.GetParamFromUncorrected(*muonTrack);
223 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
224 trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
225 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
226 trackParam.GetParamFromUncorrected(*muonTrack);
227 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0.);
228 trackParam.SetParamFor(*muonTrack); // put the new parameters in this copy of AliESDMuonTrack
231 fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
233 muonTrack->LorentzP(fV1);
235 ntrackhits = muonTrack->GetNHit();
236 fitfmin = muonTrack->GetChi2();
238 // transverse momentum
239 Float_t pt1 = fV1.Pt();
242 Float_t p1 = fV1.P();
245 Float_t rapMuon1 = fV1.Rapidity();
248 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
249 // printf(" px %f py %f pz %f NHits %d Norm.chi2 %f charge %d\n",
250 // fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge1);
252 // condition for good track (Chi2Cut and PtCut)
254 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
256 // fill histos hPtMuon and hChi2PerDof
259 hChi2PerDof->Fill(ch1);
260 hRapMuon->Fill(rapMuon1);
262 hPtMuonPlus->Fill(pt1);
263 hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
265 hPtMuonMinus->Fill(pt1);
266 hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi());
268 // loop over second track of combination
269 for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) {
271 AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2)));
273 // extrapolate to vertex if required and available
274 if (ExtrapToVertex > 0 && Vertex->GetNContributors()) {
275 trackParam.GetParamFromUncorrected(*muonTrack2);
276 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex);
277 trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
278 } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){
279 trackParam.GetParamFromUncorrected(*muonTrack2);
280 AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0.);
281 trackParam.SetParamFor(*muonTrack2); // put the new parameters in this copy of AliESDMuonTrack
284 fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum()));
286 muonTrack2->LorentzP(fV2);
288 ntrackhits = muonTrack2->GetNHit();
289 fitfmin = muonTrack2->GetChi2();
291 // transverse momentum
292 Float_t pt2 = fV2.Pt();
295 Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5);
297 // condition for good track (Chi2Cut and PtCut)
298 if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) {
300 // condition for opposite charges
301 if ((fCharge1 * fCharge2) == -1) {
305 Float_t invMass = fVtot.M();
307 // fill histos hInvMassAll and hInvMassRes
308 hInvMassAll->Fill(invMass);
309 hInvMassRes->Fill(invMass);
310 hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt());
313 ptTrig = 0x20;// mask for Hpt unlike sign pair
315 ptTrig = 0x10;// mask for Lpt unlike sign pair
317 if (esd->GetTriggerMask() & ptTrig) NbTrigger++;
318 if (invMass > massMin && invMass < massMax) {
320 if (muonTrack2->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger
323 hRapResonance->Fill(fVtot.Rapidity());
324 hPtResonance->Fill(fVtot.Pt());
327 } // if (fCharge1 * fCharge2) == -1)
328 } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax))
330 } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++)
331 } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) )
333 } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++)
335 hNumberOfTrack->Fill(nTracks);
336 // esdFile->Delete();
337 } // for (Int_t iEvent = FirstEvent;
339 // Loop over events for bg event
341 Double_t thetaPlus, phiPlus;
342 Double_t thetaMinus, phiMinus;
343 Float_t PtMinus, PtPlus;
345 for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) {
347 hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus);
348 hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus);
349 PtPlus = hPtMuonPlus->GetRandom();
350 PtMinus = hPtMuonMinus->GetRandom();
352 fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus);
353 fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus);
354 fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus);
356 fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
357 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
359 fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus);
360 fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus);
361 fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus);
363 fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2);
364 fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2);
369 // fill histos hInvMassAll and hInvMassRes
370 hInvMassBg->Fill(fVtot.M());
371 hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() );
378 cout << "EventInMass " << EventInMass << endl;
379 cout << "NbTrigger " << NbTrigger << endl;
380 cout << "EventInMass match with trigger " << EventInMassMatch << endl;