ATLAS Offline Software
Loading...
Searching...
No Matches
HIEventShapeIndex Class Reference

#include <HIEventShapeIndex.h>

Collaboration diagram for HIEventShapeIndex:

Classes

struct  range_index_t

Public Member Functions

 HIEventShapeIndex ()
 ~HIEventShapeIndex ()
unsigned int setBinning (const TH2 *h2, bool asMask)
unsigned int setBinning (HI::BinningScheme scheme)
void initializeEventShapeContainer (std::unique_ptr< xAOD::HIEventShapeContainer > &shape_container, unsigned int num_harmonics) const
unsigned int getIndex (float eta, int layer) const
unsigned int getIndexFromBin (unsigned int ebin, int layer) const
unsigned int getEtaBin (float eta, int layer) const
unsigned int getNumBins () const
xAOD::HIEventShapegetShape (float eta, int layer, xAOD::HIEventShapeContainer *shape_container) const
const xAOD::HIEventShapegetShape (float eta, int layer, const xAOD::HIEventShapeContainer *shape_container) const
xAOD::HIEventShapegetShape (float eta, int layer, std::unique_ptr< xAOD::HIEventShapeContainer > &shape_container) const
std::string print () const

Private Member Functions

unsigned int getIndex_Internal (float eta, int layer, bool etaIndex) const
std::map< int, std::vector< range_index_t > >::const_iterator getLayer (int layer) const
float roundToTenth (float d) const

Private Attributes

std::map< int, std::vector< range_index_t > > m_edges
unsigned int m_size {}

Detailed Description

Definition at line 16 of file HIEventShapeIndex.h.

Constructor & Destructor Documentation

◆ HIEventShapeIndex()

HIEventShapeIndex::HIEventShapeIndex ( )

Definition at line 15 of file HIEventShapeIndex.cxx.

15 : m_size(0)
16{
17}

◆ ~HIEventShapeIndex()

HIEventShapeIndex::~HIEventShapeIndex ( )

Definition at line 19 of file HIEventShapeIndex.cxx.

20{
21}

Member Function Documentation

◆ getEtaBin()

unsigned int HIEventShapeIndex::getEtaBin ( float eta,
int layer ) const

Definition at line 151 of file HIEventShapeIndex.cxx.

152{
153 return getIndex_Internal(eta,layer,true);
154}
Scalar eta() const
pseudorapidity method
unsigned int getIndex_Internal(float eta, int layer, bool etaIndex) const

◆ getIndex()

unsigned int HIEventShapeIndex::getIndex ( float eta,
int layer ) const

Definition at line 146 of file HIEventShapeIndex.cxx.

147{
148 return getIndex_Internal(eta,layer,false);
149}

◆ getIndex_Internal()

unsigned int HIEventShapeIndex::getIndex_Internal ( float eta,
int layer,
bool etaIndex ) const
private

Definition at line 111 of file HIEventShapeIndex.cxx.

112{
113 auto mItr=getLayer(layer);
114 const std::vector<range_index_t>& vec=mItr->second;
115 unsigned int vsize=vec.size();
116
117 if(vsize==0)
118 {
119 std::cerr << "Layer " << layer << " has no range in calorimeter" << std::endl;
120 throw std::range_error("Bad binning request");
121 return 0;
122 }
123
124 float eta_c=eta;
125 if(m_edges.size()>1)
126 {
127 float abs_eta=std::abs(eta);
128 float rmin=HICaloRange::getRange().getRangeMin(layer);
129 float rmax=HICaloRange::getRange().getRangeMax(layer);
130
131 if(abs_eta < rmin) eta_c=rmin+std::copysign(1e-2,eta);
132 else if(abs_eta > rmax) eta_c=rmax-std::copysign(1e-2,eta);
133 }
134
135 unsigned int pIndex=0;
136 for(; pIndex < vsize-1; pIndex++)
137 {
138 const range_index_t& p=vec.at(pIndex);
139 if(p(eta_c)) break;
140 }
141
142 return (etaIndex) ? pIndex : vec.at(pIndex).index;
143}
std::vector< size_t > vec
float getRangeMax(int layer) const
Definition HICaloRange.h:20
static const HICaloRange & getRange()
float getRangeMin(int layer) const
Definition HICaloRange.h:19
std::map< int, std::vector< range_index_t > >::const_iterator getLayer(int layer) const
std::map< int, std::vector< range_index_t > > m_edges
@ layer
Definition HitInfo.h:79

◆ getIndexFromBin()

unsigned int HIEventShapeIndex::getIndexFromBin ( unsigned int ebin,
int layer ) const

Definition at line 106 of file HIEventShapeIndex.cxx.

107{
108 return getLayer(layer)->second.at(ebin).index;
109}

◆ getLayer()

std::map< int, std::vector< HIEventShapeIndex::range_index_t > >::const_iterator HIEventShapeIndex::getLayer ( int layer) const
private

Definition at line 194 of file HIEventShapeIndex.cxx.

