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

KLM time calibration algorithm. More...

#include <KLMTimeAlgorithm.h>

Inheritance diagram for KLMTimeAlgorithm:
Collaboration diagram for KLMTimeAlgorithm:

Classes

struct  Event
 Event data. More...
 

Public Types

enum  ChannelCalibrationStatus {
  c_NotEnoughData = 0 ,
  c_FailedFit = 1 ,
  c_SuccessfulCalibration = 2
}
 Channel calibration status.
 
enum  EResult {
  c_OK ,
  c_Iterate ,
  c_NotEnoughData ,
  c_Failure ,
  c_Undefined
}
 The result of calibration. More...
 

Public Member Functions

 KLMTimeAlgorithm ()
 Constructor.
 
 ~KLMTimeAlgorithm ()
 Destructor.
 
void setDebug ()
 Turn on debug mode (prints histograms and output running log).
 
void setMC (bool mc)
 Set flag indicating whether the input is MC sample. More...
 
void useEvtT0 ()
 Use event T0 as the initial time point or not.
 
void setMinimalDigitNumber (int minimalDigitNumber)
 Set minimal digit number (total).
 
void setLowerLimit (int counts)
 Set the lower number of hits collected on one sigle strip. More...
 
void saveHist ()
 Save histograms to file.
 
double esti_timeShift (const KLMChannelIndex &klmChannel)
 Estimate value of calibration constant for uncalibrated channels. More...
 
std::pair< int, double > tS_upperStrip (const KLMChannelIndex &klmChannel)
 Tracing avaiable channels with increasing strip number. More...
 
std::pair< int, double > tS_lowerStrip (const KLMChannelIndex &klmChannel)
 Tracing avaiable channels with decreasing strip number. More...
 
double esti_timeRes (const KLMChannelIndex &klmChannel)
 Estimate value of calibration constant for calibrated channels. More...
 
std::pair< int, double > tR_upperStrip (const KLMChannelIndex &klmChannel)
 Tracing avaiable channels with increasing strip number. More...
 
std::pair< int, double > tR_lowerStrip (const KLMChannelIndex &klmChannel)
 Tracing avaiable channels with decreasing strip number. More...
 
std::string getPrefix () const
 Get the prefix used for getting calibration data.
 
bool checkPyExpRun (PyObject *pyObj)
 Checks that a PyObject can be successfully converted to an ExpRun type. More...
 
Calibration::ExpRun convertPyExpRun (PyObject *pyObj)
 Performs the conversion of PyObject to ExpRun. More...
 
std::string getCollectorName () const
 Alias for prefix. More...
 
void setPrefix (const std::string &prefix)
 Set the prefix used to identify datastore objects.
 
void setInputFileNames (PyObject *inputFileNames)
 Set the input file names used for this algorithm from a Python list. More...
 
PyObject * getInputFileNames ()
 Get the input file names used for this algorithm and pass them out as a Python list of unicode strings.
 
std::vector< Calibration::ExpRun > getRunListFromAllData () const
 Get the complete list of runs from inspection of collected data.
 
RunRange getRunRangeFromAllData () const
 Get the complete RunRange from inspection of collected data.
 
IntervalOfValidity getIovFromAllData () const
 Get the complete IoV from inspection of collected data.
 
void fillRunToInputFilesMap ()
 Fill the mapping of ExpRun -> Files.
 
std::string getGranularity () const
 Get the granularity of collected data.
 
EResult execute (std::vector< Calibration::ExpRun > runs={}, int iteration=0, IntervalOfValidity iov=IntervalOfValidity())
 Runs calibration over vector of runs for a given iteration. More...
 
EResult execute (PyObject *runs, int iteration=0, IntervalOfValidity iov=IntervalOfValidity())
 Runs calibration over Python list of runs. Converts to C++ and then calls the other execute() function.
 
std::list< Database::DBImportQuery > & getPayloads ()
 Get constants (in TObjects) for database update from last execution.
 
std::list< Database::DBImportQuerygetPayloadValues ()
 Get constants (in TObjects) for database update from last execution but passed by VALUE.
 
bool commit ()
 Submit constants from last calibration into database.
 
bool commit (std::list< Database::DBImportQuery > payloads)
 Submit constants from a (potentially previous) set of payloads.
 
const std::string & getDescription () const
 Get the description of the algoithm (set by developers in constructor)
 
bool loadInputJson (const std::string &jsonString)
 Load the m_inputJson variable from a string (useful from Python interface). The rturn bool indicates success or failure.
 
const std::string dumpOutputJson () const
 Dump the JSON string of the output JSON object.
 
const std::vector< Calibration::ExpRun > findPayloadBoundaries (std::vector< Calibration::ExpRun > runs, int iteration=0)
 Used to discover the ExpRun boundaries that you want the Python CAF to execute on. This is optional and only used in some.
 
template<>
std::shared_ptr< TTree > getObjectPtr (const std::string &name, const std::vector< Calibration::ExpRun > &requestedRuns)
 Specialization of getObjectPtr<TTree>.
 

Protected Member Functions

virtual EResult calibrate () override
 Run algorithm on data. More...
 
void setInputFileNames (std::vector< std::string > inputFileNames)
 Set the input file names used for this algorithm. More...
 
virtual bool isBoundaryRequired (const Calibration::ExpRun &)
 Given the current collector data, make a decision about whether or not this run should be the start of a payload boundary.
 
virtual void boundaryFindingSetup (std::vector< Calibration::ExpRun >, int)
 If you need to make some changes to your algorithm class before 'findPayloadBoundaries' is run, make them in this function.
 
virtual void boundaryFindingTearDown ()
 Put your algorithm back into a state ready for normal execution if you need to.
 
const std::vector< Calibration::ExpRun > & getRunList () const
 Get the list of runs for which calibration is called.
 
int getIteration () const
 Get current iteration.
 
std::vector< std::string > getVecInputFileNames () const
 Get the input file names used for this algorithm as a STL vector.
 
template<class T >
std::shared_ptr< T > getObjectPtr (const std::string &name, const std::vector< Calibration::ExpRun > &requestedRuns)
 Get calibration data object by name and list of runs, the Merge function will be called to generate the overall object.
 
template<class T >
std::shared_ptr< T > getObjectPtr (std::string name)
 Get calibration data object (for all runs the calibration is requested for) This function will only work during or after execute() has been called once.
 
template<>
shared_ptr< TTree > getObjectPtr (const string &name, const vector< ExpRun > &requestedRuns)
 We cheekily cast the TChain to TTree for the returned pointer so that the user never knows Hopefully this doesn't cause issues if people do low level stuff to the tree...
 
std::string getGranularityFromData () const
 Get the granularity of collected data.
 
void saveCalibration (TClonesArray *data, const std::string &name)
 Store DBArray payload with given name with default IOV.
 
void saveCalibration (TClonesArray *data, const std::string &name, const IntervalOfValidity &iov)
 Store DBArray with given name and custom IOV.
 
void saveCalibration (TObject *data)
 Store DB payload with default name and default IOV.
 
void saveCalibration (TObject *data, const IntervalOfValidity &iov)
 Store DB payload with default name and custom IOV.
 
void saveCalibration (TObject *data, const std::string &name)
 Store DB payload with given name with default IOV.
 
void saveCalibration (TObject *data, const std::string &name, const IntervalOfValidity &iov)
 Store DB payload with given name and custom IOV.
 
void updateDBObjPtrs (const unsigned int event, const int run, const int experiment)
 Updates any DBObjPtrs by calling update(event) for DBStore.
 
void setDescription (const std::string &description)
 Set algorithm description (in constructor)
 
void clearCalibrationData ()
 Clear calibration data.
 
Calibration::ExpRun getAllGranularityExpRun () const
 Returns the Exp,Run pair that means 'Everything'. Currently unused.
 
void resetInputJson ()
 Clears the m_inputJson member variable.
 
void resetOutputJson ()
 Clears the m_outputJson member variable.
 
template<class T >
void setOutputJsonValue (const std::string &key, const T &value)
 Set a key:value pair for the outputJson object, expected to used interally during calibrate()
 
template<class T >
const T getOutputJsonValue (const std::string &key) const
 Get a value using a key from the JSON output object, not sure why you would want to do this.
 
template<class T >
const T getInputJsonValue (const std::string &key) const
 Get an input JSON value using a key. The normal exceptions are raised when the key doesn't exist.
 
const nlohmann::json & getInputJsonObject () const
 Get the entire top level JSON object. We explicitly say this must be of object type so that we might pick.
 
bool inputJsonKeyExists (const std::string &key) const
 Test for a key in the input JSON object.
 

Protected Attributes

std::vector< Calibration::ExpRun > m_boundaries
 When using the boundaries functionality from isBoundaryRequired, this is used to store the boundaries. It is cleared when.
 

Private Member Functions

void setupDatabase ()
 Setup the database.
 
CalibrationAlgorithm::EResult readCalibrationData ()
 Read calibration data. More...
 
void createHistograms ()
 Create histograms. More...
 
