#include <HIEventShapeIndex.h>
|
| | 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::HIEventShape * | getShape (float eta, int layer, xAOD::HIEventShapeContainer *shape_container) const |
| const xAOD::HIEventShape * | getShape (float eta, int layer, const xAOD::HIEventShapeContainer *shape_container) const |
| xAOD::HIEventShape * | getShape (float eta, int layer, std::unique_ptr< xAOD::HIEventShapeContainer > &shape_container) const |
| std::string | print () const |
Definition at line 16 of file HIEventShapeIndex.h.
◆ HIEventShapeIndex()
| HIEventShapeIndex::HIEventShapeIndex |
( |
| ) |
|
◆ ~HIEventShapeIndex()
| HIEventShapeIndex::~HIEventShapeIndex |
( |
| ) |
|
◆ getEtaBin()
| unsigned int HIEventShapeIndex::getEtaBin |
( |
float | eta, |
|
|
int | layer ) const |
Definition at line 151 of file HIEventShapeIndex.cxx.
152{
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 |
◆ getIndex_Internal()
| unsigned int HIEventShapeIndex::getIndex_Internal |
( |
float | eta, |
|
|
int | layer, |
|
|
bool | etaIndex ) const |
|
private |
Definition at line 111 of file HIEventShapeIndex.cxx.
112{
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
126 {
127 float abs_eta=std::abs(
eta);
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 {
140 }
141
142 return (etaIndex) ? pIndex :
vec.at(pIndex).index;
143}
std::vector< size_t > vec
float getRangeMax(int layer) const
static const HICaloRange & getRange()
float getRangeMin(int layer) const
std::map< int, std::vector< range_index_t > >::const_iterator getLayer(int layer) const
std::map< int, std::vector< range_index_t > > m_edges
◆ getIndexFromBin()
| unsigned int HIEventShapeIndex::getIndexFromBin |
( |
unsigned int | ebin, |
|
|
int | layer ) const |
◆ getLayer()
Definition at line 194 of file HIEventShapeIndex.cxx.
195{
197
199 {
200
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 |
◆ getShape() [1/3]
Definition at line 170 of file HIEventShapeIndex.cxx.
171{
172 if(shape_container==nullptr)
173 {
174
175 return nullptr;
176 }
178
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]
◆ getShape() [3/3]
Definition at line 156 of file HIEventShapeIndex.cxx.
157{
158 if(shape_container==nullptr)
159 {
160
161 return nullptr;
162 }
164
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{
210 {
211
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
220 return;
221
222 }
223
225 {
226 for(auto ri : pp.second)
227 {
229 if(isEmpty)
230 {
232
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;
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 }
263}
◆ 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
29 return 0;
30 }
31
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
43 else
44 {
45 auto sItr=used_indices.find(v);
46 if(sItr==used_indices.end()) used_indices.insert(v);
47 else
48 {
49
50
51
52 continue;
53 }
54 }
55 mItr->second.push_back(
range_index_t(
h2->GetXaxis()->GetBinLowEdge(xbin),
h2->GetXaxis()->GetBinUpEdge(xbin),v));
57 }
58 }
61}
int count(std::string s, const std::string ®x)
count how many occurances of a regx are in a string
◆ setBinning() [2/2]
Definition at line 63 of file HIEventShapeIndex.cxx.
64{
67 {
69 {
72 {
73
76 mItr->second[eb].index=
index;
77 }
78 }
79 }
81 {
84 {
87 unsigned int num_eta_bins= static_cast<unsigned int>(std::round(std::abs((emax-emin)/deta)));
89 for(
unsigned int k=0;
k<num_eta_bins;
k++)
90 {
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 }
103
104}
float roundToTenth(float d) const
float getBinUpEdgeEta(unsigned int eb)
float getBinLowEdgeEta(unsigned int eb)
constexpr unsigned int numLayers()
constexpr float getBinSizeEta()
constexpr unsigned int numEtaBins()
◆ m_edges
| std::map<int,std::vector<range_index_t> > HIEventShapeIndex::m_edges |
|
private |
◆ m_size
| unsigned int HIEventShapeIndex::m_size {} |
|
private |
The documentation for this class was generated from the following files: