Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-03 04:17:41

0001 #ifndef DataFormats_Luminosity_PixelClusterCounts_h
0002 #define DataFormats_Luminosity_PixelClusterCounts_h
0003 /** \class reco::PixelClusterCounts
0004  *  
0005  * Reconstructed PixelClusterCounts object that will contain the moduleID, BX, and counts.
0006  *
0007  * \authors: Sam Higginbotham shigginb@cern.ch and Chris Palmer capalmer@cern.ch
0008  * 
0009  *
0010  */
0011 #include <algorithm>
0012 #include <string>
0013 #include <sstream>
0014 #include <iostream>
0015 #include <vector>
0016 
0017 #include "DataFormats/Luminosity/interface/LumiConstants.h"
0018 #include "DataFormats/Luminosity/interface/PixelClusterCountsInEvent.h"
0019 
0020 namespace reco {
0021   class PixelClusterCounts {
0022   public:
0023     PixelClusterCounts() : m_events(LumiConstants::numBX) {}
0024 
0025     void increment(int mD, unsigned int bxID, int count) {
0026       size_t modIndex = std::distance(m_ModID.begin(), std::find(m_ModID.begin(), m_ModID.end(), mD));
0027       if (modIndex == m_ModID.size()) {
0028         m_ModID.push_back(mD);
0029         m_counts.resize(m_counts.size() + LumiConstants::numBX, 0);
0030       }
0031       m_counts.at(LumiConstants::numBX * modIndex + bxID - 1) += count;
0032     }
0033 
0034     void incrementRoc(int rD, int count) {
0035       size_t rocIndex = std::distance(m_RocID.begin(), std::find(m_RocID.begin(), m_RocID.end(), rD));
0036 
0037       if (rocIndex == m_RocID.size()) {
0038         m_RocID.push_back(rD);
0039         m_countsRoc.push_back(0);
0040       }
0041 
0042       m_countsRoc.at(rocIndex) += count;
0043     }
0044 
0045     void eventCounter(unsigned int bxID) { m_events.at(bxID - 1)++; }
0046 
0047     void add(reco::PixelClusterCountsInEvent const& pccInEvent) {
0048       std::vector<int> const& countsInEvent = pccInEvent.counts();
0049       std::vector<int> const& rocCountsInEvent = pccInEvent.countsRoc();
0050       std::vector<int> const& modIDInEvent = pccInEvent.modID();
0051       std::vector<int> const& rocIDInEvent = pccInEvent.rocID();
0052       int bxIDInEvent = pccInEvent.bxID();
0053 
0054       for (unsigned int i = 0; i < modIDInEvent.size(); i++) {
0055         increment(modIDInEvent[i], bxIDInEvent, countsInEvent.at(i));
0056       }
0057       for (unsigned int i = 0; i < rocIDInEvent.size(); i++) {
0058         incrementRoc(rocIDInEvent[i], rocCountsInEvent.at(i));
0059       }
0060     }
0061 
0062     void merge(reco::PixelClusterCounts const& pcc) {
0063       std::vector<int> const& counts = pcc.readCounts();
0064       std::vector<int> const& countsRoc = pcc.readRocCounts();
0065       std::vector<int> const& modIDs = pcc.readModID();
0066       std::vector<int> const& rocIDs = pcc.readRocID();
0067       std::vector<int> const& events = pcc.readEvents();
0068 
0069       for (unsigned int i = 0; i < modIDs.size(); i++) {
0070         for (unsigned int bxID = 0; bxID < LumiConstants::numBX; ++bxID) {
0071           increment(modIDs[i], bxID + 1, counts.at(i * LumiConstants::numBX + bxID));
0072         }
0073       }
0074       for (unsigned int i = 0; i < rocIDs.size(); i++) {
0075         incrementRoc(rocIDs[i], countsRoc.at(i));
0076       }
0077       for (unsigned int i = 0; i < LumiConstants::numBX; ++i) {
0078         m_events[i] += events[i];
0079       }
0080     }
0081 
0082     void reset() {
0083       m_counts.clear();
0084       m_countsRoc.clear();
0085       m_ModID.clear();
0086       m_RocID.clear();
0087       m_events.clear();
0088       m_events.resize(LumiConstants::numBX, 0);
0089     }
0090 
0091     std::vector<int> const& readCounts() const { return (m_counts); }
0092     std::vector<int> const& readRocCounts() const { return (m_countsRoc); }
0093     std::vector<int> const& readEvents() const { return (m_events); }
0094     std::vector<int> const& readModID() const { return (m_ModID); }
0095     std::vector<int> const& readRocID() const { return (m_RocID); }
0096 
0097   private:
0098     std::vector<int> m_counts;
0099     std::vector<int> m_countsRoc;
0100     std::vector<int> m_events;
0101     std::vector<int> m_ModID;
0102     std::vector<int> m_RocID;
0103   };
0104 
0105 }  // namespace reco
0106 #endif