Belle II Software development
SVDdEdxCalibrationAlgorithm.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11#include <calibration/CalibrationAlgorithm.h>
12
13#include <TH2F.h>
14#include <tuple>
15
16namespace Belle2 {
25 public:
31
37
41 void setMonitoringPlots(bool value = false) { m_isMakePlots = value; }
42
46 void setNumDEdxBins(const int& value) { m_numDEdxBins = value; }
47
51 void setNumPBins(const int& value) { m_numPBins = value; }
52
56 void setDEdxCutoff(const double& value) { m_dedxCutoff = value; }
57
61 void setMinEvtsPerTree(const double& value) { m_MinEvtsPerTree = value; }
62
63 protected:
67 virtual EResult calibrate() override;
68
69 private:
71 TH2F LambdaMassFit(std::shared_ptr<TTree> preselTree);
72 std::tuple<TH2F, TH2F, TH2F> DstarMassFit(std::shared_ptr<TTree> preselTree);
73 TH2F GammaHistogram(std::shared_ptr<TTree> preselTree);
74 int m_numDEdxBins = 100;
75 int m_numPBins = 69;
76 double m_dedxCutoff = 5.e6;
78 100;
82 std::vector<double> CreatePBinningScheme()
83 {
84 std::vector<double> pbins;
85 pbins.reserve(m_numPBins + 1);
86 pbins.push_back(0.0);
87 pbins.push_back(0.05);
88
89 for (int iBin = 2; iBin <= m_numPBins; iBin++) {
90 if (iBin <= 19)
91 pbins.push_back(0.025 + 0.025 * iBin);
92 else if (iBin <= 59)
93 pbins.push_back(pbins.at(19) + 0.05 * (iBin - 19));
94 else
95 pbins.push_back(pbins.at(59) + 0.3 * (iBin - 59));
96 }
97
98 return pbins;
99 }
100 };
102} // namespace Belle2
Base class for calibration algorithms.
EResult
The result of calibration.
Class implementing the SVD dEdx calibration algorithm.
std::tuple< TH2F, TH2F, TH2F > DstarMassFit(std::shared_ptr< TTree > preselTree)
produce histograms for K/pi(/mu)
void setMinEvtsPerTree(const double &value)
set the upper edge of the dEdx binning for the payloads
int m_numPBins
the number of momentum bins for the payloads
void setNumDEdxBins(const int &value)
set the number of dEdx bins for the payloads
std::vector< double > CreatePBinningScheme()
build the binning scheme
void setDEdxCutoff(const double &value)
set the upper edge of the dEdx binning for the payloads
bool m_isMakePlots
produce plots for monitoring
int m_MinEvtsPerTree
number of events in TTree below which we don't try to fit
int m_numDEdxBins
the number of dEdx bins for the payloads
void setNumPBins(const int &value)
set the number of momentum bins for the payloads
TH2F GammaHistogram(std::shared_ptr< TTree > preselTree)
produce histograms for e
virtual EResult calibrate() override
run algorithm on data
TH2F LambdaMassFit(std::shared_ptr< TTree > preselTree)
produce histograms for protons
double m_dedxCutoff
the upper edge of the dEdx binning for the payloads
void setMonitoringPlots(bool value=false)
function to enable plotting
Abstract base class for different kinds of events.