Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:00

0001 /** \class MCTrackMatcher
0002  *
0003  * \author Luca Lista, INFN
0004  *
0005  *
0006  */
0007 #include "DataFormats/Common/interface/Association.h"
0008 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0009 #include "FWCore/Framework/interface/global/EDProducer.h"
0010 #include "FWCore/Utilities/interface/InputTag.h"
0011 #include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociator.h"
0012 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0013 
0014 namespace edm {
0015   class ParameterSet;
0016 }
0017 
0018 using namespace edm;
0019 using namespace std;
0020 using namespace reco;
0021 
0022 class MCTrackMatcher : public edm::global::EDProducer<> {
0023 public:
0024   /// constructor
0025   MCTrackMatcher(const edm::ParameterSet &);
0026 
0027 private:
0028   void produce(edm::StreamID, edm::Event &evt, const edm::EventSetup &es) const override;
0029   edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator> associator_;
0030   edm::EDGetTokenT<edm::View<reco::Track>> tracks_;
0031   edm::EDGetTokenT<GenParticleCollection> genParticles_;
0032   edm::EDGetTokenT<std::vector<int>> genParticleInts_;
0033   edm::EDGetTokenT<TrackingParticleCollection> trackingParticles_;
0034   bool throwOnMissingTPCollection_;
0035   typedef edm::Association<reco::GenParticleCollection> GenParticleMatch;
0036 };
0037 
0038 #include "DataFormats/Common/interface/Handle.h"
0039 #include "FWCore/Framework/interface/Event.h"
0040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0041 #include "FWCore/Utilities/interface/EDMException.h"
0042 
0043 MCTrackMatcher::MCTrackMatcher(const ParameterSet &p)
0044     : associator_(consumes<reco::TrackToTrackingParticleAssociator>(p.getParameter<string>("associator"))),
0045       tracks_(consumes<edm::View<reco::Track>>(p.getParameter<InputTag>("tracks"))),
0046       genParticles_(consumes<GenParticleCollection>(p.getParameter<InputTag>("genParticles"))),
0047       genParticleInts_(consumes<std::vector<int>>(p.getParameter<InputTag>("genParticles"))),
0048       trackingParticles_(consumes<TrackingParticleCollection>(p.getParameter<InputTag>("trackingParticles"))),
0049       throwOnMissingTPCollection_(p.getParameter<bool>("throwOnMissingTPCollection")) {
0050   produces<GenParticleMatch>();
0051 }
0052 
0053 void MCTrackMatcher::produce(edm::StreamID, Event &evt, const EventSetup &es) const {
0054   Handle<reco::TrackToTrackingParticleAssociator> assoc;
0055   evt.getByToken(associator_, assoc);
0056   const reco::TrackToTrackingParticleAssociator *associator = assoc.product();
0057   Handle<View<Track>> tracks;
0058   evt.getByToken(tracks_, tracks);
0059   Handle<TrackingParticleCollection> trackingParticles;
0060   auto trackingParticlesFound = evt.getByToken(trackingParticles_, trackingParticles);
0061   Handle<vector<int>> barCodes;
0062   evt.getByToken(genParticleInts_, barCodes);
0063   Handle<GenParticleCollection> genParticles;
0064   evt.getByToken(genParticles_, genParticles);
0065   unique_ptr<GenParticleMatch> match(new GenParticleMatch(GenParticleRefProd(genParticles)));
0066   if (not throwOnMissingTPCollection_ and not trackingParticlesFound) {
0067     evt.put(std::move(match));
0068     return;
0069   }
0070   RecoToSimCollection associations = associator->associateRecoToSim(tracks, trackingParticles);
0071   GenParticleMatch::Filler filler(*match);
0072   size_t n = tracks->size();
0073   vector<int> indices(n, -1);
0074   for (size_t i = 0; i < n; ++i) {
0075     RefToBase<Track> track(tracks, i);
0076     RecoToSimCollection::const_iterator f = associations.find(track);
0077     if (f != associations.end()) {
0078       TrackingParticleRef tp = f->val.front().first;
0079       TrackingParticle::genp_iterator j, b = tp->genParticle_begin(), e = tp->genParticle_end();
0080       for (j = b; j != e; ++j) {
0081         const reco::GenParticle *p = j->get();
0082         if (p->status() == 1) {
0083           indices[i] = j->key();
0084           break;
0085         }
0086       }
0087     }
0088   }
0089   filler.insert(tracks, indices.begin(), indices.end());
0090   filler.fill();
0091   evt.put(std::move(match));
0092 }
0093 
0094 #include "FWCore/Framework/interface/MakerMacros.h"
0095 
0096 DEFINE_FWK_MODULE(MCTrackMatcher);