ATLAS Offline Software
Loading...
Searching...
No Matches
Root::TElectronLikelihoodTool Class Reference

#include <TElectronLikelihoodTool.h>

Inheritance diagram for Root::TElectronLikelihoodTool:
Collaboration diagram for Root::TElectronLikelihoodTool:

Public Member Functions

 TElectronLikelihoodTool (const char *name="TElectronLikelihoodTool")
 Standard constructor.
StatusCode initialize ()
 Initialize this class.
const asg::AcceptInfogetAcceptInfo () const
 accesss to the accept info object
asg::AcceptData accept (LikeEnum::LHAcceptVars_t &vars_struct) const
 The main accept method: the actual cuts are applied here.
asg::AcceptData accept (double likelihood, double eta, double eT, int nSiHitsPlusDeadSensors, int nPixHitsPlusDeadSensors, bool passBLayerRequirement, uint8_t ambiguityBit, double d0, double deltaEta, double deltaphires, double wstot, double EoverP, double ip) const
asg::AcceptData accept () const
 Return dummy accept with only info.
double calculate (LikeEnum::LHCalcVars_t &vars_struct) const
double calculate (double eta, double eT, double f3, double rHad, double rHad1, double Reta, double w2, double f1, double eratio, double deltaEta, double d0, double d0sigma, double rphi, double deltaPoverP, double deltaphires, double TRT_PID, double ip) const
void setPDFFileName (const std::string &val)
 Add an input file that holds the PDFs.
void setVariableNames (const std::string &val)
 Define the variable names.
int loadVarHistograms (const std::string &vstr, TFile *pdfFile, unsigned int varIndex)
 Load the variable histograms from the pdf file.
void setBinning (const std::string &val)
 Define the binning.
unsigned int getBitmask (void) const
void setBitmask (unsigned int val)
void setLevel (MSG::Level lvl)
 Change the current logging level.
Functions providing the same interface as AthMessaging
bool msgLvl (const MSG::Level lvl) const
 Test the output level of the object.
MsgStream & msg () const
 The standard message stream.
MsgStream & msg (const MSG::Level lvl) const
 The standard message stream.

Public Attributes

std::vector< int > m_cutBL
 cut min on b-layer hits
std::vector< int > m_cutPi
 cut min on pixel hits
std::vector< int > m_cutSi
 cut min on precision hits
std::vector< double > m_cutA0
 cut max on track d0 bit
std::vector< double > m_cutDeltaEta
 do cut on delta eta bit
std::vector< double > m_cutDeltaPhiRes
 do cut on delta phi bit
std::vector< int > m_cutAmbiguity
 do cut on ambiguity bit
bool m_doRemoveF3AtHighEt
 do remove f3 variable from likelihood at high Et (>80 GeV)
bool m_doRemoveTRTPIDAtHighEt
 do remove TRTPID variable from likelihood at high Et (>80 GeV)
bool m_doSmoothBinInterpolation
 do smooth interpolation between bins
bool m_useOneExtraHighETLHBin
 use one extra bin for high ET LH
double m_highETBinThreshold
 ET threshold for using high ET cuts and bin.
std::vector< double > m_cutWstotAtHighET
 do cut on wstot above HighETBinThreshold bit
std::vector< double > m_cutEoverPAtHighET
 do cut on EoverP above HighETBinThreshold bit
bool m_doPileupTransform
 do pileup-dependent transform on discriminant value
bool m_doCentralityTransform
 do centrality-dependent transform on discriminant value
std::vector< double > m_cutLikelihood
 cut on likelihood output
std::vector< double > m_cutLikelihoodPileupCorrection
 pileup correction factor for cut on likelihood output
std::vector< double > m_cutLikelihood4GeV
 cut on likelihood output, 4 GeV bin
std::vector< double > m_cutLikelihoodPileupCorrection4GeV
 pileup correction factor for cut on likelihood output, 4 GeV bin
std::vector< double > m_discHardCutForPileupTransform
 reference disc for very hard cut; used by pileup transform
std::vector< double > m_discHardCutSlopeForPileupTransform
 reference slope on disc for very hard cut; used by pileup transform
std::vector< double > m_discHardCutQuadForPileupTransform
 reference quadratic apr on disc for very hard cut; used by centrality transform
std::vector< double > m_discLooseForPileupTransform
 reference disc for a pileup independent loose menu; used by pileup transform
std::vector< double > m_discHardCutForPileupTransform4GeV
 reference disc for very hard cut; used by pileup transform - 4-7 GeV
std::vector< double > m_discHardCutSlopeForPileupTransform4GeV
 reference slope on disc for very hard cut; used by pileup transform
std::vector< double > m_discHardCutQuadForPileupTransform4GeV
 reference quadratic par on disc for very hard cut; used by centrality transform - 4-7 GeV
std::vector< double > m_discLooseForPileupTransform4GeV
 reference disc for a pileup independent loose menu; used by pileup transform - 4-7 GeV
double m_discMaxForPileupTransform
 max discriminant for which pileup transform is to be used
double m_pileupMaxForPileupTransform
 max nvtx or mu to be used in pileup transform
std::string m_variableNames
 variables to use in the LH
std::string m_pdfFileName
 Name of the pdf file.

Private Member Functions

double evaluateLikelihood (const std::vector< double > &varVector, double et, double eta, double ip=0) const
double evaluateLikelihood (const std::vector< float > &varVector, double et, double eta, double ip=0) const
unsigned int getLikelihoodBitmask (const std::string &vars) const
double InterpolateCuts (const std::vector< double > &cuts, const std::vector< double > &cuts_4gev, double et, double eta) const
double InterpolatePdfs (unsigned int s_or_b, unsigned int ipbin, double et, double eta, int bin, unsigned int var) const
double TransformLikelihoodOutput (double ps, double pb, double ip, double et, double eta) const
 Apply a transform to zoom into the LH output peaks.
unsigned int getLikelihoodEtDiscBin (double eT, const bool isLHbinning) const
 Fine Et binning. Used for the likelihood discriminant cuts.
void initMessaging () const
 Initialize our message level and MessageSvc.

Static Private Member Functions

static unsigned int getLikelihoodEtaBin (double eta)
 Eta binning for pdfs and discriminant cuts.
static unsigned int getLikelihoodEtHistBin (double eT)
 Coarse Et binning. Used for the likelihood pdfs.
static unsigned int getIpBin (double ip)
static std::string getBinName (int etbin, int etabin, int ipbin, const std::string &iptype)

Private Attributes

std::string m_name
 tool name
asg::AcceptInfo m_acceptInfo
 Accept info.
unsigned int m_variableBitMask
 The bitmask corresponding to the variables in the likelihood.
std::string m_ipBinning
 Deprecated.
int m_cutPosition_kinematic
 The position of the kinematic cut bit in the AcceptInfo return object.
int m_cutPosition_NSilicon
 The position of the NSilicon cut bit in the AcceptInfo return object.
int m_cutPosition_NPixel
 The position of the NPixel cut bit in the AcceptInfo return object.
int m_cutPosition_NBlayer
 The position of the NBlayer cut bit in the AcceptInfo return object.
int m_cutPosition_ambiguity
 The position of the ambiguity cut bit in the AcceptInfo return object.
int m_cutPosition_LH
 The position of the likelihood cut bit in the AcceptInfo return object.
int m_cutPositionTrackA0
 The position of the d0 cut bit in the AcceptInfo return object.
int m_cutPositionTrackMatchEta
 The position of the deltaeta cut bit in the AcceptInfo return object.
int m_cutPositionTrackMatchPhiRes
 The position of the deltaphi cut bit in the AcceptInfo return object.
int m_cutPositionWstotAtHighET
 The position of the high ET wstot cut bit in the AcceptInfo return object.
int m_cutPositionEoverPAtHighET
 The position of the high ET EoverP cut bit in the AcceptInfo return object.
std::unique_ptr< EGSelectors::SafeTH1m_fPDFbins [2][IP_BINS][s_fnEtBinsHist][s_fnEtaBins][s_fnVariables]
std::string m_nm
 Message source name.
boost::thread_specific_ptr< MsgStream > m_msg_tls
 MsgStream instance (a std::cout like with print-out levels)
std::atomic< IMessageSvc * > m_imsg { nullptr }
 MessageSvc pointer.
std::atomic< MSG::Level > m_lvl { MSG::NIL }
 Current logging level.
std::atomic_flag m_initialized ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
 Messaging initialized (initMessaging)

Static Private Attributes

static constexpr double s_fIpBounds [IP_BINS+1] = { 0., 500. }
static constexpr unsigned int s_fnEtBinsHist = 7
static constexpr unsigned int s_fnDiscEtBins = 9
static constexpr unsigned int s_fnDiscEtBinsOneExtra = 10
static constexpr unsigned int s_fnEtaBins = 10
static constexpr unsigned int s_fnVariables = 13
static const std::string s_fVariables [s_fnVariables]

Detailed Description

Definition at line 78 of file TElectronLikelihoodTool.h.

Constructor & Destructor Documentation

◆ TElectronLikelihoodTool()

Root::TElectronLikelihoodTool::TElectronLikelihoodTool ( const char * name = "TElectronLikelihoodTool")

Standard constructor.

Author : Kurt Brendlinger kurb@.nosp@m.sas..nosp@m.upenn.nosp@m..edu Please see TElectronLikelihoodTool.h for usage.

Definition at line 27 of file TElectronLikelihoodTool.cxx.

28 : asg::AsgMessaging(std::string(name))
34 , m_doPileupTransform(false)
38 , m_variableNames("")
39 , m_pdfFileName("")
40 , m_name(name)
42 , m_ipBinning("")
54{
55}
int m_cutPosition_NPixel
The position of the NPixel cut bit in the AcceptInfo return object.
int m_cutPositionTrackMatchEta
The position of the deltaeta cut bit in the AcceptInfo return object.
bool m_useOneExtraHighETLHBin
use one extra bin for high ET LH
double m_highETBinThreshold
ET threshold for using high ET cuts and bin.
int m_cutPosition_ambiguity
The position of the ambiguity cut bit in the AcceptInfo return object.
int m_cutPositionTrackA0
The position of the d0 cut bit in the AcceptInfo return object.
int m_cutPosition_NBlayer
The position of the NBlayer cut bit in the AcceptInfo return object.
std::string m_variableNames
variables to use in the LH
bool m_doRemoveTRTPIDAtHighEt
do remove TRTPID variable from likelihood at high Et (>80 GeV)
int m_cutPosition_kinematic
The position of the kinematic cut bit in the AcceptInfo return object.
std::string m_pdfFileName
Name of the pdf file.
double m_discMaxForPileupTransform
max discriminant for which pileup transform is to be used
int m_cutPosition_LH
The position of the likelihood cut bit in the AcceptInfo return object.
int m_cutPositionWstotAtHighET
The position of the high ET wstot cut bit in the AcceptInfo return object.
bool m_doRemoveF3AtHighEt
do remove f3 variable from likelihood at high Et (>80 GeV)
int m_cutPositionEoverPAtHighET
The position of the high ET EoverP cut bit in the AcceptInfo return object.
unsigned int m_variableBitMask
The bitmask corresponding to the variables in the likelihood.
bool m_doCentralityTransform
do centrality-dependent transform on discriminant value
double m_pileupMaxForPileupTransform
max nvtx or mu to be used in pileup transform
int m_cutPosition_NSilicon
The position of the NSilicon cut bit in the AcceptInfo return object.
int m_cutPositionTrackMatchPhiRes
The position of the deltaphi cut bit in the AcceptInfo return object.
bool m_doPileupTransform
do pileup-dependent transform on discriminant value
bool m_doSmoothBinInterpolation
do smooth interpolation between bins

