Belle II Software development
ECLDatabaseImporter Class Reference

ECL database importer. More...

#include <ECLDatabaseImporter.h>

Public Member Functions

 ECLDatabaseImporter (std::vector< std::string > inputFileNames, const std::string &m_name)
 Constructor.
 
virtual ~ECLDatabaseImporter ()
 Destructor.
 
void importDigitEnergyCalibration ()
 Import ECL energy calibration constants to the database.
 
void importDigitTimeCalibration ()
 Import ECL time calibration constants to the database.
 
void importShowerShapesSecondMomentCorrections ()
 Import ECL shower shape corrections to second moment to the database.
 
void importShowerCorrectorLeakageCorrections ()
 Import ECL leakage corrections to showers.
 
void importShowerEnergyCorrectionTemporary ()
 Import ECL corrections to showers energies.
 
void importTrackClusterMatchingThresholds ()
 Import threshold values for track ECL cluster matching to the database.
 
void importTrackClusterMatchingParameterizations ()
 Import parameterizations for the RMS between tracks and ECL clusters to the database.
 

Private Member Functions

template<class rootClass >
rootClass getRootObjectFromFile (TFile *file, const std::string &rootObjName) const
 Extract a root object rootObjName from a root file file.
 

Private Attributes

std::vector< std::string > m_inputFileNames
 Input file name.
 
std::string m_name
 Database object (output) file name.
 
std::string m_bgFractionBinNumBranchName = "bgFractionBinNum"
 Branch names for shower corrector leakage correction root file.
 
std::string m_regNumBranchName = "regNum"
 regNum branch name
 
std::string m_phiBinNumBranchName = "phiBinNum"
 phiBinNum branch name
 
std::string m_thetaBinNumBranchName = "thetaBinNum"
 thetaBinNum branch name
 
std::string m_energyBinNumBranchName = "energyBinNum"
 energyBinNum branch name
 
std::string m_correctionFactorBranchName = "correctionFactor"
 correctionFactor branch name
 
const int m_numAvgRecEnEntries = 15
 Number of entries in avgRecEn array.
 
std::string m_avgRecEnBranchName = "avgRecEn"
 avgRecEn branch name
 
std::string m_lReg1ThetaBranchName = "lReg1Theta"
 lReg1Theta branch name
 
std::string m_hReg1ThetaBranchName = "hReg1Theta"
 hReg1Theta branch name
 
std::string m_lReg2ThetaBranchName = "lReg2Theta"
 lReg2Theta branch name
 
std::string m_hReg2ThetaBranchName = "hReg2Theta"
 hReg2Theta branch name
 
std::string m_lReg3ThetaBranchName = "lReg3Theta"
 lReg3Theta branch name
 
std::string m_hReg3ThetaBranchName = "hReg3Theta"
 hReg3Theta branch name
 
std::string m_numOfBfBinsBranchName = "numOfBfBins"
 numOfBfBins branch name
 
std::string m_numOfEnergyBinsBranchName = "numOfEnergyBins"
 numOfEnergyBins branch name
 
std::string m_numOfPhiBinsBranchName = "numOfPhiBins"
 numOfPhiBins branch name
 
std::string m_numOfReg1ThetaBinsBranchName = "numOfReg1ThetaBins"
 numOfReg1ThetaBins branch name
 
std::string m_numOfReg2ThetaBinsBranchName = "numOfReg2ThetaBins"
 numOfReg2ThetaBins branch name

 
std::string m_numOfReg3ThetaBinsBranchName = "numOfReg3ThetaBins"
 numOfReg3ThetaBins branch name

 
std::string m_phiPeriodicityBranchName = "phiPeriodicity"
 phiPeriodicity branch name

 

Detailed Description

ECL database importer.

This module writes data from e.g. a ROOT histogram to database.

Definition at line 26 of file ECLDatabaseImporter.h.

Constructor & Destructor Documentation

◆ ECLDatabaseImporter()

ECLDatabaseImporter ( std::vector< std::string >  inputFileNames,
const std::string &  m_name 
)

Constructor.

Definition at line 44 of file ECLDatabaseImporter.cc.

45{
46 //input file names
47 for (auto& inputFileName : inputFileNames)
48 m_inputFileNames.push_back(inputFileName);
49
50 //output file name
51 m_name = name;
52}
std::vector< std::string > m_inputFileNames
Input file name.
std::string m_name
Database object (output) file name.

◆ ~ECLDatabaseImporter()

virtual ~ECLDatabaseImporter ( )
inlinevirtual

Destructor.

Definition at line 38 of file ECLDatabaseImporter.h.

