CMS 3D CMS Logo

Public Types | Public Member Functions | Private Member Functions | Private Attributes

TrackAssociatorByHits Class Reference

#include <TrackAssociatorByHits.h>

Inheritance diagram for TrackAssociatorByHits:
TrackAssociatorBase

List of all members.

Public Types

enum  SimToRecoDenomType { denomnone, denomsim, denomreco }

Public Member Functions

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
reco::RecoToSimCollectionSeed associateRecoToSim (edm::Handle< edm::View< TrajectorySeed > > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
reco::SimToRecoCollectionSeed associateSimToReco (edm::Handle< edm::View< TrajectorySeed > > &, edm::Handle< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
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
template<typename iter >
int getDoubleCount (iter, iter, TrackerHitAssociator *, TrackingParticleCollection::const_iterator) const
template<typename iter >
void getMatchedIds (std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, int &, iter, iter, TrackerHitAssociator *) const
int getShared (std::vector< SimHitIdpr > &, std::vector< SimHitIdpr > &, TrackingParticleCollection::const_iterator) const
 TrackAssociatorByHits (const edm::ParameterSet &)
 ~TrackAssociatorByHits ()

Private Member Functions

const TrackingRecHitgetHitPtr (edm::OwnVector< TrackingRecHit >::const_iterator iter) const
const TrackingRecHitgetHitPtr (trackingRecHit_iterator iter) const
int LayerFromDetid (const DetId &) const

Private Attributes

const bool AbsoluteNumberOfHits
const edm::ParameterSetconf_
const double cut_RecoToSim
const double purity_SimToReco
const double quality_SimToReco
SimToRecoDenomType SimToRecoDenominator
const bool ThreeHitTracksAreSpecial
const bool UseGrouped
const bool UsePixels
const bool UseSplitting

Detailed Description

Definition at line 21 of file TrackAssociatorByHits.h.


Member Enumeration Documentation

Enumerator:
denomnone 
denomsim 
denomreco 

Definition at line 25 of file TrackAssociatorByHits.h.


Constructor & Destructor Documentation

TrackAssociatorByHits::TrackAssociatorByHits ( const edm::ParameterSet conf) [explicit]

Definition at line 38 of file TrackAssociatorByHits.cc.

References conf_, denomnone, denomreco, denomsim, Exception, edm::ParameterSet::getParameter(), SimToRecoDenominator, and tmp.

                                                                         :  
  conf_(conf),
  AbsoluteNumberOfHits(conf_.getParameter<bool>("AbsoluteNumberOfHits")),
  SimToRecoDenominator(denomnone),
  quality_SimToReco(conf_.getParameter<double>("Quality_SimToReco")),
  purity_SimToReco(conf_.getParameter<double>("Purity_SimToReco")),
  cut_RecoToSim(conf_.getParameter<double>("Cut_RecoToSim")),
  UsePixels(conf_.getParameter<bool>("UsePixels")),
  UseGrouped(conf_.getParameter<bool>("UseGrouped")),
  UseSplitting(conf_.getParameter<bool>("UseSplitting")),
  ThreeHitTracksAreSpecial(conf_.getParameter<bool>("ThreeHitTracksAreSpecial"))
{
  std::string tmp = conf_.getParameter<string>("SimToRecoDenominator");
  if (tmp=="sim") {
    SimToRecoDenominator = denomsim;
  } else if (tmp == "reco") {
    SimToRecoDenominator = denomreco;
  } 

  if (SimToRecoDenominator == denomnone) {
    throw cms::Exception("TrackAssociatorByHits") << "SimToRecoDenominator not specified as sim or reco";
  }

}
TrackAssociatorByHits::~TrackAssociatorByHits ( )

Definition at line 65 of file TrackAssociatorByHits.cc.

{
  //do cleanup here
}

Member Function Documentation

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

Association Reco To Sim with Collections.

Implements TrackAssociatorBase.

Definition at line 75 of file TrackAssociatorByHits.cc.

References abs, AbsoluteNumberOfHits, edm::RefToBaseVector< T >::begin(), conf_, cut_RecoToSim, edm::RefToBaseVector< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), edm::RefVector< C, T, F >::size(), lumiQTWidget::t, and ThreeHitTracksAreSpecial.

Referenced by associateRecoToSim().

                                                                            {

  //edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateRecoToSim - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
  int nshared = 0;
  float quality=0;//fraction or absolute number of shared hits
  std::vector< SimHitIdpr> SimTrackIds;
  std::vector< SimHitIdpr> matchedIds; 
  RecoToSimCollection  outputCollection;
  
  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
  
  TrackingParticleCollection tPC;
  if (TPCollectionH.size()!=0)  tPC = *const_cast<TrackingParticleCollection*>(TPCollectionH.product());

  //get the ID of the recotrack  by hits 
  int tindex=0;
  for (edm::RefToBaseVector<reco::Track>::const_iterator track=tC.begin(); track!=tC.end(); track++, tindex++){
    matchedIds.clear();
    int ri=0;//valid rechits
    //LogTrace("TrackAssociator") << "\nNEW TRACK - track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found(); 
    getMatchedIds<trackingRecHit_iterator>(matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate);

    //LogTrace("TrackAssociator") << "MATCHED IDS LIST BEGIN" ;
    //for(size_t j=0; j<matchedIds.size(); j++){
    //  LogTrace("TrackAssociator") << "matchedIds[j].first=" << matchedIds[j].first;
    //}
    //LogTrace("TrackAssociator") << "MATCHED IDS LIST END" ;
    //LogTrace("TrackAssociator") << "#matched ids=" << matchedIds.size() << " #tps=" << tPC.size();

    //save id for the track
    std::vector<SimHitIdpr> idcachev;
    if(!matchedIds.empty()){

      int tpindex =0;
      for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
        //int nsimhit = t->trackPSimHit(DetId::Tracker).size(); 
        //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
        idcachev.clear();
        nshared = getShared(matchedIds, idcachev, t);

        //if electron subtract double counting
        if (abs(t->pdgId())==11&&(t->g4Track_end()-t->g4Track_begin())>1){
          nshared-=getDoubleCount<trackingRecHit_iterator>((*track)->recHitsBegin(), (*track)->recHitsEnd(), associate, t);
        }

        if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
        else if(ri!=0) quality = (static_cast<double>(nshared)/static_cast<double>(ri));
        else quality = 0;
        //cut on the fraction
        //float purity = 1.0*nshared/ri;
        if(quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3)){
          //if a track has just 3 hits we require that all 3 hits are shared
          outputCollection.insert(tC[tindex],
                                  std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),
                                                 quality));
//          LogTrace("TrackAssociator") << "reco::Track number " << tindex  << " with #hits=" << ri <<" pt=" << (*track)->pt() 
//                                      << " associated to TP (pdgId, nb segments, p) = " 
//                                      << (*t).pdgId() << " " << (*t).g4Tracks().size() 
//                                      << " " << (*t).momentum() << " #hits=" << nsimhit
//                                      << " TP index=" << tpindex<< " with quality =" << quality;
        } else {
//          LogTrace("TrackAssociator") <<"reco::Track number " << tindex << " with #hits=" << ri 
//                                      << " NOT associated with quality =" << quality;
        }
      }//TP loop
    }
  }
  //LogTrace("TrackAssociator") << "% of Assoc Tracks=" << ((double)outputCollection.size())/((double)tC.size());
  delete associate;
  outputCollection.post_insert();
  return outputCollection;
}
RecoToSimCollectionSeed TrackAssociatorByHits::associateRecoToSim ( edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
edm::Handle< TrackingParticleCollection > &  TPCollectionH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const

Definition at line 332 of file TrackAssociatorByHits.cc.

References abs, AbsoluteNumberOfHits, edm::View< T >::begin(), conf_, cut_RecoToSim, edm::View< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LogTrace, edm::AssociationMap< Tag >::post_insert(), edm::Handle< T >::product(), edm::AssociationMap< Tag >::size(), lumiQTWidget::t, and ThreeHitTracksAreSpecial.

                                                                            {

  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateRecoToSim - #seeds="
                                      <<seedCollectionH->size()<<" #TPs="<<TPCollectionH->size();
  int nshared = 0;
  float quality=0;//fraction or absolute number of shared hits
  std::vector< SimHitIdpr> SimTrackIds;
  std::vector< SimHitIdpr> matchedIds; 
  RecoToSimCollectionSeed  outputCollection;
  
  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
  
  const TrackingParticleCollection tPC   = *(TPCollectionH.product());

  const edm::View<TrajectorySeed> sC = *(seedCollectionH.product()); 
  
  //get the ID of the recotrack  by hits 
  int tindex=0;
  for (edm::View<TrajectorySeed>::const_iterator seed=sC.begin(); seed!=sC.end(); seed++, tindex++) {
    matchedIds.clear();
    int ri=0;//valid rechits
    int nsimhit = seed->recHits().second-seed->recHits().first;
    LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << nsimhit;
    getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(matchedIds, SimTrackIds, ri, seed->recHits().first, seed->recHits().second, associate );

    //save id for the track
    std::vector<SimHitIdpr> idcachev;
    if(!matchedIds.empty()){

      int tpindex =0;
      for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
        LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
        idcachev.clear();
        nshared = getShared(matchedIds, idcachev, t);

        //if electron subtract double counting
        if (abs(t->pdgId())==11&&(t->g4Track_end()-t->g4Track_begin())>1){
          nshared-=getDoubleCount<edm::OwnVector<TrackingRecHit>::const_iterator>(seed->recHits().first, seed->recHits().second, associate, t);
        }
        
        if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
        else if(ri!=0) quality = (static_cast<double>(nshared)/static_cast<double>(ri));
        else quality = 0;
        //cut on the fraction
        if(quality > cut_RecoToSim && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3) ){
          //if a track has just 3 hits we require that all 3 hits are shared
          outputCollection.insert(edm::RefToBase<TrajectorySeed>(seedCollectionH,tindex), 
                                  std::make_pair(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex),quality));
          LogTrace("TrackAssociator") << "Seed number " << tindex << " with #hits=" << ri
                                      << "associated to TP (pdgId, nb segments, p) = " 
                                      << (*t).pdgId() << " " << (*t).g4Tracks().size() 
                                      << " " << (*t).momentum() <<" number " << tpindex << " with quality =" << quality;
        } else {
          LogTrace("TrackAssociator") <<"Seed number " << tindex << " with #hits=" << ri << " NOT associated with quality =" << quality;
        }
      }//TP loop
    }
  }
  LogTrace("TrackAssociator") << "% of Assoc Seeds=" << ((double)outputCollection.size())/((double)seedCollectionH->size());
  delete associate;
  outputCollection.post_insert();
  return outputCollection;
}
reco::RecoToSimCollection TrackAssociatorByHits::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 44 of file TrackAssociatorByHits.h.

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

                                                                                      {
    return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event,setup);
  }
