CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/RecoLocalTracker/SiStripRecHitConverter/interface/SiStripRecHitMatcher.h

Go to the documentation of this file.
00001 #ifndef RECOLOCALTRACKER_SISTRIPCLUSTERIZER_SISTRIPRECHITMATCH_H
00002 #define RECOLOCALTRACKER_SISTRIPCLUSTERIZER_SISTRIPRECHITMATCH_H
00003 
00004 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00005 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00006 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00007 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00008 #include "DataFormats/DetId/interface/DetId.h"
00009 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00010 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00011 
00012 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00013 
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 
00016 class GluedGeomDet;
00017 
00018 #include <cfloat>
00019 
00020 
00021 #include <boost/function.hpp>
00022 
00023 class SiStripRecHitMatcher {
00024 public:
00025   
00026   // may become a template argument
00027   typedef SiStripMatchedRecHit2DCollectionNew::FastFiller CollectorMatched;
00028 
00029   typedef SiStripRecHit2DCollectionNew::DetSet::const_iterator RecHitIterator;
00030   typedef std::vector<const SiStripRecHit2D *>              SimpleHitCollection;
00031   typedef SimpleHitCollection::const_iterator               SimpleHitIterator;
00032 
00033   typedef boost::function<void(SiStripMatchedRecHit2D const&)>    Collector;
00034 
00035 
00036   typedef std::pair<LocalPoint,LocalPoint>                  StripPosition; 
00037 
00038   SiStripRecHitMatcher(const edm::ParameterSet& conf);
00039   SiStripRecHitMatcher(const double theScale);
00040   
00041 
00042   // optimized matching iteration (the algo is the same, just recoded)
00043   template<typename MonoIterator, typename StereoIterator,  typename CollectorHelper>
00044   void doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend,
00045                    StereoIterator seconditer, StereoIterator seconditerend,
00046                    const GluedGeomDet* gluedDet,  LocalVector trdir, 
00047                    CollectorHelper & collectorHelper) const;
00048   
00049   
00050   
00051   SiStripMatchedRecHit2D * match(const SiStripRecHit2D *monoRH, 
00052                                  const SiStripRecHit2D *stereoRH,
00053                                  const GluedGeomDet* gluedDet,
00054                                  LocalVector trackdirection) const;
00055   
00056   SiStripMatchedRecHit2D*  match(const SiStripMatchedRecHit2D *originalRH, 
00057                                  const GluedGeomDet* gluedDet,
00058                                  LocalVector trackdirection) const;
00059   
00060   edm::OwnVector<SiStripMatchedRecHit2D> 
00061   match( const SiStripRecHit2D *monoRH,
00062          RecHitIterator begin, RecHitIterator end, 
00063          const GluedGeomDet* gluedDet) const {
00064     return match(monoRH,begin, end, gluedDet,LocalVector(0.,0.,0.));
00065   }
00066   
00067   edm::OwnVector<SiStripMatchedRecHit2D> 
00068   match( const SiStripRecHit2D *monoRH,
00069          RecHitIterator begin, RecHitIterator end,
00070          const GluedGeomDet* gluedDet,
00071          LocalVector trackdirection) const;
00072   
00073   edm::OwnVector<SiStripMatchedRecHit2D> 
00074   match( const SiStripRecHit2D *monoRH,
00075          SimpleHitIterator begin, SimpleHitIterator end,
00076          const GluedGeomDet* gluedDet,
00077          LocalVector trackdirection) const;
00078   
00079   void
00080   match( const SiStripRecHit2D *monoRH,
00081          RecHitIterator begin, RecHitIterator end,
00082          CollectorMatched & collector,
00083          const GluedGeomDet* gluedDet,
00084          LocalVector trackdirection) const;
00085   
00086   void
00087   match( const SiStripRecHit2D *monoRH,
00088          SimpleHitIterator begin, SimpleHitIterator end,
00089          CollectorMatched & collector,
00090          const GluedGeomDet* gluedDet,
00091          LocalVector trackdirection) const;
00092   
00093   
00094   
00095   
00096   // project strip coordinates on Glueddet
00097   
00098   StripPosition project(const GeomDetUnit *det,const GluedGeomDet* glueddet,StripPosition strip,LocalVector trackdirection) const;
00099   
00100   
00101   //private:
00102   
00103   
00104   void
00105   match( const SiStripRecHit2D *monoRH,
00106          SimpleHitIterator begin, SimpleHitIterator end,
00107          edm::OwnVector<SiStripMatchedRecHit2D> & collector, 
00108          const GluedGeomDet* gluedDet,
00109          LocalVector trackdirection) const;
00110   
00111   
00112   void
00113   match( const SiStripRecHit2D *monoRH,
00114          SimpleHitIterator begin, SimpleHitIterator end,
00115          std::vector<SiStripMatchedRecHit2D*> & collector, 
00116          const GluedGeomDet* gluedDet,
00117          LocalVector trackdirection) const;
00118   
00119   
00121   
00122   void
00123   match( const SiStripRecHit2D *monoRH,
00124          SimpleHitIterator begin, SimpleHitIterator end,
00125          Collector & collector, 
00126          const GluedGeomDet* gluedDet,
00127          LocalVector trackdirection) const;
00128   
00129   
00130   float scale_;
00131   
00132 };
00133 
00134 
00135 
00136 
00137 #include "DataFormats/GeometryVector/interface/Basic3DVector.h"
00138 #include "DataFormats/Math/interface/SSEVec.h"
00139 #ifdef USE_SSEVECT
00140 #define DOUBLE_MATCH
00141 #endif
00142 
00143 #ifdef DOUBLE_MATCH
00144 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00145 #include "TrackingTools/TransientTrackingRecHit/interface/HelpertRecHit2DLocalPos.h"
00146 
00147 namespace matcherDetails {  
00148   
00149   struct StereoInfo {
00150     mathSSE::Vec2D c1vec;
00151     const SiStripRecHit2D * secondHit;
00152     double sigmap22;
00153     double m10, m11;
00154   };
00155   
00156 }
00157 
00158 template<typename MonoIterator, typename StereoIterator,  typename CollectorHelper>
00159 void SiStripRecHitMatcher::doubleMatch(MonoIterator monoRHiter, MonoIterator monoRHend,
00160                                        StereoIterator seconditer, StereoIterator seconditerend,
00161                                        const GluedGeomDet* gluedDet,    LocalVector trdir, 
00162                                        CollectorHelper & collectorHelper) const{
00163   
00164   using  matcherDetails::StereoInfo;  
00165   using  mathSSE::Vec3F;
00166   using  mathSSE::Vec2D;
00167   using  mathSSE::Vec3D;
00168   using  mathSSE::Rot3F;
00169   typedef  GloballyPositioned<float> ToGlobal;
00170   typedef  typename GloballyPositioned<float>::ToLocal ToLocal;
00171   
00172   // hits in both mono and stero
00173   // match
00174   bool notk = trdir.mag2()<FLT_MIN;
00175   // FIXME we shall find a faster approximation for trdir: not useful to compute it each time for each strip
00176   
00177   // stripdet = mono
00178   // partnerstripdet = stereo
00179   const GeomDetUnit* stripdet = gluedDet->monoDet();
00180   const GeomDetUnit* partnerstripdet = gluedDet->stereoDet();
00181   const StripTopology& topol=(const StripTopology&)stripdet->topology();
00182   const StripTopology& partnertopol=(const StripTopology&)partnerstripdet->topology();
00183   
00184   // toGlobal is fast,  toLocal is slow
00185   ToGlobal const & stripDetTrans =  stripdet->surface();
00186   ToGlobal const & partnerStripDetTrans = partnerstripdet->surface();
00187   ToLocal          gluedDetInvTrans(gluedDet->surface());
00188   
00189   
00190   
00191   StereoInfo cache[std::distance(seconditer,seconditerend)];
00192   //iterate on stereo rechits
00193   // fill cache with relevant info
00194   int  cacheSize=0;
00195   for (;seconditer!=seconditerend; ++seconditer){
00196     
00197     const SiStripRecHit2D & secondHit = CollectorHelper::stereoHit(seconditer);
00198     
00199     double sigmap22 =secondHit.sigmaPitch();
00200     if (sigmap22<0) {
00201       LocalError tmpError( secondHit.localPositionErrorFast());
00202       HelpertRecHit2DLocalPos::updateWithAPE(tmpError, *partnerstripdet);
00203       MeasurementError errorstereoRH=partnertopol.measurementError(secondHit.localPositionFast(),tmpError);
00204       
00205       double pitch=partnertopol.localPitch(secondHit.localPositionFast());
00206       secondHit.setSigmaPitch(sigmap22=errorstereoRH.uu()*pitch*pitch);
00207     }
00208     
00209     
00210     double STEREOpointX=partnertopol.measurementPosition( secondHit.localPositionFast()).x();
00211     MeasurementPoint STEREOpointini(STEREOpointX,-0.5);
00212     MeasurementPoint STEREOpointend(STEREOpointX,0.5);
00213     
00214     LocalPoint locp1 = partnertopol.localPosition(STEREOpointini);
00215     LocalPoint locp2 = partnertopol.localPosition(STEREOpointend);
00216     
00217     GlobalPoint globalpointini=partnerStripDetTrans.toGlobal(locp1);
00218     GlobalPoint globalpointend=partnerStripDetTrans.toGlobal(locp2);
00219     
00220     // position of the initial and final point of the strip in glued local coordinates
00221     LocalPoint positiononGluedini=gluedDetInvTrans.toLocal(globalpointini);
00222     LocalPoint positiononGluedend=gluedDetInvTrans.toLocal(globalpointend); 
00223     
00224    // in case of no track hypothesis assume a track from the origin through the center of the strip
00225     if(notk){
00226       LocalPoint lcenterofstrip=secondHit.localPositionFast();
00227       GlobalPoint gcenterofstrip= partnerStripDetTrans.toGlobal(lcenterofstrip);
00228       GlobalVector gtrackdirection=gcenterofstrip-GlobalPoint(0,0,0);
00229       trdir=gluedDetInvTrans.toLocal(gtrackdirection);
00230     }
00231   
00232 
00233     Vec3F offset = trdir.basicVector().v * positiononGluedini.basicVector().v.get1<2>()/trdir.basicVector().v.get1<2>();
00234     
00235     
00236     Vec3F ret1 = positiononGluedini.basicVector().v - offset;
00237     Vec3F ret2 = positiononGluedend.basicVector().v - offset;
00238     
00239     double m10=-(ret2.arr[1] - ret1.arr[1]); 
00240     double m11=  ret2.arr[0] - ret1.arr[0];
00241     
00242     Vec2D c1vec; c1vec.set1(m11*ret1.arr[1] + m10 * ret1.arr[0]);
00243     
00244     // store
00245     StereoInfo info = {c1vec,&secondHit,sigmap22,m10,m11};
00246     cache[cacheSize++] = info;
00247   }
00248   
00249   
00250   
00251   for (;monoRHiter != monoRHend; ++monoRHiter) {
00252     
00253     SiStripRecHit2D const & monoRH = CollectorHelper::monoHit(monoRHiter);
00254     
00255     // position of the initial and final point of the strip (RPHI cluster) in local strip coordinates
00256     double RPHIpointX = topol.measurementPosition(monoRH.localPositionFast()).x();
00257     MeasurementPoint RPHIpointini(RPHIpointX,-0.5);
00258     MeasurementPoint RPHIpointend(RPHIpointX,0.5);
00259     
00260     // position of the initial and final point of the strip in local coordinates (mono det)
00261     //StripPosition stripmono=StripPosition(topol.localPosition(RPHIpointini),topol.localPosition(RPHIpointend));
00262     LocalPoint locp1o = topol.localPosition(RPHIpointini);
00263     LocalPoint locp2o = topol.localPosition(RPHIpointend);
00264     
00265     
00266     // in case of no track hypothesis assume a track from the origin through the center of the strip
00267     if(notk){
00268       LocalPoint lcenterofstrip=monoRH.localPositionFast();
00269       GlobalPoint gcenterofstrip= stripDetTrans.toGlobal(lcenterofstrip);
00270       GlobalVector gtrackdirection=gcenterofstrip-GlobalPoint(0,0,0);
00271       trdir=gluedDetInvTrans.toLocal(gtrackdirection);
00272     }
00273     
00274     
00275     //project mono hit on glued det
00276     //StripPosition projectedstripmono=project(stripdet,gluedDet,stripmono,trackdirection);
00277     
00278     
00279     GlobalPoint globalpointini=stripDetTrans.toGlobal(locp1o);
00280     GlobalPoint globalpointend=stripDetTrans.toGlobal(locp2o);
00281     
00282     // position of the initial and final point of the strip in glued local coordinates
00283     LocalPoint positiononGluedini=gluedDetInvTrans.toLocal(globalpointini);
00284     LocalPoint positiononGluedend=gluedDetInvTrans.toLocal(globalpointend); 
00285     
00286     Vec3F offset = trdir.basicVector().v * positiononGluedini.basicVector().v.get1<2>()/trdir.basicVector().v.get1<2>();
00287     
00288     
00289     Vec3F projini= positiononGluedini.basicVector().v - offset;
00290     Vec3F projend = positiononGluedend.basicVector().v -offset;
00291     
00292     // ret1o = ret1o + (trdir * (ret1o.getSimd(2) / trdirz));
00293     // ret2o = ret2o + (trdir * (ret2o.getSimd(2) / trdirz));
00294     
00295     double m00 = -(projend.arr[1] - projini.arr[1]);//-(projectedstripmono.second.y()-projectedstripmono.first.y()); 
00296     double m01 =  (projend.arr[0] - projini.arr[0]); // (projectedstripmono.second.x()-projectedstripmono.first.x());
00297     double c0  =  m01*projini.arr[1] + m00*projini.arr[0];//m01*projectedstripmono.first.y()   + m00*projectedstripmono.first.x();
00298     
00299     Vec2D c0vec(c0,c0);
00300     Vec2D minv00(-m01, m00);
00301     
00302     //error calculation (the part that depends on mono RH only)
00303     double c1 = -m00;
00304     double s1 = -m01;
00305     double l1 = 1./(c1*c1+s1*s1);
00306     
00307     // FIXME: here for test...
00308     double sigmap12 = monoRH.sigmaPitch();
00309     if (sigmap12<0) {
00310       
00311       LocalError tmpError(monoRH.localPositionErrorFast());
00312       HelpertRecHit2DLocalPos::updateWithAPE(tmpError,*stripdet);
00313       MeasurementError errormonoRH=topol.measurementError(monoRH.localPositionFast(),tmpError);
00314       
00315       double pitch=topol.localPitch(monoRH.localPositionFast());
00316       monoRH.setSigmaPitch(sigmap12=errormonoRH.uu()*pitch*pitch);
00317     }
00318 
00319     //float code
00320     float fc1(c1), fs1(s1);
00321     Vec3F scc1(fs1, fc1, fc1, 0.f);
00322     Vec3F ssc1(fs1, fs1, fc1, 0.f);
00323     Vec3F l1vec; l1vec.set1(l1);
00324     const Vec3F cslsimd = scc1 * ssc1 * l1vec;
00325     Vec3F sigmap12simd; sigmap12simd.set1(sigmap12);
00326     
00327     for (int i=0; i!=cacheSize; ++i) {
00328       StereoInfo const si = cache[i];
00329       
00330       // match 
00331       Vec2D minv10(si.m11, -si.m10);
00332       Vec2D mult; mult.set1(1./(m00*si.m11 - m01*si.m10));
00333       Vec2D resultmatmul = mult * (minv10 * c0vec + minv00 * si.c1vec);
00334       
00335       Local2DPoint position(resultmatmul.arr[0], resultmatmul.arr[1]);
00336       
00337       // LocalError tempError (100,0,100);
00338       if (!((gluedDet->surface()).bounds().inside(position,10.f*scale_))) continue;                                                       
00339       
00340       double c2 = -si.m10;
00341       double s2 = -si.m11;
00342       double l2 = 1./(c2*c2+s2*s2);
00343       
00344       double diff=(c1*s2-c2*s1);
00345       double invdet2 = 1./(diff*diff*l1*l2);
00346       
00347       float fc2(c2), fs2(s2), fid2(invdet2);    
00348       Vec3F invdet2simd(fid2, -fid2, fid2, 0.f);
00349       Vec3F ccssimd(fs2, fc2, fc2, 0.f);
00350       Vec3F csssimd(fs2, fs2, fc2, 0.f);
00351       Vec3F l2simd; l2simd.set1(l2);
00352       Vec3F sigmap22simd; sigmap22simd.set1(si.sigmap22);
00353       Vec3F result = invdet2simd * (sigmap22simd * cslsimd + sigmap12simd * ccssimd * csssimd * l2simd);
00354       
00355       
00356       LocalError error(result.arr[0], result.arr[1], result.arr[2]);
00357       
00358       
00359       if((gluedDet->surface()).bounds().inside(position,error,scale_)){ //if it is inside the gluedet bonds
00360         
00361         //Change NSigmaInside in the configuration file to accept more hits
00362         //...and add it to the Rechit collection 
00363         
00364         collectorHelper.collector()(SiStripMatchedRecHit2D(LocalPoint(position), error,gluedDet->geographicalId() ,
00365                                                            &monoRH,si.secondHit));
00366       }
00367       
00368     } // loop on cache info
00369     
00370     collectorHelper.closure(monoRHiter);
00371   } // loop on mono hit
00372   
00373 }
00374 
00375 #endif //DOUBLE_MATCH
00376 
00377 
00378 
00379 #endif  //  RECOLOCALTRACKER_SISTRIPCLUSTERIZER_SISTRIPRECHITMATCH_H
00380 
00381 
00382