CMS 3D CMS Logo

HLT2L1TkMuonL1TkMuonMuRefDR.cc
Go to the documentation of this file.
1 #include <cmath>
2 
12 
14 
15 //
16 // constructors and destructor
17 //
19  : HLTFilter(iConfig),
20  originTag1_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag1")),
21  originTag2_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag2")),
22  inputTag1_(iConfig.getParameter<edm::InputTag>("inputTag1")),
23  inputTag2_(iConfig.getParameter<edm::InputTag>("inputTag2")),
24  inputToken1_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag1_)),
25  inputToken2_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag2_)),
26  minDR_(iConfig.getParameter<double>("MinDR")),
27  min_N_(iConfig.getParameter<int>("MinN")),
28  same_(inputTag1_.encode() == inputTag2_.encode()) {} // same collections to be compared?
29 
31 
35  std::vector<edm::InputTag> originTag1(1, edm::InputTag("hltOriginal1"));
36  std::vector<edm::InputTag> originTag2(1, edm::InputTag("hltOriginal2"));
37  desc.add<std::vector<edm::InputTag>>("originTag1", originTag1);
38  desc.add<std::vector<edm::InputTag>>("originTag2", originTag2);
39  desc.add<edm::InputTag>("inputTag1", edm::InputTag("hltFiltered1"));
40  desc.add<edm::InputTag>("inputTag2", edm::InputTag("hltFiltered2"));
41  desc.add<double>("MinDR", -1.0);
42  desc.add<int>("MinN", 1);
43 
44  descriptions.add("hlt2L1TkMuonL1TkMuonMuRefDR", desc);
45 }
46 
48  std::vector<l1t::TrackerMuonRef>& coll1,
49  std::vector<l1t::TrackerMuonRef>& coll2,
50  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
52  if (iEvent.getByToken(inputToken1_, handle1) and iEvent.getByToken(inputToken2_, handle2)) {
53  // get hold of pre-filtered object collections
56  const trigger::size_type n1(coll1.size());
57  const trigger::size_type n2(coll2.size());
58 
59  if (saveTags()) {
60  edm::InputTag tagOld;
61  for (unsigned int i = 0; i < originTag1_.size(); ++i) {
62  filterproduct.addCollectionTag(originTag1_[i]);
63  }
64  tagOld = edm::InputTag();
65  for (trigger::size_type i1 = 0; i1 != n1; ++i1) {
66  const edm::ProductID pid(coll1[i1].id());
67  const std::string& label(iEvent.getProvenance(pid).moduleLabel());
68  const std::string& instance(iEvent.getProvenance(pid).productInstanceName());
69  const std::string& process(iEvent.getProvenance(pid).processName());
71  if (tagOld.encode() != tagNew.encode()) {
72  filterproduct.addCollectionTag(tagNew);
73  tagOld = tagNew;
74  }
75  }
76  for (unsigned int i = 0; i < originTag2_.size(); ++i) {
77  filterproduct.addCollectionTag(originTag2_[i]);
78  }
79  tagOld = edm::InputTag();
80  for (trigger::size_type i2 = 0; i2 != n2; ++i2) {
81  const edm::ProductID pid(coll2[i2].id());
82  const std::string& label(iEvent.getProvenance(pid).moduleLabel());
83  const std::string& instance(iEvent.getProvenance(pid).productInstanceName());
84  const std::string& process(iEvent.getProvenance(pid).processName());
86  if (tagOld.encode() != tagNew.encode()) {
87  filterproduct.addCollectionTag(tagNew);
88  tagOld = tagNew;
89  }
90  }
91  }
92 
93  return true;
94  } else
95  return false;
96 }
97 
100  l1t::TrackerMuonRef& r2) const {
101  if (minDR_ < 0.)
102  return true;
103 
104  return (reco::deltaR2(r1->phEta(), r1->phPhi(), r2->phEta(), r2->phPhi()) > minDR_ * minDR_);
105 }
106 
107 // ------------ method called to produce the data ------------
109  const edm::EventSetup& iSetup,
110  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
111  // All HLT filters must create and fill an HLT filter object,
112  // recording any reconstructed physics objects satisfying (or not)
113  // this HLT filter, and place it in the Event.
114  bool accept(false);
115 
116  std::vector<l1t::TrackerMuonRef> coll1;
117  std::vector<l1t::TrackerMuonRef> coll2;
118 
119  if (getCollections(iEvent, coll1, coll2, filterproduct)) {
120  int n(0);
123 
124  for (unsigned int i1 = 0; i1 != coll1.size(); i1++) {
125  r1 = coll1[i1];
126  unsigned int I(0);
127  if (same_) {
128  I = i1 + 1;
129  }
130  for (unsigned int i2 = I; i2 != coll2.size(); i2++) {
131  r2 = coll2[i2];
132 
133  if (!computeDR(iEvent, r1, r2))
134  continue;
135 
136  n++;
139  }
140  }
141 
142  accept = accept || (n >= min_N_);
143  }
144 
145  return accept;
146 }
147 
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static PFTauRenderPlugin instance
std::string encode() const
Definition: InputTag.cc:159
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken1_
const std::vector< edm::InputTag > originTag1_
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
uint16_t size_type
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
char const * label
const std::vector< edm::InputTag > originTag2_
int iEvent
Definition: GenABIO.cc:224
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
HLT2L1TkMuonL1TkMuonMuRefDR(const edm::ParameterSet &)
bool getCollections(edm::Event &iEvent, std::vector< l1t::TrackerMuonRef > &coll1, std::vector< l1t::TrackerMuonRef > &coll2, trigger::TriggerFilterObjectWithRefs &filterproduct) const
const std::complex< double > I
Definition: I.h:8
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool saveTags() const
Definition: HLTFilter.h:46
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:25
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def encode(args, files)
bool computeDR(edm::Event &iEvent, l1t::TrackerMuonRef &c1, l1t::TrackerMuonRef &c2) const
HLT enums.
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_