reco::SimToRecoCollection TrackAssociatorByHits::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 52 of file TrackAssociatorByHits.h.

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

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

Association Sim To Reco with Collections.

Implements TrackAssociatorBase.

Definition at line 153 of file TrackAssociatorByHits.cc.

References AbsoluteNumberOfHits, edm::RefToBaseVector< T >::begin(), conf_, denomreco, denomsim, edm::RefToBaseVector< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LayerFromDetid(), SiStripDetId::partnerDetId(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, edm::AssociationMap< Tag >::post_insert(), edm::RefVector< C, T, F >::product(), purity_SimToReco, quality_SimToReco, DetId::rawId(), SimToRecoDenominator, edm::RefVector< C, T, F >::size(), DetId::subdetId(), lumiQTWidget::t, SiStripDetId::TEC, ThreeHitTracksAreSpecial, SiStripDetId::TIB, SiStripDetId::TID, SiStripDetId::TOB, DetId::Tracker, UseGrouped, UsePixels, and UseSplitting.

Referenced by associateSimToReco().

                                                                            {
//  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateSimToReco - #tracks="<<tC.size()<<" #TPs="<<TPCollectionH.size();
  float quality=0;//fraction or absolute number of shared hits
  int nshared = 0;
  std::vector< SimHitIdpr> SimTrackIds;
  std::vector< SimHitIdpr> matchedIds; 
  SimToRecoCollection  outputCollection;

  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
  
  TrackingParticleCollection tPC;
  if (TPCollectionH.size()!=0)  tPC = *const_cast<TrackingParticleCollection*>(TPCollectionH.product());

  //for (TrackingParticleCollection::const_iterator t = tPC.begin(); t != tPC.end(); ++t) {
  //  LogTrace("TrackAssociator") << "NEW TP DUMP";
  //  for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin();g4T !=  t -> g4Track_end(); ++g4T) {
  //    LogTrace("TrackAssociator") << "(*g4T).trackId()=" <<(*g4T).trackId() ;
  //  }
  //}
  
  //get the ID of the recotrack  by hits 
  int tindex=0;
  for (edm::RefToBaseVector<reco::Track>::const_iterator track=tC.begin(); track!=tC.end(); track++, tindex++){
    //LogTrace("TrackAssociator") << "\nNEW TRACK - hits of track number " << tindex <<" with pt =" << (*track)->pt() << " # valid=" << (*track)->found(); 
    int ri=0;//valid rechits
    getMatchedIds<trackingRecHit_iterator>(matchedIds, SimTrackIds, ri, (*track)->recHitsBegin(), (*track)->recHitsEnd(), associate);

    //save id for the track
    std::vector<SimHitIdpr> idcachev;
    if(!matchedIds.empty()){
        
      int tpindex =0;
      for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
        idcachev.clear();
        std::vector<PSimHit> trackerPSimHit( t->trackPSimHit(DetId::Tracker) );
        //int nsimhit = trackerPSimHit.size();
        float totsimhit = 0; 
        std::vector<PSimHit> tphits;
        //LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;

        nshared = getShared(matchedIds, idcachev, t);

        //for(std::vector<PSimHit>::const_iterator TPhit = t->trackerPSimHit_begin(); TPhit != t->trackerPSimHit_end(); TPhit++){
        //  unsigned int detid = TPhit->detUnitId();
        //  DetId detId = DetId(TPhit->detUnitId());
        //  LogTrace("TrackAssociator") <<  " hit trackId= " << TPhit->trackId() << " det ID = " << detid 
        //                                    << " SUBDET = " << detId.subdetId() << " layer = " << LayerFromDetid(detId); 
        //}

        if (nshared!=0) {//do not waste time recounting when it is not needed!!!!

          //count the TP simhit
          //LogTrace("TrackAssociator") << "recounting of tp hits";
          for(std::vector<PSimHit>::const_iterator TPhit = trackerPSimHit.begin(); TPhit != trackerPSimHit.end(); TPhit++){
            DetId dId = DetId(TPhit->detUnitId());
          
            unsigned int subdetId = static_cast<unsigned int>(dId.subdetId());
            if (!UsePixels && (subdetId==PixelSubdetector::PixelBarrel || subdetId==PixelSubdetector::PixelEndcap) )
              continue;

            //unsigned int dRawId = dId.rawId();
            SiStripDetId* stripDetId = 0;
            if (subdetId==SiStripDetId::TIB||subdetId==SiStripDetId::TOB||
                subdetId==SiStripDetId::TID||subdetId==SiStripDetId::TEC)
              stripDetId= new SiStripDetId(dId);
            //LogTrace("TrackAssociator") << "consider hit SUBDET = " << subdetId
            //                            << " layer = " << LayerFromDetid(dId) 
            //                            << " id = " << dId.rawId();
            bool newhit = true;
            for(std::vector<PSimHit>::const_iterator TPhitOK = tphits.begin(); TPhitOK != tphits.end(); TPhitOK++){
              DetId dIdOK = DetId(TPhitOK->detUnitId());
              //unsigned int dRawIdOK = dIdOK.rawId();
              //LogTrace("TrackAssociator") << "\t\tcompare with SUBDET = " << dIdOK.subdetId()
              //                            << " layer = " << LayerFromDetid(dIdOK)
              //                            << " id = " << dIdOK.rawId();
              //no grouped, no splitting
              if (!UseGrouped && !UseSplitting)
                if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
                    dId.subdetId()==dIdOK.subdetId()) newhit = false;
              //no grouped, splitting
              if (!UseGrouped && UseSplitting)
                if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
                    dId.subdetId()==dIdOK.subdetId() &&
                    (stripDetId==0 || stripDetId->partnerDetId()!=dIdOK.rawId()))
                  newhit = false;
              //grouped, no splitting
              if (UseGrouped && !UseSplitting)
                if (LayerFromDetid(dId)==LayerFromDetid(dIdOK) &&
                    dId.subdetId()==dIdOK.subdetId() &&
                    stripDetId!=0 && stripDetId->partnerDetId()==dIdOK.rawId())
                  newhit = false;
              //grouped, splitting
              if (UseGrouped && UseSplitting)
                newhit = true;
            }
            if (newhit) {
              //LogTrace("TrackAssociator") << "\t\tok";
              tphits.push_back(*TPhit);
            }
            //else LogTrace("TrackAssociator") << "\t\tno";
            delete stripDetId;
          }
          totsimhit = tphits.size();
        }

        if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
        else if(SimToRecoDenominator == denomsim && totsimhit!=0) quality = ((double) nshared)/((double)totsimhit);
        else if(SimToRecoDenominator == denomreco && ri!=0) quality = ((double) nshared)/((double)ri);
        else quality = 0;
        //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit << " re-counted = " << totsimhit 
        //<< " nshared = " << nshared << " nrechit = " << ri;
        
        float purity = 1.0*nshared/ri;
        if (quality>quality_SimToReco && !(ThreeHitTracksAreSpecial && totsimhit==3 && nshared<3) && (AbsoluteNumberOfHits||(purity>purity_SimToReco))) {
          //if a track has just 3 hits we require that all 3 hits are shared
          outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex), 
                                  std::make_pair(tC[tindex],quality));
//          LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit 
//                                      << " re-counted = "  << totsimhit << " nshared = " << nshared 
//                                      << " associated to track number " << tindex << " with pt=" << (*track)->pt() 
//                                      << " with hit quality =" << quality ;
        } else {
//          LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit 
//                                      << " re-counted = "  << totsimhit << " nshared = " << nshared 
//                                      << " NOT associated with quality =" << quality;
        }
      }
    }
  }
  //LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH.size());
  delete associate;
  outputCollection.post_insert();
  return outputCollection;
}
SimToRecoCollectionSeed TrackAssociatorByHits::associateSimToReco ( edm::Handle< edm::View< TrajectorySeed > > &  seedCollectionH,
edm::Handle< TrackingParticleCollection > &  TPCollectionH,
const edm::Event event = 0,
const edm::EventSetup setup = 0 
) const

