ATLAS Offline Software
Loading...
Searching...
No Matches
TForwardElectronLikelihoodTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include "TROOT.h"
7#include "TSystem.h"
8#include <cmath>
9#include <format>
10
11const double Root::TForwardElectronLikelihoodTool::fIpBounds[IP_FBINS + 1] = {
12 0.,
13 500.
14};
15// These are the variables availalble in the likelihood.
16const std::string
18 "el_secondLambda", "el_lateral", "el_longitudinal", "el_centerLambda",
19 "el_fracMax", "el_secondR", "el_significance", "el_secondDensity"
20 };
21
34
35StatusCode
37{
38 ATH_MSG_DEBUG("TForwardElectronLikelihoodTool initialize.");
39 // use an int as a StatusCode
40 StatusCode sc(StatusCode::SUCCESS);
41
42 // Check that all needed variables are setup
43 if (m_pdfFileName.empty()) {
45 "You need to specify the input PDF file name before you call "
46 "initialize() with setPDFFileName('your/file/name.root') ");
47 sc = StatusCode::FAILURE;
48 }
49
50 unsigned int number_of_expected_bin = s_fnDiscEtBins * s_fnEtaBins;
51 if (m_cutLikelihood.size() != number_of_expected_bin) {
52 ATH_MSG_ERROR("Configuration issue : CutLikelihood expected size "
53 << number_of_expected_bin << " input size "
54 << m_cutLikelihood.size());
56 "Could NOT initialize! Please fix the errors mentioned above...");
57 sc = StatusCode::FAILURE;
58 return sc;
59 }
60
61 // --------------------------------------------------------------------------
62 // Register the cuts and check that the registration worked:
63
64 // Cut position for the kineatic pre-selection, separately for eta/Et
66 m_acceptInfo.addCut("kinematicEta", "pass kinematic eta");
68 sc = StatusCode::FAILURE;
70 m_acceptInfo.addCut("kinematicEt", "pass kinematic Et ");
72 sc = StatusCode::FAILURE;
73
74 // Cut position for the likelihood selection - DO NOT CHANGE ORDER!
75 m_cutPosition_LH = m_acceptInfo.addCut("passLH", "pass Likelihood");
76 if (m_cutPosition_LH < 0) {
77 sc = StatusCode::FAILURE;
78 }
79
80 // Check that we got everything OK
81 if (sc == StatusCode::FAILURE) {
83 "! Something went wrong with the setup of the return objects...");
84 return sc;
85 }
86
87 // Get the correct bit mask for the current likelihood operating point
89
90 //----------File/Histo operation------------------------------------
91 // Load the ROOT file containing the PDFs
92 TString tmpString(m_pdfFileName);
93 gSystem->ExpandPathName(tmpString);
94 std::string fname(tmpString.Data());
95 m_pdfFile = TFile::Open(fname.c_str(), "READ");
96 // Check that we could load the ROOT file
97 if (!m_pdfFile) {
98 ATH_MSG_ERROR(" No ROOT file found here: " << m_pdfFileName);
99 return StatusCode::FAILURE;
100 }
101
102 // Load the histograms
103 for (unsigned int varIndex = 0; varIndex < s_fnVariables; varIndex++) {
104 const std::string& vstr = fVariables[varIndex];
105 // Skip the loading of PDFs for variables we don't care about for this
106 // operating point. If the string is empty (which is true in the default
107 // 2012 case), load all of them.
108 if (m_variableNames.find(vstr) == std::string::npos &&
109 !m_variableNames.empty()) {
110 continue;
111 }
112 loadVarHistograms(vstr, varIndex);
113 }
114
115 // TFile close does not free the memory
116 m_pdfFile->Close();
117 // We need the destructor to be called
118 delete m_pdfFile;
119 //----------End File/Histo operation------------------------------------
120
121 ATH_MSG_DEBUG("Initialization complete for a LH tool with these specs:"
122 << "\n - PdfFileName : "
124 << "\n - Variable bitmask : "
126 << "\n - VariableNames : "
127 << m_variableNames);
128 return sc;
129}
130
131int
133 unsigned int varIndex)
134{
135 for (unsigned int s_or_b = 0; s_or_b < 2; ++s_or_b) {
136 for (unsigned int ip = 0; ip < IP_FBINS; ++ip) {
137 for (unsigned int et = 0; et < s_fnEtBinsHist; ++et) {
138 for (unsigned int eta = 0; eta < s_fnEtaBins; ++eta) {
139
140 std::string sig_bkg = (s_or_b == 0) ? "sig" : "bkg";
141 // For the fwd LH PDFs and optimisation have matching eta bin
142 // boundaries unsigned int eta_tmp = (eta > 0) ? eta-1 : eta ;
143 unsigned int eta_tmp = eta;
144 unsigned int et_tmp = et;
145 std::string binname = getBinName(et_tmp, eta_tmp, ip, m_ipBinning);
146
147 const std::string pdfdir = std::format("{}/{}", vstr, sig_bkg);
148 std::string pdf = std::format("{}_{}_smoothed_hist_from_KDE_{}",
149 vstr, sig_bkg, binname);
150
151 if (!m_pdfFile->GetListOfKeys()->Contains(vstr.c_str())) {
152 ATH_MSG_INFO("Warning: skipping variable "
153 << vstr << " because the folder does not exist.");
154 return 1;
155 }
156 if (!((TDirectory*)m_pdfFile->Get(vstr.c_str()))
157 ->GetListOfKeys()
158 ->Contains(sig_bkg.c_str())) {
159 ATH_MSG_INFO("Warning: skipping variable "
160 << vstr << " because the folder does not exist.");
161 return 1;
162 }
163 // Use the first Et bin given in the root file for all Et ranges below
164 if (et == 0 && !((TDirectory*)m_pdfFile->Get(pdfdir.c_str()))
165 ->GetListOfKeys()
166 ->Contains(pdf.c_str())) {
167 ATH_MSG_INFO("using lowest GeV bin in place of all below.");
168 binname = getBinName(et_tmp + 1, eta_tmp, ip, m_ipBinning);
169 pdf = std::format("{}_{}_smoothed_hist_from_KDE_{}",
170 vstr, sig_bkg, binname);
171 }
172 if (((TDirectory*)m_pdfFile->Get(pdfdir.c_str()))
173 ->GetListOfKeys()
174 ->Contains(pdf.c_str())) {
175 TH1F* hist =
176 (TH1F*)(((TDirectory*)m_pdfFile->Get(pdfdir.c_str()))->Get(pdf.c_str()));
177 m_fPDFbins[s_or_b][ip][et][eta][varIndex] =
178 std::make_unique<EGSelectors::SafeTH1>(hist);
179 } else {
180 ATH_MSG_INFO("Warning: Object " << pdf << " does not exist.");
181 ATH_MSG_INFO("Skipping all other histograms with this variable.");
182 return 1;
183 }
184 }
185 }
186 }
187 }
188 return 1;
189}
190
193 double eta,
194 double eT,
195 double ip) const
196{
198
199 vars.likelihood = likelihood;
200 vars.eta = eta;
201 vars.eT = eT;
202 vars.ip = ip;
203
204 return accept(vars);
205}
206
207// This method calculates if the current electron passes the requested
208// likelihood cut
211 LikeEnumForward::LHAcceptVars_t& vars_struct) const
212{
213 // Setup return accept with AcceptInfo
214 asg::AcceptData acceptData(&m_acceptInfo);
215
216 // Set up the individual cuts
217 bool passKineEta(true);
218 bool passKineEt(true);
219 bool passLH(true);
220
221 if (std::fabs(vars_struct.eta) < 2.5) {
222 ATH_MSG_DEBUG("This forward electron has"
223 << vars_struct.eta
224 << ", which is fabs(eta)<2.5 Returning False.");
225 passKineEta = false;
226 }
227 // Return if the kinematic requirements are not fulfilled
228 acceptData.setCutResult(m_cutPosition_kinematicEta, passKineEta);
229 if (!passKineEta) {
230 return acceptData;
231 }
232
233 unsigned int etbin = getLikelihoodEtHistBin(vars_struct.eT);
234 unsigned int etabin = getLikelihoodEtaBin(vars_struct.eta);
235
236 if (etbin >= s_fnDiscEtBins) {
237 ATH_MSG_WARNING("Cannot evaluate likelihood for Et "
238 << vars_struct.eT << ". Returning false..");
239 passKineEt = false;
240 }
241
242 // Return if the kinematic requirements are not fulfilled
243 acceptData.setCutResult(m_cutPosition_kinematicEt, passKineEt);
244 if (!passKineEt) {
245 return acceptData;
246 }
247
248 double cutDiscriminant = -999.;
249 unsigned int ibin_combined = etbin * s_fnEtaBins + etabin;
250
251 if (!m_cutLikelihood.empty()) {
252
253 cutDiscriminant = m_cutLikelihood[ibin_combined];
254 // If doPileupCorrection, then correct the discriminant itself instead of
255 // the cut value
257 cutDiscriminant +=
258 vars_struct.ip * m_cutLikelihoodPileupCorrectionA[ibin_combined] +
260 }
261 }
262
263 // Determine if the calculated likelihood value passes the cut
264 ATH_MSG_DEBUG("Likelihood macro: Discriminant: ");
265 if (vars_struct.likelihood < cutDiscriminant) {
266 ATH_MSG_DEBUG("Likelihood macro: Disciminant Cut Failed.");
267 passLH = false;
268 }
269 acceptData.setCutResult(m_cutPosition_LH, passLH);
270 return acceptData;
271}
272
273double
275 double eT,
276 double secondLambda,
277 double lateral,
278 double longitudinal,
279 double centerLambda,
280 double fracMax,
281 double secondR,
282 double significance,
283 double secondDensity,
284 double ip) const
285{
287 vars.eta = eta;
288 vars.eT = eT;
289 vars.secondLambda = secondLambda;
290 vars.lateral = lateral;
291 vars.longitudinal = longitudinal;
292 vars.centerLambda = centerLambda;
293 vars.fracMax = fracMax;
294 vars.secondR = secondR;
295 vars.significance = significance;
296 vars.secondDensity = secondDensity;
297 vars.ip = ip;
298 return calculate(vars);
299}
300
301// The main public method to actually calculate the likelihood value
302double
304 LikeEnumForward::LHCalcVars_t& vars_struct) const
305{
306 // Reset the results to defaul values
307 double result = -999;
308 double arr[] = { vars_struct.secondLambda, vars_struct.lateral,
309 vars_struct.longitudinal, vars_struct.centerLambda,
310 vars_struct.fracMax, vars_struct.secondR,
311 vars_struct.significance, vars_struct.secondDensity };
312 std::vector<double> vec(arr, arr + sizeof(arr) / sizeof(double));
313
314 // Calculate the actual likelihood value and fill the return object
316 vec, vars_struct.eT, vars_struct.eta, vars_struct.ip);
317 return result;
318}
319
320double
322 const std::vector<float>& varVector,
323 double et,
324 double eta,
325 double ip) const
326{
327 std::vector<double> vec;
328 for (unsigned int var = 0; var < s_fnVariables; var++) {
329 vec.push_back(varVector[var]);
330 }
331 return evaluateLikelihood(vec, et, eta, ip);
332}
333
334double
336 const std::vector<double>& varVector,
337 double et,
338 double eta,
339 double ip) const
340{
341 // const double GeV = 1000;
342 unsigned int etbin = getLikelihoodEtHistBin(et);
343 unsigned int etabin = getLikelihoodEtaBin(eta);
344 unsigned int ipbin = getIpBin(ip);
345
346 if (etbin >= s_fnEtBinsHist) {
347 ATH_MSG_WARNING("skipping etbin " << etbin << ", et " << et);
348 return -999.;
349 }
350 if (etabin >= s_fnEtaBins) {
351 ATH_MSG_WARNING("skipping etabin " << etabin << ", eta " << eta);
352 return -999.;
353 }
354
355 if (varVector.size() != s_fnVariables) {
356 ATH_MSG_WARNING("Error! Variable vector size mismatch! Check your vector!");
357 return -999.;
358 }
359
360 double SigmaS = 1.;
361 double SigmaB = 1.;
362
363 for (unsigned int var = 0; var < s_fnVariables; var++) {
364
365 std::string varstr = fVariables[var];
366
367 // Skip variables that are masked off (not used) in the likelihood
368 if (!(m_variableBitMask & (0x1 << var))) {
369 continue;
370 }
371
372 for (unsigned int s_or_b = 0; s_or_b < 2; s_or_b++) {
373
374 int bin =
375 m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->FindBin(varVector[var]);
376 double prob = 0;
377
378 double integral =
379 double(m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->Integral());
380 if (integral == 0) { // currently, the crack always has integral == 0
382 "Error! PDF integral == 0!"); // changed it to debug message since we
383 // intend to run the tool in the
384 // derivations
385 return -1.35;
386 }
387
388 prob = double(m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->GetBinContent(
389 bin)) /
390 integral;
391
392 if (s_or_b == 0)
393 SigmaS *= prob;
394 else if (s_or_b == 1)
395 SigmaB *= prob;
396 }
397 }
398
399 return TransformLikelihoodOutput(SigmaS, SigmaB);
400}
401
402// --------------------------------------------
403double
405 double pb) const
406{
407 // returns transformed or non-transformed output
408 // (Taken mostly from TMVA likelihood code)
409 double fEpsilon = 1e-99;
410 // If both signal and bkg are 0, we want it to fail.
411 if (ps < fEpsilon)
412 ps = 0;
413 if (pb < fEpsilon)
414 pb = fEpsilon;
415 double disc = ps / double(ps + pb);
416
417 if (disc >= 1.0)
418 disc = 1. - 1.e-15;
419 else if (disc <= 0.0)
420 disc = fEpsilon;
421
422 double tau = 15.0;
423 disc = -std::log(1.0 / disc - 1.0) * (1. / double(tau));
424
425 ATH_MSG_DEBUG("disc is " << disc);
426 return disc;
427}
428
429// Gets the IP bin
430unsigned int
432{
433 for (unsigned int ipBin = 0; ipBin < IP_FBINS; ++ipBin) {
434 if (ip < fIpBounds[ipBin + 1])
435 return ipBin;
436 }
437 return 0;
438}
439
440// Gets the Eta bin given the eta . Binning uses uper bound
441unsigned int
443{
444 const unsigned int nEtaBins = s_fnEtaBins;
445 const double etaBins[nEtaBins] = { 2.6, 2.7, 2.8, 2.9, 3.0,
446 3.1, 3.16, 3.35, 3.6, 4.9 };
447 for (unsigned int etaBin = 0; etaBin < nEtaBins; ++etaBin) {
448 if (std::fabs(eta) < etaBins[etaBin])
449 return etaBin;
450 }
451 return (nEtaBins - 1);
452}
453// Gets the histogram Et bin given the et (MeV). Binning uses upper bound
454unsigned int
456{
457 const double GeV = 1000;
458 const unsigned int nEtBins = s_fnDiscEtBins;
459 const double eTBins[nEtBins] = {
460 30 * GeV, 40 * GeV, 50 * GeV, 7000 * GeV
461 }; // removed 20 GeV bin since we're upper bound oriented
462 for (unsigned int eTBin = 0; eTBin < nEtBins; ++eTBin) {
463 if (eT < eTBins[eTBin]) {
464 return eTBin;
465 }
466 }
467 return (nEtBins - 1); // Return the last bin if > the last bin.
468}
469
470// Gets the bin name. Given the HISTOGRAM bin naming used in the input file
471// which seems to be lower bound oriented ...
472std::string
474 int etbin,
475 int etabin,
476 int ipbin,
477 const std::string& iptype)
478{
479 double eta_bounds[s_fnEtaBins] = { 2.5, 2.6, 2.7, 2.8, 2.9,
480 3.0, 3.1, 3.16, 3.35, 3.6 };
481 int et_bounds[s_fnEtBinsHist] = { 20, 30, 40, 50 };
482 if (!iptype.empty()) {
483 return std::format("{}{}et{:02}eta{:.2f}",
484 iptype, int(fIpBounds[ipbin]),
485 et_bounds[etbin], eta_bounds[etabin]);
486 } else {
487 return std::format("et{}eta{:.2f}",
488 et_bounds[etbin], eta_bounds[etabin]);
489 }
490}
491
492unsigned int
494 const std::string& vars) const
495{
496 unsigned int mask = 0x0;
497 ATH_MSG_DEBUG("Variables to be used: ");
498 for (unsigned int var = 0; var < s_fnVariables; var++) {
499 if (vars.find(fVariables[var]) != std::string::npos) {
501 mask = mask | 0x1 << var;
502 }
503 }
504 ATH_MSG_DEBUG("mask: " << mask);
505 return mask;
506}
Scalar eta() const
pseudorapidity method
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::vector< size_t > vec
static Double_t sc
asg::AcceptData accept(LikeEnumForward::LHAcceptVars_t &vars_struct) const
The main accept method: the actual cuts are applied here.
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)
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
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.
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...
void setCutResult(const std::string &cutName, bool cutResult)
Set the result of a cut, based on the cut name (safer)
Definition AcceptData.h:134
AsgMessaging(const std::string &name)
Constructor with a name.
STL class.
double integral(TH1 *h)
Definition computils.cxx:59
STL namespace.
Extra patterns decribing particle interation process.