38{};

Member Function Documentation

◆ getRootObjectFromFile()

rootClass getRootObjectFromFile ( TFile *  file,
const std::string &  rootObjName 
) const
inlineprivate

Extract a root object rootObjName from a root file file.

The file is assumed to be valid (pointer valid and not zombie). If rootObjName doesn't exist in file, do B2FATAL.

Definition at line 88 of file ECLDatabaseImporter.h.

89 {
90 rootClass rootObj = (rootClass)file->Get(rootObjName.data());
91 if (!rootObj) {
92 std::string filename = file->GetName();
93 delete file;
94 B2FATAL("Could not find " << rootObjName << " in " << filename);
95 }
96 return rootObj;
97 }

◆ importDigitEnergyCalibration()

void importDigitEnergyCalibration ( )

Import ECL energy calibration constants to the database.

Definition at line 54 of file ECLDatabaseImporter.cc.

55{
56
57 TClonesArray digitCalibrationConstants("Belle2::ECLDigitEnergyConstants");
58
59 TH1F* energy = 0;
60 TH1F* amplitude = 0;
61 int nFiles = 0;
62
63 for (const string& inputFile : m_inputFileNames) {
64
65 TFile* f = TFile::Open(inputFile.c_str(), "READ");
66
67 TIter next(f->GetListOfKeys());
68 TKey* key;
69 while ((key = (TKey*) next())) {
70
71 string histconstants = key->GetName();
72
73 if (histconstants.compare("energy") == 0) {
74 energy = (TH1F*)f->Get(histconstants.c_str());
75 } else if (histconstants.compare("amplitude") == 0) {
76 amplitude = (TH1F*)f->Get(histconstants.c_str());
77 }
78
79 else { B2FATAL("Key name does not match any of the following: energy, amplitude!"); }
80 }
81
82 nFiles++;
83 }
84
85 if (nFiles != 1) { B2FATAL("Sorry, you must only import one file at a time for now!"); }
86
87 if (!amplitude) B2FATAL("None of the input files contains the histogram called 'amplitude'!");
88 if (!energy) B2FATAL("None of the input files contains the histogram called 'energy'!");
89
90 // loop over the histogram to fill the TClonesArray
91 int cell = 0;
92 for (int bin = 1; bin <= amplitude->GetNbinsX(); ++bin) {
93 float amplitudeval = amplitude->GetBinContent(bin);
94 float energyval = energy->GetBinContent(bin);
95 new (digitCalibrationConstants[cell]) ECLDigitEnergyConstants(bin, amplitudeval, energyval);
96 cell++;
97 }
98
99 IntervalOfValidity iov(0, 0, -1, -1); // IOV (0,0,-1,-1) is valid for all runs and experiments
100 Database::Instance().storeData(m_name, &digitCalibrationConstants, iov);
101}
Energy calibration constants per digit.
A class that describes the interval of experiments/runs for which an object in the database is valid.
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:42
bool storeData(const std::string &name, TObject *object, const IntervalOfValidity &iov)
Store an object in the database.
Definition: Database.cc:141

◆ importDigitTimeCalibration()

void importDigitTimeCalibration ( )

Import ECL time calibration constants to the database.

Definition at line 103 of file ECLDatabaseImporter.cc.

104{
105
106 TClonesArray digitCalibrationConstants("Belle2::ECLDigitTimeConstants");
107
108 TH1F* offset = 0;
109 int nFiles = 0;
110
111 for (const string& inputFile : m_inputFileNames) {
112
113 TFile* f = TFile::Open(inputFile.c_str(), "READ");
114
115 TIter next(f->GetListOfKeys());
116 TKey* key;
117 while ((key = (TKey*) next())) {
118
119 string histconstants = key->GetName();
120
121 if (histconstants.compare("constantB") == 0) {
122 offset = (TH1F*)f->Get(histconstants.c_str());
123 } else { B2FATAL("Key name does not match any of the following: constantC!"); }
124 }
125
126 nFiles++;
127 }
128
129 if (nFiles != 1) { B2FATAL("Sorry, you must only import one file at a time for now!"); }
130
131 if (!offset) B2FATAL("None of the input files contains the histogram called 'constantB'!");
132
133 // loop over the histogram to fill the TClonesArray
134 int cell = 0;
135 for (int bin = 1; bin <= offset->GetNbinsX(); ++bin) {
136 float offsetval = offset->GetBinContent(bin);
137 new (digitCalibrationConstants[cell]) ECLDigitTimeConstants(bin, offsetval);
138 cell++;
139 }
140
141 IntervalOfValidity iov(0, 0, -1, -1); // IOV (0,0,-1,-1) is valid for all runs and experiments
142 Database::Instance().storeData(m_name, &digitCalibrationConstants, iov);
143}
Time and time resolution calibration constants per digit.

