Belle II Software  release-05-02-19
TOPCommonT0BFAlgorithm.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * *
6  * Author: The Belle II Collaboration *
7  * Contributors: Marko Staric *
8  * *
9  * This software is provided "as is" without any warranty. *
10  **************************************************************************/
11 
12 #pragma once
13 #include <calibration/CalibrationAlgorithm.h>
14 #include <TH1F.h>
15 
16 namespace Belle2 {
21  namespace TOP {
22 
27  class TOPCommonT0BFAlgorithm : public CalibrationAlgorithm {
28  public:
29 
32 
34  virtual ~TOPCommonT0BFAlgorithm() {}
35 
39  void setMinEntries(int minEntries) {m_minEntries = minEntries;}
40 
47  void setFitInitializers(double sigmaCore, double sigmaTail, double tailFract)
48  {
49  m_sigmaCoreInit = sigmaCore;
50  m_sigmaTailInit = sigmaTail;
51  m_tailFractInit = tailFract;
52  }
53 
59  void setCutoffEntries(int cutoffEntries) {m_cutoffEntries = cutoffEntries;}
60 
65  void setMinError(double minError) {m_minError = minError;}
66 
67  private:
68 
72  virtual EResult calibrate() final;
73 
77  std::shared_ptr<TH1F> getHistogram();
78 
84  int fitSingleGaus(std::shared_ptr<TH1F> h);
85 
91  int fitDoubleGaus(std::shared_ptr<TH1F> h);
92 
93  // algorithm parameters
94  int m_minEntries = 10;
95  double m_sigmaCoreInit = 0.060;
96  double m_sigmaTailInit = 0.120;
97  double m_tailFractInit = 0.20;
98  int m_cutoffEntries = 100;
99  double m_minError = 0.020;
101  // ntuple variables
102  int m_expNo = 0;
103  int m_runNo = 0;
104  int m_runLast = 0;
105  float m_fittedOffset = 0;
106  float m_offset = 0;
107  float m_offsetError = 0;
108  float m_sigmaCore = 0;
109  float m_sigmaTail = 0;
110  float m_tailFract = 0;
111  float m_chi2 = 0;
112  float m_integral = 0;
113  int m_numEvents = 0;
114  int m_fitStatus = -1;
116  // other
117  double m_bunchTimeSep = 0;
119  };
120 
121  } // end namespace TOP
123 } // end namespace Belle2
Belle2::TOP::TOPCommonT0BFAlgorithm::fitDoubleGaus
int fitDoubleGaus(std::shared_ptr< TH1F > h)
Fit double gaus w/ same mean + constant.
Definition: TOPCommonT0BFAlgorithm.cc:170
Belle2::TOP::TOPCommonT0BFAlgorithm::m_tailFract
float m_tailFract
tail fraction (set to 0 for single gauss fit)
Definition: TOPCommonT0BFAlgorithm.h:119
Belle2::TOP::TOPCommonT0BFAlgorithm::m_offsetError
float m_offsetError
error on fitted offset
Definition: TOPCommonT0BFAlgorithm.h:116
Belle2::TOP::TOPCommonT0BFAlgorithm::m_chi2
float m_chi2
normalized chi^2
Definition: TOPCommonT0BFAlgorithm.h:120
Belle2::TOP::TOPCommonT0BFAlgorithm::calibrate
virtual EResult calibrate() final
algorithm implementation
Definition: TOPCommonT0BFAlgorithm.cc:38
Belle2::TOP::TOPCommonT0BFAlgorithm::setCutoffEntries
void setCutoffEntries(int cutoffEntries)
Sets cutoff on the number of histogram entries for steering btw.
Definition: TOPCommonT0BFAlgorithm.h:68
Belle2::TOP::TOPCommonT0BFAlgorithm::m_numEvents
int m_numEvents
number of all events in the histogram
Definition: TOPCommonT0BFAlgorithm.h:122
Belle2::TOP::TOPCommonT0BFAlgorithm::getHistogram
std::shared_ptr< TH1F > getHistogram()
Returns histogram to be fitted.
Definition: TOPCommonT0BFAlgorithm.cc:115
Belle2::TOP::TOPCommonT0BFAlgorithm::setMinEntries
void setMinEntries(int minEntries)
Sets minimal number of histogram entries to perform a fit.
Definition: TOPCommonT0BFAlgorithm.h:48
Belle2::TOP::TOPCommonT0BFAlgorithm::m_sigmaCore
float m_sigmaCore
core gaussian sigma
Definition: TOPCommonT0BFAlgorithm.h:117
Belle2::TOP::TOPCommonT0BFAlgorithm::m_fitStatus
int m_fitStatus
fit status
Definition: TOPCommonT0BFAlgorithm.h:123
Belle2::TOP::TOPCommonT0BFAlgorithm::setFitInitializers
void setFitInitializers(double sigmaCore, double sigmaTail, double tailFract)
Sets values for the initialization of several fit parameters.
Definition: TOPCommonT0BFAlgorithm.h:56
Belle2::TOP::TOPCommonT0BFAlgorithm::m_sigmaTail
float m_sigmaTail
tail gaussian sigma (set to 0 for single gauss fit)
Definition: TOPCommonT0BFAlgorithm.h:118
Belle2::TOP::TOPCommonT0BFAlgorithm::m_fittedOffset
float m_fittedOffset
fitted offset
Definition: TOPCommonT0BFAlgorithm.h:114
Belle2::TOP::TOPCommonT0BFAlgorithm::m_expNo
int m_expNo
experiment number
Definition: TOPCommonT0BFAlgorithm.h:111
Belle2::TOP::TOPCommonT0BFAlgorithm::m_runLast
int m_runLast
last run number
Definition: TOPCommonT0BFAlgorithm.h:113
Belle2::TOP::TOPCommonT0BFAlgorithm::m_minError
double m_minError
minimal commonT0 uncertainty [ns] to declare c_OK
Definition: TOPCommonT0BFAlgorithm.h:108
Belle2::TOP::TOPCommonT0BFAlgorithm::m_integral
float m_integral
number of fitted signal events
Definition: TOPCommonT0BFAlgorithm.h:121
Belle2::TOP::TOPCommonT0BFAlgorithm::m_runNo
int m_runNo
first run number
Definition: TOPCommonT0BFAlgorithm.h:112
Belle2::TOP::TOPCommonT0BFAlgorithm::m_offset
float m_offset
wrap-around of fitted offset (= common T0)
Definition: TOPCommonT0BFAlgorithm.h:115
Belle2::TOP::TOPCommonT0BFAlgorithm::setMinError
void setMinError(double minError)
Sets minimal fitted offset uncertainty to declare this calibration as c_OK.
Definition: TOPCommonT0BFAlgorithm.h:74
Belle2::TOP::TOPCommonT0BFAlgorithm::m_cutoffEntries
int m_cutoffEntries
cutoff entries for single/double gaussian fit
Definition: TOPCommonT0BFAlgorithm.h:107
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TOP::TOPCommonT0BFAlgorithm::m_sigmaTailInit
double m_sigmaTailInit
tail gaussian sigma [ns]
Definition: TOPCommonT0BFAlgorithm.h:105
Belle2::TOP::TOPCommonT0BFAlgorithm::TOPCommonT0BFAlgorithm
TOPCommonT0BFAlgorithm()
Constructor.
Definition: TOPCommonT0BFAlgorithm.cc:30
Belle2::TOP::TOPCommonT0BFAlgorithm::m_tailFractInit
double m_tailFractInit
fraction of tail gaussian
Definition: TOPCommonT0BFAlgorithm.h:106
Belle2::TOP::TOPCommonT0BFAlgorithm::m_minEntries
int m_minEntries
minimal number of histogram entries to perform fit
Definition: TOPCommonT0BFAlgorithm.h:103
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::TOP::TOPCommonT0BFAlgorithm::~TOPCommonT0BFAlgorithm
virtual ~TOPCommonT0BFAlgorithm()
Destructor.
Definition: TOPCommonT0BFAlgorithm.h:43
Belle2::TOP::TOPCommonT0BFAlgorithm::fitSingleGaus
int fitSingleGaus(std::shared_ptr< TH1F > h)
Fit single gaus + constant.
Definition: TOPCommonT0BFAlgorithm.cc:136
Belle2::TOP::TOPCommonT0BFAlgorithm::m_bunchTimeSep
double m_bunchTimeSep
bunch separation in time [ns]
Definition: TOPCommonT0BFAlgorithm.h:126
Belle2::TOP::TOPCommonT0BFAlgorithm::m_sigmaCoreInit
double m_sigmaCoreInit
core gaussian sigma [ns]
Definition: TOPCommonT0BFAlgorithm.h:104