CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TkMeasurementDetSet.h
Go to the documentation of this file.
1 #ifndef StMeasurementDetSet_H
2 #define StMeasurementDetSet_H
3 
4 #include<vector>
9 
15 
18 
20 
21 /* Struct of arrays supporting "members of Tk...MeasurementDet
22  * implemented with vectors, to be optimized...
23  */
25 public:
26 
29 
30  typedef std::vector<SiStripCluster>::const_iterator const_iterator;
31 
34 
35  enum QualityFlags { BadModules=1, // for everybody
36  /* Strips: */ BadAPVFibers=2, BadStrips=4, MaskBad128StripBlocks=8 };
37 
38 
39  struct BadStripCuts {
42  maxBad(pset.getParameter<uint32_t>("maxBad")),
43  maxConsecutiveBad(pset.getParameter<uint32_t>("maxConsecutiveBad")) {}
45  };
46 
47  struct BadStripBlock {
48  short first;
49  short last;
50  BadStripBlock(const SiStripBadStrip::data &data) : first(data.firstStrip), last(data.firstStrip+data.range-1) { }
51  };
52 
53 
56  bool regional):
57  theMatcher(matcher), theCPE(cpe), regional_(regional){}
58 
59 
60  void init(std::vector<TkStripMeasurementDet> & stripDets);
61 
62 
63  std::vector<bool> & clusterToSkip() const { return theStripsToSkip; }
64 
66 
67  void update(int i,const StripDetset & detSet ) {
68  detSet_[i] = detSet;
69  empty_[i] = false;
70  }
71 
72  void update(int i, std::vector<SiStripCluster>::const_iterator begin ,std::vector<SiStripCluster>::const_iterator end) {
73  clusterI_[2*i] = begin - regionalHandle_->begin_record();
74  clusterI_[2*i+1] = end - regionalHandle_->begin_record();
75 
76  empty_[i] = false;
77  activeThisEvent_[i] = true;
78  }
79 
80 
81  const SiStripRecHitMatcher* matcher() const { return theMatcher;}
83 
84 
85  int nDet() const { return id_.size();}
86 
87  unsigned int id(int i) const { return id_[i]; }
88  unsigned char subId(int i) const { return subId_[i];}
89 
90  int find(unsigned int jd, int i=0) {
91  return std::lower_bound(id_.begin()+i,id_.end(),jd)-id_.begin();
92  }
93 
94  bool isRegional() const { return regional_;}
95  bool empty(int i) const { return empty_[i];}
96  bool isActive(int i) const { return activeThisEvent_[i] && activeThisPeriod_[i]; }
97 
98  void setEmpty(int i) {empty_[i] = true; activeThisEvent_[i] = true; }
99 
100  void setEmpty() {
101  std::fill(empty_.begin(),empty_.end(),true);
102  std::fill(activeThisEvent_.begin(), activeThisEvent_.end(),true);
103  }
104 
107  void setActive(int i, bool active) { activeThisPeriod_[i] = active; activeThisEvent_[i] = true; if (!active) empty_[i] = true; }
110  void setActiveThisEvent(int i, bool active) { activeThisEvent_[i] = active; if (!active) empty_[i] = true; }
111 
112 
114  StripDetset & detSet(int i) { return detSet_[i];}
115 
117  unsigned int beginClusterI(int i) const {return clusterI_[2*i];}
118  unsigned int endClusterI(int i) const {return clusterI_[2*i+1];}
119 
120  int totalStrips(int i) const { return totalStrips_[i];}
121 
122 
123 
124  void setMaskBad128StripBlocks(bool maskThem) { maskBad128StripBlocks_ = maskThem; }
126 
127 
129  bool hasAny128StripBad(int i) const { return hasAny128StripBad_[i];}
130 
131  std::vector<BadStripBlock> & getBadStripBlocks(int i) { return badStripBlocks_[i]; }
132  std::vector<BadStripBlock> const & badStripBlocks(int i) const {return badStripBlocks_[i]; }
133 
134 
135  bool isMasked(int i, const SiStripCluster &cluster) const {
136  int offset = nbad128*i;
137  if ( bad128Strip_[offset+( cluster.firstStrip() >> 7)] ) {
138  if ( bad128Strip_[offset+( (cluster.firstStrip()+cluster.amplitudes().size()) >> 7)] ||
139  bad128Strip_[offset+( static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7)] ) {
140  return true;
141  }
142  } else {
143  if ( bad128Strip_[offset+( (cluster.firstStrip()+cluster.amplitudes().size()) >> 7)] &&
144  bad128Strip_[offset+( static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7)] ) {
145  return true;
146  }
147  }
148  return false;
149  }
150 
151 
152  void set128StripStatus(int i, bool good, int idx=-1);
153 
154  void initializeStripStatus(const SiStripQuality *quality, int qualityFlags, int qualityDebugFlags, edm::ParameterSet cutPset);
155 
156 private:
157 
159 
160  // globals
163 
166 
167  mutable std::vector<bool> theStripsToSkip;
168 
169  bool regional_;
170 
173 
174  // members of TkStripMeasurementDet
175  std::vector<unsigned int> id_;
176  std::vector<unsigned char> subId_;
177 
178  std::vector<int> totalStrips_;
179 
180  static const int nbad128 = 6;
181  std::vector<bool> bad128Strip_;
182  std::vector<bool> hasAny128StripBad_;
183 
184  std::vector<std::vector<BadStripBlock>> badStripBlocks_;
185 
186 
187  std::vector<bool> empty_;
188 
190 
191  // full reco
192  std::vector<StripDetset> detSet_;
193 
194  // --- regional unpacking
195 
196  // begin,end "pairs"
197  std::vector<unsigned int> clusterI_;
198 
199 
200 };
201 
202 
203 #endif // StMeasurementDetSet_H
int i
Definition: DBlmapReader.cc:9
string fill
Definition: lumiContext.py:319
StripDetset & detSet(int i)
const SiStripRecHitMatcher * matcher() const
const StripClusterParameterEstimator * theCPE
std::vector< bool > activeThisPeriod_
std::vector< BadStripBlock > & getBadStripBlocks(int i)
edm::Handle< edm::LazyGetter< SiStripCluster > > & regionalHandle()
std::vector< unsigned int > id_
void update(int i, const StripDetset &detSet)
std::vector< bool > & clusterToSkip() const
bool isMasked(int i, const SiStripCluster &cluster) const
StripDetset::const_iterator new_const_iterator
uint16_t firstStrip() const
data_type const * const_iterator
Definition: DetSetNew.h:25
void set128StripStatus(int i, bool good, int idx=-1)
BadStripBlock(const SiStripBadStrip::data &data)
std::vector< bool > bad128Strip_
const StripClusterParameterEstimator * stripCPE() const
unsigned int id(int i) const
edmNew::DetSet< SiStripCluster > StripDetset
BadStripCuts badStripCuts_[4]
std::vector< bool > hasAny128StripBad_
std::vector< std::vector< BadStripBlock > > badStripBlocks_
void initializeStripStatus(const SiStripQuality *quality, int qualityFlags, int qualityDebugFlags, edm::ParameterSet cutPset)
edm::LazyGetter< SiStripCluster > LazyGetter
const SiStripRecHitMatcher * theMatcher
edm::Handle< edm::LazyGetter< SiStripCluster > > regionalHandle_
void setActiveThisEvent(int i, bool active)
Turn on/off the module for reconstruction for one events. This per-event flag is cleared by any call ...
unsigned int endClusterI(int i) const
std::vector< unsigned int > clusterI_
unsigned int beginClusterI(int i) const
std::vector< unsigned char > subId_
StMeasurementDetSet(const SiStripRecHitMatcher *matcher, const StripClusterParameterEstimator *cpe, bool regional)
float barycenter() const
#define end
Definition: vmac.h:38
unsigned int offset(bool)
int find(unsigned int jd, int i=0)
std::vector< SiStripCluster >::const_iterator const_iterator
void setMaskBad128StripBlocks(bool maskThem)
std::vector< bool > empty_
bool maskBad128StripBlocks() const
unsigned char subId(int i) const
bool isActive(int i) const
edm::Handle< edmNew::DetSetVector< SiStripCluster > > handle_
int totalStrips(int i) const
std::vector< BadStripBlock > const & badStripBlocks(int i) const
bool empty(int i) const
void setLazyGetter(edm::Handle< LazyGetter > const &lg)
std::vector< int > totalStrips_
void setActive(int i, bool active)
Turn on/off the module for reconstruction, for the full run or lumi (using info from DB...
#define begin
Definition: vmac.h:31
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
static const int nbad128
std::vector< StripDetset > detSet_
std::vector< bool > activeThisEvent_
void update(int i, std::vector< SiStripCluster >::const_iterator begin, std::vector< SiStripCluster >::const_iterator end)
edm::Handle< edmNew::DetSetVector< SiStripCluster > > & handle()
std::vector< bool > theStripsToSkip
BadStripCuts(const edm::ParameterSet &pset)
BadStripCuts & badStripCuts(int i)
void init(std::vector< TkStripMeasurementDet > &stripDets)
bool hasAny128StripBad(int i) const
const std::vector< uint8_t > & amplitudes() const
edm::RefGetter< SiStripCluster > RefGetter