◆ importShowerCorrectorLeakageCorrections()

void importShowerCorrectorLeakageCorrections ( )

Import ECL leakage corrections to showers.

Definition at line 145 of file ECLDatabaseImporter.cc.

146{
147 if (m_inputFileNames.size() > 1)
148 B2FATAL("Sorry, you must only import one file at a time for now!");
149
150 //Open file
151 TFile* inputFile = new TFile(m_inputFileNames[0].data(), "READ");
152
153 if (!inputFile || inputFile->IsZombie())
154 B2FATAL("Could not open file " << m_inputFileNames[0]);
155
156 //Get trees
157 TTree* correctionTree = getRootObjectFromFile<TTree*>(inputFile, "ParameterNtuple");
158 TTree* helperTree = getRootObjectFromFile<TTree*>(inputFile, "ConstantNtuple");
159
160 //----------------------------------------------------------------------------------------------
161 //Fill ParameterNtuple vectors
162 //----------------------------------------------------------------------------------------------
163
164 int bgFractionBinNum;
165 int regNum;
166 int phiBinNum;
167 int thetaBinNum;
168 int energyBinNum;
169 float correctionFactor;
170
171 //Set Branch Addresses
172 correctionTree->SetBranchAddress(m_bgFractionBinNumBranchName.c_str(), &bgFractionBinNum);
173 correctionTree->SetBranchAddress(m_regNumBranchName.c_str(), &regNum);
174 correctionTree->SetBranchAddress(m_phiBinNumBranchName.c_str(), &phiBinNum);
175 correctionTree->SetBranchAddress(m_thetaBinNumBranchName.c_str(), &thetaBinNum);
176 correctionTree->SetBranchAddress(m_energyBinNumBranchName.c_str(), &energyBinNum);
177 correctionTree->SetBranchAddress(m_correctionFactorBranchName.c_str(), &correctionFactor);
178
179 //Fill vectors
180 std::vector<int> m_bgFractionBinNum;
181 std::vector<int> m_regNum;
182 std::vector<int> m_phiBinNum;
183 std::vector<int> m_thetaBinNum;
184 std::vector<int> m_energyBinNum;
185 std::vector<float> m_correctionFactor;
186
187 for (long iEntry = 0; iEntry < correctionTree->GetEntries(); ++iEntry) {
188 correctionTree->GetEntry(iEntry);
189
190 m_bgFractionBinNum.push_back(bgFractionBinNum);
191 m_regNum.push_back(regNum);
192 m_phiBinNum.push_back(phiBinNum);
193 m_thetaBinNum.push_back(thetaBinNum);
194 m_energyBinNum.push_back(energyBinNum);
195 m_correctionFactor.push_back(correctionFactor);
196 }
197
198 //----------------------------------------------------------------------------------------------
199 //Fill ConstantNtuple vectors
200 //----------------------------------------------------------------------------------------------
201
202
203 float lReg1Theta;
204 float hReg1Theta;
205 float lReg2Theta;
206 float hReg2Theta;
207 float lReg3Theta;
208 float hReg3Theta;
209 int numOfBfBins;
210 int numOfEnergyBins;
211 int numOfPhiBins;
212 int numOfReg1ThetaBins;
213 int numOfReg2ThetaBins;
214 int numOfReg3ThetaBins;
215 int phiPeriodicity;
216
217 //Ugly hack to circumvent 'stack usage might be unbounded [-Wstack-usage=]' compiler warning that's caused by the use of c-type arrays.
218 //This is not for the faint of heart
219
220 //because root GetEntry fills the whole internal array of the vector without changing it's size, we must ensure that it's the right size.
221 std::vector<float> avgRecEns(m_numAvgRecEnEntries);
222 helperTree->SetBranchAddress(m_avgRecEnBranchName.c_str(), avgRecEns.data()); //Read c-style array right into internal vector array.
223
224 helperTree->SetBranchAddress(m_lReg1ThetaBranchName.c_str(), &lReg1Theta);
225 helperTree->SetBranchAddress(m_hReg1ThetaBranchName.c_str(), &hReg1Theta);
226 helperTree->SetBranchAddress(m_lReg2ThetaBranchName.c_str(), &lReg2Theta);
227 helperTree->SetBranchAddress(m_hReg2ThetaBranchName.c_str(), &hReg2Theta);
228 helperTree->SetBranchAddress(m_lReg3ThetaBranchName.c_str(), &lReg3Theta);
229 helperTree->SetBranchAddress(m_hReg3ThetaBranchName.c_str(), &hReg3Theta);
230 helperTree->SetBranchAddress(m_numOfBfBinsBranchName.c_str(), &numOfBfBins);
231 helperTree->SetBranchAddress(m_numOfEnergyBinsBranchName.c_str(), &numOfEnergyBins);
232 helperTree->SetBranchAddress(m_numOfPhiBinsBranchName.c_str(), &numOfPhiBins);
233 helperTree->SetBranchAddress(m_numOfReg1ThetaBinsBranchName.c_str(), &numOfReg1ThetaBins);
234 helperTree->SetBranchAddress(m_numOfReg2ThetaBinsBranchName.c_str(), &numOfReg2ThetaBins);
235 helperTree->SetBranchAddress(m_numOfReg3ThetaBinsBranchName.c_str(), &numOfReg3ThetaBins);
236 helperTree->SetBranchAddress(m_phiPeriodicityBranchName.c_str(), &phiPeriodicity);
237
238 //Fill vectors
239 std::vector<float> m_avgRecEn;
240 std::vector<float> m_lReg1Theta;
241 std::vector<float> m_hReg1Theta;
242 std::vector<float> m_lReg2Theta;
243 std::vector<float> m_hReg2Theta;
244 std::vector<float> m_lReg3Theta;
245 std::vector<float> m_hReg3Theta;
246 std::vector<int> m_numOfBfBins;
247 std::vector<int> m_numOfEnergyBins;
248 std::vector<int> m_numOfPhiBins;
249 std::vector<int> m_numOfReg1ThetaBins;
250 std::vector<int> m_numOfReg2ThetaBins;
251 std::vector<int> m_numOfReg3ThetaBins;
252 std::vector<int> m_phiPeriodicity;
253
254 for (long iEntry = 0; iEntry < helperTree->GetEntries(); ++iEntry) {
255 helperTree->GetEntry(iEntry);
256 for (unsigned int iIdx = 0; iIdx < avgRecEns.size(); ++iIdx) m_avgRecEn.push_back(avgRecEns[iIdx]);
257
258 m_lReg1Theta.push_back(lReg1Theta);
259 m_hReg1Theta.push_back(hReg1Theta);
260 m_lReg2Theta.push_back(lReg2Theta);
261 m_hReg2Theta.push_back(hReg2Theta);
262 m_lReg3Theta.push_back(lReg3Theta);
263 m_hReg3Theta.push_back(hReg3Theta);
264 m_numOfBfBins.push_back(numOfBfBins);
265 m_numOfEnergyBins.push_back(numOfEnergyBins);
266 m_numOfPhiBins.push_back(numOfPhiBins);
267 m_numOfReg1ThetaBins.push_back(numOfReg1ThetaBins);
268 m_numOfReg2ThetaBins.push_back(numOfReg2ThetaBins);
269 m_numOfReg3ThetaBins.push_back(numOfReg3ThetaBins);
270 m_phiPeriodicity.push_back(phiPeriodicity);
271 }
272
273 //----------------------------------------------------------------------------------------------
274
275 //Construct DB object
276 DBImportObjPtr<ECLShowerCorrectorLeakageCorrection> dbPtr("ecl_shower_corrector_leakage_corrections");
277 dbPtr.construct(m_bgFractionBinNum,
278 m_regNum,
279 m_phiBinNum,
280 m_thetaBinNum,
281 m_energyBinNum,
282 m_correctionFactor,
283 m_avgRecEn,
284 m_lReg1Theta,
285 m_hReg1Theta,
286 m_lReg2Theta,
287 m_hReg2Theta,
288 m_lReg3Theta,
289 m_hReg3Theta,
290 m_numOfBfBins,
291 m_numOfEnergyBins,
292 m_numOfPhiBins,
293 m_numOfReg1ThetaBins,
294 m_numOfReg2ThetaBins,
295 m_numOfReg3ThetaBins,
296 m_phiPeriodicity);
297
298 //Create IOV object
299 int startExp = 0;
300 int startRun = 0;
301 int endExp = -1;
302 int endRun = -1;
303 IntervalOfValidity iov(startExp, startRun, endExp, endRun);
304
305 //Import into local db
306 dbPtr.import(iov);
307
308 delete inputFile;
309
310}
Class for importing a single object to the database.
std::string m_phiBinNumBranchName
phiBinNum branch name
std::string m_numOfPhiBinsBranchName
numOfPhiBins branch name
std::string m_avgRecEnBranchName
avgRecEn branch name
std::string m_numOfReg3ThetaBinsBranchName
numOfReg3ThetaBins branch name
std::string m_numOfReg2ThetaBinsBranchName
numOfReg2ThetaBins branch name
const int m_numAvgRecEnEntries
Number of entries in avgRecEn array.
std::string m_hReg2ThetaBranchName
hReg2Theta branch name
std::string m_hReg3ThetaBranchName
hReg3Theta branch name
std::string m_thetaBinNumBranchName
thetaBinNum branch name
std::string m_lReg2ThetaBranchName
lReg2Theta branch name
std::string m_bgFractionBinNumBranchName
Branch names for shower corrector leakage correction root file.
std::string m_numOfEnergyBinsBranchName
numOfEnergyBins branch name
std::string m_lReg3ThetaBranchName
lReg3Theta branch name
std::string m_correctionFactorBranchName
correctionFactor branch name
std::string m_numOfBfBinsBranchName
numOfBfBins branch name
std::string m_regNumBranchName
regNum branch name
std::string m_hReg1ThetaBranchName
hReg1Theta branch name
std::string m_energyBinNumBranchName
energyBinNum branch name
std::string m_phiPeriodicityBranchName
phiPeriodicity branch name
std::string m_lReg1ThetaBranchName
lReg1Theta branch name
std::string m_numOfReg1ThetaBinsBranchName
numOfReg1ThetaBins branch name

