CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTracker/IterativeTracking/src/QualityFilter.cc

Go to the documentation of this file.
00001 #include "RecoTracker/IterativeTracking/interface/QualityFilter.h"
00002 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00003 #include "DataFormats/TrackReco/interface/TrackExtraFwd.h"
00004 #include "FWCore/Utilities/interface/InputTag.h"
00005 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00006 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00007 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00008 //
00009 // class decleration
00010 //
00011 
00012 using namespace edm;
00013 using namespace reco;
00014 using namespace std;
00015 
00016 
00017 QualityFilter::QualityFilter(const edm::ParameterSet& iConfig)
00018 {
00019   copyExtras_ = iConfig.getUntrackedParameter<bool>("copyExtras", false);
00020 
00021   produces<reco::TrackCollection>();
00022   if (copyExtras_) {
00023       produces<TrackingRecHitCollection>();
00024       produces<reco::TrackExtraCollection>();
00025   }
00026   produces<std::vector<Trajectory> >();
00027   produces<TrajTrackAssociationCollection>();
00028 
00029   tkTag  = iConfig.getParameter<edm::InputTag>("recTracks");
00030   trackQuality_=TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"));
00031 }
00032 
00033 
00034 QualityFilter::~QualityFilter()
00035 {
00036  
00037    // do anything here that needs to be done at desctruction time
00038    // (e.g. close files, deallocate resources etc.)
00039 
00040 }
00041 
00042 
00043 //
00044 // member functions
00045 //
00046 
00047 // ------------ method called to produce the data  ------------
00048 void
00049 QualityFilter::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00050 {
00051 
00052   
00053   auto_ptr<TrackCollection> selTracks(new TrackCollection);
00054   auto_ptr<TrackingRecHitCollection> selHits(copyExtras_ ? new TrackingRecHitCollection() : 0);
00055   auto_ptr<TrackExtraCollection> selTrackExtras(copyExtras_ ? new TrackExtraCollection() : 0);
00056   auto_ptr<vector<Trajectory> > outputTJ(new vector<Trajectory> );
00057   auto_ptr<TrajTrackAssociationCollection> trajTrackMap( new TrajTrackAssociationCollection() );
00058   
00059   
00060   TrackExtraRefProd rTrackExtras; 
00061   TrackingRecHitRefProd rHits;
00062   if (copyExtras_) {
00063     rTrackExtras = iEvent.getRefBeforePut<TrackExtraCollection>();
00064     rHits = iEvent.getRefBeforePut<TrackingRecHitCollection>();
00065   }
00066 
00067   Handle<std::vector<Trajectory> > TrajectoryCollection;
00068   Handle<TrajTrackAssociationCollection> assoMap;
00069   
00070   Handle<TrackCollection> tkCollection;  
00071   iEvent.getByLabel( tkTag, tkCollection);
00072   const reco::TrackCollection*  tC = tkCollection.product();
00073   TrackCollection::const_iterator itxc;
00074   TrackCollection::const_iterator firstTrack = tC->begin();
00075   TrackCollection::const_iterator lastTrack = tC->end();
00076 
00077     
00078   iEvent.getByLabel(tkTag,TrajectoryCollection);
00079   iEvent.getByLabel(tkTag,assoMap);
00080 
00081 
00082 
00083   TrajTrackAssociationCollection::const_iterator it = assoMap->begin();
00084   TrajTrackAssociationCollection::const_iterator lastAssoc = assoMap->end();
00085   for( ; it != lastAssoc; ++it ) {
00086     const Ref<vector<Trajectory> > traj = it->key;
00087     const reco::TrackRef itc = it->val;
00088     bool goodTk = (itc->quality(trackQuality_));
00089  
00090     
00091     if (goodTk){
00092       Track track =(*itc);
00093       //tracks and trajectories
00094       selTracks->push_back( track );
00095       outputTJ->push_back( *traj );
00096       if (copyExtras_) {
00097           //TRACKING HITS
00098           trackingRecHit_iterator irhit   =(*itc).recHitsBegin();
00099           trackingRecHit_iterator lasthit =(*itc).recHitsEnd();
00100           for (; irhit!=lasthit; ++irhit) {
00101             selHits->push_back((*irhit)->clone() );
00102           }
00103       }
00104           
00105     }
00106 
00107   }
00108 
00109   unsigned nTracks = selTracks->size();
00110   if (copyExtras_) {
00111       //PUT TRACKING HITS IN THE EVENT
00112       OrphanHandle<TrackingRecHitCollection> theRecoHits = iEvent.put(selHits );
00113   
00114       //PUT TRACK EXTRA IN THE EVENT
00115       selTrackExtras->reserve(nTracks);
00116       unsigned hits=0;
00117 
00118       for ( unsigned index = 0; index<nTracks; ++index ) { 
00119 
00120         reco::Track& aTrack = selTracks->at(index);
00121         TrackExtra aTrackExtra(aTrack.outerPosition(),
00122                                aTrack.outerMomentum(),
00123                                aTrack.outerOk(),
00124                                aTrack.innerPosition(),
00125                                aTrack.innerMomentum(),
00126                                aTrack.innerOk(),
00127                                aTrack.outerStateCovariance(),
00128                                aTrack.outerDetId(),
00129                                aTrack.innerStateCovariance(),
00130                                aTrack.innerDetId(),
00131                                aTrack.seedDirection(),
00132                                aTrack.seedRef());
00133             
00134 
00135         //unsigned nHits = aTrack.numberOfValidHits();
00136         unsigned nHits = aTrack.recHitsSize();
00137         for ( unsigned int ih=0; ih<nHits; ++ih) {
00138           aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
00139         }
00140         selTrackExtras->push_back(aTrackExtra);
00141       }
00142 
00143       //CORRECT REF TO TRACK
00144       OrphanHandle<TrackExtraCollection> theRecoTrackExtras = iEvent.put(selTrackExtras); 
00145       for ( unsigned index = 0; index<nTracks; ++index ) { 
00146         const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
00147         (selTracks->at(index)).setExtra(theTrackExtraRef);
00148       }
00149   } // END IF COPY EXTRAS
00150 
00151   //TRACKS AND TRAJECTORIES
00152   OrphanHandle<TrackCollection> theRecoTracks = iEvent.put(selTracks);
00153   OrphanHandle<vector<Trajectory> > theRecoTrajectories = iEvent.put(outputTJ);
00154 
00155   //TRACKS<->TRAJECTORIES MAP 
00156   nTracks = theRecoTracks->size();
00157   for ( unsigned index = 0; index<nTracks; ++index ) { 
00158     Ref<vector<Trajectory> > trajRef( theRecoTrajectories, index );
00159     Ref<TrackCollection>    tkRef( theRecoTracks, index );
00160     trajTrackMap->insert(trajRef,tkRef);
00161   }
00162   //MAP IN THE EVENT
00163   iEvent.put( trajTrackMap );
00164 }
00165 
00166 // ------------ method called once each job just before starting event loop  ------------
00167 void 
00168 QualityFilter::beginRun(edm::Run & run, const edm::EventSetup&)
00169 {
00170 }
00171 
00172 // ------------ method called once each job just after ending the event loop  ------------
00173 void 
00174 QualityFilter::endJob() {
00175 }
00176 
00177