test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripRecHitMatcher.h
Go to the documentation of this file.
1 #ifndef RECOLOCALTRACKER_SISTRIPCLUSTERIZER_SISTRIPRECHITMATCH_H
2 #define RECOLOCALTRACKER_SISTRIPCLUSTERIZER_SISTRIPRECHITMATCH_H
3 
11 
13 
15 
17 
18 
19 class GluedGeomDet;
20 
21 #include <cfloat>
22 #include <memory>
23 
24 #include <boost/function.hpp>
25 
27 public:
28 
29  // may become a template argument
31 
33  typedef std::vector<const SiStripRecHit2D *> SimpleHitCollection;
34  typedef SimpleHitCollection::const_iterator SimpleHitIterator;
35 
36  typedef boost::function<void(SiStripMatchedRecHit2D const&)> Collector;
37 
38 
39  typedef std::pair<LocalPoint,LocalPoint> StripPosition;
40 
42  SiStripRecHitMatcher(const double theScale);
43 
44  bool preFilter() const { return preFilter_;}
45 
46 
47  inline
48  static float sigmaPitch(LocalPoint const& pos, LocalError const & err,
49  const StripTopology& topol) {
50  MeasurementError error = topol.measurementError(pos,err);
51  auto pitch=topol.localPitch(pos);
52  return error.uu()*pitch*pitch;
53  }
54 
55 
56  // optimized matching iteration (the algo is the same, just recoded)
57  template<typename MonoIterator, typename StereoIterator, typename CollectorHelper>
58  void doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend,
59  StereoIterator seconditer, StereoIterator seconditerend,
60  const GluedGeomDet* gluedDet, LocalVector trdir,
61  CollectorHelper & collectorHelper) const;
62 
63 
64  //match a single hit
65  std::unique_ptr<SiStripMatchedRecHit2D> match(const SiStripRecHit2D *monoRH,
66  const SiStripRecHit2D *stereoRH,
67  const GluedGeomDet* gluedDet,
68  LocalVector trackdirection, bool force=false) const;
69 
70 
71 // this is the one used by the RecHitConverter
72  void
73  match( const SiStripRecHit2D *monoRH,
75  CollectorMatched & collector,
76  const GluedGeomDet* gluedDet,
77  LocalVector trackdirection) const;
78 
79 
80  void
81  match( const SiStripRecHit2D *monoRH,
83  CollectorMatched & collector,
84  const GluedGeomDet* gluedDet,
85  LocalVector trackdirection) const;
86 
87 
88 
89 
90  // project strip coordinates on Glueddet
91 
92  StripPosition project(const GeomDetUnit *det,const GluedGeomDet* glueddet,StripPosition strip,LocalVector trackdirection) const;
93 
94 
95  // needed by the obsolete version still in use on some architectures
96  void
97  match( const SiStripRecHit2D *monoRH,
100  const GluedGeomDet* gluedDet,
101  LocalVector trackdirection) const;
102 
103 
104 
105  void
106  match( const SiStripRecHit2D *monoRH,
108  std::vector<SiStripMatchedRecHit2D*> & collector,
109  const GluedGeomDet* gluedDet,
110  LocalVector trackdirection) const;
111 
112 
113 
115  void
116  match( const SiStripRecHit2D *monoRH,
118  Collector & collector,
119  const GluedGeomDet* gluedDet,
120  LocalVector trackdirection) const;
121 
122 
123  float scale_;
124  bool preFilter_=false;
125 };
126 
127 
128 
129 
131 #if defined(USE_SSEVECT) || defined(USE_EXTVECT)
132 #define DOUBLE_MATCH
133 #endif
134 
135 #ifdef DOUBLE_MATCH
136 
137 #if defined(USE_SSEVECT)
138  using mathSSE::Vec3F;
139  using mathSSE::Vec2D;
140  using mathSSE::Vec3D;
141  using mathSSE::Rot3F;
142 #endif
143 
146 
147 namespace matcherDetails {
148 
149  struct StereoInfo {
150  Vec2D c1vec;
151  const SiStripRecHit2D * secondHit;
152  float sigmap22;
153  double m10, m11;
154  };
155 
156 }
157 
158 template<typename MonoIterator, typename StereoIterator, typename CollectorHelper>
159 void SiStripRecHitMatcher::doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend,
160  StereoIterator seconditer, StereoIterator seconditerend,
161  const GluedGeomDet* gluedDet, LocalVector trdir,
162  CollectorHelper & collectorHelper) const{
163 
164  using matcherDetails::StereoInfo;
165 
166  typedef GloballyPositioned<float> ToGlobal;
167  typedef typename GloballyPositioned<float>::ToLocal ToLocal;
168 
169  // hits in both mono and stero
170  // match
171  bool notk = trdir.mag2()<float(FLT_MIN);
172  // FIXME we shall find a faster approximation for trdir: not useful to compute it each time for each strip
173 
174  // stripdet = mono
175  // partnerstripdet = stereo
176  const GeomDetUnit* stripdet = gluedDet->monoDet();
177  const GeomDetUnit* partnerstripdet = gluedDet->stereoDet();
178  const StripTopology& topol=(const StripTopology&)stripdet->topology();
179  const StripTopology& partnertopol=(const StripTopology&)partnerstripdet->topology();
180 
181  // toGlobal is fast, toLocal is slow
182  ToGlobal const & stripDetTrans = stripdet->surface();
183  ToGlobal const & partnerStripDetTrans = partnerstripdet->surface();
184  ToLocal gluedDetInvTrans(gluedDet->surface());
185 
186 
187 
188  std::vector<StereoInfo> cache(std::distance(seconditer,seconditerend));
189  //iterate on stereo rechits
190  // fill cache with relevant info
191  int cacheSize=0;
192  for (;seconditer!=seconditerend; ++seconditer){
193 
194  const SiStripRecHit2D & secondHit = CollectorHelper::stereoHit(seconditer);
195 
196  float sigmap22 = sigmaPitch(secondHit.localPosition(),secondHit.localPositionError(),partnertopol);
197  // assert (sigmap22>=0);
198 
199  auto STEREOpointX=partnertopol.measurementPosition( secondHit.localPositionFast()).x();
200  MeasurementPoint STEREOpointini(STEREOpointX,-0.5f);
201  MeasurementPoint STEREOpointend(STEREOpointX,0.5f);
202 
203  LocalPoint locp1 = partnertopol.localPosition(STEREOpointini);
204  LocalPoint locp2 = partnertopol.localPosition(STEREOpointend);
205 
206  GlobalPoint globalpointini=partnerStripDetTrans.toGlobal(locp1);
207  GlobalPoint globalpointend=partnerStripDetTrans.toGlobal(locp2);
208 
209  // position of the initial and final point of the strip in glued local coordinates
210  LocalPoint positiononGluedini=gluedDetInvTrans.toLocal(globalpointini);
211  LocalPoint positiononGluedend=gluedDetInvTrans.toLocal(globalpointend);
212 
213  // in case of no track hypothesis assume a track from the origin through the center of the strip
214  if(notk){
215  LocalPoint lcenterofstrip=secondHit.localPositionFast();
216  GlobalPoint gcenterofstrip= partnerStripDetTrans.toGlobal(lcenterofstrip);
217  GlobalVector gtrackdirection=gcenterofstrip-GlobalPoint(0,0,0);
218  trdir=gluedDetInvTrans.toLocal(gtrackdirection);
219  }
220 
221 
222  Vec3F offset = trdir.basicVector().v * positiononGluedini.z()/trdir.z();
223 
224 
225  Vec3F ret1 = positiononGluedini.basicVector().v - offset;
226  Vec3F ret2 = positiononGluedend.basicVector().v - offset;
227 
228  double m10=-(ret2[1] - ret1[1]);
229  double m11= ret2[0] - ret1[0];
230  double dd = m11*ret1[1] + m10 * ret1[0];
231 
232  Vec2D c1vec{dd,dd};
233 
234  // store
235  StereoInfo info = {c1vec,&secondHit,sigmap22,m10,m11};
236  cache[cacheSize++] = info;
237  }
238 
239 
240 
241  for (;monoRHiter != monoRHend; ++monoRHiter) {
242 
243  SiStripRecHit2D const & monoRH = CollectorHelper::monoHit(monoRHiter);
244 
245  // position of the initial and final point of the strip (RPHI cluster) in local strip coordinates
246  auto RPHIpointX = topol.measurementPosition(monoRH.localPositionFast()).x();
247  MeasurementPoint RPHIpointini(RPHIpointX,-0.5f);
248  MeasurementPoint RPHIpointend(RPHIpointX,0.5f);
249 
250  // position of the initial and final point of the strip in local coordinates (mono det)
251  //StripPosition stripmono=StripPosition(topol.localPosition(RPHIpointini),topol.localPosition(RPHIpointend));
252  LocalPoint locp1o = topol.localPosition(RPHIpointini);
253  LocalPoint locp2o = topol.localPosition(RPHIpointend);
254 
255 
256  // in case of no track hypothesis assume a track from the origin through the center of the strip
257  if(notk){
258  LocalPoint lcenterofstrip=monoRH.localPositionFast();
259  GlobalPoint gcenterofstrip= stripDetTrans.toGlobal(lcenterofstrip);
260  GlobalVector gtrackdirection=gcenterofstrip-GlobalPoint(0,0,0);
261  trdir=gluedDetInvTrans.toLocal(gtrackdirection);
262  }
263 
264 
265  //project mono hit on glued det
266  //StripPosition projectedstripmono=project(stripdet,gluedDet,stripmono,trackdirection);
267 
268 
269  GlobalPoint globalpointini=stripDetTrans.toGlobal(locp1o);
270  GlobalPoint globalpointend=stripDetTrans.toGlobal(locp2o);
271 
272  // position of the initial and final point of the strip in glued local coordinates
273  LocalPoint positiononGluedini=gluedDetInvTrans.toLocal(globalpointini);
274  LocalPoint positiononGluedend=gluedDetInvTrans.toLocal(globalpointend);
275 
276  Vec3F offset = trdir.basicVector().v * positiononGluedini.z()/trdir.z();
277 
278 
279  Vec3F projini= positiononGluedini.basicVector().v - offset;
280  Vec3F projend = positiononGluedend.basicVector().v -offset;
281 
282  // ret1o = ret1o + (trdir * (ret1o.getSimd(2) / trdirz));
283  // ret2o = ret2o + (trdir * (ret2o.getSimd(2) / trdirz));
284 
285  double m00 = -(projend[1] - projini[1]); //-(projectedstripmono.second.y()-projectedstripmono.first.y());
286  double m01 = (projend[0] - projini[0]); // (projectedstripmono.second.x()-projectedstripmono.first.x());
287  double c0 = m01*projini[1] + m00*projini[0]; //m01*projectedstripmono.first.y() + m00*projectedstripmono.first.x();
288 
289  Vec2D c0vec{c0,c0};
290  Vec2D minv00{-m01, m00};
291 
292  //error calculation (the part that depends on mono RH only)
293  double c1 = -m00;
294  double s1 = -m01;
295  double l1 = 1./(c1*c1+s1*s1);
296 
297  float sigmap12 = sigmaPitch(monoRH.localPosition(), monoRH.localPositionError(),topol);
298  // float sigmap12 = monoRH.sigmaPitch();
299  // assert(sigmap12>=0);
300 
301  //float code
302  float fc1(c1), fs1(s1);
303  Vec3F scc1{fs1, fc1, fc1, 0.f};
304  Vec3F ssc1{fs1, fs1, fc1, 0.f};
305  const Vec3F cslsimd = scc1 * ssc1 * float(l1);
306 
307  for (int i=0; i!=cacheSize; ++i) {
308  StereoInfo const si = cache[i];
309 
310  // match
311  Vec2D minv10{si.m11, -si.m10};
312  double mult = 1./(m00*si.m11 - m01*si.m10);
313  Vec2D resultmatmul = mult * (minv10 * c0vec + minv00 * si.c1vec);
314 
315  Local2DPoint position(resultmatmul[0], resultmatmul[1]);
316 
317  // LocalError tempError (100,0,100);
318  if (!((gluedDet->surface()).bounds().inside(position,10.f*scale_))) continue;
319 
320  double c2 = -si.m10;
321  double s2 = -si.m11;
322  double l2 = 1./(c2*c2+s2*s2);
323 
324  double diff=(c1*s2-c2*s1);
325  double invdet2 = 1./(diff*diff*l1*l2);
326 
327  float fc2(c2), fs2(s2), fid2(invdet2);
328  Vec3F invdet2simd{fid2, -fid2, fid2, 0.f};
329  Vec3F ccssimd{fs2, fc2, fc2, 0.f};
330  Vec3F csssimd{fs2, fs2, fc2, 0.f};
331  Vec3F result = invdet2simd * (si.sigmap22 * cslsimd + sigmap12 * ccssimd * csssimd * float(l2));
332 
333 
334  LocalError error(result[0], result[1], result[2]);
335 
336 
337  if((gluedDet->surface()).bounds().inside(position,error,scale_)){ //if it is inside the gluedet bonds
338 
339  //Change NSigmaInside in the configuration file to accept more hits
340  //...and add it to the Rechit collection
341 
342  collectorHelper.collector()(SiStripMatchedRecHit2D(LocalPoint(position), error,
343  *gluedDet,&monoRH,si.secondHit));
344  }
345 
346  } // loop on cache info
347 
348  collectorHelper.closure(monoRHiter);
349  } // loop on mono hit
350 
351 }
352 
353 #endif //DOUBLE_MATCH
354 
355 
356 
357 #endif // RECOLOCALTRACKER_SISTRIPCLUSTERIZER_SISTRIPRECHITMATCH_H
358 
359 
360 
int i
Definition: DBlmapReader.cc:9
T mag2() const
Definition: PV3DBase.h:66
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
static const TGPicture * info(bool iBackgroundIsBlack)
Vec4< float > Vec3F
Definition: SSEVec.h:501
const GeomDetUnit * monoDet() const
Definition: GluedGeomDet.h:20
virtual const Topology & topology() const
Definition: GeomDet.cc:86
StripPosition project(const GeomDetUnit *det, const GluedGeomDet *glueddet, StripPosition strip, LocalVector trackdirection) const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Vec4< float > Vec3F
Definition: ExtVec.h:174
static float sigmaPitch(LocalPoint const &pos, LocalError const &err, const StripTopology &topol)
list diff
Definition: mps_update.py:85
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
tuple s2
Definition: indexGen.py:106
tuple result
Definition: mps_fire.py:83
virtual float localPitch(const LocalPoint &) const =0
virtual MeasurementError measurementError(const LocalPoint &, const LocalError &) const =0
std::pair< LocalPoint, LocalPoint > StripPosition
SiStripRecHitMatcher(const edm::ParameterSet &conf)
T z() const
Definition: PV3DBase.h:64
float uu() const
double f[11][100]
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
Vec4< double > Vec3D
Definition: SSEVec.h:503
#define end
Definition: vmac.h:37
Rot3< float > Rot3F
Definition: SSERot.h:73
void doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend, StereoIterator seconditer, StereoIterator seconditerend, const GluedGeomDet *gluedDet, LocalVector trdir, CollectorHelper &collectorHelper) const
SiStripMatchedRecHit2DCollectionNew::FastFiller CollectorMatched
std::vector< const SiStripRecHit2D * > SimpleHitCollection
std::unique_ptr< SiStripMatchedRecHit2D > match(const SiStripRecHit2D *monoRH, const SiStripRecHit2D *stereoRH, const GluedGeomDet *gluedDet, LocalVector trackdirection, bool force=false) const
SiStripRecHit2DCollectionNew::DetSet::const_iterator RecHitIterator
#define begin
Definition: vmac.h:30
SimpleHitCollection::const_iterator SimpleHitIterator
static int position[264][3]
Definition: ReadPGInfo.cc:509
Vec2< double > Vec2D
Definition: ExtVec.h:175
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
virtual LocalPoint localPosition(float strip) const =0
const LocalPoint & localPositionFast() const
boost::function< void(SiStripMatchedRecHit2D const &)> Collector
virtual LocalError localPositionError() const final
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
virtual LocalPoint localPosition() const final
const GeomDetUnit * stereoDet() const
Definition: GluedGeomDet.h:21
Vec2< double > Vec2D
Definition: SSEVec.h:502