CMS 3D CMS Logo

CandPtrProjector.cc
Go to the documentation of this file.
11 
13 public:
14 
15  explicit CandPtrProjector(edm::ParameterSet const& iConfig);
16  void produce(edm::StreamID, edm::Event&, edm::EventSetup const&) const override;
17  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
18 
19 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 {
30  produces<edm::PtrVector<reco::Candidate>>();
31 }
32 
33 void
35 {
36  using namespace edm;
38  iEvent.getByToken(vetoSrcToken_, vetoes);
39 
40  auto result = std::make_unique<PtrVector<reco::Candidate>>();
41  std::set<reco::CandidatePtr> vetoedPtrs;
42  for (auto const& veto : *vetoes) {
43  auto const n = veto.numberOfSourceCandidatePtrs();
44  for (size_t j {}; j<n; ++j) {
45  vetoedPtrs.insert(veto.sourceCandidatePtr(j));
46  }
47  }
48 
50  iEvent.getByToken(candSrcToken_, cands);
51  for (size_t i {}; i<cands->size(); ++i) {
52  auto const c = cands->ptrAt(i);
53  if (vetoedPtrs.find(c)==vetoedPtrs.cend()) {
54  bool addcand = true;
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  descriptions.addWithDefaultLabel(desc);
74 }
75 
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:137
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::View< reco::Candidate > > vetoSrcToken_
CandPtrProjector(edm::ParameterSet const &iConfig)
int iEvent
Definition: GenABIO.cc:230
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:510