![]() |
![]() |
#include <Alignment/TrackHitFilter/src/TrackHitFilter.cc>
Public Member Functions | |
TrackHitFilter (const edm::ParameterSet &) | |
~TrackHitFilter () | |
Protected Attributes | |
bool | rejectBadClusterPixelHits |
bool | rejectBadMods |
bool | rejectBadStoNHits |
std::vector< unsigned int > | theBadMods |
std::string | theCMNSubtractionMode |
std::string | theHitSel |
unsigned int | theMinHits |
double | thePixelClusterthreshold |
edm::InputTag | theSrc |
double | theStoNthreshold |
Private Member Functions | |
virtual void | beginJob () |
virtual void | endJob () |
bool | keepThisHit (DetId id, int type, int layer, const TrackingRecHit *, const edm::Event &, const edm::EventSetup &) |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Description: Selects some track hits for refitting input
Implementation: <Notes on="" implementation>="">
Definition at line 36 of file TrackHitFilter.h.
TrackHitFilter::TrackHitFilter | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 46 of file TrackHitFilter.cc.
: theSrc( iConfig.getParameter<edm::InputTag>( "src" ) ), theHitSel( iConfig.getParameter<std::string>( "hitSelection" ) ), theMinHits( iConfig.getParameter<unsigned int>( "minHitsForRefit" ) ), rejectBadMods( iConfig.getParameter<bool>( "rejectBadMods" ) ), theBadMods( iConfig.getParameter<std::vector<unsigned int> >( "theBadModules" ) ), rejectBadStoNHits( iConfig.getParameter<bool>( "rejectBadStoNHits" ) ), theCMNSubtractionMode(iConfig.getUntrackedParameter<std::string>( "CMNSubtractionMode" ,"Median") ) , theStoNthreshold( iConfig.getParameter<double>( "theStoNthreshold" ) ), rejectBadClusterPixelHits( iConfig.getParameter<bool>( "rejectBadClusterPixelHits" ) ), thePixelClusterthreshold( iConfig.getParameter<double>( "thePixelClusterthreshold" ) ) { //register your products, and/or set an "alias" label produces<reco::TrackCollection>(); produces<reco::TrackExtraCollection>(); produces<TrackingRecHitCollection>(); }
TrackHitFilter::~TrackHitFilter | ( | ) |
Definition at line 67 of file TrackHitFilter.cc.
{ }
void TrackHitFilter::beginJob | ( | void | ) | [private, virtual] |
void TrackHitFilter::endJob | ( | void | ) | [private, virtual] |
bool TrackHitFilter::keepThisHit | ( | DetId | id, |
int | type, | ||
int | layer, | ||
const TrackingRecHit * | therechit, | ||
const edm::Event & | iEvent, | ||
const edm::EventSetup & | iSetup | ||
) | [private] |
Definition at line 187 of file TrackHitFilter.cc.
References abs, SiPixelCluster::charge(), SiStripRecHit2D::cluster(), SiPixelRecHit::cluster(), LogDebug, SiStripMatchedRecHit2D::monoHit(), ProjectedSiStripRecHit2D::originalHit(), rejectBadClusterPixelHits, rejectBadMods, rejectBadStoNHits, SiStripClusterInfo::signalOverNoise(), SiPixelCluster::size(), SiStripMatchedRecHit2D::stereoHit(), theBadMods, theHitSel, thePixelClusterthreshold, and theStoNthreshold.
Referenced by produce().
{ bool keepthishit = true; if ( theHitSel == "PixelOnly" ) { if (abs(type)<1 || abs(type)>2) keepthishit = false; } else if ( theHitSel == "PixelBarrelOnly" ) { if (abs(type)!=1) keepthishit = false; } else if ( theHitSel == "PixelAndDSStripBarrelOnly" ) { if (!((abs(type)==1) || ((abs(type)==3 || abs(type)==5) && layer<2))) keepthishit = false; } else if ( theHitSel == "SiStripOnly" ) { if (abs(type)>=1 && abs(type)<=2) keepthishit = false; } else if ( theHitSel == "TOBOnly" ) { if (abs(type)!=5) keepthishit = false; } else if ( theHitSel == "TOBandTIBOnly" ) { if ( !(abs(type)==3 || abs(type)==5) ) keepthishit = false; } else if ( theHitSel == "TOBandTIBl4Only" ) { if (!(abs(type)==5 || (abs(type)==3 && layer>=4))) keepthishit = false; } else if ( theHitSel == "TOBandTIBl34Only" ) { if (!(abs(type)==5 || (abs(type)==3 && layer>=3))) keepthishit = false; } else if ( theHitSel == "TOBandTIBl234Only" ) { if (!(abs(type)==5 || (abs(type)==3 && layer>=2))) keepthishit = false; } else if ( theHitSel == "noPXE" ) { if (abs(type)==2) keepthishit = false; } if (rejectBadMods) { for( std::vector<unsigned int>::const_iterator iMod = theBadMods.begin(); iMod != theBadMods.end(); ++iMod ) { if (id.rawId() == *iMod) { keepthishit = false; LogDebug("HitFilter") << "TrackHitFilter **** Rejected a hit on module " << id.rawId(); } } } //********* Rejects pixel hits with 'bad' cluster if(rejectBadClusterPixelHits && (abs(type)==1 || abs(type)==2)){ const SiPixelRecHit* pixelhit = dynamic_cast<const SiPixelRecHit*>(therechit); const SiPixelCluster* pixelcluster = &*(pixelhit->cluster()); if(pixelcluster->size()==1 && pixelcluster->charge()<thePixelClusterthreshold){keepthishit = false;} } // **************** // Reject hits with bad S/N if (rejectBadStoNHits && (abs(type)>2) ) { //apply it only to Strip hits /*** RC ****/ // const uint32_t& recHitDetId = id.rawId(); const SiStripMatchedRecHit2D* matchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(therechit); const SiStripRecHit2D* hit = dynamic_cast<const SiStripRecHit2D*>(therechit); const ProjectedSiStripRecHit2D* unmatchedhit = dynamic_cast<const ProjectedSiStripRecHit2D*>(therechit); if (matchedhit) { bool keepmonohit=true; bool keepstereohit=true; const SiStripRecHit2D* monohit=matchedhit->monoHit(); const SiStripCluster* monocluster = &*(monohit->cluster()); SiStripClusterInfo monoclusterInfo = SiStripClusterInfo( *monocluster, iSetup); if (monoclusterInfo.signalOverNoise() < theStoNthreshold ) keepmonohit = false; const SiStripRecHit2D* stereohit=matchedhit->stereoHit(); const SiStripCluster* stereocluster = &*(stereohit->cluster()); SiStripClusterInfo stereoclusterInfo = SiStripClusterInfo(*stereocluster, iSetup); if (stereoclusterInfo.signalOverNoise() < theStoNthreshold )keepstereohit = false; if (!keepmonohit || !keepstereohit) keepthishit = false; } else if (hit) { const SiStripCluster* cluster = &*(hit->cluster()); SiStripClusterInfo clusterInfo = SiStripClusterInfo(*cluster, iSetup); if (clusterInfo.signalOverNoise() < theStoNthreshold )keepthishit = false; } else if (unmatchedhit) { const SiStripRecHit2D &orighit = unmatchedhit->originalHit(); const SiStripCluster* origcluster = &*(orighit.cluster()); SiStripClusterInfo clusterInfo = SiStripClusterInfo(*origcluster, iSetup); if (clusterInfo.signalOverNoise() < theStoNthreshold ) keepthishit = false; } } // end reject bad S/N return keepthishit; }
void TrackHitFilter::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 71 of file TrackHitFilter.cc.
References reco::TrackExtraBase::add(), cond::rpcobgas::detid, TrackingRecHit::geographicalId(), edm::Event::getByLabel(), edm::Event::getRefBeforePut(), i, iEvent, reco::Track::innerDetId(), reco::Track::innerMomentum(), reco::Track::innerOk(), reco::Track::innerPosition(), reco::Track::innerStateCovariance(), keepThisHit(), LogDebug, reco::Track::outerDetId(), reco::Track::outerMomentum(), reco::Track::outerOk(), reco::Track::outerPosition(), reco::Track::outerStateCovariance(), edm::Event::put(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), reco::Track::seedDirection(), reco::Track::setExtra(), reco::TrackBase::setHitPattern(), theMinHits, theSrc, and TrackerAlignableId::typeAndLayerFromDetId().
{ //Read track collection from the Event edm::Handle<reco::TrackCollection> trackAllHits; iEvent.getByLabel(theSrc,trackAllHits); // Create empty Track, TrackExtra and TrackingRecHits collections std::auto_ptr<reco::TrackCollection> trackSelectedHits( new reco::TrackCollection() ); std::auto_ptr<reco::TrackExtraCollection> txSelectedHits( new reco::TrackExtraCollection() ); std::auto_ptr<TrackingRecHitCollection> trhSelectedHits( new TrackingRecHitCollection() ); //get references TrackingRecHitRefProd rHits = iEvent.getRefBeforePut<TrackingRecHitCollection>(); reco::TrackExtraRefProd rTrackExtras = iEvent.getRefBeforePut<reco::TrackExtraCollection>(); reco::TrackRefProd rTracks = iEvent.getRefBeforePut<reco::TrackCollection>(); edm::Ref<reco::TrackExtraCollection>::key_type idx = 0; edm::Ref<reco::TrackExtraCollection>::key_type hidx = 0; TrackerAlignableId* TkMap = new TrackerAlignableId(); LogDebug("HitFilter")<< trackAllHits->size() << " track(s) found in the event with label " << theSrc<<std::endl; unsigned int nTr = 0; std::vector<unsigned int> accHits; std::vector<unsigned int> allHits; // first iteration : count hits for( reco::TrackCollection::const_iterator iTrack = trackAllHits->begin(); iTrack != trackAllHits->end(); ++iTrack ) { const reco::Track * trk = &(*iTrack); unsigned int allhits = 0; unsigned int acchits = 0; for (trackingRecHit_iterator iHit = trk->recHitsBegin(); iHit != trk->recHitsEnd(); iHit++) { allhits++; TrackingRecHit* hit = (*iHit)->clone(); //std::cout<<"See if it crashes. Geographical Id of the hit: "<<hit->geographicalId()<<std::endl; std::pair<int,int> typeAndLay = TkMap->typeAndLayerFromDetId( hit->geographicalId() ); int type = typeAndLay.first; int layer = typeAndLay.second; DetId detid= hit->geographicalId(); if (keepThisHit( detid, type, layer, hit, iEvent,iSetup )) acchits++; }//end loop on hits if (!nTr) { LogDebug("HitFilter") << "TrackHitFilter **** In first track " << acchits << " RecHits retained out of " << allhits; nTr++; } allHits.push_back(allhits); accHits.push_back(acchits); }//end loop on tracks nTr = 0; // second iteration : store tracks for( reco::TrackCollection::const_iterator iTrack2 = trackAllHits->begin(); iTrack2 != trackAllHits->end(); ++iTrack2 ) { if (accHits.at(nTr) >= theMinHits) { const reco::Track * trk = &(*iTrack2); reco::Track * myTrk = new reco::Track(*trk); PropagationDirection seedDir = trk->seedDirection(); myTrk->setExtra( reco::TrackExtraRef( rTrackExtras, idx ++ ) ); reco::TrackExtra * tx = new reco::TrackExtra( trk->outerPosition(), trk->outerMomentum(), trk->outerOk(), trk->innerPosition(), trk->innerMomentum(), trk->innerOk(), trk->outerStateCovariance(), trk->outerDetId(), trk->innerStateCovariance(), trk->innerDetId() , seedDir ) ; unsigned int i = 0; for (trackingRecHit_iterator iHit = trk->recHitsBegin(); iHit != trk->recHitsEnd(); iHit++) { TrackingRecHit* hit = (*iHit)->clone(); std::pair<int,int> typeAndLay = TkMap->typeAndLayerFromDetId( hit->geographicalId() ); int type = typeAndLay.first; int layer = typeAndLay.second; if (keepThisHit( hit->geographicalId(), type, layer, hit, iEvent, iSetup )) { myTrk->setHitPattern( * hit, i ++ ); trhSelectedHits->push_back( hit ); tx->add( TrackingRecHitRef( rHits, hidx ++ ) ); } } trackSelectedHits->push_back( *myTrk ); txSelectedHits->push_back( *tx ); } nTr++; } iEvent.put( trackSelectedHits ); iEvent.put( txSelectedHits ); iEvent.put( trhSelectedHits ); // std::cout<<"Finished this call to <TrackHitFilter::produce>"<<std::endl; }
bool TrackHitFilter::rejectBadClusterPixelHits [protected] |
Definition at line 59 of file TrackHitFilter.h.
Referenced by keepThisHit().
bool TrackHitFilter::rejectBadMods [protected] |
Definition at line 52 of file TrackHitFilter.h.
Referenced by keepThisHit().
bool TrackHitFilter::rejectBadStoNHits [protected] |
Definition at line 55 of file TrackHitFilter.h.
Referenced by keepThisHit().
std::vector<unsigned int> TrackHitFilter::theBadMods [protected] |
Definition at line 53 of file TrackHitFilter.h.
Referenced by keepThisHit().
std::string TrackHitFilter::theCMNSubtractionMode [protected] |
Definition at line 56 of file TrackHitFilter.h.
std::string TrackHitFilter::theHitSel [protected] |
Definition at line 50 of file TrackHitFilter.h.
Referenced by keepThisHit().
unsigned int TrackHitFilter::theMinHits [protected] |
Definition at line 51 of file TrackHitFilter.h.
Referenced by produce().
double TrackHitFilter::thePixelClusterthreshold [protected] |
Definition at line 60 of file TrackHitFilter.h.
Referenced by keepThisHit().
edm::InputTag TrackHitFilter::theSrc [protected] |
Definition at line 49 of file TrackHitFilter.h.
Referenced by produce().
double TrackHitFilter::theStoNthreshold [protected] |
Definition at line 57 of file TrackHitFilter.h.
Referenced by keepThisHit().