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