void fillTimeDistanceProfiles (TProfile *profileRpcPhi, TProfile *profileRpcZ, TProfile *profileBKLMScintillatorPhi, TProfile *profileBKLMScintillatorZ, TProfile *profileEKLMScintillatorPlane1, TProfile *profileEKLMScintillatorPlane2, bool fill2dHistograms)
 Fill profiles of time versus distance. More...
 
void timeDistance2dFit (const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channels, double &delay, double &delayError)
 Two-dimensional fit for individual channels. More...
 
std::string getExpRunString (Calibration::ExpRun &expRun) const
 Gets the "exp.run" string repr. of (exp,run)
 
std::string getFullObjectPath (const std::string &name, Calibration::ExpRun expRun) const
 constructs the full TDirectory + Key name of an object in a TFile based on its name and exprun
 

Private Attributes

std::map< KLMChannelNumber, std::vector< struct Event > > m_evts
 Container of hit information. More...
 
std::map< KLMChannelNumber, int > m_cFlag
 Calibration flag if the channel has enough hits collected and fitted OK.
 
std::map< KLMChannelNumber, double > m_timeShift
 Shift values of each channel.
 
std::map< KLMChannelNumber, double > m_timeRes
 Resolution values of each channel.
 
std::map< KLMChannelNumber, double > m_time_channel
 Time distribution central value of each channel.
 
std::map< KLMChannelNumber, double > m_etime_channel
 Time distribution central value Error of each channel.
 
std::map< KLMChannelNumber, double > m_ctime_channel
 Calibrated time distribution central value of each channel.
 
std::map< KLMChannelNumber, double > mc_etime_channel
 Calibrated time distribution central value Error of each channel.
 
double m_LowerTimeBoundaryRPC = -10.0
 Lower time boundary for RPC.
 
double m_UpperTimeBoundaryRPC = 10.0
 Upper time boundary for RPC.
 
double m_LowerTimeBoundaryScintilltorsBKLM = 20.0
 Lower time boundary for BKLM scintillators.
 
double m_UpperTimeBoundaryScintilltorsBKLM = 70.0
 Upper time boundary for BKLM scintillators.
 
double m_LowerTimeBoundaryScintilltorsEKLM = 20.0
 Lower time boundary for EKLM scintillators.
 
double m_UpperTimeBoundaryScintilltorsEKLM = 70.0
 Upper time boundary for BKLM scintillators.
 
double m_LowerTimeBoundaryCalibratedRPC = -40.0
 Lower time boundary for RPC (calibrated data).
 
double m_UpperTimeBoundaryCalibratedRPC = 40.0
 Upper time boundary for RPC (calibrated data).
 
double m_LowerTimeBoundaryCalibratedScintilltorsBKLM = -40.0
 Lower time boundary for BKLM scintillators (calibrated data).
 
double m_UpperTimeBoundaryCalibratedScintilltorsBKLM = 40.0
 Upper time boundary for BKLM scintillators (calibrated data).
 
double m_LowerTimeBoundaryCalibratedScintilltorsEKLM = -40.0
 Lower time boundary for EKLM scintillators (calibrated data).
 
double m_UpperTimeBoundaryCalibratedScintilltorsEKLM = 40.0
 Upper time boundary for BKLM scintillators (calibrated data).
 
double m_time_channelAvg_rpc = 0.0
 Central value of the global time distribution (BKLM RPC part).
 
double m_etime_channelAvg_rpc = 0.0
 Central value error of the global time distribution (BKLM RPC part).
 
double m_time_channelAvg_scint = 0.0
 Central value of the global time distribution (BKLM scintillator part).
 
double m_etime_channelAvg_scint = 0.0
 Central value error of the global time distribution (BKLM scintillator part).
 
double m_time_channelAvg_scint_end = 0.0
 Central value of the global time distribution (EKLM scintillator part).
 
double m_etime_channelAvg_scint_end = 0.0
 Central value error of the global time distribution (EKLM scintillator part).
 
double m_ctime_channelAvg_rpc = 0.0
 Calibrated central value of the global time distribution (BKLM RPC part).
 
double mc_etime_channelAvg_rpc = 0.0
 Calibrated central value error of the global time distribution (BKLM RPC part).
 
double m_ctime_channelAvg_scint = 0.0
 Calibrated central value of the global time distribution (BKLM scintillator part).
 
double mc_etime_channelAvg_scint = 0.0
 Calibrated central value error of the global time distribution (BKLM scintillator part).
 
double m_ctime_channelAvg_scint_end = 0.0
 Calibrated central value of the global time distribution (EKLM scintillator part).
 
double mc_etime_channelAvg_scint_end = 0.0
 Calibrated central value error of the global time distribution (EKLM scintillator part).
 
int m_MinimalDigitNumber = 100000000
 Minimal digit number (total).
 
int m_lower_limit_counts = 50
 Lower limit of hits collected for on single channel.
 
const KLMElementNumbersm_ElementNumbers
 Element numbers.
 
const bklm::GeometryParm_BKLMGeometry = nullptr
 BKLM geometry data.
 
const EKLM::GeometryDatam_EKLMGeometry = nullptr
 EKLM geometry data.
 
KLMChannelIndex m_klmChannels
 KLM ChannelIndex object.
 
ROOT::Math::MinimizerOptions m_minimizerOptions
 Minimization options.
 
KLMTimeConstantsm_timeConstants = nullptr
 DBObject of time cost on some parts of the detector.
 
KLMTimeCableDelaym_timeCableDelay = nullptr
 DBObject of the calibration constant of each channel due to cable decay.
 
KLMTimeResolutionm_timeResolution = nullptr
 DBObject of time resolution.
 
bool m_debug = false
 Debug mode.
 
bool m_mc = false
 MC or data.
 
bool m_useEventT0 = true
 Whether to use event T0 from CDC.
 
TH1I * h_calibrated = nullptr
 Calibration statistics for each channel.
 
TH1I * hc_calibrated = nullptr
 Calibration statistics for each channel.
 
TH1F * h_diff = nullptr
 Distance between global and local position.
 
TGraphErrors * gre_time_channel_rpc = nullptr
 BKLM RPC.
 
TGraphErrors * gre_time_channel_scint = nullptr
 BKLM Scintillator.
 
TGraphErrors * gre_time_channel_scint_end = nullptr
 EKLM.
 
TGraphErrors * gre_ctime_channel_rpc = nullptr
 BKLM RPC.
 
TGraphErrors * gre_ctime_channel_scint = nullptr
 BKLM Scintillator.
 
TGraphErrors * gre_ctime_channel_scint_end = nullptr
 EKLM.
 
TGraph * gr_timeShift_channel_rpc = nullptr
 BKLM RPC.
 
TGraph * gr_timeShift_channel_scint = nullptr
 BKLM scintillator.
 
TGraph * gr_timeShift_channel_scint_end = nullptr
 EKLM.
 
TGraph * gr_timeRes_channel_rpc = nullptr
 BKLM RPC.
 
TGraph * gr_timeRes_channel_scint = nullptr
 BKLM scintillator.
 
TGraph * gr_timeRes_channel_scint_end = nullptr
 EKLM.
 
TProfile * m_ProfileRpcPhi = nullptr
 For BKLM RPC phi plane.
 
TProfile * m_ProfileRpcZ = nullptr
 For BKLM RPC z plane.
 
TProfile * m_ProfileBKLMScintillatorPhi = nullptr
 For BKLM scintillator phi plane.
 
TProfile * m_ProfileBKLMScintillatorZ = nullptr
 For BKLM scintillator z plane.
 
TProfile * m_ProfileEKLMScintillatorPlane1 = nullptr
 For EKLM scintillator plane1.
 
TProfile * m_ProfileEKLMScintillatorPlane2 = nullptr
 For EKLM scintillator plane2.
 
TProfile * m_Profile2RpcPhi = nullptr
 For BKLM RPC phi plane.
 
TProfile * m_Profile2RpcZ = nullptr
 For BKLM RPC z plane.
 
TProfile * m_Profile2BKLMScintillatorPhi = nullptr
 For BKLM scintillator phi plane.
 
TProfile * m_Profile2BKLMScintillatorZ = nullptr
 For BKLM scintillator z plane.
 
TProfile * m_Profile2EKLMScintillatorPlane1 = nullptr
 For EKLM scintillator plane1.
 
TProfile * m_Profile2EKLMScintillatorPlane2 = nullptr
 For EKLM scintillator plane2.
 
TH1F * h_time_rpc_tc = nullptr
 BKLM RPC part.
 
TH1F * h_time_scint_tc = nullptr
 BKLM scintillator part.
 
TH1F * h_time_scint_tc_end = nullptr
 EKLM part.
 
TH1F * h_time_rpc = nullptr
 BKLM RPC part.
 
TH1F * h_time_scint = nullptr
 BKLM scintillator part.
 
TH1F * h_time_scint_end = nullptr
 EKLM part.
 
TH1F * hc_time_rpc = nullptr
 BKLM RPC part.
 
TH1F * hc_time_scint = nullptr
 BKLM scintillator part.
 
TH1F * hc_time_scint_end = nullptr
 EKLM part.
 
