CMS 3D CMS Logo

TrackAssociatorByChi2 Class Reference

Class that performs the association of reco::Tracks and TrackingParticles evaluating the chi2 of reco tracks parameters and sim tracks parameters. More...

#include <SimTracker/TrackAssociation/interface/TrackAssociatorByChi2.h>

Inheritance diagram for TrackAssociatorByChi2:

TrackAssociatorBase

List of all members.

Public Types

typedef std::map< double,
SimTrack
Chi2SimMap
typedef std::pair< reco::Track,
Chi2SimMap
RecoToSimPair
typedef std::vector
< RecoToSimPair
RecoToSimPairAssociation

Public Member Functions

reco::RecoToSimCollection associateRecoToSim (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0) const
 compare reco to sim the handle of reco::Track and TrackingParticle collections
reco::RecoToSimCollection associateRecoToSim (edm::RefToBaseVector< reco::Track > &, edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0) const
 Association Reco To Sim with Collections.
double associateRecoToSim (reco::TrackCollection::const_iterator, TrackingParticleCollection::const_iterator, const reco::BeamSpot &) const
 compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2
reco::SimToRecoCollection associateSimToReco (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0) const
 compare reco to sim the handle of reco::Track and TrackingParticle collections
reco::SimToRecoCollection associateSimToReco (edm::RefToBaseVector< reco::Track > &, edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0) const
 Association Sim To Reco with Collections.
RecoToSimPairAssociation compareTracksParam (const reco::TrackCollection &, const edm::SimTrackContainer &, const edm::SimVertexContainer &, const reco::BeamSpot &) const
 compare collections reco to sim
double compareTracksParam (reco::TrackCollection::const_iterator, edm::SimTrackContainer::const_iterator, const math::XYZTLorentzVectorD, GlobalVector, reco::TrackBase::CovarianceMatrix, const reco::BeamSpot &) const
 compare reco::TrackCollection and edm::SimTrackContainer iterators: returns the chi2
std::pair< bool,
reco::TrackBase::ParameterVector
parametersAtClosestApproach (Basic3DVector< double >, Basic3DVector< double >, float, const reco::BeamSpot &) const
 propagate the track parameters of TrackinParticle from production vertex to the point of closest approach to the beam line.
 TrackAssociatorByChi2 (const edm::ESHandle< MagneticField > mF, double chi2Cut, bool onlyDiag, edm::InputTag beamspotSrc)
 Constructor with magnetic field, double, bool and InputTag.
 TrackAssociatorByChi2 (const edm::ESHandle< MagneticField > mF, edm::ParameterSet conf)
 Constructor with PSet.
 ~TrackAssociatorByChi2 ()
 Destructor.

Private Attributes

edm::InputTag bsSrc
double chi2cut
bool onlyDiagonal
edm::ESHandle< MagneticFieldtheMF


Detailed Description

Class that performs the association of reco::Tracks and TrackingParticles evaluating the chi2 of reco tracks parameters and sim tracks parameters.

The cut can be tuned from the config file: see data/TrackAssociatorByChi2.cfi. Note that the Association Map is filled with -ch2 and not chi2 because it is ordered using std::greater: the track with the lowest association chi2 will be the first in the output map.It is possible to use only diagonal terms (associator by pulls) seeting onlyDiagonal = true in the PSet

Date
2008/02/22 13:56:24
Revision
1.24
Author:
cerati, magni

Definition at line 30 of file TrackAssociatorByChi2.h.


Member Typedef Documentation

typedef std::map<double, SimTrack> TrackAssociatorByChi2::Chi2SimMap

Definition at line 33 of file TrackAssociatorByChi2.h.

typedef std::pair< reco::Track, Chi2SimMap> TrackAssociatorByChi2::RecoToSimPair

Definition at line 34 of file TrackAssociatorByChi2.h.

typedef std::vector< RecoToSimPair > TrackAssociatorByChi2::RecoToSimPairAssociation

Definition at line 35 of file TrackAssociatorByChi2.h.


Constructor & Destructor Documentation

TrackAssociatorByChi2::TrackAssociatorByChi2 ( const edm::ESHandle< MagneticField mF,
edm::ParameterSet  conf 
) [inline]

