Belle II Software  release-05-02-19
EKLMAlignmentAlongStripsAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2015 Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Kirill Chilikin *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/eklm/calibration/EKLMAlignmentAlongStripsAlgorithm.h>
13 
14 /* KLM headers. */
15 #include <klm/eklm/geometry/GeometryData.h>
16 
17 /* ROOT headers. */
18 #include <TFile.h>
19 #include <TTree.h>
20 
21 using namespace Belle2;
22 
23 static bool compareSegmentSignificance(
24  const std::pair<int, double>& first, const std::pair<int, double>& second)
25 {
26  return first.second > second.second;
27 }
28 
30  CalibrationAlgorithm("EKLMAlignmentAlongStripsCollector"),
31  m_OutputFile("")
32 {
33 }
34 
36 {
37 }
38 
40 {
41  /* cppcheck-suppress variableScope */
42  int i, n, nSections, nSectors, nLayers, nDetectorLayers, nPlanes, nSegments;
43  /* cppcheck-suppress variableScope */
44  int segment, segmentGlobal;
45  std::vector<std::pair<int, double> > segmentSignificance;
46  std::vector<std::pair<int, double> >::iterator it;
47  int iSection, iLayer, iSector, iPlane, iSegment;
48  int**** *nHits, *** *averageHits;
49  double nHitsSegment, nHitsAverage, nSigma;
50  bool found;
51  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
53  struct Event* event = nullptr;
54  std::shared_ptr<TTree> t_in = getObjectPtr<TTree>("calibration_data");
55  t_in->SetBranchAddress("event", &event);
56  nSections = geoDat->getNSections();
57  nSectors = geoDat->getNSectors();
58  nLayers = geoDat->getNLayers();
59  nPlanes = geoDat->getNPlanes();
60  nSegments = geoDat->getNSegments();
61  nHits = new int**** [nSections];
62  averageHits = new int** *[nSections];
63  for (iSection = 0; iSection < nSections; iSection++) {
64  nHits[iSection] = new int** *[nLayers];
65  averageHits[iSection] = new int** [nLayers];
66  for (iLayer = 0; iLayer < nLayers; iLayer++) {
67  nHits[iSection][iLayer] = new int** [nSectors];
68  averageHits[iSection][iLayer] = new int* [nPlanes];
69  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
70  averageHits[iSection][iLayer][iPlane] = new int[nSegments];
71  for (iSegment = 0; iSegment < nSegments; iSegment++) {
72  averageHits[iSection][iLayer][iPlane][iSegment] = 0;
73  }
74  }
75  for (iSector = 0; iSector < nSectors; iSector++) {
76  nHits[iSection][iLayer][iSector] = new int* [nPlanes];
77  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
78  nHits[iSection][iLayer][iSector][iPlane] = new int[nSegments];
79  for (iSegment = 0; iSegment < nSegments; iSegment++) {
80  nHits[iSection][iLayer][iSector][iPlane][iSegment] = 0;
81  }
82  }
83  }
84  }
85  }
86  n = t_in->GetEntries();
87  for (i = 0; i < n; i++) {
88  t_in->GetEntry(i);
89  segment = (event->strip - 1) / elementNumbers->getNStripsSegment();
90  nHits[event->section - 1][event->layer - 1][event->sector - 1]
91  [event->plane - 1][segment]++;
92  averageHits[event->section - 1][event->layer - 1]
93  [getAveragedPlane(event->sector, event->plane)][segment]++;
94  }
95  for (iSection = 0; iSection < nSections; iSection++) {
96  for (iLayer = 0; iLayer < nLayers; iLayer++) {
97  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
98  for (iSegment = 0; iSegment < nSegments; iSegment++) {
99  averageHits[iSection][iLayer][iPlane][iSegment] /= nSectors;
100  }
101  }
102  }
103  }
104  for (iSection = 0; iSection < nSections; iSection++) {
105  nDetectorLayers = geoDat->getNDetectorLayers(iSection + 1);
106  for (iLayer = 0; iLayer < nDetectorLayers; iLayer++) {
107  for (iSector = 0; iSector < nSectors; iSector++) {
108  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
109  for (iSegment = 0; iSegment < nSegments; iSegment++) {
110  segmentGlobal = elementNumbers->segmentNumber(
111  iSection + 1, iLayer + 1, iSector + 1, iPlane + 1,
112  iSegment + 1);
113  nHitsSegment = nHits[iSection][iLayer][iSector][iPlane][iSegment];
114  nHitsAverage = averageHits[iSection][iLayer]
115  [getAveragedPlane(iSector + 1, iPlane + 1)]
116  [iSegment];
117  nSigma = (nHitsSegment - nHitsAverage) /
118  sqrt(nHitsSegment + nHitsAverage -
119  2.0 / nSectors * sqrt(nHitsSegment * nHitsAverage));
120  segmentSignificance.push_back(
121  std::pair<int, double>(segmentGlobal, nSigma));
122  }
123  }
124  }
125  }
126  }
127  sort(segmentSignificance.begin(), segmentSignificance.end(),
128  compareSegmentSignificance);
129  printf("Checking for significantly shifted segments:\n");
130  found = false;
131  for (it = segmentSignificance.begin(); it != segmentSignificance.end();
132  ++it) {
133  if (it->second < 3.0)
134  break;
135  elementNumbers->segmentNumberToElementNumbers(
136  it->first, &iSection, &iLayer, &iSector, &iPlane, &iSegment);
137  printf("Segment %d (section %d, layer %d, sector %d, plane %d, segment %d):"
138  " %.1f sigma\n", it->first, iSection, iLayer, iSector, iPlane,
139  iSegment, it->second);
140  found = true;
141  }
142  if (!found)
143  printf("none found.\n");
144  for (iSection = 0; iSection < nSections; iSection++) {
145  for (iLayer = 0; iLayer < nLayers; iLayer++) {
146  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
147  delete[] averageHits[iSection][iLayer][iPlane];
148  }
149  for (iSector = 0; iSector < nSectors; iSector++) {
150  for (iPlane = 0; iPlane < nPlanes; iPlane++) {
151  delete[] nHits[iSection][iLayer][iSector][iPlane];
152  }
153  delete[] nHits[iSection][iLayer][iSector];
154  }
155  delete[] nHits[iSection][iLayer];
156  delete[] averageHits[iSection][iLayer];
157  }
158  delete[] nHits[iSection];
159  delete[] averageHits[iSection];
160  }
161  delete[] nHits;
162  delete[] averageHits;
163  if (m_OutputFile != "") {
164  TFile* f_out = new TFile(m_OutputFile.c_str(), "recreate");
165  TTree* t_out = new TTree("tree", "");
166  t_out->Branch("event", event);
167  n = t_in->GetEntries();
168  for (i = 0; i < n; i++) {
169  t_in->GetEntry(i);
170  t_out->Fill();
171  }
172  f_out->cd();
173  t_out->Write();
174  delete t_out;
175  delete f_out;
176  }
178 }
179 
181 {
182  m_OutputFile = outputFile;
183 }
184 
186 getAveragedPlane(int sector, int plane) const
187 {
188  if (sector == 1 || sector == 4)
189  return plane - 1;
190  return 2 - plane;
191 }
192 
Belle2::EKLMAlignmentAlongStripsAlgorithm::getAveragedPlane
int getAveragedPlane(int sector, int plane) const
Get plane number for average values (0-based).
Definition: EKLMAlignmentAlongStripsAlgorithm.cc:186
Belle2::EKLMElementNumbers::segmentNumberToElementNumbers
void segmentNumberToElementNumbers(int segmentGlobal, int *section, int *layer, int *sector, int *plane, int *segment) const
Get element numbers by segment global number.
Definition: EKLMElementNumbers.cc:199
Belle2::EKLMGeometry::getNSegments
int getNSegments() const
Get number of segments.
Definition: EKLMGeometry.h:1725
Belle2::EKLMGeometry::getNLayers
int getNLayers() const
Get number of layers.
Definition: EKLMGeometry.h:1695
Belle2::EKLMAlignmentAlongStripsAlgorithm::Event
Event: time, distance from hit to SiPM.
Definition: EKLMAlignmentAlongStripsAlgorithm.h:40
Belle2::EKLMElementNumbers
EKLM element numbers.
Definition: EKLMElementNumbers.h:34
Belle2::EKLMGeometry::getNPlanes
int getNPlanes() const
Get number of planes.
Definition: EKLMGeometry.h:1717
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::EKLMAlignmentAlongStripsAlgorithm::~EKLMAlignmentAlongStripsAlgorithm
~EKLMAlignmentAlongStripsAlgorithm()
Destructor.
Definition: EKLMAlignmentAlongStripsAlgorithm.cc:35
Belle2::EKLMAlignmentAlongStripsAlgorithm::setOutputFile
void setOutputFile(const char *outputFile)
Set output file name.
Definition: EKLMAlignmentAlongStripsAlgorithm.cc:180
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::EKLMElementNumbers::segmentNumber
int segmentNumber(int section, int layer, int sector, int plane, int segment) const
Get segment number.
Definition: EKLMElementNumbers.cc:191
Belle2::EKLMElementNumbers::Instance
static const EKLMElementNumbers & Instance()
Instantiation.
Definition: EKLMElementNumbers.cc:20
Belle2::EKLMAlignmentAlongStripsAlgorithm::calibrate
CalibrationAlgorithm::EResult calibrate() override
Calibration.
Definition: EKLMAlignmentAlongStripsAlgorithm.cc:39
Belle2::EKLMGeometry::getNDetectorLayers
int getNDetectorLayers(int section) const
Get number of detector layers.
Definition: EKLMGeometry.cc:295
Belle2::EKLM::GeometryData
EKLM geometry data.
Definition: GeometryData.h:40
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::EKLMElementNumbers::getNStripsSegment
static constexpr int getNStripsSegment()
Get number of strips in a segment.
Definition: EKLMElementNumbers.h:401
Belle2::EKLMAlignmentAlongStripsAlgorithm::m_OutputFile
std::string m_OutputFile
Output file name.
Definition: EKLMAlignmentAlongStripsAlgorithm.h:85
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::EKLMGeometry::getNSections
int getNSections() const
Get number of sections.
Definition: EKLMGeometry.h:1687
Belle2::EKLMGeometry::getNSectors
int getNSectors() const
Get number of sectors.
Definition: EKLMGeometry.h:1709
Belle2::EKLMAlignmentAlongStripsAlgorithm::EKLMAlignmentAlongStripsAlgorithm
EKLMAlignmentAlongStripsAlgorithm()
Constructor.
Definition: EKLMAlignmentAlongStripsAlgorithm.cc:29
Belle2::EKLM::GeometryData::Instance
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:35