Member Function Documentation

◆ accept() [1/3]

asg::AcceptData Root::TElectronLikelihoodTool::accept ( ) const
inline

Return dummy accept with only info.

Definition at line 109 of file TElectronLikelihoodTool.h.

109{ return asg::AcceptData(&m_acceptInfo); }
asg::AcceptInfo m_acceptInfo
Accept info.

◆ accept() [2/3]

asg::AcceptData Root::TElectronLikelihoodTool::accept ( double likelihood,
double eta,
double eT,
int nSiHitsPlusDeadSensors,
int nPixHitsPlusDeadSensors,
bool passBLayerRequirement,
uint8_t ambiguityBit,
double d0,
double deltaEta,
double deltaphires,
double wstot,
double EoverP,
double ip ) const

Definition at line 400 of file TElectronLikelihoodTool.cxx.

413{
414 LikeEnum::LHAcceptVars_t vars{};
415 vars.likelihood = likelihood;
416 vars.eta = eta;
417 vars.eT = eT;
418 vars.nSiHitsPlusDeadSensors = nSiHitsPlusDeadSensors;
421 vars.ambiguityBit = ambiguityBit;
422 vars.d0 = d0;
423 vars.deltaEta = deltaEta;
424 vars.deltaphires = deltaphires;
425 vars.wstot = wstot;
426 vars.EoverP = EoverP;
427 vars.ip = ip;
428
429 return accept(vars);
430}
Scalar eta() const
pseudorapidity method
asg::AcceptData accept() const
Return dummy accept with only info.
bool passBLayerRequirement(const xAOD::TrackParticle &tp)
return true if effective number of BL hits + outliers is at least one
float eT(const U &p)
Accessor utility function for getting the value of Tranverse energy.
double deltaEta(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
Definition P4Helpers.h:66
setEt setPhi setE277 setWeta2 setEta1 setE2tsts1 wstot

◆ accept() [3/3]

asg::AcceptData Root::TElectronLikelihoodTool::accept ( LikeEnum::LHAcceptVars_t & vars_struct) const

The main accept method: the actual cuts are applied here.

Definition at line 435 of file TElectronLikelihoodTool.cxx.

437{
438 // Setup return accept with AcceptInfo
439 asg::AcceptData acceptData(&m_acceptInfo);
440
441 // Set up the individual cuts
442 bool passKine(true);
443 bool passNSilicon(true);
444 bool passNPixel(true);
445 bool passNBlayer(true);
446 bool passAmbiguity(true);
447 bool passLH(true);
448 bool passTrackA0(true);
449 bool passDeltaEta(true);
450 bool passDeltaPhiRes(true);
451 bool passWstotAtHighET(true);
452 bool passEoverPAtHighET(true);
453
454 if (std::abs(vars_struct.eta) > 2.47) {
455 ATH_MSG_DEBUG("This electron is std::abs(eta)>2.47 Returning False.");
456 passKine = false;
457 }
458
459 unsigned int etbinLH = getLikelihoodEtDiscBin(vars_struct.eT, true);
460 unsigned int etbinOther = getLikelihoodEtDiscBin(vars_struct.eT, false);
461 unsigned int etabin = getLikelihoodEtaBin(vars_struct.eta);
462
463 // sanity
464 if (etbinLH >= s_fnDiscEtBinsOneExtra) {
465 ATH_MSG_WARNING("Cannot evaluate likelihood for Et "
466 << vars_struct.eT << ". Returning false..");
467 passKine = false;
468 }
469 // sanity
470 if (etbinOther >= s_fnDiscEtBins) {
471 ATH_MSG_WARNING("Cannot evaluate likelihood for Et "
472 << vars_struct.eT << ". Returning false..");
473 passKine = false;
474 }
475
476 // Return if the kinematic requirements are not fulfilled
477 acceptData.setCutResult(m_cutPosition_kinematic, passKine);
478 if (!passKine) {
479 return acceptData;
480 }
481
482 // ambiguity bit
483 if (!m_cutAmbiguity.empty()) {
486 m_cutAmbiguity[etabin])) {
487 ATH_MSG_DEBUG("Likelihood macro: ambiguity Bit Failed.");
488 passAmbiguity = false;
489 }
490 }
491
492 // blayer cut
493 if (!m_cutBL.empty()) {
494 if (m_cutBL[etabin] == 1 && !vars_struct.passBLayerRequirement) {
495 ATH_MSG_DEBUG("Likelihood macro: Blayer cut failed.");
496 passNBlayer = false;
497 }
498 }
499 // pixel cut
500 if (!m_cutPi.empty()) {
501 if (vars_struct.nPixHitsPlusDeadSensors < m_cutPi[etabin]) {
502 ATH_MSG_DEBUG("Likelihood macro: Pixels Failed.");
503 passNPixel = false;
504 }
505 }
506 // silicon cut
507 if (!m_cutSi.empty()) {
508 if (vars_struct.nSiHitsPlusDeadSensors < m_cutSi[etabin]) {
509 ATH_MSG_DEBUG("Likelihood macro: Silicon Failed.");
510 passNSilicon = false;
511 }
512 }
513
514 double cutDiscriminant;
515 unsigned int ibin_combinedLH =
516 etbinLH * s_fnEtaBins + etabin; // Must change if number of eta bins
517 // changes!. Also starts from 7-10 GeV bin.
518 unsigned int ibin_combinedOther =
519 etbinOther * s_fnEtaBins +
520 etabin; // Must change if number of eta bins changes!. Also
521 // starts from 7-10 GeV bin.
522
523 if (!m_cutLikelihood.empty()) {
524 // To protect against a binning mismatch, which should never happen
525 if (ibin_combinedLH >= m_cutLikelihood.size()) {
526 ATH_MSG_ERROR("The desired eta/pt bin "
527 << ibin_combinedLH
528 << " is outside of the range specified by the input"
529 << m_cutLikelihood.size() << "This should never happen!");
530 return acceptData;
531 }
532
534 cutDiscriminant = InterpolateCuts(
535 m_cutLikelihood, m_cutLikelihood4GeV, vars_struct.eT, vars_struct.eta);
538 cutDiscriminant +=
541 vars_struct.eT,
542 vars_struct.eta);
543 } else {
544 if (vars_struct.eT > 7000. || m_cutLikelihood4GeV.empty()) {
545 cutDiscriminant = m_cutLikelihood[ibin_combinedLH];
546 // If doPileupTransform, then correct the discriminant itself instead of
547 // the cut value
549 cutDiscriminant +=
550 vars_struct.ip * m_cutLikelihoodPileupCorrection[ibin_combinedLH];
551 }
552 } else {
553 cutDiscriminant = m_cutLikelihood4GeV[etabin];
554 if (!m_doPileupTransform &&
556 cutDiscriminant +=
557 vars_struct.ip * m_cutLikelihoodPileupCorrection4GeV[etabin];
558 }
559 }
560
561 // Determine if the calculated likelihood value passes the cut
562 ATH_MSG_DEBUG("Likelihood macro: Discriminant: ");
563 if (vars_struct.likelihood < cutDiscriminant) {
564 ATH_MSG_DEBUG("Likelihood macro: Disciminant Cut Failed.");
565 passLH = false;
566 }
567 }
568
569 // d0 cut
570 if (!m_cutA0.empty()) {
571 if (std::abs(vars_struct.d0) > m_cutA0[ibin_combinedOther]) {
572 ATH_MSG_DEBUG("Likelihood macro: D0 Failed.");
573 passTrackA0 = false;
574 }
575 }
576
577 // deltaEta cut
578 if (!m_cutDeltaEta.empty()) {
579 if (std::abs(vars_struct.deltaEta) > m_cutDeltaEta[ibin_combinedOther]) {
580 ATH_MSG_DEBUG("Likelihood macro: deltaEta Failed.");
581 passDeltaEta = false;
582 }
583 }
584
585 // deltaPhiRes cut
586 if (!m_cutDeltaPhiRes.empty()) {
587 if (std::abs(vars_struct.deltaphires) >
588 m_cutDeltaPhiRes[ibin_combinedOther]) {
589 ATH_MSG_DEBUG("Likelihood macro: deltaphires Failed.");
590 passDeltaPhiRes = false;
591 }
592 }
593
594 // Only do this above HighETBinThreshold [in GeV]
595 if (vars_struct.eT > m_highETBinThreshold * 1000) {
596 // wstot cut
597 if (!m_cutWstotAtHighET.empty()) {
598 if (std::abs(vars_struct.wstot) > m_cutWstotAtHighET[etabin]) {
599 ATH_MSG_DEBUG("Likelihood macro: wstot Failed.");
600 passWstotAtHighET = false;
601 }
602 }
603
604 // EoverP cut
605 if (!m_cutEoverPAtHighET.empty()) {
606 if (std::abs(vars_struct.EoverP) > m_cutEoverPAtHighET[etabin]) {
607 ATH_MSG_DEBUG("Likelihood macro: EoverP Failed.");
608 passEoverPAtHighET = false;
609 }
610 }
611 }
612
613 // Set the individual cut bits in the return object
614 acceptData.setCutResult(m_cutPosition_NSilicon, passNSilicon);
615 acceptData.setCutResult(m_cutPosition_NPixel, passNPixel);
616 acceptData.setCutResult(m_cutPosition_NBlayer, passNBlayer);
617 acceptData.setCutResult(m_cutPosition_ambiguity, passAmbiguity);
618 acceptData.setCutResult(m_cutPosition_LH, passLH);
619 acceptData.setCutResult(m_cutPositionTrackA0, passTrackA0);
620 acceptData.setCutResult(m_cutPositionTrackMatchEta, passDeltaEta);
621 acceptData.setCutResult(m_cutPositionTrackMatchPhiRes, passDeltaPhiRes);
622 acceptData.setCutResult(m_cutPositionWstotAtHighET, passWstotAtHighET);
623 acceptData.setCutResult(m_cutPositionEoverPAtHighET, passEoverPAtHighET);
624 return acceptData;
625}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static constexpr unsigned int s_fnDiscEtBinsOneExtra
static constexpr unsigned int s_fnDiscEtBins
std::vector< int > m_cutSi
cut min on precision hits
double InterpolateCuts(const std::vector< double > &cuts, const std::vector< double > &cuts_4gev, double et, double eta) const
std::vector< int > m_cutBL
cut min on b-layer hits
std::vector< int > m_cutAmbiguity
do cut on ambiguity bit
std::vector< int > m_cutPi
cut min on pixel hits
std::vector< double > m_cutLikelihoodPileupCorrection4GeV
pileup correction factor for cut on likelihood output, 4 GeV bin
std::vector< double > m_cutA0
cut max on track d0 bit
std::vector< double > m_cutEoverPAtHighET
do cut on EoverP above HighETBinThreshold bit
std::vector< double > m_cutLikelihood4GeV
cut on likelihood output, 4 GeV bin
std::vector< double > m_cutWstotAtHighET
do cut on wstot above HighETBinThreshold bit
static unsigned int getLikelihoodEtaBin(double eta)
Eta binning for pdfs and discriminant cuts.
std::vector< double > m_cutDeltaPhiRes
do cut on delta phi bit
std::vector< double > m_cutLikelihoodPileupCorrection
pileup correction factor for cut on likelihood output
std::vector< double > m_cutLikelihood
cut on likelihood output
unsigned int getLikelihoodEtDiscBin(double eT, const bool isLHbinning) const
Fine Et binning. Used for the likelihood discriminant cuts.
static constexpr unsigned int s_fnEtaBins
std::vector< double > m_cutDeltaEta
do cut on delta eta bit
bool passAmbiguity(xAOD::AmbiguityTool::AmbiguityType type, const uint16_t criterion)
return true if the ambiguity type is one of several that are stored in a bitmask