◆ importShowerEnergyCorrectionTemporary()

void importShowerEnergyCorrectionTemporary ( )

Import ECL corrections to showers energies.

Temperary - there will be additional improvements done to these corrections in the future. Input file should be .txt file and have the format: generated-energy bkg-scale-factor theta-min theta-max corr-factor. The numbers should be seperated by spaces. For each line, the theta value used is the average of theta-min and theta-max

Definition at line 354 of file ECLDatabaseImporter.cc.

355{
356 if (m_inputFileNames.size() > 1)
357 B2FATAL("Sorry, you must only import one file at a time for now!");
358
359 //Expect a root file
360 std::filesystem::path path(m_inputFileNames[0]);
361 if (path.extension() != ".root")
362 B2FATAL("Expecting a .root file. Aborting");
363
364 TFile* inputFile = new TFile(m_inputFileNames[0].data(), "READ");
365
366 TGraph2D* theta_geo_graph = getRootObjectFromFile<TGraph2D*>(inputFile, "LeakageCorrections_theta_geometry");
367 TGraph2D* phi_geo_graph = getRootObjectFromFile<TGraph2D*>(inputFile, "LeakageCorrections_phi_geometry");
368 TGraph2D* theta_en_graph = getRootObjectFromFile<TGraph2D*>(inputFile, "LeakageCorrections_theta_energy");
369 TGraph2D* phi_en_graph = getRootObjectFromFile<TGraph2D*>(inputFile, "LeakageCorrections_phi_energy");
370 TH1F* bg_histo = getRootObjectFromFile<TH1F*>(inputFile, "LeakageCorrections_background_fraction");
371
372 double bkgFactor = bg_histo->GetBinContent(1);
373
374 double thetaMin = theta_en_graph->GetXmin();
375 double thetaMax = theta_en_graph->GetXmax();
376 double phiMin = phi_en_graph->GetXmin();
377 double phiMax = phi_en_graph->GetXmax();
378
379 double energyMin = theta_en_graph->GetYmin();
380 double energyMax = theta_en_graph->GetYmax();
381
382
383 B2DEBUG(28, "Leakage DBobjects angle boundaries: thetaMin=" << thetaMin << " thetaMax=" << thetaMax << " phiMin= " << phiMin <<
384 " phiMax= " << phiMax << " enmin=" << energyMin <<
385 " enmax=" << energyMax);
386
387 // Import to DB
388 int startExp = 0;
389 int startRun = 0;
390 int endExp = -1;
391 int endRun = -1;
392 IntervalOfValidity iov(startExp, startRun, endExp, endRun);
393
394 if (std::abs(bkgFactor - 1.0) < 1e-9) { //bkgFactor == 1 -> phase 2 backgrounds
395
396 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_theta_geo("ECLLeakageCorrection_thetaGeometry_phase2");
397 dbPtr_theta_geo.construct(*theta_geo_graph, thetaMin, thetaMax, energyMin, energyMax);
398 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_phi_geo("ECLLeakageCorrection_phiGeometry_phase2");
399 dbPtr_phi_geo.construct(*phi_geo_graph, phiMin, phiMax, energyMin, energyMax);
400 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_theta_en("ECLLeakageCorrection_thetaEnergy_phase2");
401 dbPtr_theta_en.construct(*theta_en_graph, thetaMin, thetaMax, energyMin, energyMax);
402 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_phi_en("ECLLeakageCorrection_phiEnergy_phase2");
403 dbPtr_phi_en.construct(*phi_en_graph, phiMin, phiMax, energyMin, energyMax);
404
405 //Import into local db
406 dbPtr_theta_geo.import(iov);
407 dbPtr_phi_geo.import(iov);
408 dbPtr_theta_en.import(iov);
409 dbPtr_phi_en.import(iov);
410 }
411 /*else (because currently phase_2 and phase_3 are same payload*/ if (std::abs(bkgFactor - 1.0) < 1e-9) {
412
413 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_theta_geo("ECLLeakageCorrection_thetaGeometry_phase3");
414 dbPtr_theta_geo.construct(*theta_geo_graph, thetaMin, thetaMax, energyMin, energyMax);
415 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_phi_geo("ECLLeakageCorrection_phiGeometry_phase3");
416 dbPtr_phi_geo.construct(*phi_geo_graph, phiMin, phiMax, energyMin, energyMax);
417 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_theta_en("ECLLeakageCorrection_thetaEnergy_phase3");
418 dbPtr_theta_en.construct(*theta_en_graph, thetaMin, thetaMax, energyMin, energyMax);
419 DBImportObjPtr<ECLShowerEnergyCorrectionTemporary> dbPtr_phi_en("ECLLeakageCorrection_phiEnergy_phase3");
420 dbPtr_phi_en.construct(*phi_en_graph, phiMin, phiMax, energyMin, energyMax);
421
422 //Import into local db
423 dbPtr_theta_geo.import(iov);
424 dbPtr_phi_geo.import(iov);
425 dbPtr_theta_en.import(iov);
426 dbPtr_phi_en.import(iov);
427 }
428
429
430}