Definition at line 401 of file TrackAssociatorByHits.cc.

References AbsoluteNumberOfHits, edm::View< T >::begin(), conf_, edm::View< T >::end(), getShared(), edm::AssociationMap< Tag >::insert(), LogTrace, edm::AssociationMap< Tag >::post_insert(), edm::Handle< T >::product(), quality_SimToReco, edm::AssociationMap< Tag >::size(), lumiQTWidget::t, ThreeHitTracksAreSpecial, and DetId::Tracker.

                                                                            {

  edm::LogVerbatim("TrackAssociator") << "Starting TrackAssociatorByHits::associateSimToReco - #seeds="
                                      <<seedCollectionH->size()<<" #TPs="<<TPCollectionH->size();
  float quality=0;//fraction or absolute number of shared hits
  int nshared = 0;
  std::vector< SimHitIdpr> SimTrackIds;
  std::vector< SimHitIdpr> matchedIds; 
  SimToRecoCollectionSeed  outputCollection;

  TrackerHitAssociator * associate = new TrackerHitAssociator(*e, conf_);
  
  TrackingParticleCollection tPC =*const_cast<TrackingParticleCollection*>(TPCollectionH.product());

  const edm::View<TrajectorySeed> sC = *(seedCollectionH.product()); 

  //get the ID of the recotrack  by hits 
  int tindex=0;
  for (edm::View<TrajectorySeed>::const_iterator seed=sC.begin(); seed!=sC.end(); seed++, tindex++) {
    int ri=0;//valid rechits
    LogTrace("TrackAssociator") << "\nNEW SEED - seed number " << tindex << " # valid=" << seed->recHits().second-seed->recHits().first;
    getMatchedIds<edm::OwnVector<TrackingRecHit>::const_iterator>(matchedIds, SimTrackIds, ri, seed->recHits().first, seed->recHits().second, associate );

    //save id for the track
    std::vector<SimHitIdpr> idcachev;
    if(!matchedIds.empty()){
      int tpindex =0;
      for (TrackingParticleCollection::iterator t = tPC.begin(); t != tPC.end(); ++t, ++tpindex) {
        idcachev.clear();
        int nsimhit = t->trackPSimHit(DetId::Tracker).size(); 
        LogTrace("TrackAssociator") << "TP number " << tpindex << " pdgId=" << t->pdgId() << " with number of PSimHits: "  << nsimhit;
        nshared = getShared(matchedIds, idcachev, t);
        
        if (AbsoluteNumberOfHits) quality = static_cast<double>(nshared);
        else if(ri!=0) quality = ((double) nshared)/((double)ri);
        else quality = 0;
        //LogTrace("TrackAssociator") << "Final count: nhit(TP) = " << nsimhit 
        //<< " nshared = " << nshared 
        //<< " nrechit = " << ri;
        if(quality > quality_SimToReco && !(ThreeHitTracksAreSpecial && ri==3 && nshared<3) ){
          outputCollection.insert(edm::Ref<TrackingParticleCollection>(TPCollectionH, tpindex), 
                                  std::make_pair(edm::RefToBase<TrajectorySeed>(seedCollectionH,tindex), quality));
          LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit
                                      << " associated to seed number " << tindex << " with #hits=" << ri 
                                      << " with hit quality =" << quality ;
        }
        else {
          LogTrace("TrackAssociator") << "TrackingParticle number " << tpindex << " with #hits=" << nsimhit << " NOT associated with quality =" << quality;
        }
      }
    }
  }
  LogTrace("TrackAssociator") << "% of Assoc TPs=" << ((double)outputCollection.size())/((double)TPCollectionH->size());
  delete associate;
  outputCollection.post_insert();
  return outputCollection;
}
template<typename iter >
int TrackAssociatorByHits::getDoubleCount ( iter  begin,
iter  end,
TrackerHitAssociator associate,
TrackingParticleCollection::const_iterator  t 
) const

