CMS 3D CMS Logo

SiStripClusterInfo.cc
Go to the documentation of this file.
2 
3 #include <cmath>
4 
6  : siStripNoisesToken_{iC.esConsumes<SiStripNoises, SiStripNoisesRcd>()},
7  siStripGainToken_{iC.esConsumes<SiStripGain, SiStripGainRcd>()},
8  siStripQualityToken_{iC.esConsumes<SiStripQuality, SiStripQualityRcd>(edm::ESInputTag("", qualityLabel))} {}
9 
14 }
15 
16 void SiStripClusterInfo::setCluster(const SiStripCluster& cluster, int detId) {
18  detId_ = detId;
19 }
20 
21 std::pair<uint16_t, uint16_t> SiStripClusterInfo::chargeLR() const {
22  std::vector<uint8_t>::const_iterator begin(stripCharges().begin()), end(stripCharges().end()), max;
23  max = max_element(begin, end);
24  return std::make_pair(accumulate(begin, max, uint16_t(0)), accumulate(max + 1, end, uint16_t(0)));
25 }
26 
28  float q(0), x1(0), x2(0);
29  for (auto begin(stripCharges().begin()), end(stripCharges().end()), it(begin); it != end; ++it) {
30  unsigned i = it - begin;
31  q += (*it);
32  x1 += (*it) * (i + 0.5);
33  x2 += (*it) * (i * i + i + 1. / 3);
34  }
35  return (x2 - x1 * x1 / q) / q;
36 }
37 
41 
42  std::vector<float> results;
43  results.reserve(width());
44  for (size_t i = 0, e = width(); i < e; i++) {
45  results.push_back(siStripNoises_->getNoise(firstStrip() + i, detNoiseRange) /
46  siStripGain_->getStripGain(firstStrip() + i, detGainRange));
47  }
48  return results;
49 }
50 
51 std::vector<float> SiStripClusterInfo::stripNoises() const {
53 
54  std::vector<float> noises;
55  noises.reserve(width());
56  for (size_t i = 0; i < width(); i++) {
57  noises.push_back(siStripNoises_->getNoise(firstStrip() + i, detNoiseRange));
58  }
59  return noises;
60 }
61 
62 std::vector<float> SiStripClusterInfo::stripGains() const {
64 
65  std::vector<float> gains;
66  gains.reserve(width());
67  for (size_t i = 0; i < width(); i++) {
68  gains.push_back(siStripGain_->getStripGain(firstStrip() + i, detGainRange));
69  }
70  return gains;
71 }
72 
73 std::vector<bool> SiStripClusterInfo::stripQualitiesBad() const {
74  std::vector<bool> isBad;
75  isBad.reserve(width());
76  for (int i = 0; i < width(); i++) {
77  isBad.push_back(siStripQuality_->IsStripBad(detId_, firstStrip() + i));
78  }
79  return isBad;
80 }
81 
82 float SiStripClusterInfo::calculate_noise(const std::vector<float>& noise) const {
83  float noiseSumInQuadrature = 0;
84  int numberStripsOverThreshold = 0;
85  for (int i = 0; i < width(); i++) {
86  if (stripCharges()[i] != 0) {
87  noiseSumInQuadrature += noise.at(i) * noise.at(i);
88  numberStripsOverThreshold++;
89  }
90  }
91  return std::sqrt(noiseSumInQuadrature / numberStripsOverThreshold);
92 }
93 
95  std::vector<bool> stripBad = stripQualitiesBad();
96  return IsApvBad() || IsFiberBad() || IsModuleBad() ||
97  accumulate(stripBad.begin(), stripBad.end(), false, std::logical_or<bool>());
98 }
99 
101  return siStripQuality_->IsApvBad(detId_, firstStrip() / 128) ||
103 }
104 
106  return siStripQuality_->IsFiberBad(detId_, firstStrip() / 256) ||
108 }
109 
111 
bool IsApvBad(uint32_t detid, short apvNb) 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
bool IsModuleBad(uint32_t detid) const
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:75
auto stripCharges() const -> decltype(cluster() ->amplitudes())
static float getNoise(uint16_t strip, const Range &range)
Definition: SiStripNoises.h:73
uint32_t detId() const
static float getStripGain(const uint16_t &strip, const SiStripApvGain::Range &range)
Definition: SiStripGain.h:77
T sqrt(T t)
Definition: SSEVec.h:19
const SiStripCluster * cluster() const
std::pair< ContainerIterator, ContainerIterator > Range
SiStripClusterInfo(edm::ConsumesCollector &&, const std::string &qualityLabel="")
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const SiStripNoises * siStripNoises_
std::vector< bool > stripQualitiesBad() const
std::vector< float > stripNoises() const
edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > siStripNoisesToken_
const SiStripCluster * cluster_ptr
std::vector< float > stripNoisesRescaledByGain() const
uint16_t firstStrip() const
constexpr float gains[NGAINS]
Definition: EcalConstants.h:11
const SiStripGain * siStripGain_
std::pair< uint16_t, uint16_t > chargeLR() const
bool IsStripBad(uint32_t detid, short strip) const
const SiStripQuality * siStripQuality_
uint16_t width() const
const Range getRange(const uint32_t detID) const
edm::ESGetToken< SiStripGain, SiStripGainRcd > siStripGainToken_
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:47
void initEvent(const edm::EventSetup &iSetup)
bool IsModuleUsable(uint32_t detid) const
bool IsFiberBad(uint32_t detid, short fiberNb) const