Belle II Software  release-06-01-15
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...
 
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
 

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 ecah 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.
 
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 = -20.0
 Lower time boundary for BKLM scintillators (calibrated data).
 
double m_UpperTimeBoundaryCalibratedScintilltorsBKLM = 20.0
 Upper time boundary for BKLM scintillators (calibrated data).
 
double m_LowerTimeBoundaryCalibratedScintilltorsEKLM = -20.0
 Lower time boundary for EKLM scintillators (calibrated data).
 
double m_UpperTimeBoundaryCalibratedScintilltorsEKLM = 20.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).
 
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.
 
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.
 
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.
 
TGraph * gr_timeShift_channel_rpc = nullptr
 BKLM RPC.
 
TGraph * gr_timeShift_channel_scint = nullptr
 BKLM scintillator.
 
TGraph * gr_timeShift_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)
 

Detailed Description

KLM time calibration algorithm.

Definition at line 40 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 731 of file KLMTimeAlgorithm.cc.

732 {
733  int channelId;
734  gROOT->SetBatch(kTRUE);
735  setupDatabase();
738 
739  fcn_gaus = new TF1("fcn_gaus", "gaus");
740  fcn_land = new TF1("fcn_land", "landau");
741  fcn_pol1 = new TF1("fcn_pol1", "pol1");
742  fcn_const = new TF1("fcn_const", "pol0");
743 
745  if (result != CalibrationAlgorithm::c_OK)
746  return result;
747 
748  /* Choose non-existing file name. */
749  std::string name = "time_calibration.root";
750  int i = 1;
751  while (1) {
752  struct stat buffer;
753  if (stat(name.c_str(), &buffer) != 0)
754  break;
755  name = "time_calibration_" + std::to_string(i) + ".root";
756  i = i + 1;
757  /* Overflow. */
758  if (i < 0)
759  break;
760  }
761  m_outFile = new TFile(name.c_str(), "recreate");
763 
764  std::vector<struct Event>::iterator it;
765  std::vector<struct Event> eventsChannel;
766 
767  eventsChannel.clear();
768  m_cFlag.clear();
769  m_minimizerOptions.SetDefaultStrategy(2);
770 
771  /* Sort channels by number of events. */
772  std::vector< std::pair<KLMChannelNumber, unsigned int> > channelsBKLM;
773  std::vector< std::pair<KLMChannelNumber, unsigned int> > channelsEKLM;
774  KLMChannelIndex klmChannels;
775  for (KLMChannelIndex& klmChannel : klmChannels) {
776  KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
777  m_cFlag[channel] = ChannelCalibrationStatus::c_NotEnoughData;
778  if (m_evts.find(channel) == m_evts.end())
779  continue;
780  int nEvents = m_evts[channel].size();
781  if (nEvents < m_lower_limit_counts) {
782  B2WARNING("Not enough calibration data collected."
783  << LogVar("channel", channel)
784  << LogVar("number of digit", nEvents));
785  continue;
786  }
787  m_cFlag[channel] = ChannelCalibrationStatus::c_FailedFit;
788  if (klmChannel.getSubdetector() == KLMElementNumbers::c_BKLM &&
789  klmChannel.getLayer() < BKLMElementNumbers::c_FirstRPCLayer) {
790  channelsBKLM.push_back(
791  std::pair<KLMChannelNumber, unsigned int>(channel, nEvents));
792  }
793  if (klmChannel.getSubdetector() == KLMElementNumbers::c_EKLM) {
794  channelsEKLM.push_back(
795  std::pair<KLMChannelNumber, unsigned int>(channel, nEvents));
796  }
797  }
798  std::sort(channelsBKLM.begin(), channelsBKLM.end(), compareEventNumber);
799  std::sort(channelsEKLM.begin(), channelsEKLM.end(), compareEventNumber);
800 
801  /* Two-dimensional fit for the channel with the maximal number of events. */
802  double delayBKLM, delayBKLMError;
803  double delayEKLM, delayEKLMError;
804  timeDistance2dFit(channelsBKLM, delayBKLM, delayBKLMError);
805  timeDistance2dFit(channelsEKLM, delayEKLM, delayEKLMError);
806 
807  /**********************************
808  * First loop
809  * Estimation of effective light speed for Scintillators and RPCs, separately.
810  **********************************/
811  B2INFO("Effective light speed Estimation.");
812  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
813  channelId = klmChannel.getKLMChannelNumber();
814  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
815  continue;
816 
817  eventsChannel = m_evts[channelId];
818  int iSub = klmChannel.getSubdetector();
819 
820  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
821  TVector3 diffD = TVector3(it->diffDistX, it->diffDistY, it->diffDistZ);
822  h_diff->Fill(diffD.Mag());
823  double timeHit = it->time();
824  if (m_useEventT0)
825  timeHit = timeHit - it->t0;
826  if (iSub == KLMElementNumbers::c_BKLM) {
827  int iF = klmChannel.getSection();
828  int iS = klmChannel.getSector() - 1;
829  int iL = klmChannel.getLayer() - 1;
830  int iP = klmChannel.getPlane();
831  int iC = klmChannel.getStrip() - 1;
832  h_timeFSLPC_tc[iF][iS][iL][iP][iC]->Fill(timeHit);
833  if (iL > 1) {
834  h_time_rpc_tc->Fill(timeHit);
835  } else {
836  h_time_scint_tc->Fill(timeHit);
837  }
838  } else {
839  int iF = klmChannel.getSection() - 1;
840  int iS = klmChannel.getSector() - 1;
841  int iL = klmChannel.getLayer() - 1;
842  int iP = klmChannel.getPlane() - 1;
843  int iC = klmChannel.getStrip() - 1;
844  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC]->Fill(timeHit);
845  h_time_scint_tc_end->Fill(timeHit);
846  }
847  }
848  }
849  B2INFO("Effective light speed Estimation! Hists and Graph filling done.");
850 
851  m_timeShift.clear();
852 
853  double tmpMean_rpc_global = h_time_rpc_tc->GetMean();
854  double tmpMean_scint_global = h_time_scint_tc->GetMean();
855  double tmpMean_scint_global_end = h_time_scint_tc_end->GetMean();
856 
857  B2INFO("Global Mean for Raw." << LogVar("RPC", tmpMean_rpc_global) << LogVar("Scint BKLM",
858  tmpMean_scint_global) << LogVar("Scint EKLM", tmpMean_scint_global_end));
859 
860  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
861  channelId = klmChannel.getKLMChannelNumber();
862  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
863  continue;
864 
865  int iSub = klmChannel.getSubdetector();
866  if (iSub == KLMElementNumbers::c_BKLM) {
867  int iF = klmChannel.getSection();
868  int iS = klmChannel.getSector() - 1;
869  int iL = klmChannel.getLayer() - 1;
870  int iP = klmChannel.getPlane();
871  int iC = klmChannel.getStrip() - 1;
872  h_timeFSLPC_tc[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
873  double tmpMean_channel = fcn_gaus->GetParameter(1);
874  if (iL > 1) {
875  m_timeShift[channelId] = tmpMean_channel - tmpMean_rpc_global;
876  } else {
877  m_timeShift[channelId] = tmpMean_channel - tmpMean_scint_global;
878  }
879  } else {
880  int iF = klmChannel.getSection() - 1;
881  int iS = klmChannel.getSector() - 1;
882  int iL = klmChannel.getLayer() - 1;
883  int iP = klmChannel.getPlane() - 1;
884  int iC = klmChannel.getStrip() - 1;
885  h_timeFSLPC_tc_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
886  double tmpMean_channel = fcn_gaus->GetParameter(1);
887  m_timeShift[channelId] = tmpMean_channel - tmpMean_scint_global_end;
888  }
889  }
890 
891  delete h_time_scint_tc;
892  delete h_time_scint_tc_end;
893  delete h_time_rpc_tc;
894  B2INFO("Effective Light m_timeShift obtained. done.");
895 
900 
901  B2INFO("Effective light speed fitting.");
902  m_ProfileRpcPhi->Fit("fcn_pol1", "EMQ");
903  double delayRPCPhi = fcn_pol1->GetParameter(1);
904  double e_slope_rpc_phi = fcn_pol1->GetParError(1);
905 
906  m_ProfileRpcZ->Fit("fcn_pol1", "EMQ");
907  double delayRPCZ = fcn_pol1->GetParameter(1);
908  double e_slope_rpc_z = fcn_pol1->GetParError(1);
909 
910  m_ProfileBKLMScintillatorPhi->Fit("fcn_pol1", "EMQ");
911  double slope_scint_phi = fcn_pol1->GetParameter(1);
912  double e_slope_scint_phi = fcn_pol1->GetParError(1);
913 
914  m_ProfileBKLMScintillatorZ->Fit("fcn_pol1", "EMQ");
915  double slope_scint_z = fcn_pol1->GetParameter(1);
916  double e_slope_scint_z = fcn_pol1->GetParError(1);
917 
918  m_ProfileEKLMScintillatorPlane1->Fit("fcn_pol1", "EMQ");
919  double slope_scint_plane1_end = fcn_pol1->GetParameter(1);
920  double e_slope_scint_plane1_end = fcn_pol1->GetParError(1);
921 
922  m_ProfileEKLMScintillatorPlane2->Fit("fcn_pol1", "EMQ");
923  double slope_scint_plane2_end = fcn_pol1->GetParameter(1);
924  double e_slope_scint_plane2_end = fcn_pol1->GetParError(1);
925 
926  TString logStr_phi, logStr_z;
927  logStr_phi = Form("%.4f +/- %.4f ns/cm", delayRPCPhi, e_slope_rpc_phi);
928  logStr_z = Form("%.4f +/- %.4f ns/cm", delayRPCZ, e_slope_rpc_z);
929  B2INFO("Delay in RPCs:"
930  << LogVar("Fitted Value (phi readout) ", logStr_phi.Data())
931  << LogVar("Fitted Value (z readout) ", logStr_z.Data()));
932  logStr_phi = Form("%.4f +/- %.4f ns/cm", slope_scint_phi, e_slope_scint_phi);
933  logStr_z = Form("%.4f +/- %.4f ns/cm", slope_scint_z, e_slope_scint_z);
934  B2INFO("Delay in BKLM scintillators:"
935  << LogVar("Fitted Value (phi readout) ", logStr_phi.Data())
936  << LogVar("Fitted Value (z readout) ", logStr_z.Data()));
937  logStr_phi = Form("%.4f +/- %.4f ns/cm", slope_scint_plane1_end,
938  e_slope_scint_plane1_end);
939  logStr_z = Form("%.4f +/- %.4f ns/cm", slope_scint_plane2_end,
940  e_slope_scint_plane2_end);
941  B2INFO("Delay in EKLM scintillators:"
942  << LogVar("Fitted Value (plane1 readout) ", logStr_phi.Data())
943  << LogVar("Fitted Value (plane2 readout) ", logStr_z.Data()));
944 
945  logStr_z = Form("%.4f +/- %.4f ns/cm", delayBKLM, delayBKLMError);
946  B2INFO("Delay in BKLM scintillators:"
947  << LogVar("Fitted Value (2d fit) ", logStr_z.Data()));
948  logStr_z = Form("%.4f +/- %.4f ns/cm", delayEKLM, delayEKLMError);
949  B2INFO("Delay in EKLM scintillators:"
950  << LogVar("Fitted Value (2d fit) ", logStr_z.Data()));
951 
952  // Default Effective light speed in current Database
953  //delayEKLM = 0.5 * (slope_scint_plane1_end + slope_scint_plane2_end);
954  //delayBKLM = 0.5 * (slope_scint_phi + slope_scint_z);
955 
960 
962  B2INFO("Time distribution filling begins.");
963  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
964  channelId = klmChannel.getKLMChannelNumber();
965  int iSub = klmChannel.getSubdetector();
966 
967  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
968  continue;
969  eventsChannel = m_evts[channelId];
970 
971  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
972  double timeHit = it->time();
973  if (m_useEventT0)
974  timeHit = timeHit - it->t0;
975  if (iSub == KLMElementNumbers::c_BKLM) {
976  int iF = klmChannel.getSection();
977  int iS = klmChannel.getSector() - 1;
978  int iL = klmChannel.getLayer() - 1;
979  int iP = klmChannel.getPlane();
980  int iC = klmChannel.getStrip() - 1;
981  if (iL > 1) {
982  double propgationT;
984  propgationT = it->dist * delayRPCZ;
985  else
986  propgationT = it->dist * delayRPCPhi;
987  double time = timeHit - propgationT;
988  h_time_rpc->Fill(time);
989  h_timeF_rpc[iF]->Fill(time);
990  h_timeFS_rpc[iF][iS]->Fill(time);
991  h_timeFSL[iF][iS][iL]->Fill(time);
992  h_timeFSLP[iF][iS][iL][iP]->Fill(time);
993  h_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
994  h2_timeF_rpc[iF]->Fill(iS, time);
995  h2_timeFS[iF][iS]->Fill(iL, time);
996  h2_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
997  } else {
998  double propgationT = it->dist * delayBKLM;
999  double time = timeHit - propgationT;
1000  h_time_scint->Fill(time);
1001  h_timeF_scint[iF]->Fill(time);
1002  h_timeFS_scint[iF][iS]->Fill(time);
1003  h_timeFSL[iF][iS][iL]->Fill(time);
1004  h_timeFSLP[iF][iS][iL][iP]->Fill(time);
1005  h_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1006  h2_timeF_scint[iF]->Fill(iS, time);
1007  h2_timeFS[iF][iS]->Fill(iL, time);
1008  h2_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1009  }
1010  } else {
1011  int iF = klmChannel.getSection() - 1;
1012  int iS = klmChannel.getSector() - 1;
1013  int iL = klmChannel.getLayer() - 1;
1014  int iP = klmChannel.getPlane() - 1;
1015  int iC = klmChannel.getStrip() - 1;
1016  double propgationT = it->dist * delayEKLM;
1017  double time = timeHit - propgationT;
1018  h_time_scint_end->Fill(time);
1019  h_timeF_scint_end[iF]->Fill(time);
1020  h_timeFS_scint_end[iF][iS]->Fill(time);
1021  h_timeFSL_end[iF][iS][iL]->Fill(time);
1022  h_timeFSLP_end[iF][iS][iL][iP]->Fill(time);
1023  h_timeFSLPC_end[iF][iS][iL][iP][iC]->Fill(time);
1024  h2_timeF_scint_end[iF]->Fill(iS, time);
1025  h2_timeFS_end[iF][iS]->Fill(iL, time);
1026  h2_timeFSLP_end[iF][iS][iL][iP]->Fill(iC, time);
1027  }
1028  }
1029  }
1030 
1031  B2INFO("Orignal filling done.");
1032 
1033  int iChannel_rpc = 0;
1034  int iChannel = 0;
1035  int iChannel_end = 0;
1036  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1037  channelId = klmChannel.getKLMChannelNumber();
1038  if (m_cFlag[channelId] == ChannelCalibrationStatus::c_NotEnoughData)
1039  continue;
1040  int iSub = klmChannel.getSubdetector();
1041 
1042  if (iSub == KLMElementNumbers::c_BKLM) {
1043  int iF = klmChannel.getSection();
1044  int iS = klmChannel.getSector() - 1;
1045  int iL = klmChannel.getLayer() - 1;
1046  int iP = klmChannel.getPlane();
1047  int iC = klmChannel.getStrip() - 1;
1048 
1049  TFitResultPtr r = h_timeFSLPC[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1050  if (int(r) != 0)
1051  continue;
1052  if (int(r) == 0)
1053  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1054  m_time_channel[channelId] = fcn_gaus->GetParameter(1);
1055  m_etime_channel[channelId] = fcn_gaus->GetParError(1);
1056  if (iL > 1) {
1057  gre_time_channel_rpc->SetPoint(iChannel, channelId, m_time_channel[channelId]);
1058  gre_time_channel_rpc->SetPointError(iChannel, 0., m_etime_channel[channelId]);
1059  iChannel++;
1060  } else {
1061  gre_time_channel_scint->SetPoint(iChannel_rpc, channelId, m_time_channel[channelId]);
1062  gre_time_channel_scint->SetPointError(iChannel_rpc, 0., m_etime_channel[channelId]);
1063  iChannel_rpc++;
1064  }
1065  } else {
1066  int iF = klmChannel.getSection() - 1;
1067  int iS = klmChannel.getSector() - 1;
1068  int iL = klmChannel.getLayer() - 1;
1069  int iP = klmChannel.getPlane() - 1;
1070  int iC = klmChannel.getStrip() - 1;
1071 
1072  TFitResultPtr r = h_timeFSLPC_end[iF][iS][iL][iP][iC]->Fit(fcn_gaus, "LESQ");
1073  if (int(r) != 0)
1074  continue;
1075  if (int(r) == 0)
1076  m_cFlag[channelId] = ChannelCalibrationStatus::c_SuccessfulCalibration;
1077  m_time_channel[channelId] = fcn_gaus->GetParameter(1);
1078  m_etime_channel[channelId] = fcn_gaus->GetParError(1);
1079  gre_time_channel_scint_end->SetPoint(iChannel_end, channelId, m_time_channel[channelId]);
1080  gre_time_channel_scint_end->SetPointError(iChannel_end, 0., m_etime_channel[channelId]);
1081  iChannel_end++;
1082  }
1083  }
1084 
1085  gre_time_channel_scint->Fit("fcn_const", "EMQ");
1086  m_time_channelAvg_scint = fcn_const->GetParameter(0);
1087  m_etime_channelAvg_scint = fcn_const->GetParError(0);
1088 
1089  gre_time_channel_scint_end->Fit("fcn_const", "EMQ");
1090  m_time_channelAvg_scint_end = fcn_const->GetParameter(0);
1091  m_etime_channelAvg_scint_end = fcn_const->GetParError(0);
1092 
1093  gre_time_channel_rpc->Fit("fcn_const", "EMQ");
1094  m_time_channelAvg_rpc = fcn_const->GetParameter(0);
1095  m_etime_channelAvg_rpc = fcn_const->GetParError(0);
1096 
1097  B2INFO("Channel's time distribution fitting done.");
1098  B2DEBUG(20, LogVar("Average time (RPC)", m_time_channelAvg_rpc)
1099  << LogVar("Average time (BKLM scintillators)", m_time_channelAvg_scint)
1100  << LogVar("Average time (EKLM scintillators)", m_time_channelAvg_scint_end));
1101 
1102  B2INFO("Calibrated channel's time distribution filling begins.");
1103 
1104  m_timeShift.clear();
1105  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1106  channelId = klmChannel.getKLMChannelNumber();
1107  h_calibrated->Fill(m_cFlag[channelId]);
1108  if (m_time_channel.find(channelId) == m_time_channel.end())
1109  continue;
1110  double timeShift = m_time_channel[channelId];
1111  int iSub = klmChannel.getSubdetector();
1112  if (iSub == KLMElementNumbers::c_BKLM) {
1113  int iL = klmChannel.getLayer() - 1;
1114  if (iL > 1)
1115  m_timeShift[channelId] = timeShift;
1116  else
1117  m_timeShift[channelId] = timeShift;
1118  } else {
1119  m_timeShift[channelId] = timeShift;
1120  }
1121  m_timeCableDelay->setTimeDelay(channelId, m_timeShift[channelId]);
1122  }
1123 
1124  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1125  channelId = klmChannel.getKLMChannelNumber();
1126  if (m_timeShift.find(channelId) != m_timeShift.end())
1127  continue;
1128  m_timeShift[channelId] = esti_timeShift(klmChannel);
1129  m_timeCableDelay->setTimeDelay(channelId, m_timeShift[channelId]);
1130  B2DEBUG(20, "Uncalibrated Estimation " << LogVar("Channel", channelId) << LogVar("Estimated value", m_timeShift[channelId]));
1131  }
1132 
1133  iChannel_rpc = 0;
1134  iChannel = 0;
1135  iChannel_end = 0;
1136  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1137  channelId = klmChannel.getKLMChannelNumber();
1138  if (m_timeShift.find(channelId) == m_timeShift.end()) {
1139  B2ERROR("!!! Not All Channels Calibration Constant Set. Error Happended on " << LogVar("Channel", channelId));
1140  continue;
1141  }
1142  int iSub = klmChannel.getSubdetector();
1143  if (iSub == KLMElementNumbers::c_BKLM) {
1144  int iL = klmChannel.getLayer();
1145  if (iL > 2) {
1146  gr_timeShift_channel_rpc->SetPoint(iChannel_rpc, channelId, m_timeShift[channelId]);
1147  iChannel_rpc++;
1148  } else {
1149  gr_timeShift_channel_scint->SetPoint(iChannel, channelId, m_timeShift[channelId]);
1150  iChannel++;
1151  }
1152  } else {
1153  gr_timeShift_channel_scint_end->SetPoint(iChannel_end, channelId, m_timeShift[channelId]);
1154  iChannel_end++;
1155  }
1156  }
1157 
1162  for (KLMChannelIndex klmChannel = m_klmChannels.begin(); klmChannel != m_klmChannels.end(); ++klmChannel) {
1163  channelId = klmChannel.getKLMChannelNumber();
1164  int iSub = klmChannel.getSubdetector();
1165  eventsChannel = m_evts[channelId];
1166  for (it = eventsChannel.begin(); it != eventsChannel.end(); ++it) {
1167  double timeHit = it->time();
1168  if (m_useEventT0)
1169  timeHit = timeHit - it->t0;
1170  if (iSub == KLMElementNumbers::c_BKLM) {
1171  int iF = klmChannel.getSection();
1172  int iS = klmChannel.getSector() - 1;
1173  int iL = klmChannel.getLayer() - 1;
1174  int iP = klmChannel.getPlane();
1175  int iC = klmChannel.getStrip() - 1;
1176  if (iL > 1) {
1177  double propgationT;
1178  if (iP == BKLMElementNumbers::c_ZPlane)
1179  propgationT = it->dist * delayRPCZ;
1180  else
1181  propgationT = it->dist * delayRPCPhi;
1182  double time = timeHit - propgationT - m_timeShift[channelId];
1183  hc_time_rpc->Fill(time);
1184  hc_timeF_rpc[iF]->Fill(time);
1185  hc_timeFS_rpc[iF][iS]->Fill(time);
1186  hc_timeFSL[iF][iS][iL]->Fill(time);
1187  hc_timeFSLP[iF][iS][iL][iP]->Fill(time);
1188  hc_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1189  h2c_timeF_rpc[iF]->Fill(iS, time);
1190  h2c_timeFS[iF][iS]->Fill(iL, time);
1191  h2c_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1192  } else {
1193  double propgationT = it->dist * delayBKLM;
1194  double time = timeHit - propgationT - m_timeShift[channelId];
1195  hc_time_scint->Fill(time);
1196  hc_timeF_scint[iF]->Fill(time);
1197  hc_timeFS_scint[iF][iS]->Fill(time);
1198  hc_timeFSL[iF][iS][iL]->Fill(time);
1199  hc_timeFSLP[iF][iS][iL][iP]->Fill(time);
1200  hc_timeFSLPC[iF][iS][iL][iP][iC]->Fill(time);
1201  h2c_timeF_scint[iF]->Fill(iS, time);
1202  h2c_timeFS[iF][iS]->Fill(iL, time);
1203  h2c_timeFSLP[iF][iS][iL][iP]->Fill(iC, time);
1204  }
1205  } else {
1206  int iF = klmChannel.getSection() - 1;
1207  int iS = klmChannel.getSector() - 1;
1208  int iL = klmChannel.getLayer() - 1;
1209  int iP = klmChannel.getPlane() - 1;
1210  int iC = klmChannel.getStrip() - 1;
1211  double propgationT = it->dist * delayEKLM;
1212  double time = timeHit - propgationT - m_timeShift[channelId];
1213  hc_time_scint_end->Fill(time);
1214  hc_timeF_scint_end[iF]->Fill(time);
1215  hc_timeFS_scint_end[iF][iS]->Fill(time);
1216  hc_timeFSL_end[iF][iS][iL]->Fill(time);
1217  hc_timeFSLP_end[iF][iS][iL][iP]->Fill(time);
1218  hc_timeFSLPC_end[iF][iS][iL][iP][iC]->Fill(time);
1219  h2c_timeF_scint_end[iF]->Fill(iS, time);
1220  h2c_timeFS_end[iF][iS]->Fill(iL, time);
1221  h2c_timeFSLP_end[iF][iS][iL][iP]->Fill(iC, time);
1222  }
1223  }
1224  }
1225 
1226  saveHist();
1227 
1228  delete fcn_const;
1229  m_evts.clear();
1230  m_timeShift.clear();
1231  m_cFlag.clear();
1232 
1233  saveCalibration(m_timeCableDelay, "KLMTimeCableDelay");
1234  saveCalibration(m_timeConstants, "KLMTimeConstants");
1236 }
@ 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.
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.
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.
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 ecah 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.
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.
CalibrationAlgorithm::EResult readCalibrationData()
Read calibration data.
TGraph * gr_timeShift_channel_scint_end
EKLM.
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.
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.
void saveHist()
Save histograms to file.
TH2F * h2c_timeFSLP[2][8][15][2]
BKLM 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.
void timeDistance2dFit(const std::vector< std::pair< KLMChannelNumber, unsigned int > > &channels, double &delay, double &delayError)
Two-dimensional fit for individual channels.
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.
TH2F * h2c_timeFSLP_end[2][4][14][2]
EKLM 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.
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.
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.
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.
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 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_timeShift()

