CMS 3D CMS Logo

SiStripClusterInfo.cc
Go to the documentation of this file.
2 
9 #include <cmath>
10 
11 
12 
14  const edm::EventSetup& setup,
15  const int detId,
16  const std::string & quality)
17  : cluster_ptr(&cluster),
18  es(setup),
19  qualityLabel(quality),
20  detId_(detId) {
24 }
25 
26 std::pair<uint16_t,uint16_t > SiStripClusterInfo::
27 chargeLR() const {
28  std::vector<uint8_t>::const_iterator
29  begin( stripCharges().begin() ),
30  end( stripCharges().end() ),
31  max; max = max_element(begin,end);
32  return std::make_pair( accumulate(begin, max, uint16_t(0) ),
33  accumulate(max+1, end, uint16_t(0) ) );
34 }
35 
36 
38 variance() const {
39  float q(0), x1(0), x2(0);
40  for(auto
41  begin(stripCharges().begin()), end(stripCharges().end()), it(begin);
42  it!=end; ++it) {
43  unsigned i = it-begin;
44  q += (*it);
45  x1 += (*it) * (i+0.5);
46  x2 += (*it) * (i*i+i+1./3);
47  }
48  return (x2 - x1*x1/q ) / q;
49 }
50 
51 std::vector<float> SiStripClusterInfo::
55 
56  std::vector<float> results;
57  results.reserve(width());
58  for(size_t i = 0, e = width(); i < e; i++){
59  results.push_back(noiseHandle->getNoise(firstStrip()+i, detNoiseRange) / gainHandle->getStripGain( firstStrip()+i, detGainRange));
60  }
61  return results;
62 }
63 
64 std::vector<float> SiStripClusterInfo::
65 stripNoises() const {
67 
68  std::vector<float> noises;
69  noises.reserve(width());
70  for(size_t i=0; i < width(); i++){
71  noises.push_back( noiseHandle->getNoise( firstStrip()+i, detNoiseRange) );
72  }
73  return noises;
74 }
75 
76 std::vector<float> SiStripClusterInfo::
77 stripGains() const {
79 
80  std::vector<float> gains;
81  gains.reserve(width());
82  for(size_t i=0; i< width(); i++){
83  gains.push_back( gainHandle->getStripGain( firstStrip()+i, detGainRange) );
84  }
85  return gains;
86 }
87 
88 std::vector<bool> SiStripClusterInfo::
90  std::vector<bool> isBad;
91  isBad.reserve(width());
92  for(int i=0; i< width(); i++) {
93  isBad.push_back( qualityHandle->IsStripBad( detId_,
94  firstStrip()+i) );
95  }
96  return isBad;
97 }
98 
100 calculate_noise(const std::vector<float>& noise) const {
101  float noiseSumInQuadrature = 0;
102  int numberStripsOverThreshold = 0;
103  for(int i=0;i<width();i++) {
104  if(stripCharges()[i]!=0) {
105  noiseSumInQuadrature += noise.at(i) * noise.at(i);
106  numberStripsOverThreshold++;
107  }
108  }
109  return std::sqrt( noiseSumInQuadrature / numberStripsOverThreshold );
110 }
111 
112 
114 IsAnythingBad() const {
115  std::vector<bool> stripBad = stripQualitiesBad();
116  return
117  IsApvBad() ||
118  IsFiberBad() ||
119  IsModuleBad() ||
120  accumulate(stripBad.begin(), stripBad.end(),
121  false,
122  std::logical_or<bool>());
123 }
124 
126 IsApvBad() const {
127  return
129  qualityHandle->IsApvBad( detId_, (firstStrip()+width())/128 ) ;
130 }
131 
133 IsFiberBad() const {
134  return
137 }
138 
140 IsModuleBad() const {
141  return qualityHandle->IsModuleBad( detId_ );
142 }
143 
145 IsModuleUsable() const {
147 }
148 
149 std::vector<SiStripCluster> SiStripClusterInfo::
150 reclusterize(const edm::ParameterSet& conf) const {
151 
152  std::vector<SiStripCluster> clusters;
153 
154  std::vector<uint8_t> charges(stripCharges().begin(),stripCharges().end());
155  std::vector<float> gains = stripGains();
156  for(unsigned i=0; i < charges.size(); i++)
157  charges[i] = (charges[i] < 254)
158  ? static_cast<uint8_t>(charges[i] * gains[i])
159  : charges[i];
160 
161  std::auto_ptr<StripClusterizerAlgorithm>
163  algorithm->initialize(es);
164  auto const & det = algorithm->stripByStripBegin( detId_ );
165  if(det.valid()) {
167  for(unsigned i = 0; i<width(); i++)
168  algorithm->stripByStripAdd(state, firstStrip()+i, charges[i], clusters );
169  algorithm->stripByStripEnd(state, clusters );
170  }
171 
172  return clusters;
173 }
174 
bool IsFiberBad(const uint32_t &detid, const short &fiberNb) const
bool IsApvBad(const uint32_t &detid, const short &apvNb) 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
std::pair< uint16_t, uint16_t > chargeLR() const
bool IsStripBad(const uint32_t &detid, const short &strip) const
std::vector< bool > stripQualitiesBad() const
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
float variance() const
float calculate_noise(const std::vector< float > &) const
auto stripCharges() const -> decltype(cluster() ->amplitudes())
static float getNoise(uint16_t strip, const Range &range)
Definition: SiStripNoises.h:72
bool IsModuleUsable(const uint32_t &detid) const
static float getStripGain(const uint16_t &strip, const SiStripApvGain::Range &range)
Definition: SiStripGain.h:68
static std::auto_ptr< StripClusterizerAlgorithm > create(const edm::ParameterSet &)
T sqrt(T t)
Definition: SSEVec.h:18
uint16_t width() const
std::pair< ContainerIterator, ContainerIterator > Range
bool IsModuleBad(const uint32_t &detid) const
#define end
Definition: vmac.h:37
const T & get() const
Definition: EventSetup.h:56
const edm::EventSetup & es
bool IsModuleUsable() const
edm::ESHandle< SiStripNoises > noiseHandle
#define begin
Definition: vmac.h:30
std::vector< float > stripNoises() const
const Range getRange(const uint32_t detID) const
std::vector< SiStripCluster > reclusterize(const edm::ParameterSet &) const
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:48
std::vector< float > stripNoisesRescaledByGain() const
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:66