Definition at line 546 of file TrackAssociatorByHits.cc.

References TrackerHitAssociator::associateHitId(), end, spr::find(), and getHitPtr().

                                                                                            {
  int doublecount = 0 ;
  std::vector<SimHitIdpr> SimTrackIdsDC;
  //  cout<<begin-end<<endl;
  for (iter it = begin;  it != end; it++){
    int idcount = 0;
    SimTrackIdsDC.clear();
    associate->associateHitId(*getHitPtr(it), SimTrackIdsDC);
    //    cout<<SimTrackIdsDC.size()<<endl;
    if(SimTrackIdsDC.size()>1){
      //     cout<<(t->g4Track_end()-t->g4Track_begin())<<endl;
      for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin(); g4T !=  t -> g4Track_end(); ++g4T) {
        if(find(SimTrackIdsDC.begin(), SimTrackIdsDC.end(),SimHitIdpr((*g4T).trackId(), SimTrackIdsDC.begin()->second)) != SimTrackIdsDC.end() ){
          idcount++;
        }
      }
    }
    if (idcount>1) doublecount+=(idcount-1);
  }
  return doublecount;
}
const TrackingRecHit* TrackAssociatorByHits::getHitPtr ( trackingRecHit_iterator  iter) const [inline, private]

Definition at line 99 of file TrackAssociatorByHits.h.

