Belle II Software  release-08-01-10
ECLDatabaseImporter Class Reference

ECL database importer. More...

#include <ECLDatabaseImporter.h>

Collaboration diagram for ECLDatabaseImporter:

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. More...
 
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. More...
 

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. More...
 
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. More...
 
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.

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  }

◆ 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 349 of file ECLDatabaseImporter.cc.

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

Member Data Documentation

◆ 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_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.


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