ATLAS Offline Software
Loading...
Searching...
No Matches
TForwardElectronLikelihoodTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
11 details
12
13 The forward likelihood ID relies on the following list of topo-cluster
14 variables: el_secondLambda : longitudinal second moment el_lateral :
15 normalised lateral moment el_longitudinal : normalised longitudinal moment
16 el_centerLambda : shower depth
17 el_fracMax : fraction of cluster energy in the most energetic cell
18 el_secondR : transverse second moment
19 el_significance : cluster energy divided by cluster noise
20 el_secondDensity: second moment of signal density
21
22 Major revision when migrating to rel 21 in August 2018
23*/
24
25#ifndef TFORWARDELECTRONLIKELIHOODTOOL_H
26#define TFORWARDELECTRONLIKELIHOODTOOL_H
27
28#include <string> // for string
29#include <vector> // for vector
30
31// ROOT includes
32#include "TFile.h"
33#include "TH1.h"
34
35// Include the return objects and the base class
36//#include "PATCore/TAccept.h"
37//#include "PATCore/TResult.h"
38#include "PATCore/AcceptData.h"
39#include "PATCore/AcceptInfo.h"
40//#include "PATCore/TCalculatorToolBase.h"
41//#include "PATCore/TSelectorToolBase.h"
43#include "AsgTools/IAsgTool.h"
44#include "SafeTH1.h"
45#include "memory"
46
47namespace {
48const unsigned int IP_FBINS = 1;
49}
50
51namespace LikeEnumForward {
59
61{
62 double likelihood;
63 double eta;
64 double eT;
65 double ip;
66};
67
69{
70 double eta;
71 double eT;
73 double lateral;
75 double fracMax;
77 double secondR;
80 double ip;
81};
82
83} // namespace LikeEnumForward
84
85namespace Root {
87{
88
89public:
92 const char* name = "TForwardElectronLikelihoodTool");
93
94public:
96 StatusCode initialize();
97
99 const asg::AcceptInfo& getAcceptInfo() const { return m_acceptInfo; }
100
103 asg::AcceptData accept(double likelihood,
104 double eta,
105 double eT,
106 double ip) const;
107
110
111 double calculate(LikeEnumForward::LHCalcVars_t& vars_struct) const;
112 double calculate(double eta,
113 double eT,
114 double secondLambda,
115 double lateral,
116 double longitudinal,
117 double centerLambda,
118 double fracMax,
119 double secondR,
120 double significance,
121 double secondDensity,
122 double ip) const;
123
125 inline void setPDFFileName(const std::string& val) { m_pdfFileName = val; }
126
128 inline void setVariableNames(const std::string& val)
129 {
130 m_variableNames = val;
132 }
133
135 int loadVarHistograms(const std::string& vstr, unsigned int varIndex);
136
138 inline void setBinning(const std::string& val) { m_ipBinning = val; }
139
140 unsigned int getBitmask(void) const { return m_variableBitMask; }
141 inline void setBitmask(unsigned int val) { m_variableBitMask = val; };
142
143 // Private methods
144private:
145 // For every input "varVector", make sure elements of vector are
146 // in the same order as prescribed in fVariables
147 // Internal methods to calculate the LH discriminant from a set of variables
148 double evaluateLikelihood(const std::vector<double>& varVector,
149 double et,
150 double eta,
151 double ip = 0) const;
152 double evaluateLikelihood(const std::vector<float>& varVector,
153 double et,
154 double eta,
155 double ip = 0) const;
156
160 unsigned int getLikelihoodBitmask(const std::string& vars) const;
161
162 // PHILIP interpolation is unused
163 /* double InterpolateCuts(const std::vector<double>& cuts,const
164 * std::vector<double>& cuts_4gev,double et,double eta) const; */
165 /* double InterpolatePdfs(unsigned int s_or_b,unsigned int ipbin,double
166 * et,double eta,int bin,unsigned int var) const; */
167
169 // double TransformLikelihoodOutput(double ps, double pb) const;
171 // unsigned int getLikelihoodEtaBin(double eta) const;
173 // unsigned int getLikelihoodEtBin(double et) const;
174 // Pile-up binning
175 // unsigned int getIpBin(double ip) const;
176 // get the bin names as is in the input file
177 // void getBinName(char *buffer, int etbin, int etabin) const;
178
179public:
183 std::vector<double> m_cutLikelihood;
193 std::string m_variableNames;
195 std::string m_pdfFileName;
196
197 // Private methods
198private:
201 double TransformLikelihoodOutput(double ps, double pb) const;
202
204 static unsigned int getLikelihoodEtaBin(double eta) ;
205
207 static unsigned int getLikelihoodEtHistBin(double eT) ;
208
209 // Private member variables
210private:
213 unsigned int m_variableBitMask;
214
217
219 std::string m_ipBinning;
220
222 TFile* m_pdfFile{};
223
230
231 static const double fIpBounds[IP_FBINS + 1];
232 static const unsigned int s_fnEtBinsHist =
233 4; // number of hists stored for nominal LH
234 static const unsigned int s_fnDiscEtBins =
235 4; // number of discs stored for original LH
236 static const unsigned int s_fnEtaBins = 10;
237 static const unsigned int s_fnVariables = 8;
238 static const std::string fVariables[s_fnVariables];
239 // 5D array of unique_ptr to SafeTH1 // [sig(0)/bkg(1)][ip][et][eta][variable]
240 std::unique_ptr<EGSelectors::SafeTH1> m_fPDFbins[2][IP_FBINS][s_fnEtBinsHist][s_fnEtaBins][s_fnVariables];
241
242 static unsigned int getIpBin(double ip) ;
243 static std::string getBinName(int etbin,
244 int etabin,
245 int ipbin,
246 const std::string& iptype) ;
247};
248
249} // namespace Root
250
251#endif
Scalar eta() const
pseudorapidity method
static const std::string fVariables[s_fnVariables]
unsigned int m_variableBitMask
The bitmask corresponding to the variables in the likelihood.
std::vector< double > m_cutLikelihood
cut on likelihood output
int m_cutPosition_kinematicEta
The position of the kinematic cuts bit in the AcceptInfo return object, separate for eta/Et.
static std::string getBinName(int etbin, int etabin, int ipbin, const std::string &iptype)
const asg::AcceptInfo & getAcceptInfo() const
accesss to the accept info object
double TransformLikelihoodOutput(double ps, double pb) const
Apply a transform to zoom into the LH output peaks.
std::vector< double > m_cutLikelihoodPileupCorrectionB
pileup constant factor for cut on likelihood output
void setVariableNames(const std::string &val)
Define the variable names.
void setBinning(const std::string &val)
Define the binning.
unsigned int getLikelihoodBitmask(const std::string &vars) const
Mask out the variables ,out of all possible ones, that are not employed in the current configuration ...
bool m_doPileupCorrection
Apply a transform to zoom into the LH output peaks.
int m_cutPosition_LH
The position of the likelihood cut bit in the AcceptInfo return object.
int loadVarHistograms(const std::string &vstr, unsigned int varIndex)
Load the variable histograms from the pdf file.
static unsigned int getLikelihoodEtaBin(double eta)
Eta binning for pdfs and discriminant cuts.
TFile * m_pdfFile
Pointer to the opened TFile that holds the PDFs.
void setPDFFileName(const std::string &val)
Add an input file that holds the PDFs.
std::string m_variableNames
variables to use in the LH
asg::AcceptData accept() const
Return dummy accept with only info.
std::unique_ptr< EGSelectors::SafeTH1 > m_fPDFbins[2][IP_FBINS][s_fnEtBinsHist][s_fnEtaBins][s_fnVariables]
double calculate(LikeEnumForward::LHCalcVars_t &vars_struct) const
TForwardElectronLikelihoodTool(const char *name="TForwardElectronLikelihoodTool")
Standard constructor.
double evaluateLikelihood(const std::vector< double > &varVector, double et, double eta, double ip=0) const
static unsigned int getLikelihoodEtHistBin(double eT)
Coarse Et binning. Used for the likelihood and discriminant pdfs.
std::vector< double > m_cutLikelihoodPileupCorrectionA
the cut on the PU discriminant is adjusted as a function of nVtx cut + nVtx*cutA + cutB this is diffe...
Class mimicking the AthMessaging class from the offline software.
Extra patterns decribing particle interation process.