Constructor with PSet.

Definition at line 38 of file TrackAssociatorByChi2.h.

References onlyDiagonal, and theMF.

00038                                                                                   :
00039     chi2cut(conf.getParameter<double>("chi2cut")),
00040     onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
00041     bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
00042     theMF=mF;  
00043     if (onlyDiagonal)
00044       edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " <<  "\n";
00045     else 
00046       edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " <<  "\n";
00047   }

TrackAssociatorByChi2::TrackAssociatorByChi2 ( const edm::ESHandle< MagneticField mF,
double  chi2Cut,
bool  onlyDiag,
edm::InputTag  beamspotSrc 
) [inline]

Constructor with magnetic field, double, bool and InputTag.

Definition at line 50 of file TrackAssociatorByChi2.h.

References bsSrc, chi2cut, onlyDiagonal, and theMF.

00050                                                                                                                     {
00051     chi2cut=chi2Cut;
00052     onlyDiagonal=onlyDiag;
00053     theMF=mF;  
00054     bsSrc = beamspotSrc;
00055   }

TrackAssociatorByChi2::~TrackAssociatorByChi2 (  )  [inline]

Destructor.

Definition at line 58 of file TrackAssociatorByChi2.h.

00058 {}


Member Function Documentation

reco::RecoToSimCollection TrackAssociatorByChi2::associateRecoToSim ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0 
) const [inline, virtual]

compare reco to sim the handle of reco::Track and TrackingParticle collections

Reimplemented from TrackAssociatorBase.

Definition at line 94 of file TrackAssociatorByChi2.h.

References TrackAssociatorBase::associateRecoToSim(), and event().

00096                                                                                  {
00097     return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event);
00098   }

RecoToSimCollection TrackAssociatorByChi2::associateRecoToSim ( edm::RefToBaseVector< reco::Track > &  tC,
edm::RefVector< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0 
) const [virtual]

Association Reco To Sim with Collections.

Implements TrackAssociatorBase.

Definition at line 161 of file TrackAssociatorByChi2.cc.

References edm::RefToBaseVector< T >::begin(), bsSrc, chi2cut, edm::RefToBaseVector< T >::end(), edm::Event::getByLabel(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, LogTrace, onlyDiagonal, parametersAtClosestApproach(), reco::BeamSpot::position(), edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), and edm::RefVector< C, T, F >::size().

00163                                                                                         {
00164   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00165   e->getByLabel(bsSrc,recoBeamSpotHandle);
00166   reco::BeamSpot bs = *recoBeamSpotHandle;      
00167 
00168   RecoToSimCollection  outputCollection;
00169   double chi2;
00170 
00171   TrackingParticleCollection tPC;
00172   if (tPCH.size()!=0)  tPC = *tPCH.product();
00173 
00174   int tindex=0;
00175   for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
00176 
00177     LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
00178                                 << "rec::Track #"<<tindex<<" with pt=" << (*rt)->pt() <<  "\n"
00179                                 << "===========================================" << "\n";
00180  
00181     TrackBase::ParameterVector rParameters = (*rt)->parameters();
00182 
00183     TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
00184     if (onlyDiagonal){
00185       for (unsigned int i=0;i<5;i++){
00186         for (unsigned int j=0;j<5;j++){
00187           if (i!=j) recoTrackCovMatrix(i,j)=0;
00188         }
00189       }
00190     } 
00191 
00192     recoTrackCovMatrix.Invert();
00193 
00194     int tpindex =0;
00195     for (TrackingParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
00196         
00197       //skip tps with a very small pt
00198       //if (sqrt(tp->momentum().perp2())<0.5) continue;
00199       Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
00200       Basic3DVector<double> vert=(Basic3DVector<double>) tp->vertex();
00201 
00202       std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, tp->charge(), bs);
00203       if (params.first){
00204         TrackBase::ParameterVector sParameters=params.second;
00205         
00206         TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00207         
00208         chi2 = ROOT::Math::Similarity(diffParameters, recoTrackCovMatrix);
00209         chi2 /= 5;
00210         
00211         LogTrace("TrackAssociator") << "====TP index=" << tpindex << "  RECO TRACK index="<<tindex<<" WITH PT=" << (*rt)->pt() << "====\n" 
00212                                     << "qoverp sim: " << sParameters[0] << "\n" 
00213                                     << "lambda sim: " << sParameters[1] << "\n" 
00214                                     << "phi    sim: " << sParameters[2] << "\n" 
00215                                     << "dxy    sim: " << sParameters[3] << "\n" 
00216                                     << "dsz    sim: " << sParameters[4] << "\n" 
00217                                     << ": " /*<< */ << "\n" 
00218                                     << "qoverp rec: " << (*rt)->qoverp() << " err: " << (*rt)->error(0) << "\n"
00219                                     << "lambda rec: " << (*rt)->lambda() << " err: " << (*rt)->error(1) << "\n"
00220                                     << "phi    rec: " << (*rt)->phi() << " err: " << (*rt)->error(2) << "\n"
00221                                     << "dxy    rec: " << (*rt)->dxy(bs.position()) << " err: " << (*rt)->error(3) << "\n"
00222                                     << "dsz    rec: " << (*rt)->dsz(bs.position()) << " err: " << (*rt)->error(4) << "\n"
00223                                     << ": " /*<< */ << "\n" 
00224                                     << "chi2: " << chi2 << "\n";
00225         
00226         if (chi2<chi2cut) {
00227           outputCollection.insert(tC[tindex], 
00228                                   std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH, tpindex),
00229                                                  -chi2));//-chi2 because the Association Map is ordered using std::greater
00230         }
00231       }
00232     }
00233   }
00234   outputCollection.post_insert();
00235   return outputCollection;
00236 }