◆ importShowerShapesSecondMomentCorrections()

void importShowerShapesSecondMomentCorrections ( )

Import ECL shower shape corrections to second moment to the database.

Definition at line 312 of file ECLDatabaseImporter.cc.

313{
314 if (m_inputFileNames.size() > 1)
315 B2FATAL("Sorry, you must only import one file at a time for now!");
316
317 //Extract TGraphs from file
318 DBImportArray<ECLShowerShapeSecondMomentCorrection> dbArray("ecl_shower_shape_second_moment_corrections");
319 TFile* inputFile = new TFile(m_inputFileNames[0].data(), "READ");
320
321 if (!inputFile || inputFile->IsZombie())
322 B2FATAL("Could not open file " << m_inputFileNames[0]);
323
324 //N1 theta
325 TGraph* theta_N1_graph = getRootObjectFromFile<TGraph*>(inputFile, "SecondMomentCorrections_theta_N1");
326 dbArray.appendNew(ECLShower::c_nPhotons, ECLShowerShapeModule::c_thetaType, *theta_N1_graph);
327
328 //N1 phi
329 TGraph* phi_N1_graph = getRootObjectFromFile<TGraph*>(inputFile, "SecondMomentCorrections_phi_N1");
330 dbArray.appendNew(ECLShower::c_nPhotons, ECLShowerShapeModule::c_phiType, *phi_N1_graph);
331
332 //N2 theta
333 TGraph* theta_N2_graph = getRootObjectFromFile<TGraph*>(inputFile, "SecondMomentCorrections_theta_N2");
334 dbArray.appendNew(ECLShower::c_neutralHadron, ECLShowerShapeModule::c_thetaType, *theta_N2_graph);
335
336 //N2 phi
337 TGraph* phi_N2_graph = getRootObjectFromFile<TGraph*>(inputFile, "SecondMomentCorrections_phi_N2");
338 dbArray.appendNew(ECLShower::c_neutralHadron, ECLShowerShapeModule::c_phiType, *phi_N2_graph);
339
340
341 //Import to DB
342 int startExp = 0;
343 int startRun = 0;
344 int endExp = -1;
345 int endRun = -1;
346 IntervalOfValidity iov(startExp, startRun, endExp, endRun);
347
348 //Import into local db
349 dbArray.import(iov);
350
351 delete inputFile;
352}
Class for importing array of objects to the database.
Definition: DBImportArray.h:25
@ c_phiType
type of phi identifier
@ c_thetaType
type of theta identifier
@ c_neutralHadron
CR is reconstructed as a neutral hadron (N2)
Definition: ECLShower.h:44
@ c_nPhotons
CR is split into n photons (N1)
Definition: ECLShower.h:42