◆ calculate() [1/2]

double Root::TElectronLikelihoodTool::calculate ( double eta,
double eT,
double f3,
double rHad,
double rHad1,
double Reta,
double w2,
double f1,
double eratio,
double deltaEta,
double d0,
double d0sigma,
double rphi,
double deltaPoverP,
double deltaphires,
double TRT_PID,
double ip ) const

Definition at line 628 of file TElectronLikelihoodTool.cxx.

645{
646
647 LikeEnum::LHCalcVars_t vars{};
648
649 vars.eta = eta;
650 vars.eT = eT;
651 vars.f3 = f3;
652 vars.rHad = rHad;
653 vars.rHad1 = rHad1;
654 vars.Reta = Reta;
655 vars.w2 = w2;
656 vars.f1 = f1;
657 vars.eratio = eratio;
658 vars.deltaEta = deltaEta;
659 vars.d0 = d0;
660 vars.d0sigma = d0sigma;
661 vars.rphi = rphi;
662 vars.deltaPoverP = deltaPoverP;
663 vars.deltaphires = deltaphires;
664 vars.TRT_PID = TRT_PID;
665 vars.ip = ip;
666
667 return calculate(vars);
668}
double calculate(LikeEnum::LHCalcVars_t &vars_struct) const
@ f3
fraction of energy reconstructed in 3rd sampling
Definition EgammaEnums.h:55
setCharge setNTRTHiThresholdHits eratio

◆ calculate() [2/2]

double Root::TElectronLikelihoodTool::calculate ( LikeEnum::LHCalcVars_t & vars_struct) const

Definition at line 672 of file TElectronLikelihoodTool.cxx.

674{
675 // Reset the results to defaul values
676 double result = -999;
677
678 unsigned int etabin = getLikelihoodEtaBin(vars_struct.eta);
679 double rhad_corr;
680 if (etabin == 3 || etabin == 4) {
681 rhad_corr = vars_struct.rHad;
682 } else {
683 rhad_corr = vars_struct.rHad1;
684 }
685 double d0significance = vars_struct.d0sigma == 0
686 ? 0.
687 : std::abs(vars_struct.d0) / vars_struct.d0sigma;
688
689 std::vector<double> vec = {
690 d0significance, vars_struct.eratio, vars_struct.deltaEta,
691 vars_struct.f1, vars_struct.f3, vars_struct.Reta,
692 rhad_corr, vars_struct.rphi, vars_struct.d0,
693 vars_struct.w2, vars_struct.deltaPoverP, vars_struct.deltaphires,
694 vars_struct.TRT_PID
695 };
696 // Calculate the actual likelihood value and fill the return object
698 vec, vars_struct.eT, vars_struct.eta, vars_struct.ip);
699
700 return result;
701}
std::vector< size_t > vec
double evaluateLikelihood(const std::vector< double > &varVector, double et, double eta, double ip=0) const

◆ evaluateLikelihood() [1/2]

double Root::TElectronLikelihoodTool::evaluateLikelihood ( const std::vector< double > & varVector,
double et,
double eta,
double ip = 0 ) const
private

Definition at line 718 of file TElectronLikelihoodTool.cxx.

723{
724
725 const double GeV = 1000;
726 unsigned int etbin = getLikelihoodEtHistBin(et); // hist binning
727 unsigned int etabin = getLikelihoodEtaBin(eta);
728 unsigned int ipbin = getIpBin(ip);
729
730 ATH_MSG_DEBUG("et: " << et << " eta: " << eta << " etbin: " << etbin
731 << " etabin: " << etabin);
732
733 if (etbin >= s_fnEtBinsHist) {
734 ATH_MSG_WARNING("skipping etbin " << etbin << ", et " << et);
735 return -999.;
736 }
737 if (etabin >= s_fnEtaBins) {
738 ATH_MSG_WARNING("skipping etabin " << etabin << ", eta " << eta);
739 return -999.;
740 }
741
742 if (varVector.size() != s_fnVariables) {
743 ATH_MSG_WARNING("Error! Variable vector size mismatch! Check your vector!");
744 }
745
746 double SigmaS = 1.;
747 double SigmaB = 1.;
748
749 // define some string constants
750 const std::string TRT_string = "TRT";
751 const std::string el_f3_string = "el_f3";
752 const std::string el_TRT_PID_string = "el_TRT_PID";
753
754 for (unsigned int var = 0; var < s_fnVariables; var++) {
755
756 const std::string& varstr = s_fVariables[var];
757
758 // Skip variables that are masked off (not used) in the likelihood
759 if (!(m_variableBitMask & (0x1 << var))) {
760 continue;
761 }
762 // Don't use TRT for outer eta bins (2.01,2.37)
763 if (((etabin == 8) || (etabin == 9)) &&
764 (varstr.find(TRT_string) != std::string::npos)) {
765 continue;
766 }
767 // Don't use f3 for outer eta bin (2.37)
768 if ((etabin == 9) && (varstr.find(el_f3_string) != std::string::npos)) {
769 continue;
770 }
771 // Don't use f3 for high et (>80 GeV)
772 if (m_doRemoveF3AtHighEt && (et > 80 * GeV) &&
773 (varstr.find(el_f3_string) != std::string::npos)) {
774 continue;
775 }
776 // Don't use TRTPID for high et (>80 GeV)
777 if (m_doRemoveTRTPIDAtHighEt && (et > 80 * GeV) &&
778 (varstr.find(el_TRT_PID_string) != std::string::npos)) {
779 continue;
780 }
781 for (unsigned int s_or_b = 0; s_or_b < 2; s_or_b++) {
782
783 int bin =
784 m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->FindBin(varVector[var]);
785
786 double prob = 0;
788 prob = InterpolatePdfs(s_or_b, ipbin, et, eta, bin, var);
789 } else {
790 double integral =
791 double(m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->Integral());
792 if (integral == 0) {
793 ATH_MSG_WARNING("Error! PDF integral == 0!");
794 return -1.35;
795 }
796
797 prob =
798 double(
799 m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->GetBinContent(bin)) /
800 integral;
801 }
802
803 if (s_or_b == 0) {
804 SigmaS *= prob;
805 } else if (s_or_b == 1) {
806 SigmaB *= prob;
807 }
808 }
809 }
810
811 return TransformLikelihoodOutput(SigmaS, SigmaB, ip, et, eta);
812}
float et(const xAOD::jFexSRJetRoI *j)
std::unique_ptr< EGSelectors::SafeTH1 > m_fPDFbins[2][IP_BINS][s_fnEtBinsHist][s_fnEtaBins][s_fnVariables]
static constexpr unsigned int s_fnEtBinsHist
static const std::string s_fVariables[s_fnVariables]
static unsigned int getIpBin(double ip)
static unsigned int getLikelihoodEtHistBin(double eT)
Coarse Et binning. Used for the likelihood pdfs.
static constexpr unsigned int s_fnVariables
double InterpolatePdfs(unsigned int s_or_b, unsigned int ipbin, double et, double eta, int bin, unsigned int var) const
double TransformLikelihoodOutput(double ps, double pb, double ip, double et, double eta) const
Apply a transform to zoom into the LH output peaks.
double integral(TH1 *h)
Definition computils.cxx:59

◆ evaluateLikelihood() [2/2]

double Root::TElectronLikelihoodTool::evaluateLikelihood ( const std::vector< float > & varVector,
double et,
double eta,
double ip = 0 ) const
private

Definition at line 704 of file TElectronLikelihoodTool.cxx.

709{
710 std::vector<double> vec;
711 for (unsigned int var = 0; var < s_fnVariables; var++) {
712 vec.push_back(varVector[var]);
713 }
714 return evaluateLikelihood(vec, et, eta, ip);
715}

◆ getAcceptInfo()

const asg::AcceptInfo & Root::TElectronLikelihoodTool::getAcceptInfo ( ) const
inline

accesss to the accept info object

Definition at line 91 of file TElectronLikelihoodTool.h.

91{ return m_acceptInfo; }

◆ getBinName()

std::string Root::TElectronLikelihoodTool::getBinName ( int etbin,
int etabin,
int ipbin,
const std::string & iptype )
staticprivate

Definition at line 1062 of file TElectronLikelihoodTool.cxx.

1066{
1067 const double eta_bounds[9] = { 0.0, 0.6, 0.8, 1.15, 1.37, 1.52, 1.81, 2.01, 2.37 };
1068 const int et_bounds[s_fnEtBinsHist] = { 4, 7, 10, 15, 20, 30, 40 };
1069 if (!iptype.empty()) {
1070 return std::format("{}{}et{:02}eta{:.2f}",
1071 iptype, int(s_fIpBounds[ipbin]),
1072 et_bounds[etbin], eta_bounds[etabin]);
1073 } else {
1074 return std::format("et{}eta{:.2f}",
1075 et_bounds[etbin], eta_bounds[etabin]);
1076 }
1077}
static constexpr double s_fIpBounds[IP_BINS+1]

◆ getBitmask()

unsigned int Root::TElectronLikelihoodTool::getBitmask ( void ) const
inline

Definition at line 148 of file TElectronLikelihoodTool.h.

148{ return m_variableBitMask; }

◆ getIpBin()

unsigned int Root::TElectronLikelihoodTool::getIpBin ( double ip)
staticprivate

Definition at line 975 of file TElectronLikelihoodTool.cxx.

976{
977 for (unsigned int ipBin = 0; ipBin < IP_BINS; ++ipBin) {
978 if (ip < s_fIpBounds[ipBin + 1])
979 return ipBin;
980 }
981 return 0;
982}

◆ getLikelihoodBitmask()

unsigned int Root::TElectronLikelihoodTool::getLikelihoodBitmask ( const std::string & vars) const
private

Definition at line 1080 of file TElectronLikelihoodTool.cxx.

1082{
1083 unsigned int mask = 0x0;
1084 ATH_MSG_DEBUG("Variables to be used: ");
1085 for (unsigned int var = 0; var < s_fnVariables; var++) {
1086 if (vars.find(s_fVariables[var]) != std::string::npos) {
1088 mask = mask | 0x1 << var;
1089 }
1090 }
1091 ATH_MSG_DEBUG("mask: " << mask);
1092 return mask;
1093}

◆ getLikelihoodEtaBin()

unsigned int Root::TElectronLikelihoodTool::getLikelihoodEtaBin ( double eta)
staticprivate

Eta binning for pdfs and discriminant cuts.

Definition at line 987 of file TElectronLikelihoodTool.cxx.

988{
989 const unsigned int nEtaBins = s_fnEtaBins;
990 const double etaBins[nEtaBins] = { 0.1, 0.6, 0.8, 1.15, 1.37,
991 1.52, 1.81, 2.01, 2.37, 2.47 };
992
993 for (unsigned int etaBin = 0; etaBin < nEtaBins; ++etaBin) {
994 if (std::abs(eta) < etaBins[etaBin])
995 return etaBin;
996 }
997
998 return 9;
999}
constexpr unsigned nEtaBins
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap etaBin