{return &**iter;}
const TrackingRecHit* TrackAssociatorByHits::getHitPtr ( edm::OwnVector< TrackingRecHit >::const_iterator  iter) const [inline, private]

Definition at line 98 of file TrackAssociatorByHits.h.

Referenced by getDoubleCount(), and getMatchedIds().

{return &*iter;}
template<typename iter >
void TrackAssociatorByHits::getMatchedIds ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  SimTrackIds,
int &  ri,
iter  begin,
iter  end,
TrackerHitAssociator associate 
) const

Definition at line 463 of file TrackAssociatorByHits.cc.

References TrackerHitAssociator::associateHitId(), end, TrackingRecHit::geographicalId(), getHitPtr(), TrackingRecHit::isValid(), j, LogTrace, and DetId::rawId().

                                                                                  {
    matchedIds.clear();
    ri=0;//valid rechits
    for (iter it = begin;  it != end; it++){
      const TrackingRecHit *hit=getHitPtr(it);
      if (hit->isValid()){
        ri++;
        uint32_t t_detID=  hit->geographicalId().rawId();
        SimTrackIds.clear();      
        associate->associateHitId(*hit, SimTrackIds);
        //save all the id of matched simtracks
        if(!SimTrackIds.empty()){
         for(size_t j=0; j<SimTrackIds.size(); j++){
           LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << hit->isValid() 
                                       << " det id = " << t_detID << " SimId " << SimTrackIds[j].first 
                                       << " evt=" << SimTrackIds[j].second.event() 
                                       << " bc=" << SimTrackIds[j].second.bunchCrossing();  
           matchedIds.push_back(SimTrackIds[j]);                        
         }
        }
        //****to be used when the crossing frame is in the event and with flag TrackAssociatorByHitsESProducer.associateRecoTracks = false
        //std::vector<PSimHit> tmpSimHits = associate->associateHit(*getHitPtr(it));
        //for(size_t j=0; j<tmpSimHits.size(); j++) {
        //  LogTrace("TrackAssociator") << " hit # " << ri << " valid=" << getHitPtr(it)->isValid()
        //                              << " det id = " << t_detID << " SimId " << SimTrackIds[j].first
        //                              << " evt=" << SimTrackIds[j].second.event()
        //                              << " bc=" << SimTrackIds[j].second.bunchCrossing()
        //                              << " process=" << tmpSimHits[j].processType()
        //                              << " eloss=" << tmpSimHits[j].energyLoss()
        //                              << " part type=" << tmpSimHits[j].particleType()
        //                              << " pabs=" << tmpSimHits[j].pabs()
        //                              << " trId=" << tmpSimHits[j].trackId();
        //}
      }else{
        LogTrace("TrackAssociator") <<"\t\t Invalid Hit On "<<hit->geographicalId().rawId();
      }
    }//trackingRecHit loop
}
int TrackAssociatorByHits::getShared ( std::vector< SimHitIdpr > &  matchedIds,
std::vector< SimHitIdpr > &  idcachev,
TrackingParticleCollection::const_iterator  t 
) const

