ATLAS Offline Software
Loading...
Searching...
No Matches
TrigEgammaTopoHypoTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5/*************************************************************************************************
6 **
7 ** File: Trigger/TrigHypothesis/TrigEgammaHypo/src/combos/TrigEgammaTopoHypoTool.h
8 **
9 ** Description: - Hypothesis Tool: search for electron/photon pairs with
10 ** invariant mass or deltaPhi in some interval; intended for e.g Z->ee, H->gg, etc
11 **
12 **
13 ** Author: Debottam Bakshi Gupta <debottam.bakshi.gupta@cern.ch>
14 **
15 ** Created: May 15 2021
16 **
17 **************************************************************************************************/
18
21
22namespace TCU = TrigCompositeUtils;
23
24TrigEgammaTopoHypoTool::TrigEgammaTopoHypoTool(const std::string& type, const std::string& name, const IInterface* parent)
25 : ComboHypoToolBase(type, name, parent) {}
26
27
29{
30 ATH_MSG_DEBUG("AcceptAll = " << m_acceptAll );
31 ATH_MSG_DEBUG("ApplyMassCut = " << m_applyMassCut );
32 ATH_MSG_DEBUG("ApplyDPhiCut = " << m_applyDPhiCut );
33 ATH_MSG_DEBUG("LowerMassCut = " << m_lowerMassEgammaClusterCut );
34 ATH_MSG_DEBUG("UpperMassCut = " << m_upperMassEgammaClusterCut );
35
36 if ( not m_monTool.name().empty() ) {
37 ATH_CHECK( m_monTool.retrieve() );
38 ATH_MSG_DEBUG("m_monTool name: " << m_monTool);
39 }
40
41 return StatusCode::SUCCESS;
42}
43
44bool TrigEgammaTopoHypoTool::executeAlg(const std::vector<Combo::LegDecision> &combination, const EventContext& /*ctx*/) const {
45 auto massOfProcessed = Monitored::Scalar( "MassOfProcessed" , -1.0);
46 auto dphiOfProcessed = Monitored::Scalar( "DphiOfProcessed" , -99. );
47 auto massOfAccepted = Monitored::Scalar( "MassOfAccepted" , -1.0);
48 auto dphiOfAccepted = Monitored::Scalar( "DphiOfAccepted" , -99. );
49 auto monitorIt = Monitored::Group( m_monTool, massOfProcessed, dphiOfProcessed, massOfAccepted, dphiOfAccepted);
50 //retrieve the elements
51 std::vector<ElementLink<xAOD::IParticleContainer>> selected_electrons;
52
53 // Expecting to only run over chains with two legs and one electron or photon required on each leg
54 // So should always have two objects from which to form the invariant mass
55 if(combination.size() != 2){
57 "Expecting to combine exactly two electrons/photons, but instead found "
58 << combination.size() << ". Will throw a runtime error");
59 throw std::runtime_error(
60 "Expecting to combine exactly two electrons/photons, but instead found " +
61 std::to_string(combination.size()));
62 }
63
64 for (auto el: combination){
65 ATH_MSG_DEBUG("found Combination: "<<combination);
66 auto EL= el.second;
68 selected_electrons.push_back(electronLink);
69 }
70 auto electronLink1=selected_electrons[0];
71 auto electronLink2=selected_electrons[1];
72 TLorentzVector hlv1 = (*electronLink1)->p4();
73 TLorentzVector hlv2 = (*electronLink2)->p4();
74 massOfProcessed = (hlv1+hlv2).M();
75 dphiOfProcessed = hlv1.DeltaPhi(hlv2);
76
77 ATH_MSG_DEBUG("Found two Electrons/Photons with deltaPhi " <<dphiOfProcessed);
78
79 ATH_MSG_DEBUG("Found two Electrons/Photons with mass " <<massOfProcessed);
80
81 // apply the cut
82 bool pass=true;
83
84 if(m_acceptAll)
85 ATH_MSG_DEBUG("Applying no selections");
86
88 bool FailsMassCut = massOfProcessed<m_lowerMassEgammaClusterCut;
89 if (m_upperMassEgammaClusterCut>=0) FailsMassCut |= massOfProcessed>m_upperMassEgammaClusterCut;
90 if (FailsMassCut){
91 pass = false;
92 ATH_MSG_DEBUG("Combination failed mass cut: " << massOfProcessed << " not in [" << m_lowerMassEgammaClusterCut << "," << m_upperMassEgammaClusterCut << "]");
93 }
94 else
95 massOfAccepted = (hlv1+hlv2).M();
96 ATH_MSG_DEBUG( " Invariant mass " << massOfAccepted << " is within [" <<m_lowerMassEgammaClusterCut<< "," << m_upperMassEgammaClusterCut << "] This selection passed! ");
97 }
98
99 if (m_applyDPhiCut){
100 bool FailsDPhiCut = dphiOfProcessed < m_thresholdDPhiCut;
101 if (FailsDPhiCut){
102 ATH_MSG_DEBUG("Combination failed deltaPhi cut: " << dphiOfProcessed << " is below " << m_thresholdDPhiCut);
103 pass = false;
104 }
105 else {
106 dphiOfAccepted = hlv1.DeltaPhi(hlv2);
107 ATH_MSG_DEBUG( " deltaPhi " << dphiOfAccepted << " is above the threshold "<<m_thresholdDPhiCut<<" This selection passed! ");
108 }
109 }
110
111 return pass;
112
113}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
ComboHypoToolBase(const std::string &type, const std::string &name, const IInterface *parent)
Group of local monitoring quantities and retain correlation when filling histograms
Declare a monitored scalar variable.
Gaudi::Property< float > m_lowerMassEgammaClusterCut
lower inv mass cut (e,cluster)
Gaudi::Property< bool > m_applyDPhiCut
virtual StatusCode initialize() override
Gaudi::Property< float > m_thresholdDPhiCut
TrigEgammaTopoHypoTool(const std::string &type, const std::string &name, const IInterface *parent)
ToolHandle< GenericMonitoringTool > m_monTool
virtual bool executeAlg(const std::vector< Combo::LegDecision > &combination, const EventContext &ctx) const override
Only a dummy implementation exists in ComboHypoToolBase.
Gaudi::Property< bool > m_acceptAll
Gaudi::Property< bool > m_applyMassCut
Gaudi::Property< float > m_upperMassEgammaClusterCut
upper inv mass cut (e,cluster)
This module defines the arguments passed from the BATCH driver to the BATCH worker.
const std::string & featureString()
LinkInfo< T > findLink(const Decision *start, const std::string &linkName, const bool suppressMultipleLinksWarning=false)
Perform a recursive search for ElementLinks of type T and name 'linkName', starting from Decision obj...