◆ getLikelihoodEtDiscBin()

unsigned int Root::TElectronLikelihoodTool::getLikelihoodEtDiscBin ( double eT,
const bool isLHbinning ) const
private

Fine Et binning. Used for the likelihood discriminant cuts.

Definition at line 1023 of file TElectronLikelihoodTool.cxx.

1026{
1027 const double GeV = 1000;
1028
1029 if (m_useOneExtraHighETLHBin && isLHbinning) {
1030 const unsigned int nEtBins = s_fnDiscEtBinsOneExtra;
1031 const double eTBins[nEtBins] = {
1032 10 * GeV, 15 * GeV, 20 * GeV,
1033 25 * GeV, 30 * GeV, 35 * GeV,
1034 40 * GeV, 45 * GeV, m_highETBinThreshold * GeV,
1035 6000 * GeV
1036 };
1037
1038 for (unsigned int eTBin = 0; eTBin < nEtBins; ++eTBin) {
1039 if (eT < eTBins[eTBin])
1040 return eTBin;
1041 }
1042
1043 return nEtBins - 1; // Return the last bin if > the last bin.
1044 }
1045
1046 const unsigned int nEtBins = s_fnDiscEtBins;
1047 const double eTBins[nEtBins] = { 10 * GeV, 15 * GeV, 20 * GeV,
1048 25 * GeV, 30 * GeV, 35 * GeV,
1049 40 * GeV, 45 * GeV, 50 * GeV };
1050
1051 for (unsigned int eTBin = 0; eTBin < nEtBins; ++eTBin) {
1052 if (eT < eTBins[eTBin])
1053 return eTBin;
1054 }
1055
1056 return nEtBins - 1; // Return the last bin if > the last bin.
1057}

◆ getLikelihoodEtHistBin()

unsigned int Root::TElectronLikelihoodTool::getLikelihoodEtHistBin ( double eT)
staticprivate

Coarse Et binning. Used for the likelihood pdfs.

Definition at line 1003 of file TElectronLikelihoodTool.cxx.

1004{
1005 const double GeV = 1000;
1006
1007 const unsigned int nEtBins = s_fnEtBinsHist;
1008 const double eTBins[nEtBins] = { 7 * GeV, 10 * GeV, 15 * GeV, 20 * GeV,
1009 30 * GeV, 40 * GeV, 50 * GeV };
1010
1011 for (unsigned int eTBin = 0; eTBin < nEtBins; ++eTBin) {
1012 if (eT < eTBins[eTBin]) {
1013 return eTBin;
1014 }
1015 }
1016
1017 return nEtBins - 1; // Return the last bin if > the last bin.
1018}

◆ initialize()

StatusCode Root::TElectronLikelihoodTool::initialize ( )

Initialize this class.

Definition at line 58 of file TElectronLikelihoodTool.cxx.

59{
60 ATH_MSG_DEBUG("TElectronLikelihoodTool initialize.");
61
62 // use an int as a StatusCode
63 StatusCode sc(StatusCode::SUCCESS);
64
65 // Check that all needed variables are setup
66 if (m_pdfFileName.empty()) {
67 ATH_MSG_WARNING("You need to specify the input PDF file name before you "
68 "call initialize() with "
69 "setPDFFileName('your/file/name.root') ");
70 sc = StatusCode::FAILURE;
71 }
72
73 unsigned int number_of_expected_bin_combinedLH;
75 number_of_expected_bin_combinedLH = s_fnDiscEtBinsOneExtra * s_fnEtaBins;
76 else
77 number_of_expected_bin_combinedLH = s_fnDiscEtBins * s_fnEtaBins;
78 unsigned int number_of_expected_bin_combinedOther =
80
81 if (m_cutLikelihood.size() != number_of_expected_bin_combinedLH) {
82 ATH_MSG_ERROR("Configuration issue : cutLikelihood expected size "
83 << number_of_expected_bin_combinedLH << " input size "
84 << m_cutLikelihood.size());
85 sc = StatusCode::FAILURE;
86 }
87
90 number_of_expected_bin_combinedLH) {
92 "Configuration issue : DiscHardCutForPileupTransform expected size "
93 << number_of_expected_bin_combinedLH << " input size "
95 sc = StatusCode::FAILURE;
96 }
97 }
100 number_of_expected_bin_combinedLH) {
101 ATH_MSG_ERROR("Configuration issue : "
102 "DiscHardCutSlopeForPileupTransform expected size "
103 << number_of_expected_bin_combinedLH << " input size "
105 sc = StatusCode::FAILURE;
106 }
107 }
108 if (!m_discLooseForPileupTransform.empty()) {
110 number_of_expected_bin_combinedLH) {
112 "Configuration issue : DiscLooseForPileupTransform expected size "
113 << number_of_expected_bin_combinedLH << " input size "
115 sc = StatusCode::FAILURE;
116 }
117 }
118
119 // d0 cut
120 if (!m_cutA0.empty()) {
121 if (m_cutA0.size() != number_of_expected_bin_combinedOther) {
122 ATH_MSG_ERROR("Configuration issue : CutA0 expected size "
123 << number_of_expected_bin_combinedOther << " input size "
124 << m_cutA0.size());
125 sc = StatusCode::FAILURE;
126 }
127 }
128
129 // deltaEta cut
130 if (!m_cutDeltaEta.empty()) {
131 if (m_cutDeltaEta.size() != number_of_expected_bin_combinedOther) {
132 ATH_MSG_ERROR("Configuration issue : CutDeltaEta expected size "
133 << number_of_expected_bin_combinedOther << " input size "
134 << m_cutDeltaEta.size());
135 sc = StatusCode::FAILURE;
136 }
137 }
138
139 // deltaPhiRes cut
140 if (!m_cutDeltaPhiRes.empty()) {
141 if (m_cutDeltaPhiRes.size() != number_of_expected_bin_combinedOther) {
142 ATH_MSG_ERROR("Configuration issue : CutDeltaPhiRes expected size "
143 << number_of_expected_bin_combinedOther << " input size "
144 << m_cutDeltaPhiRes.size());
145 sc = StatusCode::FAILURE;
146 }
147 }
148 if (sc == StatusCode::FAILURE) {
150 "Could NOT initialize! Please fix the errors mentioned above...");
151 return sc;
152 }
153
154 // --------------------------------------------------------------------------
155 // Register the cuts and check that the registration worked:
156 // NOTE: THE ORDER IS IMPORTANT!!! Cut0 corresponds to bit 0, Cut1 to bit
157 // 1,... Cut position for the kineatic pre-selection
158 m_cutPosition_kinematic = m_acceptInfo.addCut("kinematic", "pass kinematic");
159 if (m_cutPosition_kinematic < 0) {
160 sc = StatusCode::FAILURE;
161 }
162
163 // NSilicon
164 m_cutPosition_NSilicon = m_acceptInfo.addCut("NSilicon", "pass NSilicon");
165 if (m_cutPosition_NSilicon < 0) {
166 sc = StatusCode::FAILURE;
167 }
168
169 // NPixel
170 m_cutPosition_NPixel = m_acceptInfo.addCut("NPixel", "pass NPixel");
171 if (m_cutPosition_NPixel < 0) {
172 sc = StatusCode::FAILURE;
173 }
174
175 // NBlayer
176 m_cutPosition_NBlayer = m_acceptInfo.addCut("NBlayer", "pass NBlayer");
177 if (m_cutPosition_NBlayer < 0) {
178 sc = StatusCode::FAILURE;
179 }
180
181 // Ambiguity
182 m_cutPosition_ambiguity = m_acceptInfo.addCut("ambiguity", "pass ambiguity");
183 if (m_cutPosition_ambiguity < 0) {
184 sc = StatusCode::FAILURE;
185 }
186
187 // Cut position for the likelihood selection - DO NOT CHANGE ORDER!
188 m_cutPosition_LH = m_acceptInfo.addCut("passLH", "pass Likelihood");
189 if (m_cutPosition_LH < 0) {
190 sc = StatusCode::FAILURE;
191 }
192
193 // D0
195 m_acceptInfo.addCut("TrackA0", "A0 (aka d0) wrt beam spot < Cut");
196 if (m_cutPositionTrackA0 < 0) {
197 sc = StatusCode::FAILURE;
198 }
199
200 // deltaeta
202 "TrackMatchEta", "Track match deta in 1st sampling < Cut");
204 sc = StatusCode::FAILURE;
205 }
206
207 // deltaphi
209 "TrackMatchPhiRes", "Track match dphi in 2nd sampling, rescaled < Cut");
211 sc = StatusCode::FAILURE;
212 }
213
214 // Wstot
216 "WstotAtHighET", "Above HighETBinThreshold, Wstot < Cut");
218 sc = StatusCode::FAILURE;
219 }
220
221 // EoverP
223 "EoverPAtHighET", "Above HighETBinThreshold, EoverP < Cut");
225 sc = StatusCode::FAILURE;
226 }
227
228 // Check that we got everything OK
229 if (sc == StatusCode::FAILURE) {
231 "! Something went wrong with the setup of the decision objects...");
232 return sc;
233 }
234
235 // ----------------------------------
236 // Get the correct bit mask for the current likelihood operating point
238
239 //----------File/Histo operation------------------------------------
240 // Load the ROOT file containing the PDFs
241 TString tmpString(m_pdfFileName);
242 gSystem->ExpandPathName(tmpString);
243 std::string fname(tmpString.Data());
244 auto pdfFile = std::unique_ptr<TFile>(TFile::Open(fname.c_str(), "READ"));
245 // Check that we could load the ROOT file
246 if (!pdfFile) {
247 ATH_MSG_ERROR(" No ROOT file found here: " << m_pdfFileName);
248 return StatusCode::FAILURE;
249 }
250
251 // Load the histograms
252 for (unsigned int varIndex = 0; varIndex < s_fnVariables; varIndex++) {
253 const std::string& vstr = s_fVariables[varIndex];
254 // Skip the loading of PDFs for variables we don't care about for this
255 // operating point. If the string is empty (which is true in the default
256 // 2012 case), load all of them.
257 if (m_variableNames.find(vstr) == std::string::npos &&
258 !m_variableNames.empty()) {
259 continue;
260 }
261 loadVarHistograms(vstr, pdfFile.get(), varIndex);
262 }
263
264 // TFile close does not free the memory
265 pdfFile->Close();
266 //----------End File/Histo operation------------------------------------
267
268 ATH_MSG_DEBUG("Initialization complete for a LH tool with these specs:"
269 << "\n - pdfFileName : "
271 << "\n - Variable bitmask : "
273
274 ATH_MSG_DEBUG("\n - VariableNames : "
276 << "\n - (bool)CutBL (yes/no) : "
277 << (!m_cutBL.empty() ? "yes" : "no")
278 << "\n - (bool)CutPi (yes/no) : "
279 << (!m_cutPi.empty() ? "yes" : "no")
280 << "\n - (bool)CutSi (yes/no) : "
281 << (!m_cutSi.empty() ? "yes" : "no")
282 << "\n - (bool)CutAmbiguity (yes/no) : "
283 << (!m_cutAmbiguity.empty() ? "yes" : "no")
284 << "\n - (bool)doRemoveF3AtHighEt (yes/no) : "
285 << (m_doRemoveF3AtHighEt ? "yes" : "no")
286 << "\n - (bool)doRemoveTRTPIDAtHighEt (yes/no) : "
287 << (m_doRemoveTRTPIDAtHighEt ? "yes" : "no")
288 << "\n - (bool)doSmoothBinInterpolation (yes/no) : "
289 << (m_doSmoothBinInterpolation ? "yes" : "no")
290 << "\n - (bool)useOneExtraHighETLHBin(yes/no) : "
291 << (m_useOneExtraHighETLHBin ? "yes" : "no")
292 << "\n - (double)HighETBinThreshold : "
294 << "\n - (bool)doPileupTransform (yes/no) : "
295 << (m_doPileupTransform ? "yes" : "no")
296 << "\n - (bool)doCentralityTransform (yes/no) : "
297 << (m_doCentralityTransform ? "yes" : "no")
298 << "\n - (bool)CutLikelihood (yes/no) : "
299 << (!m_cutLikelihood.empty() ? "yes" : "no")
300 << "\n - (bool)CutLikelihoodPileupCorrection (yes/no) : "
301 << (!m_cutLikelihoodPileupCorrection.empty() ? "yes" : "no")
302 << "\n - (bool)CutA0 (yes/no) : "
303 << (!m_cutA0.empty() ? "yes" : "no")
304 << "\n - (bool)CutDeltaEta (yes/no) : "
305 << (!m_cutDeltaEta.empty() ? "yes" : "no")
306 << "\n - (bool)CutDeltaPhiRes (yes/no) : "
307 << (!m_cutDeltaPhiRes.empty() ? "yes" : "no")
308 << "\n - (bool)CutWstotAtHighET (yes/no) : "
309 << (!m_cutWstotAtHighET.empty() ? "yes" : "no")
310 << "\n - (bool)CutEoverPAtHighET (yes/no) : "
311 << (!m_cutEoverPAtHighET.empty() ? "yes" : "no"));
312 return sc;
313}
static Double_t sc
std::vector< double > m_discHardCutForPileupTransform
reference disc for very hard cut; used by pileup transform
std::vector< double > m_discHardCutSlopeForPileupTransform
reference slope on disc for very hard cut; used by pileup transform
std::vector< double > m_discLooseForPileupTransform
reference disc for a pileup independent loose menu; used by pileup transform
int loadVarHistograms(const std::string &vstr, TFile *pdfFile, unsigned int varIndex)
Load the variable histograms from the pdf file.
unsigned int getLikelihoodBitmask(const std::string &vars) const
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ initMessaging()

