ATLAS Offline Software
Loading...
Searching...
No Matches
xAOD::IParticlesLookUpTable< T > Class Template Reference

2D look up table for iParticles More...

#include <IParticlesLookUpTable.h>

Collaboration diagram for xAOD::IParticlesLookUpTable< T >:

Public Member Functions

 IParticlesLookUpTable (unsigned int nbins=50, float minPt=1e-3)
 constructor taking the desired binsize
void init (const DataVector< T > &particles)
 initialize the look up table with an iParticle container
void clear ()
 clear the look up table
template<class O>
bool iParticlesInCone (float eta, float phi, float dr, O &output) const
 collect IParticles in a given cone
bool isInitialized () const
 return whether table is initialized

Private Member Functions

float phiInRange (float phi) const
 hepler function to ensure phi is within +-Pi
int phiIndex (float phi) const
 calculate phi index for a given phi
void addEntry (int i, std::vector< const T * > &output) const
 add an entry into a vector of pointers
void addEntry (int i, std::vector< ElementLink< DataVector< T > > > &output) const
 add an entry into a vector of ElementLinks

Private Attributes

const DataVector< T > * m_container
int m_nphiBins
float m_phiBinSize
 number of bins
float m_minPt
 bin size
std::vector< std::vector< int > > m_phiBinnedLookUpTable
 cut on minimum Pt for the considered particles

Static Private Attributes

static constexpr float m_2PI = 2*M_PI
 define 2*Pi

Detailed Description

template<class T>
class xAOD::IParticlesLookUpTable< T >

2D look up table for iParticles

Definition at line 17 of file IParticlesLookUpTable.h.

Constructor & Destructor Documentation

◆ IParticlesLookUpTable()

template<class T>
xAOD::IParticlesLookUpTable< T >::IParticlesLookUpTable ( unsigned int nbins = 50,
float minPt = 1e-3 )
inline

constructor taking the desired binsize

Definition at line 20 of file IParticlesLookUpTable.h.

20 :m_container(0),
25 {}
2D look up table for iParticles
const DataVector< T > * m_container
static constexpr float m_2PI
define 2*Pi
std::vector< std::vector< int > > m_phiBinnedLookUpTable
cut on minimum Pt for the considered particles

Member Function Documentation

◆ addEntry() [1/2]

template<class T>
void xAOD::IParticlesLookUpTable< T >::addEntry ( int i,
std::vector< const T * > & output ) const
private

add an entry into a vector of pointers

Definition at line 153 of file IParticlesLookUpTable.h.

153 {
154 output.push_back( (*m_container)[i]);
155 }

◆ addEntry() [2/2]

template<class T>
void xAOD::IParticlesLookUpTable< T >::addEntry ( int i,
std::vector< ElementLink< DataVector< T > > > & output ) const
private

add an entry into a vector of ElementLinks

Definition at line 158 of file IParticlesLookUpTable.h.

158 {
160 }

◆ clear()

template<class T>
void xAOD::IParticlesLookUpTable< T >::clear ( )

clear the look up table

Definition at line 73 of file IParticlesLookUpTable.h.

73 {
75 }

◆ init()

template<class T>
void xAOD::IParticlesLookUpTable< T >::init ( const DataVector< T > & particles)

initialize the look up table with an iParticle container

resize hash table

loop over iparticles and copy them into the look-up struct use hashing for phi look-up and sorting for eta look-up

Definition at line 78 of file IParticlesLookUpTable.h.

78 {
84 unsigned int size = particles.size();
85 for( unsigned int i=0; i<size;++i ){
86 if(particles[i]->pt()>m_minPt){//sanity check (mainly due to Truth)
87 int index = phiIndex(particles[i]->phi());
88 m_phiBinnedLookUpTable[index].push_back(i);
89 }
90 }
91 //std::cout << " initializing lookup " << m_container->size() << std::endl;
92 for( auto& vec : m_phiBinnedLookUpTable ) {
93 if( vec.empty() ) continue;
94 std::stable_sort(vec.begin(),vec.end(),[&](int i, int j) { return (*m_container)[i]->eta() < (*m_container)[j]->eta(); });
95 // std::cout << " new phi slice " << vec.size() << " eta " ;
96 // for( auto index : vec ){
97 // std::cout << " " << (*m_container)[index]->eta();
98 // }
99 // std::cout << std::endl;
100 }
101 }
int phiIndex(float phi) const
calculate phi index for a given phi
void stable_sort(DataModel_detail::iterator< DVL > beg, DataModel_detail::iterator< DVL > end)
Specialization of stable_sort for DataVector/List.
setRcore setEtHad setFside pt