double esti_timeShift ( const KLMChannelIndex klmChannel)

Estimate value of calibration constant for uncalibrated channels.

Parameters
[in]klmChannelKLM channel index.

Definition at line 1409 of file KLMTimeAlgorithm.cc.

◆ execute()

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.

◆ 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 586 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 160 of file KLMTimeAlgorithm.h.

161  {
162  m_lower_limit_counts = counts;
163  }

◆ 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 134 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 647 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 1466 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 1440 of file KLMTimeAlgorithm.cc.

Member Data Documentation

◆ fcn_const

TF1* fcn_const = nullptr
private

Const function.

Global time distribution fitting.

Definition at line 672 of file KLMTimeAlgorithm.h.

◆ fcn_gaus

TF1* fcn_gaus = nullptr
private

Gaussian function.

Scitillator time ditribution fitting.

Definition at line 675 of file KLMTimeAlgorithm.h.

◆ fcn_land

TF1* fcn_land = nullptr
private

Landau function.

RPC time ditribution fitting.

Definition at line 678 of file KLMTimeAlgorithm.h.

◆ fcn_pol1

TF1* fcn_pol1 = nullptr
private

Pol1 function.

Effective light speed fitting.

Definition at line 669 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 260 of file KLMTimeAlgorithm.h.


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