Belle II Software  release-05-01-25
KLMDisplacementGenerator.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2012 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/calibration/KLMDisplacementGenerator.h>
13 
14 /* KLM headers. */
15 #include <klm/eklm/geometry/AlignmentChecker.h>
16 #include <klm/eklm/geometry/GeometryData.h>
17 
18 /* Belle 2 headers. */
19 #include <framework/gearbox/Unit.h>
20 
21 /* ROOT headers. */
22 #include <TFile.h>
23 #include <TRandom.h>
24 #include <TTree.h>
25 
26 using namespace Belle2;
27 
29  m_GeoDat(&(EKLM::GeometryData::Instance())),
30  m_ElementNumbers(&(KLMElementNumbers::Instance())),
31  m_eklmElementNumbers(&(EKLMElementNumbers::Instance()))
32 {
33 }
34 
36 {
37 }
38 
40  EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment)
41 {
42  KLMAlignmentData alignmentData(0, 0, 0, 0, 0, 0);
43  int iSection, iLayer, iSector, iPlane, iSegment, segment;
44  for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
45  for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
46  iLayer++) {
47  for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
48  int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
49  alignment->setModuleAlignment(module, &alignmentData);
50  for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
51  for (iSegment = 1; iSegment <= m_GeoDat->getNSegments(); iSegment++) {
53  iSection, iLayer, iSector, iPlane, iSegment);
54  segmentAlignment->setSegmentAlignment(segment, &alignmentData);
55  }
56  }
57  }
58  }
59  }
60 }
61 
63  double deltaU, double deltaV, double deltaGamma)
64 {
65  EKLMAlignment alignment;
66  EKLMSegmentAlignment segmentAlignment;
67  KLMAlignmentData moduleAlignment(deltaU, deltaV, 0, 0, 0, deltaGamma);
68  int iSection, iLayer, iSector;
69  fillZeroDisplacements(&alignment, &segmentAlignment);
70  for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
71  for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
72  iLayer++) {
73  for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
74  int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
75  alignment.setModuleAlignment(module, &moduleAlignment);
76  }
77  }
78  }
79 }
80 
82  EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment,
83  bool displaceModule, bool displaceSegment, bool moduleSameDisplacement,
84  bool moduleZeroDeltaU, bool moduleZeroDeltaV, bool moduleZeroDeltaGamma)
85 {
86  const double moduleMaxDeltaU = 5. * Unit::cm;
87  const double moduleMinDeltaU = -5. * Unit::cm;
88  const double moduleMaxDeltaV = 5. * Unit::cm;
89  const double moduleMinDeltaV = -5. * Unit::cm;
90  const double moduleMaxDeltaGamma = 0.02 * Unit::rad;
91  const double moduleMinDeltaGamma = -0.02 * Unit::rad;
92  const double segmentMaxDeltaU = 1. * Unit::cm;
93  const double segmentMinDeltaU = -1. * Unit::cm;
94  const double segmentMaxDeltaV = 0.2 * Unit::cm;
95  const double segmentMinDeltaV = -0.2 * Unit::cm;
96  const double segmentMaxDeltaGamma = 0.003 * Unit::rad;
97  const double segmentMinDeltaGamma = -0.003 * Unit::rad;
98  KLMAlignmentData moduleAlignmentData, segmentAlignmentData, *alignmentData;
99  EKLM::AlignmentChecker alignmentChecker(false);
100  int iSection, iLayer, iSector, iPlane, iSegment, segment;
101  double d;
102  fillZeroDisplacements(alignment, segmentAlignment);
103  for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
104  for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
105  iLayer++) {
106  for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
107  if (moduleSameDisplacement) {
108  if (iSection > 1 || iLayer > 1 || iSector > 1) {
109  int module = m_ElementNumbers->moduleNumberEKLM(1, 1, 1);
110  alignmentData = const_cast<KLMAlignmentData*>(
111  alignment->getModuleAlignment(module));
112  module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
113  alignment->setModuleAlignment(module, alignmentData);
114  for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
115  for (iSegment = 1; iSegment <= m_GeoDat->getNSegments();
116  iSegment++) {
117  segment = m_eklmElementNumbers->segmentNumber(1, 1, 1, iPlane, iSegment);
118  alignmentData = const_cast<KLMAlignmentData*>(
119  segmentAlignment->getSegmentAlignment(segment));
121  iSection, iLayer, iSector, iPlane, iSegment);
122  segmentAlignment->setSegmentAlignment(segment, alignmentData);
123  }
124  }
125  continue;
126  }
127  }
128 module:
129  if (displaceModule) {
130  do {
131  if (moduleZeroDeltaU)
132  d = 0;
133  else
134  d = gRandom->Uniform(moduleMinDeltaU, moduleMaxDeltaU);
135  moduleAlignmentData.setDeltaU(d);
136  if (moduleZeroDeltaV)
137  d = 0;
138  else
139  d = gRandom->Uniform(moduleMinDeltaV, moduleMaxDeltaV);
140  moduleAlignmentData.setDeltaV(d);
141  if (moduleZeroDeltaGamma)
142  d = 0;
143  else
144  d = gRandom->Uniform(moduleMinDeltaGamma, moduleMaxDeltaGamma);
145  moduleAlignmentData.setDeltaGamma(d);
146  } while (!alignmentChecker.checkSectorAlignment(
147  iSection, iLayer, iSector, &moduleAlignmentData));
148  int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
149  alignment->setModuleAlignment(module, &moduleAlignmentData);
150  }
151  if (displaceSegment) {
152  for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
153  for (iSegment = 1; iSegment <= m_GeoDat->getNSegments();
154  iSegment++) {
155  do {
156  segmentAlignmentData.setDeltaU(
157  gRandom->Uniform(segmentMinDeltaU, segmentMaxDeltaU));
158  segmentAlignmentData.setDeltaV(
159  gRandom->Uniform(segmentMinDeltaV, segmentMaxDeltaV));
160  segmentAlignmentData.setDeltaGamma(
161  gRandom->Uniform(segmentMinDeltaGamma, segmentMaxDeltaGamma));
162  } while (!alignmentChecker.checkSegmentAlignment(
163  iSection, iLayer, iSector, iPlane, iSegment,
164  &moduleAlignmentData, &segmentAlignmentData, false));
166  iSection, iLayer, iSector, iPlane, iSegment);
167  segmentAlignment->setSegmentAlignment(
168  segment, &segmentAlignmentData);
169  }
170  }
171  } else {
172  for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
173  for (iSegment = 1; iSegment <= m_GeoDat->getNSegments();
174  iSegment++) {
175  if (!alignmentChecker.checkSegmentAlignment(
176  iSection, iLayer, iSector, iPlane, iSegment,
177  &moduleAlignmentData, &segmentAlignmentData, false))
178  goto module;
179  }
180  }
181  }
182  }
183  }
184  }
185 }
186 
188  EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment,
189  const char* inputFile)
190 {
191  /* cppcheck-suppress variableScope */
192  int i, n, iSection, iLayer, iSector, iPlane, iSegment, segment, param;
193  float value;
194  TFile* f;
195  TTree* tEKLMModule, *tEKLMSegment;
196  KLMAlignmentData* alignmentData;
197  fillZeroDisplacements(alignment, segmentAlignment);
198  f = new TFile(inputFile);
199  tEKLMModule = (TTree*)f->Get("eklm_module");
200  tEKLMModule->SetBranchAddress("section", &iSection);
201  tEKLMModule->SetBranchAddress("layer", &iLayer);
202  tEKLMModule->SetBranchAddress("sector", &iSector);
203  tEKLMModule->SetBranchAddress("param", &param);
204  tEKLMModule->SetBranchAddress("value", &value);
205  tEKLMSegment = (TTree*)f->Get("eklm_segment");
206  tEKLMSegment->SetBranchAddress("section", &iSection);
207  tEKLMSegment->SetBranchAddress("layer", &iLayer);
208  tEKLMSegment->SetBranchAddress("sector", &iSector);
209  tEKLMSegment->SetBranchAddress("plane", &iPlane);
210  tEKLMSegment->SetBranchAddress("segment", &iSegment);
211  tEKLMSegment->SetBranchAddress("param", &param);
212  tEKLMSegment->SetBranchAddress("value", &value);
213  n = tEKLMModule->GetEntries();
214  for (i = 0; i < n; i++) {
215  tEKLMModule->GetEntry(i);
216  int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
217  alignmentData = const_cast<KLMAlignmentData*>(
218  alignment->getModuleAlignment(module));
219  switch (param) {
220  case 1:
221  alignmentData->setDeltaU(value);
222  break;
223  case 2:
224  alignmentData->setDeltaV(value);
225  break;
226  case 3:
227  alignmentData->setDeltaGamma(value);
228  break;
229  }
230  }
231  n = tEKLMSegment->GetEntries();
232  for (i = 0; i < n; i++) {
233  tEKLMSegment->GetEntry(i);
235  iSection, iLayer, iSector, iPlane, iSegment);
236  alignmentData = const_cast<KLMAlignmentData*>(
237  segmentAlignment->getSegmentAlignment(segment));
238  switch (param) {
239  case 1:
240  alignmentData->setDeltaV(value);
241  break;
242  case 2:
243  alignmentData->setDeltaGamma(value);
244  break;
245  }
246  }
247 }
248 
250 {
251  const int nPoints = 1000;
252  const float maxDeltaU = 5. * Unit::cm;
253  const float minDeltaU = -5. * Unit::cm;
254  const float maxDeltaV = 5. * Unit::cm;
255  const float minDeltaV = -5. * Unit::cm;
256  const float maxDeltaGamma = 0.02 * Unit::rad;
257  const float minDeltaGamma = -0.02 * Unit::rad;
258  float deltaU, deltaV, deltaGamma;
259  int i, alignmentStatus, iSection, iLayer, iSector;
260  EKLMAlignment alignment;
261  EKLMSegmentAlignment segmentAlignment;
262  KLMAlignmentData alignmentDataRandom;
263  EKLM::AlignmentChecker alignmentChecker(false);
264  TTree* t;
265  f->cd();
266  t = new TTree("module", "");
267  t->Branch("deltaU", &deltaU, "deltaU/F");
268  t->Branch("deltaV", &deltaV, "deltaV/F");
269  t->Branch("deltaGamma", &deltaGamma, "deltaGamma/F");
270  t->Branch("status", &alignmentStatus, "status/I");
271  for (i = 0; i < nPoints; i++) {
272  fillZeroDisplacements(&alignment, &segmentAlignment);
273  deltaU = gRandom->Uniform(minDeltaU, maxDeltaU);
274  deltaV = gRandom->Uniform(minDeltaV, maxDeltaV);
275  deltaGamma = gRandom->Uniform(minDeltaGamma, maxDeltaGamma);
276  alignmentDataRandom.setDeltaU(deltaU);
277  alignmentDataRandom.setDeltaV(deltaV);
278  alignmentDataRandom.setDeltaGamma(deltaGamma);
279  for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
280  for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
281  iLayer++) {
282  for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
283  int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
284  alignment.setModuleAlignment(module, &alignmentDataRandom);
285  }
286  }
287  }
288  if (alignmentChecker.checkAlignment(&alignment, &segmentAlignment))
289  alignmentStatus = 1;
290  else
291  alignmentStatus = 0;
292  t->Fill();
293  }
294  t->Write();
295 }
296 
298 {
299  const int nPoints = 1000;
300  const float maxDeltaU = 9. * Unit::cm;
301  const float minDeltaU = -4. * Unit::cm;
302  const float maxDeltaV = 0.2 * Unit::cm;
303  const float minDeltaV = -0.2 * Unit::cm;
304  const float maxDeltaGamma = 0.003 * Unit::rad;
305  const float minDeltaGamma = -0.003 * Unit::rad;
306  float deltaU, deltaV, deltaGamma;
307  int i, alignmentStatus;
308  int iSection, iLayer, iSector, jPlane, jSegment, segment;
309  EKLMAlignment alignment;
310  EKLMSegmentAlignment segmentAlignment;
311  KLMAlignmentData alignmentDataRandom;
312  EKLM::AlignmentChecker alignmentChecker(false);
313  TTree* t;
314  f->cd();
315  t = new TTree("segment", "");
316  t->Branch("plane", &jPlane, "plane/I");
317  t->Branch("segment", &jSegment, "segment/I");
318  t->Branch("deltaU", &deltaU, "deltaU/F");
319  t->Branch("deltaV", &deltaV, "deltaV/F");
320  t->Branch("deltaGamma", &deltaGamma, "deltaGamma/F");
321  t->Branch("status", &alignmentStatus, "status/I");
322  for (jPlane = 1; jPlane <= m_GeoDat->getNPlanes(); jPlane++) {
323  printf("Plane %d\n", jPlane);
324  for (jSegment = 1; jSegment <= m_GeoDat->getNSegments(); jSegment++) {
325  printf("Segment %d\n", jSegment);
326  for (i = 0; i < nPoints; i++) {
327  fillZeroDisplacements(&alignment, &segmentAlignment);
328  deltaU = gRandom->Uniform(minDeltaU, maxDeltaU);
329  deltaV = gRandom->Uniform(minDeltaV, maxDeltaV);
330  deltaGamma = gRandom->Uniform(minDeltaGamma, maxDeltaGamma);
331  alignmentDataRandom.setDeltaU(deltaU);
332  alignmentDataRandom.setDeltaV(deltaV);
333  alignmentDataRandom.setDeltaGamma(deltaGamma);
334  for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
335  for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
336  iLayer++) {
337  for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
339  iSection, iLayer, iSector, jPlane, jSegment);
340  segmentAlignment.setSegmentAlignment(
341  segment, &alignmentDataRandom);
342  }
343  }
344  }
345  if (alignmentChecker.checkAlignment(&alignment, &segmentAlignment))
346  alignmentStatus = 1;
347  else
348  alignmentStatus = 0;
349  t->Fill();
350  }
351  }
352  }
353  t->Write();
354 }
355 
357 {
358  TFile* f;
359  f = new TFile(outputFile, "recreate");
362  delete f;
363 }
364 
366  EKLMAlignment* alignment, EKLMSegmentAlignment* segmentAlignment,
367  const char* outputFile)
368 {
369  int iSection, iLayer, iSector, iPlane, iSegment, segment, param;
370  float value;
371  KLMAlignmentData* alignmentData;
372  TFile* f;
373  TTree* tEKLMModule, *tEKLMSegment;
374  f = new TFile(outputFile, "recreate");
375  tEKLMModule = new TTree("eklm_module", "");
376  tEKLMModule->Branch("section", &iSection, "section/I");
377  tEKLMModule->Branch("layer", &iLayer, "layer/I");
378  tEKLMModule->Branch("sector", &iSector, "sector/I");
379  tEKLMModule->Branch("param", &param, "param/I");
380  tEKLMModule->Branch("value", &value, "value/F");
381  tEKLMSegment = new TTree("eklm_segment", "");
382  tEKLMSegment->Branch("section", &iSection, "section/I");
383  tEKLMSegment->Branch("layer", &iLayer, "layer/I");
384  tEKLMSegment->Branch("sector", &iSector, "sector/I");
385  tEKLMSegment->Branch("plane", &iPlane, "plane/I");
386  tEKLMSegment->Branch("segment", &iSegment, "segment/I");
387  tEKLMSegment->Branch("param", &param, "param/I");
388  tEKLMSegment->Branch("value", &value, "value/F");
389  for (iSection = 1; iSection <= m_GeoDat->getNSections(); iSection++) {
390  for (iLayer = 1; iLayer <= m_GeoDat->getNDetectorLayers(iSection);
391  iLayer++) {
392  for (iSector = 1; iSector <= m_GeoDat->getNSectors(); iSector++) {
393  int module = m_ElementNumbers->moduleNumberEKLM(iSection, iSector, iLayer);
394  alignmentData = const_cast<KLMAlignmentData*>(
395  alignment->getModuleAlignment(module));
396  param = 1;
397  value = alignmentData->getDeltaU();
398  tEKLMModule->Fill();
399  /* cppcheck-suppress redundantAssignment */
400  param = 2;
401  /* cppcheck-suppress redundantAssignment */
402  value = alignmentData->getDeltaV();
403  tEKLMModule->Fill();
404  /* cppcheck-suppress redundantAssignment */
405  param = 6;
406  /* cppcheck-suppress redundantAssignment */
407  value = alignmentData->getDeltaGamma();
408  tEKLMModule->Fill();
409  for (iPlane = 1; iPlane <= m_GeoDat->getNPlanes(); iPlane++) {
410  for (iSegment = 1; iSegment <= m_GeoDat->getNSegments(); iSegment++) {
412  iSection, iLayer, iSector, iPlane, iSegment);
413  alignmentData = const_cast<KLMAlignmentData*>(
414  segmentAlignment->getSegmentAlignment(segment));
415  param = 2;
416  value = alignmentData->getDeltaV();
417  tEKLMSegment->Fill();
418  /* cppcheck-suppress redundantAssignment */
419  param = 6;
420  /* cppcheck-suppress redundantAssignment */
421  value = alignmentData->getDeltaGamma();
422  tEKLMSegment->Fill();
423  }
424  }
425  }
426  }
427  }
428  f->cd();
429  tEKLMModule->Write();
430  tEKLMSegment->Write();
431  delete tEKLMModule;
432  delete tEKLMSegment;
433  delete f;
434 }
Belle2::Unit::cm
static const double cm
Standard units with the value = 1.
Definition: Unit.h:57
Belle2::EKLMGeometry::getNSegments
int getNSegments() const
Get number of segments.
Definition: EKLMGeometry.h:1725
Belle2::KLMDisplacementGenerator::generateFixedModuleDisplacement
void generateFixedModuleDisplacement(double deltaU, double deltaV, double deltaGamma)
Generation of fixed module displacements.
Definition: KLMDisplacementGenerator.cc:62
Belle2::KLMDisplacementGenerator::~KLMDisplacementGenerator
~KLMDisplacementGenerator()
Destructor.
Definition: KLMDisplacementGenerator.cc:35
Belle2::EKLMSegmentAlignment::setSegmentAlignment
void setSegmentAlignment(uint16_t segment, KLMAlignmentData *dat)
Set segment alignment data.
Definition: EKLMSegmentAlignment.cc:24
Belle2::EKLMElementNumbers
EKLM element numbers.
Definition: EKLMElementNumbers.h:34
Belle2::KLMDisplacementGenerator::readDisplacementFromROOTFile
void readDisplacementFromROOTFile(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment, const char *inputFile)
Read displacement from ROOT file.
Definition: KLMDisplacementGenerator.cc:187
Belle2::EKLMGeometry::getNPlanes
int getNPlanes() const
Get number of planes.
Definition: EKLMGeometry.h:1717
Belle2::EKLM::AlignmentChecker::checkSegmentAlignment
bool checkSegmentAlignment(int section, int layer, int sector, int plane, int segment, const KLMAlignmentData *sectorAlignment, const KLMAlignmentData *segmentAlignment, bool calledFromSectorCheck) const
Check segment alignment.
Definition: AlignmentChecker.cc:197
Belle2::KLMDisplacementGenerator::saveDisplacement
void saveDisplacement(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment, const char *outputFile)
Save displacements to a ROOT file.
Definition: KLMDisplacementGenerator.cc:365
Belle2::KLMDisplacementGenerator::studyModuleAlignmentLimits
void studyModuleAlignmentLimits(TFile *f)
Generate random module displacements and check if they are correct (no overlaps).
Definition: KLMDisplacementGenerator.cc:249
Belle2::KLMDisplacementGenerator::studySegmentAlignmentLimits
void studySegmentAlignmentLimits(TFile *f)
Generate random segment displacements and check if they are correct (no overlaps).
Definition: KLMDisplacementGenerator.cc:297
Belle2::EKLMAlignment
Class to store EKLM alignment data in the database.
Definition: EKLMAlignment.h:40
Belle2::EKLM::AlignmentChecker::checkSectorAlignment
bool checkSectorAlignment(int section, int layer, int sector, const KLMAlignmentData *sectorAlignment) const
Check sector alignment.
Definition: AlignmentChecker.cc:84
Belle2::EKLM::AlignmentChecker
Class for EKLM alignment checking.
Definition: AlignmentChecker.h:38
Belle2::Unit::rad
static const double rad
Standard of [angle].
Definition: Unit.h:60
Belle2::KLMDisplacementGenerator::m_eklmElementNumbers
const EKLMElementNumbers * m_eklmElementNumbers
EKLM element numbers.
Definition: KLMDisplacementGenerator.h:154
Belle2::KLMDisplacementGenerator::generateRandomDisplacement
void generateRandomDisplacement(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment, bool displaceModule, bool displaceSegment, bool moduleSameDisplacement=false, bool moduleZeroDeltaU=false, bool moduleZeroDeltaV=false, bool moduleZeroDeltaGamma=false)
Generation of random displacements.
Definition: KLMDisplacementGenerator.cc:81
Belle2::KLMAlignmentData::setDeltaU
void setDeltaU(float deltaU)
Set shift in U.
Definition: KLMAlignmentData.h:96
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KLMAlignmentData::setDeltaGamma
void setDeltaGamma(float deltaGamma)
Set rotation in alpha.
Definition: KLMAlignmentData.h:181
Belle2::EKLMElementNumbers::segmentNumber
int segmentNumber(int section, int layer, int sector, int plane, int segment) const
Get segment number.
Definition: EKLMElementNumbers.cc:191
Belle2::KLMDisplacementGenerator::fillZeroDisplacements
void fillZeroDisplacements(EKLMAlignment *alignment, EKLMSegmentAlignment *segmentAlignment)
Fill EKLMAlignment with zero displacements.
Definition: KLMDisplacementGenerator.cc:39
Belle2::KLMElementNumbers::moduleNumberEKLM
uint16_t moduleNumberEKLM(int section, int sector, int layer) const
Get module number for EKLM.
Definition: KLMElementNumbers.cc:166
Belle2::EKLMGeometry::getNDetectorLayers
int getNDetectorLayers(int section) const
Get number of detector layers.
Definition: EKLMGeometry.cc:295
Belle2::EKLMSegmentAlignment
Class to store EKLM alignment data in the database.
Definition: EKLMSegmentAlignment.h:40
Belle2::KLMDisplacementGenerator::m_ElementNumbers
const KLMElementNumbers * m_ElementNumbers
Element numbers.
Definition: KLMDisplacementGenerator.h:151
Belle2::EKLMSegmentAlignment::getSegmentAlignment
const KLMAlignmentData * getSegmentAlignment(uint16_t segment) const
Get segment alignment data.
Definition: EKLMSegmentAlignment.cc:36
Belle2::KLMDisplacementGenerator::m_GeoDat
const EKLM::GeometryData * m_GeoDat
Geometry data.
Definition: KLMDisplacementGenerator.h:148
Belle2::EKLM::AlignmentChecker::checkAlignment
bool checkAlignment(const EKLMAlignment *alignment, const EKLMSegmentAlignment *segmentAlignment) const
Check alignment.
Definition: AlignmentChecker.cc:310
Belle2::KLMElementNumbers
KLM element numbers.
Definition: KLMElementNumbers.h:37
Belle2::KLMAlignmentData
KLM Alignment data.
Definition: KLMAlignmentData.h:33
Belle2::EKLMGeometry::getNSections
int getNSections() const
Get number of sections.
Definition: EKLMGeometry.h:1687
Belle2::KLMDisplacementGenerator::studyAlignmentLimits
void studyAlignmentLimits(const char *outputFile)
Generate random displacements and check if they are correct (no overlaps).
Definition: KLMDisplacementGenerator.cc:356
Belle2::EKLMGeometry::getNSectors
int getNSectors() const
Get number of sectors.
Definition: EKLMGeometry.h:1709
Belle2::KLMAlignmentData::setDeltaV
void setDeltaV(float deltaV)
Set shift in V.
Definition: KLMAlignmentData.h:113
Belle2::KLMDisplacementGenerator::KLMDisplacementGenerator
KLMDisplacementGenerator()
Constructor.
Definition: KLMDisplacementGenerator.cc:28