◆ iParticlesInCone()

template<class T>
template<class O>
bool xAOD::IParticlesLookUpTable< T >::iParticlesInCone ( float eta,
float phi,
float dr,
O & output ) const

collect IParticles in a given cone

check if initialized

comparison functions for upper and lower bound

get phi hash ranges

Definition at line 104 of file IParticlesLookUpTable.h.

104 {
105
107 if( m_phiBinnedLookUpTable.empty() ) return false;
108
110 auto compEta1 = [&](int i,double val ) { return (*m_container)[i]->eta() < val; };
111 auto compEta2 = [&](double val,int i ) { return val < (*m_container)[i]->eta(); };
112
116 //std::cout << " eta " << eta << " phi " << phi << " phi ranges " << indexMin << " " << indexMax << std::endl;
117 // special treatment for boundary region
119 if( indexMin > indexMax ) {
120 ranges.push_back( std::make_pair(0,indexMax) );
121 ranges.push_back( std::make_pair(indexMin,m_nphiBins-1) );
122 }else{
124 }
125
126 float dr2Cut = dr*dr;
127
128 // loop over ranges
129 for( auto& range : ranges ){
130 indexMin = range.first;
131 indexMax = range.second;
132 for( ; indexMin <= indexMax; ++indexMin ){
133 // get iterators for iparticles to be included
135 auto it_min = std::lower_bound (tps.begin(),tps.end(),eta-dr,compEta1 );
137 //std::cout << " new range, entries in bounds " << std::distance(it_max,it_min) << std::endl;
138 // add iparticles in cone
139 for( ;it_min!=it_max;++it_min ){
140 const T* entry = (*m_container)[*it_min];
141 float deta = eta- entry->eta();
142 float dphi = phiInRange(phi-entry->phi());
143 float dr2 = deta*deta + dphi*dphi;
144 //std::cout << " new entry: " << *it_min << " eta,phi " << entry->eta() << " phi " << entry->phi() << " dr " << sqrt(dr2) << std::endl;
145 if( dr2 < dr2Cut ) addEntry( *it_min, output );
146 }
147 }
148 }
149 return true;
150 }
float phiInRange(float phi) const
hepler function to ensure phi is within +-Pi
void addEntry(int i, std::vector< const T * > &output) const
add an entry into a vector of pointers

◆ isInitialized()

template<class T>
bool xAOD::IParticlesLookUpTable< T >::isInitialized ( ) const

return whether table is initialized

Definition at line 69 of file IParticlesLookUpTable.h.

69{ return !m_phiBinnedLookUpTable.empty(); }

◆ phiIndex()

template<class T>
int xAOD::IParticlesLookUpTable< T >::phiIndex ( float phi) const
inlineprivate

calculate phi index for a given phi

Definition at line 52 of file IParticlesLookUpTable.h.

52{ return (phi + M_PI)/m_phiBinSize; };

◆ phiInRange()

template<class T>
float xAOD::IParticlesLookUpTable< T >::phiInRange ( float phi) const
inlineprivate

hepler function to ensure phi is within +-Pi

Definition at line 45 of file IParticlesLookUpTable.h.

45 {
46 while (phi >= M_PI) phi -= m_2PI;
47 while (phi < -M_PI) phi += m_2PI;
48 return phi;
49 }

Member Data Documentation

◆ m_2PI

template<class T>
float xAOD::IParticlesLookUpTable< T >::m_2PI = 2*M_PI
staticconstexprprivate

define 2*Pi

Definition at line 42 of file IParticlesLookUpTable.h.

◆ m_container

template<class T>
const DataVector<T>* xAOD::IParticlesLookUpTable< T >::m_container
private

Definition at line 61 of file IParticlesLookUpTable.h.

◆ m_minPt

template<class T>
float xAOD::IParticlesLookUpTable< T >::m_minPt
private

bin size

Definition at line 64 of file IParticlesLookUpTable.h.

◆ m_nphiBins

template<class T>
int xAOD::IParticlesLookUpTable< T >::m_nphiBins
private

Definition at line 62 of file IParticlesLookUpTable.h.

◆ m_phiBinnedLookUpTable

template<class T>
std::vector< std::vector< int > > xAOD::IParticlesLookUpTable< T >::m_phiBinnedLookUpTable
private

cut on minimum Pt for the considered particles

Definition at line 65 of file IParticlesLookUpTable.h.

◆ m_phiBinSize

template<class T>
float xAOD::IParticlesLookUpTable< T >::m_phiBinSize
private

number of bins

Definition at line 63 of file IParticlesLookUpTable.h.


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