CMS 3D CMS Logo

CandMergerCleanOthersByDR.cc
Go to the documentation of this file.
1 //****************************************************************
2 //
3 // A simple class to combine two candidate collections into a single
4 // collection of ptrs to the candidates
5 // note: it is a std::vector as the candidates are from different
6 // collections
7 //
8 // collection 1 is added fully while collection 2 is deltaR cross
9 // cleaned against collection 2
10 //
11 // usecase: getting a unified list of e/gammas + jets for jet core
12 // regional tracking
13 //
14 //****************************************************************
15 
18 
21 
24 
29 
31 public:
34 
35  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
36 
37 private:
38  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
39 
40 private:
43  const float maxDR2ToClean_;
44 };
45 
46 namespace {
47  double pow2(double val) { return val * val; }
48 } // namespace
49 
51  : coll1Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("coll1"))),
52  coll2Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("coll2"))),
53  maxDR2ToClean_(pow2(iConfig.getParameter<double>("maxDRToClean"))) {
54  produces<std::vector<edm::Ptr<reco::Candidate>>>();
55 }
56 
57 namespace {
58  template <typename T>
59  edm::Handle<T> getHandle(const edm::Event& iEvent, const edm::EDGetTokenT<T>& token) {
61  iEvent.getByToken(token, handle);
62  return handle;
63  }
64 
65  bool hasDRMatch(const reco::Candidate& cand,
66  const std::vector<std::pair<float, float>>& etaPhisToMatch,
67  const float maxDR2) {
68  const float candEta = cand.eta();
69  const float candPhi = cand.phi();
70  for (const auto& etaPhi : etaPhisToMatch) {
71  if (reco::deltaR2(candEta, candPhi, etaPhi.first, etaPhi.second) <= maxDR2) {
72  return true;
73  }
74  }
75  return false;
76  }
77 } // namespace
78 
79 // ------------ method called to produce the data ------------
82  const edm::EventSetup& iSetup) const {
83  auto outColl = std::make_unique<std::vector<edm::Ptr<reco::Candidate>>>();
84 
85  auto coll1Handle = getHandle(iEvent, coll1Token_);
86  auto coll2Handle = getHandle(iEvent, coll2Token_);
87 
88  std::vector<std::pair<float, float>> coll1EtaPhis;
89  for (size_t objNr = 0; objNr < coll1Handle->size(); objNr++) {
90  edm::Ptr<reco::Candidate> objPtr(coll1Handle, objNr);
91  coll1EtaPhis.push_back({objPtr->eta(), objPtr->phi()}); //just to speed up the DR match
92  outColl->push_back(objPtr);
93  }
94  for (size_t objNr = 0; objNr < coll2Handle->size(); objNr++) {
95  edm::Ptr<reco::Candidate> objPtr(coll2Handle, objNr);
96  if (!hasDRMatch(*objPtr, coll1EtaPhis, maxDR2ToClean_)) {
97  outColl->push_back(objPtr);
98  }
99  }
100 
101  iEvent.put(std::move(outColl));
102 }
103 
104 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
107  desc.add<edm::InputTag>("coll1", edm::InputTag("egammasForCoreTracking"));
108  desc.add<edm::InputTag>("coll2", edm::InputTag("jetsForCoreTracking"));
109  desc.add<double>("maxDRToClean", 0.05);
110  descriptions.add("candMergerCleanOthersByDR", desc);
111 }
112 
113 //define this as a plug-in
const edm::EDGetTokenT< edm::View< reco::Candidate > > coll2Token_
const edm::EDGetTokenT< edm::View< reco::Candidate > > coll1Token_
constexpr int pow2(int x)
Definition: TauNNIdHW.h:51
int iEvent
Definition: GenABIO.cc:224
CandMergerCleanOthersByDR(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
def move(src, dest)
Definition: eostools.py:511
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity