ATLAS Offline Software
Loading...
Searching...
No Matches
TauShotFinder.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef TAURECTOOLS_TAUSHOTFINDER_H
6#define TAURECTOOLS_TAUSHOTFINDER_H
7
9
14#include "GaudiKernel/ToolHandle.h"
16
27
28class CaloCell_ID;
29
31
32public:
33
35
36 TauShotFinder(const std::string& name);
37 virtual ~TauShotFinder() = default;
38
39 virtual StatusCode initialize() override;
40 virtual StatusCode executeShotFinder(xAOD::TauJet& pTau, xAOD::CaloClusterContainer& tauShotCaloClusContainer, xAOD::PFOContainer& tauShotPFOContainer) const override;
41
42private:
43
44 struct ptSort
45 {
46 ptSort( const TauShotFinder& info );
48 bool operator()( const CaloCell* c1, const CaloCell* c2 );
49 };
50
54 StatusCode selectCells(const xAOD::TauJet& tau,
55 std::vector<const CaloCell*>& cells) const;
56
63 StatusCode selectSeedCells(const xAOD::TauJet& tau,
64 const CaloCellContainer& cellContainer,
65 std::vector<const CaloCell*>& seedCells) const;
66
68 bool isPhiNeighbour(IdentifierHash cell1Hash, IdentifierHash cell2Hash) const;
69
71 const CaloCell* getPhiNeighbour(const CaloCell& seedCell, const std::vector<const CaloCell*>& seedCells) const;
72
74 std::vector<const CaloCell*> getEtaNeighbours(const CaloCell& cell, const CaloCellContainer& cellContainer, int maxDepth) const;
75
77 void addEtaNeighbours(const CaloCell& cell,
78 const CaloCellContainer& cellContainer,
79 std::vector<const CaloCell*>& cells,
80 int depth,
81 int maxDepth,
82 bool next) const;
83
89 const CaloCell* phiNeighCell,
90 const CaloCellContainer& cellContainer,
91 xAOD::CaloClusterContainer* clusterContainer) const;
92
94 int getEtaBin(float eta) const;
95
97 int getNPhotons(float eta, float energy) const;
98
99 Gaudi::Property<int> m_nCellsInEta {this, "NCellsInEta"};
100 Gaudi::Property<std::vector<float>> m_minPtCut {this, "MinPtCut"};
101 Gaudi::Property<std::vector<float>> m_doubleShotCut {this, "AutoDoubleShotCut"};
102 Gaudi::Property<bool> m_removeElectronCells {this, "RemoveElectronCells", false};
103 Gaudi::Property<float> m_dRThreshold{this, "DRThreshold", 0.4, "dR(cell,tau) to select cells"};
104 Gaudi::Property<float> m_energyThreshold{this, "EnergyThreshold", 100., "energy threshold (in MeV) to select cells"};
105
106 SG::ReadHandleKey<CaloCellContainer> m_caloCellInputContainer{this,"Key_caloCellInputContainer", "AllCalo", "input calo cell container key"};
107 ToolHandle<IHadronicCalibrationTool> m_caloWeightTool {this, "CaloWeightTool", "H1WeightToolCSC12Generic"};
108 SG::ReadHandleKey<xAOD::CaloClusterContainer> m_removedClusterInputContainer {this,"Key_RemovedClusterInputContainer", "", "input removed cluster key"};
109
111 const CaloCell_ID* m_calo_id = nullptr;
112
113};
114
115#endif // TAURECTOOLS_TAUSHOTFINDER_H
Scalar eta() const
pseudorapidity method
Definition of CaloDetDescrManager.
Container class for CaloCell.
Helper class for offline cell identifiers.
Definition CaloCell_ID.h:34
Data object for each calorimeter readout cell.
Definition CaloCell.h:57
The base class for all tau tools.
This is a "hash" representation of an Identifier.
Property holding a SG store/key/clid from which a ReadHandle is made.
TauRecToolBase(const std::string &name)
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
ASG_TOOL_CLASS2(TauShotFinder, TauRecToolBase, ITauToolBase)
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.
virtual ~TauShotFinder()=default
Gaudi::Property< float > m_energyThreshold
Gaudi::Property< float > m_dRThreshold
std::string depth
tag string for intendation
Definition fastadd.cxx:46
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