CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/SimTracker/TrackHistory/plugins/GenTrackMatcher.cc

Go to the documentation of this file.
00001 
00011 #include "FWCore/Framework/interface/EDProducer.h"
00012 #include "DataFormats/Common/interface/Association.h"
00013 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00014 #include "SimTracker/TrackHistory/interface/TrackHistory.h"
00015 
00016 namespace edm
00017 {
00018 class ParameterSet;
00019 }
00020 
00021 class GenTrackMatcher : public edm::EDProducer
00022 {
00023 public:
00025     GenTrackMatcher( const edm::ParameterSet & );
00026 
00027 private:
00028     void produce( edm::Event& evt, const edm::EventSetup& es );
00029     TrackHistory tracer_;
00030     edm::InputTag tracks_, genParticles_;
00031     typedef edm::Association<reco::GenParticleCollection> GenParticleMatch;
00032 };
00033 
00034 #include "DataFormats/Common/interface/Handle.h"
00035 #include "FWCore/Framework/interface/ESHandle.h"
00036 #include "FWCore/Framework/interface/Event.h"
00037 #include "FWCore/Utilities/interface/EDMException.h"
00038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00039 
00040 using namespace edm;
00041 using namespace std;
00042 using namespace reco;
00043 
00044 GenTrackMatcher::GenTrackMatcher(const ParameterSet & p) :
00045         tracer_(p),
00046         tracks_(p.getUntrackedParameter<edm::InputTag>("trackProducer")),
00047         genParticles_(p.getUntrackedParameter<edm::InputTag>("genParticles"))
00048 {
00049     produces<GenParticleMatch>();
00050 }
00051 
00052 void GenTrackMatcher::produce(Event& evt, const EventSetup& es)
00053 {
00054     Handle<View<Track> > tracks;
00055     evt.getByLabel(tracks_, tracks);
00056     Handle<vector<int> > barCodes;
00057     evt.getByLabel(genParticles_, barCodes);
00058     Handle<GenParticleCollection> genParticles;
00059     evt.getByLabel(genParticles_, genParticles);
00060     auto_ptr<GenParticleMatch> match(new GenParticleMatch(GenParticleRefProd(genParticles)));
00061     GenParticleMatch::Filler filler(*match);
00062     size_t n = tracks->size();
00063     vector<int> indices(n,-1);
00064     tracer_.newEvent(evt, es);
00065     for (size_t i = 0; i < n; ++ i )
00066     {
00067         RefToBase<Track> track(tracks, i);
00068         if (tracer_.evaluate(track))
00069         {
00070             const HepMC::GenParticle * particle = tracer_.genParticle();
00071             if (particle)
00072             {
00073                 int barCode = particle->barcode();
00074                 vector<int>::const_iterator b = barCodes->begin(), e = barCodes->end(), f = find( b, e, barCode );
00075                 if (f == e) throw edm::Exception(errors::InvalidReference)
00076                     << "found matching particle with barcode" << *f
00077                     << " which has not been found in " << genParticles_;
00078                 indices[i] = *f;
00079             }
00080         }
00081     }
00082     filler.insert(tracks, indices.begin(), indices.end());
00083     filler.fill();
00084     evt.put(match);
00085 }
00086 
00087 #include "FWCore/Framework/interface/MakerMacros.h"
00088 
00089 DEFINE_FWK_MODULE( GenTrackMatcher );
00090