Belle II Software development
SVDdEdxValidationAlgorithm.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#include <TTree.h>
13#include <map>
14
15namespace Belle2 {
24 public:
30
36
40 void setMonitoringPlots(bool value = false) { m_isMakePlots = value; }
41
45 void setNumROCpoints(const unsigned int& value) { m_NumROCpoints = value; }
46
50 void setMinROCMomentum(const double& value) { m_MomLowROC = value; }
51
55 void setMaxROCMomentum(const double& value) { m_MomHighROC = value; }
56
60 void setNumEffBins(const unsigned int& value) { m_NumEffBins = value; }
61
65 void setMaxEffMomentum(const double& value) { m_MomHighEff = value; }
66
70 void setMinEvtsPerTree(const double& value) { m_MinEvtsPerTree = value; }
71
72 protected:
76 virtual EResult calibrate() override;
77
78 private:
83 void PlotEfficiencyPlots(const TString& PIDDetectorsName, TTree* SignalTree, TString SignalWeightName, TString SignalVarName,
84 TString SignalVarNameFull, TTree* FakeTree, TString FakeWeightName, TString FakeVarName, TString FakeVarNameFull,
85 TString PIDVarName, TString PIDCut, unsigned int nbins, double MomLow, double MomHigh);
86
90 void PlotROCCurve(TTree* SignalTree, TString SignalWeightName, TString SignalVarName, TString SignalVarNameFull, TTree* FakeTree,
91 TString FakeWeightName, TString FakeVarName, TString FakeVarNameFull, TString PIDVarName);
92
93 TTree* LambdaMassFit(std::shared_ptr<TTree> preselTree);
94 TTree* DstarMassFit(std::shared_ptr<TTree> preselTree);
96 int m_MinEvtsPerTree = 100;
97 unsigned int m_NumROCpoints = 250;
98 double m_MomLowROC = 0.;
99 double m_MomHighROC = 7.;
100 unsigned int m_NumEffBins = 30;
101 double m_MomHighEff = 2.5;
103 };
105} // namespace Belle2
Base class for calibration algorithms.
EResult
The result of calibration.
Class implementing the SVD dEdx calibration algorithm.
double m_MomLowROC
lower edge of the momentum interval considered for the ROC curve
void setNumROCpoints(const unsigned int &value)
set the number of points for ROC curve plotting
void setMinEvtsPerTree(const double &value)
set the upper edge of the dEdx binning for the payloads
void setMaxROCMomentum(const double &value)
set the lower edge of the momentum range for ROC curve plotting
void setNumEffBins(const unsigned int &value)
set the number of bins for the efficiency scan
void PlotROCCurve(TTree *SignalTree, TString SignalWeightName, TString SignalVarName, TString SignalVarNameFull, TTree *FakeTree, TString FakeWeightName, TString FakeVarName, TString FakeVarNameFull, TString PIDVarName)
a generic function to produce ROC curves
unsigned int m_NumEffBins
number of bins for the efficiency/fake rate plot
bool m_isMakePlots
produce plots for monitoring of the fit quality
void setMinROCMomentum(const double &value)
set the lower edge of the momentum range for ROC curve plotting
int m_MinEvtsPerTree
number of events in TTree below which we don't try to fit
TTree * LambdaMassFit(std::shared_ptr< TTree > preselTree)
produce histograms for protons
double m_MomHighROC
upper edge of the momentum interval considered for the ROC curve
void setMaxEffMomentum(const double &value)
set the upper edge of the momentum range for the efficiency scan
double m_MomHighEff
upper edge of the momentum interval for the efficiency/fake rate plot
unsigned int m_NumROCpoints
number of points for the ROC curve plotting
TTree * DstarMassFit(std::shared_ptr< TTree > preselTree)
produce histograms for K/pi(/mu)
virtual EResult calibrate() override
run algorithm on data
void PlotEfficiencyPlots(const TString &PIDDetectorsName, TTree *SignalTree, TString SignalWeightName, TString SignalVarName, TString SignalVarNameFull, TTree *FakeTree, TString FakeWeightName, TString FakeVarName, TString FakeVarNameFull, TString PIDVarName, TString PIDCut, unsigned int nbins, double MomLow, double MomHigh)
a generic function to produce efficiency plots
void setMonitoringPlots(bool value=false)
function to enable plotting
Abstract base class for different kinds of events.