Belle II Software  release-05-02-19
SVDChargeSharingAnalysisModule.cc
1 /*************************************************************************
2 * BASF2 (Belle Analysis Framework 2) *
3 * Copyright(C) 2018 - Belle II Collaboration *
4 * *
5 * Author: The Belle II Collaboration *
6 * Contributors: Mateusz Kaleta *
7 * *
8 * This software is provided "as is" without any warranty. *
9 **************************************************************************/
10 
11 #include <svd/modules/svdChargeSharing/SVDChargeSharingAnalysisModule.h>
12 
13 #include <TDirectory.h>
14 #include <TCollection.h>
15 #include <reconstruction/dataobjects/VXDDedxTrack.h>
16 #include <vxd/geometry/GeoCache.h>
17 #include <mdst/dataobjects/HitPatternVXD.h>
18 
19 #include <boost/foreach.hpp>
20 
21 #include <TCanvas.h>
22 #include <TSystem.h>
23 #include <TLegend.h>
24 #include <TROOT.h>
25 #include <TStyle.h>
26 
27 
28 /* --------------- WARNING ---------------------------------------------- *
29 If you have more complex parameter types in your class then simple int,
30 double or std::vector of those you might need to uncomment the following
31 include directive to avoid an undefined reference on compilation.
32 * ---------------------------------------------------------------------- */
33 // #include <framework/core/ModuleParam.templateDetails.h>
34 
35 
36 
37 using namespace Belle2;
38 
39 //-----------------------------------------------------------------
40 // Register the Module
41 //-----------------------------------------------------------------
42 REG_MODULE(SVDChargeSharingAnalysis)
43 
44 //-----------------------------------------------------------------
45 // Implementation
46 //-----------------------------------------------------------------
47 
49  : Module()
50 {
51  // Set module properties
52  setDescription("Module for monitoring DSSD cluster charge deposition");
53 
54  // Parameter definitions
55  addParam("outputDirName", m_outputDirName, "Name of the output directory.");
56  addParam("outputRootFileName", m_outputRootFileName, "Name of output rootfile.");
57  addParam("useTrackInfo", m_useTrackInfo, "True if using clusters related to tracks in the analysis", bool(true));
58 }
59 
61 {
62 }
63 
65 {
66  m_outputRootFile = new TFile((m_outputDirName + "/" + m_outputRootFileName).c_str(), "RECREATE");
67 
68  //StoreArrays
69  m_svdClusters.isRequired();
70  m_Tracks.isOptional();
71 
72  // Tracks
73  m_histoList_Tracks = new TList;
74  h_nTracks = createHistogram1D("N_tracks", "Number of tracks", 10, 0., 10., "N tracks", "Entries", m_histoList_Tracks);
75  h_TracksPvalue = createHistogram1D("Tracks_Pvalue", "P value of tracks", 10, 0., 1., "PValue", "Entries",
77  h_TracksMomentum = createHistogram1D("Tracks_momentum", "Tracks momentum", 500, 0., 5., "Momentum[GeV/c]", "Entries",
79  h_TracksnSVDhits = createHistogram1D("SVD_hits", "SVD hits", 20, 0., 20., "# SVD hits", "Entries", m_histoList_Tracks);
80 
81  // clCharge, clSize, clSNR etc.
82  TString NameOfHisto;
83  TString TitleOfHisto;
84 
85  for (int i = 0; i < m_nSensorTypes; i++) {
86  m_histoList_clCharge[i] = new TList;
87  m_histoList_clChargeVsMomentum[i] = new TList;
89  m_histoList_clChargeVsSNR[i] = new TList;
90 
91  m_histoList_clSize[i] = new TList;
92  m_histoList_clSizeVsMomentum[i] = new TList;
93  m_histoList_clSizeVsIncidentAngle[i] = new TList;
94  m_histoList_clSizeVsSNR[i] = new TList;
95 
96  m_histoList_clSNR[i] = new TList;
97  m_histoList_clSNRVsMomentum[i] = new TList;
98  m_histoList_clSNRVsIncidentAngle[i] = new TList;
99 
100  TString nameSensorType = "";
101  nameSensorType += m_nameSensorTypes[i];
102 
103  for (int j = 0; j < m_nSides; j++) {
104  TString nameSide = "";
105  nameSide = (j > 0 ? "U" : "V");
106 
107  // Cluster size
108  NameOfHisto = "DSSDClusterSize_" + nameSensorType + "Side_" + nameSide;
109  TitleOfHisto = "DSSD Cluster Size: " + nameSensorType + ", Side " + nameSide;
110  h_clSize[i][j] = createHistogram1D(NameOfHisto, TitleOfHisto, 4, 1, 5, "Cluster Size", "Entries", m_histoList_clSize[i]);
111 
112  NameOfHisto = "DSSDClusterSizeVsMomentum_" + nameSensorType + "Side_" + nameSide;
113  TitleOfHisto = "DSSD Cluster Size vs. Track Momentum: " + nameSensorType + ", Side " + nameSide;
114  h_clSizeVsMomentum[i][j] = createHistogram2D(NameOfHisto, TitleOfHisto, 500, 0., 5., "Track Momentum [GeV/c]", 4, 1, 5,
115  "Cluster Size", m_histoList_clSizeVsMomentum[i]);
116 
117  NameOfHisto = "DSSDClusterSizeVsIncidentAngle_" + nameSensorType + "Side_" + nameSide;
118  TitleOfHisto = "DSSD Cluster Size vs. Incident Angle: " + nameSensorType + ", Side " + nameSide;
119  h_clSizeVsIncidentAngle[i][j] = createHistogram2D(NameOfHisto, TitleOfHisto, 0, 90., 90., "Incident Angle [deg]", 4, 1, 5,
120  "Cluster Size", m_histoList_clSizeVsIncidentAngle[i]);
121 
122  NameOfHisto = "DSSDClusterSizeVsSNR_" + nameSensorType + "Side_" + nameSide;
123  TitleOfHisto = "DSSD Cluster Size vs. SNR: " + nameSensorType + ", Side " + nameSide;
124  h_clSizeVsSNR[i][j] = createHistogram2D(NameOfHisto, TitleOfHisto, 100, 0., 100., "Cluster SNR", 4, 1, 5,
125  "Cluster Size", m_histoList_clSizeVsSNR[i]);
126 
127  for (int k = 0; k < m_nClSizes; k++) {
128  TString nameclSize = "";
129  if (k < 2) nameclSize += k + 1;
130  else nameclSize += ">= 3";
131 
132  // Cluster charge
133  NameOfHisto = "DSSDClusterCharge_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
134  TitleOfHisto = "DSSD Cluster Charge: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
135  h_clCharge[i][j][k] = createHistogram1D(NameOfHisto, TitleOfHisto, m_nBins, m_minCharge, m_maxCharge,
136  "Cluster Charge [e]", "Entries / (375 e)", m_histoList_clCharge[i]);
137 
138  NameOfHisto = "DSSDClusterChargeVsMomentum_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
139  TitleOfHisto = "DSSD Cluster Charge vs. Track Momentum: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
140  h_clChargeVsMomentum[i][j][k] = createHistogram2D(NameOfHisto, TitleOfHisto, 500, 0, 5., "Track Momentum[GeV/c]", m_nBins,
141  m_minCharge, m_maxCharge, "Cluster Charge [e]",
143 
144  NameOfHisto = "DSSDClusterChargeVsIncidentAngle_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
145  TitleOfHisto = "DSSD Cluster Charge vs. Incident Angle: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
146  h_clChargeVsIncidentAngle[i][j][k] = createHistogram2D(NameOfHisto, TitleOfHisto, 90, 0., 90., "Incident Angle [deg]", m_nBins,
147  m_minCharge, m_maxCharge, "Cluster Charge [e]",
149 
150  NameOfHisto = "DSSDClusterChargeVsSNR_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
151  TitleOfHisto = "DSSD Cluster Charge vs. SNR: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
152  h_clChargeVsSNR[i][j][k] = createHistogram2D(NameOfHisto, TitleOfHisto, 100, 0., 100., "Cluster SNR", m_nBins, m_minCharge,
153  m_maxCharge, "Cluster Charge [e]", m_histoList_clChargeVsSNR[i]);
154 
155  // cluster SNR
156  NameOfHisto = "DSSDClusterSNR_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
157  TitleOfHisto = "DSSD Cluster SNR: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
158  h_clSNR[i][j][k] = createHistogram1D(NameOfHisto, TitleOfHisto, 100, 0., 100., "Cluster SNR", "Entries",
159  m_histoList_clSNR[i]);
160 
161  NameOfHisto = "DSSDClusterSNRVsMomentum_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
162  TitleOfHisto = "DSSD Cluster SNR vs. Track Momentum: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
163  h_clSNRVsMomentum[i][j][k] = createHistogram2D(NameOfHisto, TitleOfHisto, 500, 0, 5., "Track Momentum[GeV/c]",
164  100, 0., 100., "Cluster SNR", m_histoList_clSNRVsMomentum[i]);
165 
166  NameOfHisto = "DSSDClusterSNRVsIncAngle_" + nameSensorType + "Side_" + nameSide + "clSize_" + nameclSize;
167  TitleOfHisto = "DSSD Cluster SNR vs. Incident Angle: " + nameSensorType + ", Side " + nameSide + ", clSize " + nameclSize;
168  h_clSNRVsIncidentAngle[i][j][k] = createHistogram2D(NameOfHisto, TitleOfHisto, 90, 0., 90., "Incident Angle [deg]",
169  100, 0., 100., "Cluster SNR", m_histoList_clSNRVsIncidentAngle[i]);
170  } //cl sizes
171  } //sides
172  } //sensor types
173 } // initialize
174 
176 {
177  if (m_useTrackInfo == false) {
178  B2ERROR(" you must use track information! change the module parameter to TRUE");
179  return;
180  }
181 
182  BOOST_FOREACH(Track & track, m_Tracks) {
183 
184  h_nTracks->Fill(m_Tracks.getEntries());
185  // Obtaining track momentum, P value & SVD hits, track hypothesis made for pions(or electrons in case of TB)
186  const TrackFitResult* tfr = NULL;
187  tfr = track.getTrackFitResult(Const::pion);
188 
189  if (tfr) {
190  h_TracksPvalue->Fill(tfr->getPValue());
191  h_TracksMomentum->Fill(tfr->getMomentum().Mag());
192  h_TracksnSVDhits->Fill((tfr->getHitPatternVXD()).getNSVDHits());
193  } // trf
194 
195  RelationVector<RecoTrack> theRC = DataStore::getRelationsWithObj<RecoTrack>(&track);
196  RelationVector<SVDCluster> svdClustersTrack = DataStore::getRelationsWithObj<SVDCluster>(theRC[0]);
197  RelationVector<VXDDedxTrack> VXDdedxTr = DataStore::getRelationsWithObj<VXDDedxTrack>(&track);
198  // if (theRC.size() == 0 || svdClustersTrack.size() == 0 || VXDdedxTr.size() == 0) continue;
199 
200  // double trkMom = VXDdedxTr[0]->getMomentum();
201 
202  for (int cl = 0; cl < static_cast<int>(svdClustersTrack.size()); cl++) {
203  VxdID cl_VxdID = svdClustersTrack[cl]->getSensorID();
204  int layer = cl_VxdID.getLayerNumber() - 3;
205  // int ladder = cl_VxdID.getLadderNumber() - 1;
206  int sensor = cl_VxdID.getSensorNumber() - 1;
207  int side = static_cast<int>(svdClustersTrack[cl]->isUCluster());
208  float clCharge = svdClustersTrack[cl]->getCharge();
209  float clSNR = svdClustersTrack[cl]->getSNR();
210  unsigned short clSize = svdClustersTrack[cl]->getSize();
211  unsigned short clSizeIndex;
212  unsigned short sensTypeIndex;
213  VXD::SensorInfoBase info = VXD::GeoCache::get(cl_VxdID);
214  // double sensThickness = info.getThickness();
215  // double pathLength = VXDdedxTr[0]->getDx(cl); // path length in the silicon, for a given sensor
216  // double incAngle = 180. / M_PI * acos(sensThickness / pathLength);
217 
218  // clSizeIndex = 0,1 for clSize 1,2 respectively and =2 for clSize >=3
219  clSizeIndex = clSize >= 3 ? 2 : clSize - 1;
220 
221  if (layer == 0) { // L3 small rectangular
222  sensTypeIndex = 0;
223  } else if ((sensor != 0) && (sensor != (m_nSensorsOnLayer[layer] - 1))) { // L456 large rectangular (origami)
224  sensTypeIndex = 1;
225  } else if (sensor == m_nSensorsOnLayer[layer] - 1) { // L456 large rectangular (BWD)
226  sensTypeIndex = 2;
227  } else { // L456 trapezoidal (FWD)
228  sensTypeIndex = 3;
229  }
230 
231  // cl charge histograms
232  h_clCharge[sensTypeIndex][side][clSizeIndex]->Fill(clCharge);
233  // h_clChargeVsMomentum[sensTypeIndex][side][clSizeIndex]->Fill(trkMom, clCharge);
234  // h_clChargeVsIncidentAngle[sensTypeIndex][side][clSizeIndex]->Fill(incAngle, clCharge);
235  h_clChargeVsSNR[sensTypeIndex][side][clSizeIndex]->Fill(clSNR, clCharge);
236 
237  //cl size
238  h_clSize[sensTypeIndex][side]->Fill(clSize);
239  // h_clSizeVsMomentum[sensTypeIndex][side]->Fill(trkMom, clSize);
240  // h_clSizeVsIncidentAngle[sensTypeIndex][side]->Fill(incAngle, clSize);
241  h_clSizeVsSNR[sensTypeIndex][side]->Fill(clSNR, clSize);
242 
243  //cl SNR
244  h_clSNR[sensTypeIndex][side][clSizeIndex]->Fill(clSNR);
245  // h_clSNRVsMomentum[sensTypeIndex][side][clSizeIndex]->Fill(trkMom, clSize);
246  // h_clSNRVsIncidentAngle[sensTypeIndex][side][clSizeIndex]->Fill(incAngle, clSize);
247 
248  } //clusters
249  } //tracks
250 
251 } //event()
252 
253 
255 {
256  gROOT->SetBatch(1);
257  // comparison plots
258  TCanvas* c;
259  for (int i = 0; i < m_nSensorTypes; i++) {
260  for (int j = 0; j < m_nSides; j++) {
261  std::string nameSensorType = "";
262  nameSensorType += m_nameSensorTypes[i];
263  std::string nameSide = "";
264  nameSide = (j > 0 ? "U" : "V");
265 
266  c = comparisonPlot(h_clCharge[i][j][0], h_clCharge[i][j][1], h_clCharge[i][j][2]);
267  c->SaveAs((m_outputDirName + "/clChargeComparison_" + nameSensorType + nameSide + ".png").c_str());
268  if (c) {
269  c->Close();
270  gSystem->ProcessEvents();
271  }
272 
273  c = comparisonPlot(h_clSNR[i][j][0], h_clSNR[i][j][1], h_clSNR[i][j][2]);
274  c->SaveAs((m_outputDirName + "/clSNRComparison_" + nameSensorType + nameSide + ".png").c_str());
275  if (c) {
276  c->Close();
277  gSystem->ProcessEvents();
278  }
279  }
280  }
281  // save to .root file
282  if (m_outputRootFile != NULL) {
283  m_outputRootFile->cd();
284  TDirectory* oldDir = gDirectory;
285 
286  TDirectory* dir_clCharge = oldDir->mkdir("clCharge");
287  TDirectory* dir_clChargeVsSNR = oldDir->mkdir("clChargeVsSNR");
288  TDirectory* dir_clSize = oldDir->mkdir("clSize");
289  TDirectory* dir_clSNR = oldDir->mkdir("clSNR");
290 
291  for (int i = 0; i < m_nSensorTypes; i++) {
292  // cluster charge
293  dir_clCharge->cd();
294  TDirectory* dir_clChargeSt = dir_clCharge->mkdir(m_nameSensorTypes[i].c_str());
295  dir_clChargeSt->cd();
296  TIter nextH_clCharge(m_histoList_clCharge[i]);
297  TObject* obj;
298 
299  while ((obj = dynamic_cast<TH1F*>(nextH_clCharge()))) {
300  obj->Write();
301  }
302 
303  // cluster charge vs. SNR
304  dir_clChargeVsSNR->cd();
305  TDirectory* dir_clChargeVsSNRSt = dir_clChargeVsSNR->mkdir(m_nameSensorTypes[i].c_str());
306  dir_clChargeVsSNRSt->cd();
307  TIter nextH_clChargeVsSNR(m_histoList_clChargeVsSNR[i]);
308  while ((obj = dynamic_cast<TH2F*>(nextH_clChargeVsSNR()))) {
309  obj->Write();
310  }
311 
312  // cluster size
313  dir_clSize->cd();
314  TDirectory* dir_clSizeSt = dir_clSize->mkdir(m_nameSensorTypes[i].c_str());
315  dir_clSizeSt->cd();
316  TIter nextH_clSize(m_histoList_clSize[i]);
317  while ((obj = dynamic_cast<TH1F*>(nextH_clSize()))) {
318  obj->Write();
319  }
320 
321  // cluster SNR
322  dir_clSNR->cd();
323  TDirectory* dir_clSNRSt = dir_clSNR->mkdir(m_nameSensorTypes[i].c_str());
324  dir_clSNRSt->cd();
325  TIter nextH_clSNR(m_histoList_clSNR[i]);
326  while ((obj = dynamic_cast<TH1F*>(nextH_clSNR()))) {
327  obj->Write();
328  }
329  }
330  m_outputRootFile->Close();
331  }
332 } //terminate
333 
334 TH1F* SVDChargeSharingAnalysisModule::createHistogram1D(const char* name, const char* title,
335  Int_t nbins, Double_t min, Double_t max,
336  const char* xtitle, const char* ytitle, TList* histoList)
337 {
338  TH1F* h = new TH1F(name, title, nbins, min, max);
339  h->GetXaxis()->SetTitle(xtitle);
340  h->GetYaxis()->SetTitle(ytitle);
341  if (histoList) {
342  histoList->Add(h);
343  }
344  return h;
345 }
346 
347 TH2F* SVDChargeSharingAnalysisModule::createHistogram2D(const char* name, const char* title,
348  Int_t nbinsX, Double_t minX, Double_t maxX,
349  const char* titleX,
350  Int_t nbinsY, Double_t minY, Double_t maxY,
351  const char* titleY, TList* histoList)
352 {
353 
354  TH2F* h = new TH2F(name, title, nbinsX, minX, maxX, nbinsY, minY, maxY);
355 
356  h->GetXaxis()->SetTitle(titleX);
357  h->GetYaxis()->SetTitle(titleY);
358 
359  if (histoList)
360  histoList->Add(h);
361 
362  return h;
363 }
364 
365 TCanvas* SVDChargeSharingAnalysisModule::comparisonPlot(TH1F* h1, TH1F* h2, TH1F* h3)
366 {
367  TCanvas* c = new TCanvas("c", "c");
368  c->cd();
369 
370  gStyle->SetOptStat(0);
371  gStyle->SetOptTitle(0);
372  gStyle->SetCanvasPreferGL(1);
373 
374  h1->SetLineColor(kRed);
375  h1->SetFillColorAlpha(kRed, 0.35);
376 
377  h2->SetLineColor(kBlue);
378  h2->SetFillColorAlpha(kBlue, 0.35);
379 
380  h3->SetLineColor(kGreen);
381  h3->SetFillColorAlpha(kGreen, 0.35);
382 
383  h2->GetYaxis()->SetTitleOffset(1.4);
384  h2->Draw();
385  h1->Draw("][sames");
386  h3->Draw("][sames");
387 
388  auto legend = new TLegend(0.7, 0.7, 0.9, 0.9);
389  legend->SetTextAlign(22);
390  legend->AddEntry(h1, "clSize 1", "f");
391  legend->AddEntry(h2, "clSize 2", "f");
392  legend->AddEntry(h3, "clSize >= 3", "f");
393  legend->Draw();
394 
395  return c;
396 }
Belle2::RelationVector::size
size_t size() const
Get number of relations.
Definition: RelationVector.h:98
Belle2::SVDChargeSharingAnalysisModule::m_outputRootFileName
std::string m_outputRootFileName
root output file name.
Definition: SVDChargeSharingAnalysisModule.h:67
Belle2::SVDChargeSharingAnalysisModule::m_outputRootFile
TFile * m_outputRootFile
root ouput file pointer.
Definition: SVDChargeSharingAnalysisModule.h:65
Belle2::SVDChargeSharingAnalysisModule::h_clChargeVsMomentum
TH2F * h_clChargeVsMomentum[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster charge vs.
Definition: SVDChargeSharingAnalysisModule.h:101
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSNR
TList * m_histoList_clSNR[m_nSensorTypes]
SVDCluster SNR histogram list.
Definition: SVDChargeSharingAnalysisModule.h:130
Belle2::SVDChargeSharingAnalysisModule::m_nClSizes
static const int m_nClSizes
Distinction between clSize=1, clSize=2 & clSize>=3.
Definition: SVDChargeSharingAnalysisModule.h:81
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::TrackFitResult::getMomentum
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
Definition: TrackFitResult.h:116
Belle2::VXD::GeoCache::get
static const SensorInfoBase & get(Belle2::VxdID id)
Return a reference to the SensorInfo of a given SensorID.
Definition: GeoCache.h:141
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::SVDChargeSharingAnalysisModule::h_clCharge
TH1F * h_clCharge[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster charge.
Definition: SVDChargeSharingAnalysisModule.h:100
Belle2::SVDChargeSharingAnalysisModule::h_TracksnSVDhits
TH1F * h_TracksnSVDhits
Number of SVDhits for a track.
Definition: SVDChargeSharingAnalysisModule.h:97
Belle2::TrackFitResult::getPValue
double getPValue() const
Getter for Chi2 Probability of the track fit.
Definition: TrackFitResult.h:163
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSize
TList * m_histoList_clSize[m_nSensorTypes]
SVDCluster size histogram list.
Definition: SVDChargeSharingAnalysisModule.h:125
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clChargeVsMomentum
TList * m_histoList_clChargeVsMomentum[m_nSensorTypes]
SVDCluster charge vs.
Definition: SVDChargeSharingAnalysisModule.h:121
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSNRVsMomentum
TList * m_histoList_clSNRVsMomentum[m_nSensorTypes]
SVDCluster SNR vs.
Definition: SVDChargeSharingAnalysisModule.h:131
Belle2::SVDChargeSharingAnalysisModule::h_clSNRVsIncidentAngle
TH2F * h_clSNRVsIncidentAngle[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster SNR vs.
Definition: SVDChargeSharingAnalysisModule.h:114
Belle2::VXD::SensorInfoBase
Base class to provide Sensor Information for PXD and SVD.
Definition: SensorInfoBase.h:40
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clChargeVsSNR
TList * m_histoList_clChargeVsSNR[m_nSensorTypes]
SVDCluster charge vs.
Definition: SVDChargeSharingAnalysisModule.h:123
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSizeVsIncidentAngle
TList * m_histoList_clSizeVsIncidentAngle[m_nSensorTypes]
SVDCluster size vs.
Definition: SVDChargeSharingAnalysisModule.h:127
Belle2::SVDChargeSharingAnalysisModule::h_clSize
TH1F * h_clSize[m_nSensorTypes][m_nSides]
SVDCluster size.
Definition: SVDChargeSharingAnalysisModule.h:106
Belle2::SVDChargeSharingAnalysisModule::m_Tracks
StoreArray< Track > m_Tracks
SVD Tracks strore array.
Definition: SVDChargeSharingAnalysisModule.h:74
Belle2::SVDChargeSharingAnalysisModule::h_clSizeVsIncidentAngle
TH2F * h_clSizeVsIncidentAngle[m_nSensorTypes][m_nSides]
SVDCluster size vs.
Definition: SVDChargeSharingAnalysisModule.h:108
Belle2::SVDChargeSharingAnalysisModule
Module for monitoring DSSD cluster charge deposition in regard of capacitive charge sharing between a...
Definition: SVDChargeSharingAnalysisModule.h:42
Belle2::TrackFitResult
Values of the result of a track fit with a given particle hypothesis.
Definition: TrackFitResult.h:59
Belle2::SVDChargeSharingAnalysisModule::m_nSensorTypes
static const int m_nSensorTypes
L3 small rect,.
Definition: SVDChargeSharingAnalysisModule.h:82
Belle2::SVDChargeSharingAnalysisModule::event
virtual void event() override
Core method, called for each processed event.
Definition: SVDChargeSharingAnalysisModule.cc:175
Belle2::SVDChargeSharingAnalysisModule::createHistogram1D
TH1F * createHistogram1D(const char *name, const char *title, Int_t nbins, Double_t min, Double_t max, const char *xtitle, const char *ytitle, TList *histoList=NULL)
1D histogram creator.
Definition: SVDChargeSharingAnalysisModule.cc:334
Belle2::SVDChargeSharingAnalysisModule::terminate
virtual void terminate() override
Summary method, called after processing of all events.
Definition: SVDChargeSharingAnalysisModule.cc:254
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::SVDChargeSharingAnalysisModule::m_outputDirName
std::string m_outputDirName
output directory.
Definition: SVDChargeSharingAnalysisModule.h:66
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
Belle2::SVDChargeSharingAnalysisModule::h_clChargeVsIncidentAngle
TH2F * h_clChargeVsIncidentAngle[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster charge vs.
Definition: SVDChargeSharingAnalysisModule.h:102
Belle2::SVDChargeSharingAnalysisModule::~SVDChargeSharingAnalysisModule
virtual ~SVDChargeSharingAnalysisModule() override
Default destructor.
Definition: SVDChargeSharingAnalysisModule.cc:60
Belle2::RelationVector
Class for type safe access to objects that are referred to in relations.
Definition: DataStore.h:38
Belle2::SVDChargeSharingAnalysisModule::m_minCharge
const Double_t m_minCharge
Minimum charge in histogram.
Definition: SVDChargeSharingAnalysisModule.h:90
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::SVDChargeSharingAnalysisModule::h_clChargeVsSNR
TH2F * h_clChargeVsSNR[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster charge vs.
Definition: SVDChargeSharingAnalysisModule.h:103
Belle2::SVDChargeSharingAnalysisModule::m_nSides
static const int m_nSides
DSSD sides.
Definition: SVDChargeSharingAnalysisModule.h:80
Belle2::SVDChargeSharingAnalysisModule::createHistogram2D
TH2F * createHistogram2D(const char *name, const char *title, Int_t nbinsX, Double_t minX, Double_t maxX, const char *titleX, Int_t nbinsY, Double_t minY, Double_t maxY, const char *titleY, TList *histoList)
2D histogram creator.
Definition: SVDChargeSharingAnalysisModule.cc:347
Belle2::SVDChargeSharingAnalysisModule::m_nSensorsOnLayer
int m_nSensorsOnLayer[m_nLayers]
Total number of DSSDs on L3,4,5,6.
Definition: SVDChargeSharingAnalysisModule.h:85
Belle2::SVDChargeSharingAnalysisModule::h_TracksMomentum
TH1F * h_TracksMomentum
Tracks momentum.
Definition: SVDChargeSharingAnalysisModule.h:96
Belle2::SVDChargeSharingAnalysisModule::h_clSizeVsSNR
TH2F * h_clSizeVsSNR[m_nSensorTypes][m_nSides]
SVDCluster size vs.
Definition: SVDChargeSharingAnalysisModule.h:109
Belle2::VxdID::getSensorNumber
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:110
Belle2::SVDChargeSharingAnalysisModule::h_nTracks
TH1F * h_nTracks
Number of tracks.
Definition: SVDChargeSharingAnalysisModule.h:94
Belle2::SVDChargeSharingAnalysisModule::m_useTrackInfo
bool m_useTrackInfo
True if using clusters related to tracks.
Definition: SVDChargeSharingAnalysisModule.h:70
Belle2::SVDChargeSharingAnalysisModule::m_maxCharge
const Double_t m_maxCharge
1ADU bin width.
Definition: SVDChargeSharingAnalysisModule.h:91
Belle2::TrackFitResult::getHitPatternVXD
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
Definition: TrackFitResult.cc:125
Belle2::SVDChargeSharingAnalysisModule::m_svdClusters
StoreArray< SVDCluster > m_svdClusters
SVD clusters store array.
Definition: SVDChargeSharingAnalysisModule.h:73
Belle2::SVDChargeSharingAnalysisModule::h_TracksPvalue
TH1F * h_TracksPvalue
Tracks P value.
Definition: SVDChargeSharingAnalysisModule.h:95
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clCharge
TList * m_histoList_clCharge[m_nSensorTypes]
SVDCluster charge histogram list.
Definition: SVDChargeSharingAnalysisModule.h:120
Belle2::SVDChargeSharingAnalysisModule::h_clSizeVsMomentum
TH2F * h_clSizeVsMomentum[m_nSensorTypes][m_nSides]
SVDCluster size vs.
Definition: SVDChargeSharingAnalysisModule.h:107
Belle2::SVDChargeSharingAnalysisModule::initialize
virtual void initialize() override
Initialize the SVDChargeSharingAnalysi module.
Definition: SVDChargeSharingAnalysisModule.cc:64
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clChargeVsIncidentAngle
TList * m_histoList_clChargeVsIncidentAngle[m_nSensorTypes]
SVDCluster charge vs.
Definition: SVDChargeSharingAnalysisModule.h:122
Belle2::Track
Class that bundles various TrackFitResults.
Definition: Track.h:35
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSizeVsSNR
TList * m_histoList_clSizeVsSNR[m_nSensorTypes]
SVDCluster size vs.
Definition: SVDChargeSharingAnalysisModule.h:128
Belle2::VxdID::getLayerNumber
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:106
Belle2::SVDChargeSharingAnalysisModule::m_nameSensorTypes
std::string m_nameSensorTypes[m_nSensorTypes]
Sensor type names.
Definition: SVDChargeSharingAnalysisModule.h:84
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSizeVsMomentum
TList * m_histoList_clSizeVsMomentum[m_nSensorTypes]
SVDCluster size vs.
Definition: SVDChargeSharingAnalysisModule.h:126
Belle2::SVDChargeSharingAnalysisModule::h_clSNR
TH1F * h_clSNR[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster SNR.
Definition: SVDChargeSharingAnalysisModule.h:112
Belle2::SVDChargeSharingAnalysisModule::m_histoList_Tracks
TList * m_histoList_Tracks
List of all histograms concerning track information: nTracks, Track P Value, momentum and SVD hits.
Definition: SVDChargeSharingAnalysisModule.h:118
Belle2::SVDChargeSharingAnalysisModule::m_nBins
const Int_t m_nBins
Number of bins.
Definition: SVDChargeSharingAnalysisModule.h:89
Belle2::SVDChargeSharingAnalysisModule::comparisonPlot
TCanvas * comparisonPlot(TH1F *h1, TH1F *h2, TH1F *h3)
Used to compare charge histograms for different cluster sizes.
Definition: SVDChargeSharingAnalysisModule.cc:365
Belle2::SVDChargeSharingAnalysisModule::m_histoList_clSNRVsIncidentAngle
TList * m_histoList_clSNRVsIncidentAngle[m_nSensorTypes]
SVDCluster SNR vs.
Definition: SVDChargeSharingAnalysisModule.h:132
Belle2::SVDChargeSharingAnalysisModule::h_clSNRVsMomentum
TH2F * h_clSNRVsMomentum[m_nSensorTypes][m_nSides][m_nClSizes]
SVDCluster SNR vs.
Definition: SVDChargeSharingAnalysisModule.h:113