TH1F * h_timeF_rpc [2] = {nullptr}
 BKLM RPC part.
 
TH1F * h_timeF_scint [2] = {nullptr}
 BKLM scintillator part.
 
TH1F * h_timeF_scint_end [2] = {nullptr}
 EKLM part.
 
TH1F * hc_timeF_rpc [2] = {nullptr}
 BKLM RPC part.
 
TH1F * hc_timeF_scint [2] = {nullptr}
 BKLM scintillator part.
 
TH1F * hc_timeF_scint_end [2] = {nullptr}
 EKLM part.
 
TH2F * h2_timeF_rpc [2] = {nullptr}
 BKLM RPC part.
 
TH2F * h2_timeF_scint [2] = {nullptr}
 BKLM scintillator part.
 
TH2F * h2_timeF_scint_end [2] = {nullptr}
 EKLM part.
 
TH2F * h2c_timeF_rpc [2] = {nullptr}
 BKLM RPC part.
 
TH2F * h2c_timeF_scint [2] = {nullptr}
 BKLM scintillator part.
 
TH2F * h2c_timeF_scint_end [2] = {nullptr}
 EKLM part.
 
TH1F * h_timeFS_rpc [2][8] = {{nullptr}}
 BKLM RPC part.
 
TH1F * h_timeFS_scint [2][8] = {{nullptr}}
 BKLM scintillator part.
 
TH1F * h_timeFS_scint_end [2][4] = {{nullptr}}
 EKLM part.
 
TH1F * hc_timeFS_rpc [2][8] = {{nullptr}}
 BKLM RPC part.
 
TH1F * hc_timeFS_scint [2][8] = {{nullptr}}
 BKLM scintillator part.
 
TH1F * hc_timeFS_scint_end [2][4] = {{nullptr}}
 EKLM part.
 
TH2F * h2_timeFS [2][8] = {{nullptr}}
 BKLM part.
 
TH2F * h2_timeFS_end [2][4] = {{nullptr}}
 EKLM part.
 
TH2F * h2c_timeFS [2][8] = {{nullptr}}
 BKLM part.
 
TH2F * h2c_timeFS_end [2][4] = {{nullptr}}
 EKLM part.
 
TH1F * h_timeFSL [2][8][15] = {{{nullptr}}}
 BKLM part.
 
TH1F * h_timeFSL_end [2][4][14] = {{{nullptr}}}
 EKLM part.
 
TH1F * hc_timeFSL [2][8][15] = {{{nullptr}}}
 BKLM part.
 
TH1F * hc_timeFSL_end [2][4][14] = {{{nullptr}}}
 EKLM part.
 
TH1F * h_timeFSLP [2][8][15][2] = {{{{nullptr}}}}
 BKLM part.
 
TH1F * h_timeFSLP_end [2][4][14][2] = {{{{nullptr}}}}
 EKLM part.
 
TH1F * hc_timeFSLP [2][8][15][2] = {{{{nullptr}}}}
 BKLM part.
 
TH1F * hc_timeFSLP_end [2][4][14][2] = {{{{nullptr}}}}
 EKLM part.
 
TH2F * h2_timeFSLP [2][8][15][2] = {{{{nullptr}}}}
 BKLM part.
 
TH2F * h2_timeFSLP_end [2][4][14][2] = {{{{nullptr}}}}
 EKLM part.
 
TH2F * h2c_timeFSLP [2][8][15][2] = {{{{nullptr}}}}
 BKLM part.
 
TH2F * h2c_timeFSLP_end [2][4][14][2] = {{{{nullptr}}}}
 EKLM part.
 
TH1F * h_timeFSLPC_tc [2][8][15][2][54] = {{{{{nullptr}}}}}
 BKLM part, used for effective light speed estimation.
 
TH1F * h_timeFSLPC [2][8][15][2][54] = {{{{{nullptr}}}}}
 BKLM part.
 
TH2F * m_HistTimeLengthBKLM [2][8][15][2][54] = {{{{{nullptr}}}}}
 Two-dimensional distributions of time versus propagation length.
 
TH1F * h_timeFSLPC_tc_end [2][4][14][2][75] = {{{{{nullptr}}}}}
 EKLM part, used for effective light speed estimation.
 
TH1F * h_timeFSLPC_end [2][4][14][2][75] = {{{{{nullptr}}}}}
 EKLM part.
 
TH2F * m_HistTimeLengthEKLM [2][4][14][2][75] = {{{{{nullptr}}}}}
 Two-dimensional distributions of time versus propagation length.
 
TH1F * hc_timeFSLPC [2][8][15][2][54] = {{{{{nullptr}}}}}
 BKLM part.
 
TH1F * hc_timeFSLPC_end [2][4][14][2][75] = {{{{{nullptr}}}}}
 EKLM part.
 
TF1 * fcn_pol1 = nullptr
 Pol1 function. More...
 
TF1 * fcn_const = nullptr
 Const function. More...
 
TF1 * fcn_gaus = nullptr
 Gaussian function. More...
 
TF1 * fcn_land = nullptr
 Landau function. More...
 
TFile * m_outFile = nullptr
 Output file.
 
std::vector< std::string > m_inputFileNames
 List of input files to the Algorithm, will initially be user defined but then gets the wildcards expanded during execute()
 
std::map< Calibration::ExpRun, std::vector< std::string > > m_runsToInputFiles
 Map of Runs to input files. Gets filled when you call getRunRangeFromAllData, gets cleared when setting input files again.
 
std::string m_granularityOfData
 Granularity of input data. This only changes when the input files change so it isn't specific to an execution.
 
ExecutionData m_data
 Data specific to a SINGLE execution of the algorithm. Gets reset at the beginning of execution.
 
std::string m_description {""}
 Description of the algorithm.
 
std::string m_prefix {""}
 The name of the TDirectory the collector objects are contained within.
 
nlohmann::json m_jsonExecutionInput = nlohmann::json::object()
 Optional input JSON object used to make decisions about how to execute the algorithm code.
 
nlohmann::json m_jsonExecutionOutput = nlohmann::json::object()
 Optional output JSON object that can be set during the execution by the underlying algorithm code.
 

Static Private Attributes

static const Calibration::ExpRun m_allExpRun = make_pair(-1, -1)
 allExpRun
 

Detailed Description

KLM time calibration algorithm.

Definition at line 41 of file KLMTimeAlgorithm.h.

Member Enumeration Documentation

◆ EResult

enum EResult
inherited

The result of calibration.

Enumerator
c_OK 

Finished successfuly =0 in Python.

c_Iterate 

Needs iteration =1 in Python.

c_NotEnoughData 

Needs more data =2 in Python.

c_Failure 

Failed =3 in Python.

c_Undefined 

Not yet known (before execution) =4 in Python.

Definition at line 40 of file CalibrationAlgorithm.h.

Member Function Documentation

◆ calibrate()

CalibrationAlgorithm::EResult calibrate ( )
overrideprotectedvirtual

Run algorithm on data.

=======================================================================================

Implements CalibrationAlgorithm.

Definition at line 740 of file KLMTimeAlgorithm.cc.

