Loading [MathJax]/extensions/tex2jax.js
ATLAS Offline Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
HIEventShapeIndex.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #ifndef HIEVENTUTILS_HIEVENTSHAPEINDEX_H
6 #define HIEVENTUTILS_HIEVENTSHAPEINDEX_H
7 
8 #include "HIEventUtils/HIEventDefs.h" //HI::BinningScheme
9 #include "xAODHIEvent/HIEventShapeContainer.h" //typedef
10 #include <memory>
11 #include <map>
12 #include <vector>
13 class TH2;
14 
15 
17 
18 public:
21 
22  //initialize indexing given object or scheme
23  //returns total number of bins
24 
25  unsigned int setBinning(const TH2* h2, bool asMask);
26  unsigned int setBinning(HI::BinningScheme scheme);
27 
28  void initializeEventShapeContainer(std::unique_ptr<xAOD::HIEventShapeContainer>& shape_container, unsigned int num_harmonics) const;
29  //can associate ptr to shape container w/ binning
30  //not required
31 
32  //lookup functions
33  unsigned int getIndex(float eta, int layer) const;
34  unsigned int getIndexFromBin(unsigned int ebin, int layer) const;
35  unsigned int getEtaBin(float eta, int layer) const;
36  unsigned int getNumBins() const {return m_size;};
37  xAOD::HIEventShape* getShape(float eta, int layer, xAOD::HIEventShapeContainer* shape_container) const;
38  const xAOD::HIEventShape* getShape(float eta, int layer, const xAOD::HIEventShapeContainer* shape_container) const;
39  xAOD::HIEventShape* getShape(float eta, int layer, std::unique_ptr<xAOD::HIEventShapeContainer>& shape_container) const;
40  //HIEventShapeIndex& operator=(const HIEventShapeIndex& in);
41  std::string print() const;
42 private:
43  struct range_index_t{
44  float eta_min{};
45  float eta_max{};
46  unsigned int index{};
47 
48  range_index_t(float emin, float emax, unsigned int ii) : eta_min(emin), eta_max(emax),index(ii) {};
49  range_index_t()=default;
50  bool operator() (float eta) const
51  {
52  if (eta > this->eta_min && eta < this->eta_max) return true;
53  if(eta==this->eta_min) return true;
54  return false;
55  };
56  bool operator< (const range_index_t& rhs) const {return this->eta_min < rhs.eta_min;};
57  };
58  std::map<int,std::vector<range_index_t> > m_edges;
59 
60  unsigned int m_size{};
61 
62  unsigned int getIndex_Internal(float eta, int layer, bool etaIndex) const;
63  std::map<int,std::vector<range_index_t> >::const_iterator getLayer(int layer) const;
64 
65  inline float roundToTenth(float d) const {return std::floor(d)+std::floor((d-std::floor(d))*10.0+0.5)/10.0;};
66 
67 
68 
69 };
70 
71 #endif
HIEventShapeIndex::getEtaBin
unsigned int getEtaBin(float eta, int layer) const
Definition: HIEventShapeIndex.cxx:151
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
index
Definition: index.py:1
HIEventShapeIndex::getIndexFromBin
unsigned int getIndexFromBin(unsigned int ebin, int layer) const
Definition: HIEventShapeIndex.cxx:106
hist_file_dump.d
d
Definition: hist_file_dump.py:143
HIEventShapeIndex::roundToTenth
float roundToTenth(float d) const
Definition: HIEventShapeIndex.h:65
HIEventShapeIndex::HIEventShapeIndex
HIEventShapeIndex()
Definition: HIEventShapeIndex.cxx:15
HIEventShapeIndex::getShape
xAOD::HIEventShape * getShape(float eta, int layer, xAOD::HIEventShapeContainer *shape_container) const
Definition: HIEventShapeIndex.cxx:156
HIEventDefs.h
HIEventShapeIndex::~HIEventShapeIndex
~HIEventShapeIndex()
Definition: HIEventShapeIndex.cxx:19
HIEventShapeIndex::range_index_t::operator()
bool operator()(float eta) const
Definition: HIEventShapeIndex.h:50
HIEventShapeIndex::initializeEventShapeContainer
void initializeEventShapeContainer(std::unique_ptr< xAOD::HIEventShapeContainer > &shape_container, unsigned int num_harmonics) const
Definition: HIEventShapeIndex.cxx:207
HIEventShapeIndex::getLayer
std::map< int, std::vector< range_index_t > >::const_iterator getLayer(int layer) const
Definition: HIEventShapeIndex.cxx:194
xAOD::HIEventShape_v2
Interface class for the HI reconstruction EDM.
Definition: HIEventShape_v2.h:32
HIEventShapeIndex::range_index_t::range_index_t
range_index_t()=default
HIEventShapeIndex::range_index_t::eta_min
float eta_min
Definition: HIEventShapeIndex.h:44
HIEventShapeIndex::print
std::string print() const
Definition: HIEventShapeIndex.cxx:247
TRT::Hit::layer
@ layer
Definition: HitInfo.h:79
HIEventShapeIndex::setBinning
unsigned int setBinning(const TH2 *h2, bool asMask)
Definition: HIEventShapeIndex.cxx:24
HIEventShapeIndex::getNumBins
unsigned int getNumBins() const
Definition: HIEventShapeIndex.h:36
HIEventShapeContainer.h
HIEventShapeIndex::range_index_t
Definition: HIEventShapeIndex.h:43
HIEventShapeIndex::range_index_t::range_index_t
range_index_t(float emin, float emax, unsigned int ii)
Definition: HIEventShapeIndex.h:48
HIEventShapeIndex::getIndex_Internal
unsigned int getIndex_Internal(float eta, int layer, bool etaIndex) const
Definition: HIEventShapeIndex.cxx:111
DataVector
Derived DataVector<T>.
Definition: DataVector.h:794
HIEventShapeIndex
Definition: HIEventShapeIndex.h:16
HIEventShapeIndex::m_size
unsigned int m_size
Definition: HIEventShapeIndex.h:60
HIEventShapeIndex::getIndex
unsigned int getIndex(float eta, int layer) const
Definition: HIEventShapeIndex.cxx:146
HIEventShapeIndex::range_index_t::operator<
bool operator<(const range_index_t &rhs) const
Definition: HIEventShapeIndex.h:56
HIEventShapeIndex::range_index_t::eta_max
float eta_max
Definition: HIEventShapeIndex.h:45
HI::BinningScheme
BinningScheme
Definition: HIEventDefs.h:16
HIEventShapeIndex::m_edges
std::map< int, std::vector< range_index_t > > m_edges
Definition: HIEventShapeIndex.h:58