195{
196 auto mItr=m_edges.find(layer);
197
198 if(mItr==m_edges.end())
199 {
200 //ATH_MSG_ERROR("Layer " << layer << " has no binning specification");
201 std::cerr << "Layer " << layer << " has no binning specification" << std::endl;
202 throw std::range_error("Bad binning request");
203 }
204 return mItr;
205}

◆ getNumBins()

unsigned int HIEventShapeIndex::getNumBins ( ) const
inline

Definition at line 36 of file HIEventShapeIndex.h.

36{return m_size;};

◆ getShape() [1/3]

const xAOD::HIEventShape * HIEventShapeIndex::getShape ( float eta,
int layer,
const xAOD::HIEventShapeContainer * shape_container ) const

Definition at line 170 of file HIEventShapeIndex.cxx.

171{
172 if(shape_container==nullptr)
173 {
174 //ATH_MSG_WARNING("No default HIEventShapeContainer to query. Returning nullptr_t.");
175 return nullptr;
176 }
177 unsigned int bin=getIndex(eta,layer);
178 //correct for round-off errors in cell eta values, just need to make sure they end up in the right bin
179 if( bin >= shape_container->size() ) bin=getIndex(eta-std::copysign(1e-2,eta),layer);
180
181 return shape_container->at(bin);
182}
const T * at(size_type n) const
Access an element, as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
unsigned int getIndex(float eta, int layer) const

◆ getShape() [2/3]

xAOD::HIEventShape * HIEventShapeIndex::getShape ( float eta,
int layer,
std::unique_ptr< xAOD::HIEventShapeContainer > & shape_container ) const

Definition at line 185 of file HIEventShapeIndex.cxx.

186{
187 unsigned int bin=getIndex(eta,layer);
188 //correct for round-off errors in cell eta values, just need to make sure they end up in the right bin
189 if( bin >= shape_container->size() ) bin=getIndex(eta-std::copysign(1e-2,eta),layer);
190
191 return shape_container->at(bin);
192}

◆ getShape() [3/3]

xAOD::HIEventShape * HIEventShapeIndex::getShape ( float eta,
int layer,
xAOD::HIEventShapeContainer * shape_container ) const

Definition at line 156 of file HIEventShapeIndex.cxx.

157{
158 if(shape_container==nullptr)
159 {
160 //ATH_MSG_WARNING("No default HIEventShapeContainer to query. Returning nullptr_t.");
161 return nullptr;
162 }
163 unsigned int bin=getIndex(eta,layer);
164 //correct for round-off errors in cell eta values, just need to make sure they end up in the right bin
165 if( bin >= shape_container->size() ) bin=getIndex(eta-std::copysign(1e-2,eta),layer);
166
167 return shape_container->at(bin);
168}

◆ initializeEventShapeContainer()

void HIEventShapeIndex::initializeEventShapeContainer ( std::unique_ptr< xAOD::HIEventShapeContainer > & shape_container,
unsigned int num_harmonics ) const

Definition at line 207 of file HIEventShapeIndex.cxx.

208{
209 if(getNumBins()==0)
210 {
211 //ATH_MSG_WARNING("Cannot initialize HIEventShapeContainer. Index has not been initialized.");
212 return;
213
214 }
215 bool isEmpty=(shape_container->size()==0);
216 if(isEmpty) shape_container->resize(getNumBins());
217 else if(shape_container->size()!=getNumBins())
218 {
219 //ATH_MSG_WARNING("Cannot initialize HIEventShapeContainer. Number of elements incompatible w/ requested binning.");
220 return;
221
222 }
223
224 for(const auto& pp : m_edges)
225 {
226 for(auto ri : pp.second)
227 {
229 if(isEmpty)
230 {
232 //shape_container->push_back(slice);
233 shape_container->operator[](ri.index)=slice;
234 }
235 else slice=shape_container->at(ri.index);
236
237 slice->setEtaMin(ri.eta_min);
238 slice->setEtaMax(ri.eta_max);
239 slice->setLayer(pp.first);
240 slice->etCos().assign(num_harmonics,0);
241 slice->etSin().assign(num_harmonics,0);
242 }
243 }
244}
unsigned int getNumBins() const
HIEventShape_v2 HIEventShape
Definition of the latest event info version.

◆ print()

std::string HIEventShapeIndex::print ( ) const

Definition at line 247 of file HIEventShapeIndex.cxx.

248{
249 std::stringstream ss;
250 ss.precision(3);
251 for(const auto& pp : m_edges)
252 {
253 for(const auto & ri : pp.second)
254 {
255 ss << std::setw(10) << pp.first
256 << std::setw(10) << ri.eta_min
257 << std::setw(10) << ri.eta_max
258 << std::setw(10) << ri.index
259 << std::endl;
260 }
261 }
262 return ss.str();
263}
static Double_t ss

◆ roundToTenth()

float HIEventShapeIndex::roundToTenth ( float d) const
inlineprivate

Definition at line 65 of file HIEventShapeIndex.h.

