ATLAS Offline Software
Loading...
Searching...
No Matches
SubjetBuilder.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
7// fastjet includes
8#include "fastjet/ClusterSequenceArea.hh"
9#include "fastjet/AreaDefinition.hh"
10
12 const std::string& name,
13 const IInterface * parent) :
14 DiTauToolBase(type, name, parent)
15{
16 declareInterface<DiTauToolBase > (this);
17}
18
19
21
22
24
25 return StatusCode::SUCCESS;
26}
27
28
30 const EventContext& /*ctx*/) const {
31
32 ATH_MSG_DEBUG("subjet builder executing...");
33
34 // get seed jet
35 xAOD::DiTauJet* pDiTau = data->xAODDiTau;
36
37 if (pDiTau == nullptr) {
38 ATH_MSG_ERROR("no di-tau candidate given");
39 return StatusCode::FAILURE;
40 }
41
42 const xAOD::Jet* pSeed = data->seed;
43 if (!pSeed) {
44 ATH_MSG_WARNING("no seed given. Di-Tau can not be reconstructed.");
45 return StatusCode::FAILURE;
46 }
47
48 // retrieve seed jet clusters
49
50 using namespace fastjet;
51 const xAOD::JetConstituentVector vConst = pSeed->getConstituents();
52
53 if (vConst.empty()) {
54 ATH_MSG_WARNING("cluster constituents could not be retrieved from seed jet");
55 return StatusCode::FAILURE;
56 }
57
58 std::vector<PseudoJet> vpjClusters;
59
60 for (const auto *cl: vConst) {
61
62 TLorentzVector temp_p4;
63 temp_p4.SetPtEtaPhiM(cl->pt(), cl->eta(), cl->phi(), cl->m());
64 PseudoJet c( temp_p4.Px(), temp_p4.Py(), temp_p4.Pz(), temp_p4.E());
65
66 vpjClusters.push_back(c);
67 }
68
69 // reconstruct subjets
70
71 // Jet and area definitions
72 JetDefinition jd = JetDefinition(antikt_algorithm, m_Rsubjet);
73 AreaDefinition area_def(active_area_explicit_ghosts,GhostedAreaSpec(SelectorAbsRapMax(4.0)));
74 ClusterSequenceArea cs(vpjClusters, jd, area_def);
75
76 // store (pt-sorted) subjets
77 std::vector<PseudoJet> vSubjets = sorted_by_pt( cs.inclusive_jets(m_ptmin) );
78 if (vSubjets.size()<=1) {
79 ATH_MSG_DEBUG("Found less than 2 subjets. Reject ditau candidate");
80 return StatusCode::FAILURE;
81 }
82
83 ATH_MSG_DEBUG("found "<< vSubjets.size() << " subjets");
84 for (const auto& subjet: vSubjets) {
85 ATH_MSG_DEBUG("pt: " << subjet.pt() << " eta: " << subjet.eta() << " phi: " << subjet.phi());
86 }
87
88 data->subjets = vSubjets;
89
90 vSubjets.clear();
91
92 return StatusCode::SUCCESS;
93}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
DiTauToolBase(const std::string &type, const std::string &name, const IInterface *parent)
virtual StatusCode execute(DiTauCandidateData *data, const EventContext &ctx) const override
Execute - called for each Ditau candidate.
Gaudi::Property< float > m_ptmin
SubjetBuilder(const std::string &type, const std::string &name, const IInterface *parent)
Constructor.
virtual ~SubjetBuilder()
Destructor.
virtual StatusCode initialize() override
Tool initializer.
Gaudi::Property< float > m_Rsubjet
A vector of jet constituents at the scale used during jet finding.
bool empty() const
true if vector is empty()
JetConstituentVector getConstituents() const
Return a vector of consituents. The object behaves like vector<const IParticle*>. See JetConstituentV...
Definition Jet_v1.cxx:147
Jet_v1 Jet
Definition of the current "jet version".
DiTauJet_v1 DiTauJet
Definition of the current version.
Definition DiTauJet.h:17