◆ importTrackClusterMatchingParameterizations()

void importTrackClusterMatchingParameterizations ( )

Import parameterizations for the RMS between tracks and ECL clusters to the database.

Definition at line 477 of file ECLDatabaseImporter.cc.

478{
479 if (m_inputFileNames.size() > 1)
480 B2FATAL("Sorry, you must only import one file at a time for now!");
481
482 // Open file
483 TFile* inputFile = new TFile(m_inputFileNames[0].data(), "READ");
484
485 if (!inputFile || inputFile->IsZombie())
486 B2FATAL("Could not open file " << m_inputFileNames[0]);
487
488 map<string, TF1> m_parametrizationFunctions;
489 vector<string> angles = {"Theta", "Phi"};
490 vector<string> regions = {"BRL", "BWD", "FWD"};
491 vector<string> hittypes = {"CROSS", "DL", "NEAR"};
492
493 for (const auto& angle : angles) {
494 for (const auto& region : regions) {
495 for (const auto& hittype : hittypes) {
496 m_parametrizationFunctions.insert(make_pair(angle + region + hittype, *(getRootObjectFromFile<TF1*>(inputFile,
497 "RMSParameterization" + angle + region + hittype))));
498 }
499 }
500 }
501
502 DBImportObjPtr<ECLTrackClusterMatchingParameterizations> dbPtr("ECLTrackClusterMatchingParameterizations");
503 dbPtr.construct(m_parametrizationFunctions);
504
505 IntervalOfValidity iov(0, 0, -1, -1);
506
507 //Import into local db
508 dbPtr.import(iov);
509
510 delete inputFile;
511}

