CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TkStripMeasurementDet.h
Go to the documentation of this file.
1 #ifndef TkStripMeasurementDet_H
2 #define TkStripMeasurementDet_H
3 
5 #include "TkMeasurementDetSet.h"
6 
20 
22 
24 
25 class TkStripMeasurementDet GCC11_FINAL : public MeasurementDet {
26 public:
27 
30 
31  typedef SiStripRecHit2D::ClusterRef SiStripClusterRef;
32 
34 
37 
38  typedef std::vector<SiStripCluster>::const_iterator const_iterator;
39 
41 
42  TkStripMeasurementDet( const GeomDet* gdet, StMeasurementDetSet & dets);
43 
44 
45  void setIndex(int i) { index_=i;}
46 
47  void update( const detset &detSet ) {
48  theDets().update(index(),detSet);
49  }
50  void update( std::vector<SiStripCluster>::const_iterator begin ,std::vector<SiStripCluster>::const_iterator end ) {
51  theDets().update(index(), begin, end);
52  }
53 
54  bool isRegional() const { return theDets().isRegional();}
55 
56  void setEmpty(){ theDets().setEmpty(index()); }
57 
58  bool isEmpty() const {return theDets().empty(index());}
59 
60  int index() const { return index_;}
61 
62  unsigned int rawId() const { return theDets().id(index()); }
63  unsigned char subId() const { return theDets().subId(index());}
64 
65 
66  const detset & theSet() const {return theDets().detSet(index());}
67  const detset & detSet() const {return theDets().detSet(index());}
68  detset & detSet() { return theDets().detSet(index());}
69  unsigned int beginClusterI() const {return theDets().beginClusterI(index());}
70  unsigned int endClusterI() const {return theDets().endClusterI(index());}
71 
72  int size() const {return endClusterI() - beginClusterI() ; }
73 
74 
76  bool isActive() const { return theDets().isActive(index()); }
77 
78  //TO BE IMPLEMENTED
79  bool hasBadComponents( const TrajectoryStateOnSurface &tsos ) const {return false;}
80 
81 
82 
83  virtual RecHitContainer recHits( const TrajectoryStateOnSurface&) const;
84  void simpleRecHits( const TrajectoryStateOnSurface& ts, std::vector<SiStripRecHit2D> &result) const ;
85 
86  virtual bool recHits( const TrajectoryStateOnSurface& stateOnThisDet, const MeasurementEstimator& est,
87  RecHitContainer & result, std::vector<float> & diffs) const;
88 
89  virtual bool measurements( const TrajectoryStateOnSurface& stateOnThisDet,
90  const MeasurementEstimator& est,
91  TempMeasurements & result) const;
92 
93  const StripGeomDetUnit& specificGeomDet() const {return static_cast<StripGeomDetUnit const &>(fastGeomDet());}
94 
95 
96  template<class ClusterRefT>
98  buildRecHit( const ClusterRefT &cluster, const TrajectoryStateOnSurface& ltp) const {
99  const GeomDetUnit& gdu( specificGeomDet());
100  LocalValues lv = cpe()->localParameters( *cluster, gdu, ltp);
101  return TSiStripRecHit2DLocalPos::build( lv.first, lv.second, &fastGeomDet(), cluster, cpe());
102  }
103 
104 
105  template<class ClusterRefT>
106  void
107  buildRecHits( const ClusterRefT& cluster, const TrajectoryStateOnSurface& ltp, const RecHitContainer& _res) const {
108  RecHitContainer res = _res;
109  const GeomDetUnit& gdu( specificGeomDet());
110  VLocalValues vlv = cpe()->localParametersV( *cluster, gdu, ltp);
111  for(VLocalValues::const_iterator it=vlv.begin();it!=vlv.end();++it)
112  res.push_back(TSiStripRecHit2DLocalPos::build( it->first, it->second, &fastGeomDet(), cluster, cpe()));
113  }
114 
115 
116  template<class ClusterRefT>
117  bool filteredRecHits( const ClusterRefT& cluster, const TrajectoryStateOnSurface& ltp, const MeasurementEstimator& est,
118  RecHitContainer & result, std::vector<float> & diffs) const {
119  if (isMasked(*cluster)) return true;
120  const GeomDetUnit& gdu( specificGeomDet());
121  if (!accept(cluster)) return true;
122  VLocalValues const & vlv = cpe()->localParametersV( *cluster, gdu, ltp);
123  bool isCompatible(false);
124  for(auto vl : vlv) {
125  auto && recHit = TSiStripRecHit2DLocalPos::build( vl.first, vl.second, &fastGeomDet(), cluster, cpe());
126  std::pair<bool,double> diffEst = est.estimate(ltp, *recHit);
127  if ( diffEst.first ) {
128  result.push_back(std::move(recHit));
129  diffs.push_back(diffEst.second);
130  isCompatible = true;
131  }
132  }
133  return isCompatible;
134  }
135 
136 
137 
140  void setActive(bool active) { theDets().setActive(index(),active);}
143  void setActiveThisEvent(bool active) { theDets().setActiveThisEvent(index(),active); }
144 
146  bool hasAllGoodChannels() const { return (!hasAny128StripBad()) && badStripBlocks().empty(); }
147 
149  void set128StripStatus(bool good, int idx=-1) {
150  theDets().set128StripStatus(index(),good,idx);
151  }
152 
154 
156  bool testStrips(float utraj, float uerr) const;
157 
159 
160  std::vector<BadStripBlock> & getBadStripBlocks() { return theDets().getBadStripBlocks(index()); }
161  std::vector<BadStripBlock> const & badStripBlocks() const { return theDets().badStripBlocks(index()); }
162 
163  bool maskBad128StripBlocks() const { return theDets().maskBad128StripBlocks();}
164 
165 
166 
167 private:
168 
169  StMeasurementDetSet & theDets() { return *theDets_;}
170  StMeasurementDetSet & theDets() const { return *theDets_;}
171 
173  int index_;
174 
175 
176 
177  edm::Handle<edmNew::DetSetVector<SiStripCluster> > const & handle() const { return theDets().handle();}
178  edm::Handle<edm::LazyGetter<SiStripCluster> > const & regionalHandle() const { return theDets().regionalHandle();}
179 
180  const StripClusterParameterEstimator* cpe() const { return theDets().stripCPE(); }
181 
182 
183  const std::vector<bool> & skipClusters() const { return theDets().clusterToSkip();}
184 
185  // --- regional unpacking
186 
187  int totalStrips() const { return theDets().totalStrips(index()); }
188  BadStripCuts const & badStripCuts() const { return theDets().badStripCuts(index());}
189 
190  bool hasAny128StripBad() const { return theDets().hasAny128StripBad(index()); }
191 
192 
193 
194 
195  inline bool isMasked(const SiStripCluster &cluster) const {
196  return theDets().isMasked(index(), cluster);
197  }
198 
199 
200  template<class ClusterRefT>
201  void buildSimpleRecHit( const ClusterRefT& cluster,
202  const TrajectoryStateOnSurface& ltp,
203  std::vector<SiStripRecHit2D>& res) const {
204  const GeomDetUnit& gdu( specificGeomDet());
205  VLocalValues const & vlv = cpe()->localParametersV( *cluster, gdu, ltp);
206  for(VLocalValues::const_iterator it=vlv.begin();it!=vlv.end();++it){
207  res.push_back(SiStripRecHit2D( it->first, it->second, rawId(), cluster));
208  }
209  }
210 
211 
212 
213 
214 public:
215  inline bool accept(SiStripClusterRef const & r) const {
216  if(skipClusters().empty()) return true;
217  if (r.key()>=skipClusters().size()){
218  edm::LogError("WrongStripMasking")<<r.key()<<" is larger than: "<<skipClusters().size()<<" no skipping done";
219  return true;
220  }
221  return (not (skipClusters()[r.key()]));
222  }
223  inline bool accept(SiStripRegionalClusterRef const &r) const{
224  if(skipClusters().empty()) return true;
225  if (r.key()>=skipClusters().size()){
226  LogDebug("TkStripMeasurementDet")<<r.key()<<" is larger than: "<<skipClusters().size()
227  <<"\n This must be a new cluster, and therefore should not be skiped most likely.";
228  return true;
229  }
230  return (not (skipClusters()[r.key()]));
231  }
232 
233 };
234 
235 #endif
#define LogDebug(id)
virtual ~TkStripMeasurementDet()
int i
Definition: DBlmapReader.cc:9
virtual HitReturnType estimate(const TrajectoryStateOnSurface &ts, const TransientTrackingRecHit &hit) const =0
bool hasBadComponents(const TrajectoryStateOnSurface &tsos) const
TransientTrackingRecHit::ConstRecHitContainer RecHitContainer
bool isRegional() const
bool hasAny128StripBad() const
unsigned int beginClusterI() const
StMeasurementDetSet::BadStripBlock BadStripBlock
data_type const * const_iterator
Definition: DetSetNew.h:25
void set128StripStatus(bool good, int idx=-1)
Sets the status of a block of 128 strips (or all blocks if idx=-1)
bool hasAllGoodChannels() const
does this module have at least one bad strip, APV or channel?
int totalStrips() const
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:22
int size() const
std::pair< LocalPoint, LocalError > LocalValues
const StripGeomDetUnit & specificGeomDet() const
const StripClusterParameterEstimator * cpe() const
int index() const
const GeomDet & fastGeomDet() const
virtual RecHitContainer recHits(const TrajectoryStateOnSurface &) const =0
edmNew::DetSet< SiStripCluster > detset
bool accept(SiStripRegionalClusterRef const &r) const
tuple result
Definition: query.py:137
void setActiveThisEvent(bool active)
Turn on/off the module for reconstruction for one events. This per-event flag is cleared by any call ...
bool isEmpty() const
StMeasurementDetSet * theDets_
edm::Handle< edm::LazyGetter< SiStripCluster > > const & regionalHandle() const
void buildSimpleRecHit(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp, std::vector< SiStripRecHit2D > &res) const
SiStripRecHit2D::ClusterRef SiStripClusterRef
#define end
Definition: vmac.h:38
edm::Handle< edmNew::DetSetVector< SiStripCluster > > const & handle() const
StMeasurementDetSet & theDets()
unsigned char subId() const
StMeasurementDetSet & theDets() const
StMeasurementDetSet::BadStripCuts BadStripCuts
std::vector< BadStripBlock > const & badStripBlocks() const
const detset & detSet() const
const detset & theSet() const
void setIndex(int i)
std::vector< BadStripBlock > & getBadStripBlocks()
bool accept(SiStripClusterRef const &r) const
void buildRecHits(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp, const RecHitContainer &_res) const
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
key_type key() const
Accessor for product key.
Definition: Ref.h:266
void update(const detset &detSet)
bool maskBad128StripBlocks() const
void update(std::vector< SiStripCluster >::const_iterator begin, std::vector< SiStripCluster >::const_iterator end)
unsigned int rawId() const
edm::LazyGetter< SiStripCluster >::value_ref SiStripRegionalClusterRef
#define begin
Definition: vmac.h:31
bool filteredRecHits(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp, const MeasurementEstimator &est, RecHitContainer &result, std::vector< float > &diffs) const
StripClusterParameterEstimator::VLocalValues VLocalValues
unsigned int endClusterI() const
bool isActive() const
Is this module active in reconstruction? It must be both &#39;setActiveThisEvent&#39; and &#39;setActive&#39;...
const std::vector< bool > & skipClusters() const
void setActive(bool active)
Turn on/off the module for reconstruction, for the full run or lumi (using info from DB...
BadStripCuts const & badStripCuts() const
std::vector< SiStripCluster >::const_iterator const_iterator
TransientTrackingRecHit::RecHitPointer buildRecHit(const ClusterRefT &cluster, const TrajectoryStateOnSurface &ltp) const
virtual bool measurements(const TrajectoryStateOnSurface &stateOnThisDet, const MeasurementEstimator &est, TempMeasurements &result) const =0
bool isMasked(const SiStripCluster &cluster) const
detset::const_iterator new_const_iterator
StripClusterParameterEstimator::LocalValues LocalValues
Unlimited (trivial) bounds.