void AthMessaging::initMessaging ( ) const
privateinherited

Initialize our message level and MessageSvc.

This method should only be called once.

Definition at line 39 of file AthMessaging.cxx.

40{
42 // If user did not set an explicit level, set a default
43 if (m_lvl == MSG::NIL) {
44 m_lvl = m_imsg ?
45 static_cast<MSG::Level>( m_imsg.load()->outputLevel(m_nm) ) :
46 MSG::INFO;
47 }
48}
std::string m_nm
Message source name.
std::atomic< IMessageSvc * > m_imsg
MessageSvc pointer.
std::atomic< MSG::Level > m_lvl
Current logging level.
IMessageSvc * getMessageSvc(bool quiet=false)

◆ InterpolateCuts()

double Root::TElectronLikelihoodTool::InterpolateCuts ( const std::vector< double > & cuts,
const std::vector< double > & cuts_4gev,
double et,
double eta ) const
private

Definition at line 1099 of file TElectronLikelihoodTool.cxx.

1104{
1105
1106
1107 //get the value of the cut
1108 unsigned int etbinLH = getLikelihoodEtDiscBin(et, true);
1109 unsigned int etabin = getLikelihoodEtaBin(eta);
1110 unsigned int ibin_combinedLH = etbinLH * s_fnEtaBins + etabin;
1111 double cut = cuts.at(ibin_combinedLH);
1112
1113 //Special low pt cuts
1114 if (!cuts_4gev.empty()) {
1115 if (et < 7000.) {
1116 cut = cuts_4gev.at(etabin);
1117 }
1118 // Below 6 GeV we do not interpolate
1119 if (et < 6000) {
1120 return cut;
1121 }
1122 } else {// No special low Et cuts
1123 // and below 8500 we do not interpolate
1124 if (et < 8500.) {
1125 return cut;
1126 }
1127 }
1128 // We interpolate until 47500 (last bin)
1129 // size of array is s_fnDiscEtBins
1130 if (et > 47500. || !(etbinLH < s_fnDiscEtBins)) {
1131 return cut;
1132 }
1133 // Interpolate
1134 double bin_width = 5000.;
1135 if (7000. < et && et < 10000.) {
1136 bin_width = 3000.;
1137 }
1138 if (et < 7000.) {
1139 bin_width = 2000.;
1140 }
1141 const double GeV = 1000;
1142 const double eTBins[s_fnDiscEtBins] = { 8.5 * GeV, 12.5 * GeV, 17.5 * GeV,
1143 22.5 * GeV, 27.5 * GeV, 32.5 * GeV,
1144 37.5 * GeV, 42.5 * GeV, 47.5 * GeV };
1145 double bin_center = eTBins[etbinLH];
1146 if (et > bin_center) {
1147 double cut_next = cut;
1148 if (etbinLH + 1 <= 8)
1149 cut_next = cuts.at((etbinLH + 1) * s_fnEtaBins + etabin);
1150 return cut + (cut_next - cut) * (et - bin_center) / (bin_width);
1151 }
1152 // or else if et < bin_center :
1153 double cut_before = cut;
1154 if (etbinLH >= 1) {
1155 cut_before = cuts.at((etbinLH - 1) * s_fnEtaBins + etabin);
1156 } else if (etbinLH == 0 && !cuts_4gev.empty()) {
1157 cut_before = cuts_4gev.at(etabin);
1158 }
1159
1160 return cut - (cut - cut_before) * (bin_center - et) / (bin_width);
1161}
cut
This script demonstrates how to call a C++ class from Python Also how to use PyROOT is shown.

◆ InterpolatePdfs()

double Root::TElectronLikelihoodTool::InterpolatePdfs ( unsigned int s_or_b,
unsigned int ipbin,
double et,
double eta,
int bin,
unsigned int var ) const
private

Definition at line 1167 of file TElectronLikelihoodTool.cxx.

1173{
1174 // histograms exist for the following bins: 4, 7, 10, 15, 20, 30, 40.
1175 // Interpolation between histograms must follow fairly closely the
1176 // interpolation scheme between cuts - so be careful!
1177 int etbin = getLikelihoodEtHistBin(et); // hist binning
1178 int etabin = getLikelihoodEtaBin(eta);
1179 double integral =
1180 double(m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->Integral());
1181 double prob =
1182 double(m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->GetBinContent(bin)) /
1183 integral;
1184
1185 int Nbins = m_fPDFbins[s_or_b][ipbin][etbin][etabin][var]->GetNbinsX();
1186 if (et > 42500.) {
1187 return prob; // interpolation stops here.
1188 }
1189 if (et < 6000.) {
1190 return prob; // interpolation stops here.
1191 }
1192 if (22500. < et && et < 27500.) {
1193 return prob; // region of non-interpolation for pdfs
1194 }
1195 if (32500. < et && et < 37500.) {
1196 return prob; // region of non-interpolation for pdfs
1197 }
1198 double bin_width = 5000.;
1199 if (7000. < et && et < 10000.) {
1200 bin_width = 3000.;
1201 }
1202 if (et < 7000.) {
1203 bin_width = 2000.;
1204 }
1205 const double GeV = 1000;
1206 const double eTHistBins[s_fnEtBinsHist] = { 6. * GeV, 8.5 * GeV,
1207 12.5 * GeV, 17.5 * GeV,
1208 22.5 * GeV, 32.5 * GeV,
1209 42.5 * GeV };
1210 double bin_center = eTHistBins[etbin];
1211 if (etbin == 4 && et >= 27500.) {
1212 bin_center = 27500.; // special: interpolate starting from 27.5 here
1213 }
1214 if (etbin == 5 && et >= 37500.) {
1215 bin_center = 37500.; // special: interpolate starting from 37.5 here
1216 }
1217 if (et > bin_center) {
1218 double prob_next = prob;
1219 if (etbin + 1 <= 6) {
1220 // account for potential histogram bin inequalities
1221 int NbinsPlus =
1222 m_fPDFbins[s_or_b][ipbin][etbin + 1][etabin][var]->GetNbinsX();
1223 int binplus = bin;
1224 if (Nbins < NbinsPlus) {
1225 binplus = int(round(bin * (Nbins / NbinsPlus)));
1226 } else if (Nbins > NbinsPlus) {
1227 binplus = int(round(bin * (NbinsPlus / Nbins)));
1228 }
1229 // do interpolation
1230 double integral_next =
1231 double(m_fPDFbins[s_or_b][ipbin][etbin + 1][etabin][var]->Integral());
1232 prob_next =
1233 double(m_fPDFbins[s_or_b][ipbin][etbin + 1][etabin][var]->GetBinContent(
1234 binplus)) /
1235 integral_next;
1236 return prob + (prob_next - prob) * (et - bin_center) / (bin_width);
1237 }
1238 }
1239 // or else if et < bin_center :
1240 double prob_before = prob;
1241 if (etbin - 1 >= 0) {
1242 // account for potential histogram bin inequalities
1243 int NbinsMinus =
1244 m_fPDFbins[s_or_b][ipbin][etbin - 1][etabin][var]->GetNbinsX();
1245 int binminus = bin;
1246 if (Nbins < NbinsMinus) {
1247 binminus = int(round(bin * (Nbins / NbinsMinus)));
1248 } else if (Nbins > NbinsMinus) {
1249 binminus = int(round(bin * (NbinsMinus / Nbins)));
1250 }
1251 double integral_before =
1252 double(m_fPDFbins[s_or_b][ipbin][etbin - 1][etabin][var]->Integral());
1253 prob_before =
1254 double(m_fPDFbins[s_or_b][ipbin][etbin - 1][etabin][var]->GetBinContent(
1255 binminus)) /
1256 integral_before;
1257 }
1258 return prob - (prob - prob_before) * (bin_center - et) / (bin_width);
1259}
float round(const float toRound, const unsigned int decimals)
Definition Mdt.cxx:27

◆ loadVarHistograms()

int Root::TElectronLikelihoodTool::loadVarHistograms ( const std::string & vstr,
TFile * pdfFile,
unsigned int varIndex )

Load the variable histograms from the pdf file.

Definition at line 316 of file TElectronLikelihoodTool.cxx.

