CMS 3D CMS Logo

SiStripClusterInfo.h
Go to the documentation of this file.
1 #ifndef SISTRIPCLUSTERIZER_SISTRIPCLUSTERINFO_H
2 #define SISTRIPCLUSTERIZER_SISTRIPCLUSTERINFO_H
3 
7 #include <algorithm>
8 #include <numeric>
9 
10 class SiStripNoises;
11 class SiStripGain;
12 class SiStripQuality;
13 
15 public:
17  const edm::EventSetup& es,
18  const int detid,
19  const std::string& qualityLabel = "");
20 
21  const SiStripCluster* cluster() const { return cluster_ptr; }
22 
23  uint32_t detId() const { return detId_; }
24  uint16_t width() const { return cluster()->amplitudes().size(); }
25  uint16_t firstStrip() const { return cluster()->firstStrip(); }
26  float baryStrip() const { return cluster()->barycenter(); }
27  uint16_t maxStrip() const { return firstStrip() + maxIndex(); }
28  float variance() const;
29 
30  auto stripCharges() const -> decltype(cluster()->amplitudes()) { return cluster()->amplitudes(); }
31  std::vector<float> stripGains() const;
32  std::vector<float> stripNoises() const;
33  std::vector<float> stripNoisesRescaledByGain() const;
34  std::vector<bool> stripQualitiesBad() const;
35 
36  uint16_t charge() const { return std::accumulate(stripCharges().begin(), stripCharges().end(), uint16_t(0)); }
37  uint8_t maxCharge() const { return *std::max_element(stripCharges().begin(), stripCharges().end()); }
38  uint16_t maxIndex() const {
39  return std::max_element(stripCharges().begin(), stripCharges().end()) - stripCharges().begin();
40  }
41  std::pair<uint16_t, uint16_t> chargeLR() const;
42 
43  float noise() const { return calculate_noise(stripNoises()); }
45 
46  float signalOverNoise() const { return charge() / noiseRescaledByGain(); }
47 
48  bool IsAnythingBad() const;
49  bool IsApvBad() const;
50  bool IsFiberBad() const;
51  bool IsModuleBad() const;
52  bool IsModuleUsable() const;
53 
54  std::vector<SiStripCluster> reclusterize(const edm::ParameterSet&) const;
55 
56 private:
57  float calculate_noise(const std::vector<float>&) const;
58 
65  uint32_t detId_;
66 };
67 
68 #endif
uint8_t maxCharge() const
edm::ESHandle< SiStripGain > gainHandle
uint16_t firstStrip() const
edm::ESHandle< SiStripQuality > qualityHandle
SiStripClusterInfo(const SiStripCluster &cluster, const edm::EventSetup &es, const int detid, const std::string &qualityLabel="")
std::vector< float > stripGains() const
bool IsAnythingBad() const
float noise() const
float noiseRescaledByGain() const
std::pair< uint16_t, uint16_t > chargeLR() const
std::vector< bool > stripQualitiesBad() const
float baryStrip() const
uint16_t firstStrip() const
float variance() const
float calculate_noise(const std::vector< float > &) const
uint16_t maxIndex() const
auto stripCharges() const -> decltype(cluster() ->amplitudes())
const SiStripCluster * cluster() const
float signalOverNoise() const
uint16_t charge() const
uint16_t width() const
float barycenter() const
#define end
Definition: vmac.h:39
const SiStripCluster * cluster_ptr
uint16_t maxStrip() const
const edm::EventSetup & es
uint32_t detId() const
bool IsModuleUsable() const
edm::ESHandle< SiStripNoises > noiseHandle
#define begin
Definition: vmac.h:32
std::vector< float > stripNoises() const
std::vector< SiStripCluster > reclusterize(const edm::ParameterSet &) const
const std::vector< uint8_t > & amplitudes() const
std::vector< float > stripNoisesRescaledByGain() const