CMS 3D CMS Logo

CandPtrProjector.cc
Go to the documentation of this file.
11 
12 
14 public:
15 
16  explicit CandPtrProjector(edm::ParameterSet const& iConfig);
17  void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
18  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
19 
20 private:
24 };
25 
27  candSrcToken_{consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("src"))},
28  vetoSrcToken_{consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("veto"))},
29  useDeltaRforFootprint_(iConfig.getParameter<bool>("useDeltaRforFootprint"))
30 {
31  produces<edm::PtrVector<reco::Candidate>>();
32 }
33 
34 void
36 {
37  using namespace edm;
39  iEvent.getByToken(vetoSrcToken_, vetoes);
40 
41  auto result = std::make_unique<PtrVector<reco::Candidate>>();
42  std::set<reco::CandidatePtr> vetoedPtrs;
43  for (auto const& veto : *vetoes) {
44  auto const n = veto.numberOfSourceCandidatePtrs();
45  for (size_t j {}; j<n; ++j) {
46  vetoedPtrs.insert(veto.sourceCandidatePtr(j));
47  }
48  }
49 
51  iEvent.getByToken(candSrcToken_, cands);
52  for (size_t i {}; i<cands->size(); ++i) {
53  auto const c = cands->ptrAt(i);
54  if (vetoedPtrs.find(c)==vetoedPtrs.cend()) {
55  bool addcand = true;
57  for (const auto& it : vetoedPtrs)
58  if (it.isNonnull() && it.isAvailable() && reco::deltaR2(*it, *c) < 0.00000025) {
59  addcand = false;
60  break;
61  }
62  if (addcand)
63  result->push_back(c);
64  }
65  }
66  iEvent.put(std::move(result));
67 }
68 
71  desc.add<edm::InputTag>("src");
72  desc.add<edm::InputTag>("veto");
73  desc.add<bool>("useDeltaRforFootprint", false);
74  descriptions.addWithDefaultLabel(desc);
75 }
76 
T getParameter(std::string const &) const
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< edm::View< reco::Candidate > > vetoSrcToken_
CandPtrProjector(edm::ParameterSet const &iConfig)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< edm::View< reco::Candidate > > candSrcToken_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
HLT enums.
def move(src, dest)
Definition: eostools.py:511