ATLAS Offline Software
Loading...
Searching...
No Matches
DiTauConstituentFinder.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
7
9 const std::string& name,
10 const IInterface * parent) :
11 DiTauToolBase(type, name, parent)
12{
13 declareInterface<DiTauToolBase > (this);
14}
15
16
18
19
21 return StatusCode::SUCCESS;
22}
23
24
26 const EventContext& /*ctx*/) const {
27
28 ATH_MSG_DEBUG("execute DiTauConstituentFinder...");
29
30 // get ditau and its seed jet
31
32 xAOD::DiTauJet* pDiTau = data->xAODDiTau;
33 if (!pDiTau) {
34 ATH_MSG_ERROR("no di-tau candidate given");
35 return StatusCode::FAILURE;
36 }
37
38 const xAOD::Jet* pSeed = data->seed;
39 if (!pSeed) {
40 ATH_MSG_WARNING("No jet seed given.");
41 return StatusCode::FAILURE;
42 }
43
44 std::vector<fastjet::PseudoJet> vSubjets = data->subjets;
45 if (vSubjets.empty()) {
46 ATH_MSG_WARNING("No subjets given. Continue without cluster information.");
47 return StatusCode::SUCCESS;
48 }
49
50 std::vector<TLorentzVector> subjetConstituentsP4;
51 // loop over seed jet constituents
52 for (const auto *const seedConst: pSeed->getConstituents()) {
53 TLorentzVector tmp_const_p4;
54 if (m_useRawConstit) {
55 // cast jet constituent to cluster object
56 const xAOD::CaloCluster* seedConstCluster = dynamic_cast<const xAOD::CaloCluster*>( seedConst->rawConstituent() );
57 if(seedConstCluster == nullptr) continue;
58 tmp_const_p4.SetPtEtaPhiM(seedConstCluster->pt(), seedConstCluster->eta(), seedConstCluster->phi(), seedConstCluster->m());
59 } else {
60 tmp_const_p4.SetPtEtaPhiM(seedConst->pt(), seedConst->eta(), seedConst->phi(), seedConst->m());
61 }
62 for (const auto& subjet : vSubjets) {
63 TLorentzVector tmp_subjet_p4;
64 tmp_subjet_p4.SetPtEtaPhiM(subjet.pt(), subjet.eta(), subjet.phi_std(), subjet.m());
65 if (tmp_const_p4.DeltaR(tmp_subjet_p4) < m_Rsubjet) {
66 subjetConstituentsP4.push_back(tmp_const_p4);
67 }
68 }
69 }
70
71 ATH_MSG_DEBUG("subjetConstituents.size()=" << subjetConstituentsP4.size());
72
73 std::vector<float> vec_f_core(vSubjets.size(), 0);
74 for (unsigned int i = 0; i < vSubjets.size(); i++) {
75 const fastjet::PseudoJet& subjet = vSubjets.at(i);
76 float ptAll = 0.;
77 float ptCore = 0.;
78 float f_core = 0.;
79
80 TLorentzVector temp_sub_p4;
81 temp_sub_p4.SetPtEtaPhiM(subjet.pt(), subjet.eta(), subjet.phi_std(), subjet.m());
82
83 for (TLorentzVector const_p4 : subjetConstituentsP4) {
84 // ATH_MSG_DEBUG("dR(const, subjet" << i << ") = " << const_p4.DeltaR(temp_sub_p4) << "const pt = " << const_p4.Pt());
85 if (const_p4.DeltaR(temp_sub_p4) < data->Rsubjet) {
86 ptAll += const_p4.Pt();
87 }
88
89 if (const_p4.DeltaR(temp_sub_p4) < data->Rcore) {
90 ptCore += const_p4.Pt();
91 }
92 }
93
94 if (ptAll != 0.)
95 f_core = ptCore/ptAll;
96 else
97 f_core = -999.;
98
99 ATH_MSG_DEBUG("subjet "<< i << ": f_cluster_core = " << f_core);
100 vec_f_core.at(i) = f_core;
101 if (not m_useRawConstit){
102 // set f_core for the subjet in the DiTauJet object in HLT only
103 pDiTau->setfCore(i, f_core);
104 }
105 }
106
107 const static SG::Accessor<std::vector<float>> mDecor("f_cluster_core");
108 mDecor(*pDiTau) = std::move(vec_f_core);
109 return StatusCode::SUCCESS;
110}
#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
Gaudi::Property< float > m_Rsubjet
DiTauConstituentFinder(const std::string &type, const std::string &name, const IInterface *parent)
Gaudi::Property< bool > m_useRawConstit
virtual StatusCode execute(DiTauCandidateData *data, const EventContext &ctx) const override
Execute - called for each Ditau candidate.
virtual ~DiTauConstituentFinder()
virtual StatusCode initialize() override
Tool initializer.
DiTauToolBase(const std::string &type, const std::string &name, const IInterface *parent)
Helper class to provide type-safe access to aux data.
virtual double pt() const
The transverse momentum ( ) of the particle (negative for negative-energy clusters)
virtual double m() const
The invariant mass of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
void setfCore(unsigned int numSubjet, float fCore)
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".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
DiTauJet_v1 DiTauJet
Definition of the current version.
Definition DiTauJet.h:17