ATLAS Offline Software
Loading...
Searching...
No Matches
LVL1::EFexEMEnergyWeightedClusterTool Class Reference

#include <EFexEMEnergyWeightedClusterTool.h>

Inheritance diagram for LVL1::EFexEMEnergyWeightedClusterTool:
Collaboration diagram for LVL1::EFexEMEnergyWeightedClusterTool:

Public Member Functions

 EFexEMEnergyWeightedClusterTool (const std::string &type, const std::string &name, const IInterface *parent)
 Name : EFexEMEnergyWeightedClusterTool.cxx PACKAGE : Trigger/TrigT1/TrigT1CaloFexPerf AUTHOR : Denis Oliveira Damazio PURPOSE : emulate the eFex EM algorithm for phase 1 L1Calo (energy weighted clustering)
void findCellsAbove_EMB2_EMEC2 (const CaloConstCellContainer *, const float &thr, std::vector< const CaloCell * > &out) const
void findCellsAround (const CaloConstCellContainer *, const CaloCell *cell, std::vector< const CaloCell * > &out, const float deta, const float dphi) const
void findCellsAround (const CaloConstCellContainer *, const float eta, const float phi, std::vector< const CaloCell * > &out, const float deta, const float dphi) 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.
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 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
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.72 (tile region)
void findCluster (const std::vector< const CaloCell * > &scells, float &etaCluster, float &phiCluster) const
 detect central cluster position
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

void findCellsAbove (const CaloConstCellContainer *, const float &thr, std::vector< const CaloCell * > &out) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

const float m_detaTT {0.125}
 member variables
const float m_dphiTT {0.15}
 dphi for the cluster to TT definition
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 23 of file EFexEMEnergyWeightedClusterTool.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ EFexEMEnergyWeightedClusterTool()

LVL1::EFexEMEnergyWeightedClusterTool::EFexEMEnergyWeightedClusterTool ( const std::string & type,
const std::string & name,
const IInterface * parent )

Name : EFexEMEnergyWeightedClusterTool.cxx PACKAGE : Trigger/TrigT1/TrigT1CaloFexPerf AUTHOR : Denis Oliveira Damazio PURPOSE : emulate the eFex EM algorithm for phase 1 L1Calo (energy weighted clustering)

Definition at line 12 of file EFexEMEnergyWeightedClusterTool.cxx.

13 : AthAlgTool(type, name, parent)
14{}
AthAlgTool()
Default constructor:

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ findCellsAbove()

void LVL1::EFexEMEnergyWeightedClusterTool::findCellsAbove ( const CaloConstCellContainer * scells,
const float & thr,
std::vector< const CaloCell * > & out ) const
private

Definition at line 18 of file EFexEMEnergyWeightedClusterTool.cxx.

19{
20 out.clear();
21 for(auto scell : *scells) {
22 if ( scell->et() < Thr ) continue;
23 out.push_back(scell);
24 }
25}

◆ findCellsAbove_EMB2_EMEC2()

void LVL1::EFexEMEnergyWeightedClusterTool::findCellsAbove_EMB2_EMEC2 ( const CaloConstCellContainer * scells,
const float & thr,
std::vector< const CaloCell * > & out ) const

Definition at line 29 of file EFexEMEnergyWeightedClusterTool.cxx.

30{
31 out.clear();
32 for(auto scell : *scells) {
33 if ( scell->et() < Thr ) continue;
34 if ( scell->caloDDE()->getSampling()==2 || scell->caloDDE()->getSampling()==6 ) {
35 out.push_back(scell);
36 }
37 }
38}

◆ findCellsAround() [1/2]

void LVL1::EFexEMEnergyWeightedClusterTool::findCellsAround ( const CaloConstCellContainer * scells,
const CaloCell * cell,
std::vector< const CaloCell * > & out,
const float deta,
const float dphi ) const

Definition at line 42 of file EFexEMEnergyWeightedClusterTool.cxx.

