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 
16 
19 
22 
25 
30 
32 public:
35 
36  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
37 
38 private:
39  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&)const override;
40 
41 private:
44  const float maxDR2ToClean_;
45 
46 };
47 
48 
49 namespace {
50  double pow2(double val){return val*val;}
51 }
52 
54  coll1Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("coll1"))),
55  coll2Token_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("coll2"))),
56  maxDR2ToClean_(pow2(iConfig.getParameter<double>("maxDRToClean")))
57 {
58  produces<std::vector<edm::Ptr<reco::Candidate> > >();
59 }
60 
61 namespace {
62  template<typename T> edm::Handle<T> getHandle(const edm::Event& iEvent,
63  const edm::EDGetTokenT<T>& token)
64  {
66  iEvent.getByToken(token,handle);
67  return handle;
68  }
69 
70  bool hasDRMatch(const reco::Candidate& cand,const std::vector<std::pair<float,float> >& etaPhisToMatch,const float maxDR2)
71  {
72  const float candEta = cand.eta();
73  const float candPhi = cand.phi();
74  for(const auto& etaPhi : etaPhisToMatch){
75  if(reco::deltaR2(candEta,candPhi,etaPhi.first,etaPhi.second)<=maxDR2){
76  return true;
77  }
78  }
79  return false;
80  }
81 }
82 
83 // ------------ method called to produce the data ------------
84 void
86 {
87  auto outColl = std::make_unique<std::vector<edm::Ptr<reco::Candidate> > >();
88 
89  auto coll1Handle = getHandle(iEvent,coll1Token_);
90  auto coll2Handle = getHandle(iEvent,coll2Token_);
91 
92  std::vector<std::pair<float,float> > coll1EtaPhis;
93  for(size_t objNr=0;objNr<coll1Handle->size();objNr++){
94  edm::Ptr<reco::Candidate> objPtr(coll1Handle,objNr);
95  coll1EtaPhis.push_back({objPtr->eta(),objPtr->phi()}); //just to speed up the DR match
96  outColl->push_back(objPtr);
97  }
98  for(size_t objNr=0;objNr<coll2Handle->size();objNr++){
99  edm::Ptr<reco::Candidate> objPtr(coll2Handle,objNr);
100  if(!hasDRMatch(*objPtr,coll1EtaPhis,maxDR2ToClean_)){
101  outColl->push_back(objPtr);
102  }
103  }
104 
105  iEvent.put(std::move(outColl));
106 }
107 
108 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
109 void
111 
113  desc.add<edm::InputTag>("coll1", edm::InputTag("egammasForCoreTracking"));
114  desc.add<edm::InputTag>("coll2", edm::InputTag("jetsForCoreTracking"));
115  desc.add<double>("maxDRToClean", 0.05);
116  descriptions.add("candMergerCleanOthersByDR", desc);
117 }
118 
119 //define this as a plug-in
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const edm::EDGetTokenT< edm::View< reco::Candidate > > coll2Token_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const edm::EDGetTokenT< edm::View< reco::Candidate > > coll1Token_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CandMergerCleanOthersByDR(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
virtual double eta() const =0
momentum pseudorapidity
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
virtual double phi() const =0
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511