Belle II Software prerelease-10-00-00a
CDCDedxValidation.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <cdc/modules/CDCDedxValidation/CDCDedxValidation.h>
10
11#include <mdst/dataobjects/Track.h>
12#include <mdst/dataobjects/ECLCluster.h>
13
14#include <TH2D.h>
15
16#include <cmath>
17
18using namespace Belle2;
19
20
21REG_MODULE(CDCDedxValidation);
22
23//-------------------------------------------------
26 fD0Window(1.0),
27 fZ0Window(1.0),
28 fnRunCounter(0),
29 fiRun(0),
30 fnBinsdedx(120),
31 fnBinsdedxLE(0.4),
32 fnBinsdedxUE(1.6)
33{
34 setDescription("Make data quality monitoring plots for CDC dE/dx");
35 addParam("outputFileName", fOutFileName, "Name for output file", std::string("CDCdEdxValidation.root"));
36 addParam("SampleType", fCollType, "Switch to hadron (false) vs bhabha files", std::string("temp"));
37 addParam("fnRuns", fnRuns, "Number of input runs");
38
39}
40
41//-----------------------------------------
43{
44
45 m_cdcDedxTracks.isRequired();
46
47 if (fCollType != "radbhabha" and fCollType != "bhabha" and fCollType != "hadron") {
48 printf("Wrong input file type (%s): choose bhabha or radbhabha or hadron \n", fCollType.data());
49 return;
50 } else {
51 fBasic = new TList(); fBasic->SetOwner(); fBasic->SetName("AllRunBasics");
52 fPRdEdx = new TList(); fPRdEdx->SetOwner(); fPRdEdx->SetName("PerRunDEdx");
53 fPRdEdxinP = new TList(); fPRdEdxinP->SetOwner(); fPRdEdxinP->SetName("PerRundEdxinP");
54 DefineHistograms("AR", 0);
55 }
56}
57
58
59//---------------------------------------
61{
62
63 StoreObjPtr<EventMetaData> eventMetaDataPtr;
64 fCurrentRunNum = eventMetaDataPtr->getRun();
65
67 fcRunGain = -1.0;
68
70
72}
73
74
75//------------------------------------
77{
78
79 //Loop over CDC Tracks
80 for (Int_t idedx = 0; idedx < m_cdcDedxTracks.getEntries(); idedx++) {
81
82 CDCDedxTrack* dedxTrack = m_cdcDedxTracks[idedx];
83 if (!dedxTrack) {
84 ((TH1D*)fBasic->FindObject("hTrkPerEvtStats"))->Fill(0.0);
85 continue;
86 }
87
88 const Track* track = dedxTrack->getRelatedFrom<Track>();
89 if (!track) {
90 ((TH1D*)fBasic->FindObject("hTrkPerEvtStats"))->Fill(1.0);
91 continue;
92 }
93
94 const TrackFitResult* mTrack = NULL;
95 if (fCollType == "bhabha" || fCollType == "radbhabha") {
96 mTrack = track->getTrackFitResultWithClosestMass(Const::electron);
97 } else {
98 mTrack = track->getTrackFitResultWithClosestMass(Const::pion);
99 }
100 if (!mTrack) {
101 ((TH1D*)fBasic->FindObject("hTrkPerEvtStats"))->Fill(2.0);
102 continue;
103 }
104
105 bool IsTrkSelected = IsSelectedTrack(mTrack);
106 if (!IsTrkSelected) {
107 ((TH1D*)fBasic->FindObject("hTrkPerEvtStats"))->Fill(3.0);
108 continue;
109 }
110
111 if (dedxTrack->getNLayerHits() <= 20) continue;
112
113 const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
114 if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
116 if (fCollType == "bhabha" || fCollType == "radbhabha") {
117 if (std::abs(fTrkEoverP - 1.0) >= 0.2)continue;
118 ((TH1D*)fBasic->FindObject(Form("hEOverP_AR")))->Fill(double(fTrkEoverP));
119 }
120 } else {
121 ((TH1D*)fBasic->FindObject("hTrkPerEvtStats"))->Fill(4.0);
122 continue;
123
124 }
125
126 ((TH1D*)fBasic->FindObject("hTrkPerEvtStats"))->Fill(5.0);
127 FillHistograms(dedxTrack, mTrack);
128 }
129}
130
131//----------------------------------------------------------------------------------------------------
133{
134
135 fcRunGain = dedxTrack->getRunGain();
136
137 Int_t TrkCharge = mTrack->getChargeSign();
138 Double_t TrkdEdxnosat = dedxTrack->getDedxNoSat();
139 Double_t TrkdEdx = dedxTrack->getDedx();
140 Double_t TrkCosTheta = dedxTrack->getCosTheta();
141 Double_t TrkMom = dedxTrack->getMomentum();
142
143 if (TrkMom >= 8.00)TrkMom = 7.999;
144 Double_t BWMom = 0.250; //in GeV
145 Int_t iMomBin = Int_t(TrkMom / BWMom);
146
147 if (fCollType == "radbhabha" or fCollType == "bhabha") {
148
149 if (TrkCharge > 0) {
150 ((TH1D*)fBasic->FindObject(Form("hP_Positron_AR")))->Fill(TrkMom);
151 ((TH1D*)fBasic->FindObject(Form("hdEdx_Positron_AR")))->Fill(TrkdEdxnosat);
152 ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", iMomBin)))->Fill(TrkdEdxnosat);
153 } else if (TrkCharge < 0) {
154 ((TH1D*)fBasic->FindObject(Form("hP_Electron_AR")))->Fill(TrkMom);
155 ((TH1D*)fBasic->FindObject(Form("hdEdx_Electron_AR")))->Fill(TrkdEdxnosat);
156 ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Elec_Pbin_AR%d", iMomBin)))->Fill(TrkdEdxnosat);
157 }
158
159 ((TH1D*)fBasic->FindObject(Form("hdEdx_AR")))->Fill(double(TrkdEdxnosat));
160 ((TH2D*)fBasic->FindObject(Form("hdEdxvsPhi_AR")))->Fill(double(mTrack->getPhi()), double(TrkdEdxnosat));
161 ((TH2D*)fBasic->FindObject(Form("hPvsdEdx_AR")))->Fill(TrkMom * TrkCharge, double(TrkdEdxnosat));
162 ((TH2D*)fBasic->FindObject(Form("hPvsCosth_AR")))->Fill(TrkMom * TrkCharge, double(TrkCosTheta));
163
164 hdEdx_PR[fiRun]->Fill(double(TrkdEdxnosat));
165
166 } else if (fCollType == "hadron") {
167
168 // double ChiE = dedxTrack->getChi(0);
169 // double ChiMu = dedxTrack->getChi(1);
170 double ChiPi = dedxTrack->getChi(2);
171 double ChiK = dedxTrack->getChi(3);
172 double ChiP = dedxTrack->getChi(4);
173 // double ChiD = dedxTrack->getChi(5);
174
175 ((TH2D*)fBasic->FindObject(Form("hPvsdEdx_hadAR")))->Fill(TrkMom, TrkdEdx);
176
177 if ((TrkMom < 0.40) && (fTrkEoverP < 0.4) && (TrkdEdx < (0.6 + 0.10 / (TrkMom * TrkMom)))
178 && (TrkdEdx > (0.4 + 0.012 / (TrkMom * TrkMom)))) {
179 ((TH1D*)fBasic->FindObject(Form("hPionChiallP")))->Fill(ChiPi);
180 if (TrkMom < 0.300)((TH1D*)fBasic->FindObject(Form("hPionChiLowP")))->Fill(ChiPi);
181 else ((TH1D*)fBasic->FindObject(Form("hPionChiHighP")))->Fill(ChiPi);
182 ((TH2D*)fBasic->FindObject(Form("hPvsdEdxPion_hadAR")))->Fill(TrkMom, TrkdEdx);
183 }
184
185 if ((TrkMom < 0.40) && (TrkdEdx > 1.35) && (TrkdEdx < (0.6 + 0.40 / (TrkMom * TrkMom)))
186 && (TrkdEdx > (0.6 + 0.10 / (TrkMom * TrkMom)))) {
187 ((TH1D*)fBasic->FindObject(Form("hKaonChiallP")))->Fill(ChiK);
188 if (TrkMom < 0.350)((TH1D*)fBasic->FindObject(Form("hKaonChiLowP")))->Fill(ChiK);
189 else ((TH1D*)fBasic->FindObject(Form("hKaonChiHighP")))->Fill(ChiK);
190 ((TH2D*)fBasic->FindObject(Form("hPvsdEdxKaon_hadAR")))->Fill(TrkMom, TrkdEdx);
191 }
192
193 if ((TrkMom < 0.80) && (TrkdEdx > 1.35) && (TrkdEdx < (0.6 + 1.20 / (TrkMom * TrkMom)))
194 && (TrkdEdx > (0.6 + 0.40 / (TrkMom * TrkMom)))) {
195 ((TH1D*)fBasic->FindObject(Form("hProtonChiallP")))->Fill(ChiP);
196 if (TrkMom < 0.600)((TH1D*)fBasic->FindObject(Form("hProtonChiLowP")))->Fill(ChiP);
197 else ((TH1D*)fBasic->FindObject(Form("hProtonChiHighP")))->Fill(ChiP);
198 ((TH2D*)fBasic->FindObject(Form("hPvsdEdxProton_hadAR")))->Fill(TrkMom, TrkdEdx);
199 }
200 }
201
202}
203
204
205
206//----------------------------------
208{
209 if (fCollType != "hadron")ExtractHistograms("PR");
210}
211
212
213//----------------------------------------------------------------------------------
214void CDCDedxValidationModule::DefineHistograms(TString level = "XR", Int_t iR = 123)
215{
216 if (level == "AR") {
217
218 TH1D* hTrkPerEvtStats = new TH1D("hTrkPerEvtStats", "Track selections", 6, -0.5, 5.5);
219 hTrkPerEvtStats->GetXaxis()->SetBinLabel(1, "no dedxobject");
220 hTrkPerEvtStats->GetXaxis()->SetBinLabel(2, "no assoc-track");
221 hTrkPerEvtStats->GetXaxis()->SetBinLabel(3, "no TrackFit");
222 hTrkPerEvtStats->GetXaxis()->SetBinLabel(4, "no pass cuts");
223 hTrkPerEvtStats->GetXaxis()->SetBinLabel(5, "no eclCluster");
224 hTrkPerEvtStats->GetXaxis()->SetBinLabel(6, "Selected");
225 hTrkPerEvtStats->SetFillColor(kRed);
226 hTrkPerEvtStats->SetFillStyle(3015);
227 hTrkPerEvtStats->SetMinimum(0);
228 fBasic->Add(hTrkPerEvtStats);
229
230 if (fCollType == "radbhabha" or fCollType == "bhabha") {
231
232 TH1D* hdEdx_AR = new TH1D("hdEdx_AR", "dE/dx (nohad sat)", fnBinsdedx, fnBinsdedxLE, fnBinsdedxUE);
233 hdEdx_AR->GetXaxis()->SetTitle(Form("dE/dx truncMean of %s tracks", fCollType.data()));
234 hdEdx_AR->GetYaxis()->SetTitle("Entries");
235 fBasic->Add(hdEdx_AR);
236
237 TH1D* hEOverP_AR = new TH1D("hEOverP_AR", "E/p distribution", 100, 0.5, 1.5);
238 hEOverP_AR->GetXaxis()->SetTitle("E/p distribution");
239 hEOverP_AR->GetYaxis()->SetTitle("Entries");
240 fBasic->Add(hEOverP_AR);
241
242 TH1D* hRunGainPR = new TH1D("hRunGainPR", "bla-bla", fnRuns, -0.5, fnRuns - 0.5);
243 hRunGainPR->SetTitle("Run gain variation vs. RunNumber;Run Numbers;dE/dx mean");
244 hRunGainPR->GetYaxis()->SetRangeUser(0.85, 1.15);
245 fBasic->Add(hRunGainPR);
246
247 TH1D* hP_Electron_AR = new TH1D("hP_Electron_AR", "bla-bla", 320, 0.0, 8.0);
248 hP_Electron_AR->SetTitle("Momentum distribution of e-; Momentum of (e-); Entries");
249 fBasic->Add(hP_Electron_AR);
250
251 TH1D* hP_Positron_AR = new TH1D("hP_Positron_AR", "bla-bla", 320, 0.0, 8.0);
252 hP_Positron_AR->SetTitle("Momentum distribution of e+;Momentum of (e+);Entries");
253 fBasic->Add(hP_Positron_AR);
254
255 TH1D* hdEdx_Electron_AR = new TH1D("hdEdx_Electron_AR", "bla-bla", fnBinsdedx, fnBinsdedxLE, fnBinsdedxUE);
256 hdEdx_Electron_AR->SetTitle("dE/dx (nohad sat) of e- ;dE/dx distribution (e-);Entries");
257 fBasic->Add(hdEdx_Electron_AR);
258
259 TH1D* hdEdx_Positron_AR = new TH1D("hdEdx_Positron_AR", "bla-bla", fnBinsdedx, fnBinsdedxLE, fnBinsdedxUE);
260 hdEdx_Positron_AR->SetTitle("dE/dx (nohad sat) of e+;dE/dx distribution (e+);Entries");
261 fBasic->Add(hdEdx_Positron_AR);
262
263 TH2D* hPvsdEdx_AR = new TH2D("hPvsdEdx_AR", "bla-bla", 320, -8.0, 8.0, 100, 0.0, 2.0);
264 hPvsdEdx_AR->SetTitle("dE/dx band plots for e+ and e-; Momentum of (e+(right)) and e-);dE/dx");
265 fBasic->Add(hPvsdEdx_AR);
266
267 TH2D* hdEdxvsPhi_AR = new TH2D("hdEdxvsPhi_AR", "dE/dx (no had sat) vs #phi", 64, -3.14, 3.14, 200, 0., 2.0);
268 hdEdxvsPhi_AR->SetTitle("dE/dx (no Had Sat) vs #phi;track #phi;dE/dx");
269 fBasic->Add(hdEdxvsPhi_AR);
270
271 TH2D* hPvsCosth_AR = new TH2D("hPvsCosth_AR", "cos(#theta) vs. p: all Runs", 2 * 48, -10., 10., 60, -1.2, 1.2);
272 hPvsCosth_AR->GetXaxis()->SetTitle(Form("Momentum of %s tracks", fCollType.data()));
273 hPvsCosth_AR->GetYaxis()->SetTitle("cos(#theta)");
274 fBasic->Add(hPvsCosth_AR);
275
276
277 TH1D* hdEdx_Posi_Pbin_AR[32], *hdEdx_Elec_Pbin_AR[32];
278 for (int ip = 0; ip < 32; ip++) {
279
280 hdEdx_Posi_Pbin_AR[ip] = new TH1D(Form("hdEdx_Posi_Pbin_AR%d", ip), Form("hdEdx_Posi_Pbin_AR%d", ip), fnBinsdedx, fnBinsdedxLE,
282 hdEdx_Posi_Pbin_AR[ip]->GetXaxis()->SetTitle("dE/dx distribution (e+)");
283 hdEdx_Posi_Pbin_AR[ip]->GetYaxis()->SetTitle("Entries");
284 hdEdx_Posi_Pbin_AR[ip]->SetTitle(Form("Momentum range %0.03f to %0.03f", ip * 0.250, (ip + 1) * 0.250));
285 fPRdEdxinP->Add(hdEdx_Posi_Pbin_AR[ip]);
286
287 hdEdx_Elec_Pbin_AR[ip] = new TH1D(Form("hdEdx_Elec_Pbin_AR%d", ip), Form("hdEdx_Elec_Pbin_AR%d", ip), fnBinsdedx, fnBinsdedxLE,
289 hdEdx_Elec_Pbin_AR[ip]->GetXaxis()->SetTitle("dE/dx distribution (e-)");
290 hdEdx_Elec_Pbin_AR[ip]->GetYaxis()->SetTitle("Entries");
291 hdEdx_Elec_Pbin_AR[ip]->SetTitle(Form("Momentum range %0.03f to %0.03f", ip * 0.250, (ip + 1) * 0.250));
292
293 }
294
295 for (int ip = 0; ip < 32; ip++) fPRdEdxinP->Add(hdEdx_Elec_Pbin_AR[ip]); //Adding later for simplicity
296
297 hdEdx_PR.reserve(fnRuns);
298 }
299
300 if (fCollType == "hadron") {
301
302 TH2D* hPvsdEdx_hadAR = new TH2D("hPvsdEdx_hadAR", "bla-bla", 500, 0.10, 15.0, 750, 0.05, 15);
303 hPvsdEdx_hadAR->SetTitle("dE/dx band plot; Momentum;dE/dx");
304 fBasic->Add(hPvsdEdx_hadAR);
305
306 TH2D* hPvsdEdxPion_hadAR = new TH2D("hPvsdEdxPion_hadAR", "bla-bla", 500, 0.10, 15.0, 750, 0.05, 15);
307 hPvsdEdxPion_hadAR->SetTitle("dE/dx band plot (Pion); Momentum;dE/dx");
308 hPvsdEdxPion_hadAR->SetMarkerColor(kRed);
309 fBasic->Add(hPvsdEdxPion_hadAR);
310
311 TH2D* hPvsdEdxKaon_hadAR = new TH2D("hPvsdEdxKaon_hadAR", "bla-bla", 500, 0.10, 15.0, 750, 0.05, 15);
312 hPvsdEdxKaon_hadAR->SetTitle("dE/dx band plot (Kaon); Momentum;dE/dx");
313 hPvsdEdxKaon_hadAR->SetMarkerColor(kGreen);
314 fBasic->Add(hPvsdEdxKaon_hadAR);
315
316 TH2D* hPvsdEdxProton_hadAR = new TH2D("hPvsdEdxProton_hadAR", "bla-bla", 500, 0.10, 15.0, 750, 0.05, 15);
317 hPvsdEdxProton_hadAR->SetTitle("dE/dx band plot (Proton); Momentum;dE/dx");
318 hPvsdEdxKaon_hadAR->SetMarkerColor(kBlue);
319 fBasic->Add(hPvsdEdxProton_hadAR);
320
321 //Pions chi values
322 TH1D* hPionChiallP = new TH1D("hPionChiallP", "bla-bla", 240, -6.0, 6.0);
323 hPionChiallP->SetTitle("Chi value (Pion);chi value; Entries");
324 fBasic->Add(hPionChiallP);
325
326 TH1D* hPionChiLowP = new TH1D("hPionChiLowP", "bla-bla", 240, -6.0, 6.0);
327 hPionChiLowP->SetTitle("Chi value (Pion), Momentum (0-300) MeV; chi value; Entries");
328 fBasic->Add(hPionChiLowP);
329
330 TH1D* hPionChiHighP = new TH1D("hPionChiHighP", "bla-bla", 240, -6.0, 6.0);
331 hPionChiHighP->SetTitle("Chi value (Pion), Momentum (300-400) MeV; chi value; Entries");
332 fBasic->Add(hPionChiHighP);
333
334 //Kaons chi values
335 TH1D* hKaonChiallP = new TH1D("hKaonChiallP", "bla-bla", 240, -6.0, 6.0);
336 hKaonChiallP->SetTitle("Chi value (Kaon);chi value; Entries");
337 fBasic->Add(hKaonChiallP);
338
339 TH1D* hKaonChiLowP = new TH1D("hKaonChiLowP", "bla-bla", 240, -6.0, 6.0);
340 hKaonChiLowP->SetTitle("Chi value (Kaon), Momentum (0-350) MeV; chi value; Entries");
341 fBasic->Add(hKaonChiLowP);
342
343 TH1D* hKaonChiHighP = new TH1D("hKaonChiHighP", "bla-bla", 240, -6.0, 6.0);
344 hKaonChiHighP->SetTitle("Chi value (Kaon), Momentum (350-800) MeV; chi value; Entries");
345 fBasic->Add(hKaonChiHighP);
346
347 //Protons chi values
348 TH1D* hProtonChiallP = new TH1D("hProtonChiallP", "bla-bla", 240, -6.0, 6.0);
349 hProtonChiallP->SetTitle("Chi value (Proton);chi value; Entries");
350 fBasic->Add(hProtonChiallP);
351
352 TH1D* hProtonChiLowP = new TH1D("hProtonChiLowP", "bla-bla", 240, -6.0, 6.0);
353 hProtonChiLowP->SetTitle("Chi value (Proton), Momentum (0-600) MeV; chi value; Entries");
354 fBasic->Add(hProtonChiLowP);
355
356 TH1D* hProtonChiHighP = new TH1D("hProtonChiHighP", "bla-bla", 240, -6.0, 6.0);
357 hProtonChiHighP->SetTitle("Chi value (Proton), Momentum (600-800) MeV; chi value; Entries");
358 fBasic->Add(hProtonChiHighP);
359 }
360 } else if (level == "PR") {
361 if (fCollType != "hadron") {
362 hdEdx_PR[iR] = new TH1D(Form("hdEdx_Run%d", fCurrentRunNum), Form("dE/dx (no had sat): Run # = %d", fCurrentRunNum), fnBinsdedx,
364 hdEdx_PR[iR]->GetXaxis()->SetTitle(Form("dE/dx trucMean of %s tracks", fCollType.data()));
365 hdEdx_PR[iR]->GetYaxis()->SetTitle("Entries");
366 fPRdEdx->Add(hdEdx_PR[iR]);
367 }
368 } else {
369 B2ERROR("Run Gain: Enter AR or PR mode only");
370 }
371}
372
373
374//-----------------------------------------------------------------------
376{
377
378 if (level == "PR") {
379
380 Double_t mean = 0., meanError = 0.;
381 Double_t sigma = 0., sigmaError = 0.;
382
383 if (hdEdx_PR[fiRun]->GetEntries() > 100) {
384 Int_t fitStatus = -1;
385 fitStatus = hdEdx_PR[fiRun]->Fit("gaus", "Q"); //Q = No printing
386 if (fitStatus == 0) {
387 TF1* fit = (TF1*)hdEdx_PR[fiRun]->GetFunction("gaus");
388 mean = fit->GetParameter(1);
389 meanError = fit->GetParError(1);
390 sigma = fit->GetParameter(2);
391 sigmaError = fit->GetParError(2);
392 hdEdx_PR[fiRun]->GetXaxis()->SetRangeUser(mean - 7 * sigma, mean + 7 * sigma);
393 fit->Clear();
394 }
395 }
396
397 TotRunN.push_back(fCurrentRunNum);
398 TotMean.push_back(mean);
399 TotMeanE.push_back(meanError);
400 TotSigma.push_back(sigma);
401 TotSigmaE.push_back(sigmaError);
402
403 if (fiRun % 10 == 0)((TH1D*)fBasic->FindObject("hRunGainPR"))->GetXaxis()->SetBinLabel(fiRun + 1, Form("%d", TotRunN.at(fiRun)));
404 ((TH1D*)fBasic->FindObject("hRunGainPR"))->SetBinContent(fiRun + 1, fcRunGain);
405 ((TH1D*)fBasic->FindObject("hRunGainPR"))->SetBinError(fiRun + 1, 0.001 * fcRunGain); // no meaning but histogramming only
406
407 } else if (level == "AR") {
408
409 const Int_t allNRun = TotMean.size();
410 ((TH1D*)fBasic->FindObject("hRunGainPR"))->GetXaxis()->SetRange(1, allNRun);
411
412 TH1D* hFitdEdxMeanPR = new TH1D("hFitdEdxMeanPR", "dE/dx(nohad-sat) #mu via fit vs. Runs", allNRun, 0, allNRun);
413 hFitdEdxMeanPR->GetYaxis()->SetRangeUser(0.90, 1.10);
414 hFitdEdxMeanPR->SetMarkerStyle(21);
415 hFitdEdxMeanPR->SetMarkerColor(kRed);
416 hFitdEdxMeanPR->SetMarkerSize(1);
417 hFitdEdxMeanPR->GetXaxis()->SetTitle("Run numbers");
418 hFitdEdxMeanPR->GetYaxis()->SetTitle("dEdx mean (fit)");
419 hFitdEdxMeanPR->GetXaxis()->LabelsOption("v");
420
421 TH1D* hFitdEdxSigmaPR = new TH1D("hFitdEdxSigmaPR", "dE/dx(nohad-sat) #sigma via fit vs. Runs", allNRun, 0, allNRun);
422 hFitdEdxSigmaPR->GetYaxis()->SetRangeUser(0, 0.30);
423 hFitdEdxSigmaPR->SetMarkerStyle(21);
424 hFitdEdxSigmaPR->SetMarkerColor(kRed);
425 hFitdEdxSigmaPR->SetMarkerSize(1);
426 hFitdEdxSigmaPR->GetXaxis()->SetTitle("Run numbers");
427 hFitdEdxSigmaPR->GetYaxis()->SetTitle("dEdx sigma (fit)");
428 hFitdEdxSigmaPR->GetXaxis()->LabelsOption("v");
429
430 for (Int_t i = 0; i < allNRun; i++) {
431
432 if (i % 10 == 0)hFitdEdxMeanPR->GetXaxis()->SetBinLabel(i + 1, Form("%d", TotRunN.at(i)));
433 hFitdEdxMeanPR->SetBinContent(i + 1, TotMean.at(i));
434 hFitdEdxMeanPR->SetBinError(i + 1, TotMeanE.at(i));
435
436 if (i % 10 == 0)hFitdEdxSigmaPR->GetXaxis()->SetBinLabel(i + 1, Form("%d", TotRunN.at(i)));
437 hFitdEdxSigmaPR->SetBinContent(i + 1, TotSigma.at(i));
438 hFitdEdxSigmaPR->SetBinError(i + 1, TotSigmaE.at(i));
439 }
440
441 fBasic->Add(hFitdEdxMeanPR);
442 fBasic->Add(hFitdEdxSigmaPR);
443
444 TH1D* hdEdxFit_allRun = (TH1D*)(fBasic->FindObject(Form("hdEdx_AR"))->Clone("hdEdxFit_allRun"));
445 if (hdEdxFit_allRun->GetEntries() > 100) {
446 hdEdxFit_allRun->Fit("gaus", "Q");
447 TF1* hGfit = (TF1*)hdEdxFit_allRun->GetFunction("gaus");
448 Double_t meanGFit = hGfit->GetParameter(1);
449 Double_t sigmaGFit = hGfit->GetParameter(2);
450 hdEdxFit_allRun->GetXaxis()->SetRangeUser(meanGFit - 7 * sigmaGFit, meanGFit + 7 * sigmaGFit);
451 hdEdxFit_allRun->SetFillColor(kYellow);
452 hdEdxFit_allRun->SetStats(kTRUE);
453 fBasic->Add(hdEdxFit_allRun);
454 }
455
456 TH1D* hdEdxMeanVsMomentum = new TH1D("hdEdxMeanVsMomentum", "dEdx-mean vs P bins (BW = 250MeV)", 64, -8.0, 8.0);
457 hdEdxMeanVsMomentum->GetXaxis()->SetTitle("Track Momentum");
458 hdEdxMeanVsMomentum->GetYaxis()->SetTitle("dEdx Mean");
459
460 TH1D* hdEdxSigmaVsMomentum = new TH1D("hdEdxSigmaVsMomentum", "dEdx-sigma vs P bins (BW = 250MeV)", 64, -8.0, 8.0);
461 hdEdxSigmaVsMomentum->GetXaxis()->SetTitle("Track Momentum");
462 hdEdxSigmaVsMomentum->GetYaxis()->SetTitle("dEdx Sigma");
463
464 for (int ip = 0; ip < 32; ip++) {
465 Int_t nTrack = ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", ip)))->GetEntries();
466 ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", ip)))->SetFillColor(kYellow);
467 Double_t iPMean = 1.0, iPSigma = 0.0;
468 if (nTrack > 100) {
469 ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", ip)))->Fit("gaus", "0");
470 iPMean = ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", ip)))->GetFunction("gaus")->GetParameter(1);
471 iPSigma = ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", ip)))->GetFunction("gaus")->GetParameter(2);
472 }
473 hdEdxMeanVsMomentum->SetBinContent(32 + ip + 1, iPMean);
474 hdEdxSigmaVsMomentum->SetBinContent(32 + ip + 1, iPSigma);
475 }
476
477 for (int ip = 0; ip < 32; ip++) {
478 Int_t nTrack = ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Posi_Pbin_AR%d", ip)))->GetEntries();
479 ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Elec_Pbin_AR%d", ip)))->SetFillColor(kYellow);
480 Double_t iPMean = 1.0, iPSigma = 0.0;
481
482 if (nTrack > 100) {
483 ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Elec_Pbin_AR%d", ip)))->Fit("gaus", "0");
484 iPMean = ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Elec_Pbin_AR%d", ip)))->GetFunction("gaus")->GetParameter(1);
485 iPSigma = ((TH1D*)fPRdEdxinP->FindObject(Form("hdEdx_Elec_Pbin_AR%d", ip)))->GetFunction("gaus")->GetParameter(2);
486 }
487 hdEdxMeanVsMomentum->SetBinContent(32 - ip, iPMean);
488 hdEdxSigmaVsMomentum->SetBinContent(32 - ip, iPSigma);
489 }
490
491 fBasic->Add(hdEdxSigmaVsMomentum);
492 fBasic->Add(hdEdxMeanVsMomentum);
493 } else {
494 B2ERROR("RunGain >> NO-REQUEST-FOUND for PR or AR level plots, exiting..");
495 }
496}
497
498//-----------------------------------------
500{
501
502 B2INFO("Terminating plots for all runs ");
503 fFileOutput = new TFile(Form("fvalidate%s", fOutFileName.data()), "RECREATE");
504 if (fCollType != "hadron")ExtractHistograms("AR");
505 fFileOutput->cd();
506 fBasic->Write("ARBasics", 1);
507 if (fCollType == "radbhabha" or fCollType == "bhabha") {
508 fPRdEdx->Write("PRdedx", 1);
509 fPRdEdxinP->Write("MeanSigmavsP", 1);
510 }
511 fFileOutput->Close();
512}
513
514//-----------------------------------------
516{
517
518 if (std::abs(mTrack->getD0()) >= fD0Window || std::abs(mTrack->getZ0()) >= fZ0Window)return kFALSE;
519 return kTRUE;
520
521}
double R
typedef autogenerated by FFTW
Debug output for CDCDedxPID module.
double getDedx() const
Get dE/dx truncated mean for this track.
int getNLayerHits() const
Return the number of layer hits for this track.
double getCosTheta() const
Return cos(theta) for this track.
double getRunGain() const
Return the run gain for this track.
double getChi(int i) const
Return the PID (chi) value.
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
double getMomentum() const
Return the track momentum valid in the CDC.
std::vector< Double_t > TotMean
Mean of dedx by Fit.
Double_t fnBinsdedxUE
up edge of dedx
void FillHistograms(CDCDedxTrack *dedxTrack, const TrackFitResult *mTrack)
Filling histograms This will fill histogram defined histograms in above function.
TFile * fFileOutput
Write final objects to file for RG.
std::string fOutFileName
name of output ROOT file
Double_t fcRunGain
existing run gain
Int_t fCurrentRunNum
current run number
Int_t fnRunCounter
Total runs used counter.
TList * fBasic
List of basic histos.
virtual void initialize() override
Initialize This is inherited from base class.
std::vector< TH1D * > hdEdx_PR
histogram array per run
virtual void event() override
function to execute event (event by event) This is inherited from base class
std::vector< Double_t > TotSigma
Sigma of dedx by Fit.
CDCDedxValidationModule()
Default constructor contain list of members with initial values.
TList * fPRdEdxinP
list per run dedx in P histos
Double_t fTrkEoverP
E/p ratio for cut.
virtual void endRun() override
function is called after each event This is inherited from base class
std::vector< Int_t > TotRunN
veector array of runs processed
virtual void terminate() override
Terminate after all data processed This is inherited from base class.
std::vector< Double_t > TotSigmaE
Sigma Error of dedx by Fit.
std::vector< Double_t > TotMeanE
Mean Error of dedx by Fit.
virtual void beginRun() override
Function to execute each run This is inherited from base class.
StoreArray< CDCDedxTrack > m_cdcDedxTracks
Data members for objects, cuts and others.
Double_t fnBinsdedxLE
low edge of dedx
TList * fPRdEdx
List of per run dedx histos.
void ExtractHistograms(TString level)
Extracting histogram and some calculation Higher level histograms are filled after each run or full p...
void DefineHistograms(TString level, Int_t iR)
Definition of histograms This contain a list of histogram for validation.
Int_t fnRuns
Number of runs ref.
Int_t fiRun
Current run number.
Int_t fnBinsdedx
nbin of dedx range
std::string fCollType
collision type
Bool_t IsSelectedTrack(const TrackFitResult *mTrack)
Track selection A clean way to implement selections on tracks (so far few only)
static const ChargedStable pion
charged pion particle
Definition Const.h:661
static const ChargedStable electron
electron particle
Definition Const.h:659
ECL cluster data.
Definition ECLCluster.h:27
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
Definition ECLCluster.h:351
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition ECLCluster.cc:23
@ c_nPhotons
CR is split into n photons (N1)
Definition ECLCluster.h:41
HistoModule()
Constructor.
Definition HistoModule.h:32
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
Values of the result of a track fit with a given particle hypothesis.
double getPhi() const
Getter for phi0 with CDF naming convention.
short getChargeSign() const
Return track charge (1 or -1).
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Class that bundles various TrackFitResults.
Definition Track.h:25
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.