65{return std::floor(d)+std::floor((d-std::floor(d))*10.0+0.5)/10.0;};

◆ setBinning() [1/2]

unsigned int HIEventShapeIndex::setBinning ( const TH2 * h2,
bool asMask )

Definition at line 24 of file HIEventShapeIndex.cxx.

25{
26 if(h2==nullptr)
27 {
28 //ATH_MSG_ERROR("Bad histogram pointer");
29 return 0;
30 }
31
32 unsigned int count=0;
33 std::set<unsigned int> used_indices;
34 for(int ybin=1; ybin<h2->GetNbinsY(); ybin++)
35 {
36 auto mItr=m_edges.emplace_hint(m_edges.end(),ybin-1,std::vector<range_index_t>());
37 for(int xbin=1; xbin<h2->GetNbinsX(); xbin++)
38 {
39 if(h2->GetBinContent(xbin,ybin) < 0) continue;
40 unsigned int v=h2->GetBinContent(xbin,ybin);
41
42 if(asMask) v=count;
43 else
44 {
45 auto sItr=used_indices.find(v);
46 if(sItr==used_indices.end()) used_indices.insert(v);
47 else
48 {
49 //ATH_MSG_ERROR("Could not set binning. initialization histogram has duplicate indices " << v);
50 //m_edges.clear();
51 //return 0;
52 continue;
53 }
54 }
55 mItr->second.push_back(range_index_t(h2->GetXaxis()->GetBinLowEdge(xbin),h2->GetXaxis()->GetBinUpEdge(xbin),v));
56 count++;
57 }
58 }
60 return count;
61}
int count(std::string s, const std::string &regx)
count how many occurances of a regx are in a string
Definition hcg.cxx:146

◆ setBinning() [2/2]

unsigned int HIEventShapeIndex::setBinning ( HI::BinningScheme scheme)

Definition at line 63 of file HIEventShapeIndex.cxx.

64{
65 unsigned int index=0;
66 if(scheme==HI::BinningScheme::TOWER)
67 {
68 for(unsigned int layer=0; layer < static_cast<int>(HI::TowerBins::numLayers()); layer++)
69 {
70 auto mItr=m_edges.emplace_hint(m_edges.end(),layer,std::vector<range_index_t>(HI::TowerBins::numEtaBins(),range_index_t()));
71 for(unsigned int eb=0; eb<HI::TowerBins::numEtaBins(); eb++,index++)
72 {
73
74 mItr->second[eb].eta_min=HI::TowerBins::getBinLowEdgeEta(eb);
75 mItr->second[eb].eta_max=HI::TowerBins::getBinUpEdgeEta(eb);
76 mItr->second[eb].index=index;
77 }
78 }
79 }
80 else if(scheme==HI::BinningScheme::COMPACT)
81 {
83 for(int layer=0; layer < static_cast<int>(HI::TowerBins::numLayers()); layer++)
84 {
85 float emax=HICaloRange::getRange().getRangeMax(layer);
86 float emin=HICaloRange::getRange().getRangeMin(layer);
87 unsigned int num_eta_bins= static_cast<unsigned int>(std::round(std::abs((emax-emin)/deta)));
88 auto mItr=m_edges.emplace_hint(m_edges.end(),layer,std::vector<range_index_t>(2*num_eta_bins,range_index_t()));
89 for(unsigned int k=0; k<num_eta_bins; k++)
90 {
91 mItr->second[k].eta_min=roundToTenth(-emax+k*deta);
92 mItr->second[k].eta_max=roundToTenth(-emax+k*deta+deta);
93 mItr->second[k].index=roundToTenth(k+index);
94 mItr->second[k+num_eta_bins].eta_min=roundToTenth(emin+k*deta);
95 mItr->second[k+num_eta_bins].eta_max=roundToTenth(emin+k*deta+deta);
96 mItr->second[k+num_eta_bins].index=roundToTenth(k+num_eta_bins+index);
97 }
98 index+=2*num_eta_bins;
99 }
100 }
102 return index;
103
104}
float roundToTenth(float d) const
str index
Definition DeMoScan.py:362
float getBinUpEdgeEta(unsigned int eb)
Definition HIEventDefs.h:39
float getBinLowEdgeEta(unsigned int eb)
Definition HIEventDefs.h:38
constexpr unsigned int numLayers()
Definition HIEventDefs.h:23
constexpr float getBinSizeEta()
Definition HIEventDefs.h:31
constexpr unsigned int numEtaBins()
Definition HIEventDefs.h:19
@ COMPACT
Definition HIEventDefs.h:16
@ TOWER
Definition HIEventDefs.h:16

Member Data Documentation

◆ m_edges

std::map<int,std::vector<range_index_t> > HIEventShapeIndex::m_edges
private

Definition at line 58 of file HIEventShapeIndex.h.

◆ m_size

unsigned int HIEventShapeIndex::m_size {}
private

Definition at line 60 of file HIEventShapeIndex.h.

60{};

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