TrackingParticleConversionRefSelector

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61
#include "FWCore/Framework/interface/global/EDProducer.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

#include "DataFormats/Common/interface/EDProductfwd.h"
#include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
#include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"

class TrackingParticleConversionRefSelector : public edm::global::EDProducer<> {
public:
  TrackingParticleConversionRefSelector(const edm::ParameterSet& iConfig);

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;

private:
  edm::EDGetTokenT<TrackingParticleCollection> tpToken_;
};

TrackingParticleConversionRefSelector::TrackingParticleConversionRefSelector(const edm::ParameterSet& iConfig)
    : tpToken_(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("src"))) {
  produces<TrackingParticleRefVector>();
}

void TrackingParticleConversionRefSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("src", edm::InputTag("mix", "MergedTrackTruth"));
  descriptions.add("trackingParticleConversionRefSelectorDefault", desc);
}

void TrackingParticleConversionRefSelector::produce(edm::StreamID,
                                                    edm::Event& iEvent,
                                                    const edm::EventSetup& iSetup) const {
  edm::Handle<TrackingParticleCollection> h_tps;
  iEvent.getByToken(tpToken_, h_tps);

  auto ret = std::make_unique<TrackingParticleRefVector>();

  // Logic is similar to Validation/RecoEgamma/plugins/PhotonValidator.cc
  // and RecoEgamma/EgammaMCTools/src/PhotonMCTruthFinder.cc,
  // but implemented purely in terms of TrackingParticles (simpler and works for pileup too)
  for (const auto& tp : *h_tps) {
    if (tp.pdgId() == 22) {
      for (const auto& vertRef : tp.decayVertices()) {
        for (const auto& tpRef : vertRef->daughterTracks()) {
          if (std::abs(tpRef->pdgId()) == 11) {
            ret->push_back(tpRef);
          }
        }
      }
    }
  }

  iEvent.put(std::move(ret));
}

DEFINE_FWK_MODULE(TrackingParticleConversionRefSelector);