ATLAS Offline Software
Loading...
Searching...
No Matches
eflowRec::EtaPhiLUT Class Reference

2D look up table for eflowRecClusters More...

#include <EtaPhiLUT.h>

Collaboration diagram for eflowRec::EtaPhiLUT:

Public Member Functions

 EtaPhiLUT (unsigned int nbins=50)
 constructor taking the desired binsize
void fill (eflowRecClusterContainer &clustersin)
std::vector< eflowRecCluster * > clustersInCone (float eta, float phi, float dr) const
 collect eflowRecClusters in a given cone

Private Attributes

unsigned int m_nphiBins
float m_phiBinSize
 number of bins
std::vector< std::vector< eflowRecCluster * > > m_phiBinnedLookUpTable
 bin size

Detailed Description

2D look up table for eflowRecClusters

Definition at line 26 of file EtaPhiLUT.h.

Constructor & Destructor Documentation

◆ EtaPhiLUT()

eflowRec::EtaPhiLUT::EtaPhiLUT ( unsigned int nbins = 50)

constructor taking the desired binsize

Definition at line 25 of file EtaPhiLUT.cxx.

25 :
26 m_nphiBins(nbins),
29 {}
std::vector< std::vector< eflowRecCluster * > > m_phiBinnedLookUpTable
bin size
Definition EtaPhiLUT.h:40
unsigned int m_nphiBins
Definition EtaPhiLUT.h:38
float m_phiBinSize
number of bins
Definition EtaPhiLUT.h:39
static constexpr float TWOPI
define 2*Pi
Definition EtaPhiLUT.cxx:13

Member Function Documentation

◆ clustersInCone()

std::vector< eflowRecCluster * > eflowRec::EtaPhiLUT::clustersInCone ( float eta,
float phi,
float dr ) const

collect eflowRecClusters in a given cone

Definition at line 43 of file EtaPhiLUT.cxx.

43 {
44 std::vector<eflowRecCluster*> result;
45
46 // Indices corresponding to phi range
47 unsigned int iPhiMin = phiIndex( phiInRange(phi-dr), m_phiBinSize );
48 unsigned int iPhiMax = phiIndex( phiInRange(phi+dr), m_phiBinSize );
49
50 // Extract index ranges to iterate over
51 std::vector< std::pair<int,int> > iPhiRanges;
52 if( iPhiMin < iPhiMax ) {
53 iPhiRanges.emplace_back(iPhiMin,iPhiMax );
54 } else { // special treatment for phi-wrapping
55 iPhiRanges.emplace_back(0,iPhiMax );
56 iPhiRanges.emplace_back(iPhiMin,m_nphiBins-1 );
57 }
58
59 float dr2Cut = dr*dr;
60
61 // loop over ranges
62 for( auto& range : iPhiRanges ){
63 unsigned int indexMin = range.first;
64 unsigned int indexMax = range.second;
65 for( ; indexMin <= indexMax; ++indexMin ){
66 const std::vector<eflowRecCluster*>& phiClusters = m_phiBinnedLookUpTable[indexMin];
67
68 // Get iterators for clusters in relevant eta range
69 auto it_min = std::lower_bound (phiClusters.begin(), phiClusters.end(), eta-dr, [] (const eflowRecCluster* cl, float etaval) {return cl->getCluster()->eta()<etaval;} );
70 auto it_max = std::upper_bound (it_min, phiClusters.end(), eta+dr, [] (float etaval, const eflowRecCluster* cl) {return etaval<cl->getCluster()->eta();} );
71
72 // Apply deltaR cut
73 for( ;it_min!=it_max;++it_min ){
74 const xAOD::CaloCluster& xcluster = *(*it_min)->getCluster();
75 float deta = eta - xcluster.eta();
76 float dphi = phiInRange(phi - xcluster.phi());
77 float dr2 = deta*deta + dphi*dphi;
78 if( dr2 < dr2Cut ) result.push_back(*it_min);
79 }
80 }
81 }
82 return result;
83 }
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
unsigned int phiIndex(float phi, float binsize)
calculate phi index for a given phi
Definition EtaPhiLUT.cxx:23
float phiInRange(float phi)
hepler function to ensure phi is within +-Pi
Definition EtaPhiLUT.cxx:16
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.

◆ fill()

void eflowRec::EtaPhiLUT::fill ( eflowRecClusterContainer & clustersin)

Definition at line 31 of file EtaPhiLUT.cxx.

32 {
33 for(eflowRecCluster* cluster : clustersin) {
34 int index = phiIndex(cluster->getCluster()->phi(), m_phiBinSize);
35 m_phiBinnedLookUpTable[index].push_back(cluster);
36 }
37
38 for( auto& vec : m_phiBinnedLookUpTable ) {
39 std::stable_sort(vec.begin(),vec.end(),[](const eflowRecCluster* c1, const eflowRecCluster* c2) { return c1->getCluster()->eta() < c2->getCluster()->eta(); });
40 }
41 }
std::vector< size_t > vec
str index
Definition DeMoScan.py:362
void stable_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > end)
Specialization of stable_sort for DataVector/List.

Member Data Documentation

◆ m_nphiBins

unsigned int eflowRec::EtaPhiLUT::m_nphiBins
private

Definition at line 38 of file EtaPhiLUT.h.

◆ m_phiBinnedLookUpTable

std::vector< std::vector<eflowRecCluster*> > eflowRec::EtaPhiLUT::m_phiBinnedLookUpTable
private

bin size

Definition at line 40 of file EtaPhiLUT.h.

◆ m_phiBinSize

float eflowRec::EtaPhiLUT::m_phiBinSize
private

number of bins

Definition at line 39 of file EtaPhiLUT.h.


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