741 {
742  int channelId;
743  gROOT->SetBatch(kTRUE);
744  setupDatabase();
748 
749  fcn_gaus = new TF1("fcn_gaus", "gaus");
750  fcn_land = new TF1("fcn_land", "landau");
751  fcn_pol1 = new TF1("fcn_pol1", "pol1");
752  fcn_const = new TF1("fcn_const", "pol0");
753 
755  if (result != CalibrationAlgorithm::c_OK)
756  return result;
757 
758  /* Choose non-existing file name. */
759  std::string name = "time_calibration.root";
760  int i = 1;
761  while (1) {
762  struct stat buffer;
763  if (stat(name.c_str(), &buffer) != 0)
764  break;
765  name = "time_calibration_" + std::to_string(i) + ".root";
766  i = i + 1;
767  /* Overflow. */
768  if (i < 0)
769  break;
770  }
771  m_outFile = new TFile(name.c_str(), "recreate");
773 
774  std::vector<struct Event>::iterator it;
775  std::vector<struct Event> eventsChannel;
776 
777  eventsChannel.clear();
778  m_cFlag.clear();
779  m_minimizerOptions.SetDefaultStrategy(2);
780 
781  /* Sort channels by number of events. */
782  std::vector< std::pair<KLMChannelNumber, unsigned int> > channelsBKLM;
783  std::vector< std::pair<KLMChannelNumber, unsigned int> > channelsEKLM;
784  KLMChannelIndex klmChannels;
785  for (KLMChannelIndex& klmChannel : klmChannels) {
786  KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
787  m_cFlag[channel] = ChannelCalibrationStatus::c_NotEnoughData;
788  if (m_evts.find(channel) == m_evts.end())
789  continue;
790  int nEvents = m_evts[channel].size();
791  if (nEvents < m_lower_limit_counts) {
792  B2WARNING("Not enough calibration data collected."
793  << LogVar("channel", channel)
794  << LogVar("number of digit", nEvents));
795  continue;
796  }
797  m_cFlag[channel] = ChannelCalibrationStatus::c_FailedFit;
798  if (klmChannel.getSubdetector() == KLMElementNumbers::c_BKLM &&
799  klmChannel.getLayer() < BKLMElementNumbers::c_FirstRPCLayer) {
800  channelsBKLM.push_back(
801  std::pair<KLMChannelNumber, unsigned int>(channel, nEvents));
802  }
803  if (klmChannel.getSubdetector() == KLMElementNumbers::c_EKLM) {
804  channelsEKLM.push_back(
805  std::pair<KLMChannelNumber, unsigned int>(channel, nEvents));
806  }
807  }
808  std::sort(channelsBKLM.begin(), channelsBKLM.end(), compareEventNumber);
809  std::sort(channelsEKLM.begin(), channelsEKLM.end(), compareEventNumber);
810 
811  /* Two-dimensional fit for the channel with the maximal number of events. */
812  double delayBKLM, delayBKLMError;
813  double delayEKLM, delayEKLMError;
814  timeDistance2dFit(channelsBKLM, delayBKLM, delayBKLMError);
815  timeDistance2dFit(channelsEKLM, delayEKLM, delayEKLMError);
816 
817  /**********************************
818  * First loop
819  * Estimation of effective light speed for Scintillators and RPCs, separately.
820  **********************************/
821  B2INFO("Effective light speed Estimation.");
822  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
823  channelId = klmChannel.getKLMChannelNumber();
824  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
825  continue;
826 
827  eventsChannel = m_evts[channelId];
828  int iSub = klmChannel.getSubdetector();
829 
830  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
831  TVector3 diffD = TVector3(it->diffDistX, it->diffDistY, it->diffDistZ);
832  h_diff->Fill(diffD.Mag());
833  double timeHit = it->time();
834  if (m_useEventT0)
835  timeHit = timeHit - it->t0;
836  if (iSub == KLMElementNumbers::c_BKLM) {
837  int iF = klmChannel.getSection();
838  int iS = klmChannel.getSector() - 1;
839  int iL = klmChannel.getLayer() - 1;
840  int iP = klmChannel.getPlane();
841  int iC = klmChannel.getStrip() - 1;
842  h_timeFSLPC_tc[iF][iS][iL][iP][iC]->Fill(timeHit);
843  if (iL > 1) {
844  h_time_rpc_tc->Fill(timeHit);
845  } else {
846  h_time_scint_tc->Fill(timeHit);
847  }
848  } else {
849  int iF = klmChannel.getSection() - 1;
850  int iS = klmChannel.getSector() - 1;
851  int iL = klmChannel.getLayer() - 1;
852  int iP = klmChannel.getPlane() - 1;
853  int iC = klmChannel.getStrip() - 1;
854  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC]->Fill(timeHit);
855  h_time_scint_tc_end->Fill(timeHit);
856  }
857  }
858  }
859  B2INFO("Effective light speed Estimation! Hists and Graph filling done.");
860 
861  m_timeShift.clear();
862 
863  double tmpMean_rpc_global = h_time_rpc_tc->GetMean();
864  double tmpMean_scint_global = h_time_scint_tc->GetMean();
865  double tmpMean_scint_global_end = h_time_scint_tc_end->GetMean();
866 
867  B2INFO("Global Mean for Raw." << LogVar("RPC", tmpMean_rpc_global) << LogVar("Scint BKLM",
868  tmpMean_scint_global) << LogVar("Scint EKLM", tmpMean_scint_global_end));
869 
870  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
871  channelId = klmChannel.getKLMChannelNumber();
872  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
873  continue;
874 
875  int iSub = klmChannel.getSubdetector();
876  if (iSub == KLMElementNumbers::c_BKLM) {
877  int iF = klmChannel.getSection();
878  int iS = klmChannel.getSector() - 1;
879  int iL = klmChannel.getLayer() - 1;
880  int iP = klmChannel.getPlane();
881  int iC = klmChannel.getStrip() - 1;
882  h_timeFSLPC_tc[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
883  double tmpMean_channel = fcn_gaus->GetParameter(1);
884  if (iL > 1) {
885  m_timeShift[channelId] = tmpMean_channel - tmpMean_rpc_global;
886  } else {
887  m_timeShift[channelId] = tmpMean_channel - tmpMean_scint_global;
888  }
889  } else {
890  int iF = klmChannel.getSection() - 1;
891  int iS = klmChannel.getSector() - 1;
892  int iL = klmChannel.getLayer() - 1;
893  int iP = klmChannel.getPlane() - 1;
894  int iC = klmChannel.getStrip() - 1;
895  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
896  double tmpMean_channel = fcn_gaus->GetParameter(1);
897  m_timeShift[channelId] = tmpMean_channel - tmpMean_scint_global_end;
898  }
899  }
900 
901  delete h_time_scint_tc;
902  delete h_time_scint_tc_end;
903  delete h_time_rpc_tc;
904  B2INFO("Effective Light m_timeShift obtained. done.");
905 
910 
911  B2INFO("Effective light speed fitting.");
912  m_ProfileRpcPhi->Fit("fcn_pol1", "EMQ");
913  double delayRPCPhi = fcn_pol1->GetParameter(1);
914  double e_slope_rpc_phi = fcn_pol1->GetParError(1);
915 
916  m_ProfileRpcZ->Fit("fcn_pol1", "EMQ");
917  double delayRPCZ = fcn_pol1->GetParameter(1);
918  double e_slope_rpc_z = fcn_pol1->GetParError(1);
919 
920  m_ProfileBKLMScintillatorPhi->Fit("fcn_pol1", "EMQ");
921  double slope_scint_phi = fcn_pol1->GetParameter(1);
922  double e_slope_scint_phi = fcn_pol1->GetParError(1);
923 
924  m_ProfileBKLMScintillatorZ->Fit("fcn_pol1", "EMQ");
925  double slope_scint_z = fcn_pol1->GetParameter(1);
926  double e_slope_scint_z = fcn_pol1->GetParError(1);
927 
928  m_ProfileEKLMScintillatorPlane1->Fit("fcn_pol1", "EMQ");
929  double slope_scint_plane1_end = fcn_pol1->GetParameter(1);
930  double e_slope_scint_plane1_end = fcn_pol1->GetParError(1);
931 
932  m_ProfileEKLMScintillatorPlane2->Fit("fcn_pol1", "EMQ");
933  double slope_scint_plane2_end = fcn_pol1->GetParameter(1);
934  double e_slope_scint_plane2_end = fcn_pol1->GetParError(1);
935 
936  TString logStr_phi, logStr_z;
937  logStr_phi = Form("%.4f +/- %.4f ns/cm", delayRPCPhi, e_slope_rpc_phi);
938  logStr_z = Form("%.4f +/- %.4f ns/cm", delayRPCZ, e_slope_rpc_z);
939  B2INFO("Delay in RPCs:"
940  << LogVar("Fitted Value (phi readout) ", logStr_phi.Data())
941  << LogVar("Fitted Value (z readout) ", logStr_z.Data()));
942  logStr_phi = Form("%.4f +/- %.4f ns/cm", slope_scint_phi, e_slope_scint_phi);
943  logStr_z = Form("%.4f +/- %.4f ns/cm", slope_scint_z, e_slope_scint_z);
944  B2INFO("Delay in BKLM scintillators:"
945  << LogVar("Fitted Value (phi readout) ", logStr_phi.Data())
946  << LogVar("Fitted Value (z readout) ", logStr_z.Data()));
947  logStr_phi = Form("%.4f +/- %.4f ns/cm", slope_scint_plane1_end,
948  e_slope_scint_plane1_end);
949  logStr_z = Form("%.4f +/- %.4f ns/cm", slope_scint_plane2_end,
950  e_slope_scint_plane2_end);
951  B2INFO("Delay in EKLM scintillators:"
952  << LogVar("Fitted Value (plane1 readout) ", logStr_phi.Data())
953  << LogVar("Fitted Value (plane2 readout) ", logStr_z.Data()));
954 
955  logStr_z = Form("%.4f +/- %.4f ns/cm", delayBKLM, delayBKLMError);
956  B2INFO("Delay in BKLM scintillators:"
957  << LogVar("Fitted Value (2d fit) ", logStr_z.Data()));
958  logStr_z = Form("%.4f +/- %.4f ns/cm", delayEKLM, delayEKLMError);
959  B2INFO("Delay in EKLM scintillators:"
960  << LogVar("Fitted Value (2d fit) ", logStr_z.Data()));
961 
962  // Default Effective light speed in current Database
963  //delayEKLM = 0.5 * (slope_scint_plane1_end + slope_scint_plane2_end);
964  //delayBKLM = 0.5 * (slope_scint_phi + slope_scint_z);
965 
970 
972  B2INFO("Time distribution filling begins.");
973  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
974  channelId = klmChannel.getKLMChannelNumber();
975  int iSub = klmChannel.getSubdetector();
976 
977  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
978  continue;
979  eventsChannel = m_evts[channelId];
980 
981  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
982  double timeHit = it->time();
983  if (m_useEventT0)
984  timeHit = timeHit - it->t0;
985  if (iSub == KLMElementNumbers::c_BKLM) {
986  int iF = klmChannel.getSection();
987  int iS = klmChannel.getSector() - 1;
988  int iL = klmChannel.getLayer() - 1;
989  int iP = klmChannel.getPlane();
990  int iC = klmChannel.getStrip() - 1;
991  if (iL > 1) {
992  double propgationT;
994  propgationT = it->dist * delayRPCZ;
995  else
996  propgationT = it->dist * delayRPCPhi;
997  double time = timeHit - propgationT;
998  h_time_rpc->Fill(time);
999  h_timeF_rpc[iF]->Fill(time);
1000  h_timeFS_rpc[iF][iS]->Fill(time);
1001  h_timeFSL[iF][iS][iL]->Fill(time);
1002  h_timeFSLP[iF][iS][iL][iP]->Fill(time);
1003  h_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1004  h2_timeF_rpc[iF]->Fill(iS, time);
1005  h2_timeFS[iF][iS]->Fill(iL, time);
1006  h2_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1007  } else {
1008  double propgationT = it->dist * delayBKLM;
1009  double time = timeHit - propgationT;
1010  h_time_scint->Fill(time);
1011  h_timeF_scint[iF]->Fill(time);
1012  h_timeFS_scint[iF][iS]->Fill(time);
1013  h_timeFSL[iF][iS][iL]->Fill(time);
1014  h_timeFSLP[iF][iS][iL][iP]->Fill(time);
1015  h_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1016  h2_timeF_scint[iF]->Fill(iS, time);
1017  h2_timeFS[iF][iS]->Fill(iL, time);
1018  h2_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1019  }
1020  } else {
1021  int iF = klmChannel.getSection() - 1;
1022  int iS = klmChannel.getSector() - 1;
1023  int iL = klmChannel.getLayer() - 1;
1024  int iP = klmChannel.getPlane() - 1;
1025  int iC = klmChannel.getStrip() - 1;
1026  double propgationT = it->dist * delayEKLM;
1027  double time = timeHit - propgationT;
1028  h_time_scint_end->Fill(time);
1029  h_timeF_scint_end[iF]->Fill(time);
1030  h_timeFS_scint_end[iF][iS]->Fill(time);
1031  h_timeFSL_end[iF][iS][iL]->Fill(time);
1032  h_timeFSLP_end[iF][iS][iL][iP]->Fill(time);
1033  h_timeFSLPC_end[iF][iS][iL][iP][iC]->Fill(time);
1034  h2_timeF_scint_end[iF]->Fill(iS, time);
1035  h2_timeFS_end[iF][iS]->Fill(iL, time);
1036  h2_timeFSLP_end[iF][iS][iL][iP]->Fill(iC, time);
1037  }
1038  }
1039  }
1040 
1041  B2INFO("Orignal filling done.");
1042 
1043  int iChannel_rpc = 0;
1044  int iChannel = 0;
1045  int iChannel_end = 0;
1046  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1047  channelId = klmChannel.getKLMChannelNumber();
1048  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
1049  continue;
1050  int iSub = klmChannel.getSubdetector();
1051 
1052  if (iSub == KLMElementNumbers::c_BKLM) {
1053  int iF = klmChannel.getSection();
1054  int iS = klmChannel.getSector() - 1;
1055  int iL = klmChannel.getLayer() - 1;
1056  int iP = klmChannel.getPlane();
1057  int iC = klmChannel.getStrip() - 1;
1058 
1059  TFitResultPtr r = h_timeFSLPC[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1060  if (int(r) != 0)
1061  continue;
1062  if (int(r) == 0)
1063  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1064  m_time_channel[channelId] = fcn_gaus->GetParameter(1);
1065  m_etime_channel[channelId] = fcn_gaus->GetParError(1);
1066  if (iL > 1) {
1067  gre_time_channel_rpc->SetPoint(iChannel_rpc, channelId, m_time_channel[channelId]);
1068  gre_time_channel_rpc->SetPointError(iChannel_rpc, 0., m_etime_channel[channelId]);
1069  iChannel_rpc++;
1070  } else {
1071  gre_time_channel_scint->SetPoint(iChannel, channelId, m_time_channel[channelId]);
1072  gre_time_channel_scint->SetPointError(iChannel, 0., m_etime_channel[channelId]);
1073  iChannel++;
1074  }
1075  } else {
1076  int iF = klmChannel.getSection() - 1;
1077  int iS = klmChannel.getSector() - 1;
1078  int iL = klmChannel.getLayer() - 1;
1079  int iP = klmChannel.getPlane() - 1;
1080  int iC = klmChannel.getStrip() - 1;
1081 
1082  TFitResultPtr r = h_timeFSLPC_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1083  if (int(r) != 0)
1084  continue;
1085  if (int(r) == 0)
1086  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1087  m_time_channel[channelId] = fcn_gaus->GetParameter(1);
1088  m_etime_channel[channelId] = fcn_gaus->GetParError(1);
1089  gre_time_channel_scint_end->SetPoint(iChannel_end, channelId, m_time_channel[channelId]);
1090  gre_time_channel_scint_end->SetPointError(iChannel_end, 0., m_etime_channel[channelId]);
1091  iChannel_end++;
1092  }
1093  }
1094 
1095  gre_time_channel_scint->Fit("fcn_const", "EMQ");
1096  m_time_channelAvg_scint = fcn_const->GetParameter(0);
1097  m_etime_channelAvg_scint = fcn_const->GetParError(0);
1098 
1099  gre_time_channel_scint_end->Fit("fcn_const", "EMQ");
1100  m_time_channelAvg_scint_end = fcn_const->GetParameter(0);
1101  m_etime_channelAvg_scint_end = fcn_const->GetParError(0);
1102 
1103  gre_time_channel_rpc->Fit("fcn_const", "EMQ");
1104  m_time_channelAvg_rpc = fcn_const->GetParameter(0);
1105  m_etime_channelAvg_rpc = fcn_const->GetParError(0);
1106 
1107  B2INFO("Channel's time distribution fitting done.");
1108  B2DEBUG(20, LogVar("Average time (RPC)", m_time_channelAvg_rpc)
1109  << LogVar("Average time (BKLM scintillators)", m_time_channelAvg_scint)
1110  << LogVar("Average time (EKLM scintillators)", m_time_channelAvg_scint_end));
1111 
1112  B2INFO("Calibrated channel's time distribution filling begins.");
1113 
1114  m_timeShift.clear();
1115  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1116  channelId = klmChannel.getKLMChannelNumber();
1117  h_calibrated->Fill(m_cFlag[channelId]);
1118  if (m_time_channel.find(channelId) == m_time_channel.end())
1119  continue;
1120  double timeShift = m_time_channel[channelId];
1121  int iSub = klmChannel.getSubdetector();
1122  if (iSub == KLMElementNumbers::c_BKLM) {
1123  int iL = klmChannel.getLayer() - 1;
1124  if (iL > 1)
1125  m_timeShift[channelId] = timeShift;
1126  else
1127  m_timeShift[channelId] = timeShift;
1128  } else {
1129  m_timeShift[channelId] = timeShift;
1130  }
1131  m_timeCableDelay->setTimeDelay(channelId, m_timeShift[channelId]);
1132  }
1133 
1134  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1135  channelId = klmChannel.getKLMChannelNumber();
1136  if (m_timeShift.find(channelId) != m_timeShift.end())
1137  continue;
1138  m_timeShift[channelId] = esti_timeShift(klmChannel);
1139  m_timeCableDelay->setTimeDelay(channelId, m_timeShift[channelId]);
1140  B2DEBUG(20, "Uncalibrated Estimation " << LogVar("Channel", channelId) << LogVar("Estimated value", m_timeShift[channelId]));
1141  }
1142 
1143  iChannel_rpc = 0;
1144  iChannel = 0;
1145  iChannel_end = 0;
1146  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1147  channelId = klmChannel.getKLMChannelNumber();
1148  if (m_timeShift.find(channelId) == m_timeShift.end()) {
1149  B2ERROR("!!! Not All Channels Calibration Constant Set. Error Happended on " << LogVar("Channel", channelId));
1150  continue;
1151  }
1152  int iSub = klmChannel.getSubdetector();
1153  if (iSub == KLMElementNumbers::c_BKLM) {
1154  int iL = klmChannel.getLayer();
1155  if (iL > 2) {
1156  gr_timeShift_channel_rpc->SetPoint(iChannel_rpc, channelId, m_timeShift[channelId]);
1157  iChannel_rpc++;
1158  } else {
1159  gr_timeShift_channel_scint->SetPoint(iChannel, channelId, m_timeShift[channelId]);
1160  iChannel++;
1161  }
1162  } else {
1163  gr_timeShift_channel_scint_end->SetPoint(iChannel_end, channelId, m_timeShift[channelId]);
1164  iChannel_end++;
1165  }
1166  }
1167 
1172 
1173  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1174  channelId = klmChannel.getKLMChannelNumber();
1175  int iSub = klmChannel.getSubdetector();
1176  eventsChannel = m_evts[channelId];
1177  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
1178  double timeHit = it->time();
1179  if (m_useEventT0)
1180  timeHit = timeHit - it->t0;
1181  if (iSub == KLMElementNumbers::c_BKLM) {
1182  int iF = klmChannel.getSection();
1183  int iS = klmChannel.getSector() - 1;
1184  int iL = klmChannel.getLayer() - 1;
1185  int iP = klmChannel.getPlane();
1186  int iC = klmChannel.getStrip() - 1;
1187  if (iL > 1) {
1188  double propgationT;
1189  if (iP == BKLMElementNumbers::c_ZPlane)
1190  propgationT = it->dist * delayRPCZ;
1191  else
1192  propgationT = it->dist * delayRPCPhi;
1193  double time = timeHit - propgationT - m_timeShift[channelId];
1194  hc_time_rpc->Fill(time);
1195  hc_timeF_rpc[iF]->Fill(time);
1196  hc_timeFS_rpc[iF][iS]->Fill(time);
1197  hc_timeFSL[iF][iS][iL]->Fill(time);
1198  hc_timeFSLP[iF][iS][iL][iP]->Fill(time);
1199  hc_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1200  h2c_timeF_rpc[iF]->Fill(iS, time);
1201  h2c_timeFS[iF][iS]->Fill(iL, time);
1202  h2c_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1203  } else {
1204  double propgationT = it->dist * delayBKLM;
1205  double time = timeHit - propgationT - m_timeShift[channelId];
1206  hc_time_scint->Fill(time);
1207  hc_timeF_scint[iF]->Fill(time);
1208  hc_timeFS_scint[iF][iS]->Fill(time);
1209  hc_timeFSL[iF][iS][iL]->Fill(time);
1210  hc_timeFSLP[iF][iS][iL][iP]->Fill(time);
1211  hc_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1212  h2c_timeF_scint[iF]->Fill(iS, time);
1213  h2c_timeFS[iF][iS]->Fill(iL, time);
1214  h2c_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1215  }
1216  } else {
1217  int iF = klmChannel.getSection() - 1;
1218  int iS = klmChannel.getSector() - 1;
1219  int iL = klmChannel.getLayer() - 1;
1220  int iP = klmChannel.getPlane() - 1;
1221  int iC = klmChannel.getStrip() - 1;
1222  double propgationT = it->dist * delayEKLM;
1223  double time = timeHit - propgationT - m_timeShift[channelId];
1224  hc_time_scint_end->Fill(time);
1225  hc_timeF_scint_end[iF]->Fill(time);
1226  hc_timeFS_scint_end[iF][iS]->Fill(time);
1227  hc_timeFSL_end[iF][iS][iL]->Fill(time);
1228  hc_timeFSLP_end[iF][iS][iL][iP]->Fill(time);
1229  hc_timeFSLPC_end[iF][iS][iL][iP][iC]->Fill(time);
1230  h2c_timeF_scint_end[iF]->Fill(iS, time);
1231  h2c_timeFS_end[iF][iS]->Fill(iL, time);
1232  h2c_timeFSLP_end[iF][iS][iL][iP]->Fill(iC, time);
1233  }
1234  }
1235  }
1236 
1237  int icChannel_rpc = 0;
1238  int icChannel = 0;
1239  int icChannel_end = 0;
1240  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1241  channelId = klmChannel.getKLMChannelNumber();
1242  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
1243  continue;
1244  int iSub = klmChannel.getSubdetector();
1245 
1246  if (iSub == KLMElementNumbers::c_BKLM) {
1247  int iF = klmChannel.getSection();
1248  int iS = klmChannel.getSector() - 1;
1249  int iL = klmChannel.getLayer() - 1;
1250  int iP = klmChannel.getPlane();
1251  int iC = klmChannel.getStrip() - 1;
1252 
1253  TFitResultPtr rc = hc_timeFSLPC[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1254  if (int(rc) != 0)
1255  continue;
1256  if (int(rc) == 0)
1257  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1258  m_ctime_channel[channelId] = fcn_gaus->GetParameter(1);
1259  mc_etime_channel[channelId] = fcn_gaus->GetParError(1);
1260  if (iL > 1) {
1261  gre_ctime_channel_rpc->SetPoint(icChannel_rpc, channelId, m_ctime_channel[channelId]);
1262  gre_ctime_channel_rpc->SetPointError(icChannel_rpc, 0., mc_etime_channel[channelId]);
1263  icChannel_rpc++;
1264  } else {
1265  gre_ctime_channel_scint->SetPoint(icChannel, channelId, m_ctime_channel[channelId]);
1266  gre_ctime_channel_scint->SetPointError(icChannel, 0., mc_etime_channel[channelId]);
1267  icChannel++;
1268  }
1269  } else {
1270  int iF = klmChannel.getSection() - 1;
1271  int iS = klmChannel.getSector() - 1;
1272  int iL = klmChannel.getLayer() - 1;
1273  int iP = klmChannel.getPlane() - 1;
1274  int iC = klmChannel.getStrip() - 1;
1275 
1276  TFitResultPtr rc = hc_timeFSLPC_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1277  if (int(rc) != 0)
1278  continue;
1279  if (int(rc) == 0)
1280  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1281  m_ctime_channel[channelId] = fcn_gaus->GetParameter(1);
1282  mc_etime_channel[channelId] = fcn_gaus->GetParError(1);
1283  gre_ctime_channel_scint_end->SetPoint(icChannel_end, channelId, m_ctime_channel[channelId]);
1284  gre_ctime_channel_scint_end->SetPointError(icChannel_end, 0., mc_etime_channel[channelId]);
1285  icChannel_end++;
1286  }
1287  }
1288 
1289  gre_ctime_channel_scint->Fit("fcn_const", "EMQ");
1290  m_ctime_channelAvg_scint = fcn_const->GetParameter(0);
1291  mc_etime_channelAvg_scint = fcn_const->GetParError(0);
1292 
1293  gre_ctime_channel_scint_end->Fit("fcn_const", "EMQ");
1294  m_ctime_channelAvg_scint_end = fcn_const->GetParameter(0);
1295  mc_etime_channelAvg_scint_end = fcn_const->GetParError(0);
1296 
1297  gre_ctime_channel_rpc->Fit("fcn_const", "EMQ");
1298  m_ctime_channelAvg_rpc = fcn_const->GetParameter(0);
1299  mc_etime_channelAvg_rpc = fcn_const->GetParError(0);
1300 
1301  B2INFO("Channel's time distribution fitting done.");
1302  B2DEBUG(20, LogVar("Average calibrated time (RPC)", m_ctime_channelAvg_rpc)
1303  << LogVar("Average calibrated time (BKLM scintillators)", m_ctime_channelAvg_scint)
1304  << LogVar("Average calibrated time (EKLM scintillators)", m_ctime_channelAvg_scint_end));
1305 
1306  B2INFO("Calibrated channel's time distribution filling begins.");
1307 
1308  m_timeRes.clear();
1309  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1310  channelId = klmChannel.getKLMChannelNumber();
1311  hc_calibrated->Fill(m_cFlag[channelId]);
1312  if (m_ctime_channel.find(channelId) == m_ctime_channel.end())
1313  continue;
1314  double timeRes = m_ctime_channel[channelId];
1315  int iSub = klmChannel.getSubdetector();
1316  if (iSub == KLMElementNumbers::c_BKLM) {
1317  int iL = klmChannel.getLayer() - 1;
1318  if (iL > 1)
1319  m_timeRes[channelId] = timeRes;
1320  else
1321  m_timeRes[channelId] = timeRes;
1322  } else {
1323  m_timeRes[channelId] = timeRes;
1324  }
1325  m_timeResolution->setTimeResolution(channelId, m_timeRes[channelId]);
1326  }
1327 
1328  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1329  channelId = klmChannel.getKLMChannelNumber();
1330  if (m_timeRes.find(channelId) != m_timeRes.end())
1331  continue;
1332  m_timeRes[channelId] = esti_timeRes(klmChannel);
1333  m_timeResolution->setTimeResolution(channelId, m_timeRes[channelId]);
1334  B2DEBUG(20, "Calibrated Estimation " << LogVar("Channel", channelId) << LogVar("Estimated value", m_timeRes[channelId]));
1335  }
1336 
1337  icChannel_rpc = 0;
1338  icChannel = 0;
1339  icChannel_end = 0;
1340  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1341  channelId = klmChannel.getKLMChannelNumber();
1342  if (m_timeRes.find(channelId) == m_timeRes.end()) {
1343  B2ERROR("!!! Not All Channels Calibration Constant Set. Error Happended on " << LogVar("Channel", channelId));
1344  continue;
1345  }
1346  int iSub = klmChannel.getSubdetector();
1347  if (iSub == KLMElementNumbers::c_BKLM) {
1348  int iL = klmChannel.getLayer();
1349  if (iL > 2) {
1350  gr_timeRes_channel_rpc->SetPoint(icChannel_rpc, channelId, m_timeRes[channelId]);
1351  icChannel_rpc++;
1352  } else {
1353  gr_timeRes_channel_scint->SetPoint(icChannel, channelId, m_timeRes[channelId]);
1354  icChannel++;
1355  }
1356  } else {
1357  gr_timeRes_channel_scint_end->SetPoint(icChannel_end, channelId, m_timeRes[channelId]);
1358  icChannel_end++;
1359  }
1360  }
1361  saveHist();
1362 
1363  delete fcn_const;
1364  m_evts.clear();
1365  m_timeShift.clear();
1366  m_timeRes.clear();
1367  m_cFlag.clear();
1368 
1369  saveCalibration(m_timeCableDelay, "KLMTimeCableDelay");
1370  saveCalibration(m_timeConstants, "KLMTimeConstants");
1371  saveCalibration(m_timeResolution, "KLMTimeResolution");
1373 }
@ c_FirstRPCLayer
First RPC layer.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
KLM channel index.
KLMChannelIndex begin()
First channel.
KLMChannelIndex & end()
Last channel.
TProfile * m_Profile2EKLMScintillatorPlane2
For EKLM scintillator plane2.
double mc_etime_channelAvg_rpc
Calibrated central value error of the global time distribution (BKLM RPC part).
TH1F * h_timeFSLPC_tc[2][8][15][2][54]
BKLM part, used for effective light speed estimation.
TH2F * h2c_timeF_scint_end[2]
EKLM part.
TH1F * h_timeFSLPC_tc_end[2][4][14][2][75]
EKLM part, used for effective light speed estimation.
KLMTimeResolution * m_timeResolution
DBObject of time resolution.
TH1F * h_time_scint_tc_end
EKLM part.
void createHistograms()
Create histograms.
TGraphErrors * gre_time_channel_scint
BKLM Scintillator.
TH1F * h_timeFSL[2][8][15]
BKLM part.
TH1F * hc_timeFSLPC_end[2][4][14][2][75]
EKLM part.
TH1F * hc_timeFSL_end[2][4][14]
EKLM part.
TH1F * h_timeFSLP_end[2][4][14][2]
EKLM part.
TGraph * gr_timeRes_channel_rpc
BKLM RPC.
TH1F * hc_timeFSLP_end[2][4][14][2]
EKLM part.
TH1F * h_timeFSLP[2][8][15][2]
BKLM part.
TH1F * hc_timeF_scint_end[2]
EKLM part.
std::map< KLMChannelNumber, double > m_timeShift
Shift values of each channel.
TH1F * h_time_scint
BKLM scintillator part.
double m_time_channelAvg_scint
Central value of the global time distribution (BKLM scintillator part).
TH1F * hc_timeFS_scint_end[2][4]
EKLM part.
double esti_timeRes(const KLMChannelIndex &klmChannel)
Estimate value of calibration constant for calibrated channels.
double m_ctime_channelAvg_rpc
Calibrated central value of the global time distribution (BKLM RPC part).
KLMTimeConstants * m_timeConstants
DBObject of time cost on some parts of the detector.
void setupDatabase()
Setup the database.
TH1F * hc_timeFS_scint[2][8]
BKLM scintillator part.
std::map< KLMChannelNumber, double > m_time_channel
Time distribution central value of each channel.
double m_ctime_channelAvg_scint_end
Calibrated central value of the global time distribution (EKLM scintillator part).
CalibrationAlgorithm::EResult readCalibrationData()
Read calibration data.
TGraph * gr_timeShift_channel_scint_end
EKLM.
TGraph * gr_timeRes_channel_scint
BKLM scintillator.
TH1F * hc_timeF_scint[2]
BKLM scintillator part.
TH1F * h_timeFS_scint[2][8]
BKLM scintillator part.
TH1F * hc_timeF_rpc[2]
BKLM RPC part.
TH2F * h2c_timeFS_end[2][4]
EKLM part.
TH1F * h_timeFSLPC_end[2][4][14][2][75]
EKLM part.
TGraphErrors * gre_ctime_channel_scint_end
EKLM.
TProfile * m_Profile2BKLMScintillatorPhi
For BKLM scintillator phi plane.
TH1F * hc_time_scint_end
EKLM part.
TGraphErrors * gre_time_channel_scint_end
EKLM.
TH2F * h2_timeFSLP[2][8][15][2]
BKLM part.
TGraph * gr_timeShift_channel_scint
BKLM scintillator.
double m_time_channelAvg_scint_end
Central value of the global time distribution (EKLM scintillator part).
TProfile * m_Profile2EKLMScintillatorPlane1
For EKLM scintillator plane1.
TH1F * hc_timeFS_rpc[2][8]
BKLM RPC part.
void fillTimeDistanceProfiles(TProfile *profileRpcPhi, TProfile *profileRpcZ, TProfile *profileBKLMScintillatorPhi, TProfile *profileBKLMScintillatorZ, TProfile *profileEKLMScintillatorPlane1, TProfile *profileEKLMScintillatorPlane2, bool fill2dHistograms)
Fill profiles of time versus distance.
TFile * m_outFile
Output file.
double esti_timeShift(const KLMChannelIndex &klmChannel)
Estimate value of calibration constant for uncalibrated channels.
TH1F * hc_timeFSLP[2][8][15][2]
BKLM part.
TGraphErrors * gre_ctime_channel_rpc
BKLM RPC.
void saveHist()
Save histograms to file.
TH2F * h2c_timeFSLP[2][8][15][2]
BKLM part.
double m_ctime_channelAvg_scint
Calibrated central value of the global time distribution (BKLM scintillator part).
TF1 * fcn_const
Const function.
TProfile * m_Profile2RpcZ
For BKLM RPC z plane.
TH1F * h_timeFSLPC[2][8][15][2][54]
BKLM part.
TH2F * h2_timeFSLP_end[2][4][14][2]
EKLM part.
TH1I * hc_calibrated
Calibration statistics for each channel.
void timeDistance2dFit(const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channels, double &delay, double &delayError)
Two-dimensional fit for individual channels.
TGraph * gr_timeRes_channel_scint_end
EKLM.
TH1I * h_calibrated
Calibration statistics for each channel.
TProfile * m_ProfileBKLMScintillatorZ
For BKLM scintillator z plane.
TH1F * h_time_rpc_tc
BKLM RPC part.
TH1F * h_time_scint_end
EKLM part.
TH2F * h2c_timeF_scint[2]
BKLM scintillator part.
TF1 * fcn_pol1
Pol1 function.
double m_etime_channelAvg_scint_end
Central value error of the global time distribution (EKLM scintillator part).
TH1F * hc_time_rpc
BKLM RPC part.
double mc_etime_channelAvg_scint
Calibrated central value error of the global time distribution (BKLM scintillator part).
TH2F * h2c_timeFSLP_end[2][4][14][2]
EKLM part.
double mc_etime_channelAvg_scint_end
Calibrated central value error of the global time distribution (EKLM scintillator part).
TH2F * h2_timeF_scint_end[2]
EKLM part.
TF1 * fcn_gaus
Gaussian function.
TH1F * h_timeF_rpc[2]
BKLM RPC part.
TProfile * m_ProfileBKLMScintillatorPhi
For BKLM scintillator phi plane.
TH1F * hc_time_scint
BKLM scintillator part.
TH2F * h2_timeFS[2][8]
BKLM part.
double m_etime_channelAvg_scint
Central value error of the global time distribution (BKLM scintillator part).
TH1F * h_timeF_scint_end[2]
EKLM part.
TProfile * m_ProfileRpcPhi
For BKLM RPC phi plane.
TGraphErrors * gre_ctime_channel_scint
BKLM Scintillator.
TProfile * m_ProfileEKLMScintillatorPlane2
For EKLM scintillator plane2.
TH1F * h_time_rpc
BKLM RPC part.
TH1F * h_timeFSL_end[2][4][14]
EKLM part.
TProfile * m_Profile2RpcPhi
For BKLM RPC phi plane.
TH1F * h_timeF_scint[2]
BKLM scintillator part.
TH1F * hc_timeFSL[2][8][15]
BKLM part.
TProfile * m_Profile2BKLMScintillatorZ
For BKLM scintillator z plane.
TH1F * h_timeFS_rpc[2][8]
BKLM RPC part.
KLMChannelIndex m_klmChannels
KLM ChannelIndex object.
TGraph * gr_timeShift_channel_rpc
BKLM RPC.
std::map< KLMChannelNumber, double > m_timeRes
Resolution values of each channel.
TH1F * h_diff
Distance between global and local position.
TH2F * h2_timeF_scint[2]
BKLM scintillator part.
TH1F * h_time_scint_tc
BKLM scintillator part.
std::map< KLMChannelNumber, double > m_ctime_channel
Calibrated time distribution central value of each channel.
bool m_useEventT0
Whether to use event T0 from CDC.
TProfile * m_ProfileEKLMScintillatorPlane1
For EKLM scintillator plane1.
TH2F * h2c_timeF_rpc[2]
BKLM RPC part.
TF1 * fcn_land
Landau function.
KLMTimeCableDelay * m_timeCableDelay
DBObject of the calibration constant of each channel due to cable decay.
TProfile * m_ProfileRpcZ
For BKLM RPC z plane.
std::map< KLMChannelNumber, double > mc_etime_channel
Calibrated time distribution central value Error of each channel.
TH1F * hc_timeFSLPC[2][8][15][2][54]
BKLM part.
ROOT::Math::MinimizerOptions m_minimizerOptions
Minimization options.
std::map< KLMChannelNumber, int > m_cFlag
Calibration flag if the channel has enough hits collected and fitted OK.
TGraphErrors * gre_time_channel_rpc
BKLM RPC.
TH2F * h2_timeFS_end[2][4]
EKLM part.
TH2F * h2_timeF_rpc[2]
BKLM RPC part.
std::map< KLMChannelNumber, std::vector< struct Event > > m_evts
Container of hit information.
TH1F * h_timeFS_scint_end[2][4]
EKLM part.
double m_time_channelAvg_rpc
Central value of the global time distribution (BKLM RPC part).
TH2F * h2c_timeFS[2][8]
BKLM part.
double m_etime_channelAvg_rpc
Central value error of the global time distribution (BKLM RPC part).
std::map< KLMChannelNumber, double > m_etime_channel
Time distribution central value Error of each channel.
int m_lower_limit_counts
Lower limit of hits collected for on single channel.
Class to store BKLM delay time coused by cable in the database.
void setTimeDelay(KLMChannelNumber channel, float timeDelay)
Set the time calibration constant.
Class to store KLM constants related to time.
@ c_BKLM
BKLM scintillator.
@ c_EKLM
EKLM scintillator.
void setDelay(float delay, int cType)
Set effective light speed of scintillators.
Class to store KLM time resolution in the database.
void setTimeResolution(KLMChannelNumber channel, float resolution)
Set time resolution.
Class to store variables with their name which were sent to the logging service.
uint16_t KLMChannelNumber
Channel number.

