ATLAS Offline Software
Loading...
Searching...
No Matches
AsgElectronSelectorTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Dear emacs, this is -*-c++-*-
6
7#ifndef __ASGELECTRONSELECTORTOOL__
8#define __ASGELECTRONSELECTORTOOL__
9
10// This include is needed at the top before any includes regarding Eigen
11// since it includes Eigen in a specific way which causes compilation errors
12// if not included before Eigen
14
15// Atlas includes
16#include "AsgTools/AsgTool.h"
19#include <Eigen/Dense>
20
21class EventContext;
22
24
26 virtual public IAsgElectronLikelihoodTool
27{
29
30public:
32 AsgElectronSelectorTool( const std::string& myname );
33
34
37public:
39 virtual StatusCode initialize() override;
40
41 // Main methods for IAsgSelectorTool interface
46 virtual const asg::AcceptInfo& getAcceptInfo() const override;
47
49 asg::AcceptData accept( const xAOD::IParticle* part ) const override;
50 asg::AcceptData accept( const EventContext& ctx, const xAOD::IParticle* part ) const override;
51
53 asg::AcceptData accept( const EventContext& ctx, const xAOD::Electron* eg ) const override {
54 return accept (ctx, eg, -99); // mu = -99 as input will force accept to grab the pileup variable from the xAOD object
55 }
56
58 asg::AcceptData accept( const EventContext& ctx, const xAOD::Egamma* eg ) const override {
59 return accept (ctx, eg, -99); // mu = -99 as input will force accept to grab the pileup variable from the xAOD object
60 }
61
63 asg::AcceptData accept( const EventContext& ctx, const xAOD::Electron* eg, double mu ) const override;
64
66 asg::AcceptData accept( const EventContext& ctx, const xAOD::Egamma* eg, double mu ) const override;
67
68 // Main methods for IAsgCalculatorTool interface
69public:
71 double calculate( const xAOD::IParticle* part ) const;
72 double calculate( const EventContext &ctx, const xAOD::IParticle* part ) const override;
73
75 double calculate( const EventContext &ctx, const xAOD::Electron* eg ) const override {
76 return calculate (ctx, eg, -99); // mu = -99 as input will force accept to grab the pileup variable from the xAOD object
77 }
78
80 double calculate( const EventContext &ctx, const xAOD::Egamma* eg ) const override {
81 return calculate (ctx, eg, -99); // mu = -99 as input will force accept to grab the pileup variable from the xAOD object
82 }
83
85 double calculate( const EventContext &ctx, const xAOD::Electron* eg, double mu ) const override;
86
88 double calculate( const EventContext &ctx, const xAOD::Egamma* eg, double mu ) const override;
89
91 double getDiscriminant(std::vector<float>& mvaOutputs, const xAOD::Electron* egu ) const;
92
94 std::vector<float> calculateMultipleOutputs( const EventContext &ctx, const xAOD::Electron *eg, double mu = -99) const override;
95
96 virtual std::string getOperatingPointName() const override;
97
98 // Private methods
99private:
102
104 bool isForwardElectron( const xAOD::Egamma* eg, const float eta ) const;
105
107 double transformMLOutput( float score ) const;
108
110 double combineOutputs(const std::vector<float>& mvaScores, double eta) const;
111 static double combineOutputsCF(const std::vector<float>& mvaScores) ;
112
114 static unsigned int getDiscEtaBin( double eta ) ;
115
117 static unsigned int getDiscEtBin( double et ) ;
118
119 // NOTE that this will only perform the cut interpolation up to ~45 GeV, so
120 // no smoothing is done above this for the high ET MVA binning yet
122 static double interpolateCuts( const std::vector<double>& cuts, double et, double eta ) ;
123
124
125
126 // Private member variables
127private:
128
130 std::string m_workingPoint;
131
133 std::string m_configFile;
134
136 std::unique_ptr<const ElectronDNNCalculator> m_mvaTool;
137
139 std::string m_modelFileName;
140
143
145 std::vector<std::string> m_variables;
146
148 std::vector<int> m_enum_variables;
149
152
154
162 std::vector<double> m_fractions;
163
165 std::vector<int> m_cutAmbiguity;
167 std::vector<int> m_cutBL;
169 std::vector<int> m_cutPi;
171 std::vector<int> m_cutSCT;
175 std::vector<double> m_cutSelector;
176 std::vector<double> m_cutSelectorCF;
177
178
191
193 std::vector<float> m_defaultVector;
194
196 static const unsigned int s_fnDiscEtBins = 10;
198 static const unsigned int s_fnDiscEtaBins = 10;
199
200
201}; // End: class definition
202
203#endif
Scalar eta() const
pseudorapidity method
#define ASG_TOOL_CLASS2(CLASSNAME, INT1, INT2)
std::vector< std::string > m_variables
Variables used in the MVA Tool.
virtual std::string getOperatingPointName() const override
Get the name of the current operating point.
bool m_doSmoothBinInterpolation
do smooth interpolation between bins
static const unsigned int s_fnDiscEtBins
number of discrimintants vs Et
std::vector< double > m_fractions
Fractions to combine the output nodes of a multiclass model into one discriminant.
asg::AcceptData accept(const EventContext &ctx, const xAOD::Electron *eg) const override
The main accept method: the actual cuts are applied here.
virtual const asg::AcceptInfo & getAcceptInfo() const override
Method to get the plain AcceptInfo.
std::string m_workingPoint
Working Point.
std::vector< int > m_enum_variables
Enum version of used variables.
static unsigned int getDiscEtaBin(double eta)
Gets the Discriminant Eta bin [0,s_fnDiscEtaBins-1] given the eta.
bool m_CFReject
Run CF rejection or not.
double getDiscriminant(std::vector< float > &mvaOutputs, const xAOD::Electron *egu) const
Computes discrimiant value from mva output based on whether multiclass is true or false.
int m_cutPosition_kinematic
The position of the kinematic cut bit in the AcceptInfo return object.
std::unique_ptr< const ElectronDNNCalculator > m_mvaTool
Pointer to the class that calculates the MVA score.
bool m_multiClass
Multiclass model or not.
static unsigned int getDiscEtBin(double et)
Gets the Descriminant Et bin the et (MeV) [0,s_fnDiscEtBins-1].
virtual StatusCode initialize() override
Gaudi Service Interface method implementations.
std::vector< int > m_cutBL
cut min on b-layer hits
std::string m_modelFileName
The input file name that holds the model.
double calculate(const EventContext &ctx, const xAOD::Egamma *eg) const override
The main result method: the actual mva score is calculated here.
bool m_cfSignal
Use the CF output node in the numerator or the denominator.
std::vector< int > m_cutPi
cut min on pixel hits
bool m_skipDeltaPoverP
Flag for skip the use of deltaPoverP in dnn calculation (like at HLT)
double combineOutputs(const std::vector< float > &mvaScores, double eta) const
Combines the six output nodes of a multiclass model into one discriminant.
virtual ~AsgElectronSelectorTool()
Standard destructor.
static double interpolateCuts(const std::vector< double > &cuts, double et, double eta)
Interpolates cut values along pt.
static const unsigned int s_fnDiscEtaBins
number of discriminants vs |eta|
asg::AcceptInfo m_acceptMVA
Accept info.
double calculate(const EventContext &ctx, const xAOD::Electron *eg) const override
The main result method: the actual mva score is calculated here.
std::vector< float > calculateMultipleOutputs(const EventContext &ctx, const xAOD::Electron *eg, double mu=-99) const override
The result method for multiple outputs: can return multiple outputs of the MVA.
int m_cutPosition_NPixel
The position of the NPixel cut bit in the AcceptInfo return object.
std::vector< double > m_cutSelectorCF
int m_cutPosition_NSilicon
The position of the NSilicon cut bit in the AcceptInfo return object.
std::vector< int > m_cutSCT
cut min on precision hits
std::vector< int > m_cutAmbiguity
do cut on ambiguity bit
asg::AcceptData accept(const EventContext &ctx, const xAOD::Egamma *eg) const override
The main accept method: the actual cuts are applied here.
std::vector< double > m_cutSelector
cut on mva output
std::string m_quantileFileName
The input file name that holds the QuantileTransformer.
bool isForwardElectron(const xAOD::Egamma *eg, const float eta) const
check for FwdElectron
AsgElectronSelectorTool(const std::string &myname)
Standard constructor.
int m_cutPosition_NBlayer
The position of the NBlayer cut bit in the AcceptInfo return object.
int m_cutPosition_MVA
The position of the MVA cut bit in the AcceptInfo return object.
double transformMLOutput(float score) const
Applies a logit transformation to the score returned by the underlying MVA tool.
std::vector< float > m_defaultVector
Default vector to return if calculation fails.
double calculate(const xAOD::IParticle *part) const
The main result method: the actual mva score is calculated here.
asg::AcceptData accept(const xAOD::IParticle *part) const override
The main accept method: using the generic interface.
static double combineOutputsCF(const std::vector< float > &mvaScores)
int m_cutPosition_ambiguity
The position of the ambiguity cut bit in the AcceptInfo return object.
std::string m_configFile
The input config file.
Used by AsgElectronSelectorTool to calculate the score of a python trained DNN using lwtnn as interfa...
Interface to tool to select electrons.
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
Class providing the definition of the 4-vector interface.
Egamma_v1 Egamma
Definition of the current "egamma version".
Definition Egamma.h:17
Electron_v1 Electron
Definition of the current "egamma version".
Extra patterns decribing particle interation process.