319{
320 for (unsigned int s_or_b = 0; s_or_b < 2; ++s_or_b) {
321 for (unsigned int ip = 0; ip < IP_BINS; ++ip) {
322 for (unsigned int et = 0; et < s_fnEtBinsHist; ++et) {
323 for (unsigned int eta = 0; eta < s_fnEtaBins; ++eta) {
324
325 std::string sig_bkg = (s_or_b == 0) ? "sig" : "bkg";
326 // Because eta bins in the root file don't match up exactly with cut
327 // menu definitions, the second eta bin is an exact copy of the first,
328 // and all subsequent eta bins are pushed back by one.
329 unsigned int eta_tmp = (eta > 0) ? eta - 1 : eta;
330 // The 7-10 GeV, crack bin uses the 10-15 Gev pdfs. WE DO NOT DO THIS
331 // ANYMORE! unsigned int et_tmp = (eta == 5 && et == 1) ? 1 : et;
332 unsigned int et_tmp = et;
333 std::string binname = getBinName(et_tmp, eta_tmp, ip, m_ipBinning);
334
335 if (((std::string(binname).find("2.37") != std::string::npos)) &&
336 (vstr.find("el_f3") != std::string::npos)) {
337 continue;
338 }
339
340 if (((std::string(binname).find("2.01") != std::string::npos) ||
341 (std::string(binname).find("2.37") != std::string::npos)) &&
342 (vstr.find("TRT") != std::string::npos)) {
343 continue;
344 }
345
346 const std::string pdfdir = std::format("{}/{}", vstr, sig_bkg);
347
348 std::string pdf = std::format("{}_{}_smoothed_hist_from_KDE_{}",
349 vstr, sig_bkg, binname);
350
351 std::string pdf_newname = std::format("{}_{}_{}_LHtool_copy_{}",
352 m_name, vstr, sig_bkg, binname);
353
354 if (!pdfFile->GetListOfKeys()->Contains(vstr.c_str())) {
355 ATH_MSG_INFO("Warning: skipping variable "
356 << vstr << " because the folder does not exist.");
357 return 1;
358 }
359 if (!((TDirectory*)pdfFile->Get(vstr.c_str()))
360 ->GetListOfKeys()
361 ->Contains(sig_bkg.c_str())) {
362 ATH_MSG_INFO("Warning: skipping variable "
363 << vstr << " because the folder does not exist.");
364 return 1;
365 }
366
367 // If the 0th et bin (4-7 GeV) histogram does not exist in the root
368 // file, then just use the 7-10 GeV bin histogram. This should
369 // preserve backward compatibility
370 if (et == 0 && !((TDirectory*)pdfFile->Get(pdfdir.c_str()))
371 ->GetListOfKeys()
372 ->Contains(pdf.c_str())) {
373 binname = getBinName(et_tmp + 1, eta_tmp, ip, m_ipBinning);
374
375 pdf = std::format("{}_{}_smoothed_hist_from_KDE_{}",
376 vstr, sig_bkg, binname);
377
378 pdf_newname = std::format("{}_{}_{}_LHtool_copy4GeV_{}",
379 m_name, vstr, sig_bkg, binname);
380 }
381 if (((TDirectory*)pdfFile->Get(pdfdir.c_str()))
382 ->GetListOfKeys()
383 ->Contains(pdf.c_str())) {
384 TH1F* hist = (TH1F*)(((TDirectory*)pdfFile->Get(pdfdir.c_str()))->Get(pdf.c_str()));
385 m_fPDFbins[s_or_b][ip][et][eta][varIndex] =
386 std::make_unique<EGSelectors::SafeTH1>(hist);
387 } else {
388 ATH_MSG_INFO("Warning: Object " << pdf << " does not exist.");
389 ATH_MSG_INFO("Skipping all other histograms with this variable.");
390 return 1;
391 }
392 }
393 }
394 }
395 }
396 return 1;
397}
#define ATH_MSG_INFO(x)
static std::string getBinName(int etbin, int etabin, int ipbin, const std::string &iptype)
std::string find(const std::string &s)
return a remapped string
Definition hcg.cxx:138
TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)

◆ msg() [1/2]

MsgStream & asg::AsgMessaging::msg ( ) const
inherited

The standard message stream.

Returns
A reference to the default message stream of this object.

Definition at line 49 of file AsgMessaging.cxx.

49 {
50#ifndef XAOD_STANDALONE
51 return ::AthMessaging::msg();
52#else // not XAOD_STANDALONE
53 return m_msg;
54#endif // not XAOD_STANDALONE
55 }

◆ msg() [2/2]

MsgStream & asg::AsgMessaging::msg ( const MSG::Level lvl) const
inherited

The standard message stream.

Parameters
lvlThe message level to set the stream to
Returns
A reference to the default message stream, set to level "lvl"

Definition at line 57 of file AsgMessaging.cxx.

57 {
58#ifndef XAOD_STANDALONE
59 return ::AthMessaging::msg( lvl );
60#else // not XAOD_STANDALONE
61 m_msg << lvl;
62 return m_msg;
63#endif // not XAOD_STANDALONE
64 }

◆ msgLvl()

bool asg::AsgMessaging::msgLvl ( const MSG::Level lvl) const
inherited

Test the output level of the object.

Parameters
lvlThe message level to test against
Returns
boolean Indicting if messages at given level will be printed
true If messages at level "lvl" will be printed

Definition at line 41 of file AsgMessaging.cxx.

41 {
42#ifndef XAOD_STANDALONE
43 return ::AthMessaging::msgLvl( lvl );
44#else // not XAOD_STANDALONE
45 return m_msg.msgLevel( lvl );
46#endif // not XAOD_STANDALONE
47 }

◆ setBinning()

void Root::TElectronLikelihoodTool::setBinning ( const std::string & val)
inline

Define the binning.

Definition at line 146 of file TElectronLikelihoodTool.h.

◆ setBitmask()

void Root::TElectronLikelihoodTool::setBitmask ( unsigned int val)
inline

Definition at line 149 of file TElectronLikelihoodTool.h.

◆ setLevel()

void AthMessaging::setLevel ( MSG::Level lvl)
inherited

Change the current logging level.

Use this rather than msg().setLevel() for proper operation with MT.

Definition at line 28 of file AthMessaging.cxx.

29{
30 m_lvl = lvl;
31}

◆ setPDFFileName()

void Root::TElectronLikelihoodTool::setPDFFileName ( const std::string & val)
inline

Add an input file that holds the PDFs.

Definition at line 131 of file TElectronLikelihoodTool.h.

131{ m_pdfFileName = val; }

◆ setVariableNames()

void Root::TElectronLikelihoodTool::setVariableNames ( const std::string & val)
inline

Define the variable names.

Definition at line 134 of file TElectronLikelihoodTool.h.

◆ TransformLikelihoodOutput()

double Root::TElectronLikelihoodTool::TransformLikelihoodOutput ( double ps,
double pb,
double ip,
double et,
double eta ) const
private

Apply a transform to zoom into the LH output peaks.

Optionally do pileup correction too

Definition at line 816 of file TElectronLikelihoodTool.cxx.

821{
822 // returns transformed or non-transformed output
823 // (Taken mostly from TMVA likelihood code)
824 double fEpsilon = 1e-99;
825 // If both signal and bkg are 0, we want it to fail.
826 if (ps < fEpsilon)
827 ps = 0;
828 if (pb < fEpsilon)
829 pb = fEpsilon;
830 double disc = ps / double(ps + pb);
831
832 if (disc >= 1.0) {
833 disc = 1. - 1.e-15;
834 } else if (disc <= 0.0) {
835 disc = fEpsilon;
836 }
837
838 double tau = 15.0;
839 disc = -std::log(1.0 / disc - 1.0) * (1. / double(tau));
840
841 // Linearly transform the discriminant as a function of pileup, rather than
842 // the old scheme of changing the cut value based on pileup. This is simpler
843 // for the tuning, as well as ensuring subsets / making discriminants more
844 // transparent. In the HI case, a quadratic centrality transform is applied
845 // instead.
847
848 // The variables used by the transform:
849 //
850 // - hard_cut_ref = extremely tight discriminant cut as reference to ensure
851 // pileup correction for looser menus is less than for tighter menus.
852 // - loose_ref = a loose menu with no pileup correction. Any tighter
853 // menu with same inputs will necessarily have pileup dependence built in
854 // - disc_max = max disc value for which pileup correction is desired.
855 // - pileup_max = max nvtx or mu for calculating the transform. Any larger
856 // pileup values will use this maximum value in the transform.
857
862 "Vectors needed for pileup-dependent transform not correctly filled! "
863 "Skipping the transform.");
864 return disc;
865 }
866
869 ATH_MSG_WARNING("Vectors needed for centrality-dependent transform not "
870 "correctly filled! "
871 "Skipping the transform.");
872 return disc;
873 }
874
875 unsigned int etabin = getLikelihoodEtaBin(eta);
876
877 double disc_hard_cut_ref = 0;
878 double disc_hard_cut_ref_slope = 0;
879 double disc_hard_cut_ref_quad =
880 0; // only used for heavy ion implementation of the LH
881 double disc_loose_ref = 0;
882 double disc_max = m_discMaxForPileupTransform;
883 double pileup_max = m_pileupMaxForPileupTransform;
884
888 et,
889 eta);
890 disc_hard_cut_ref_slope =
893 et,
894 eta);
896 disc_hard_cut_ref_quad =
899 et,
900 eta);
903 et,
904 eta);
905 } else {
906 // default situation, in the case where 4-7 GeV bin is not defined
907 if (et > 7000. || m_discHardCutForPileupTransform4GeV.empty()) {
908 unsigned int etfinebinLH = getLikelihoodEtDiscBin(et, true);
909 unsigned int ibin_combined = etfinebinLH * s_fnEtaBins + etabin;
910 disc_hard_cut_ref = m_discHardCutForPileupTransform[ibin_combined];
911 disc_hard_cut_ref_slope =
914 disc_hard_cut_ref_quad =
916 disc_loose_ref = m_discLooseForPileupTransform[ibin_combined];
917 } else {
921 ATH_MSG_WARNING("Vectors needed for pileup-dependent transform not "
922 "correctly filled for 4-7 GeV "
923 "bin! Skipping the transform.");
924 return disc;
925 }
928 ATH_MSG_WARNING("Vectors needed for centrality-dependent transform "
929 "not correctly filled for 4-7 "
930 "GeV bin! Skipping the transform.");
931 return disc;
932 }
933 disc_hard_cut_ref = m_discHardCutForPileupTransform4GeV[etabin];
934 disc_hard_cut_ref_slope =
937 disc_hard_cut_ref_quad =
939 disc_loose_ref = m_discLooseForPileupTransform4GeV[etabin];
940 }
941 }
942
943 double ip_for_corr =
944 std::min(ip, pileup_max); // turn off correction for values > pileup_max
945 double disc_hard_cut_ref_prime =
946 disc_hard_cut_ref + disc_hard_cut_ref_slope * ip_for_corr +
947 disc_hard_cut_ref_quad * ip_for_corr * ip_for_corr;
948
949 if (disc <= disc_loose_ref) {
950 // Below threshold for applying pileup correction
951 } else if (disc <= disc_hard_cut_ref_prime) {
952 // Between the loose and hard cut reference points for pileup correction
953 double denom = double(disc_hard_cut_ref_prime - disc_loose_ref);
954 if (denom < 0.001)
955 denom = 0.001;
956 disc = disc_loose_ref + (disc - disc_loose_ref) *
957 (disc_hard_cut_ref - disc_loose_ref) / denom;
958 } else if (disc_hard_cut_ref_prime < disc && disc <= disc_max) {
959 // Between the hard cut and max reference points for pileup correction
960 double denom = double(disc_max - disc_hard_cut_ref_prime);
961 if (denom < 0.001)
962 denom = 0.001;
963 disc = disc_hard_cut_ref + (disc - disc_hard_cut_ref_prime) *
964 (disc_max - disc_hard_cut_ref) / denom;
965 }
966 }
967
968 ATH_MSG_DEBUG("disc is " << disc);
969 return disc;
970}
std::vector< double > m_discHardCutQuadForPileupTransform
reference quadratic apr on disc for very hard cut; used by centrality transform
std::vector< double > m_discLooseForPileupTransform4GeV
reference disc for a pileup independent loose menu; used by pileup transform - 4-7 GeV
std::vector< double > m_discHardCutQuadForPileupTransform4GeV
reference quadratic par on disc for very hard cut; used by centrality transform - 4-7 GeV
std::vector< double > m_discHardCutForPileupTransform4GeV
reference disc for very hard cut; used by pileup transform - 4-7 GeV
std::vector< double > m_discHardCutSlopeForPileupTransform4GeV
reference slope on disc for very hard cut; used by pileup transform