double TrackAssociatorByChi2::associateRecoToSim ( reco::TrackCollection::const_iterator  ,
TrackingParticleCollection::const_iterator  ,
const reco::BeamSpot  
) const

compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2

Referenced by KVFTest::analyze().

reco::SimToRecoCollection TrackAssociatorByChi2::associateSimToReco ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0 
) const [inline, virtual]

compare reco to sim the handle of reco::Track and TrackingParticle collections

Reimplemented from TrackAssociatorBase.

Definition at line 101 of file TrackAssociatorByChi2.h.

References TrackAssociatorBase::associateSimToReco(), and event().

00103                                                                                  {
00104     return TrackAssociatorBase::associateSimToReco(tCH,tPCH,event);
00105   }  

SimToRecoCollection TrackAssociatorByChi2::associateSimToReco ( edm::RefToBaseVector< reco::Track > &  tC,
edm::RefVector< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0 
) const [virtual]

Association Sim To Reco with Collections.

Implements TrackAssociatorBase.

Definition at line 240 of file TrackAssociatorByChi2.cc.

References edm::RefToBaseVector< T >::begin(), bsSrc, chi2cut, edm::RefToBaseVector< T >::end(), edm::Event::getByLabel(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, LogTrace, onlyDiagonal, parametersAtClosestApproach(), reco::BeamSpot::position(), edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), edm::RefVector< C, T, F >::size(), and funct::sqrt().

