Belle II Software development
HadronPrep.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/calibration/CDCdEdx/HadronPrep.h>
10using namespace Belle2;
11
13{
14 m_bgBins = 10;
15 m_bgMin = 2.85;
16 m_bgMax = 40;
17 m_cosBins = 8;
18 m_cosMin = 0.0;
19 m_cosMax = 0.95;
20 m_cut = 0.5;
21}
22
23HadronPrep::HadronPrep(int bgbins, double lowerbg, double upperbg, int cosbins, double lowercos, double uppercos, double cut)
24{
25 m_bgBins = bgbins;
26 m_bgMin = lowerbg;
27 m_bgMax = upperbg;
28 m_cosBins = cosbins;
29 m_cosMin = lowercos;
30 m_cosMax = uppercos;
31 m_cut = cut;
32}
33
34void HadronPrep::prepareSample(std::shared_ptr<TTree> hadron, TFile*& outfile, const std::string& suffix,
35 const std::string& pdg, bool ismakePlots, bool correct)
36{
37
38 std::map<int, std::vector<TH1F*>> hdedx_bgcosth;
39 std::vector<TH2F*> hdedxvscosth_bg;
40
41 //define histograms
42 for (int i = 0; i < m_bgBins; ++i) {
43
44 std::string title = Form("%s_dedxvscos_bg_%d", pdg.data(), i);
45 hdedxvscosth_bg.push_back(new TH2F(title.data(), Form("%s;costh;dEdx", title.data()),
46 440, -1.1, 1.1, 2600, -1.0, 25));
47
48 defineHisto(hdedx_bgcosth[i], Form("dedx_bg_%d_cos", i), pdg.data());
49 }
50
51
52 // Create some containers to calculate averages
53 for (int i = 0; i < m_bgBins; ++i) {
54 for (int j = 0; j < m_cosBins; ++j) {
55 m_sumcos[i].push_back(0.);
56 m_sumbg[i].push_back(0.);
57 m_sumsize[i].push_back(0.);
58 m_means[i].push_back(0.);
59 m_errors[i].push_back(0.);
60
61 }
62 }
63
64 // // --------------------------------------------------
65 // // LOOP OVER EVENTS AND FILL CONTAINERS
66 // // --------------------------------------------------
67 // // Fill the histograms to be fitted
68
69 double dedxnosat; // dE/dx w/o HS correction (use to get new parameters)
70 double p; // track momentum
71 double costh; // cosine of track polar angle
72 double timereso;
73 int nhit; // number of hits on this track
74
75 hadron->SetBranchAddress("dedxnosat", &dedxnosat); //must be without any HS correction
76 hadron->SetBranchAddress("p", &p);
77 hadron->SetBranchAddress("costh", &costh);
78 hadron->SetBranchAddress("timereso", &timereso);
79 hadron->SetBranchAddress("nhits", &nhit);
80
81 HadronBgPrep prep;
82 double mass = prep.getParticleMass(pdg);
83 if (mass == 0.0) B2FATAL("Mass of particle " << pdg.data() << " is zero");
84
85 CDCDedxHadSat had;
86
87 if (correct) {
88 std::string infile = "sat-pars.fit.txt";
89 std::ifstream parfile(infile);
90
91 if (!parfile.fail()) {
92 B2INFO("new parameters are using");
93 had.setParameters(infile);
94 } else {
95 B2INFO("default parameters are using");
96 had.setParameters();
97 }
98 }
99
100 for (unsigned int index = 0; index < hadron->GetEntries(); ++index) {
101
102 hadron->GetEntry(index);
103
104 double bg; // track beta-gamma
105 bg = fabs(p) / mass;
106
107 double uncosth = costh;
108 costh = fabs(costh);
109
110 // clean up bad events and restrict the momentum range
111 if (nhit < 0 || nhit > 100)continue;
112 if (dedxnosat <= 0)continue;
113 if (costh != costh)continue;
114 if (bg <= m_bgMin || bg >= m_bgMax)continue;
115 if (costh <= m_cosMin || costh >= m_cosMax)continue;
116
117 if (pdg == "proton") if ((dedxnosat - 0.45) * abs(p) * abs(p) < m_cut)continue;
118
119 int bgBin = static_cast<int>((bg - m_bgMin) / (m_bgMax - m_bgMin) * m_bgBins);
120 bgBin = std::max(0, std::min(bgBin, m_bgBins - 1));
121
122 int cosBin = static_cast<int>((costh - m_cosMin) / (m_cosMax - m_cosMin) * m_cosBins);
123 cosBin = std::max(0, std::min(cosBin, m_cosBins - 1));
124
125 double dedx_new = dedxnosat;
126
127 if (correct) dedx_new = had.D2I(costh, had.I2D(costh, 1.0) * dedxnosat);
128
129 hdedx_bgcosth[bgBin][cosBin]->Fill(dedx_new);
130 hdedxvscosth_bg[bgBin]->Fill(uncosth, dedx_new);
131
132 m_sumcos[bgBin][cosBin] += costh;
133 m_sumbg[bgBin][cosBin] += bg;
134 m_sumsize[bgBin][cosBin] += 1;
135
136 }// end of event loop
137
138 // // --------------------------------------------------
139 // // FIT IN BINS OF BETA-GAMMA AND COS(THETA)
140 // // --------------------------------------------------
141 // // fit the histograms with Gaussian functions
142 // // and extract the means and errors
143
144 setPars(outfile, hdedx_bgcosth, pdg.data());
145
146 if (ismakePlots) {
147
148 plotDist(hdedx_bgcosth, Form("fits_dedx_inbg_%s", suffix.data()), pdg.data());
149 plotDist(hdedxvscosth_bg, Form("dist_dedx_vscos_scat_%s", suffix.data()), pdg.data());
150 plotGraph(suffix.data(), pdg.data());
151 }
152
153 hdedxvscosth_bg.clear();
154 hdedx_bgcosth.clear();
155 clear();
156
157}
158
159//------------------------------------
160void HadronPrep::defineHisto(std::vector<TH1F*>& htemp, const std::string& svar, const std::string& pdg)
161{
162 for (int j = 0; j < m_cosBins; ++j) {
163 // initialize the histograms
164 std::string title = Form("%s_%s_%d", pdg.data(), svar.data(), j);
165 double dedxMax = 4.0;
166 int cosbins = int(50 * dedxMax);
167 if (pdg == "pion") {
168 dedxMax = 2.0; cosbins = int(100 * dedxMax);
169 }
170 if (pdg == "proton") dedxMax = 20.0;
171 if (svar.substr(0, svar.find("_")) == "dedx")
172 htemp.push_back(new TH1F(title.data(), title.data(), cosbins, 0, dedxMax));
173 else
174 htemp.push_back(new TH1F(title.data(), title.data(), cosbins, -30, 30));
175
176 }
177}
178
179//------------------------------------
180void HadronPrep::plotDist(std::map<int, std::vector<TH1F*>>& hist, const std::string& sname, const std::string& pdg)
181{
182
183 TCanvas ctmp(Form("cdcdedx_%s_%s", sname.data(), pdg.data()), "", 1200, 1200);
184 ctmp.Divide(2, 2);
185 ctmp.SetBatch(kTRUE);
186
187 std::stringstream psname;
188 psname << Form("plots/HadronSat/%s_%s.pdf[", sname.data(), pdg.data());
189 ctmp.Print(psname.str().c_str());
190 psname.str("");
191 psname << Form("plots/HadronSat/%s_%s.pdf", sname.data(), pdg.data());
192
193 for (int i = 0 ; i < m_bgBins; ++i) {
194 for (int j = 0; j < m_cosBins; ++j) {
195 hist[i][j]->SetFillColor(kYellow - 9);
196
197 ctmp.cd(j % 4 + 1);
198 hist[i][j]->Draw();
199
200 if (m_cosBins <= 4) {
201 if (j + 1 == m_cosBins) ctmp.Print(psname.str().c_str());
202 } else if ((j + 1) % 4 == 0 || (j + 1 == m_cosBins)) {
203 ctmp.Print(psname.str().c_str());
204 ctmp.Clear("D");
205 }
206 }
207 }
208 psname.str("");
209 psname << Form("plots/HadronSat/%s_%s.pdf]", sname.data(), pdg.data());
210 ctmp.Print(psname.str().c_str());
211}
212
213//------------------------------------
214void HadronPrep::plotDist(std::vector<TH2F*>& hist, const std::string& sname, const std::string& pdg)
215{
216
217 TCanvas ctmp(Form("cdcdedx_%s_%s", sname.data(), pdg.data()), "", 1200, 1200);
218 ctmp.Divide(2, 2);
219 ctmp.SetBatch(kTRUE);
220
221 std::stringstream psname;
222 psname << Form("plots/HadronSat/%s_%s.pdf[", sname.data(), pdg.data());
223 ctmp.Print(psname.str().c_str());
224 psname.str("");
225 psname << Form("plots/HadronSat/%s_%s.pdf", sname.data(), pdg.data());
226 double bgstep = (m_bgMax - m_bgMin) / m_bgBins;
227 double dmax = 2.5;
228 if (pdg == "kaon")dmax = 3.5;
229 else if (pdg == "proton")dmax = 15.5;
230
231 for (int i = 0 ; i < m_bgBins; ++i) {
232 ctmp.cd(i % 4 + 1);
233 double frombg = m_bgMin + i * bgstep;
234 double tobg = frombg + bgstep;
235
236 hist[i]->SetTitle(Form("%s, pF=(%0.02f,%0.02f)", hist[i]->GetTitle(), frombg, tobg));
237 hist[i]->GetYaxis()->SetRangeUser(-0.5, dmax);
238 hist[i]->DrawCopy("colz");
239 if ((i + 1) % 4 == 0 || (i + 1) == m_bgBins) {
240 ctmp.Print(psname.str().c_str());
241 ctmp.Clear("D");
242 }
243 }
244 psname.str("");
245 psname << Form("plots/HadronSat/%s_%s.pdf]", sname.data(), pdg.data());
246 ctmp.Print(psname.str().c_str());
247}
248
249//----------------------------------------
250void HadronPrep::setPars(TFile*& outfile, std::map<int, std::vector<TH1F*>>& hdedx_bgcosth, const std::string& pdg)
251{
252 outfile->cd();
253 TTree* satTree = new TTree(Form("%s", pdg.data()), "dE/dx means and errors");
254 double satbg; // beta-gamma value for this bin
255 double satcosth; // cos(theta) value for this bin
256 double satdedx; // mean dE/dx value for this bin
257 double satdedxerr; // error on ^
258 double satdedxwidth;// width of ^ distribution
259 double satbg_avg; // average beta-gamma value for this sample
260 double satcosth_avg; // average cos(theta) value for this sample
261
262 satTree->Branch("bg", &satbg, "bg/D");
263 satTree->Branch("costh", &satcosth, "costh/D");
264
265 //dEdx related variables
266 satTree->Branch("dedx", &satdedx, "dedx/D");
267 satTree->Branch("dedxerr", &satdedxerr, "dedxerr/D");
268 satTree->Branch("dedxwidth", &satdedxwidth, "dedxwidth/D");
269
270 // Other variables
271 satTree->Branch("bg_avg", &satbg_avg, "bg_avg/D");
272 satTree->Branch("costh_avg", &satcosth_avg, "costh_avg/D");
273
274 double bgstep = (m_bgMax - m_bgMin) / m_bgBins;
275 double cosstep = (m_cosMax - m_cosMin) / m_cosBins;
276
277 for (int i = 0; i < m_bgBins; ++i) {
278 for (int j = 0; j < m_cosBins; ++j) {
279
280 // fill some details for this bin
281 satbg = m_bgMin + 0.5 * bgstep + i * bgstep;
282 satcosth = m_cosMin + 0.5 * cosstep + j * cosstep;
283
284 if (m_sumsize[i][j] > 0) {
285 satbg_avg = m_sumbg[i][j] / m_sumsize[i][j];
286 satcosth_avg = m_sumcos[i][j] / m_sumsize[i][j];
287 } else {
288 satbg_avg = satcosth_avg = 0.0;
289 }
290
291 //1. -------------------------
292 // fit the dE/dx distribution in bins of beta-gamma and cosine
293 HadronBgPrep prep;
294 gstatus stats;
295 prep.fit(hdedx_bgcosth[i][j], pdg.data(), stats);
296 if (stats == OK) {
297 satdedx = m_means[i][j] = hdedx_bgcosth[i][j]->GetFunction("gaus")->GetParameter(1);
298 satdedxerr = m_errors[i][j] = hdedx_bgcosth[i][j]->GetFunction("gaus")->GetParError(1);
299 satdedxwidth = hdedx_bgcosth[i][j]->GetFunction("gaus")->GetParameter(2);
300 } else { satdedx = 0.0; satdedxerr = 0.0; satdedxwidth = 0.0;}
301
302
303 //Be careful to not set for every particle
304 if (satdedx > m_dedxmax) m_dedxmax = satdedx;
305 if (satdedx < m_dedxmin) m_dedxmin = satdedx;
306
307 // fill the tree for this bin
308 satTree->Fill();
309 }
310 }
311 outfile->Write();
312}
313
314//------------------------------------
315void HadronPrep::plotGraph(const std::string& sname, const std::string& pdg)
316{
317 std::vector<double> cbcenters(m_cosBins), cberrors(m_cosBins);
318 double cosstep = (m_cosMax - m_cosMin) / m_cosBins;
319 double bgstep = (m_bgMax - m_bgMin) / m_bgBins;
320
321 for (int j = 0; j < m_cosBins; ++j) {
322 cbcenters[j] = m_cosMin + 0.5 * cosstep + j * cosstep;
323 cberrors[j] = 0;
324 }
325
326 // Plot the dE/dx means vs. cos(theta) for validation
327 TCanvas ctmp6(Form("cdcdedx_costh_%s_%s", sname.data(), pdg.data()), "", 500, 500);
328 ctmp6.SetGridy(1);
329
330 TH1F base("base", "bla-bla", 25, 0, 1.0);
331 base.SetTitle(Form("%s: dE/dx fit means vs. cos(#theta);cos(#theta);dE/dx-fit mean", pdg.data()));
332 base.GetXaxis()->SetTitleOffset(1.2);
333
334 double margin = 0.1;
335 if (pdg == "proton") margin = 0.95;
336
337 base.SetMaximum(m_dedxmax + margin);
338 base.SetMinimum(m_dedxmin - margin);
339 base.SetStats(0);
340 base.DrawCopy("");
341
342 TLegend legend(0.75, 0.75, 0.95, 0.9);
343
344 std::vector<TGraphErrors*> gdedx_costh(m_bgBins);
345
346 legend.SetBorderSize(0);
347
348 for (int i = 0; i < m_bgBins; ++i) {
349 std::vector<double> mean_d(m_cosBins);
350 std::vector<double> error_d(m_cosBins);
351
352 for (int j = 0; j < m_cosBins; ++j) {
353 mean_d [j] = m_means[i][j];
354 error_d[j] = m_errors[i][j];
355 }
356 gdedx_costh[i] = new TGraphErrors(m_cosBins, cbcenters.data(), mean_d.data(), cberrors.data(), error_d.data());
357
358 gdedx_costh[i]->SetMarkerSize(0.9);
359 gdedx_costh[i]->SetMarkerColor(50 + i * 3);
360 gdedx_costh[i]->SetLineColor(i + 1);
361 gdedx_costh[i]->SetLineWidth(1);
362 gdedx_costh[i]->Draw("same");
363 legend.AddEntry(gdedx_costh[i], Form("bg = (%0.03f - %0.03f)", m_bgMin + i * bgstep, m_bgMin + (i + 1) * bgstep), "lep");
364 }
365
366 legend.Draw("same");
367
368 std::stringstream psname;
369 psname.str("");
370 psname << Form("plots/HadronSat/gr_dedx_vscosth_%s_%s.pdf", sname.data(), pdg.data());
371 ctmp6.SaveAs(psname.str().c_str());
372
373 for (int i = 0; i < m_bgBins; ++i)
374 delete gdedx_costh[i];
375}
376
377
379{
380 m_sumcos.clear();
381 m_sumbg.clear();
382 m_sumsize.clear();
383 m_means.clear();
384 m_errors.clear();
385}
Class to hold the hadron saturation functions.
double I2D(double cosTheta, double I) const
hadron saturation parameterization part 2
double D2I(double cosTheta, double D) const
hadron saturation parameterization part 1
void setParameters()
set the parameters
Class to prepare sample for fitting in beta gamma bins.
void fit(TH1F *&hist, const std::string &pdg, gstatus &status)
function to fit the histograms
double getParticleMass(const std::string &particle)
function to get the particle mass
void prepareSample(std::shared_ptr< TTree > hadron, TFile *&outfile, const std::string &suffix, const std::string &pdg, bool ismakePlots, bool correct)
function to prepare sample for monitoring plots, bg curve fitting and sigma vs ionz fitting
Definition HadronPrep.cc:34
double m_dedxmax
variables to set maximum dedx mean
Definition HadronPrep.h:104
int m_cosBins
bins for dedx histogram
Definition HadronPrep.h:111
std::map< int, std::vector< double > > m_errors
error variable
Definition HadronPrep.h:100
double m_cosMax
max range of dedx
Definition HadronPrep.h:113
std::map< int, std::vector< int > > m_sumsize
variables for size
Definition HadronPrep.h:102
HadronPrep()
Constructor: Sets the description, the properties and the parameters of the algorithm.
Definition HadronPrep.cc:12
double m_dedxmin
variables to set minimum dedx mean
Definition HadronPrep.h:105
void plotDist(std::map< int, std::vector< TH1F * > > &hist, const std::string &sname, const std::string &pdg)
function to plot the map of histograms
int m_bgBins
bins for dedx histogram
Definition HadronPrep.h:107
void setPars(TFile *&outfile, std::map< int, std::vector< TH1F * > > &hdedx_bgcosth, const std::string &pdg)
function to fill the parameters like mean and reso in the tree
double m_bgMax
max range of dedx
Definition HadronPrep.h:109
void plotGraph(const std::string &sname, const std::string &pdg)
function to make graph dedx vs cos in different bg bins
void defineHisto(std::vector< TH1F * > &htemp, const std::string &title, const std::string &pdg)
function to define histograms
double m_bgMin
min range of dedx
Definition HadronPrep.h:108
double m_cosMin
min range of dedx
Definition HadronPrep.h:112
double m_cut
cut to clean protons
Definition HadronPrep.h:115
void clear()
function to clear the variables
std::map< int, std::vector< double > > m_means
mean variable
Definition HadronPrep.h:99
std::map< int, std::vector< double > > m_sumbg
variables to add bg values
Definition HadronPrep.h:98
std::map< int, std::vector< double > > m_sumcos
variables to add cos values
Definition HadronPrep.h:97
Abstract base class for different kinds of events.