File indexing completed on 2024-04-06 12:29:40
0001
0002
0003
0004 #ifndef SimDataFormats_CaloAnalysis_MtdSimCluster_h
0005 #define SimDataFormats_CaloAnalysis_MtdSimCluster_h
0006
0007 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0008 #include "SimDataFormats/CaloAnalysis/interface/SimCluster.h"
0009 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0010 #include <vector>
0011
0012 class MtdSimCluster : public SimCluster {
0013 friend std::ostream &operator<<(std::ostream &s, MtdSimCluster const &tp);
0014
0015 public:
0016 MtdSimCluster();
0017 MtdSimCluster(const SimTrack &simtrk);
0018 MtdSimCluster(EncodedEventId eventID, uint32_t particleID);
0019
0020
0021 ~MtdSimCluster();
0022
0023
0024 void addHitTime(float time) {
0025 times_.emplace_back(time);
0026 ++nsimhits_;
0027 }
0028
0029
0030 void addHitAndFraction(uint64_t hit, float fraction) {
0031 mtdHits_.emplace_back(hit);
0032 fractions_.emplace_back(fraction);
0033 }
0034
0035
0036 void addHitPosition(LocalPoint pos) { positions_.emplace_back(pos); }
0037
0038
0039 std::vector<std::pair<uint64_t, float>> hits_and_fractions() const {
0040 assert(mtdHits_.size() == fractions_.size());
0041 std::vector<std::pair<uint64_t, float>> result;
0042 result.reserve(mtdHits_.size());
0043 for (size_t i = 0; i < mtdHits_.size(); ++i) {
0044 result.emplace_back(mtdHits_[i], fractions_[i]);
0045 }
0046 return result;
0047 }
0048
0049
0050 std::vector<std::pair<uint64_t, float>> hits_and_energies() const {
0051 assert(mtdHits_.size() == energies_.size());
0052 std::vector<std::pair<uint64_t, float>> result;
0053 result.reserve(mtdHits_.size());
0054 for (size_t i = 0; i < mtdHits_.size(); ++i) {
0055 result.emplace_back(mtdHits_[i], energies_[i]);
0056 }
0057 return result;
0058 }
0059
0060
0061 void clearHitsAndFractions() {
0062 std::vector<uint64_t>().swap(mtdHits_);
0063 std::vector<float>().swap(fractions_);
0064 }
0065
0066
0067 std::vector<std::pair<uint64_t, float>> hits_and_times() const {
0068 assert(mtdHits_.size() == times_.size());
0069 std::vector<std::pair<uint64_t, float>> result;
0070 result.reserve(mtdHits_.size());
0071 for (size_t i = 0; i < mtdHits_.size(); ++i) {
0072 result.emplace_back(mtdHits_[i], times_[i]);
0073 }
0074 return result;
0075 }
0076
0077
0078 std::vector<std::pair<uint64_t, LocalPoint>> hits_and_positions() const {
0079 assert(mtdHits_.size() == times_.size());
0080 std::vector<std::pair<uint64_t, LocalPoint>> result;
0081 result.reserve(mtdHits_.size());
0082 for (size_t i = 0; i < mtdHits_.size(); ++i) {
0083 result.emplace_back(mtdHits_[i], positions_[i]);
0084 }
0085 return result;
0086 }
0087
0088
0089 std::vector<std::pair<uint32_t, std::pair<uint8_t, uint8_t>>> detIds_and_rows() const {
0090 std::vector<std::pair<uint32_t, std::pair<uint8_t, uint8_t>>> result;
0091 result.reserve(mtdHits_.size());
0092 for (size_t i = 0; i < mtdHits_.size(); ++i) {
0093 result.emplace_back(
0094 mtdHits_[i] >> 32,
0095 std::pair<uint8_t, uint8_t>(static_cast<uint8_t>(mtdHits_[i] >> 16), static_cast<uint8_t>(mtdHits_[i])));
0096 }
0097 return result;
0098 }
0099
0100
0101 void clearHitsTime() { std::vector<float>().swap(times_); }
0102
0103
0104 void clearHitsPosition() { std::vector<LocalPoint>().swap(positions_); }
0105
0106 void clear() {
0107 clearHitsAndFractions();
0108 clearHitsEnergy();
0109 clearHitsTime();
0110 clearHitsPosition();
0111 }
0112
0113
0114 void addSimHit(const PSimHit &hit) {
0115 simhit_energy_ += hit.energyLoss();
0116 ++nsimhits_;
0117 }
0118
0119 void setTrackIdOffset(unsigned int offset) { idOffset_ = offset; }
0120
0121 unsigned int trackIdOffset() const { return idOffset_; }
0122
0123 protected:
0124 std::vector<uint64_t> mtdHits_;
0125 std::vector<float> times_;
0126 std::vector<LocalPoint> positions_;
0127 unsigned int idOffset_{0};
0128 };
0129
0130 #endif