44{
45 out.clear();
46 if ( !cell ) return;
47 float etacell = cell->eta();
48 float phicell = cell->phi();
49 for(auto scell : *scells) {
50 if ( std::abs( scell->eta() - etacell) > detaSize ) continue;
51 float dphi = std::abs( scell->phi() - phicell);
52 dphi = std::abs( M_PI - dphi );
53 dphi = std::abs( M_PI - dphi );
54 if ( std::abs( dphi ) > dphiSize ) continue;
55 out.push_back(scell);
56 }
57}
#define M_PI

◆ findCellsAround() [2/2]

void LVL1::EFexEMEnergyWeightedClusterTool::findCellsAround ( const CaloConstCellContainer * scells,
const float eta,
const float phi,
std::vector< const CaloCell * > & out,
const float deta,
const float dphi ) const

Definition at line 61 of file EFexEMEnergyWeightedClusterTool.cxx.

63{
64 out.clear();
65 for(auto scell : *scells) {
66 if ( std::abs( scell->eta() - etacell) > detaSize ) continue;
67 float dphi = std::abs( scell->phi() - phicell);
68 dphi = std::abs( M_PI - dphi );
69 dphi = std::abs( M_PI - dphi );
70 if ( std::abs( dphi ) > dphiSize ) continue;
71 out.push_back(scell);
72 }
73}

◆ findCluster()

void LVL1::EFexEMEnergyWeightedClusterTool::findCluster ( const std::vector< const CaloCell * > & scells,
float & etaCluster,
float & phiCluster ) const

detect central cluster position

Definition at line 161 of file EFexEMEnergyWeightedClusterTool.cxx.

162{
163 etaCluster=0.0;
164 phiCluster=0.0;
165 double etaClusterD=0.0;
166 double phiClusterD=0.0;
167 double energyCluster=0.0;
168 bool cross_phi_bound=false;
169 int last_sign=0;
170 for(auto scell : scells){
171 if ( std::abs( scell->phi() ) < 2.7 ) continue;
172 int layer = scell->caloDDE()->getSampling();
173 if ( ( layer != 2 ) && ( layer != 6 ) ) continue;
174 int cell_sign = ( scell->phi() >=0 ? 1 : -1 );
175 if ( ( last_sign!=0 ) && ( last_sign != cell_sign ) ) cross_phi_bound = true;
176 last_sign = cell_sign;
177 }
178 for(auto scell : scells){
179 int layer = scell->caloDDE()->getSampling();
180 if ( ( layer != 2 ) && ( layer != 6 ) ) continue;
181 double scelleta = scell->eta();
182 double scellphi = scell->phi();
183 double scellet = scell->et();
184 etaClusterD+= (scellet * scelleta);
185 if (cross_phi_bound && scellphi < 0 ) scellphi += 2 * M_PI;
186 phiClusterD+= (scellet * scellphi);
187 energyCluster+= (scellet) ;
188 }
189 if ( energyCluster > 0.1 ) {
190 etaClusterD/=energyCluster;
191 phiClusterD/=energyCluster;
192 etaCluster = (float)etaClusterD;
193 phiCluster = (float)phiClusterD;
194 if ( phiCluster > M_PI ) phiCluster-=2*M_PI;
195 } else {
196 etaCluster=-999.0;
197 phiCluster=-999.0;
198 }
199}
@ layer
Definition HitInfo.h:79

◆ findTTsAround()

void LVL1::EFexEMEnergyWeightedClusterTool::findTTsAround ( const xAOD::TriggerTowerContainer * scells,
const float eta,
const float phi,
std::vector< const xAOD::TriggerTower * > & out ) const

finds TTs around a seed cell.

These TTs will be part of the cluster. This helps to cover the part related to TileCall

Definition at line 77 of file EFexEMEnergyWeightedClusterTool.cxx.

79{
80 out.clear();
81 for(auto scell : *scells) {
82 if ( std::abs( scell->eta() - etacell) > m_detaTT ) continue;
83 float dphi = std::abs( scell->phi() - phicell);
84 dphi = std::abs( M_PI - dphi );
85 dphi = std::abs( M_PI - dphi );
86 if ( std::abs( dphi ) > m_dphiTT ) continue;
87 out.push_back(scell);
88 }
89}
const float m_dphiTT
dphi for the cluster to TT definition

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ isCellEmMaximum()

