]>
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> |
21939432 | 14 | #include <TROOT.h> |
5473f16a | 15 | |
16 | // STEER includes | |
17 | #include "AliRun.h" | |
0ff94351 | 18 | #include "AliLog.h" |
5473f16a | 19 | #include "AliRunLoader.h" |
20 | #include "AliHeader.h" | |
21 | #include "AliLoader.h" | |
22 | #include "AliStack.h" | |
0ff94351 | 23 | #include "AliMagFMaps.h" |
1c2bdf00 | 24 | #include "AliESDEvent.h" |
995a61aa | 25 | #include "AliESDVertex.h" |
0ff94351 | 26 | #include "AliTracker.h" |
103e6575 | 27 | #include "AliESDMuonTrack.h" |
5473f16a | 28 | |
29 | // MUON includes | |
211c52eb | 30 | #include "AliMUONTrackParam.h" |
37827b29 | 31 | #include "AliMUONTrackExtrap.h" |
103e6575 | 32 | #include "AliMUONESDInterface.h" |
5473f16a | 33 | #endif |
5473f16a | 34 | |
e54bf126 | 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 | |
5473f16a | 50 | |
690d2205 | 51 | Bool_t MUONmassPlot(char* filename = "generated/galice.root", Int_t ExtrapToVertex = -1, char* geoFilename = "geometry.root", |
8cde4af5 | 52 | Int_t FirstEvent = 0, Int_t LastEvent = 10000, char* esdFileName = "AliESDs.root", Int_t ResType = 553, |
5473f16a | 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 | { | |
e54bf126 | 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 | ||
5473f16a | 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 | ||
5473f16a | 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.); | |
ac3c5325 | 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.); | |
5473f16a | 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.); | |
ac3c5325 | 93 | TH1F *hInvMassBg = new TH1F("hInvMassBg", "Mu+Mu- invariant mass BG(GeV/c2)", 480, 0., 12.); |
f57d136a | 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; | |
6ab68e5f | 97 | TH1F *hPrimaryVertex = new TH1F("hPrimaryVertex","SPD reconstructed Z vertex",150,-15,15); |
5473f16a | 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.); | |
ac3c5325 | 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.); | |
5473f16a | 111 | |
112 | ||
113 | // settings | |
114 | Int_t EventInMass = 0; | |
6678cd54 | 115 | Int_t EventInMassMatch = 0; |
116 | Int_t NbTrigger = 0; | |
117 | ||
5473f16a | 118 | Float_t muonMass = 0.105658389; |
119 | // Float_t UpsilonMass = 9.46037; | |
120 | // Float_t JPsiMass = 3.097; | |
121 | ||
22ccc301 | 122 | Int_t fCharge1, fCharge2; |
5473f16a | 123 | Double_t fPxRec1, fPyRec1, fPzRec1, fE1; |
124 | Double_t fPxRec2, fPyRec2, fPzRec2, fE2; | |
5473f16a | 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; | |
690d2205 | 131 | Double_t errXVtx=0; |
132 | Double_t errYVtx=0; | |
5473f16a | 133 | |
134 | TLorentzVector fV1, fV2, fVtot; | |
cb75342e | 135 | |
8cde4af5 | 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 | ||
0ff94351 | 145 | // set mag field |
146 | // waiting for mag field in CDB | |
147 | printf("Loading field map...\n"); | |
b97b210c | 148 | AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG); |
0ff94351 | 149 | AliTracker::SetFieldMap(field, kFALSE); |
cb75342e | 150 | |
5473f16a | 151 | // open run loader and load gAlice, kinematics and header |
152 | AliRunLoader* runLoader = AliRunLoader::Open(filename); | |
153 | if (!runLoader) { | |
8cde4af5 | 154 | Error("MUONmass_ESD", "getting run loader from file %s failed", filename); |
5473f16a | 155 | return kFALSE; |
156 | } | |
6b092dfc | 157 | /* |
158 | runLoader->LoadgAlice(); | |
5473f16a | 159 | if (!gAlice) { |
160 | Error("MUONmass_ESD", "no galice object found"); | |
161 | return kFALSE; | |
162 | } | |
6b092dfc | 163 | */ |
5473f16a | 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 | ||
1c2bdf00 | 172 | AliESDEvent* esd = new AliESDEvent(); |
48d6b312 | 173 | TTree* tree = (TTree*) esdFile->Get("esdTree"); |
174 | if (!tree) { | |
175 | Error("CheckESD", "no ESD tree found"); | |
176 | return kFALSE; | |
177 | } | |
8b4a9b89 | 178 | // tree->SetBranchAddress("ESD", &esd); |
179 | esd->ReadFromTree(tree); | |
f9ebb3bd | 180 | |
48d6b312 | 181 | |
5473f16a | 182 | runLoader->LoadHeader(); |
183 | nevents = runLoader->GetNumberOfEvents(); | |
211c52eb | 184 | |
185 | AliMUONTrackParam trackParam; | |
0ff94351 | 186 | |
5473f16a | 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); | |
48d6b312 | 192 | |
5473f16a | 193 | // get the event summary data |
48d6b312 | 194 | tree->GetEvent(iEvent); |
5473f16a | 195 | if (!esd) { |
48d6b312 | 196 | Error("CheckESD", "no ESD object found for event %d", iEvent); |
5473f16a | 197 | return kFALSE; |
198 | } | |
199 | ||
f9ebb3bd | 200 | // get the SPD reconstructed vertex (vertexer) and fill the histogram |
43939bd8 | 201 | AliESDVertex* Vertex = (AliESDVertex*) esd->GetVertex(); |
43939bd8 | 202 | if (Vertex->GetNContributors()) { |
211c52eb | 203 | fZVertex = Vertex->GetZv(); |
204 | fYVertex = Vertex->GetYv(); | |
205 | fXVertex = Vertex->GetXv(); | |
690d2205 | 206 | errXVtx = Vertex->GetXRes(); |
207 | errYVtx = Vertex->GetYRes(); | |
211c52eb | 208 | } |
f9ebb3bd | 209 | hPrimaryVertex->Fill(fZVertex); |
210 | ||
48d6b312 | 211 | Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; |
5473f16a | 212 | |
213 | // printf("\n Nb of events analysed: %d\r",iEvent); | |
48d6b312 | 214 | // cout << " number of tracks: " << nTracks <<endl; |
5473f16a | 215 | |
37827b29 | 216 | // set the magnetic field for track extrapolations |
217 | AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap()); | |
5473f16a | 218 | // loop over all reconstructed tracks (also first track of combination) |
219 | for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) { | |
220 | ||
b1fea02e | 221 | // skip ghosts |
222 | if (!esd->GetMuonTrack(iTrack)->ContainTrackerData()) continue; | |
223 | ||
8cde4af5 | 224 | AliESDMuonTrack* muonTrack = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack))); |
5473f16a | 225 | |
8cde4af5 | 226 | // extrapolate to vertex if required and available |
227 | if (ExtrapToVertex > 0 && Vertex->GetNContributors()) { | |
103e6575 | 228 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam); |
690d2205 | 229 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx); |
103e6575 | 230 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack |
8cde4af5 | 231 | } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){ |
103e6575 | 232 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack, trackParam); |
690d2205 | 233 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.); |
103e6575 | 234 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack); // put the new parameters in this copy of AliESDMuonTrack |
211c52eb | 235 | } |
8cde4af5 | 236 | |
22ccc301 | 237 | fCharge1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum())); |
238 | ||
239 | muonTrack->LorentzP(fV1); | |
8cde4af5 | 240 | |
5473f16a | 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", | |
22ccc301 | 256 | // fPxRec1, fPyRec1, fPzRec1, ntrackhits, ch1, fCharge1); |
5473f16a | 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); | |
22ccc301 | 267 | if (fCharge1 > 0) { |
ac3c5325 | 268 | hPtMuonPlus->Fill(pt1); |
22ccc301 | 269 | hThetaPhiPlus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi()); |
ac3c5325 | 270 | } else { |
271 | hPtMuonMinus->Fill(pt1); | |
22ccc301 | 272 | hThetaPhiMinus->Fill(fV1.Phi()*180./TMath::Pi(),fV1.Theta()*180./TMath::Pi()); |
ac3c5325 | 273 | } |
5473f16a | 274 | // loop over second track of combination |
275 | for (Int_t iTrack2 = iTrack + 1; iTrack2 < nTracks; iTrack2++) { | |
276 | ||
b1fea02e | 277 | // skip ghosts |
278 | if (!esd->GetMuonTrack(iTrack2)->ContainTrackerData()) continue; | |
279 | ||
8cde4af5 | 280 | AliESDMuonTrack* muonTrack2 = new AliESDMuonTrack(*(esd->GetMuonTrack(iTrack2))); |
281 | ||
282 | // extrapolate to vertex if required and available | |
283 | if (ExtrapToVertex > 0 && Vertex->GetNContributors()) { | |
103e6575 | 284 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam); |
690d2205 | 285 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, fXVertex, fYVertex, fZVertex, errXVtx, errYVtx); |
103e6575 | 286 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack |
8cde4af5 | 287 | } else if ((ExtrapToVertex > 0 && !Vertex->GetNContributors()) || ExtrapToVertex == 0){ |
103e6575 | 288 | AliMUONESDInterface::GetParamAtFirstCluster(*muonTrack2, trackParam); |
690d2205 | 289 | AliMUONTrackExtrap::ExtrapToVertex(&trackParam, 0., 0., 0., 0., 0.); |
103e6575 | 290 | AliMUONESDInterface::SetParamAtVertex(trackParam, *muonTrack2); // put the new parameters in this copy of AliESDMuonTrack |
211c52eb | 291 | } |
8cde4af5 | 292 | |
8cde4af5 | 293 | fCharge2 = Int_t(TMath::Sign(1.,muonTrack2->GetInverseBendingMomentum())); |
5473f16a | 294 | |
22ccc301 | 295 | muonTrack2->LorentzP(fV2); |
5473f16a | 296 | |
8cde4af5 | 297 | ntrackhits = muonTrack2->GetNHit(); |
298 | fitfmin = muonTrack2->GetChi2(); | |
5473f16a | 299 | |
300 | // transverse momentum | |
301 | Float_t pt2 = fV2.Pt(); | |
302 | ||
303 | // chi2 per d.o.f. | |
304 | Float_t ch2 = fitfmin / (2.0 * ntrackhits - 5); | |
305 | ||
306 | // condition for good track (Chi2Cut and PtCut) | |
307 | if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) { | |
308 | ||
309 | // condition for opposite charges | |
22ccc301 | 310 | if ((fCharge1 * fCharge2) == -1) { |
5473f16a | 311 | |
312 | // invariant mass | |
313 | fVtot = fV1 + fV2; | |
314 | Float_t invMass = fVtot.M(); | |
315 | ||
316 | // fill histos hInvMassAll and hInvMassRes | |
317 | hInvMassAll->Fill(invMass); | |
318 | hInvMassRes->Fill(invMass); | |
ac3c5325 | 319 | hInvMassAll_vs_Pt->Fill(invMass,fVtot.Pt()); |
6678cd54 | 320 | Int_t ptTrig; |
321 | if (ResType == 553) | |
8d4fefab | 322 | ptTrig = 0x20;// mask for Hpt unlike sign pair |
6678cd54 | 323 | else |
8d4fefab | 324 | ptTrig = 0x10;// mask for Lpt unlike sign pair |
6678cd54 | 325 | |
f57d136a | 326 | if (esd->GetTriggerMask() & ptTrig) NbTrigger++; |
5473f16a | 327 | if (invMass > massMin && invMass < massMax) { |
328 | EventInMass++; | |
8cde4af5 | 329 | if (muonTrack2->GetMatchTrigger() && (esd->GetTriggerMask() & ptTrig))// match with trigger |
6678cd54 | 330 | EventInMassMatch++; |
331 | ||
5473f16a | 332 | hRapResonance->Fill(fVtot.Rapidity()); |
333 | hPtResonance->Fill(fVtot.Pt()); | |
334 | } | |
335 | ||
22ccc301 | 336 | } // if (fCharge1 * fCharge2) == -1) |
5473f16a | 337 | } // if ((ch2 < Chi2Cut) && (pt2 > PtCutMin) && (pt2 < PtCutMax)) |
8cde4af5 | 338 | delete muonTrack2; |
5473f16a | 339 | } // for (Int_t iTrack2 = iTrack + 1; iTrack2 < iTrack; iTrack2++) |
340 | } // if (ch1 < Chi2Cut) && (pt1 > PtCutMin)&& (pt1 < PtCutMax) ) | |
8cde4af5 | 341 | delete muonTrack; |
5473f16a | 342 | } // for (Int_t iTrack = 0; iTrack < nrectracks; iTrack++) |
343 | ||
344 | hNumberOfTrack->Fill(nTracks); | |
48d6b312 | 345 | // esdFile->Delete(); |
5473f16a | 346 | } // for (Int_t iEvent = FirstEvent; |
347 | ||
ac3c5325 | 348 | // Loop over events for bg event |
349 | ||
350 | Double_t thetaPlus, phiPlus; | |
351 | Double_t thetaMinus, phiMinus; | |
352 | Float_t PtMinus, PtPlus; | |
353 | ||
354 | for (Int_t iEvent = 0; iEvent < hInvMassAll->Integral(); iEvent++) { | |
355 | ||
356 | hThetaPhiPlus->GetRandom2(phiPlus, thetaPlus); | |
357 | hThetaPhiMinus->GetRandom2(phiMinus,thetaMinus); | |
358 | PtPlus = hPtMuonPlus->GetRandom(); | |
359 | PtMinus = hPtMuonMinus->GetRandom(); | |
360 | ||
361 | fPxRec1 = PtPlus * TMath::Cos(TMath::Pi()/180.*phiPlus); | |
362 | fPyRec1 = PtPlus * TMath::Sin(TMath::Pi()/180.*phiPlus); | |
363 | fPzRec1 = PtPlus / TMath::Tan(TMath::Pi()/180.*thetaPlus); | |
364 | ||
365 | fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1); | |
366 | fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1); | |
367 | ||
368 | fPxRec2 = PtMinus * TMath::Cos(TMath::Pi()/180.*phiMinus); | |
369 | fPyRec2 = PtMinus * TMath::Sin(TMath::Pi()/180.*phiMinus); | |
370 | fPzRec2 = PtMinus / TMath::Tan(TMath::Pi()/180.*thetaMinus); | |
371 | ||
372 | fE2 = TMath::Sqrt(muonMass * muonMass + fPxRec2 * fPxRec2 + fPyRec2 * fPyRec2 + fPzRec2 * fPzRec2); | |
373 | fV2.SetPxPyPzE(fPxRec2, fPyRec2, fPzRec2, fE2); | |
374 | ||
375 | // invariant mass | |
376 | fVtot = fV1 + fV2; | |
377 | ||
378 | // fill histos hInvMassAll and hInvMassRes | |
379 | hInvMassBg->Fill(fVtot.M()); | |
380 | hInvMassBgk_vs_Pt->Fill( fVtot.M(), fVtot.Pt() ); | |
381 | } | |
382 | ||
5473f16a | 383 | histoFile->Write(); |
384 | histoFile->Close(); | |
385 | ||
6678cd54 | 386 | cout << endl; |
5473f16a | 387 | cout << "EventInMass " << EventInMass << endl; |
6678cd54 | 388 | cout << "NbTrigger " << NbTrigger << endl; |
389 | cout << "EventInMass match with trigger " << EventInMassMatch << endl; | |
5473f16a | 390 | |
391 | return kTRUE; | |
392 | } | |
393 |