CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

TrackAssociatorByChi2 Class Reference

#include <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::GenToRecoCollection associateGenToReco (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 Association Sim To Reco with Collections (Gen Particle version)
virtual reco::GenToRecoCollection associateGenToReco (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< reco::GenParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 compare reco to sim the handle of reco::Track and GenParticle collections
reco::RecoToGenCollection associateRecoToGen (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 Association Sim To Reco with Collections (Gen Particle version)
virtual reco::RecoToGenCollection associateRecoToGen (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< reco::GenParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 compare reco to sim the handle of reco::Track and GenParticle collections
reco::RecoToSimCollection associateRecoToSim (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 Association Reco To Sim with Collections.
reco::RecoToSimCollection associateRecoToSim (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 compare reco to sim the handle of reco::Track and TrackingParticle 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 (const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 Association Sim To Reco with Collections.
reco::SimToRecoCollection associateSimToReco (edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
 compare reco to sim the handle of reco::Track and TrackingParticle 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
double getChi2 (reco::TrackBase::ParameterVector &rParameters, reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, Basic3DVector< double > &momAtVtx, Basic3DVector< double > &vert, int &charge, const reco::BeamSpot &) const
 basic method where chi2 is computed
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, edm::ParameterSet conf)
 Constructor with PSet.
 TrackAssociatorByChi2 (const edm::ESHandle< MagneticField > mF, double chi2Cut, bool onlyDiag, edm::InputTag beamspotSrc)
 Constructor with magnetic field, double, bool and InputTag.
 ~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:
2012/12/03 10:49:23
Revision:
1.29
Author:
cerati, magni

Definition at line 42 of file TrackAssociatorByChi2.h.


Member Typedef Documentation

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

Definition at line 45 of file TrackAssociatorByChi2.h.

Definition at line 46 of file TrackAssociatorByChi2.h.

Definition at line 47 of file TrackAssociatorByChi2.h.


Constructor & Destructor Documentation

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

Constructor with PSet.

Definition at line 50 of file TrackAssociatorByChi2.h.

References onlyDiagonal, and theMF.

                                                                                  :
    chi2cut(conf.getParameter<double>("chi2cut")),
    onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
    bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
    theMF=mF;  
    if (onlyDiagonal)
      edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " <<  "\n";
    else 
      edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " <<  "\n";
  }
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 62 of file TrackAssociatorByChi2.h.

References bsSrc, chi2cut, onlyDiagonal, and theMF.

                                                                                                                    {
    chi2cut=chi2Cut;
    onlyDiagonal=onlyDiag;
    theMF=mF;  
    bsSrc = beamspotSrc;
  }
TrackAssociatorByChi2::~TrackAssociatorByChi2 ( ) [inline]

Destructor.

Definition at line 70 of file TrackAssociatorByChi2.h.

{}

Member Function Documentation

GenToRecoCollection TrackAssociatorByChi2::associateGenToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< reco::GenParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const

Association Sim To Reco with Collections (Gen Particle version)

Definition at line 345 of file TrackAssociatorByChi2.cc.

References edm::RefToBaseVector< T >::begin(), edm::RefToBaseVector< T >::end(), edm::Event::getByLabel(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), edm::RefVector< C, T, F >::size(), and mathSSE::sqrt().

Referenced by MultiTrackValidatorGenPs::analyze(), and associateGenToReco().

                                                                                                 {

  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
  e->getByLabel(bsSrc,recoBeamSpotHandle);
  reco::BeamSpot bs = *recoBeamSpotHandle;      

  GenToRecoCollection  outputCollection;

  GenParticleCollection tPC;
  if (tPCH.size()!=0)  tPC = *tPCH.product();

  int tpindex =0;
  for (GenParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
    
    //skip tps with a very small pt
    //if (sqrt(tp->momentum().perp2())<0.5) continue;
    int charge = tp->charge();
    if (charge==0) continue;
    
    LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
                                << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt(tp->momentum().perp2()) << "\n"
                                << "===========================================" << "\n";
    
    Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
    Basic3DVector<double> vert(tp->vertex().x(),tp->vertex().y(),tp->vertex().z());
      
    int tindex=0;
    for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
      
      TrackBase::ParameterVector rParameters = (*rt)->parameters();      
      TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
      if (onlyDiagonal) {
        for (unsigned int i=0;i<5;i++){
          for (unsigned int j=0;j<5;j++){
            if (i!=j) recoTrackCovMatrix(i,j)=0;
          }
        }
      }
      recoTrackCovMatrix.Invert();
      
      double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
      
      if (chi2<chi2cut) {
        outputCollection.insert(edm::Ref<GenParticleCollection>(tPCH, tpindex),
                                std::make_pair(tC[tindex],
                                               -chi2));//-chi2 because the Association Map is ordered using std::greater
      }
    }
  }
  outputCollection.post_insert();
  return outputCollection;
}
virtual reco::GenToRecoCollection TrackAssociatorByChi2::associateGenToReco ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< reco::GenParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const [inline, virtual]

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

Definition at line 159 of file TrackAssociatorByChi2.h.

References associateGenToReco(), event(), edm::HandleBase::id(), j, edm::RefVector< C, T, F >::push_back(), edm::RefToBaseVector< T >::push_back(), and HcalObjRepresent::setup().

                                                                                              {
    edm::RefToBaseVector<reco::Track> tc(tCH);
    for (unsigned int j=0; j<tCH->size();j++)
      tc.push_back(edm::RefToBase<reco::Track>(tCH,j));

    edm::RefVector<reco::GenParticleCollection> tpc(tPCH.id());
    for (unsigned int j=0; j<tPCH->size();j++)
      tpc.push_back(edm::Ref<reco::GenParticleCollection>(tPCH,j));

    return associateGenToReco(tc,tpc,event,setup);
  }  
RecoToGenCollection TrackAssociatorByChi2::associateRecoToGen ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< reco::GenParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const

Association Sim To Reco with Collections (Gen Particle version)

Definition at line 288 of file TrackAssociatorByChi2.cc.

References edm::RefToBaseVector< T >::begin(), edm::RefToBaseVector< T >::end(), edm::Event::getByLabel(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), and edm::RefVector< C, T, F >::size().

Referenced by MultiTrackValidatorGenPs::analyze(), and associateRecoToGen().

                                                                                                {
  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
  e->getByLabel(bsSrc,recoBeamSpotHandle);
  reco::BeamSpot bs = *recoBeamSpotHandle;      

  RecoToGenCollection  outputCollection;

  GenParticleCollection tPC;
  if (tPCH.size()!=0)  tPC = *tPCH.product();

  int tindex=0;
  for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){

    LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
                                << "rec::Track #"<<tindex<<" with pt=" << (*rt)->pt() <<  "\n"
                                << "===========================================" << "\n";
 
    TrackBase::ParameterVector rParameters = (*rt)->parameters();

    TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
    if (onlyDiagonal){
      for (unsigned int i=0;i<5;i++){
        for (unsigned int j=0;j<5;j++){
          if (i!=j) recoTrackCovMatrix(i,j)=0;
        }
      }
    } 

    recoTrackCovMatrix.Invert();

    int tpindex =0;
    for (GenParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
        
      //skip tps with a very small pt
      //if (sqrt(tp->momentum().perp2())<0.5) continue;
      int charge = tp->charge();
      if (charge==0) continue;
      Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
      Basic3DVector<double> vert=(Basic3DVector<double>) tp->vertex();

      double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
      
      if (chi2<chi2cut) {
        outputCollection.insert(tC[tindex], 
                                std::make_pair(edm::Ref<GenParticleCollection>(tPCH, tpindex),
                                               -chi2));//-chi2 because the Association Map is ordered using std::greater
      }
    }
  }
  outputCollection.post_insert();
  return outputCollection;
}
virtual reco::RecoToGenCollection TrackAssociatorByChi2::associateRecoToGen ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< reco::GenParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const [inline, virtual]

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

Definition at line 143 of file TrackAssociatorByChi2.h.

References associateRecoToGen(), event(), edm::HandleBase::id(), j, edm::RefVector< C, T, F >::push_back(), edm::RefToBaseVector< T >::push_back(), and HcalObjRepresent::setup().

                                                                                              {
    edm::RefToBaseVector<reco::Track> tc(tCH);
    for (unsigned int j=0; j<tCH->size();j++)
      tc.push_back(edm::RefToBase<reco::Track>(tCH,j));

    edm::RefVector<reco::GenParticleCollection> tpc(tPCH.id());
    for (unsigned int j=0; j<tPCH->size();j++)
      tpc.push_back(edm::Ref<reco::GenParticleCollection>(tPCH,j));

    return associateRecoToGen(tc,tpc,event,setup);
  }
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(), and associateRecoToSim().

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

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

Reimplemented from TrackAssociatorBase.

Definition at line 116 of file TrackAssociatorByChi2.h.

References associateRecoToSim(), event(), and HcalObjRepresent::setup().

                                                                                      {
    return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event,setup);
  }
RecoToSimCollection TrackAssociatorByChi2::associateRecoToSim ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const [virtual]

Association Reco To Sim with Collections.

Implements TrackAssociatorBase.

Definition at line 173 of file TrackAssociatorByChi2.cc.

References edm::RefToBaseVector< T >::begin(), edm::RefToBaseVector< T >::end(), edm::Event::getByLabel(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), and edm::RefVector< C, T, F >::size().

                                                                                                {
  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
  e->getByLabel(bsSrc,recoBeamSpotHandle);
  reco::BeamSpot bs = *recoBeamSpotHandle;      

  RecoToSimCollection  outputCollection;

  TrackingParticleCollection tPC;
  if (tPCH.size()!=0)  tPC = *tPCH.product();

  int tindex=0;
  for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){

    LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
                                << "rec::Track #"<<tindex<<" with pt=" << (*rt)->pt() <<  "\n"
                                << "===========================================" << "\n";
 
    TrackBase::ParameterVector rParameters = (*rt)->parameters();

    TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
    if (onlyDiagonal){
      for (unsigned int i=0;i<5;i++){
        for (unsigned int j=0;j<5;j++){
          if (i!=j) recoTrackCovMatrix(i,j)=0;
        }
      }
    } 

    recoTrackCovMatrix.Invert();

    int tpindex =0;
    for (TrackingParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
        
      //skip tps with a very small pt
      //if (sqrt(tp->momentum().perp2())<0.5) continue;
      int charge = tp->charge();
      if (charge==0) continue;
      Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
      Basic3DVector<double> vert=(Basic3DVector<double>) tp->vertex();

      double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
      
      if (chi2<chi2cut) {
        outputCollection.insert(tC[tindex], 
                                std::make_pair(edm::Ref<TrackingParticleCollection>(tPCH, tpindex),
                                               -chi2));//-chi2 because the Association Map is ordered using std::greater
      }
    }
  }
  outputCollection.post_insert();
  return outputCollection;
}
reco::SimToRecoCollection TrackAssociatorByChi2::associateSimToReco ( edm::Handle< edm::View< reco::Track > > &  tCH,
edm::Handle< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const [inline, virtual]

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

Reimplemented from TrackAssociatorBase.

Definition at line 124 of file TrackAssociatorByChi2.h.

References associateSimToReco(), event(), and HcalObjRepresent::setup().

                                                                                      {
    return TrackAssociatorBase::associateSimToReco(tCH,tPCH,event,setup);
  }  
SimToRecoCollection TrackAssociatorByChi2::associateSimToReco ( const edm::RefToBaseVector< reco::Track > &  tC,
const edm::RefVector< TrackingParticleCollection > &  tPCH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const [virtual]

Association Sim To Reco with Collections.

Implements TrackAssociatorBase.

Definition at line 230 of file TrackAssociatorByChi2.cc.

References edm::RefToBaseVector< T >::begin(), edm::RefToBaseVector< T >::end(), edm::Event::getByLabel(), i, edm::AssociationMap< Tag >::insert(), j, LogDebug, edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), edm::RefVector< C, T, F >::size(), and mathSSE::sqrt().

Referenced by associateSimToReco().

                                                                                                 {
  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
  e->getByLabel(bsSrc,recoBeamSpotHandle);
  reco::BeamSpot bs = *recoBeamSpotHandle;      

  SimToRecoCollection  outputCollection;

  TrackingParticleCollection tPC;
  if (tPCH.size()!=0)  tPC = *tPCH.product();

  int tpindex =0;
  for (TrackingParticleCollection::const_iterator tp=tPC.begin(); tp!=tPC.end(); tp++, ++tpindex){
    
    //skip tps with a very small pt
    //if (sqrt(tp->momentum().perp2())<0.5) continue;
    int charge = tp->charge();
    if (charge==0) continue;
    
    LogDebug("TrackAssociator") << "=========LOOKING FOR ASSOCIATION===========" << "\n"
                                << "TrackingParticle #"<<tpindex<<" with pt=" << sqrt(tp->momentum().perp2()) << "\n"
                                << "===========================================" << "\n";
    
    Basic3DVector<double> momAtVtx(tp->momentum().x(),tp->momentum().y(),tp->momentum().z());
    Basic3DVector<double> vert(tp->vertex().x(),tp->vertex().y(),tp->vertex().z());
      
    int tindex=0;
    for (RefToBaseVector<reco::Track>::const_iterator rt=tC.begin(); rt!=tC.end(); rt++, tindex++){
      
      TrackBase::ParameterVector rParameters = (*rt)->parameters();      
      TrackBase::CovarianceMatrix recoTrackCovMatrix = (*rt)->covariance();
      if (onlyDiagonal) {
        for (unsigned int i=0;i<5;i++){
          for (unsigned int j=0;j<5;j++){
            if (i!=j) recoTrackCovMatrix(i,j)=0;
          }
        }
      }
      recoTrackCovMatrix.Invert();
      
      double chi2 = getChi2(rParameters,recoTrackCovMatrix,momAtVtx,vert,charge,bs);
      
      if (chi2<chi2cut) {
        outputCollection.insert(edm::Ref<TrackingParticleCollection>(tPCH, tpindex),
                                std::make_pair(tC[tindex],
                                               -chi2));//-chi2 because the Association Map is ordered using std::greater
      }
    }
  }
  outputCollection.post_insert();
  return outputCollection;
}
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 41 of file TrackAssociatorByChi2.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi, f, i, and j.

                                                                       {
  
  RecoToSimPairAssociation outputVec;

  for (TrackCollection::const_iterator track=rtColl.begin(); track!=rtColl.end(); track++){
     Chi2SimMap outMap;

    TrackBase::ParameterVector rParameters = track->parameters();

    TrackBase::CovarianceMatrix recoTrackCovMatrix = track->covariance();
    if (onlyDiagonal){
      for (unsigned int i=0;i<5;i++){
        for (unsigned int j=0;j<5;j++){
          if (i!=j) recoTrackCovMatrix(i,j)=0;
        }
      }
    }
    recoTrackCovMatrix.Invert();

    for (SimTrackContainer::const_iterator st=stColl.begin(); st!=stColl.end(); st++){

      Basic3DVector<double> momAtVtx(st->momentum().x(),st->momentum().y(),st->momentum().z());
      Basic3DVector<double> vert = (Basic3DVector<double>)  svColl[st->vertIndex()].position();

      std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, st->charge(), bs);
      if (params.first){
        TrackBase::ParameterVector sParameters = params.second;
      
        TrackBase::ParameterVector diffParameters = rParameters - sParameters;
        diffParameters[2] = reco::deltaPhi(diffParameters[2],0.f);
        double chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
        chi2/=5;
        if (chi2<chi2cut) outMap[chi2]=*st;
      }
    }
    outputVec.push_back(RecoToSimPair(*track,outMap));
  }
  return outputVec;
}
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

double TrackAssociatorByChi2::getChi2 ( reco::TrackBase::ParameterVector rParameters,
reco::TrackBase::CovarianceMatrix recoTrackCovMatrix,
Basic3DVector< double > &  momAtVtx,
Basic3DVector< double > &  vert,
int &  charge,
const reco::BeamSpot bs 
) const

basic method where chi2 is computed

Definition at line 84 of file TrackAssociatorByChi2.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi, f, LogDebug, and funct::sin().

                                                                   {
  
  double chi2;
  
  std::pair<bool,reco::TrackBase::ParameterVector> params = parametersAtClosestApproach(vert, momAtVtx, charge, bs);
  if (params.first){
    TrackBase::ParameterVector sParameters=params.second;
    
    TrackBase::ParameterVector diffParameters = rParameters - sParameters;
    diffParameters[2] = reco::deltaPhi(diffParameters[2],0.f);
    chi2 = ROOT::Math::Dot(diffParameters * recoTrackCovMatrix, diffParameters);
    chi2 /= 5;
    
    LogDebug("TrackAssociator") << "====NEW RECO TRACK WITH PT=" << sin(rParameters[1])*float(charge)/rParameters[0] << "====\n" 
                                << "qoverp sim: " << sParameters[0] << "\n" 
                                << "lambda sim: " << sParameters[1] << "\n" 
                                << "phi    sim: " << sParameters[2] << "\n" 
                                << "dxy    sim: " << sParameters[3] << "\n" 
                                << "dsz    sim: " << sParameters[4] << "\n" 
                                << ": " /*<< */ << "\n" 
                                << "qoverp rec: " << rParameters[0] << "\n" 
                                << "lambda rec: " << rParameters[1] << "\n" 
                                << "phi    rec: " << rParameters[2] << "\n" 
                                << "dxy    rec: " << rParameters[3] << "\n" 
                                << "dsz    rec: " << rParameters[4] << "\n" 
                                << ": " /*<< */ << "\n" 
                                << "chi2: " << chi2 << "\n";
    
    return chi2;  
  } else {
    return 10000000000.;
  }
}
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 145 of file TrackAssociatorByChi2.cc.

References FreeTrajectoryState::charge(), funct::cos(), Geom::halfPi(), PV3DBase< T, PVType, FrameType >::mag(), FreeTrajectoryState::momentum(), AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), FreeTrajectoryState::position(), funct::sin(), 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 VertexFitterResult::fill().

                                                                            {
  
  TrackBase::ParameterVector sParameters;
  try {
    FreeTrajectoryState ftsAtProduction(GlobalPoint(vertex.x(),vertex.y(),vertex.z()),
                                        GlobalVector(momAtVtx.x(),momAtVtx.y(),momAtVtx.z()),
                                        TrackCharge(charge),
                                        theMF.product());
    TSCBLBuilderNoMaterial tscblBuilder;
    TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,bs);//as in TrackProducerAlgorithm
    
    GlobalPoint v = tsAtClosestApproach.trackStateAtPCA().position();
    GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
    sParameters[0] = tsAtClosestApproach.trackStateAtPCA().charge()/p.mag();
    sParameters[1] = Geom::halfPi() - p.theta();
    sParameters[2] = p.phi();
    sParameters[3] = (-v.x()*sin(p.phi())+v.y()*cos(p.phi()));
    sParameters[4] = v.z()*p.perp()/p.mag() - (v.x()*p.x()+v.y()*p.y())/p.perp() * p.z()/p.mag();
    
    return pair<bool,TrackBase::ParameterVector>(true,sParameters);
  } catch ( ... ) {
    return pair<bool,TrackBase::ParameterVector>(false,sParameters);
  }
}

Member Data Documentation

Definition at line 179 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().

Definition at line 177 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().

Definition at line 178 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().

Definition at line 176 of file TrackAssociatorByChi2.h.

Referenced by TrackAssociatorByChi2().