Member Data Documentation

◆ ATLAS_THREAD_SAFE

std::atomic_flag m_initialized AthMessaging::ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
mutableprivateinherited

Messaging initialized (initMessaging)

Definition at line 141 of file AthMessaging.h.

◆ m_acceptInfo

asg::AcceptInfo Root::TElectronLikelihoodTool::m_acceptInfo
private

Accept info.

Definition at line 277 of file TElectronLikelihoodTool.h.

◆ m_cutA0

std::vector<double> Root::TElectronLikelihoodTool::m_cutA0

cut max on track d0 bit

Definition at line 187 of file TElectronLikelihoodTool.h.

◆ m_cutAmbiguity

std::vector<int> Root::TElectronLikelihoodTool::m_cutAmbiguity

do cut on ambiguity bit

Definition at line 193 of file TElectronLikelihoodTool.h.

◆ m_cutBL

std::vector<int> Root::TElectronLikelihoodTool::m_cutBL

cut min on b-layer hits

Definition at line 181 of file TElectronLikelihoodTool.h.

◆ m_cutDeltaEta

std::vector<double> Root::TElectronLikelihoodTool::m_cutDeltaEta

do cut on delta eta bit

Definition at line 189 of file TElectronLikelihoodTool.h.

◆ m_cutDeltaPhiRes

std::vector<double> Root::TElectronLikelihoodTool::m_cutDeltaPhiRes

do cut on delta phi bit

Definition at line 191 of file TElectronLikelihoodTool.h.

◆ m_cutEoverPAtHighET

std::vector<double> Root::TElectronLikelihoodTool::m_cutEoverPAtHighET

do cut on EoverP above HighETBinThreshold bit

Definition at line 207 of file TElectronLikelihoodTool.h.

◆ m_cutLikelihood

std::vector<double> Root::TElectronLikelihoodTool::m_cutLikelihood

cut on likelihood output

Definition at line 213 of file TElectronLikelihoodTool.h.

◆ m_cutLikelihood4GeV

std::vector<double> Root::TElectronLikelihoodTool::m_cutLikelihood4GeV

cut on likelihood output, 4 GeV bin

Definition at line 217 of file TElectronLikelihoodTool.h.

◆ m_cutLikelihoodPileupCorrection

std::vector<double> Root::TElectronLikelihoodTool::m_cutLikelihoodPileupCorrection

pileup correction factor for cut on likelihood output

Definition at line 215 of file TElectronLikelihoodTool.h.

◆ m_cutLikelihoodPileupCorrection4GeV

std::vector<double> Root::TElectronLikelihoodTool::m_cutLikelihoodPileupCorrection4GeV

pileup correction factor for cut on likelihood output, 4 GeV bin

Definition at line 219 of file TElectronLikelihoodTool.h.

◆ m_cutPi

std::vector<int> Root::TElectronLikelihoodTool::m_cutPi

cut min on pixel hits

Definition at line 183 of file TElectronLikelihoodTool.h.

◆ m_cutPosition_ambiguity

int Root::TElectronLikelihoodTool::m_cutPosition_ambiguity
private

The position of the ambiguity cut bit in the AcceptInfo return object.

Definition at line 299 of file TElectronLikelihoodTool.h.

◆ m_cutPosition_kinematic

int Root::TElectronLikelihoodTool::m_cutPosition_kinematic
private

The position of the kinematic cut bit in the AcceptInfo return object.

Definition at line 287 of file TElectronLikelihoodTool.h.

◆ m_cutPosition_LH

int Root::TElectronLikelihoodTool::m_cutPosition_LH
private

The position of the likelihood cut bit in the AcceptInfo return object.

Definition at line 302 of file TElectronLikelihoodTool.h.

◆ m_cutPosition_NBlayer

int Root::TElectronLikelihoodTool::m_cutPosition_NBlayer
private

The position of the NBlayer cut bit in the AcceptInfo return object.

Definition at line 296 of file TElectronLikelihoodTool.h.

◆ m_cutPosition_NPixel

int Root::TElectronLikelihoodTool::m_cutPosition_NPixel
private

The position of the NPixel cut bit in the AcceptInfo return object.

Definition at line 293 of file TElectronLikelihoodTool.h.

◆ m_cutPosition_NSilicon

int Root::TElectronLikelihoodTool::m_cutPosition_NSilicon
private

The position of the NSilicon cut bit in the AcceptInfo return object.

Definition at line 290 of file TElectronLikelihoodTool.h.

◆ m_cutPositionEoverPAtHighET

int Root::TElectronLikelihoodTool::m_cutPositionEoverPAtHighET
private

The position of the high ET EoverP cut bit in the AcceptInfo return object.

Definition at line 319 of file TElectronLikelihoodTool.h.

◆ m_cutPositionTrackA0

int Root::TElectronLikelihoodTool::m_cutPositionTrackA0
private

The position of the d0 cut bit in the AcceptInfo return object.

Definition at line 305 of file TElectronLikelihoodTool.h.

◆ m_cutPositionTrackMatchEta

int Root::TElectronLikelihoodTool::m_cutPositionTrackMatchEta
private

The position of the deltaeta cut bit in the AcceptInfo return object.

Definition at line 308 of file TElectronLikelihoodTool.h.

◆ m_cutPositionTrackMatchPhiRes

int Root::TElectronLikelihoodTool::m_cutPositionTrackMatchPhiRes
private

The position of the deltaphi cut bit in the AcceptInfo return object.

Definition at line 311 of file TElectronLikelihoodTool.h.

◆ m_cutPositionWstotAtHighET

int Root::TElectronLikelihoodTool::m_cutPositionWstotAtHighET
private

The position of the high ET wstot cut bit in the AcceptInfo return object.

Definition at line 315 of file TElectronLikelihoodTool.h.

◆ m_cutSi

std::vector<int> Root::TElectronLikelihoodTool::m_cutSi

cut min on precision hits

Definition at line 185 of file TElectronLikelihoodTool.h.

◆ m_cutWstotAtHighET

std::vector<double> Root::TElectronLikelihoodTool::m_cutWstotAtHighET

do cut on wstot above HighETBinThreshold bit

Definition at line 205 of file TElectronLikelihoodTool.h.

◆ m_discHardCutForPileupTransform

std::vector<double> Root::TElectronLikelihoodTool::m_discHardCutForPileupTransform

reference disc for very hard cut; used by pileup transform

Definition at line 221 of file TElectronLikelihoodTool.h.

◆ m_discHardCutForPileupTransform4GeV

std::vector<double> Root::TElectronLikelihoodTool::m_discHardCutForPileupTransform4GeV

reference disc for very hard cut; used by pileup transform - 4-7 GeV

Definition at line 233 of file TElectronLikelihoodTool.h.

◆ m_discHardCutQuadForPileupTransform

std::vector<double> Root::TElectronLikelihoodTool::m_discHardCutQuadForPileupTransform

reference quadratic apr on disc for very hard cut; used by centrality transform

Definition at line 227 of file TElectronLikelihoodTool.h.

◆ m_discHardCutQuadForPileupTransform4GeV

std::vector<double> Root::TElectronLikelihoodTool::m_discHardCutQuadForPileupTransform4GeV

reference quadratic par on disc for very hard cut; used by centrality transform - 4-7 GeV

Definition at line 239 of file TElectronLikelihoodTool.h.

◆ m_discHardCutSlopeForPileupTransform

std::vector<double> Root::TElectronLikelihoodTool::m_discHardCutSlopeForPileupTransform

reference slope on disc for very hard cut; used by pileup transform

Definition at line 224 of file TElectronLikelihoodTool.h.

◆ m_discHardCutSlopeForPileupTransform4GeV

std::vector<double> Root::TElectronLikelihoodTool::m_discHardCutSlopeForPileupTransform4GeV

reference slope on disc for very hard cut; used by pileup transform

  • 4-7 GeV

Definition at line 236 of file TElectronLikelihoodTool.h.

◆ m_discLooseForPileupTransform

std::vector<double> Root::TElectronLikelihoodTool::m_discLooseForPileupTransform

reference disc for a pileup independent loose menu; used by pileup transform

Definition at line 230 of file TElectronLikelihoodTool.h.

◆ m_discLooseForPileupTransform4GeV

std::vector<double> Root::TElectronLikelihoodTool::m_discLooseForPileupTransform4GeV

reference disc for a pileup independent loose menu; used by pileup transform - 4-7 GeV

Definition at line 242 of file TElectronLikelihoodTool.h.

◆ m_discMaxForPileupTransform

double Root::TElectronLikelihoodTool::m_discMaxForPileupTransform

max discriminant for which pileup transform is to be used

Definition at line 244 of file TElectronLikelihoodTool.h.

◆ m_doCentralityTransform

bool Root::TElectronLikelihoodTool::m_doCentralityTransform

do centrality-dependent transform on discriminant value

Definition at line 211 of file TElectronLikelihoodTool.h.

◆ m_doPileupTransform

bool Root::TElectronLikelihoodTool::m_doPileupTransform

do pileup-dependent transform on discriminant value

Definition at line 209 of file TElectronLikelihoodTool.h.

◆ m_doRemoveF3AtHighEt

bool Root::TElectronLikelihoodTool::m_doRemoveF3AtHighEt

do remove f3 variable from likelihood at high Et (>80 GeV)

Definition at line 195 of file TElectronLikelihoodTool.h.

◆ m_doRemoveTRTPIDAtHighEt

bool Root::TElectronLikelihoodTool::m_doRemoveTRTPIDAtHighEt

do remove TRTPID variable from likelihood at high Et (>80 GeV)

Definition at line 197 of file TElectronLikelihoodTool.h.

◆ m_doSmoothBinInterpolation

bool Root::TElectronLikelihoodTool::m_doSmoothBinInterpolation

do smooth interpolation between bins

Definition at line 199 of file TElectronLikelihoodTool.h.

◆ m_fPDFbins

std::unique_ptr<EGSelectors::SafeTH1> Root::TElectronLikelihoodTool::m_fPDFbins[2][IP_BINS][s_fnEtBinsHist][s_fnEtaBins][s_fnVariables]
private

Definition at line 335 of file TElectronLikelihoodTool.h.

◆ m_highETBinThreshold

double Root::TElectronLikelihoodTool::m_highETBinThreshold

ET threshold for using high ET cuts and bin.

Definition at line 203 of file TElectronLikelihoodTool.h.

◆ m_imsg

std::atomic<IMessageSvc*> AthMessaging::m_imsg { nullptr }
mutableprivateinherited

MessageSvc pointer.

Definition at line 135 of file AthMessaging.h.

135{ nullptr };

◆ m_ipBinning

std::string Root::TElectronLikelihoodTool::m_ipBinning
private

Deprecated.

Definition at line 284 of file TElectronLikelihoodTool.h.

◆ m_lvl

std::atomic<MSG::Level> AthMessaging::m_lvl { MSG::NIL }
mutableprivateinherited

Current logging level.

Definition at line 138 of file AthMessaging.h.

138{ MSG::NIL };

◆ m_msg_tls

boost::thread_specific_ptr<MsgStream> AthMessaging::m_msg_tls
mutableprivateinherited

MsgStream instance (a std::cout like with print-out levels)