00242                                                                                          {
00243   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00244   e->getByLabel(bsSrc,recoBeamSpotHandle);
00245   reco::BeamSpot bs = *recoBeamSpotHandle;      
00246 
00247   SimToRecoCollection  outputCollection;
00248   double chi2;
00249 
00250   TrackingParticleCollection tPC;
00251   if (tPCH.size()!=0)  tPC = *tPCH.product();
00252 
00253   int tpindex =0;
00254   for (TrackingParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
00255     
00256     //skip tps with a very small pt
00257     //if (sqrt(tp->momentum().perp2())<0.5) continue;
00258     
00259     LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
00260                                 << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt(tp->momentum().perp2()) << "\n"
00261                                 << "===========================================" << "\n";
00262     
00263     Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
00264     Basic3DVector<double> vert(tp->vertex().x(),tp->vertex().y(),tp->vertex().z());
00265     
00266     std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, tp->charge(), bs);
00267     if (params.first){
00268       TrackBase::ParameterVector sParameters=params.second;
00269       
00270       int tindex=0;
00271       for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
00272         
00273         TrackBase::ParameterVector rParameters = (*rt)->parameters();
00274 
00275         TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
00276         if (onlyDiagonal) {
00277           for (unsigned int i=0;i<5;i++){
00278             for (unsigned int j=0;j<5;j++){
00279               if (i!=j) recoTrackCovMatrix(i,j)=0;
00280             }
00281           }
00282         }
00283         
00284         recoTrackCovMatrix.Invert();
00285         
00286         TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00287         
00288         chi2 = ROOT::Math::Similarity(recoTrackCovMatrix, diffParameters);
00289         chi2 /= 5;
00290         
00291         LogTrace("TrackAssociator") << "====TP index=" << tpindex << "  RECO TRACK index="<<tindex<<" WITH PT=" << (*rt)->pt() << "====\n" 
00292                                     << "qoverp sim: " << sParameters[0] << "\n" 
00293                                     << "lambda sim: " << sParameters[1] << "\n" 
00294                                     << "phi    sim: " << sParameters[2] << "\n" 
00295                                     << "dxy    sim: " << sParameters[3] << "\n" 
00296                                     << "dsz    sim: " << sParameters[4] << "\n" 
00297                                     << ": " /*<< */ << "\n" 
00298                                     << "qoverp rec: " << (*rt)->qoverp() << " err: " << (*rt)->error(0) << "\n"
00299                                     << "lambda rec: " << (*rt)->lambda() << " err: " << (*rt)->error(1) << "\n"
00300                                     << "phi    rec: " << (*rt)->phi() << " err: " << (*rt)->error(2) << "\n"
00301                                     << "dxy    rec: " << (*rt)->dxy(bs.position()) << " err: " << (*rt)->error(3) << "\n"
00302                                     << "dsz    rec: " << (*rt)->dsz(bs.position()) << " err: " << (*rt)->error(4) << "\n"
00303                                     << ": " /*<< */ << "\n" 
00304                                     << "chi2: " << chi2 << "\n";
00305         
00306         if (chi2<chi2cut) {
00307           outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH, tpindex),
00308                                   std::make_pair(tC[tindex],
00309                                                  -chi2));//-chi2 because the Association Map is ordered using std::greater
00310         }
00311       }
00312     }
00313   }
00314   outputCollection.post_insert();
00315   return outputCollection;
00316 }

TrackAssociatorByChi2::RecoToSimPairAssociation TrackAssociatorByChi2::compareTracksParam ( const reco::TrackCollection rtColl,
const edm::SimTrackContainer stColl,
const edm::SimVertexContainer svColl,
const reco::BeamSpot bs 
) const

compare collections reco to sim

Definition at line 39 of file TrackAssociatorByChi2.cc.

References chi2cut, i, j, onlyDiagonal, parametersAtClosestApproach(), and track.

00042                                                                        {
00043   
00044   RecoToSimPairAssociation outputVec;
00045 
00046   for (TrackCollection::const_iterator track=rtColl.begin(); track!=rtColl.end(); track++){
00047      Chi2SimMap outMap;
00048 
00049     TrackBase::ParameterVector rParameters = track->parameters();
00050 
00051     TrackBase::CovarianceMatrix recoTrackCovMatrix = track->covariance();
00052     if (onlyDiagonal){
00053       for (unsigned int i=0;i<5;i++){
00054         for (unsigned int j=0;j<5;j++){
00055           if (i!=j) recoTrackCovMatrix(i,j)=0;
00056         }
00057       }
00058     }
00059     recoTrackCovMatrix.Invert();
00060 
00061     for (SimTrackContainer::const_iterator st=stColl.begin(); st!=stColl.end(); st++){
00062 
00063       Basic3DVector<double> momAtVtx(st->momentum().x(),st->momentum().y(),st->momentum().z());
00064       Basic3DVector<double> vert = (Basic3DVector<double>)  svColl[st->vertIndex()].position();
00065 
00066       std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, st->charge(), bs);
00067       if (params.first){
00068         TrackBase::ParameterVector sParameters = params.second;
00069       
00070         TrackBase::ParameterVector diffParameters = rParameters - sParameters;
00071         double chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
00072         chi2/=5;
00073         if (chi2<chi2cut) outMap[chi2]=*st;
00074       }
00075     }
00076     outputVec.push_back(RecoToSimPair(*track,outMap));
00077   }
00078   return outputVec;
00079 }

