00001
00002
00003
00004 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripRecHitMatcher.h"
00005 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00006 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00007 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00008 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00009
00010 #include "TrackingTools/TransientTrackingRecHit/interface/HelpertRecHit2DLocalPos.h"
00011 #include<boost/bind.hpp>
00012
00013 #include <DataFormats/TrackingRecHit/interface/AlignmentPositionError.h>
00014
00015
00016
00017
00018 SiStripRecHitMatcher::SiStripRecHitMatcher(const edm::ParameterSet& conf){
00019 scale_=conf.getParameter<double>("NSigmaInside");
00020 }
00021
00022 SiStripRecHitMatcher::SiStripRecHitMatcher(const double theScale){
00023 scale_=theScale;
00024 }
00025
00026
00027
00028 SiStripMatchedRecHit2D * SiStripRecHitMatcher::match(const SiStripRecHit2D *monoRH,
00029 const SiStripRecHit2D *stereoRH,
00030 const GluedGeomDet* gluedDet,
00031 LocalVector trackdirection) const{
00032 SimpleHitCollection stereoHits(1,stereoRH);
00033 std::vector<SiStripMatchedRecHit2D*> collection;
00034 match( monoRH,
00035 stereoHits.begin(), stereoHits.end(),
00036 collection,
00037 gluedDet,trackdirection);
00038
00039 return collection.empty() ? (SiStripMatchedRecHit2D*)(0) : collection.front();
00040 }
00041
00042
00043 SiStripMatchedRecHit2D* SiStripRecHitMatcher::match(const SiStripMatchedRecHit2D *origRH,
00044 const GluedGeomDet* gluedDet,
00045 LocalVector trackdirection) const{
00046
00047 throw "SiStripRecHitMatcher::match(const SiStripMatchedRecHit2D *,..) is obsoltete since 5.2.0";
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064 return nullptr;
00065 }
00066
00067
00068 edm::OwnVector<SiStripMatchedRecHit2D>
00069 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00070 RecHitIterator begin, RecHitIterator end,
00071 const GluedGeomDet* gluedDet,
00072 LocalVector trackdirection) const
00073 {
00074 SimpleHitCollection stereoHits;
00075 stereoHits.reserve(end-begin);
00076
00077 for (RecHitIterator i=begin; i != end; ++i) {
00078 stereoHits.push_back( &(*i));
00079 }
00080 return match( monoRH,
00081 stereoHits.begin(), stereoHits.end(),
00082 gluedDet,trackdirection);
00083 }
00084
00085 edm::OwnVector<SiStripMatchedRecHit2D>
00086 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00087 SimpleHitIterator begin, SimpleHitIterator end,
00088 const GluedGeomDet* gluedDet,
00089 LocalVector trackdirection) const {
00090 edm::OwnVector<SiStripMatchedRecHit2D> collector;
00091 collector.reserve(end-begin);
00092 match(monoRH,begin,end,collector,gluedDet,trackdirection);
00093 return collector;
00094 }
00095
00096
00097 void
00098 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00099 SimpleHitIterator begin, SimpleHitIterator end,
00100 edm::OwnVector<SiStripMatchedRecHit2D> & collector,
00101 const GluedGeomDet* gluedDet,
00102 LocalVector trackdirection) const {
00103
00104 std::vector<SiStripMatchedRecHit2D*> result;
00105 result.reserve(end-begin);
00106 match(monoRH,begin,end,result,gluedDet,trackdirection);
00107 for (std::vector<SiStripMatchedRecHit2D*>::iterator p=result.begin(); p!=result.end();
00108 p++) collector.push_back(*p);
00109 }
00110
00111 namespace {
00112
00113 inline void pb(std::vector<SiStripMatchedRecHit2D*> & v, SiStripMatchedRecHit2D* h) {
00114 v.push_back(h);
00115 }
00116 }
00117
00118 void
00119 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00120 SimpleHitIterator begin, SimpleHitIterator end,
00121 std::vector<SiStripMatchedRecHit2D*> & collector,
00122 const GluedGeomDet* gluedDet,
00123 LocalVector trackdirection) const {
00124 Collector result(boost::bind(&pb,boost::ref(collector),
00125 boost::bind(&SiStripMatchedRecHit2D::clone,_1)));
00126 match(monoRH,begin,end,result,gluedDet,trackdirection);
00127 }
00128
00129
00130
00131 void
00132 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00133 RecHitIterator begin, RecHitIterator end,
00134 CollectorMatched & collector,
00135 const GluedGeomDet* gluedDet,
00136 LocalVector trackdirection) const {
00137
00138
00139 SimpleHitCollection stereoHits;
00140 stereoHits.reserve(end-begin);
00141 for (RecHitIterator i=begin; i != end; ++i)
00142 stereoHits.push_back( &(*i));
00143
00144 return match( monoRH,
00145 stereoHits.begin(), stereoHits.end(),
00146 collector,
00147 gluedDet,trackdirection);
00148 }
00149
00150 void
00151 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00152 SimpleHitIterator begin, SimpleHitIterator end,
00153 CollectorMatched & collector,
00154 const GluedGeomDet* gluedDet,
00155 LocalVector trackdirection) const {
00156
00157 Collector result(boost::bind(&CollectorMatched::push_back,boost::ref(collector),_1));
00158 match(monoRH,begin,end,result,gluedDet,trackdirection);
00159
00160 }
00161
00162
00163
00164 void
00165 SiStripRecHitMatcher::match( const SiStripRecHit2D *monoRH,
00166 SimpleHitIterator begin, SimpleHitIterator end,
00167 Collector & collector,
00168 const GluedGeomDet* gluedDet,
00169 LocalVector trackdirection) const {
00170
00171
00172 const GeomDetUnit* stripdet = gluedDet->monoDet();
00173 const GeomDetUnit* partnerstripdet = gluedDet->stereoDet();
00174 const StripTopology& topol=(const StripTopology&)stripdet->topology();
00175
00176
00177 double RPHIpointX = topol.measurementPosition(monoRH->localPositionFast()).x();
00178 MeasurementPoint RPHIpointini(RPHIpointX,-0.5);
00179 MeasurementPoint RPHIpointend(RPHIpointX,0.5);
00180
00181
00182 StripPosition stripmono=StripPosition(topol.localPosition(RPHIpointini),topol.localPosition(RPHIpointend));
00183
00184 if(trackdirection.mag2()<FLT_MIN){
00185 LocalPoint lcenterofstrip=monoRH->localPositionFast();
00186 GlobalPoint gcenterofstrip=(stripdet->surface()).toGlobal(lcenterofstrip);
00187 GlobalVector gtrackdirection=gcenterofstrip-GlobalPoint(0,0,0);
00188 trackdirection=(gluedDet->surface()).toLocal(gtrackdirection);
00189 }
00190
00191
00192 StripPosition projectedstripmono=project(stripdet,gluedDet,stripmono,trackdirection);
00193 const StripTopology& partnertopol=(const StripTopology&)partnerstripdet->topology();
00194
00195 double m00 = -(projectedstripmono.second.y()-projectedstripmono.first.y());
00196 double m01 = (projectedstripmono.second.x()-projectedstripmono.first.x());
00197 double c0 = m01*projectedstripmono.first.y() + m00*projectedstripmono.first.x();
00198
00199
00200
00201
00202
00203
00204
00205
00206 double c1 = -m00;
00207 double s1 = -m01;
00208 double l1 = 1./(c1*c1+s1*s1);
00209
00210
00211
00212 double sigmap12 = monoRH->sigmaPitch();
00213 if (sigmap12<0) {
00214
00215
00216
00217
00218
00219
00220
00221
00222 LocalError tmpError(monoRH->localPositionErrorFast());
00223 HelpertRecHit2DLocalPos::updateWithAPE(tmpError,*stripdet);
00224 MeasurementError errormonoRH=topol.measurementError(monoRH->localPositionFast(),tmpError);
00225
00226
00227
00228
00229
00230
00231
00232 double pitch=topol.localPitch(monoRH->localPositionFast());
00233 monoRH->setSigmaPitch(sigmap12=errormonoRH.uu()*pitch*pitch);
00234 }
00235
00236 SimpleHitIterator seconditer;
00237
00238 for(seconditer=begin;seconditer!=end;++seconditer){
00239
00240
00241 double STEREOpointX=partnertopol.measurementPosition((*seconditer)->localPositionFast()).x();
00242 MeasurementPoint STEREOpointini(STEREOpointX,-0.5);
00243 MeasurementPoint STEREOpointend(STEREOpointX,0.5);
00244
00245
00246 StripPosition stripstereo(partnertopol.localPosition(STEREOpointini),partnertopol.localPosition(STEREOpointend));
00247
00248
00249 StripPosition projectedstripstereo=project(partnerstripdet,gluedDet,stripstereo,trackdirection);
00250
00251
00252 double m10=-(projectedstripstereo.second.y()-projectedstripstereo.first.y());
00253 double m11=(projectedstripstereo.second.x()-projectedstripstereo.first.x());
00254
00255
00256
00257 AlgebraicMatrix22 m; AlgebraicVector2 c;
00258 m(0,0)=m00;
00259 m(0,1)=m01;
00260 m(1,0)=m10;
00261 m(1,1)=m11;
00262 c(0)=c0;
00263 c(1)=m11*projectedstripstereo.first.y()+m10*projectedstripstereo.first.x();
00264 m.Invert();
00265 AlgebraicVector2 solution = m * c;
00266 LocalPoint position(solution(0),solution(1));
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286 LocalError tempError (100,0,100);
00287 if (!((gluedDet->surface()).bounds().inside(position,tempError,scale_))) continue;
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297 double c2 = -m10;
00298 double s2 = -m11;
00299 double l2 = 1./(c2*c2+s2*s2);
00300
00301
00302
00303 double sigmap22 = (*seconditer)->sigmaPitch();
00304 if (sigmap22<0) {
00305
00306 LocalError tmpError((*seconditer)->localPositionErrorFast());
00307 HelpertRecHit2DLocalPos::updateWithAPE(tmpError, *partnerstripdet);
00308 MeasurementError errorstereoRH=partnertopol.measurementError((*seconditer)->localPositionFast(),tmpError);
00309
00310 double pitch=partnertopol.localPitch((*seconditer)->localPositionFast());
00311 (*seconditer)->setSigmaPitch(sigmap22=errorstereoRH.uu()*pitch*pitch);
00312 }
00313
00314 double diff=(c1*s2-c2*s1);
00315 double invdet2=1/(diff*diff*l1*l2);
00316 float xx= invdet2*(sigmap12*s2*s2*l2+sigmap22*s1*s1*l1);
00317 float xy=-invdet2*(sigmap12*c2*s2*l2+sigmap22*c1*s1*l1);
00318 float yy= invdet2*(sigmap12*c2*c2*l2+sigmap22*c1*c1*l1);
00319 LocalError error(xx,xy,yy);
00320
00321 if((gluedDet->surface()).bounds().inside(position,error,scale_)){
00322
00323
00324
00325 const SiStripRecHit2D* secondHit = *seconditer;
00326 collector(SiStripMatchedRecHit2D(position, error,gluedDet->geographicalId() ,
00327 monoRH,secondHit));
00328 }
00329 }
00330 }
00331
00332
00333 SiStripRecHitMatcher::StripPosition SiStripRecHitMatcher::project(const GeomDetUnit *det,const GluedGeomDet* glueddet,StripPosition strip,LocalVector trackdirection)const
00334 {
00335
00336 GlobalPoint globalpointini=(det->surface()).toGlobal(strip.first);
00337 GlobalPoint globalpointend=(det->surface()).toGlobal(strip.second);
00338
00339
00340 LocalPoint positiononGluedini=(glueddet->surface()).toLocal(globalpointini);
00341 LocalPoint positiononGluedend=(glueddet->surface()).toLocal(globalpointend);
00342
00343
00344
00345 float scale=-positiononGluedini.z()/trackdirection.z();
00346
00347 LocalPoint projpositiononGluedini= positiononGluedini + scale*trackdirection;
00348 LocalPoint projpositiononGluedend= positiononGluedend + scale*trackdirection;
00349
00350 return StripPosition(projpositiononGluedini,projpositiononGluedend);
00351 }