Definition at line 511 of file TrackAssociatorByHits.cc.

References prof2calltree::count, spr::find(), and j.

Referenced by associateRecoToSim(), and associateSimToReco().

                                                                                       {
  int nshared = 0;
  if (t->trackPSimHit().size()==0) return nshared;//should use trackerPSimHit but is not const

  for(size_t j=0; j<matchedIds.size(); j++){
    //LogTrace("TrackAssociator") << "now matchedId=" << matchedIds[j].first;
    if(find(idcachev.begin(), idcachev.end(),matchedIds[j]) == idcachev.end() ){
      //only the first time we see this ID 
      idcachev.push_back(matchedIds[j]);
      
      for (TrackingParticle::g4t_iterator g4T = t -> g4Track_begin(); g4T !=  t -> g4Track_end(); ++g4T) {
//      LogTrace("TrackAssociator") << " TP   (ID, Ev, BC) = " << (*g4T).trackId() 
//                                  << ", " << t->eventId().event() << ", "<< t->eventId().bunchCrossing()
//                                  << " Match(ID, Ev, BC) = " <<  matchedIds[j].first
//                                  << ", " << matchedIds[j].second.event() << ", "
//                                  << matchedIds[j].second.bunchCrossing() ;
                                    //<< "\t G4  Track Momentum " << (*g4T).momentum() 
                                    //<< " \t reco Track Momentum " << track->momentum();             
        if((*g4T).trackId() == matchedIds[j].first && t->eventId() == matchedIds[j].second){
                int countedhits = std::count(matchedIds.begin(), matchedIds.end(), matchedIds[j]);
                nshared += countedhits;

//              LogTrace("TrackAssociator") << "hits shared by this segment : " << countedhits;
//              LogTrace("TrackAssociator") << "hits shared so far : " << nshared;
        }
      }//g4Tracks loop
    }
  }
  return nshared;
}
int TrackAssociatorByHits::LayerFromDetid ( const DetId detId) const [private]