bool LVL1::EFexEMEnergyWeightedClusterTool::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.

This is to make sure we have a real local EM maximum

Definition at line 93 of file EFexEMEnergyWeightedClusterTool.cxx.

94{
95 if ( !cell ) return false;
96 int samp = cell->caloDDE()->getSampling();
97 if ( (samp >= 8) && (samp!=21) ) return false; // include FCAL0 EM
98 float cellpt = 1.0001*cell->et(); //make sure you don't get thecell itself
99 for(auto scell : scells){
100 int samp1 = scell->caloDDE()->getSampling();
101 if ( ( samp1 >= 8 ) && (samp1!=21) ) continue;
102 if ( scell->ID() == cell->ID() ) continue;
103 if ( scell->et() > cellpt ) return false;
104 }
105 return true;
106}

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sumEmCells()

float LVL1::EFexEMEnergyWeightedClusterTool::sumEmCells ( const std::vector< const CaloCell * > & scells) const

sum all cells from the vector that are in the EM calorimeter part

Definition at line 110 of file EFexEMEnergyWeightedClusterTool.cxx.

111{
112 float totalSum = 0.0;
113 for(auto scell : scells) {
114 int samp1 = scell->caloDDE()->getSampling();
115 if ( (samp1<8) || (samp1==21) ) totalSum+= scell->energy();
116 }
117 return totalSum;
118}

◆ sumEmCells2nd()

float LVL1::EFexEMEnergyWeightedClusterTool::sumEmCells2nd ( const std::vector< const CaloCell * > & scells) const

sum all cells from the vector that are in the EM calorimeter part (only 2nd layer)

Definition at line 122 of file EFexEMEnergyWeightedClusterTool.cxx.

123{
124 float totalSum = 0.0;
125 for(auto scell : scells) {
126 if ( (scell->caloDDE()->getSampling()==2) ||(scell->caloDDE()->getSampling()==6) ) {
127 totalSum+= scell->energy();
128 }
129 }
130 return totalSum;
131}

◆ sumHadCells()

float LVL1::EFexEMEnergyWeightedClusterTool::sumHadCells ( const std::vector< const CaloCell * > & scells) const

sum all cells from the vector that are in the HAD calorimeter part

Definition at line 135 of file EFexEMEnergyWeightedClusterTool.cxx.

136{
137 float totalSum = 0.0;
138 for(auto scell : scells){
139 if ( (scell->caloDDE()->getSampling() <8) || ( scell->caloDDE()->getSampling()>=22) ) continue;
140 //totalSum+= (scell->et())*TMath::CosH(scell->eta());
141 totalSum+= (scell->energy());
142 }
143 return totalSum;
144}

◆ sumHadTTs()

float LVL1::EFexEMEnergyWeightedClusterTool::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.72 (tile region)

Definition at line 148 of file EFexEMEnergyWeightedClusterTool.cxx.

149{
150 float totalSum = 0.0;
151 for(auto scell : scells){
152 if ( std::abs( scell->eta() ) > 1.5 ) continue;
153 if ( scell->sampling() == 0 ) continue;
154 totalSum+= (scell->pt())*TMath::CosH(scell->eta());
155 }
156 return totalSum * 1e3; // express in MeV
157}

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_detaTT

const float LVL1::EFexEMEnergyWeightedClusterTool::m_detaTT {0.125}
private

member variables

deta for the cluster to TT definition

Definition at line 67 of file EFexEMEnergyWeightedClusterTool.h.

67{0.125};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_dphiTT

const float LVL1::EFexEMEnergyWeightedClusterTool::m_dphiTT {0.15}
private

dphi for the cluster to TT definition

Definition at line 68 of file EFexEMEnergyWeightedClusterTool.h.

68{0.15};

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: