CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TkStripMeasurementDet.cc
Go to the documentation of this file.
10 
11 #include <typeinfo>
14 #include<cassert>
15 
17  MeasurementDet (gdet), index_(-1), theDetConditions(&conditions)
18  {
19  if (dynamic_cast<const StripGeomDetUnit*>(gdet) == 0) {
20  throw MeasurementDetException( "TkStripMeasurementDet constructed with a GeomDet which is not a StripGeomDetUnit");
21  }
22  }
23 
24 
25 // fast check if the det contains any useful cluster
27  if unlikely( (!isActive(data)) || isEmpty(data.stripData())) return true;
28 
29  const detset & detSet = data.stripData().detSet(index());
30  for ( auto ci = detSet.begin(); ci != detSet.end(); ++ ci ) {
31  if (isMasked(*ci)) continue;
32  SiStripClusterRef cluster = edmNew::makeRefTo( data.stripData().handle(), ci );
33  if (accept(cluster, data.stripClustersToSkip()))
34  return false;
35  }
36  return true;
37 }
38 
39 
42 {
44  if unlikely( (!isActive(data)) || isEmpty(data.stripData())) return result;
45  const detset & detSet = data.stripData().detSet(index());
46  result.reserve(detSet.size());
47  for ( new_const_iterator ci = detSet.begin(); ci != detSet.end(); ++ ci ) {
48  if (isMasked(*ci)) continue;
49  // for ( ClusterIterator ci=theClusterRange.first; ci != theClusterRange.second; ci++) {
50  SiStripClusterRef cluster = edmNew::makeRefTo( data.stripData().handle(), ci );
51  if (accept(cluster, data.stripClustersToSkip()))
52  result.push_back( buildRecHit( cluster, ts));
53  else LogDebug("TkStripMeasurementDet")<<"skipping this str from last iteration on"<<rawId()<<" key: "<<cluster.key();
54  }
55  return result;
56 
57 }
58 
59 
60 // FIXME need to be merged with simpleRecHits
62  const TrajectoryStateOnSurface& stateOnThisDet,
63  const MeasurementEstimator& est, const MeasurementTrackerEvent & data) const {
64  if unlikely( (!isActive(data)) || isEmpty(data.stripData())) return false;
65  auto oldSize = result.size();
66 
67  float utraj = specificGeomDet().specificTopology().measurementPosition( stateOnThisDet.localPosition()).x();
68  const detset & detSet = data.stripData().detSet(index());
69  auto rightCluster =
70  std::find_if( detSet.begin(), detSet.end(), [utraj](const SiStripCluster& hit) { return hit.barycenter() > utraj; });
71 
72 
73  std::vector<SiStripRecHit2D> tmp;
74  if ( rightCluster != detSet.begin()) {
75  // there are hits on the left of the utraj
76  auto leftCluster = rightCluster;
77  while ( --leftCluster >= detSet.begin()) {
78  SiStripClusterRef clusterref = edmNew::makeRefTo( data.stripData().handle(), leftCluster );
79  bool isCompatible = filteredRecHits(clusterref, stateOnThisDet, est, data.stripClustersToSkip(), tmp);
80  if(!isCompatible) break; // exit loop on first incompatible hit
81  for (auto && h: tmp) result.push_back(new SiStripRecHit2D(std::move(h))); tmp.clear();
82  }
83  }
84  for ( ; rightCluster != detSet.end(); rightCluster++) {
85  SiStripClusterRef clusterref = edmNew::makeRefTo( data.stripData().handle(), rightCluster );
86  bool isCompatible = filteredRecHits(clusterref, stateOnThisDet, est, data.stripClustersToSkip(), tmp);
87  if(!isCompatible) break; // exit loop on first incompatible hit
88  for (auto && h: tmp) result.push_back(new SiStripRecHit2D(std::move(h))); tmp.clear();
89  }
90 
91  return result.size()>oldSize;
92 }
93 
94 
95 
96 
99  std::vector<SiStripRecHit2D> &result) const {
100  if unlikely( (!isActive(data)) || isEmpty(data.stripData())) return false;
101 
102  auto oldSize = result.size();
103 
104  float utraj = specificGeomDet().specificTopology().measurementPosition( stateOnThisDet.localPosition()).x();
105  const detset & detSet = data.stripData().detSet(index());
106  auto rightCluster =
107  std::find_if( detSet.begin(), detSet.end(), [utraj](const SiStripCluster& hit) { return hit.barycenter() > utraj; });
108 
109  if ( rightCluster != detSet.begin()) {
110  // there are hits on the left of the utraj
111  auto leftCluster = rightCluster;
112  while ( --leftCluster >= detSet.begin()) {
113  SiStripClusterRef clusterref = edmNew::makeRefTo( data.stripData().handle(), leftCluster );
114  bool isCompatible = filteredRecHits(clusterref, stateOnThisDet, est, data.stripClustersToSkip(), result);
115  if(!isCompatible) break; // exit loop on first incompatible hit
116  }
117  }
118  for ( ; rightCluster != detSet.end(); rightCluster++) {
119  SiStripClusterRef clusterref = edmNew::makeRefTo( data.stripData().handle(), rightCluster );
120  bool isCompatible = filteredRecHits(clusterref, stateOnThisDet, est, data.stripClustersToSkip(), result);
121  if(!isCompatible) break; // exit loop on first incompatible hit
122  }
123 
124  return result.size()>oldSize;
125 }
126 
127 
128 
129 bool
131  RecHitContainer & result, std::vector<float> & diffs ) const {
132  if unlikely( (!isActive(data)) || isEmpty(data.stripData())) return false;
133 
134  auto oldSize = result.size();
135 
136  float utraj = specificGeomDet().specificTopology().measurementPosition( stateOnThisDet.localPosition()).x();
137 
138  const detset & detSet = data.stripData().detSet(index());
139  auto rightCluster =
140  std::find_if( detSet.begin(), detSet.end(), [utraj](const SiStripCluster& hit) { return hit.barycenter() > utraj; });
141 
142  if ( rightCluster != detSet.begin()) {
143  // there are hits on the left of the utraj
144  auto leftCluster = rightCluster;
145  while ( --leftCluster >= detSet.begin()) {
146  SiStripClusterRef clusterref = edmNew::makeRefTo( data.stripData().handle(), leftCluster );
147  bool isCompatible = filteredRecHits(clusterref, stateOnThisDet, est, data.stripClustersToSkip(), result, diffs);
148  if(!isCompatible) break; // exit loop on first incompatible hit
149  }
150  }
151  for ( ; rightCluster != detSet.end(); rightCluster++) {
152  SiStripClusterRef clusterref = edmNew::makeRefTo( data.stripData().handle(), rightCluster );
153  bool isCompatible = filteredRecHits(clusterref, stateOnThisDet, est, data.stripClustersToSkip(), result,diffs);
154  if(!isCompatible) break; // exit loop on first incompatible hit
155  }
156 
157  return result.size()>oldSize;
158 }
159 
162  TempMeasurements & result) const {
163 
164  if (!isActive(data)) {
165  LogDebug("TkStripMeasurementDet")<<" found an inactive module "<<rawId();
166  result.add(theInactiveHit, 0.F);
167  return true;
168  }
169 
170  if (!isEmpty(data.stripData())){
171  LogDebug("TkStripMeasurementDet")<<" found hit on this module "<<rawId();
173  std::vector<float> diffs;
174  if (recHits(stateOnThisDet,est,data,result.hits,result.distances)) return true;
175  }
176 
177 
178  // create a TrajectoryMeasurement with an invalid RecHit and zero estimate
179 
180  if (!stateOnThisDet.hasError()) {
181  result.add(theMissingHit, 0.F);
182  return false;
183  }
184 
185  float utraj = specificGeomDet().specificTopology().measurementPosition( stateOnThisDet.localPosition()).x();
186  float uerr= sqrt(specificGeomDet().specificTopology().measurementError(stateOnThisDet.localPosition(),stateOnThisDet.localError().positionError()).uu());
187  if (testStrips(utraj,uerr)) {
188  //LogDebug("TkStripMeasurementDet") << " DetID " << id_ << " empty after search, but active ";
189  result.add(theMissingHit, 0.F);
190  return false;
191  }
192 
193  //LogDebug("TkStripMeasurementDet") << " DetID " << id_ << " empty after search, and inactive ";
194  result.add(theInactiveHit, 0.F);
195  return true;
196 
197 }
198 
199 
200 
201 
202 
203 
204 
205 void
207 {
208  if (isEmpty(data.stripData()) || !isActive(data)) return;
209 
210  const detset & detSet = data.stripData().detSet(index());
211  result.reserve(detSet.size());
212  for ( new_const_iterator ci = detSet.begin(); ci != detSet.end(); ++ ci ) {
213  if (isMasked(*ci)) continue;
214  // for ( ClusterIterator ci=theClusterRange.first; ci != theClusterRange.second; ci++) {
215  SiStripClusterRef cluster = edmNew::makeRefTo( data.stripData().handle(), ci );
216  if (accept(cluster, data.stripClustersToSkip()))
217  buildSimpleRecHit( cluster, ts,result);
218  else LogDebug("TkStripMeasurementDet")<<"skipping this str from last iteration on"<<rawId()<<" key: "<<cluster.key();
219  }
220 }
221 
222 
223 
224 std::tuple<TkStripRecHitIter,TkStripRecHitIter>
226 {
227  if (isEmpty(data.stripData()) || !isActive(data)) return std::tuple<TkStripRecHitIter,TkStripRecHitIter>();
228  const detset & detSet = data.stripData().detSet(index());
229  return std::make_tuple(TkStripRecHitIter(detSet.begin(),detSet.end(),*this,ts,data),
230  TkStripRecHitIter(detSet.end(),detSet.end(),*this,ts,data)
231  );
232 }
233 
235  while (!hi.empty()) {
236  auto ci = hi.clusterI;
237  auto const & data = *hi.data;
238  if (isMasked(*ci)) continue;
239  SiStripClusterRef cluster = edmNew::makeRefTo( data.stripData().handle(), ci );
240  if (accept(cluster, data.stripClustersToSkip())) return;
241  ++hi.clusterI;
242  }
243 }
244 
246  const GeomDetUnit& gdu( specificGeomDet());
247  auto ci = hi.clusterI;
248  auto const & data = *hi.data;
249  auto const & ltp = *hi.tsos;
250 
251  SiStripClusterRef cluster = edmNew::makeRefTo( data.stripData().handle(), ci );
252  LocalValues lv = cpe()->localParameters( *cluster, gdu, ltp);
253  return SiStripRecHit2D(lv.first,lv.second, gdu, cluster);
254 }
255 
256 bool
257 TkStripMeasurementDet::testStrips(float utraj, float uerr) const {
258  int16_t start = (int16_t) std::max<float>(utraj - 3.f*uerr, 0);
259  int16_t end = (int16_t) std::min<float>(utraj + 3.f*uerr, totalStrips());
260 
261  if (start >= end) { // which means either end <=0 or start >= totalStrips_
262  /* LogDebug("TkStripMeasurementDet") << "Testing module " << id_ <<","<<
263  " U = " << utraj << " +/- " << uerr <<
264  "; Range [" << (utraj - 3*uerr) << ", " << (utraj + 3*uerr) << "] " <<
265  ": YOU'RE COMPLETELY OFF THE MODULE."; */
266  //return false;
267  return true; // Wolfgang thinks this way is better
268  // and solves some problems with grouped ckf
269  }
270 
271  typedef std::vector<BadStripBlock>::const_iterator BSBIT;
272 
273  int16_t bad = 0, largestBadBlock = 0;
274  for (BSBIT bsbc = badStripBlocks().begin(), bsbe = badStripBlocks().end(); bsbc != bsbe; ++bsbc) {
275  if (bsbc->last < start) continue;
276  if (bsbc->first > end) break;
277  int16_t thisBad = std::min(bsbc->last, end) - std::max(bsbc->first, start);
278  if (thisBad > largestBadBlock) largestBadBlock = thisBad;
279  bad += thisBad;
280  }
281 
282  bool ok = (bad < (end-start)) &&
283  (uint16_t(bad) <= badStripCuts().maxBad) &&
284  (uint16_t(largestBadBlock) <= badStripCuts().maxConsecutiveBad);
285 
286 // if (bad) {
287 // edm::LogWarning("TkStripMeasurementDet") << "Testing module " << id_ <<" (subdet: "<< SiStripDetId(id_).subdetId() << ")" <<
288 // " U = " << utraj << " +/- " << uerr <<
289 // "; Range [" << (utraj - 3*uerr) << ", " << (utraj + 3*uerr) << "] " <<
290 // "= [" << start << "," << end << "]" <<
291 // " total strips:" << (end-start) << ", good:" << (end-start-bad) << ", bad:" << bad << ", largestBadBlock: " << largestBadBlock <<
292 // ". " << (ok ? "OK" : "NO");
293 // }
294  return ok;
295 }
296 
#define LogDebug(id)
virtual LocalValues localParameters(const SiStripCluster &, const GeomDetUnit &) const
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
void buildSimpleRecHit(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp, std::vector< SiStripRecHit2D > &res) const
const StripDetset & detSet(int i) const
const TrajectoryStateOnSurface * tsos
std::tuple< TkStripRecHitIter, TkStripRecHitIter > hitRange(const TrajectoryStateOnSurface &, const MeasurementTrackerEvent &data) const
bool empty(const MeasurementTrackerEvent &data) const
StripClusterParameterEstimator::LocalValues LocalValues
TkStripMeasurementDet(const GeomDet *gdet, StMeasurementConditionSet &conditionSet)
BadStripCuts const & badStripCuts() const
void simpleRecHits(const TrajectoryStateOnSurface &ts, const MeasurementTrackerEvent &data, std::vector< SiStripRecHit2D > &result) const
bool accept(SiStripClusterRef const &r, const std::vector< bool > &skipClusters) const
const StMeasurementDetSet & stripData() const
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
const std::vector< bool > & stripClustersToSkip() const
LocalError positionError() const
bool isMasked(const SiStripCluster &cluster) const
SiStripRecHit2D hit(TkStripRecHitIter const &hi) const
int bad(Items const &cont)
TrackingRecHit::ConstRecHitPointer theMissingHit
const MeasurementTrackerEvent * data
#define unlikely(x)
virtual RecHitContainer recHits(const TrajectoryStateOnSurface &, const MeasurementTrackerEvent &data) const
virtual bool measurements(const TrajectoryStateOnSurface &stateOnThisDet, const MeasurementEstimator &est, const MeasurementTrackerEvent &data, TempMeasurements &result) const
bool filteredRecHits(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp, const MeasurementEstimator &est, const std::vector< bool > &skipClusters, RecHitContainer &result, std::vector< float > &diffs) const
T sqrt(T t)
Definition: SSEVec.h:48
void add(ConstRecHitPointer const &h, float d)
std::vector< BadStripBlock > const & badStripBlocks() const
tuple result
Definition: query.py:137
double f[11][100]
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
#define end
Definition: vmac.h:37
T min(T a, T b)
Definition: MathUtil.h:58
const LocalTrajectoryError & localError() const
void advance(TkStripRecHitIter &hi) const
TrackingRecHit::RecHitPointer buildRecHit(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp) const
TrackingRecHit::ConstRecHitPointer theInactiveHit
key_type key() const
Accessor for product key.
Definition: Ref.h:266
bool isActive(const MeasurementTrackerEvent &data) const
Is this module active in reconstruction? It must be both &#39;setActiveThisEvent&#39; and &#39;setActive&#39;...
const StripGeomDetUnit & specificGeomDet() const
bool isEmpty(const StMeasurementDetSet &theDets) const
new_const_iterator clusterI
const StripClusterParameterEstimator * cpe() const
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
TrackingRecHit::ConstRecHitContainer RecHitContainer
std::vector< BaseTrackerRecHit * > SimpleHitContainer
const detset & detSet(const StMeasurementDetSet &theDets) const
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
iterator end()
Definition: DetSetNew.h:70
detset::const_iterator new_const_iterator
bool testStrips(float utraj, float uerr) const
return true if there are &#39;enough&#39; good strips in the utraj +/- 3 uerr range.
unsigned int rawId() const
edm::Handle< edmNew::DetSetVector< SiStripCluster > > & handle()
Definition: DDAxes.h:10
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
size_type size() const
Definition: DetSetNew.h:86
iterator begin()
Definition: DetSetNew.h:67