◆ importTrackClusterMatchingThresholds()

void importTrackClusterMatchingThresholds ( )

Import threshold values for track ECL cluster matching to the database.

Definition at line 432 of file ECLDatabaseImporter.cc.

433{
434 if (m_inputFileNames.size() > 1)
435 B2FATAL("Sorry, you must only import one file at a time for now!");
436
437 //Expect a txt file
438 std::filesystem::path path(m_inputFileNames[0]);
439 if (path.extension() != ".txt")
440 B2FATAL("Expecting a .txt file. Aborting");
441
442 vector<pair<double, double>> m_matchingThresholdPairsFWD;
443 vector<pair<double, double>> m_matchingThresholdPairsBWD;
444 vector<pair<double, pair<double, double>>> m_matchingThresholdPairsBRL;
445 pair<double, double> m_matchingThresholdPair;
446 pair<double, pair<double, double>> m_thetaMatchingThresholdPair;
447 double pt, threshold, thetalimit;
448 string eclregion;
449
450 ifstream infile(m_inputFileNames[0]);
451 string line;
452 while (getline(infile, line)) {
453 istringstream iss(line);
454 iss >> eclregion;
455 if (eclregion == "FWD" || eclregion == "BWD") {
456 iss >> pt >> threshold;
457 m_matchingThresholdPair = make_pair(pt, threshold);
458 if (eclregion == "FWD") m_matchingThresholdPairsFWD.push_back(m_matchingThresholdPair);
459 else m_matchingThresholdPairsBWD.push_back(m_matchingThresholdPair);
460 } else if (eclregion == "BRL") {
461 iss >> thetalimit >> pt >> threshold;
462 m_matchingThresholdPair = make_pair(pt, threshold);
463 m_thetaMatchingThresholdPair = make_pair(thetalimit, m_matchingThresholdPair);
464 m_matchingThresholdPairsBRL.push_back(m_thetaMatchingThresholdPair);
465 }
466 }
467
468 DBImportObjPtr<ECLTrackClusterMatchingThresholds> dbPtr("ECLTrackClusterMatchingThresholds");
469 dbPtr.construct(m_matchingThresholdPairsFWD, m_matchingThresholdPairsBWD, m_matchingThresholdPairsBRL);
470
471 IntervalOfValidity iov(0, 0, -1, -1);
472
473 //Import into local db
474 dbPtr.import(iov);
475}

Member Data Documentation

◆ m_avgRecEnBranchName

std::string m_avgRecEnBranchName = "avgRecEn"
private

avgRecEn branch name

Definition at line 113 of file ECLDatabaseImporter.h.

◆ m_bgFractionBinNumBranchName

std::string m_bgFractionBinNumBranchName = "bgFractionBinNum"
private

Branch names for shower corrector leakage correction root file.

