ATLAS Offline Software
Loading...
Searching...
No Matches
CellFinder.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
8CellFinder::CellFinder(const std::string& type,
9 const std::string& name,
10 const IInterface * parent) :
11 DiTauToolBase(type, name, parent)
12{
13 declareInterface<DiTauToolBase > (this);
14}
15
16
17CellFinder::~CellFinder() = default;
18
19
21
22 return StatusCode::SUCCESS;
23}
24
25
27 const EventContext& /*ctx*/) const {
28
29 ATH_MSG_DEBUG("execute CellFinder...");
30
31 // get ditau and its seed jet
32
33 xAOD::DiTauJet* pDiTau = data->xAODDiTau;
34 if (!pDiTau) {
35 ATH_MSG_ERROR("no di-tau candidate given");
36 return StatusCode::FAILURE;
37 }
38
39 const xAOD::Jet* pSeed = data->seed;
40 if (!pSeed) {
41 ATH_MSG_WARNING("No jet seed given.");
42 return StatusCode::FAILURE;
43 }
44
45 std::vector<fastjet::PseudoJet> vSubjets = data->subjets;
46 if (vSubjets.empty()) {
47 ATH_MSG_WARNING("No subjets given. Continue without cell information.");
48 return StatusCode::SUCCESS;
49 }
50
51 // get clusters linked to the seed jet. Loop over clusters to get linked cells
52
53 std::bitset<200000> cellSeen;
54 std::vector<const CaloCell*> subjetCells;
55
56 // loop over seed jet constituents
57 for (const auto *const seedConst: pSeed->getConstituents()) {
58 // cast jet constituent to cluster object
59 const xAOD::CaloCluster* cluster = dynamic_cast<const xAOD::CaloCluster*>( seedConst->rawConstituent() );
60
61 // loop over cells which are linked to the cluster
62 for (const auto *const cc : *(cluster->getCellLinks())) {
63 // skip if pt<0 or cell already encountered
64 if (cc->pt() < 0) continue;
65 if (cellSeen.test(cc->caloDDE()->calo_hash())) continue;
66 // register cell hash as already seen
67 cellSeen.set(cc->caloDDE()->calo_hash());
68
69 TLorentzVector temp_cc_p4;
70 temp_cc_p4.SetPtEtaPhiM(cc->pt(), cc->eta(), cc->phi(), cc->m());
71
72 // check if cell is in one of the subjets cones
73 for (const auto& subjet : vSubjets) {
74 TLorentzVector temp_sub_p4;
75 temp_sub_p4.SetPtEtaPhiM(subjet.pt(), subjet.eta(), subjet.phi_std(), subjet.m());
76 if (temp_cc_p4.DeltaR(temp_sub_p4) < m_Rsubjet) {
77 subjetCells.push_back(cc);
78 }
79 }
80 }
81 }
82
83 ATH_MSG_DEBUG("subjetCells.size()=" << subjetCells.size());
84 data->subjetCells = subjetCells;
85
86 // write f_core
87 float f_core;
88 for (unsigned int i = 0; i < vSubjets.size(); i++) {
89 const fastjet::PseudoJet& subjet = vSubjets.at(i);
90 float ptAll = 0.;
91 float ptCore = 0.;
92
93 TLorentzVector temp_sub_p4;
94 temp_sub_p4.SetPtEtaPhiM(subjet.pt(), subjet.eta(), subjet.phi_std(), subjet.m());
95
96 for (const auto& cc : data->subjetCells) {
97
98 TLorentzVector temp_cc_p4;
99 temp_cc_p4.SetPtEtaPhiM(cc->pt(), cc->eta(), cc->phi(), cc->m());
100
101 if (temp_cc_p4.DeltaR(temp_sub_p4) < data->Rsubjet) {
102 ptAll += cc->pt();
103 }
104
105 if (temp_cc_p4.DeltaR(temp_sub_p4) < data->Rcore) {
106 ptCore += cc->pt();
107 }
108 }
109
110 if (ptAll != 0.)
111 f_core = ptCore/ptAll;
112 else
113 f_core = -999.;
114
115 ATH_MSG_DEBUG("subjet "<< i << ": f_core=" << f_core);
116 pDiTau->setfCore(i, f_core);
117 }
118
119 return StatusCode::SUCCESS;
120}
#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
Definition CellFinder.h:32
virtual StatusCode initialize() override
Tool initializer.
virtual ~CellFinder()
CellFinder(const std::string &type, const std::string &name, const IInterface *parent)
Definition CellFinder.cxx:8
virtual StatusCode execute(DiTauCandidateData *data, const EventContext &ctx) const override
Execute - called for each Ditau candidate.
DiTauToolBase(const std::string &type, const std::string &name, const IInterface *parent)
const CaloClusterCellLink * getCellLinks() const
Get a pointer to the CaloClusterCellLink object (const version)
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