ATLAS Offline Software
Loading...
Searching...
No Matches
TauShotFinder.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#ifndef XAOD_ANALYSIS
6
7#include "TauShotFinder.h"
10
14
15
16TauShotFinder::TauShotFinder(const std::string& name) :
17 TauRecToolBase(name) {}
18
20
21 ATH_CHECK(m_caloWeightTool.retrieve());
24 ATH_CHECK(detStore()->retrieve (m_calo_id, "CaloCell_ID"));
25 ATH_MSG_INFO("Find TauShot in context: " << (inEleRM() ? "`EleRM`" : "`Standard`") << ", with Electron cell removal Flag: " << m_removeElectronCells);
26 return StatusCode::SUCCESS;
27}
28
29
30
32 xAOD::PFOContainer& shotPFOContainer) const {
33
34 // Any tau needs to have shot PFO vectors. Set empty vectors before nTrack cut
35 std::vector<ElementLink<xAOD::PFOContainer>> empty;
37
38 // Only run on 0-5 prong taus
40 return StatusCode::SUCCESS;
41 }
42
44 if (!caloCellInHandle.isValid()) {
45 ATH_MSG_ERROR ("Could not retrieve HiveDataObj with key " << caloCellInHandle.key());
46 return StatusCode::FAILURE;
47 }
48 const CaloCellContainer *cellContainer = caloCellInHandle.cptr();
49
50 // Select seed cells:
51 // -- dR < 0.4, EM1, pt > 100
52 // -- largest pt among the neighbours in eta direction
53 // -- no other seed cell as neighbour in eta direction
54 std::vector<const CaloCell*> seedCells;
55 ATH_CHECK(selectSeedCells(tau, *cellContainer, seedCells));
56 ATH_MSG_DEBUG("seedCells.size() = " << seedCells.size());
57
58 // Construt shot by merging neighbour cells in phi direction
59 while (!seedCells.empty()) {
60 // Find the neighbour in phi direction, and choose the one with highest pt
61 const CaloCell* cell = seedCells.front();
62 const CaloCell* phiNeigCell = getPhiNeighbour(*cell, seedCells);
63
64 // Construct shot PFO candidate
65 xAOD::PFO* shot = new xAOD::PFO();
66 shotPFOContainer.push_back(shot);
67
68 // -- Construct the shot cluster
69 xAOD::CaloCluster* shotCluster = createShotCluster(cell, phiNeigCell, *cellContainer, &shotClusterContainer);
70
72 clusElementLink.toContainedElement( shotClusterContainer, shotCluster );
73 shot->setClusterLink( clusElementLink );
74
75 // -- Calculate the four momentum
76 // TODO: simplify the calculation
77 if (phiNeigCell) {
78 // interpolate position
79 double wtCell_phiNeigCell = m_caloWeightTool->wtCell(phiNeigCell);
80 double wtCell_cell = m_caloWeightTool->wtCell(cell);
81 double dPhi = TVector2::Phi_mpi_pi( phiNeigCell->phi() - cell->phi());
82 double ratio = phiNeigCell->pt()*wtCell_phiNeigCell/(cell->pt()*wtCell_cell + phiNeigCell->pt()*wtCell_phiNeigCell);
83 float phi = cell->phi()+dPhi*ratio;
84 float pt = cell->pt()*wtCell_cell+phiNeigCell->pt()*wtCell_phiNeigCell;
85
86 shot->setP4( static_cast<float>(pt), static_cast<float>(cell->eta()), static_cast<float>(phi), static_cast<float>(cell->m()));
87 }
88 else {
89 shot->setP4( static_cast<float>(cell->pt()*m_caloWeightTool->wtCell(cell)), static_cast<float>(cell->eta()), static_cast<float>(cell->phi()), static_cast<float>(cell->m()));
90 }
91
92 // -- Set the Attribute
93 shot->setBDTPi0Score(-9999.);
94 shot->setCharge(0);
95 shot->setCenterMag(0.0);
96
97 shot->setAttribute<int>(xAOD::PFODetails::PFOAttributes::tauShots_nCellsInEta, m_nCellsInEta);
98
99 const IdentifierHash seedHash = cell->caloDDE()->calo_hash();
100 shot->setAttribute<int>(xAOD::PFODetails::PFOAttributes::tauShots_seedHash, seedHash);
101
102 std::vector<std::vector<const CaloCell*>> cellBlock = TauShotVariableHelpers::getCellBlock(*shot, m_calo_id);
103
104 float pt1 = TauShotVariableHelpers::ptWindow(cellBlock, 1, m_caloWeightTool);
105 shot->setAttribute<float>(xAOD::PFODetails::PFOAttributes::tauShots_pt1, pt1);
106
107 float pt3 = TauShotVariableHelpers::ptWindow(cellBlock, 3, m_caloWeightTool);
108 shot->setAttribute<float>(xAOD::PFODetails::PFOAttributes::tauShots_pt3, pt3);
109
110 float pt5 = TauShotVariableHelpers::ptWindow(cellBlock, 5, m_caloWeightTool);
111 shot->setAttribute<float>(xAOD::PFODetails::PFOAttributes::tauShots_pt5, pt5);
112
113 int nPhotons = getNPhotons(cell->eta(), pt1);
114 shot->setAttribute<int>(xAOD::PFODetails::PFOAttributes::tauShots_nPhotons, nPhotons);
115
116 // Add Element link to the shot PFO container
117 ElementLink<xAOD::PFOContainer> PFOElementLink;
118 PFOElementLink.toContainedElement(shotPFOContainer, shot);
119 tau.addShotPFOLink(PFOElementLink);
120
121 // Remove used cells from list
122 auto cellIndex = std::find(seedCells.begin(), seedCells.end(), cell);
123 seedCells.erase(cellIndex);
124 if (phiNeigCell) {
125 cellIndex = std::find(seedCells.begin(), seedCells.end(), phiNeigCell);
126 seedCells.erase(cellIndex);
127 }
128 } // Loop over seed cells
129
130 return StatusCode::SUCCESS;
131}
132
133
134
136 float absEta=std::abs(eta);
137
138 if (absEta < 0.80) {
139 return 0; // Central Barrel
140 }
141 if (absEta<1.39) {
142 return 1; // Outer Barrel
143 }
144 if (absEta<1.51) {
145 return 2; // Crack region
146 }
147 if (absEta<1.80) {
148 return 3; // Endcap, fine granularity
149 }
150 return 4; // Endcap, coarse granularity
151}
152
153
154
155int TauShotFinder::getNPhotons(float eta, float energy) const {
156 int etaBin = getEtaBin(eta);
157
158 // No photons in crack region
159 if(etaBin==2) return 0;
160
161 const std::vector<float>& minPtCut = m_minPtCut.value();
162 const std::vector<float>& doubleShotCut = m_doubleShotCut.value();
163 ATH_MSG_DEBUG("etaBin = " << etaBin << ", energy = " << energy);
164 ATH_MSG_DEBUG("MinPtCut: " << minPtCut.at(etaBin) << "DoubleShotCut: " << doubleShotCut.at(etaBin));
165
166 if (energy < minPtCut.at(etaBin)) return 0;
167 if (energy > doubleShotCut.at(etaBin)) return 2;
168 return 1;
169}
170
171
172
174 std::vector<const CaloCell*>& cells) const {
175 // if in EleRM tau reco, do electron cell removal
176 std::vector<const CaloCell*> removed_cells;
179 if (!removedClustersHandle.isValid()){
180 ATH_MSG_ERROR ("Could not retrieve HiveDataObj with key " << removedClustersHandle.key());
181 return StatusCode::FAILURE;
182 }
183 const xAOD::CaloClusterContainer *removed_clusters_cont = removedClustersHandle.cptr();
184
185 for (auto cluster : *removed_clusters_cont){
186 for(auto cell_it = cluster->cell_cbegin(); cell_it != cluster->cell_cend(); cell_it++){
187 removed_cells.push_back(*cell_it);
188 }
189 }
190 }
191
192 // retrieve EM1 cells within dR=0.4, pre-selected by TauPi0CreateROI
193 static const SG::ConstAccessor<std::vector<const CaloCell*>> acc_shotCells("shotCells");
194 std::vector<const CaloCell*> shotCells = acc_shotCells(tau);
195
196 for (const CaloCell* cell : shotCells) {
197 // Require cells above threshold (100 MeV by default)
198 // FIXME: cells are not corrected to point at tau vertex
199 if (cell->pt() * m_caloWeightTool->wtCell(cell) < m_energyThreshold) continue;
200 // if in EleRM, check the clusters do not include electron activities
201 if (m_removeElectronCells && inEleRM() && std::find(removed_cells.cbegin(), removed_cells.cend(), cell) != removed_cells.cend()) continue;
202
203 cells.push_back(cell);
204 }
205 return StatusCode::SUCCESS;
206}
207
208
209
211 const CaloCellContainer& cellContainer,
212 std::vector<const CaloCell*>& seedCells) const {
213
214 // Apply pre-selection of the cells
215 assert(seedCells.empty());
216 std::vector<const CaloCell*> cells;
217 ATH_CHECK(selectCells(tau, cells));
218 std::sort(cells.begin(),cells.end(),ptSort(*this));
219
220 std::set<IdentifierHash> seedCellHashes;
221
222 // Loop the pt sorted cells, and select the seed cells
223 for (const CaloCell* cell: cells) {
224 const IdentifierHash cellHash = cell->caloDDE()->calo_hash();
225
226 std::vector<IdentifierHash> nextEtaHashes;
227 m_calo_id->get_neighbours(cellHash, LArNeighbours::nextInEta, nextEtaHashes);
228 std::vector<IdentifierHash> prevEtaHashes;
229 m_calo_id->get_neighbours(cellHash, LArNeighbours::prevInEta, prevEtaHashes);
230
231 std::vector<IdentifierHash> neighHashes = nextEtaHashes;
232 neighHashes.insert(neighHashes.end(),prevEtaHashes.begin(),prevEtaHashes.end());
233
234 // Check whether it is a seed cell
235 bool status = true;
236 for (const IdentifierHash& neighHash : neighHashes) {
237 // Seed cells must not have seed cells as neighbours
238 // TODO: maybe this requirement can be removed
239 if (seedCellHashes.find(neighHash) != seedCellHashes.end()) {
240 status = false;
241 break;
242 }
243
244 // Pt of seed cells must be larger than neighbours'
245 const CaloCell* neighCell = cellContainer.findCell(neighHash);
246 if (!neighCell) continue;
247 if (neighCell->pt() * m_caloWeightTool->wtCell(neighCell) >= cell->pt() * m_caloWeightTool->wtCell(cell)) {
248 status = false;
249 break;
250 }
251 } // End of the loop of neighbour cells
252
253 if (!status) continue;
254
255 seedCells.push_back(cell);
256 seedCellHashes.insert(cellHash);
257 } // End of the loop of cells
258
259 return StatusCode::SUCCESS;
260}
261
262
263
265 std::vector<IdentifierHash> neigHashes;
266
267 // Next cell in phi direction
268 m_calo_id->get_neighbours(cell1Hash,LArNeighbours::nextInPhi,neigHashes);
269 if (neigHashes.size() > 1) {
270 ATH_MSG_DEBUG(cell1Hash << " has " << neigHashes.size() << " neighbours in the next phi direction !");
271 }
272 if (std::find(neigHashes.begin(), neigHashes.end(), cell2Hash) != neigHashes.end()) {
273 return true;
274 }
275
276 // Previous cell in phi direction
277 m_calo_id->get_neighbours(cell1Hash,LArNeighbours::prevInPhi,neigHashes);
278 if (neigHashes.size() > 1) {
279 ATH_MSG_DEBUG(cell1Hash << " has " << neigHashes.size() << " neighbours in the previous phi direction !");
280 }
281 return std::find(neigHashes.begin(), neigHashes.end(), cell2Hash) != neigHashes.end();
282}
283
284
285
287 const std::vector<const CaloCell*>& seedCells) const {
288
289 const IdentifierHash seedHash = seedCell.caloDDE()->calo_hash();
290
291 // Obtain the neighbour cells in the phi direction
292 std::vector<const CaloCell*> neighCells;
293 for (const CaloCell* neighCell : seedCells) {
294 if (neighCell == &seedCell) continue;
295
296 IdentifierHash neighHash = neighCell->caloDDE()->calo_hash();
297 if (this->isPhiNeighbour(seedHash, neighHash)) {
298 neighCells.push_back(neighCell);
299 }
300 }
301 std::sort(neighCells.begin(),neighCells.end(),ptSort(*this));
302
303 // Select the one with largest pt
304 const CaloCell* phiNeigCell = nullptr;
305 if (!neighCells.empty()) {
306 phiNeigCell = neighCells[0];
307 }
308
309 return phiNeigCell;
310}
311
312
313
314std::vector<const CaloCell*> TauShotFinder::getEtaNeighbours(const CaloCell& cell,
315 const CaloCellContainer& cellContainer,
316 int maxDepth) const {
317 std::vector<const CaloCell*> cells;
318
319 // Add neighbours in next eta direction
320 this->addEtaNeighbours(cell, cellContainer, cells, 0, maxDepth, true);
321 // Add neighbours in previous eta direction
322 this->addEtaNeighbours(cell, cellContainer, cells, 0, maxDepth, false);
323
324 return cells;
325}
326
327
328
330 const CaloCellContainer& cellContainer,
331 std::vector<const CaloCell*>& cells,
332 int depth,
333 int maxDepth,
334 bool next) const {
335 ++depth;
336
337 if (depth > maxDepth) return;
338
339 const IdentifierHash cellHash = cell.caloDDE()->calo_hash();
340
341 std::vector<IdentifierHash> neigHashes;
342 if (next) {
343 m_calo_id->get_neighbours(cellHash,LArNeighbours::nextInEta,neigHashes);
344 }
345 else {
346 m_calo_id->get_neighbours(cellHash,LArNeighbours::prevInEta,neigHashes);
347 }
348
349 for (const IdentifierHash& hash : neigHashes) {
350 const CaloCell* newCell = cellContainer.findCell(hash);
351
352 if (!newCell) continue;
353
354 cells.push_back(newCell);
355 this->addEtaNeighbours(*newCell, cellContainer, cells, depth, maxDepth, next);
356
357 if (neigHashes.size() > 1) {
358 ATH_MSG_DEBUG(cellHash << " has " << neigHashes.size() << " neighbours in the eta direction !");
359 break;
360 }
361 }
362}
363
364
365
367 const CaloCell* phiNeigCell,
368 const CaloCellContainer& cellContainer,
369 xAOD::CaloClusterContainer* clusterContainer) const {
370
371 xAOD::CaloCluster* shotCluster = CaloClusterStoreHelper::makeCluster(clusterContainer,&cellContainer);
372
373 int maxDepth = (m_nCellsInEta - 1) / 2;
374
375 std::vector<const CaloCell*> windowNeighbours = this->getEtaNeighbours(*cell, cellContainer, maxDepth);
376 if (phiNeigCell) {
377 std::vector<const CaloCell*> mergeCells = this->getEtaNeighbours(*phiNeigCell, cellContainer, maxDepth);
378 windowNeighbours.push_back(phiNeigCell);
379 windowNeighbours.insert(windowNeighbours.end(), mergeCells.begin(), mergeCells.end());
380 }
381
382 shotCluster->getOwnCellLinks()->reserve(windowNeighbours.size()+1);
383 const IdentifierHash seedHash = cell->caloDDE()->calo_hash();
384 shotCluster->addCell(cellContainer.findIndex(seedHash), 1.);
385
386 for (const CaloCell* cell : windowNeighbours) {
387 shotCluster->addCell(cellContainer.findIndex(cell->caloDDE()->calo_hash()), 1.0);
388 }
389
390 CaloClusterKineHelper::calculateKine(shotCluster,true,true);
391
392 return shotCluster;
393}
394
395
396
398bool TauShotFinder::ptSort::operator()( const CaloCell* cell1, const CaloCell* cell2 ){
399 double pt1 = cell1->pt()*m_info.m_caloWeightTool->wtCell(cell1);
400 double pt2 = cell2->pt()*m_info.m_caloWeightTool->wtCell(cell2);
401 return pt1 > pt2;
402}
403
404#endif
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
static const Attributes_t empty
const ServiceHandle< StoreGateSvc > & detStore() const
Container class for CaloCell.
const CaloCell * findCell(const IdentifierHash theHash) const
fast find method given identifier hash.
int findIndex(const IdentifierHash theHash) const
Return index of the cell with a given hash.
Data object for each calorimeter readout cell.
Definition CaloCell.h:57
virtual double phi() const override final
get phi (through CaloDetDescrElement)
Definition CaloCell.h:375
const CaloDetDescrElement * caloDDE() const
get pointer to CaloDetDescrElement (data member)
Definition CaloCell.h:321
static void calculateKine(xAOD::CaloCluster *clu, const bool useweight=true, const bool updateLayers=true, const bool useGPUCriteria=false)
Helper class to calculate cluster kinematics based on cells.
static std::unique_ptr< xAOD::CaloCluster > makeCluster(const CaloCellContainer *cellCont)
Creates a valid CaloCluster with a private Aux-Store and CellLink container.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
This is a "hash" representation of an Identifier.
virtual double pt() const
transverse momentum
Helper class to provide constant type-safe access to aux data.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
virtual const std::string & key() const override final
Return the StoreGate ID for the referenced object.
TauRecToolBase(const std::string &name)
bool inEleRM() const
StatusCode selectCells(const xAOD::TauJet &tau, std::vector< const CaloCell * > &cells) const
Apply preselection of the cells Cells within dR < 0.4, in EM1, and pt > 100 MeV are selected.
virtual StatusCode executeShotFinder(xAOD::TauJet &pTau, xAOD::CaloClusterContainer &tauShotCaloClusContainer, xAOD::PFOContainer &tauShotPFOContainer) const override
void addEtaNeighbours(const CaloCell &cell, const CaloCellContainer &cellContainer, std::vector< const CaloCell * > &cells, int depth, int maxDepth, bool next) const
Get neighbour cells in the eta direction.
std::vector< const CaloCell * > getEtaNeighbours(const CaloCell &cell, const CaloCellContainer &cellContainer, int maxDepth) const
Get neighbour cells in the eta direction.
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_removedClusterInputContainer
StatusCode selectSeedCells(const xAOD::TauJet &tau, const CaloCellContainer &cellContainer, std::vector< const CaloCell * > &seedCells) const
Select the seed cells used to construct the shot Cells must sastisfy:
Gaudi::Property< int > m_nCellsInEta
xAOD::CaloCluster * createShotCluster(const CaloCell *cell, const CaloCell *phiNeighCell, const CaloCellContainer &cellContainer, xAOD::CaloClusterContainer *clusterContainer) const
Create the shot cluster Shot cluster contains 5x1 cells from the seed cell and hottestneighbour cell ...
SG::ReadHandleKey< CaloCellContainer > m_caloCellInputContainer
const CaloCell_ID * m_calo_id
calo cell navigation
ToolHandle< IHadronicCalibrationTool > m_caloWeightTool
Gaudi::Property< bool > m_removeElectronCells
TauShotFinder(const std::string &name)
int getNPhotons(float eta, float energy) const
Get NPhotons in shot.
const CaloCell * getPhiNeighbour(const CaloCell &seedCell, const std::vector< const CaloCell * > &seedCells) const
Get the hottest neighbour cell in the phi direction.
Gaudi::Property< std::vector< float > > m_minPtCut
Gaudi::Property< std::vector< float > > m_doubleShotCut
virtual StatusCode initialize() override
Tool initializer.
bool isPhiNeighbour(IdentifierHash cell1Hash, IdentifierHash cell2Hash) const
Check whether two cells are neighbours in the phi direction.
int getEtaBin(float eta) const
Get eta bin.
Gaudi::Property< float > m_energyThreshold
CaloClusterCellLink * getOwnCellLinks()
Get a pointer to the owned CaloClusterCellLink object (non-const version)
bool addCell(const unsigned index, const double weight)
Method to add a cell to the cluster (Beware: Kinematics not updated!)
void setAttribute(PFODetails::PFOAttributes AttributeType, const T &anAttribute)
Set a PFO Variable via enum - overwrite is allowed.
bool setClusterLink(const ElementLink< xAOD::CaloClusterContainer > &theCluster)
Set a cluster constituent - does NOT append to existing container.
Definition PFO_v1.cxx:549
void setCenterMag(float CenterMag)
set CenterMag moment needed for vertex correction
void setP4(const FourMom_t &vec)
set the 4-vec
Definition PFO_v1.cxx:107
void setBDTPi0Score(float BDTPi0Score)
set BDT Score used to classify clusters as Pi0 like or not
void setCharge(float charge)
set charge of PFO
void setShotPFOLinks(const PFOLinks_t &shotPFOs)
void addShotPFOLink(const ElementLink< PFOContainer > &pfo)
add a shot PFO to the tau
std::string depth
tag string for intendation
Definition fastadd.cxx:46
std::vector< std::vector< const CaloCell * > > getCellBlock(const xAOD::PFO &shot, const CaloCell_ID *calo_id)
Get cell block with (currently) 2 x 5 cells in correct order for variable calculations.
float ptWindow(const std::vector< std::vector< const CaloCell * > > &shotCells, int windowSize, const ToolHandle< IHadronicCalibrationTool > &caloWeightTool)
pt in a window of (currently) 2 x windowSize cells
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
bool doPi0andShots(const xAOD::TauJet &tau)
Determines whether pi0s and shots should be built for a tau candidate.
PFO_v1 PFO
Definition of the current "pfo version".
Definition PFO.h:17
PFOContainer_v1 PFOContainer
Definition of the current "pfo container version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TauJet_v3 TauJet
Definition of the current "tau version".
CaloClusterContainer_v1 CaloClusterContainer
Define the latest version of the calorimeter cluster container.
ptSort(const TauShotFinder &info)
bool operator()(const CaloCell *c1, const CaloCell *c2)
const TauShotFinder & m_info