bgFractionBinNum branch name

Definition at line 104 of file ECLDatabaseImporter.h.

◆ m_correctionFactorBranchName

std::string m_correctionFactorBranchName = "correctionFactor"
private

correctionFactor branch name

Definition at line 109 of file ECLDatabaseImporter.h.

◆ m_energyBinNumBranchName

std::string m_energyBinNumBranchName = "energyBinNum"
private

energyBinNum branch name

Definition at line 108 of file ECLDatabaseImporter.h.

◆ m_hReg1ThetaBranchName

std::string m_hReg1ThetaBranchName = "hReg1Theta"
private

hReg1Theta branch name

Definition at line 115 of file ECLDatabaseImporter.h.

◆ m_hReg2ThetaBranchName

std::string m_hReg2ThetaBranchName = "hReg2Theta"
private

hReg2Theta branch name

Definition at line 117 of file ECLDatabaseImporter.h.

◆ m_hReg3ThetaBranchName

std::string m_hReg3ThetaBranchName = "hReg3Theta"
private

hReg3Theta branch name

Definition at line 119 of file ECLDatabaseImporter.h.

◆ m_inputFileNames

std::vector<std::string> m_inputFileNames
private

Input file name.

Definition at line 82 of file ECLDatabaseImporter.h.

◆ m_lReg1ThetaBranchName

std::string m_lReg1ThetaBranchName = "lReg1Theta"
private

lReg1Theta branch name

Definition at line 114 of file ECLDatabaseImporter.h.

◆ m_lReg2ThetaBranchName

std::string m_lReg2ThetaBranchName = "lReg2Theta"
private

lReg2Theta branch name

Definition at line 116 of file ECLDatabaseImporter.h.

◆ m_lReg3ThetaBranchName

std::string m_lReg3ThetaBranchName = "lReg3Theta"
private

lReg3Theta branch name

Definition at line 118 of file ECLDatabaseImporter.h.

◆ m_name

std::string m_name
private

Database object (output) file name.

Definition at line 83 of file ECLDatabaseImporter.h.

◆ m_numAvgRecEnEntries

const int m_numAvgRecEnEntries = 15
private

Number of entries in avgRecEn array.

If this is wrong bad things will happen

Definition at line 112 of file ECLDatabaseImporter.h.

◆ m_numOfBfBinsBranchName

std::string m_numOfBfBinsBranchName = "numOfBfBins"
private

numOfBfBins branch name

Definition at line 120 of file ECLDatabaseImporter.h.

◆ m_numOfEnergyBinsBranchName

std::string m_numOfEnergyBinsBranchName = "numOfEnergyBins"
private

numOfEnergyBins branch name

Definition at line 121 of file ECLDatabaseImporter.h.

◆ m_numOfPhiBinsBranchName

std::string m_numOfPhiBinsBranchName = "numOfPhiBins"
private

numOfPhiBins branch name

Definition at line 122 of file ECLDatabaseImporter.h.

◆ m_numOfReg1ThetaBinsBranchName

std::string m_numOfReg1ThetaBinsBranchName = "numOfReg1ThetaBins"
private

numOfReg1ThetaBins branch name

Definition at line 123 of file ECLDatabaseImporter.h.

◆ m_numOfReg2ThetaBinsBranchName

std::string m_numOfReg2ThetaBinsBranchName = "numOfReg2ThetaBins"
private

numOfReg2ThetaBins branch name

Definition at line 124 of file ECLDatabaseImporter.h.

◆ m_numOfReg3ThetaBinsBranchName

std::string m_numOfReg3ThetaBinsBranchName = "numOfReg3ThetaBins"
private

numOfReg3ThetaBins branch name

Definition at line 125 of file ECLDatabaseImporter.h.

◆ m_phiBinNumBranchName

std::string m_phiBinNumBranchName = "phiBinNum"
private

phiBinNum branch name

Definition at line 106 of file ECLDatabaseImporter.h.

◆ m_phiPeriodicityBranchName

std::string m_phiPeriodicityBranchName = "phiPeriodicity"
private

phiPeriodicity branch name

Definition at line 126 of file ECLDatabaseImporter.h.

◆ m_regNumBranchName

std::string m_regNumBranchName = "regNum"
private

regNum branch name

Definition at line 105 of file ECLDatabaseImporter.h.

◆ m_thetaBinNumBranchName

std::string m_thetaBinNumBranchName = "thetaBinNum"
private

thetaBinNum branch name

Definition at line 107 of file ECLDatabaseImporter.h.


The documentation for this class was generated from the following files: