Belle II Software  release-05-01-25
KLMDatabaseImporter.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Kirill Chilikin, Giacomo De Pietro *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/calibration/KLMDatabaseImporter.h>
13 
14 /* Belle 2 headers. */
15 #include <framework/database/DBImportObjPtr.h>
16 #include <framework/database/IntervalOfValidity.h>
17 #include <framework/logging/Logger.h>
18 
19 /* ROOT headers. */
20 #include <TFile.h>
21 #include <TTree.h>
22 
23 /* C++ headers. */
24 #include <string>
25 
26 using namespace Belle2;
27 
29 {
30 }
31 
33 {
34 }
35 
36 void KLMDatabaseImporter::setIOV(int experimentLow, int runLow,
37  int experimentHigh, int runHigh)
38 {
39  m_ExperimentLow = experimentLow;
40  m_RunLow = runLow;
41  m_ExperimentHigh = experimentHigh;
42  m_RunHigh = runHigh;
43 }
44 
46  const KLMChannelStatus* channelStatus)
47 {
48  DBImportObjPtr<KLMChannelStatus> channelStatusImport;
49  channelStatusImport.construct(*channelStatus);
52  channelStatusImport.import(iov);
53 }
54 
56  const KLMScintillatorDigitizationParameters* digitizationParameters)
57 {
59  digPar.construct(*digitizationParameters);
62  digPar.import(iov);
63 }
64 
66  const KLMTimeConversion* timeConversion)
67 {
68  DBImportObjPtr<KLMTimeConversion> timeConversionImport;
69  timeConversionImport.construct(*timeConversion);
72  timeConversionImport.import(iov);
73 }
74 
76 {
77  DBImportObjPtr<KLMTimeWindow> timeWindowImport;
78  timeWindowImport.construct(*timeWindow);
81  timeWindowImport.import(iov);
82 }
83 
85  KLMStripEfficiency* stripEfficiency, std::string fileName)
86 {
87  TFile* file = TFile::Open(fileName.c_str(), "r");
88  if (!file) {
89  B2ERROR("KLMDatabaseImporter: calibration file " << fileName << " *** failed to open");
90  } else {
91  TTree* tree = (TTree*)file->Get("tree");
92  if (!tree) {
93  B2ERROR("KLMDatabaseImporter: calibration file " << fileName << " *** no tree named 'tree' found");
94  file->Close();
95  } else {
96  int isBarrel = 0;
97  tree->SetBranchAddress("isBarrel", &isBarrel);
98  int section = 0;
99  tree->SetBranchAddress("isForward", &section);
100  int sector = 0;
101  tree->SetBranchAddress("sector", &sector);
102  int layer = 0;
103  tree->SetBranchAddress("layer", &layer);
104  int plane = 0;
105  tree->SetBranchAddress("plane", &plane);
106  int strip = 0;
107  tree->SetBranchAddress("strip", &strip);
108  float efficiency = 1.;
109  tree->SetBranchAddress("efficiency", &efficiency);
110  float efficiencyError = 0.;
111  tree->SetBranchAddress("efficiencyError", &efficiencyError);
112 
113  for (int i = 0; i < tree->GetEntries(); i++) {
114  tree->GetEntry(i);
115  if (isBarrel)
116  stripEfficiency->setBarrelEfficiency(section, sector, layer, plane, strip, efficiency, efficiencyError);
117  else
118  stripEfficiency->setEndcapEfficiency(section, sector, layer, plane, strip, efficiency, efficiencyError);
119  }
120  }
121  file->Close();
122  }
123 }
124 
126  const KLMStripEfficiency* stripEfficiency)
127 {
128  DBImportObjPtr<KLMStripEfficiency> stripEfficiencyImport;
129  stripEfficiencyImport.construct(*stripEfficiency);
132  stripEfficiencyImport.import(iov);
133 }
134 
136  const BKLMAlignment* bklmAlignment, bool displacement)
137 {
138  std::string payloadName;
139  if (displacement)
140  payloadName = "BKLMDisplacement";
141  else
142  payloadName = "BKLMAlignment";
143  DBImportObjPtr<BKLMAlignment> bklmAlignmentImport(payloadName);
144  bklmAlignmentImport.construct(*bklmAlignment);
147  bklmAlignmentImport.import(iov);
148 }
149 
151  const EKLMAlignment* eklmAlignment, bool displacement)
152 {
153  std::string payloadName;
154  if (displacement)
155  payloadName = "EKLMDisplacement";
156  else
157  payloadName = "EKLMAlignment";
158  DBImportObjPtr<EKLMAlignment> eklmAlignmentImport(payloadName);
159  eklmAlignmentImport.construct(*eklmAlignment);
162  eklmAlignmentImport.import(iov);
163 }
164 
166  const EKLMSegmentAlignment* eklmSegmentAlignment, bool displacement)
167 {
168  std::string payloadName;
169  if (displacement)
170  payloadName = "EKLMSegmentDisplacement";
171  else
172  payloadName = "EKLMSegmentAlignment";
173  DBImportObjPtr<EKLMSegmentAlignment> eklmSegmentAlignmentImport(payloadName);
174  eklmSegmentAlignmentImport.construct(*eklmSegmentAlignment);
177  eklmSegmentAlignmentImport.import(iov);
178 }
179 
181  const BKLMAlignment* bklmAlignment, const EKLMAlignment* eklmAlignment,
182  const EKLMSegmentAlignment* eklmSegmentAlignment, bool displacement)
183 {
184  importBKLMAlignment(bklmAlignment, displacement);
185  importEKLMAlignment(eklmAlignment, displacement);
186  importEKLMSegmentAlignment(eklmSegmentAlignment, displacement);
187 }
Belle2::IntervalOfValidity
A class that describes the interval of experiments/runs for which an object in the database is valid.
Definition: IntervalOfValidity.h:35
Belle2::KLMDatabaseImporter::setIOV
void setIOV(int experimentLow, int runLow, int experimentHigh, int runHigh)
Set interval of validity.
Definition: KLMDatabaseImporter.cc:36
Belle2::KLMDatabaseImporter::m_ExperimentLow
int m_ExperimentLow
Low experiment.
Definition: KLMDatabaseImporter.h:140
Belle2::KLMDatabaseImporter::importEKLMAlignment
void importEKLMAlignment(const EKLMAlignment *eklmAlignment, bool displacement=false)
Import EKLM alignment.
Definition: KLMDatabaseImporter.cc:150
Belle2::KLMChannelStatus
KLM channel status.
Definition: KLMChannelStatus.h:37
Belle2::KLMDatabaseImporter::importStripEfficiency
void importStripEfficiency(const KLMStripEfficiency *stripEfficiency)
Import strip efficiencies.
Definition: KLMDatabaseImporter.cc:125
Belle2::KLMDatabaseImporter::importTimeConversion
void importTimeConversion(const KLMTimeConversion *timeConversion)
Import time conversion parameters.
Definition: KLMDatabaseImporter.cc:65
Belle2::DBImportObjPtr::construct
void construct(Args &&... params)
Construct an object of type T in this DBImportObjPtr using the provided constructor arguments.
Definition: DBImportObjPtr.h:57
Belle2::DBImportBase::import
bool import(const IntervalOfValidity &iov)
Import the object to database.
Definition: DBImportBase.cc:38
Belle2::KLMTimeWindow
DBObject containing KLM time window parameters used in KLMReconstructor module.
Definition: KLMTimeWindow.h:34
Belle2::KLMDatabaseImporter::m_ExperimentHigh
int m_ExperimentHigh
High experiment.
Definition: KLMDatabaseImporter.h:146
Belle2::EKLMAlignment
Class to store EKLM alignment data in the database.
Definition: EKLMAlignment.h:40
Belle2::KLMStripEfficiency
DBObject used to store the efficiencies of KLM strips.
Definition: KLMStripEfficiency.h:41
Belle2::KLMDatabaseImporter::importBKLMAlignment
void importBKLMAlignment(const BKLMAlignment *bklmAlignment, bool displacement=false)
Import BKLM alignment.
Definition: KLMDatabaseImporter.cc:135
Belle2::KLMDatabaseImporter::importTimeWindow
void importTimeWindow(KLMTimeWindow *timeWindow)
Import KLM time window parameters.
Definition: KLMDatabaseImporter.cc:75
Belle2::KLMDatabaseImporter::KLMDatabaseImporter
KLMDatabaseImporter()
Constructor.
Definition: KLMDatabaseImporter.cc:28
Belle2::KLMDatabaseImporter::m_RunHigh
int m_RunHigh
High run.
Definition: KLMDatabaseImporter.h:149
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KLMScintillatorDigitizationParameters
Class to store KLM scintillator simulation parameters in the database.
Definition: KLMScintillatorDigitizationParameters.h:33
Belle2::KLMTimeConversion
KLM time conversion.
Definition: KLMTimeConversion.h:33
Belle2::DBImportObjPtr
Class for importing a single object to the database.
Definition: DBImportObjPtr.h:33
Belle2::KLMDatabaseImporter::importScintillatorDigitizationParameters
void importScintillatorDigitizationParameters(const KLMScintillatorDigitizationParameters *digitizationParameters)
Import scintillator simulation parameters.
Definition: KLMDatabaseImporter.cc:55
Belle2::KLMDatabaseImporter::loadStripEfficiency
void loadStripEfficiency(KLMStripEfficiency *stripEfficiency, std::string fileName)
Load strip efficiencies.
Definition: KLMDatabaseImporter.cc:84
Belle2::KLMDatabaseImporter::m_RunLow
int m_RunLow
Low run.
Definition: KLMDatabaseImporter.h:143
Belle2::EKLMSegmentAlignment
Class to store EKLM alignment data in the database.
Definition: EKLMSegmentAlignment.h:40
Belle2::KLMDatabaseImporter::importChannelStatus
void importChannelStatus(const KLMChannelStatus *channelStatus)
Import channel status.
Definition: KLMDatabaseImporter.cc:45
Belle2::BKLMAlignment
Class to store BKLM alignment data in the database.
Definition: BKLMAlignment.h:40
Belle2::KLMDatabaseImporter::importAlignment
void importAlignment(const BKLMAlignment *bklmAlignment, const EKLMAlignment *eklmAlignment, const EKLMSegmentAlignment *eklmSegmentAlignment, bool displacement=false)
Import alignment.
Definition: KLMDatabaseImporter.cc:180
Belle2::KLMDatabaseImporter::~KLMDatabaseImporter
~KLMDatabaseImporter()
Destructor.
Definition: KLMDatabaseImporter.cc:32
Belle2::KLMDatabaseImporter::importEKLMSegmentAlignment
void importEKLMSegmentAlignment(const EKLMSegmentAlignment *eklmSegmentAlignment, bool displacement=false)
Import EKLM segment alignment.
Definition: KLMDatabaseImporter.cc:165