Belle II Software  release-05-02-19
TOPTimeBaseCalibratorModule.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Marko Staric, Hiromichi Kichimi and Xiaolong Wang *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #include <framework/core/Module.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <top/dataobjects/TOPDigit.h>
16 
17 #include <string>
18 #include <vector>
19 #include <TMatrixDSym.h>
20 #include <TH1F.h>
21 #include <TH2F.h>
22 
23 namespace Belle2 {
32  struct Hit {
33  double time = 0;
34  double timeErr = 0;
35  int pulseHeight = 0;
36  double pulseWidth = 0;
41  Hit(double t, double et, int height, double width)
42  {
43  time = t;
44  timeErr = et;
45  pulseHeight = height;
46  pulseWidth = width;
47  }
48  };
49 
53  struct TwoTimes {
54  float t1 = 0;
55  float t2 = 0;
56  float sigma = 0;
57  bool good = false;
62  TwoTimes(double time1, double time2, double sig)
63  {
64  if (time1 < time2) {
65  t1 = time1;
66  t2 = time2;
67  } else {
68  t1 = time2;
69  t2 = time1;
70  }
71  sigma = sig;
72  good = true;
73  }
74  };
75 
79  class TOPTimeBaseCalibratorModule : public Module {
80 
81  public:
82 
87 
92 
97  virtual void initialize() override;
98 
103  virtual void beginRun() override;
104 
108  virtual void event() override;
109 
114  virtual void endRun() override;
115 
120  virtual void terminate() override;
121 
122  private:
123 
127  enum {c_NumChannels = 512,
128  c_NumBoardstacks = 4,
129  c_NumScrodChannels = c_NumChannels / c_NumBoardstacks,
130  c_WindowSize = 64,
132  };
133 
144  bool calibrateChannel(std::vector<TwoTimes>& ntuple,
145  unsigned scrodID, unsigned scrodChannel,
146  TH1F& Hchi2, TH1F& Hndf, TH1F& HDeltaT);
147 
159  bool matrixInversion(const std::vector<TwoTimes>& ntuple,
160  unsigned scrodID, unsigned scrodChannel,
161  double meanTimeDifference,
162  TH1F& Hchi2, TH1F& Hndf, TH1F& HDeltaT);
163 
176  bool iterativeTBC(const std::vector<TwoTimes>& ntuple,
177  unsigned scrodID, unsigned scrodChannel,
178  double meanTimeDifference,
179  TH1F& Hchi2, TH1F& Hndf, TH1F& HDeltaT);
180 
187  void Iteration(const std::vector<TwoTimes>& ntuple, std::vector<double>& xval);
188 
194  double Chisq(const std::vector<TwoTimes>& ntuple, const std::vector<double>& xxval);
195 
204  void saveAsHistogram(const std::vector<double>& vec,
205  const std::string& name,
206  const std::string& title,
207  const std::string& xTitle = "",
208  const std::string& yTitle = "") const;
209 
219  void saveAsHistogram(const std::vector<double>& vec,
220  const std::vector<double>& err,
221  const std::string& name,
222  const std::string& title,
223  const std::string& xTitle = "",
224  const std::string& yTitle = "") const;
225 
232  void saveAsHistogram(const TMatrixDSym& M,
233  const std::string& name,
234  const std::string& title) const;
235 
236 
237  int m_moduleID = 0;
238  double m_minTimeDiff = 0;
239  double m_maxTimeDiff = 0;
240  unsigned m_minHits = 0;
241  unsigned m_numIterations = 100;
242  std::string m_directoryName;
243  unsigned m_method = 0;
244  bool m_useFallingEdge = false;
246  std::vector<TwoTimes> m_ntuples[c_NumChannels];
247  double m_syncTimeBase = 0;
249  TH2F m_goodHits;
253  // parameters for iTBC
254  int m_good = 0;
255  double m_dt_min = 20.0;
256  double m_dt_max = 24.0;
257  double m_dev_step = 0.001;
258  double m_xstep = 0.020;
259  double m_dchi2dxv = 0.0;
260  double m_change_xstep = 0.015;
261  double m_new_xstep = 2.0 * m_xstep;
262  double m_min_binwidth = 0.05;
263  double m_max_binwidth = 2.0;
264  double m_dchi2_min = 0.2;
265  unsigned m_conv_iter = 100;
266  double m_deltamin = 0.01 * 0.01;
267  double m_DtSigma = 0.0424;
270  };
271 
273 } // Belle2 namespace
274 
Belle2::TOPTimeBaseCalibratorModule::m_calPulseFirst
TH2F m_calPulseFirst
pulse height versus width of the first calibration pulse
Definition: TOPTimeBaseCalibratorModule.h:258
Belle2::TOPTimeBaseCalibratorModule::m_dchi2dxv
double m_dchi2dxv
rms of 255 dchi2/dxval values
Definition: TOPTimeBaseCalibratorModule.h:267
Belle2::TOPTimeBaseCalibratorModule::Iteration
void Iteration(const std::vector< TwoTimes > &ntuple, std::vector< double > &xval)
Iteration function called by iterativeTBC()
Definition: TOPTimeBaseCalibratorModule.cc:620
Belle2::TOPTimeBaseCalibratorModule::m_calPulseSecond
TH2F m_calPulseSecond
pulse height versus width of the second calibration pulse
Definition: TOPTimeBaseCalibratorModule.h:259
Belle2::TOPTimeBaseCalibratorModule::Chisq
double Chisq(const std::vector< TwoTimes > &ntuple, const std::vector< double > &xxval)
Return the chisqure of one set of TBC constants (xval) in iTBC calculaton.
Definition: TOPTimeBaseCalibratorModule.cc:666
Belle2::TOPTimeBaseCalibratorModule::c_WindowSize
@ c_WindowSize
samples per window
Definition: TOPTimeBaseCalibratorModule.h:138
Belle2::Hit::timeErr
double timeErr
raw time uncertainty [samples]
Definition: TOPTimeBaseCalibratorModule.h:42
Belle2::TOPTimeBaseCalibratorModule::c_TimeAxisSize
@ c_TimeAxisSize
number of samples to calibrate
Definition: TOPTimeBaseCalibratorModule.h:139
Belle2::TOPTimeBaseCalibratorModule::m_maxTimeDiff
double m_maxTimeDiff
upper bound for time difference [samples]
Definition: TOPTimeBaseCalibratorModule.h:247
Belle2::TwoTimes::good
bool good
flag
Definition: TOPTimeBaseCalibratorModule.h:65
Belle2::TOPTimeBaseCalibratorModule::event
virtual void event() override
Event processor.
Definition: TOPTimeBaseCalibratorModule.cc:146
Belle2::TOPTimeBaseCalibratorModule::m_moduleID
int m_moduleID
slot ID
Definition: TOPTimeBaseCalibratorModule.h:245
Belle2::Hit::pulseHeight
int pulseHeight
pulse height [ADC counts]
Definition: TOPTimeBaseCalibratorModule.h:43
Belle2::TOPTimeBaseCalibratorModule::iterativeTBC
bool iterativeTBC(const std::vector< TwoTimes > &ntuple, unsigned scrodID, unsigned scrodChannel, double meanTimeDifference, TH1F &Hchi2, TH1F &Hndf, TH1F &HDeltaT)
Method by iteration of chi2 minimization.
Definition: TOPTimeBaseCalibratorModule.cc:535
Belle2::TOPTimeBaseCalibratorModule::terminate
virtual void terminate() override
Termination action.
Definition: TOPTimeBaseCalibratorModule.cc:219
Belle2::TOPTimeBaseCalibratorModule::c_NumBoardstacks
@ c_NumBoardstacks
number of boardstacks per module
Definition: TOPTimeBaseCalibratorModule.h:136
Belle2::TOPTimeBaseCalibratorModule::m_max_binwidth
double m_max_binwidth
maximum time interval of one sample
Definition: TOPTimeBaseCalibratorModule.h:271
Belle2::TOPTimeBaseCalibratorModule::m_dchi2_min
double m_dchi2_min
quit if chisq increase in iteratons is larger than this value.
Definition: TOPTimeBaseCalibratorModule.h:272
Belle2::Hit::time
double time
raw time [samples]
Definition: TOPTimeBaseCalibratorModule.h:41
Belle2::TOPTimeBaseCalibratorModule::m_goodHits
TH2F m_goodHits
pulse height versus width of all good hits
Definition: TOPTimeBaseCalibratorModule.h:257
Belle2::TOPTimeBaseCalibratorModule::beginRun
virtual void beginRun() override
Called when entering a new run.
Definition: TOPTimeBaseCalibratorModule.cc:142
Belle2::TOPTimeBaseCalibratorModule::c_NumChannels
@ c_NumChannels
number of channels per module
Definition: TOPTimeBaseCalibratorModule.h:135
Belle2::Hit::Hit
Hit(double t, double et, int height, double width)
Full constructor.
Definition: TOPTimeBaseCalibratorModule.h:49
Belle2::TOPTimeBaseCalibratorModule::m_min_binwidth
double m_min_binwidth
minimum time interval of one sample
Definition: TOPTimeBaseCalibratorModule.h:270
Belle2::TOPTimeBaseCalibratorModule::m_change_xstep
double m_change_xstep
update m_xstep if m_dchi2dxv < m_change_step
Definition: TOPTimeBaseCalibratorModule.h:268
Belle2::TOPTimeBaseCalibratorModule::m_sigm2_exp
double m_sigm2_exp
(sigma_0(dT))**2 for nomarlization of chisq = sum{dT^2/sigma^2}
Definition: TOPTimeBaseCalibratorModule.h:276
Belle2::TOPTimeBaseCalibratorModule::~TOPTimeBaseCalibratorModule
virtual ~TOPTimeBaseCalibratorModule()
Destructor.
Definition: TOPTimeBaseCalibratorModule.cc:101
Belle2::TOPTimeBaseCalibratorModule::endRun
virtual void endRun() override
End-of-run action.
Definition: TOPTimeBaseCalibratorModule.cc:215
Belle2::TOPTimeBaseCalibratorModule::m_minHits
unsigned m_minHits
minimal required hits per channel
Definition: TOPTimeBaseCalibratorModule.h:248
Belle2::TwoTimes
Structure to hold calpulse raw times expressed in samples since sample 0 of window 0.
Definition: TOPTimeBaseCalibratorModule.h:61
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TOPTimeBaseCalibratorModule::m_ntuples
std::vector< TwoTimes > m_ntuples[c_NumChannels]
channel wise data
Definition: TOPTimeBaseCalibratorModule.h:254
Belle2::TOPTimeBaseCalibratorModule::m_dev_step
double m_dev_step
a step size to calculate the value of d(chisq)/dxval
Definition: TOPTimeBaseCalibratorModule.h:265
Belle2::TOPTimeBaseCalibratorModule::saveAsHistogram
void saveAsHistogram(const std::vector< double > &vec, const std::string &name, const std::string &title, const std::string &xTitle="", const std::string &yTitle="") const
Save vector to histogram and write it out.
Definition: TOPTimeBaseCalibratorModule.cc:708
Belle2::TwoTimes::sigma
float sigma
uncertainty of time difference (r.m.s.)
Definition: TOPTimeBaseCalibratorModule.h:64
Belle2::TOPTimeBaseCalibratorModule::m_useFallingEdge
bool m_useFallingEdge
if true, use falling edge instead of rising
Definition: TOPTimeBaseCalibratorModule.h:252
Belle2::TOPTimeBaseCalibratorModule::matrixInversion
bool matrixInversion(const std::vector< TwoTimes > &ntuple, unsigned scrodID, unsigned scrodChannel, double meanTimeDifference, TH1F &Hchi2, TH1F &Hndf, TH1F &HDeltaT)
Method by matrix inversion.
Definition: TOPTimeBaseCalibratorModule.cc:397
Belle2::TOPTimeBaseCalibratorModule::m_DtSigma
double m_DtSigma
a reference resolution of sigma_0(dT)=42.4 ps from Run3524
Definition: TOPTimeBaseCalibratorModule.h:275
Belle2::TOPTimeBaseCalibratorModule::m_minTimeDiff
double m_minTimeDiff
lower bound for time difference [samples]
Definition: TOPTimeBaseCalibratorModule.h:246
Belle2::TwoTimes::t2
float t2
time of the second pulse [samples]
Definition: TOPTimeBaseCalibratorModule.h:63
Belle2::TOPTimeBaseCalibratorModule::m_new_xstep
double m_new_xstep
m_xstep = m_new_xstep if m_dchi2dxv < m_change_step
Definition: TOPTimeBaseCalibratorModule.h:269
Belle2::TwoTimes::t1
float t1
time of the first pulse [samples]
Definition: TOPTimeBaseCalibratorModule.h:62
Belle2::TOPTimeBaseCalibratorModule::m_dt_max
double m_dt_max
maximum Delta T of raw calpulse
Definition: TOPTimeBaseCalibratorModule.h:264
Belle2::TOPTimeBaseCalibratorModule::m_good
int m_good
good events used for chisq cal.
Definition: TOPTimeBaseCalibratorModule.h:262
Belle2::TOPTimeBaseCalibratorModule::m_syncTimeBase
double m_syncTimeBase
synchronization time (two ASIC windows)
Definition: TOPTimeBaseCalibratorModule.h:255
Belle2::TOPTimeBaseCalibratorModule::m_conv_iter
unsigned m_conv_iter
Number of iteration with chisq changes less than deltamin.
Definition: TOPTimeBaseCalibratorModule.h:273
Belle2::TOPTimeBaseCalibratorModule::initialize
virtual void initialize() override
Initialize the Module.
Definition: TOPTimeBaseCalibratorModule.cc:105
Belle2::TOPTimeBaseCalibratorModule::calibrateChannel
bool calibrateChannel(std::vector< TwoTimes > &ntuple, unsigned scrodID, unsigned scrodChannel, TH1F &Hchi2, TH1F &Hndf, TH1F &HDeltaT)
calibrate single channel
Definition: TOPTimeBaseCalibratorModule.cc:310
Belle2::TOPTimeBaseCalibratorModule::m_xstep
double m_xstep
unit for an interation of delta(X_s)
Definition: TOPTimeBaseCalibratorModule.h:266
Belle2::TOPTimeBaseCalibratorModule::m_digits
StoreArray< TOPDigit > m_digits
collection of digits
Definition: TOPTimeBaseCalibratorModule.h:256
Belle2::TOPTimeBaseCalibratorModule::TOPTimeBaseCalibratorModule
TOPTimeBaseCalibratorModule()
Constructor.
Definition: TOPTimeBaseCalibratorModule.cc:54
Belle2::TOPTimeBaseCalibratorModule::m_directoryName
std::string m_directoryName
directory name for the output root files
Definition: TOPTimeBaseCalibratorModule.h:250
Belle2::TOPTimeBaseCalibratorModule::m_dt_min
double m_dt_min
minimum Delta T of raw calpulse
Definition: TOPTimeBaseCalibratorModule.h:263
Belle2::TOPTimeBaseCalibratorModule::m_method
unsigned m_method
method to use
Definition: TOPTimeBaseCalibratorModule.h:251
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::TOPTimeBaseCalibratorModule::m_numIterations
unsigned m_numIterations
Number of Iterations of iTBC.
Definition: TOPTimeBaseCalibratorModule.h:249
Belle2::Hit::pulseWidth
double pulseWidth
pulse width [ns]
Definition: TOPTimeBaseCalibratorModule.h:44
Belle2::TwoTimes::TwoTimes
TwoTimes(double time1, double time2, double sig)
Full constructor.
Definition: TOPTimeBaseCalibratorModule.h:70
Belle2::TOPTimeBaseCalibratorModule::m_deltamin
double m_deltamin
minumum chisq change in an iteration.
Definition: TOPTimeBaseCalibratorModule.h:274