◆ checkPyExpRun()

bool checkPyExpRun ( PyObject *  pyObj)
inherited

Checks that a PyObject can be successfully converted to an ExpRun type.

Checks if the PyObject can be converted to ExpRun.

Definition at line 28 of file CalibrationAlgorithm.cc.

◆ convertPyExpRun()

ExpRun convertPyExpRun ( PyObject *  pyObj)
inherited

Performs the conversion of PyObject to ExpRun.

Converts the PyObject to an ExpRun. We've preoviously checked the object so this assumes a lot about the PyObject.

Definition at line 70 of file CalibrationAlgorithm.cc.

◆ createHistograms()

void createHistograms ( )
private

Create histograms.

Hist declaration Global time distribution

Definition at line 174 of file KLMTimeAlgorithm.cc.

◆ esti_timeRes()

double esti_timeRes ( const KLMChannelIndex klmChannel)

Estimate value of calibration constant for calibrated channels.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1633 of file KLMTimeAlgorithm.cc.

◆ esti_timeShift()

double esti_timeShift ( const KLMChannelIndex klmChannel)

Estimate value of calibration constant for uncalibrated channels.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1553 of file KLMTimeAlgorithm.cc.

◆ execute()

CalibrationAlgorithm::EResult execute ( std::vector< Calibration::ExpRun >  runs = {},
int  iteration = 0,
IntervalOfValidity  iov = IntervalOfValidity() 
)
inherited

