Belle II Software development
create_dedx_PDFs.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 <framework/gearbox/Const.h>
10#include <reconstruction/dataobjects/DedxConstants.h>
11
12#include <TFile.h>
13#include <TTree.h>
14#include <TCanvas.h>
15#include <TH2F.h>
16#include <TApplication.h>
17#include <TROOT.h>
18#include <TStyle.h>
19#include <TSystem.h>
20
21#include <iostream>
22
23using namespace Belle2;
24
26const int num_dedx_bins = 100;
27
28
34int main(int argc, char* argv[])
35{
36 const int num_pdg_codes = Const::ChargedStable::c_SetSize;
37
38 if (argc < 2 or argc > 3) {
39 std::cerr << "Usage: " << argv[0] << " INPUT_FILE [OUTPUT_FILE]\n\n";
40 std::cerr <<
41 "Generates PDFs for the DedxPID module, requires an input .root file created by running DedxPID with enableDebug=True.\n";
42 std::cerr << "If OUTPUT_FILE is not given, histograms will be drawn instead of being saved.\n";
43 return 1;
44 }
45
46 const TString input_filename(argv[1]);
47 TString output_filename;
48 if (argc == 3)
49 output_filename = TString(argv[2]);
50
51 //TApplication eats command line arguments, don't use argc/argv after that
52 TApplication app("noname", &argc, argv);
53 gROOT->SetStyle("Plain");
54 gStyle->SetPalette(1);
55
56 TFile* file = new TFile(input_filename, "READ");
57 TTree* tree = dynamic_cast<TTree*>(file->Get("tree"));
58 if (!tree) {
59 std::cerr << "'tree' not found, aborting.\n";
60 return 1;
61 }
62
63 TCanvas* c1 = 0;
64
65
66 const bool save_graphs = !output_filename.IsNull();
67 TFile* pdffile = 0;
68 if (save_graphs) {
69 pdffile = new TFile(output_filename, "RECREATE");
70 } else {
71 c1 = new TCanvas("c1", input_filename, 600, 900);
72 c1->Divide(num_pdg_codes, Dedx::c_num_detectors);
73 }
74
75 //create momentum binning (larger bins for high momentum)
76 const int num_p_bins = 69;
77 double pbins[num_p_bins + 1];
78 pbins[0] = 0.0; pbins[1] = 0.05;
79 std::cout << ": " << pbins[0] << "\n";
80 std::cout << ": " << pbins[1] << "\n";
81 for (int bin = 2; bin <= num_p_bins; bin++) {
82 if (bin <= 19)
83 pbins[bin] = 0.025 + 0.025 * bin;
84 else if (bin <= 59)
85 pbins[bin] = pbins[19] + 0.05 * (bin - 19);
86 else
87 pbins[bin] = pbins[59] + 0.3 * (bin - 59);
88 std::cout << ": " << pbins[bin] << "\n";
89 }
90
91 for (int i = 0; i < 2; i++) { //normal/truncated
92 const bool use_truncated_mean = (i == 1);
93 const char* suffix = use_truncated_mean ? "_trunc" : "";
94
95 for (int iPart = 0; iPart < num_pdg_codes; iPart++) {
96 TH2F* hists[Dedx::c_num_detectors]; //one for each Detector
97 const int pdg_code = Const::chargedStableSet.at(iPart).getPDGCode();
98 //const int num_dedx_bins = use_truncated_mean ? (num_dedx_bins_hits) : num_dedx_bins_hits;
99
100 // first prepare the histograms for the VXD
101 for (int detector = 0; detector <= Dedx::c_SVD; detector++) {
102 std::cout << "i " << i << ", d" << detector << "\n";
103 const char* flayer_var = "VXDDedxTracks.dedxLayer";
104 TString flayer_selection;
105 double dedx_cutoff = 0;
106 switch (Dedx::Detector(detector)) {
107 case Dedx::c_PXD:
108 flayer_selection = TString::Format("%s < 0 && %s >= -2", flayer_var, flayer_var);
109 dedx_cutoff = 10e3;
110 break;
111 case Dedx::c_SVD:
112 flayer_selection = TString::Format("%s < -2", flayer_var);
113 dedx_cutoff = 5e6;
114 break;
115 case Dedx::c_CDC:
116 break;
117 }
118
119 const TString histname = TString::Format("hist_d%i_%i%s", detector, pdg_code, suffix);
120 const char* varname = use_truncated_mean ? "VXDDedxTracks.m_dedxAvgTruncated" : "VXDDedxTracks.dedx";
121
122 hists[detector] = new TH2F(histname.Data(), histname.Data(),
123 num_p_bins, pbins,
124 num_dedx_bins, 0, dedx_cutoff);
125 hists[detector]->Sumw2(); //save weights (important for fitting)
126 if (use_truncated_mean)
127 tree->Project(histname.Data(),
128 TString::Format("%s[][%i]:VXDDedxTracks.m_p", varname, detector),
129 TString::Format("%s < %g && abs(VXDDedxTracks.m_pdg) == %i && %s", varname, dedx_cutoff, pdg_code, flayer_selection.Data()));
130 else
131 tree->Project(histname.Data(),
132 TString::Format("%s:VXDDedxTracks.m_p", varname),
133 TString::Format("%s < %g && abs(VXDDedxTracks.m_pdg) == %i && %s", varname, dedx_cutoff, pdg_code, flayer_selection.Data()));
134 } //detector type
135
136 // now add the CDC histograms
137 int detector = (int)Dedx::c_CDC;
138 std::cout << "i " << i << ", d" << detector << "\n";
139 double dedx_cutoff = 0;
140 dedx_cutoff = 3.0;
141
142 const TString histname = TString::Format("hist_d%i_%i%s", detector, pdg_code, suffix);
143 const char* varname = use_truncated_mean ? "CDCDedxTracks.m_dedxAvgTruncated" : "CDCDedxTracks.m_lDedx";
144
145 hists[detector] = new TH2F(histname.Data(), histname.Data(),
146 num_p_bins, pbins,
147 num_dedx_bins, 0, dedx_cutoff);
148 hists[detector]->Sumw2(); //save weights (important for fitting)
149 if (use_truncated_mean)
150 tree->Project(histname.Data(),
151 TString::Format("%s[][%i]:CDCDedxTracks.m_pCDC", varname, detector),
152 TString::Format("CDCDedxTracks.m_lNHitsUsed > 15 && %s < %g && abs(CDCDedxTracks.m_pdg) == %i", varname, dedx_cutoff, pdg_code));
153 else
154 tree->Project(histname.Data(),
155 TString::Format("%s:CDCDedxTracks.m_pCDC", varname),
156 TString::Format("CDCDedxTracks.m_lNHitsUsed > 15 && %s < %g && abs(CDCDedxTracks.m_pdg) == %i", varname, dedx_cutoff, pdg_code));
157
158 //{{{ for each momentum bin, normalize pdf (disable if you want to keep the originals, e.g. for fitting)
159 if (true) {
160 for (int d = 0; d < Dedx::c_num_detectors; d++) {
161 for (int pbin = 0; pbin <= num_p_bins + 1; pbin++) {
162 // get number of entries in this pbin
163 double integral = 0;
164 for (int dedxbin = 0; dedxbin <= num_dedx_bins + 1; dedxbin++) {
165 integral += hists[d]->GetBinContent(pbin, dedxbin);
166 }
167
168 if (integral == 0)
169 continue; //nothing to do
170
171 // normalize this pbin to 1
172 const double normal_width = -(hists[d]->GetYaxis()->GetBinLowEdge(1) - hists[d]->GetYaxis()->GetBinUpEdge(
173 num_dedx_bins)) / num_dedx_bins;
174 for (int dedxbin = 0; dedxbin <= num_dedx_bins + 1; dedxbin++) {
175 hists[d]->SetBinContent(pbin, dedxbin, hists[d]->GetBinContent(pbin,
176 dedxbin) / integral * hists[d]->GetYaxis()->GetBinWidth(dedxbin) / normal_width);
177 //std::cout << d << ", " << dedxbin << ": " << hists[d]->GetYaxis()->GetBinWidth(dedxbin)/normal_width << "\n";
178 }
179 }
180 }
181 }
182 //}}}
183
184 for (int d = 0; d < Dedx::c_num_detectors; d++) {
185 if (save_graphs) {
186 hists[d]->Write();
187 } else if (i == 0) {
188 TPad* pad = static_cast<TPad*>(c1->cd(num_pdg_codes * d + iPart + 1));
189 pad->SetLogz();
190 hists[d]->Draw("colz");
191 c1->Update();
192 }
193 }
194 } //particle type
195 } //normal / trunc
196
197 if (save_graphs) {
198 pdffile->Close();
199 delete pdffile;
200
201 std::cout << "Saved histograms into '" << output_filename.Data() << "'.\n";
202
203 //no output, so we'll just exit when we're done
204 gSystem->Exit(0);
205 }
206
207 app.Run(true);
208 return 0;
209}
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
Definition: Const.h:615
const ParticleType & at(unsigned int index) const
Return particle at given index, or end() if out of range.
Definition: Const.h:549
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:618
Abstract base class for different kinds of events.