Definition at line 291 of file TrackAssociatorByHits.cc.

References PXFDetId::disk(), TIBDetId::layer(), TOBDetId::layer(), PXBDetId::layer(), align::tib::layerNumber(), LogTrace, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, DetId::rawId(), DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, StripSubdetector::TOB, TIDDetId::wheel(), and TECDetId::wheel().

Referenced by associateSimToReco().

{
  int layerNumber=0;
  unsigned int subdetId = static_cast<unsigned int>(detId.subdetId()); 
  if ( subdetId == StripSubdetector::TIB) 
    { 
      TIBDetId tibid(detId.rawId()); 
      layerNumber = tibid.layer();
    }
  else if ( subdetId ==  StripSubdetector::TOB )
    { 
      TOBDetId tobid(detId.rawId()); 
      layerNumber = tobid.layer();
    }
  else if ( subdetId ==  StripSubdetector::TID) 
    { 
      TIDDetId tidid(detId.rawId());
      layerNumber = tidid.wheel();
    }
  else if ( subdetId ==  StripSubdetector::TEC )
    { 
      TECDetId tecid(detId.rawId()); 
      layerNumber = tecid.wheel(); 
    }
  else if ( subdetId ==  PixelSubdetector::PixelBarrel ) 
    { 
      PXBDetId pxbid(detId.rawId()); 
      layerNumber = pxbid.layer();  
    }
  else if ( subdetId ==  PixelSubdetector::PixelEndcap ) 
    { 
      PXFDetId pxfid(detId.rawId()); 
      layerNumber = pxfid.disk();  
    }
  else
    LogTrace("TrackAssociator") << "Unknown subdetid: " <<  subdetId;
  
  return layerNumber;
} 

Member Data Documentation

Definition at line 87 of file TrackAssociatorByHits.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const double TrackAssociatorByHits::cut_RecoToSim [private]

Definition at line 91 of file TrackAssociatorByHits.h.

Referenced by associateRecoToSim().

Definition at line 90 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

Definition at line 89 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

Definition at line 88 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco(), and TrackAssociatorByHits().

Definition at line 95 of file TrackAssociatorByHits.h.

Referenced by associateRecoToSim(), and associateSimToReco().

const bool TrackAssociatorByHits::UseGrouped [private]

Definition at line 93 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

const bool TrackAssociatorByHits::UsePixels [private]

Definition at line 92 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().

Definition at line 94 of file TrackAssociatorByHits.h.

Referenced by associateSimToReco().