Belle II Software  release-05-01-25
KLMCalibrationChecker.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2020 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Giacomo De Pietro *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own include. */
12 #include <klm/calibration/KLMCalibrationChecker.h>
13 
14 /* KLM headers. */
15 #include <klm/dataobjects/bklm/BKLMElementNumbers.h>
16 #include <klm/dataobjects/eklm/EKLMElementNumbers.h>
17 #include <klm/dataobjects/KLMChannelIndex.h>
18 #include <klm/dbobjects/bklm/BKLMAlignment.h>
19 #include <klm/dbobjects/eklm/EKLMAlignment.h>
20 #include <klm/dbobjects/eklm/EKLMSegmentAlignment.h>
21 #include <klm/dbobjects/KLMStripEfficiency.h>
22 
23 /* Belle II headers. */
24 #include <framework/database/Database.h>
25 #include <framework/database/DBStore.h>
26 #include <framework/database/Configuration.h>
27 #include <framework/datastore/DataStore.h>
28 
29 /* ROOT include. */
30 #include <TCanvas.h>
31 #include <TFile.h>
32 #include <TH1F.h>
33 #include <TString.h>
34 #include <TTree.h>
35 
36 using namespace Belle2;
37 
39  m_experiment(0),
40  m_run(0),
41  m_ElementNumbers(&(KLMElementNumbers::Instance()))
42 {
43 }
44 
46 {
47 }
48 
49 void KLMCalibrationChecker::setExperimentRun(int experiment, int run)
50 {
51  m_experiment = experiment;
52  m_run = run;
53  if (m_EventMetaData.isValid()) {
54  m_EventMetaData->setExperiment(experiment);
55  m_EventMetaData->setRun(run);
56  }
57 }
58 
60 {
61  /* Mimic a module initialization. */
62  StoreObjPtr<EventMetaData> eventMetaData;
64  m_EventMetaData.registerInDataStore();
66  if (!m_EventMetaData.isValid())
67  m_EventMetaData.construct(1, m_run, m_experiment);
68  /* Database instance and configuration. */
69  DBStore& dbStore = DBStore::Instance();
70  dbStore.update();
71  dbStore.updateEvent();
72  auto& dbConfiguration = Conditions::Configuration::getInstance();
73  if ((m_testingPayloadName != "") and (m_GlobalTagName == ""))
74  dbConfiguration.prependTestingPayloadLocation(m_testingPayloadName);
75  else if ((m_testingPayloadName == "") and (m_GlobalTagName != ""))
76  dbConfiguration.prependGlobalTag(m_GlobalTagName);
77  else
78  B2FATAL("Setting both testing payload and Global Tag or setting no one of them.");
79 }
80 
82 {
83  /* Reset both DataStore and Database. */
85  Database::Instance().reset(false);
86  DBStore::Instance().reset(false);
87 }
88 
90 {
91  /* Initialize the database. */
93  /* Now we can read the payload. */
94  DBObjPtr<BKLMAlignment> bklmAlignment;
95  DBObjPtr<EKLMAlignment> eklmAlignment;
96  DBObjPtr<EKLMSegmentAlignment> eklmSegmentAlignment;
97  DBObjPtr<BKLMAlignment> bklmAlignmentErrors("BKLMAlignment_ERRORS");
98  DBObjPtr<EKLMAlignment> eklmAlignmentErrors("EKLMAlignment_ERRORS");
99  DBObjPtr<EKLMSegmentAlignment> eklmSegmentAlignmentErrors("EKLMSegmentAlignment_ERRORS");
100  DBObjPtr<BKLMAlignment> bklmAlignmentCorrections("BKLMAlignment_CORRECTIONS");
101  DBObjPtr<EKLMAlignment> eklmAlignmentCorrections("EKLMAlignment_CORRECTIONS");
102  DBObjPtr<EKLMSegmentAlignment> eklmSegmentAlignmentCorrections("EKLMSegmentAlignment_CORRECTIONS");
103  if (!bklmAlignment.isValid() ||
104  !eklmAlignment.isValid() ||
105  !eklmSegmentAlignment.isValid() ||
106  !bklmAlignmentErrors.isValid() ||
107  !eklmAlignmentErrors.isValid() ||
108  !eklmSegmentAlignmentErrors.isValid() ||
109  !bklmAlignmentCorrections.isValid() ||
110  !eklmAlignmentCorrections.isValid() ||
111  !eklmSegmentAlignmentCorrections.isValid())
112  B2FATAL("Alignment data are not valid.");
113  if (m_GlobalTagName != "") {
114  printPayloadInformation(bklmAlignment);
115  printPayloadInformation(eklmAlignment);
116  printPayloadInformation(eklmSegmentAlignment);
117  printPayloadInformation(bklmAlignmentErrors);
118  printPayloadInformation(eklmAlignmentErrors);
119  printPayloadInformation(eklmSegmentAlignmentErrors);
120  printPayloadInformation(bklmAlignmentCorrections);
121  printPayloadInformation(eklmAlignmentCorrections);
122  printPayloadInformation(eklmSegmentAlignmentCorrections);
123  }
124  /* Create trees with alignment results. */
125  int section, sector, layer, plane, segment, param;
126  float value, correction, error;
127  TFile* alignmentResults = new TFile(m_AlignmentResultsFile.c_str(),
128  "recreate");
129  TTree* bklmModuleTree = new TTree("bklm_module",
130  "BKLM module alignment data.");
131  bklmModuleTree->Branch("experiment", &m_experiment, "experiment/I");
132  bklmModuleTree->Branch("run", &m_run, "run/I");
133  bklmModuleTree->Branch("section", &section, "section/I");
134  bklmModuleTree->Branch("sector", &sector, "sector/I");
135  bklmModuleTree->Branch("layer", &layer, "layer/I");
136  bklmModuleTree->Branch("param", &param, "param/I");
137  bklmModuleTree->Branch("value", &value, "value/F");
138  bklmModuleTree->Branch("correction", &correction, "correction/F");
139  bklmModuleTree->Branch("error", &error, "error/F");
140  TTree* eklmModuleTree = new TTree("eklm_module",
141  "EKLM module alignment data.");
142  eklmModuleTree->Branch("experiment", &m_experiment, "experiment/I");
143  eklmModuleTree->Branch("run", &m_run, "run/I");
144  eklmModuleTree->Branch("section", &section, "section/I");
145  eklmModuleTree->Branch("sector", &sector, "sector/I");
146  eklmModuleTree->Branch("layer", &layer, "layer/I");
147  eklmModuleTree->Branch("param", &param, "param/I");
148  eklmModuleTree->Branch("value", &value, "value/F");
149  eklmModuleTree->Branch("correction", &correction, "correction/F");
150  eklmModuleTree->Branch("error", &error, "error/F");
151  TTree* eklmSegmentTree = new TTree("eklm_segment",
152  "EKLM segment alignment data.");
153  eklmSegmentTree->Branch("experiment", &m_experiment, "experiment/I");
154  eklmSegmentTree->Branch("run", &m_run, "run/I");
155  eklmSegmentTree->Branch("section", &section, "section/I");
156  eklmSegmentTree->Branch("sector", &sector, "sector/I");
157  eklmSegmentTree->Branch("layer", &layer, "layer/I");
158  eklmSegmentTree->Branch("plane", &plane, "plane/I");
159  eklmSegmentTree->Branch("segment", &segment, "segment/I");
160  eklmSegmentTree->Branch("param", &param, "param/I");
161  eklmSegmentTree->Branch("value", &value, "value/F");
162  eklmSegmentTree->Branch("correction", &correction, "correction/F");
163  eklmSegmentTree->Branch("error", &error, "error/F");
164  const KLMAlignmentData* alignment, *alignmentError, *alignmentCorrection;
165  KLMAlignmentData zeroAlignment(0, 0, 0, 0, 0, 0);
167  for (KLMChannelIndex& klmModule : klmModules) {
168  uint16_t module = klmModule.getKLMModuleNumber();
169  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM) {
170  alignment = bklmAlignment->getModuleAlignment(module);
171  alignmentError = bklmAlignmentErrors->getModuleAlignment(module);
172  alignmentCorrection =
173  bklmAlignmentCorrections->getModuleAlignment(module);
174  } else {
175  alignment = eklmAlignment->getModuleAlignment(module);
176  alignmentError = eklmAlignmentErrors->getModuleAlignment(module);
177  alignmentCorrection =
178  eklmAlignmentCorrections->getModuleAlignment(module);
179  }
180  if (alignment == nullptr)
181  B2FATAL("Incomplete KLM alignment data.");
182  if ((alignmentError == nullptr) && (alignmentCorrection == nullptr)) {
183  B2WARNING("Alignment is not determined for KLM module."
184  << LogVar("Module", module));
185  alignmentError = &zeroAlignment;
186  alignmentCorrection = &zeroAlignment;
187  } else if ((alignmentError == nullptr) ||
188  (alignmentCorrection == nullptr)) {
189  B2FATAL("Inconsistent undtermined parameters.");
190  }
191  section = klmModule.getSection();
192  sector = klmModule.getSector();
193  layer = klmModule.getLayer();
195  value = alignment->getDeltaU();
196  error = alignmentError->getDeltaU();
197  correction = alignmentCorrection->getDeltaU();
198  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
199  bklmModuleTree->Fill();
200  else
201  eklmModuleTree->Fill();
202  /* cppcheck-suppress redundantAssignment */
204  /* cppcheck-suppress redundantAssignment */
205  value = alignment->getDeltaV();
206  /* cppcheck-suppress redundantAssignment */
207  error = alignmentError->getDeltaV();
208  /* cppcheck-suppress redundantAssignment */
209  correction = alignmentCorrection->getDeltaV();
210  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
211  bklmModuleTree->Fill();
212  else
213  eklmModuleTree->Fill();
214  /* cppcheck-suppress redundantAssignment */
216  /* cppcheck-suppress redundantAssignment */
217  value = alignment->getDeltaW();
218  /* cppcheck-suppress redundantAssignment */
219  error = alignmentError->getDeltaW();
220  /* cppcheck-suppress redundantAssignment */
221  correction = alignmentCorrection->getDeltaW();
222  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
223  bklmModuleTree->Fill();
224  else
225  eklmModuleTree->Fill();
226  /* cppcheck-suppress redundantAssignment */
228  /* cppcheck-suppress redundantAssignment */
229  value = alignment->getDeltaAlpha();
230  /* cppcheck-suppress redundantAssignment */
231  error = alignmentError->getDeltaAlpha();
232  /* cppcheck-suppress redundantAssignment */
233  correction = alignmentCorrection->getDeltaAlpha();
234  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
235  bklmModuleTree->Fill();
236  else
237  eklmModuleTree->Fill();
238  /* cppcheck-suppress redundantAssignment */
240  /* cppcheck-suppress redundantAssignment */
241  value = alignment->getDeltaBeta();
242  /* cppcheck-suppress redundantAssignment */
243  error = alignmentError->getDeltaBeta();
244  /* cppcheck-suppress redundantAssignment */
245  correction = alignmentCorrection->getDeltaBeta();
246  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
247  bklmModuleTree->Fill();
248  else
249  eklmModuleTree->Fill();
250  /* cppcheck-suppress redundantAssignment */
252  /* cppcheck-suppress redundantAssignment */
253  value = alignment->getDeltaGamma();
254  /* cppcheck-suppress redundantAssignment */
255  error = alignmentError->getDeltaGamma();
256  /* cppcheck-suppress redundantAssignment */
257  correction = alignmentCorrection->getDeltaGamma();
258  if (klmModule.getSubdetector() == KLMElementNumbers::c_BKLM)
259  bklmModuleTree->Fill();
260  else
261  eklmModuleTree->Fill();
262  }
265  eklmSegment = eklmSegments.beginEKLM();
266  eklmSegment.useEKLMSegments();
267  for (; eklmSegment != eklmSegments.endEKLM(); ++eklmSegment) {
268  int eklmSegmentNumber = eklmSegment.getEKLMSegmentNumber();
269  alignment = eklmSegmentAlignment->getSegmentAlignment(eklmSegmentNumber);
270  alignmentError =
271  eklmSegmentAlignmentErrors->getSegmentAlignment(eklmSegmentNumber);
272  alignmentCorrection =
273  eklmSegmentAlignmentCorrections->getSegmentAlignment(eklmSegmentNumber);
274  if (alignment == nullptr)
275  B2FATAL("Incomplete KLM alignment data.");
276  if ((alignmentError == nullptr) && (alignmentCorrection == nullptr)) {
277  /*
278  * The segment alignment is not determined now.
279  * TODO: If it will be turned on in the future, add a warning here.
280  */
281  alignmentError = &zeroAlignment;
282  alignmentCorrection = &zeroAlignment;
283  } else if ((alignmentError == nullptr) ||
284  (alignmentCorrection == nullptr)) {
285  B2FATAL("Inconsistent undtermined parameters.");
286  }
287  section = eklmSegment.getSection();
288  sector = eklmSegment.getSector();
289  layer = eklmSegment.getLayer();
290  plane = eklmSegment.getPlane();
291  segment = eklmSegment.getStrip();
293  value = alignment->getDeltaU();
294  error = alignmentError->getDeltaU();
295  correction = alignmentCorrection->getDeltaU();
296  eklmSegmentTree->Fill();
297  /* cppcheck-suppress redundantAssignment */
299  /* cppcheck-suppress redundantAssignment */
300  value = alignment->getDeltaV();
301  /* cppcheck-suppress redundantAssignment */
302  error = alignmentError->getDeltaV();
303  /* cppcheck-suppress redundantAssignment */
304  correction = alignmentCorrection->getDeltaV();
305  eklmSegmentTree->Fill();
306  /* cppcheck-suppress redundantAssignment */
308  /* cppcheck-suppress redundantAssignment */
309  value = alignment->getDeltaW();
310  /* cppcheck-suppress redundantAssignment */
311  error = alignmentError->getDeltaW();
312  /* cppcheck-suppress redundantAssignment */
313  correction = alignmentCorrection->getDeltaW();
314  eklmSegmentTree->Fill();
315  /* cppcheck-suppress redundantAssignment */
317  /* cppcheck-suppress redundantAssignment */
318  value = alignment->getDeltaAlpha();
319  /* cppcheck-suppress redundantAssignment */
320  error = alignmentError->getDeltaAlpha();
321  /* cppcheck-suppress redundantAssignment */
322  correction = alignmentCorrection->getDeltaAlpha();
323  eklmSegmentTree->Fill();
324  /* cppcheck-suppress redundantAssignment */
326  /* cppcheck-suppress redundantAssignment */
327  value = alignment->getDeltaBeta();
328  /* cppcheck-suppress redundantAssignment */
329  error = alignmentError->getDeltaBeta();
330  /* cppcheck-suppress redundantAssignment */
331  correction = alignmentCorrection->getDeltaBeta();
332  eklmSegmentTree->Fill();
333  /* cppcheck-suppress redundantAssignment */
335  /* cppcheck-suppress redundantAssignment */
336  value = alignment->getDeltaGamma();
337  /* cppcheck-suppress redundantAssignment */
338  error = alignmentError->getDeltaGamma();
339  /* cppcheck-suppress redundantAssignment */
340  correction = alignmentCorrection->getDeltaGamma();
341  eklmSegmentTree->Fill();
342  }
343  bklmModuleTree->Write();
344  eklmModuleTree->Write();
345  eklmSegmentTree->Write();
346  delete bklmModuleTree;
347  delete eklmModuleTree;
348  delete eklmSegmentTree;
349  delete alignmentResults;
350  /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
351  resetDatabase();
352 }
353 
355 {
356  /* Initialize the database. */
358  /* Now we can read the payload. */
359  DBObjPtr<KLMStripEfficiency> stripEfficiency;
360  if (!stripEfficiency.isValid())
361  B2FATAL("Strip efficiency data are not valid.");
362  if (m_GlobalTagName != "")
363  printPayloadInformation(stripEfficiency);
364  /* Create trees with strip efficiency measurement results. */
365  int subdetector, section, sector, layer, plane;
366  float efficiency, error;
367  TFile* stripEfficiencyResults =
368  new TFile(m_StripEfficiencyResultsFile.c_str(), "recreate");
369  TTree* efficiencyTree = new TTree("efficiency", "KLM strip efficiency data.");
370  efficiencyTree->Branch("experiment", &m_experiment, "experiment/I");
371  efficiencyTree->Branch("run", &m_run, "run/I");
372  efficiencyTree->Branch("subdetector", &subdetector, "subdetector/I");
373  efficiencyTree->Branch("section", &section, "section/I");
374  efficiencyTree->Branch("sector", &sector, "sector/I");
375  efficiencyTree->Branch("layer", &layer, "layer/I");
376  efficiencyTree->Branch("plane", &plane, "plane/I");
377  efficiencyTree->Branch("efficiency", &efficiency, "efficiency/F");
378  efficiencyTree->Branch("error", &error, "error/F");
380  for (KLMChannelIndex& klmPlane : klmPlanes) {
381  subdetector = klmPlane.getSubdetector();
382  section = klmPlane.getSection();
383  sector = klmPlane.getSector();
384  layer = klmPlane.getLayer();
385  plane = klmPlane.getPlane();
386  uint16_t channel = m_ElementNumbers->channelNumber(
387  subdetector, section, sector, layer, plane, 1);
388  efficiency = stripEfficiency->getEfficiency(channel);
389  error = stripEfficiency->getEfficiencyError(channel);
390  efficiencyTree->Fill();
391  }
392  efficiencyTree->Write();
393  delete efficiencyTree;
394  delete stripEfficiencyResults;
395  /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
396  resetDatabase();
397 }
398 
400 {
401  /* Initialize the database. */
403  /* Now we can read the payload. */
404  DBObjPtr<KLMStripEfficiency> stripEfficiency;
405  if (not stripEfficiency.isValid())
406  B2FATAL("Strip efficiency data are not valid.");
407  if (m_GlobalTagName != "")
408  printPayloadInformation(stripEfficiency);
409  /* Finally, loop over KLM sectors to check the efficiency. */
411  TCanvas* canvas = new TCanvas();
412  for (KLMChannelIndex& klmSector : klmSectors) {
413  int subdetector = klmSector.getSubdetector();
414  int section = klmSector.getSection();
415  int sector = klmSector.getSector();
416  /* Setup the histogram. */
417  TH1F* hist = new TH1F("plane_histogram", "", 30, 0.5, 30.5);
418  hist->GetYaxis()->SetTitle("Efficiency");
419  hist->SetMinimum(0.4);
420  hist->SetMaximum(1.);
421  hist->SetMarkerStyle(20);
422  hist->SetMarkerSize(0.5);
423  TString title;
424  if (subdetector == KLMElementNumbers::c_BKLM) {
426  title.Form("BKLM backward sector %d", sector);
427  else
428  title.Form("BKLM forward sector %d", sector);
429  hist->SetTitle(title.Data());
430  hist->GetXaxis()->SetTitle("(Layer - 1) * 2 + plane + 1");
431  for (int layer = 1; layer <= BKLMElementNumbers::getMaximalLayerNumber(); layer++) {
432  for (int plane = 0; plane <= BKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
433  int bin = (layer - 1) * 2 + plane + 1;
434  float efficiency = stripEfficiency->getBarrelEfficiency(section, sector, layer, plane, 2);
435  float efficiencyError = stripEfficiency->getBarrelEfficiencyError(section, sector, layer, plane, 2);
436  hist->SetBinContent(bin, efficiency);
437  hist->SetBinError(bin, efficiencyError);
438  }
439  }
440  } else {
441  if (section == EKLMElementNumbers::c_BackwardSection) {
442  hist->SetBins(24, 0.5, 24.5);
443  title.Form("EKLM backward sector %d", sector);
444  } else {
445  hist->SetBins(28, 0.5, 28.5);
446  title.Form("EKLM forward sector %d", sector);
447  }
448  hist->SetTitle(title.Data());
449  hist->GetXaxis()->SetTitle("(Layer - 1) * 2 + plane");
450  const EKLMElementNumbers* elementNumbersEKLM = &(EKLMElementNumbers::Instance());
451  for (int layer = 1; layer <= elementNumbersEKLM->getMaximalDetectorLayerNumber(section); layer++) {
452  for (int plane = 1; plane <= EKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
453  int bin = (layer - 1) * 2 + plane;
454  float efficiency = stripEfficiency->getEndcapEfficiency(section, sector, layer, plane, 2);
455  float efficiencyError = stripEfficiency->getEndcapEfficiencyError(section, sector, layer, plane, 2);
456  hist->SetBinContent(bin, efficiency);
457  hist->SetBinError(bin, efficiencyError);
458  }
459  }
460  }
461  hist->Draw("e");
462  TString name;
463  name.Form("efficiency_subdetector_%d_section_%d_sector_%d.pdf", subdetector, section, sector);
464  canvas->Print(name.Data());
465  canvas->Update();
466  delete hist;
467  }
468  delete canvas;
469  /* Reset the database. Needed to avoid mess if we call this method multiple times with different GTs. */
470  resetDatabase();
471 }
Belle2::KLMChannelIndex::endEKLM
KLMChannelIndex & endEKLM()
Last channel for EKLM.
Definition: KLMChannelIndex.cc:211
Belle2::KLMCalibrationChecker::checkStripEfficiency
void checkStripEfficiency()
Check strip efficiency.
Definition: KLMCalibrationChecker.cc:354
Belle2::KLMChannelIndex::getSector
int getSector() const
Get sector.
Definition: KLMChannelIndex.h:132
Belle2::KLMCalibrationChecker::initializeDatabase
void initializeDatabase()
Initialize the database.
Definition: KLMCalibrationChecker.cc:59
Belle2::KLMCalibrationChecker::checkAlignment
void checkAlignment()
Check alignment.
Definition: KLMCalibrationChecker.cc:89
Belle2::KLMAlignmentData::c_DeltaV
@ c_DeltaV
Shift in V (EKLM: local Y).
Definition: KLMAlignmentData.h:46
Belle2::EKLMElementNumbers
EKLM element numbers.
Definition: EKLMElementNumbers.h:34
Belle2::KLMCalibrationChecker::~KLMCalibrationChecker
~KLMCalibrationChecker()
Destructor.
Definition: KLMCalibrationChecker.cc:45
Belle2::KLMChannelIndex::getStrip
int getStrip() const
Get strip.
Definition: KLMChannelIndex.h:156
Belle2::DataStore::Instance
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
Belle2::KLMCalibrationChecker::m_StripEfficiencyResultsFile
std::string m_StripEfficiencyResultsFile
Output file for alignment results.
Definition: KLMCalibrationChecker.h:149
Belle2::KLMChannelIndex::c_IndexLevelLayer
@ c_IndexLevelLayer
Layer.
Definition: KLMChannelIndex.h:52
Belle2::KLMChannelIndex::getSection
int getSection() const
Get section.
Definition: KLMChannelIndex.h:124
Belle2::KLMChannelIndex::beginEKLM
KLMChannelIndex beginEKLM()
First channel for EKLM.
Definition: KLMChannelIndex.cc:205
Belle2::KLMCalibrationChecker::createStripEfficiencyHistograms
void createStripEfficiencyHistograms()
Create strip efficiency histograms.
Definition: KLMCalibrationChecker.cc:399
Belle2::KLMAlignmentData::c_DeltaW
@ c_DeltaW
Shift in W.
Definition: KLMAlignmentData.h:49
Belle2::DataStore::setInitializeActive
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
Belle2::KLMAlignmentData::getDeltaV
float getDeltaV() const
Get shift in V.
Definition: KLMAlignmentData.h:104
Belle2::KLMChannelIndex::getLayer
int getLayer() const
Get layer.
Definition: KLMChannelIndex.h:140
Belle2::KLMCalibrationChecker::m_GlobalTagName
std::string m_GlobalTagName
Global Tag name.
Definition: KLMCalibrationChecker.h:143
Belle2::KLMElementNumbers::channelNumber
uint16_t channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
Definition: KLMElementNumbers.cc:37
Belle2::DBStore::reset
void reset(bool keepEntries=false)
Invalidate all payloads.
Definition: DBStore.cc:185
Belle2::EKLMElementNumbers::c_BackwardSection
@ c_BackwardSection
Backward.
Definition: EKLMElementNumbers.h:44
Belle2::KLMAlignmentData::getDeltaAlpha
float getDeltaAlpha() const
Get rotation in alpha.
Definition: KLMAlignmentData.h:138
Belle2::KLMCalibrationChecker::m_EventMetaData
StoreObjPtr< EventMetaData > m_EventMetaData
Event metadata.
Definition: KLMCalibrationChecker.h:155
Belle2::KLMCalibrationChecker::m_run
int m_run
Run number.
Definition: KLMCalibrationChecker.h:137
Belle2::KLMChannelIndex::c_IndexLevelSector
@ c_IndexLevelSector
Sector.
Definition: KLMChannelIndex.h:49
Belle2::Database::reset
static void reset(bool keepConfig=false)
Reset the database instance.
Definition: Database.cc:62
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::DataStore::reset
void reset(EDurability durability)
Frees memory occupied by data store items and removes all objects from the map.
Definition: DataStore.cc:86
Belle2::KLMAlignmentData::getDeltaU
float getDeltaU() const
Get shift in U.
Definition: KLMAlignmentData.h:87
Belle2::KLMCalibrationChecker::printPayloadInformation
void printPayloadInformation(DBObjPtr< T > &dbObject)
Print payload information.
Definition: KLMCalibrationChecker.h:124
Belle2::KLMCalibrationChecker::m_ElementNumbers
const KLMElementNumbers * m_ElementNumbers
Element numbers.
Definition: KLMCalibrationChecker.h:152
Belle2::KLMChannelIndex::c_IndexLevelPlane
@ c_IndexLevelPlane
Plane.
Definition: KLMChannelIndex.h:55
Belle2::KLMCalibrationChecker::m_experiment
int m_experiment
Experiment number.
Definition: KLMCalibrationChecker.h:134
Belle2::BKLMElementNumbers::c_BackwardSection
@ c_BackwardSection
Backward.
Definition: BKLMElementNumbers.h:44
Belle2::BKLMElementNumbers::getMaximalLayerNumber
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
Definition: BKLMElementNumbers.h:251
Belle2::KLMCalibrationChecker::setExperimentRun
void setExperimentRun(int experiment, int run)
Set experiment and run numbers.
Definition: KLMCalibrationChecker.cc:49
Belle2::KLMAlignmentData::getDeltaBeta
float getDeltaBeta() const
Get rotation in alpha.
Definition: KLMAlignmentData.h:155
Belle2::DBStore
Singleton class to cache database objects.
Definition: DBStore.h:42
Belle2::Conditions::Configuration::getInstance
static Configuration & getInstance()
Get a reference to the instance which will be used when the Database is initialized.
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::EKLMElementNumbers::Instance
static const EKLMElementNumbers & Instance()
Instantiation.
Definition: EKLMElementNumbers.cc:20
Belle2::DBStore::Instance
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:36
Belle2::KLMCalibrationChecker::KLMCalibrationChecker
KLMCalibrationChecker()
Constructor.
Definition: KLMCalibrationChecker.cc:38
Belle2::KLMAlignmentData::getDeltaGamma
float getDeltaGamma() const
Get rotation in alpha.
Definition: KLMAlignmentData.h:172
Belle2::KLMElementNumbers::c_BKLM
@ c_BKLM
BKLM.
Definition: KLMElementNumbers.h:47
Belle2::KLMChannelIndex::useEKLMSegments
void useEKLMSegments(bool useSegments=true)
Iterate over EKLM segments instead of strips.
Definition: KLMChannelIndex.cc:117
Belle2::KLMCalibrationChecker::m_testingPayloadName
std::string m_testingPayloadName
Testing payload location.
Definition: KLMCalibrationChecker.h:140
Belle2::KLMChannelIndex::getEKLMSegmentNumber
int getEKLMSegmentNumber() const
Get EKLM segment number.
Definition: KLMChannelIndex.cc:183
Belle2::KLMChannelIndex::getPlane
int getPlane() const
Get plane.
Definition: KLMChannelIndex.h:148
Belle2::EKLMElementNumbers::getMaximalPlaneNumber
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.
Definition: EKLMElementNumbers.h:313
Belle2::DBStore::updateEvent
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:150
Belle2::DBAccessorBase::isValid
bool isValid() const
Check whether a valid object was obtained from the database.
Definition: DBAccessorBase.h:75
Belle2::Database::Instance
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:54
Belle2::BKLMElementNumbers::getMaximalPlaneNumber
static constexpr int getMaximalPlaneNumber()
Get maximal plane number (0-based).
Definition: BKLMElementNumbers.h:259
Belle2::KLMChannelIndex
KLM channel index.
Definition: KLMChannelIndex.h:33
Belle2::KLMCalibrationChecker::resetDatabase
void resetDatabase()
Reset the database.
Definition: KLMCalibrationChecker.cc:81
Belle2::KLMElementNumbers
KLM element numbers.
Definition: KLMElementNumbers.h:37
Belle2::KLMAlignmentData
KLM Alignment data.
Definition: KLMAlignmentData.h:33
Belle2::KLMAlignmentData::c_DeltaAlpha
@ c_DeltaAlpha
Rotation in alpha.
Definition: KLMAlignmentData.h:52
Belle2::KLMAlignmentData::c_DeltaGamma
@ c_DeltaGamma
Rotation in gamma (EKLM: rotation in local plane).
Definition: KLMAlignmentData.h:58
Belle2::KLMAlignmentData::getDeltaW
float getDeltaW() const
Get shift in W.
Definition: KLMAlignmentData.h:121
Belle2::KLMCalibrationChecker::m_AlignmentResultsFile
std::string m_AlignmentResultsFile
Output file for alignment results.
Definition: KLMCalibrationChecker.h:146
Belle2::KLMAlignmentData::c_DeltaBeta
@ c_DeltaBeta
Rotation in beta.
Definition: KLMAlignmentData.h:55
Belle2::DBStore::update
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:87
Belle2::KLMAlignmentData::c_DeltaU
@ c_DeltaU
Shift in U (EKLM: local X).
Definition: KLMAlignmentData.h:43
Belle2::EKLMElementNumbers::getMaximalDetectorLayerNumber
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer number.
Definition: EKLMElementNumbers.cc:279
Belle2::KLMChannelIndex::c_IndexLevelStrip
@ c_IndexLevelStrip
Strip.
Definition: KLMChannelIndex.h:58