CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/SimTracker/TrackAssociation/plugins/TrackMCQuality.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    TrackMCQuality
00004 // Class:      TrackMCQuality
00005 // 
00013 //
00014 // Original Author:  Jean-Roch Vlimant
00015 //         Created:  Fri Mar 27 15:19:03 CET 2009
00016 // $Id: TrackMCQuality.cc,v 1.2 2010/09/09 15:42:23 sunanda Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027 
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030 
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "FWCore/Framework/interface/ESHandle.h"
00033 
00034 #include "SimTracker/TrackAssociation/interface/TrackAssociatorBase.h"
00035 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
00036 
00037 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00038 #include "DataFormats/TrackReco/interface/Track.h"
00039 
00040 //
00041 // class decleration
00042 //
00043 
00044 class TrackMCQuality : public edm::EDProducer {
00045    public:
00046       explicit TrackMCQuality(const edm::ParameterSet&);
00047       ~TrackMCQuality();
00048 
00049    private:
00050       virtual void beginJob() ;
00051       virtual void produce(edm::Event&, const edm::EventSetup&);
00052       virtual void endJob() ;
00053       
00054       // ----------member data ---------------------------
00055 
00056   edm::ESHandle<TrackAssociatorBase> theAssociator;
00057   edm::InputTag label_tr;
00058   edm::InputTag label_tp;
00059   std::string associator;
00060 };
00061 
00062 //
00063 // constants, enums and typedefs
00064 //
00065 
00066 
00067 //
00068 // static data member definitions
00069 //
00070 
00071 //
00072 // constructors and destructor
00073 //
00074 TrackMCQuality::TrackMCQuality(const edm::ParameterSet& pset):
00075   label_tr(pset.getParameter< edm::InputTag >("label_tr")),
00076   label_tp(pset.getParameter< edm::InputTag >("label_tp")),
00077   associator(pset.getParameter< std::string >("associator"))
00078 {
00079   
00080   produces<reco::TrackCollection>();
00081 }
00082 
00083 
00084 TrackMCQuality::~TrackMCQuality()
00085 {
00086 }
00087 
00088 
00089 //
00090 // member functions
00091 //
00092 
00093 // ------------ method called to produce the data  ------------
00094 void
00095 TrackMCQuality::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00096 {
00097 
00098   iSetup.get<TrackAssociatorRecord>().get(associator,theAssociator);
00099 
00100 
00101    using namespace edm;
00102    Handle<TrackingParticleCollection>  TPCollection ;
00103    iEvent.getByLabel(label_tp, TPCollection);
00104      
00105    Handle<edm::View<reco::Track> > trackCollection;
00106    iEvent.getByLabel (label_tr, trackCollection );
00107 
00108    reco::RecoToSimCollection recSimColl=theAssociator->associateRecoToSim(trackCollection,
00109                                                                           TPCollection,
00110                                                                           &iEvent);
00111    
00112    //then loop the track collection
00113    std::auto_ptr<reco::TrackCollection> outTracks(new reco::TrackCollection(trackCollection->size()));
00114    
00115    for (unsigned int iT=0;iT!=trackCollection->size();++iT){
00116      edm::RefToBase<reco::Track> track(trackCollection, iT);
00117      bool matched=false;
00118      //find it in the map
00119      if (recSimColl.find(track)!=recSimColl.end()){
00120        // you can get the data if you want
00121        std::vector<std::pair<TrackingParticleRef, double> > tp= recSimColl[track];
00122        matched=true;
00123      }
00124      else{
00125        matched=false;
00126      }     
00127 
00128      //copy the track into the new container
00129      (*outTracks)[iT] = reco::Track(*track);
00130      if (matched){
00131        (*outTracks)[iT].setQuality(reco::TrackBase::qualitySize); //is not assigned to any quality. use it as a fake/matched flag
00132      }
00133    }
00134    
00135    iEvent.put(outTracks);
00136 }
00137 
00138 // ------------ method called once each job just before starting event loop  ------------
00139 void 
00140 TrackMCQuality::beginJob()
00141 {
00142 }
00143 
00144 // ------------ method called once each job just after ending the event loop  ------------
00145 void 
00146 TrackMCQuality::endJob() {
00147 }
00148 
00149 //define this as a plug-in
00150 DEFINE_FWK_MODULE(TrackMCQuality);