CMS 3D CMS Logo

PixelClusterCounts.h
Go to the documentation of this file.
1 #ifndef DataFormats_Luminosity_PixelClusterCounts_h
2 #define DataFormats_Luminosity_PixelClusterCounts_h
3 
11 #include <algorithm>
12 #include <string>
13 #include <sstream>
14 #include <iostream>
15 #include <vector>
16 
19 
20 namespace reco {
22  public:
24 
25  void increment(int mD, unsigned int bxID, int count) {
26  size_t modIndex = std::distance(m_ModID.begin(), std::find(m_ModID.begin(), m_ModID.end(), mD));
27  if (modIndex == m_ModID.size()) {
28  m_ModID.push_back(mD);
29  m_counts.resize(m_counts.size() + LumiConstants::numBX, 0);
30  }
31  m_counts.at(LumiConstants::numBX * modIndex + bxID - 1) += count;
32  }
33 
34  void incrementRoc(int rD, unsigned int bxID, int count) {
35  size_t rocIndex = std::distance(m_RocID.begin(), std::find(m_RocID.begin(), m_RocID.end(), rD));
36  if (rocIndex == m_RocID.size()) {
37  m_RocID.push_back(rD);
38  m_countsRoc.resize(m_countsRoc.size() + LumiConstants::numBX, 0);
39  }
40  m_countsRoc.at(LumiConstants::numBX * rocIndex + bxID - 1) += count;
41  }
42 
43  void eventCounter(unsigned int bxID) { m_events.at(bxID - 1)++; }
44 
45  void add(reco::PixelClusterCountsInEvent const& pccInEvent) {
46  std::vector<int> const& countsInEvent = pccInEvent.counts();
47  std::vector<int> const& rocCountsInEvent = pccInEvent.countsRoc();
48  std::vector<int> const& modIDInEvent = pccInEvent.modID();
49  std::vector<int> const& rocIDInEvent = pccInEvent.rocID();
50  int bxIDInEvent = pccInEvent.bxID();
51  for (unsigned int i = 0; i < modIDInEvent.size(); i++) {
52  increment(modIDInEvent[i], bxIDInEvent, countsInEvent.at(i));
53  }
54  for (unsigned int i = 0; i < rocIDInEvent.size(); i++) {
55  incrementRoc(rocIDInEvent[i], bxIDInEvent, rocCountsInEvent.at(i));
56  }
57  }
58 
59  void merge(reco::PixelClusterCounts const& pcc) {
60  std::vector<int> const& counts = pcc.readCounts();
61  std::vector<int> const& countsRoc = pcc.readRocCounts();
62  std::vector<int> const& modIDs = pcc.readModID();
63  std::vector<int> const& rocIDs = pcc.readRocID();
64  std::vector<int> const& events = pcc.readEvents();
65  for (unsigned int i = 0; i < modIDs.size(); i++) {
66  for (unsigned int bxID = 0; bxID < LumiConstants::numBX; ++bxID) {
67  increment(modIDs[i], bxID + 1, counts.at(i * LumiConstants::numBX + bxID));
68  }
69  }
70  for (unsigned int i = 0; i < rocIDs.size(); i++) {
71  for (unsigned int bxID = 0; bxID < LumiConstants::numBX; ++bxID) {
72  incrementRoc(rocIDs[i], bxID + 1, countsRoc.at(i * LumiConstants::numBX + bxID));
73  }
74  }
75  for (unsigned int i = 0; i < LumiConstants::numBX; ++i) {
76  m_events[i] += events[i];
77  }
78  }
79 
80  void reset() {
81  m_counts.clear();
82  m_countsRoc.clear();
83  m_ModID.clear();
84  m_RocID.clear();
85  m_events.clear();
86  m_events.resize(LumiConstants::numBX, 0);
87  }
88 
89  std::vector<int> const& readCounts() const { return (m_counts); }
90  std::vector<int> const& readRocCounts() const { return (m_countsRoc); }
91  std::vector<int> const& readEvents() const { return (m_events); }
92  std::vector<int> const& readModID() const { return (m_ModID); }
93  std::vector<int> const& readRocID() const { return (m_RocID); }
94 
95  private:
96  std::vector<int> m_counts;
97  std::vector<int> m_countsRoc;
98  std::vector<int> m_events;
99  std::vector<int> m_ModID;
100  std::vector<int> m_RocID;
101  };
102 
103 } // namespace reco
104 #endif
std::vector< int > const & rocID() const
std::vector< int > const & readModID() const
unsigned int const & bxID() const
static const unsigned int numBX
Definition: LumiConstants.h:8
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void incrementRoc(int rD, unsigned int bxID, int count)
void add(reco::PixelClusterCountsInEvent const &pccInEvent)
std::vector< int > const & readEvents() const
std::vector< int > m_RocID
std::vector< int > m_events
std::vector< int > const & readRocID() const
void merge(reco::PixelClusterCounts const &pcc)
std::vector< int > const & readCounts() const
void increment(int mD, unsigned int bxID, int count)
std::vector< int > const & readRocCounts() const
std::vector< int > const & modID() const
fixed size matrix
std::vector< int > m_counts
void eventCounter(unsigned int bxID)
std::vector< int > const & countsRoc() const
std::vector< int > m_countsRoc
std::vector< int > m_ModID
int events
std::vector< int > const & counts() const