Runs calibration over vector of runs for a given iteration.

You can also specify the IoV to save the database payload as. By default the Algorithm will create an IoV from your requested ExpRuns, or from the overall ExpRuns of the input data if you haven't specified ExpRuns in this function.

No checks are performed to make sure that a IoV you specify matches the data you ran over, it simply labels the IoV to commit to the database later.

Definition at line 114 of file CalibrationAlgorithm.cc.

◆ fillTimeDistanceProfiles()

void fillTimeDistanceProfiles ( TProfile *  profileRpcPhi,
TProfile *  profileRpcZ,
TProfile *  profileBKLMScintillatorPhi,
TProfile *  profileBKLMScintillatorZ,
TProfile *  profileEKLMScintillatorPlane1,
TProfile *  profileEKLMScintillatorPlane2,
bool  fill2dHistograms 
)
private

Fill profiles of time versus distance.

Parameters
[out]profileRpcPhiBKLM RPC phi plane.
[out]profileRpcZBKLM RPC z plane.
[out]profileBKLMScintillatorPhiBKLM scintillator phi plane.
[out]profileBKLMScintillatorZBKLM scintillator z plane.
[out]profileEKLMScintillatorPlane1EKLM scintillator plane1.
[out]profileEKLMScintillatorPlane2EKLM scintillator plane2.
[in]fill2dHistogramsWhether to fill 2d histograms.

Definition at line 595 of file KLMTimeAlgorithm.cc.

◆ getCollectorName()

std::string getCollectorName ( ) const
inlineinherited

Alias for prefix.

For convenience and less writing, we say developers to set this to default collector module name in constructor of base class. One can however use the dublets of collector+algorithm multiple times with different settings. To bind these together correctly, the prefix has to be set the same for algo and collector. So we call the setter setPrefix rather than setModuleName or whatever. This getter will work out of the box for default cases -> return the name of module you have to add to your path to collect data for this algorihtm.

Definition at line 164 of file CalibrationAlgorithm.h.

◆ readCalibrationData()

CalibrationAlgorithm::EResult readCalibrationData ( )
private

Read calibration data.

Returns
CalibrationAlgorithm::c_OK if the amount of data is sufficient, CalibrationAlgorithm::c_NotEnoughData otherwise.

Definition at line 141 of file KLMTimeAlgorithm.cc.

◆ setInputFileNames() [1/2]

void setInputFileNames ( PyObject *  inputFileNames)
inherited

Set the input file names used for this algorithm from a Python list.

Set the input file names used for this algorithm and resolve the wildcards.

Definition at line 166 of file CalibrationAlgorithm.cc.

◆ setInputFileNames() [2/2]

void setInputFileNames ( std::vector< std::string >  inputFileNames)
protectedinherited

