CMS 3D CMS Logo

SiStripClusterInfo.h
Go to the documentation of this file.
1 #ifndef SISTRIPCLUSTERIZER_SISTRIPCLUSTERINFO_H
2 #define SISTRIPCLUSTERIZER_SISTRIPCLUSTERINFO_H
3 
14 
15 #include <algorithm>
16 #include <cstdint>
17 #include <numeric>
18 #include <utility>
19 #include <vector>
20 
22 public:
24 
25  void initEvent(const edm::EventSetup& iSetup);
26  void setCluster(const SiStripCluster& cluster, int detId);
27 
28  const SiStripCluster* cluster() const { return cluster_ptr; }
29 
30  uint32_t detId() const { return detId_; }
31  uint16_t width() const { return cluster()->amplitudes().size(); }
32  uint16_t firstStrip() const { return cluster()->firstStrip(); }
33  float baryStrip() const { return cluster()->barycenter(); }
34  uint16_t maxStrip() const { return firstStrip() + maxIndex(); }
35  float variance() const;
36 
37  auto stripCharges() const -> decltype(cluster()->amplitudes()) { return cluster()->amplitudes(); }
38  std::vector<float> stripGains() const;
39  std::vector<float> stripNoises() const;
40  std::vector<float> stripNoisesRescaledByGain() const;
41  std::vector<bool> stripQualitiesBad() const;
42 
43  uint16_t charge() const { return std::accumulate(stripCharges().begin(), stripCharges().end(), uint16_t(0)); }
44  uint8_t maxCharge() const { return *std::max_element(stripCharges().begin(), stripCharges().end()); }
45  uint16_t maxIndex() const {
46  return std::max_element(stripCharges().begin(), stripCharges().end()) - stripCharges().begin();
47  }
48  std::pair<uint16_t, uint16_t> chargeLR() const;
49 
50  float noise() const { return calculate_noise(stripNoises()); }
52 
53  float signalOverNoise() const { return charge() / noiseRescaledByGain(); }
54 
55  bool IsAnythingBad() const;
56  bool IsApvBad() const;
57  bool IsFiberBad() const;
58  bool IsModuleBad() const;
59  bool IsModuleUsable() const;
60 
61  const SiStripGain* siStripGain() const { return siStripGain_; }
62  const SiStripQuality* siStripQuality() const { return siStripQuality_; }
63 
64 private:
65  float calculate_noise(const std::vector<float>&) const;
66 
67  const SiStripCluster* cluster_ptr = nullptr;
68 
72 
73  const SiStripNoises* siStripNoises_ = nullptr;
74  const SiStripGain* siStripGain_ = nullptr;
75  const SiStripQuality* siStripQuality_ = nullptr;
76 
77  uint32_t detId_ = 0;
78 };
79 
80 #endif
uint16_t firstStrip() const
void setCluster(const SiStripCluster &cluster, int detId)
float calculate_noise(const std::vector< float > &) const
bool IsAnythingBad() const
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > siStripQualityToken_
std::vector< float > stripGains() const
const SiStripGain * siStripGain() const
auto stripCharges() const -> decltype(cluster() ->amplitudes())
float signalOverNoise() const
uint32_t detId() const
SiStripCluster const & amplitudes() const
auto size() const
uint8_t maxCharge() const
const SiStripCluster * cluster() const
SiStripClusterInfo(edm::ConsumesCollector &&, const std::string &qualityLabel="")
const SiStripQuality * siStripQuality() const
const SiStripNoises * siStripNoises_
uint16_t maxStrip() const
std::vector< bool > stripQualitiesBad() const
std::vector< float > stripNoises() const
edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > siStripNoisesToken_
float noiseRescaledByGain() const
const SiStripCluster * cluster_ptr
std::vector< float > stripNoisesRescaledByGain() const
uint16_t firstStrip() const
uint16_t maxIndex() const
const SiStripGain * siStripGain_
std::pair< uint16_t, uint16_t > chargeLR() const
const SiStripQuality * siStripQuality_
uint16_t width() const
edm::ESGetToken< SiStripGain, SiStripGainRcd > siStripGainToken_
float barycenter() const
uint16_t charge() const
void initEvent(const edm::EventSetup &iSetup)
float baryStrip() const