Definition at line 132 of file AthMessaging.h.

◆ m_name

std::string Root::TElectronLikelihoodTool::m_name
private

tool name

Definition at line 274 of file TElectronLikelihoodTool.h.

◆ m_nm

std::string AthMessaging::m_nm
privateinherited

Message source name.

Definition at line 129 of file AthMessaging.h.

◆ m_pdfFileName

std::string Root::TElectronLikelihoodTool::m_pdfFileName

Name of the pdf file.

Definition at line 250 of file TElectronLikelihoodTool.h.

◆ m_pileupMaxForPileupTransform

double Root::TElectronLikelihoodTool::m_pileupMaxForPileupTransform

max nvtx or mu to be used in pileup transform

Definition at line 246 of file TElectronLikelihoodTool.h.

◆ m_useOneExtraHighETLHBin

bool Root::TElectronLikelihoodTool::m_useOneExtraHighETLHBin

use one extra bin for high ET LH

Definition at line 201 of file TElectronLikelihoodTool.h.

◆ m_variableBitMask

unsigned int Root::TElectronLikelihoodTool::m_variableBitMask
private

The bitmask corresponding to the variables in the likelihood.

For internal use.

Definition at line 281 of file TElectronLikelihoodTool.h.

◆ m_variableNames

std::string Root::TElectronLikelihoodTool::m_variableNames

variables to use in the LH

Definition at line 248 of file TElectronLikelihoodTool.h.

◆ s_fIpBounds

double Root::TElectronLikelihoodTool::s_fIpBounds[IP_BINS+1] = { 0., 500. }
staticconstexprprivate

Definition at line 321 of file TElectronLikelihoodTool.h.

321{ 0., 500. };

◆ s_fnDiscEtBins

unsigned int Root::TElectronLikelihoodTool::s_fnDiscEtBins = 9
staticconstexprprivate

Definition at line 328 of file TElectronLikelihoodTool.h.

◆ s_fnDiscEtBinsOneExtra

unsigned int Root::TElectronLikelihoodTool::s_fnDiscEtBinsOneExtra = 10
staticconstexprprivate

Definition at line 331 of file TElectronLikelihoodTool.h.

◆ s_fnEtaBins

unsigned int Root::TElectronLikelihoodTool::s_fnEtaBins = 10
staticconstexprprivate

Definition at line 332 of file TElectronLikelihoodTool.h.

◆ s_fnEtBinsHist

unsigned int Root::TElectronLikelihoodTool::s_fnEtBinsHist = 7
staticconstexprprivate

Definition at line 325 of file TElectronLikelihoodTool.h.

◆ s_fnVariables

unsigned int Root::TElectronLikelihoodTool::s_fnVariables = 13
staticconstexprprivate

Definition at line 333 of file TElectronLikelihoodTool.h.

◆ s_fVariables

const std::string Root::TElectronLikelihoodTool::s_fVariables
staticprivate
Initial value:
= {
"el_d0significance",
"el_eratio",
"el_deltaeta1",
"el_f1",
"el_f3",
"el_reta",
"el_rhad",
"el_rphi",
"el_trackd0pvunbiased",
"el_weta2",
"el_DeltaPoverP",
"el_deltaphiRescaled",
"el_TRT_PID"
}

Definition at line 1264 of file TElectronLikelihoodTool.h.

1296 {
1297constexpr unsigned int IP_BINS = 1;
1298}
1299namespace LikeEnum {
1300
1301struct LHAcceptVars_t
1302{
1303 double likelihood;
1304 double eta;
1305 double eT;
1306 int nSiHitsPlusDeadSensors;
1309 uint8_t ambiguityBit;
1310 double d0;
1311 double deltaEta;
1312 double deltaphires;
1313 double wstot;
1314 double EoverP;
1315 double ip;
1316};
1317
1318struct LHCalcVars_t
1319{
1320 double eta;
1321 double eT;
1322 double f3;
1323 double rHad;
1324 double rHad1;
1325 double Reta;
1326 double w2;
1327 double f1;
1328 double eratio;
1329 double deltaEta;
1330 double d0;
1331 double d0sigma;
1332 double rphi;
1333 double deltaPoverP;
1334 double deltaphires;
1335 double TRT_PID;
1336 double ip;
1337};
1338}
1339
1340namespace Root {
1341class TElectronLikelihoodTool : public asg::AsgMessaging
1342{
1343
1344public:
1346 TElectronLikelihoodTool(const char* name = "TElectronLikelihoodTool");
1347
1348 // Main methods
1349public:
1352
1354 const asg::AcceptInfo& getAcceptInfo() const { return m_acceptInfo; }
1355
1357 asg::AcceptData accept(LikeEnum::LHAcceptVars_t& vars_struct) const;
1358 asg::AcceptData accept(double likelihood,
1359 double eta,
1360 double eT,
1361 int nSiHitsPlusDeadSensors,
1362 int nPixHitsPlusDeadSensors,
1363 bool passBLayerRequirement,
1364 uint8_t ambiguityBit,
1365 double d0,
1366 double deltaEta,
1367 double deltaphires,
1368 double wstot,
1369 double EoverP,
1370 double ip) const;
1372 asg::AcceptData accept() const { return asg::AcceptData(&m_acceptInfo); }
1373
1374 double calculate(LikeEnum::LHCalcVars_t& vars_struct) const;
1375 double calculate(double eta,
1376 double eT,
1377 double f3,
1378 double rHad,
1379 double rHad1,
1380 double Reta,
1381 double w2,
1382 double f1,
1383 double eratio,
1384 double deltaEta,
1385 double d0,
1386 double d0sigma,
1387 double rphi,
1388 double deltaPoverP,
1389 double deltaphires,
1390 double TRT_PID,
1391 double ip) const;
1392
1394 inline void setPDFFileName(const std::string& val) { m_pdfFileName = val; }
1395
1397 inline void setVariableNames(const std::string& val)
1398 {
1401 }
1402
1404 int loadVarHistograms(const std::string& vstr,
1405 TFile* pdfFile,
1406 unsigned int varIndex);
1407
1409 inline void setBinning(const std::string& val) { m_ipBinning = val; }
1410
1411 unsigned int getBitmask(void) const { return m_variableBitMask; }
1412 inline void setBitmask(unsigned int val) { m_variableBitMask = val; };
1413
1414 // Private methods
1415private:
1416 // For every input "varVector", make sure elements of vector are
1417 // in the same order as prescribed in fVariables
1418
1419 double evaluateLikelihood(const std::vector<double>& varVector,
1420 double et,
1421 double eta,
1422 double ip = 0) const;
1423
1424 double evaluateLikelihood(const std::vector<float>& varVector,
1425 double et,
1426 double eta,
1427 double ip = 0) const;
1428
1429 unsigned int getLikelihoodBitmask(const std::string& vars) const;
1430
1431 double InterpolateCuts(const std::vector<double>& cuts,
1432 const std::vector<double>& cuts_4gev,
1433 double et,
1434 double eta) const;
1435 double InterpolatePdfs(unsigned int s_or_b,
1436 unsigned int ipbin,
1437 double et,
1438 double eta,
1439 int bin,
1440 unsigned int var) const;
1441
1442public:
1444 std::vector<int> m_cutBL;
1446 std::vector<int> m_cutPi;
1448 std::vector<int> m_cutSi;
1450 std::vector<double> m_cutA0;
1452 std::vector<double> m_cutDeltaEta;
1453 // /** @brief do cut on delta phi bit*/
1454 std::vector<double> m_cutDeltaPhiRes;
1456 std::vector<int> m_cutAmbiguity;
1466 double m_highETBinThreshold;
1467 // /** @brief do cut on wstot above HighETBinThreshold bit*/
1468 std::vector<double> m_cutWstotAtHighET;
1469 // /** @brief do cut on EoverP above HighETBinThreshold bit*/
1470 std::vector<double> m_cutEoverPAtHighET;
1476 std::vector<double> m_cutLikelihood;
1478 std::vector<double> m_cutLikelihoodPileupCorrection;
1480 std::vector<double> m_cutLikelihood4GeV;
1482 std::vector<double> m_cutLikelihoodPileupCorrection4GeV;
1484 std::vector<double> m_discHardCutForPileupTransform;
1487 std::vector<double> m_discHardCutSlopeForPileupTransform;
1490 std::vector<double> m_discHardCutQuadForPileupTransform;
1493 std::vector<double> m_discLooseForPileupTransform;
1496 std::vector<double> m_discHardCutForPileupTransform4GeV;
1499 std::vector<double> m_discHardCutSlopeForPileupTransform4GeV;
1502 std::vector<double> m_discHardCutQuadForPileupTransform4GeV;
1505 std::vector<double> m_discLooseForPileupTransform4GeV;
1511 std::string m_variableNames;
1513 std::string m_pdfFileName;
1514
1515 // Private methods
1516private:
1519 double TransformLikelihoodOutput(double ps,
1520 double pb,
1521 double ip,
1522 double et,
1523 double eta) const;
1524
1526 static unsigned int getLikelihoodEtaBin(double eta) ;
1527
1529 static unsigned int getLikelihoodEtHistBin(double eT) ;
1530
1532 unsigned int getLikelihoodEtDiscBin(double eT, const bool isLHbinning) const;
1533
1534 // Private member variables
1535private:
1537 std::string m_name;
1538
1540 asg::AcceptInfo m_acceptInfo;
1541
1544 unsigned int m_variableBitMask;
1545
1547 std::string m_ipBinning;
1548
1551
1554
1557
1560
1563
1565 int m_cutPosition_LH;
1566
1569
1572
1575
1579
1583
1584 static constexpr double s_fIpBounds[IP_BINS + 1] = { 0., 500. };
1585
1586 // number of hists stored for original LH, including 4GeV bin (for backwards
1587 // compatibility)
1588 static constexpr unsigned int s_fnEtBinsHist = 7;
1589 // number of discs stored for original LH, excluding 4GeV bin (for
1590 // backwards compatibility)
1591 static constexpr unsigned int s_fnDiscEtBins = 9;
1592 // number of discs stored for original LH plus one for
1593 // HighETBinThreshold (useOneExtraHighETLHBin), excluding 4GeV bin
1594 static constexpr unsigned int s_fnDiscEtBinsOneExtra = 10;
1595 static constexpr unsigned int s_fnEtaBins = 10;
1596 static constexpr unsigned int s_fnVariables = 13;
1597 // 5D array of unique_ptr to SafeTH1 // [sig(0)/bkg(1)][ip][et][eta][variable]
1598 std::unique_ptr<EGSelectors::SafeTH1> m_fPDFbins[2][IP_BINS][s_fnEtBinsHist][s_fnEtaBins][s_fnVariables];
1599 static const std::string s_fVariables[s_fnVariables];
1600
1601 static unsigned int getIpBin(double ip) ;
1602 static std::string getBinName(int etbin,
1603 int etabin,
1604 int ipbin,
1605 const std::string& iptype) ;
1606};
1607
1608} // End: namespace Root
1609
1610#endif
void setBinning(const std::string &val)
Define the binning.
StatusCode initialize()
Initialize this class.
const asg::AcceptInfo & getAcceptInfo() const
accesss to the accept info object
TElectronLikelihoodTool(const char *name="TElectronLikelihoodTool")
Standard constructor.
void setPDFFileName(const std::string &val)
Add an input file that holds the PDFs.
void setVariableNames(const std::string &val)
Define the variable names.

The documentation for this class was generated from the following files: