CMS 3D CMS Logo

CandPtrProjector.cc
Go to the documentation of this file.
11 
13 public:
14  explicit CandPtrProjector(edm::ParameterSet const& iConfig);
15  void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
16  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
17 
18 private:
23 };
24 
26  : candSrcToken_{consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("src"))},
27  vetoSrcToken_{consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("veto"))},
28  useDeltaRforFootprint_(iConfig.getParameter<bool>("useDeltaRforFootprint")),
29  extendVetoBySingleSourcePtr_(iConfig.getParameter<bool>("extendVetoBySingleSourcePtr")) {
30  produces<edm::PtrVector<reco::Candidate>>();
31 }
32 
34  using namespace edm;
36  iEvent.getByToken(vetoSrcToken_, vetoes);
37 
38  auto result = std::make_unique<PtrVector<reco::Candidate>>();
39  std::set<reco::CandidatePtr> vetoedPtrs;
40  for (auto const& veto : *vetoes) {
41  auto const n = veto.numberOfSourceCandidatePtrs();
42  for (size_t j{}; j < n; ++j) {
43  vetoedPtrs.insert(veto.sourceCandidatePtr(j));
44  }
45  }
46 
48  iEvent.getByToken(candSrcToken_, cands);
49  for (size_t i{}; i < cands->size(); ++i) {
50  auto const c = cands->ptrAt(i);
51  if (vetoedPtrs.find(c) == vetoedPtrs.cend()) {
52  bool addcand = true;
53  if ((extendVetoBySingleSourcePtr_) && (c->numberOfSourceCandidatePtrs() == 1))
54  addcand = (vetoedPtrs.find(c->sourceCandidatePtr(0)) == vetoedPtrs.cend());
56  for (const auto& it : vetoedPtrs)
57  if (it.isNonnull() && it.isAvailable() && reco::deltaR2(*it, *c) < 0.00000025) {
58  addcand = false;
59  break;
60  }
61  if (addcand)
62  result->push_back(c);
63  }
64  }
65  iEvent.put(std::move(result));
66 }
67 
70  desc.add<edm::InputTag>("src");
71  desc.add<edm::InputTag>("veto");
72  desc.add<bool>("useDeltaRforFootprint", false);
73  desc.add<bool>("extendVetoBySingleSourcePtr", true);
74  descriptions.addWithDefaultLabel(desc);
75 }
76 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< edm::View< reco::Candidate > > vetoSrcToken_
CandPtrProjector(edm::ParameterSet const &iConfig)
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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.
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
def move(src, dest)
Definition: eostools.py:511