ATLAS Offline Software
Loading...
Searching...
No Matches
EFexEMEnergyWeightedClusterTool.h
Go to the documentation of this file.
1// Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
2
9
10#ifndef LVL1__EFEXEMENERGYWEIGHTEDCLUSTERTOOL
11#define LVL1__EFEXEMENERGYWEIGHTEDCLUSTERTOOL
12
16
17#include <string>
18#include <vector>
19
20
21namespace LVL1
22{
24 {
25 public:
26 EFexEMEnergyWeightedClusterTool(const std::string& type, const std::string& name, const IInterface* parent);
27
28 // find the cells at EMB2 and EMEC2 above a given threshold in an event to be used as seeds to build electron/tau candidates
29 void findCellsAbove_EMB2_EMEC2(const CaloConstCellContainer*, const float& thr, std::vector<const CaloCell*>& out) const;
30
31 // finds cells around a seed cell. These cells will be part of the cluster
32 void findCellsAround(const CaloConstCellContainer*, const CaloCell* cell, std::vector<const CaloCell*>& out,
33 const float deta, const float dphi) const;
34
35 void findCellsAround(const CaloConstCellContainer*, const float eta, const float phi, std::vector<const CaloCell*>& out,
36 const float deta, const float dphi) const ;
37
38
40 void findTTsAround(const xAOD::TriggerTowerContainer*, const float eta, const float phi, std::vector<const xAOD::TriggerTower*>& out) const ;
41
43 bool isCellEmMaximum(const std::vector<const CaloCell*>& scells, const CaloCell* cell) const ;
45 float sumEmCells(const std::vector<const CaloCell*>& scells) const ;
47 float sumEmCells2nd(const std::vector<const CaloCell*>& scells) const ;
49 float sumHadCells(const std::vector<const CaloCell*>& scells) const ;
51 float sumHadTTs(const std::vector<const xAOD::TriggerTower*>& scells) const ;
52
54 void findCluster(const std::vector<const CaloCell*>& scells, float &etaCluster, float &phiCluster) const;
55
56
57
58 private:
59
60 // find the cells above a given threshold in an event to be used as seeds to build electron candidates
61 void findCellsAbove(const CaloConstCellContainer*, const float& thr, std::vector<const CaloCell*>& out) const;
62
63
67 const float m_detaTT{0.125};
68 const float m_dphiTT{0.15};
69 };
70}
71
72#endif
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
CaloCellContainer that can accept const cell pointers.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Data object for each calorimeter readout cell.
Definition CaloCell.h:57
CaloCellContainer that can accept const cell pointers.
bool isCellEmMaximum(const std::vector< const CaloCell * > &scells, const CaloCell *cell) const
checks if a give (seed) cell is the highest in a vector of cells.
float sumEmCells(const std::vector< const CaloCell * > &scells) const
sum all cells from the vector that are in the EM calorimeter part
float sumHadTTs(const std::vector< const xAOD::TriggerTower * > &scells) const
sum all TTs from the vector that are in the HAD calorimeter part, but only for |eta|<1....
void findCellsAround(const CaloConstCellContainer *, const CaloCell *cell, std::vector< const CaloCell * > &out, const float deta, const float dphi) const
void findCluster(const std::vector< const CaloCell * > &scells, float &etaCluster, float &phiCluster) const
detect central cluster position
void findCellsAbove(const CaloConstCellContainer *, const float &thr, std::vector< const CaloCell * > &out) const
const float m_dphiTT
dphi for the cluster to TT definition
float sumEmCells2nd(const std::vector< const CaloCell * > &scells) const
sum all cells from the vector that are in the EM calorimeter part (only 2nd layer)
float sumHadCells(const std::vector< const CaloCell * > &scells) const
sum all cells from the vector that are in the HAD calorimeter part
EFexEMEnergyWeightedClusterTool(const std::string &type, const std::string &name, const IInterface *parent)
Name : EFexEMEnergyWeightedClusterTool.cxx PACKAGE : Trigger/TrigT1/TrigT1CaloFexPerf AUTHOR : Denis ...
void findCellsAbove_EMB2_EMEC2(const CaloConstCellContainer *, const float &thr, std::vector< const CaloCell * > &out) const
void findTTsAround(const xAOD::TriggerTowerContainer *, const float eta, const float phi, std::vector< const xAOD::TriggerTower * > &out) const
finds TTs around a seed cell.
eFexTowerBuilder creates xAOD::eFexTowerContainer from supercells (LATOME) and triggerTowers (TREX) i...
TriggerTowerContainer_v2 TriggerTowerContainer
Define the latest version of the TriggerTower container.