CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoEgamma/EgammaElectronAlgos/src/PixelHitMatcher.cc

Go to the documentation of this file.
00001 
00002 #include "RecoEgamma/EgammaElectronAlgos/interface/PixelHitMatcher.h"
00003 #include "RecoEgamma/EgammaElectronAlgos/interface/PixelMatchNextLayers.h"
00004 #include "RecoEgamma/EgammaElectronAlgos/interface/ElectronUtilities.h"
00005 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00006 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00007 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00008 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00009 #include "DataFormats/DetId/interface/DetId.h"
00010 #include "DataFormats/GeometryCommonDetAlgo/interface/PerpendicularBoundPlaneBuilder.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 
00013 #include <typeinfo>
00014 
00015 using namespace reco ;
00016 using namespace std ;
00017 
00018 PixelHitMatcher::PixelHitMatcher
00019  ( float phi1min, float phi1max,
00020    float phi2minB, float phi2maxB, float phi2minF, float phi2maxF,
00021    float z2minB, float z2maxB, float r2minF, float r2maxF,
00022    float rMinI, float rMaxI, bool searchInTIDTEC)
00023  : //zmin1 and zmax1 are dummy at this moment, set from beamspot later
00024    meas1stBLayer(phi1min,phi1max,0.,0.), meas2ndBLayer(phi2minB,phi2maxB,z2minB,z2maxB),
00025    meas1stFLayer(phi1min,phi1max,0.,0.), meas2ndFLayer(phi2minF,phi2maxF,r2minF,r2maxF),
00026    startLayers(),
00027    prop1stLayer(0), prop2ndLayer(0),theGeometricSearchTracker(0),theLayerMeasurements(0),vertex_(0.),
00028    searchInTIDTEC_(searchInTIDTEC), useRecoVertex_(false)
00029  {
00030   meas1stFLayer.setRRangeI(rMinI,rMaxI) ;
00031   meas2ndFLayer.setRRangeI(rMinI,rMaxI) ;
00032  }
00033 
00034 PixelHitMatcher::~PixelHitMatcher()
00035  {
00036   delete prop1stLayer ;
00037   delete prop2ndLayer ;
00038   delete theLayerMeasurements ;
00039  }
00040 
00041 void PixelHitMatcher::set1stLayer( float dummyphi1min, float dummyphi1max )
00042  {
00043   meas1stBLayer.setPhiRange(dummyphi1min,dummyphi1max) ;
00044   meas1stFLayer.setPhiRange(dummyphi1min,dummyphi1max) ;
00045  }
00046 
00047 void PixelHitMatcher::set1stLayerZRange( float zmin1, float zmax1 )
00048  {
00049   meas1stBLayer.setZRange(zmin1,zmax1) ;
00050   meas1stFLayer.setRRange(zmin1,zmax1) ;
00051  }
00052 
00053 void PixelHitMatcher::set2ndLayer( float dummyphi2minB, float dummyphi2maxB, float dummyphi2minF, float dummyphi2maxF )
00054  {
00055   meas2ndBLayer.setPhiRange(dummyphi2minB,dummyphi2maxB) ;
00056   meas2ndFLayer.setPhiRange(dummyphi2minF,dummyphi2maxF) ;
00057  }
00058 
00059 void PixelHitMatcher::setUseRecoVertex( bool val )
00060  { useRecoVertex_ = val ; }
00061 
00062 void PixelHitMatcher::setES
00063  ( const MagneticField * magField,
00064    const MeasurementTracker * theMeasurementTracker,
00065    const TrackerGeometry * trackerGeometry )
00066  {
00067   if (theMeasurementTracker)
00068    {
00069     theGeometricSearchTracker=theMeasurementTracker->geometricSearchTracker() ;
00070     startLayers.setup(theGeometricSearchTracker) ;
00071     if (theLayerMeasurements ) delete theLayerMeasurements ;
00072     theLayerMeasurements = new LayerMeasurements(theMeasurementTracker) ;
00073    }
00074 
00075   theMagField = magField ;
00076   theTrackerGeometry = trackerGeometry ;
00077   float mass=.000511 ; // electron propagation
00078   if (prop1stLayer) delete prop1stLayer ;
00079   prop1stLayer = new PropagatorWithMaterial(oppositeToMomentum,mass,theMagField) ;
00080   if (prop2ndLayer) delete prop2ndLayer ;
00081   prop2ndLayer = new PropagatorWithMaterial(alongMomentum,mass,theMagField) ;
00082  }
00083 
00084 vector<CLHEP::Hep3Vector> PixelHitMatcher::predicted1Hits()
00085  { return pred1Meas ; }
00086 
00087 vector<CLHEP::Hep3Vector> PixelHitMatcher::predicted2Hits()
00088  { return pred2Meas ; }
00089 
00090 float PixelHitMatcher::getVertex()
00091  { return vertex_ ; }
00092 
00093 //CLHEP::Hep3Vector point_to_vector( const GlobalPoint & p )
00094 // { return CLHEP::Hep3Vector(p.x(),p.y(),p.z()) ; }
00095 
00096 std::vector<SeedWithInfo>
00097 PixelHitMatcher::compatibleSeeds
00098  ( TrajectorySeedCollection * seeds, const GlobalPoint & xmeas,
00099    const GlobalPoint & vprim, float energy, float fcharge )
00100  {
00101   int charge = int(fcharge) ;
00102 
00103   FreeTrajectoryState fts = myFTS(theMagField,xmeas, vprim, energy, charge);
00104   PerpendicularBoundPlaneBuilder bpb;
00105   TrajectoryStateOnSurface tsos(fts, *bpb(fts.position(), fts.momentum()));
00106 
00107   std::vector<SeedWithInfo> result ;
00108   mapTsos_.clear() ;
00109   mapTsos2_.clear() ;
00110   mapTsos_.reserve(seeds->size()) ;
00111   mapTsos2_.reserve(seeds->size()) ;
00112 
00113   for (unsigned int i=0;i<seeds->size();++i)
00114    {
00115     if ((*seeds)[i].nHits()>9)
00116      {
00117       edm::LogWarning("GsfElectronAlgo|UnexpectedSeed") <<"We cannot deal with seeds having more than 9 hits." ;
00118       continue ;
00119      }
00120     TrajectorySeed::range rhits=(*seeds)[i].recHits() ;
00121 
00122     // build all possible pairs
00123     unsigned char rank1, rank2, hitsMask ;
00124     TrajectorySeed::const_iterator it1, it2 ;
00125     for ( rank1=0, it1=rhits.first ; it1!=rhits.second ; rank1++, it1++ )
00126      {
00127       for ( rank2=rank1+1, it2=it1+1 ; it2!=rhits.second ; rank2++, it2++ )
00128        {
00129         //TrajectorySeed::range r(it1,it2) ;
00130 
00131         // first Hit
00132         TrajectorySeed::const_iterator it=it1 ;
00133         if (!(*it).isValid()) continue;
00134         DetId id=(*it).geographicalId();
00135         const GeomDet *geomdet=theTrackerGeometry->idToDet((*it).geographicalId());
00136         LocalPoint lp=(*it).localPosition() ;
00137         GlobalPoint hitPos=geomdet->surface().toGlobal(lp) ;
00138 
00139         TrajectoryStateOnSurface tsos1;
00140         bool found = false;
00141         std::vector<std::pair<const GeomDet *, TrajectoryStateOnSurface> >::iterator itTsos ;
00142         for (itTsos=mapTsos_.begin();itTsos!=mapTsos_.end();++itTsos)
00143          {
00144           if ((*itTsos).first==geomdet)
00145            { found=true ; break ; }
00146          }
00147         if (!found)
00148          {
00149           tsos1 = prop1stLayer->propagate(tsos,geomdet->surface()) ;
00150           mapTsos_.push_back(std::pair<const GeomDet *, TrajectoryStateOnSurface>(geomdet,tsos1));
00151          }
00152         else
00153          { tsos1=(*itTsos).second ; }
00154 
00155         if (tsos1.isValid())
00156          {
00157           std::pair<bool,double> est;
00158           if (id.subdetId()%2==1) est=meas1stBLayer.estimate(vprim, tsos1,hitPos);
00159           else est=meas1stFLayer.estimate(vprim, tsos1,hitPos);
00160           if (!est.first) continue ;
00161 
00162           if (std::abs(normalized_phi(hitPos.phi()-xmeas.phi()))>2.5) continue ;
00163           EleRelPointPair pp1(hitPos,tsos1.globalParameters().position(),vprim) ;
00164           int subDet1 = id.subdetId() ;
00165           float dRz1 = (subDet1%2==1)?pp1.dZ():pp1.dPerp() ;
00166           float dPhi1 = pp1.dPhi() ;
00167 
00168           // now second Hit
00169           //CC@@
00170           //it++;
00171           it=it2 ;
00172           if (!(*it).isValid()) continue ;
00173 
00174           DetId id2=(*it).geographicalId();
00175           const GeomDet *geomdet2=theTrackerGeometry->idToDet((*it).geographicalId());
00176           TrajectoryStateOnSurface tsos2;
00177 
00178           double zVertex;
00179           if (!useRecoVertex_) // we don't know the z vertex position, get it from linear extrapolation
00180            {
00181             // compute the z vertex from the cluster point and the found pixel hit
00182             double pxHit1z = hitPos.z();
00183             double pxHit1x = hitPos.x();
00184             double pxHit1y = hitPos.y();
00185             double r1diff = (pxHit1x-vprim.x())*(pxHit1x-vprim.x()) + (pxHit1y-vprim.y())*(pxHit1y-vprim.y()) ;
00186             r1diff=sqrt(r1diff) ;
00187             double r2diff = (xmeas.x()-pxHit1x)*(xmeas.x()-pxHit1x) + (xmeas.y()-pxHit1y)*(xmeas.y()-pxHit1y) ;
00188             r2diff=sqrt(r2diff);
00189             zVertex = pxHit1z - r1diff*(xmeas.z()-pxHit1z)/r2diff;
00190            }
00191           else // here use rather the reco vertex z position
00192            { zVertex = vprim.z() ; }
00193 
00194           GlobalPoint vertex(vprim.x(),vprim.y(),zVertex) ;
00195           FreeTrajectoryState fts2 = myFTS(theMagField,hitPos,vertex,energy, charge) ;
00196 
00197           found = false;
00198           std::vector<std::pair< std::pair<const GeomDet *,GlobalPoint>, TrajectoryStateOnSurface> >::iterator itTsos2 ;
00199           for (itTsos2=mapTsos2_.begin();itTsos2!=mapTsos2_.end();++itTsos2)
00200            {
00201             if (((*itTsos2).first).first==geomdet2 &&
00202                 (((*itTsos2).first).second).x()==hitPos.x() &&
00203                 (((*itTsos2).first).second).y()==hitPos.y() &&
00204                 (((*itTsos2).first).second).z()==hitPos.z()  )
00205              {
00206               found=true;
00207               break;
00208              }
00209            }
00210           if (!found)
00211            {
00212             tsos2 = prop2ndLayer->propagate(fts2,geomdet2->surface()) ;
00213             std::pair<const GeomDet *,GlobalPoint> pair(geomdet2,hitPos);
00214             mapTsos2_.push_back(std::pair<std::pair<const GeomDet *,GlobalPoint>, TrajectoryStateOnSurface> (pair,tsos2));
00215            }
00216           else
00217            { tsos2=(*itTsos2).second ; }
00218 
00219           if (tsos2.isValid())
00220            {
00221             LocalPoint lp2=(*it).localPosition() ;
00222             GlobalPoint hitPos2=geomdet2->surface().toGlobal(lp2) ;
00223             std::pair<bool,double> est2 ;
00224             if (id2.subdetId()%2==1) est2=meas2ndBLayer.estimate(vertex, tsos2,hitPos2) ;
00225             else est2=meas2ndFLayer.estimate(vertex, tsos2,hitPos2) ;
00226             if (est2.first)
00227              {
00228               EleRelPointPair pp2(hitPos2,tsos2.globalParameters().position(),vertex) ;
00229               int subDet2 = id2.subdetId() ;
00230               float dRz2 = (subDet2%2==1)?pp2.dZ():pp2.dPerp() ;
00231               float dPhi2 = pp2.dPhi() ;
00232               hitsMask = (1<<rank1)|(1<<rank2) ;
00233               result.push_back(SeedWithInfo((*seeds)[i],hitsMask,subDet2,dRz2,dPhi2,subDet1,dRz1,dPhi1)) ;
00234              }
00235            }
00236          } // end tsos1 is valid
00237        } // end loop on second seed hit
00238      } // end loop on first seed hit
00239    } // end loop on seeds
00240 
00241   mapTsos_.clear() ;
00242   mapTsos2_.clear() ;
00243 
00244   return result ;
00245  }
00246 
00247 //========================= OBSOLETE ? =========================
00248 
00249 vector< pair< RecHitWithDist, PixelHitMatcher::ConstRecHitPointer > >
00250 PixelHitMatcher::compatibleHits
00251  ( const GlobalPoint & xmeas,
00252    const GlobalPoint & vprim,
00253    float energy, float fcharge,
00254    const TrackerTopology *tTopo)
00255  {
00256   float SCl_phi = xmeas.phi();
00257 
00258   int charge = int(fcharge);
00259   // return all compatible RecHit pairs (vector< TSiPixelRecHit>)
00260   vector<pair<RecHitWithDist, ConstRecHitPointer> > result;
00261   LogDebug("") << "[PixelHitMatcher::compatibleHits] entering .. ";
00262 
00263   vector<TrajectoryMeasurement> validMeasurements;
00264   vector<TrajectoryMeasurement> invalidMeasurements;
00265 
00266   typedef vector<TrajectoryMeasurement>::const_iterator aMeas;
00267 
00268   pred1Meas.clear();
00269   pred2Meas.clear();
00270 
00271   typedef vector<BarrelDetLayer*>::const_iterator BarrelLayerIterator;
00272   BarrelLayerIterator firstLayer = startLayers.firstBLayer();
00273 
00274   FreeTrajectoryState fts =myFTS(theMagField,xmeas, vprim,
00275                                  energy, charge);
00276 
00277   PerpendicularBoundPlaneBuilder bpb;
00278   TrajectoryStateOnSurface tsos(fts, *bpb(fts.position(), fts.momentum()));
00279 
00280   if (tsos.isValid()) {
00281     vector<TrajectoryMeasurement> pixelMeasurements =
00282       theLayerMeasurements->measurements(**firstLayer,tsos,
00283                                          *prop1stLayer, meas1stBLayer);
00284 
00285     LogDebug("") <<"[PixelHitMatcher::compatibleHits] nbr of hits compatible with extrapolation to first layer: " << pixelMeasurements.size();
00286     for (aMeas m=pixelMeasurements.begin(); m!=pixelMeasurements.end(); m++){
00287      if (m->recHit()->isValid()) {
00288        float localDphi = normalized_phi(SCl_phi-m->forwardPredictedState().globalPosition().phi()) ;
00289        if(std::abs(localDphi)>2.5)continue;
00290         CLHEP::Hep3Vector prediction(m->forwardPredictedState().globalPosition().x(),
00291                               m->forwardPredictedState().globalPosition().y(),
00292                               m->forwardPredictedState().globalPosition().z());
00293         LogDebug("") << "[PixelHitMatcher::compatibleHits] compatible hit position " << m->recHit()->globalPosition();
00294         LogDebug("") << "[PixelHitMatcher::compatibleHits] predicted position " << m->forwardPredictedState().globalPosition();
00295         pred1Meas.push_back( prediction);
00296 
00297         validMeasurements.push_back(*m);
00298 
00299         LogDebug("") <<"[PixelHitMatcher::compatibleHits] Found a rechit in layer ";
00300         const BarrelDetLayer *bdetl = dynamic_cast<const BarrelDetLayer *>(*firstLayer);
00301         if (bdetl) {
00302           LogDebug("") <<" with radius "<<bdetl->specificSurface().radius();
00303         }
00304         else  LogDebug("") <<"Could not downcast!!";
00305      }
00306     }
00307 
00308 
00309     // check if there are compatible 1st hits in the second layer
00310     firstLayer++;
00311 
00312     vector<TrajectoryMeasurement> pixel2Measurements =
00313       theLayerMeasurements->measurements(**firstLayer,tsos,
00314                                          *prop1stLayer, meas1stBLayer);
00315 
00316     for (aMeas m=pixel2Measurements.begin(); m!=pixel2Measurements.end(); m++){
00317       if (m->recHit()->isValid()) {
00318         float localDphi = normalized_phi(SCl_phi-m->forwardPredictedState().globalPosition().phi()) ;
00319         if(std::abs(localDphi)>2.5)continue;
00320         CLHEP::Hep3Vector prediction(m->forwardPredictedState().globalPosition().x(),
00321                               m->forwardPredictedState().globalPosition().y(),
00322                               m->forwardPredictedState().globalPosition().z());
00323         pred1Meas.push_back( prediction);
00324         LogDebug("")  << "[PixelHitMatcher::compatibleHits] compatible hit position " << m->recHit()->globalPosition() << endl;
00325         LogDebug("") << "[PixelHitMatcher::compatibleHits] predicted position " << m->forwardPredictedState().globalPosition() << endl;
00326 
00327         validMeasurements.push_back(*m);
00328         LogDebug("") <<"[PixelHitMatcher::compatibleHits] Found a rechit in layer ";
00329         const BarrelDetLayer *bdetl = dynamic_cast<const BarrelDetLayer *>(*firstLayer);
00330         if (bdetl) {
00331           LogDebug("") <<" with radius "<<bdetl->specificSurface().radius();
00332         }
00333         else  LogDebug("") <<"Could not downcast!!";
00334       }
00335 
00336     }
00337   }
00338 
00339 
00340   // check if there are compatible 1st hits the forward disks
00341   typedef vector<ForwardDetLayer*>::const_iterator ForwardLayerIterator;
00342   ForwardLayerIterator flayer;
00343 
00344   TrajectoryStateOnSurface tsosfwd(fts, *bpb(fts.position(), fts.momentum()));
00345   if (tsosfwd.isValid()) {
00346 
00347     for (int i=0; i<2; i++) {
00348       i == 0 ? flayer = startLayers.pos1stFLayer() : flayer = startLayers.neg1stFLayer();
00349 
00350       if (i==0 && xmeas.z() < -100. ) continue;
00351       if (i==1 && xmeas.z() > 100. ) continue;
00352 
00353       vector<TrajectoryMeasurement> pixelMeasurements =
00354         theLayerMeasurements->measurements(**flayer, tsosfwd,
00355                                            *prop1stLayer, meas1stFLayer);
00356 
00357       for (aMeas m=pixelMeasurements.begin(); m!=pixelMeasurements.end(); m++){
00358         if (m->recHit()->isValid()) {
00359           float localDphi = normalized_phi(SCl_phi-m->forwardPredictedState().globalPosition().phi());
00360           if(std::abs(localDphi)>2.5)continue;
00361           CLHEP::Hep3Vector prediction(m->forwardPredictedState().globalPosition().x(),
00362                                 m->forwardPredictedState().globalPosition().y(),
00363                                 m->forwardPredictedState().globalPosition().z());
00364           pred1Meas.push_back( prediction);
00365 
00366           validMeasurements.push_back(*m);
00367         }
00368       }
00369 
00370       //check if there are compatible 1st hits the outer forward disks
00371       if (searchInTIDTEC_) {
00372         flayer++;
00373 
00374         vector<TrajectoryMeasurement> pixel2Measurements =
00375           theLayerMeasurements->measurements(**flayer, tsosfwd,
00376                                              *prop1stLayer, meas1stFLayer);
00377 
00378         for (aMeas m=pixel2Measurements.begin(); m!=pixel2Measurements.end(); m++){
00379           if (m->recHit()->isValid()) {
00380             float localDphi = normalized_phi(SCl_phi-m->forwardPredictedState().globalPosition().phi()) ;
00381             if(std::abs(localDphi)>2.5)continue;
00382             CLHEP::Hep3Vector prediction(m->forwardPredictedState().globalPosition().x(),
00383                                   m->forwardPredictedState().globalPosition().y(),
00384                                   m->forwardPredictedState().globalPosition().z());
00385             pred1Meas.push_back( prediction);
00386 
00387             validMeasurements.push_back(*m);
00388           }
00389           //    else{std::cout<<" hit non valid "<<std::endl; }
00390         }  //end 1st hit in outer f disk
00391       }
00392     }
00393   }
00394 
00395   // now we have the vector of all valid measurements of the first point
00396   for (unsigned i=0; i<validMeasurements.size(); i++){
00397 
00398     const DetLayer * newLayer = theGeometricSearchTracker->detLayer(validMeasurements[i].recHit()->det()->geographicalId());
00399 
00400     double zVertex ;
00401     if (!useRecoVertex_)
00402      {
00403       // we don't know the z vertex position, get it from linear extrapolation
00404       // compute the z vertex from the cluster point and the found pixel hit
00405       double pxHit1z = validMeasurements[i].recHit()->det()->surface().toGlobal(
00406         validMeasurements[i].recHit()->localPosition()).z();
00407       double pxHit1x = validMeasurements[i].recHit()->det()->surface().toGlobal(
00408         validMeasurements[i].recHit()->localPosition()).x();
00409       double pxHit1y = validMeasurements[i].recHit()->det()->surface().toGlobal(
00410         validMeasurements[i].recHit()->localPosition()).y();
00411       double r1diff = (pxHit1x-vprim.x())*(pxHit1x-vprim.x()) + (pxHit1y-vprim.y())*(pxHit1y-vprim.y());
00412       r1diff=sqrt(r1diff);
00413       double r2diff = (xmeas.x()-pxHit1x)*(xmeas.x()-pxHit1x) + (xmeas.y()-pxHit1y)*(xmeas.y()-pxHit1y);
00414       r2diff=sqrt(r2diff);
00415       zVertex = pxHit1z - r1diff*(xmeas.z()-pxHit1z)/r2diff;
00416      }
00417     else
00418      {
00419       // here we use the reco vertex z position
00420       zVertex = vprim.z();
00421      }
00422 
00423     if (i==0)
00424      { vertex_ = zVertex; }
00425 
00426     GlobalPoint vertexPred(vprim.x(),vprim.y(),zVertex) ;
00427     GlobalPoint hitPos( validMeasurements[i].recHit()->det()->surface().toGlobal( validMeasurements[i].recHit()->localPosition() ) ) ;
00428 
00429     FreeTrajectoryState secondFTS=myFTS(theMagField,hitPos,vertexPred,energy, charge);
00430 
00431     PixelMatchNextLayers secondHit(theLayerMeasurements, newLayer, secondFTS,
00432                                    prop2ndLayer, &meas2ndBLayer,&meas2ndFLayer,
00433                                    tTopo,searchInTIDTEC_);
00434     vector<CLHEP::Hep3Vector> predictions = secondHit.predictionInNextLayers();
00435 
00436     for (unsigned it = 0; it < predictions.size(); it++) pred2Meas.push_back(predictions[it]);
00437 
00438     // we may get more than one valid second measurements here even for single electrons:
00439     // two hits from the same layer/disk (detector overlap) or from the loop over the
00440     // next layers in EPMatchLoopNextLayers. Take only the 1st hit.
00441 
00442     if(!secondHit.measurementsInNextLayers().empty()){
00443       for(unsigned int shit=0; shit<secondHit.measurementsInNextLayers().size(); shit++)
00444         {
00445           float dphi = normalized_phi(pred1Meas[i].phi()-validMeasurements[i].recHit()->globalPosition().phi()) ;
00446           if (std::abs(dphi)<2.5)
00447             {
00448               ConstRecHitPointer pxrh = validMeasurements[i].recHit();
00449               RecHitWithDist rh(pxrh,dphi);
00450 
00451               //  pxrh = secondHit.measurementsInNextLayers()[0].recHit();
00452               pxrh = secondHit.measurementsInNextLayers()[shit].recHit();
00453 
00454               pair<RecHitWithDist,ConstRecHitPointer> compatiblePair = pair<RecHitWithDist,ConstRecHitPointer>(rh,pxrh) ;
00455               result.push_back(compatiblePair);
00456               break;
00457             }
00458         }
00459     }
00460   }
00461   return result;
00462 }
00463 
00464