double TrackAssociatorByChi2::compareTracksParam ( reco::TrackCollection::const_iterator  ,
edm::SimTrackContainer::const_iterator  ,
const math::XYZTLorentzVectorD  ,
GlobalVector  ,
reco::TrackBase::CovarianceMatrix  ,
const reco::BeamSpot  
) const

compare reco::TrackCollection and edm::SimTrackContainer iterators: returns the chi2

pair< bool, TrackBase::ParameterVector > TrackAssociatorByChi2::parametersAtClosestApproach ( Basic3DVector< double >  vertex,
Basic3DVector< double >  momAtVtx,
float  charge,
const reco::BeamSpot bs 
) const

propagate the track parameters of TrackinParticle from production vertex to the point of closest approach to the beam line.

Definition at line 133 of file TrackAssociatorByChi2.cc.

References FreeTrajectoryState::charge(), funct::cos(), Geom::halfPi(), PV3DBase< T, PVType, FrameType >::mag(), FreeTrajectoryState::momentum(), p, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), FreeTrajectoryState::position(), edm::ESHandle< T >::product(), funct::sin(), theMF, PV3DBase< T, PVType, FrameType >::theta(), TrajectoryStateClosestToBeamLine::trackStateAtPCA(), v, PV3DBase< T, PVType, FrameType >::x(), Basic3DVector< T >::x(), PV3DBase< T, PVType, FrameType >::y(), Basic3DVector< T >::y(), PV3DBase< T, PVType, FrameType >::z(), and Basic3DVector< T >::z().

Referenced by associateRecoToSim(), associateSimToReco(), compareTracksParam(), and VertexFitterResult::fill().

00136                                                                             {
00137   
00138   TrackBase::ParameterVector sParameters;
00139   try {
00140     FreeTrajectoryState ftsAtProduction(GlobalPoint(vertex.x(),vertex.y(),vertex.z()),
00141                                         GlobalVector(momAtVtx.x(),momAtVtx.y(),momAtVtx.z()),
00142                                         TrackCharge(charge),
00143                                         theMF.product());
00144     TrajectoryStateClosestToBeamLineBuilder tscblBuilder;
00145     TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,bs);//as in TrackProducerAlgorithm
00146     
00147     GlobalPoint v = tsAtClosestApproach.trackStateAtPCA().position();
00148     GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
00149     sParameters[0] = tsAtClosestApproach.trackStateAtPCA().charge()/p.mag();
00150     sParameters[1] = Geom::halfPi() - p.theta();
00151     sParameters[2] = p.phi();
00152     sParameters[3] = (-v.x()*sin(p.phi())+v.y()*cos(p.phi()));
00153     sParameters[4] = v.z()*p.perp()/p.mag() - (v.x()*p.x()+v.y()*p.y())/p.perp() * p.z()/p.mag();
00154     
00155     return make_pair<bool,TrackBase::ParameterVector>(true,sParameters);
00156   } catch ( ... ) {
00157     return make_pair<bool,TrackBase::ParameterVector>(false,sParameters);
00158   }
00159 }


Member Data Documentation

edm::InputTag TrackAssociatorByChi2::bsSrc [private]

Definition at line 110 of file TrackAssociatorByChi2.h.

Referenced by associateRecoToSim(), associateSimToReco(), and TrackAssociatorByChi2().

double TrackAssociatorByChi2::chi2cut [private]

Definition at line 108 of file TrackAssociatorByChi2.h.

Referenced by associateRecoToSim(), associateSimToReco(), compareTracksParam(), and TrackAssociatorByChi2().

bool TrackAssociatorByChi2::onlyDiagonal [private]

Definition at line 109 of file TrackAssociatorByChi2.h.

Referenced by associateRecoToSim(), associateSimToReco(), compareTracksParam(), and TrackAssociatorByChi2().

edm::ESHandle<MagneticField> TrackAssociatorByChi2::theMF [private]

Definition at line 107 of file TrackAssociatorByChi2.h.

Referenced by parametersAtClosestApproach(), and TrackAssociatorByChi2().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:33:58 2009 for CMSSW by  doxygen 1.5.4