Set the input file names used for this algorithm.

Set the input file names used for this algorithm and resolve the wildcards.

Definition at line 194 of file CalibrationAlgorithm.cc.

◆ setLowerLimit()

void setLowerLimit ( int  counts)
inline

Set the lower number of hits collected on one sigle strip.

If the hit number is lower than the limit, the strip will not be calibrated and set the average value of the calibration constant.

Definition at line 161 of file KLMTimeAlgorithm.h.

162  {
163  m_lower_limit_counts = counts;
164  }

◆ setMC()

void setMC ( bool  mc)
inline

Set flag indicating whether the input is MC sample.

The histogram ranges are different for data and MC. This setting cannot be determined automatically, because the collector output does not contain metadata.

Definition at line 135 of file KLMTimeAlgorithm.h.

◆ timeDistance2dFit()

void timeDistance2dFit ( const std::vector< std::pair< KLMChannelNumber, unsigned int > > &  channels,
double &  delay,
double &  delayError 
)
private

Two-dimensional fit for individual channels.

Parameters
[in]channelsChannels.
[out]delayDelay (ns / cm).
[out]delayErrorDelay error.

Definition at line 656 of file KLMTimeAlgorithm.cc.

◆ tR_lowerStrip()

std::pair< int, double > tR_lowerStrip ( const KLMChannelIndex klmChannel)

Tracing avaiable channels with decreasing strip number.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1690 of file KLMTimeAlgorithm.cc.

◆ tR_upperStrip()

std::pair< int, double > tR_upperStrip ( const KLMChannelIndex klmChannel)

Tracing avaiable channels with increasing strip number.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1664 of file KLMTimeAlgorithm.cc.

◆ tS_lowerStrip()

std::pair< int, double > tS_lowerStrip ( const KLMChannelIndex klmChannel)

Tracing avaiable channels with decreasing strip number.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1610 of file KLMTimeAlgorithm.cc.

◆ tS_upperStrip()

std::pair< int, double > tS_upperStrip ( const KLMChannelIndex klmChannel)

Tracing avaiable channels with increasing strip number.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1584 of file KLMTimeAlgorithm.cc.

Member Data Documentation

◆ fcn_const

TF1* fcn_const = nullptr
private

Const function.

Global time distribution fitting.

Definition at line 744 of file KLMTimeAlgorithm.h.

◆ fcn_gaus

TF1* fcn_gaus = nullptr
private

Gaussian function.

Scitillator time ditribution fitting.

Definition at line 747 of file KLMTimeAlgorithm.h.

◆ fcn_land

TF1* fcn_land = nullptr
private

Landau function.

RPC time ditribution fitting.

Definition at line 750 of file KLMTimeAlgorithm.h.

◆ fcn_pol1

TF1* fcn_pol1 = nullptr
private

Pol1 function.

Effective light speed fitting.

Definition at line 741 of file KLMTimeAlgorithm.h.

◆ m_evts

std::map<KLMChannelNumber, std::vector<struct Event> > m_evts
private

Container of hit information.

the global element number of the strip is used as the key.

Definition at line 279 of file KLMTimeAlgorithm.h.


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