CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HLT2L1TkMuonL1TkMuonMuRefDR.cc
Go to the documentation of this file.
1 #include <cmath>
2 
14 
16 
17 //
18 // constructors and destructor
19 //
21  : HLTFilter(iConfig),
22  originTag1_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag1")),
23  originTag2_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag2")),
24  inputTag1_(iConfig.getParameter<edm::InputTag>("inputTag1")),
25  inputTag2_(iConfig.getParameter<edm::InputTag>("inputTag2")),
26  inputToken1_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag1_)),
27  inputToken2_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag2_)),
28  minDR_(iConfig.getParameter<double>("MinDR")),
29  min_N_(iConfig.getParameter<int>("MinN")),
30  same_(inputTag1_.encode() == inputTag2_.encode()) {} // same collections to be compared?
31 
33 
37  std::vector<edm::InputTag> originTag1(1, edm::InputTag("hltOriginal1"));
38  std::vector<edm::InputTag> originTag2(1, edm::InputTag("hltOriginal2"));
39  desc.add<std::vector<edm::InputTag>>("originTag1", originTag1);
40  desc.add<std::vector<edm::InputTag>>("originTag2", originTag2);
41  desc.add<edm::InputTag>("inputTag1", edm::InputTag("hltFiltered1"));
42  desc.add<edm::InputTag>("inputTag2", edm::InputTag("hltFiltered2"));
43  desc.add<double>("MinDR", -1.0);
44  desc.add<int>("MinN", 1);
45 
46  descriptions.add("hlt2L1TkMuonL1TkMuonMuRefDR", desc);
47 }
48 
50  std::vector<l1t::TkMuonRef>& coll1,
51  std::vector<l1t::TkMuonRef>& coll2,
52  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
54  if (iEvent.getByToken(inputToken1_, handle1) and iEvent.getByToken(inputToken2_, handle2)) {
55  // get hold of pre-filtered object collections
56  handle1->getObjects(trigger::TriggerObjectType::TriggerL1TkMu, coll1);
57  handle2->getObjects(trigger::TriggerObjectType::TriggerL1TkMu, coll2);
58  const trigger::size_type n1(coll1.size());
59  const trigger::size_type n2(coll2.size());
60 
61  if (saveTags()) {
62  edm::InputTag tagOld;
63  for (unsigned int i = 0; i < originTag1_.size(); ++i) {
64  filterproduct.addCollectionTag(originTag1_[i]);
65  }
66  tagOld = edm::InputTag();
67  for (trigger::size_type i1 = 0; i1 != n1; ++i1) {
68  const edm::ProductID pid(coll1[i1].id());
69  const std::string& label(iEvent.getProvenance(pid).moduleLabel());
71  const std::string& process(iEvent.getProvenance(pid).processName());
73  if (tagOld.encode() != tagNew.encode()) {
74  filterproduct.addCollectionTag(tagNew);
75  tagOld = tagNew;
76  }
77  }
78  for (unsigned int i = 0; i < originTag2_.size(); ++i) {
79  filterproduct.addCollectionTag(originTag2_[i]);
80  }
81  tagOld = edm::InputTag();
82  for (trigger::size_type i2 = 0; i2 != n2; ++i2) {
83  const edm::ProductID pid(coll2[i2].id());
84  const std::string& label(iEvent.getProvenance(pid).moduleLabel());
86  const std::string& process(iEvent.getProvenance(pid).processName());
88  if (tagOld.encode() != tagNew.encode()) {
89  filterproduct.addCollectionTag(tagNew);
90  tagOld = tagNew;
91  }
92  }
93  }
94 
95  return true;
96  } else
97  return false;
98 }
99 
100 std::pair<float, float> HLT2L1TkMuonL1TkMuonMuRefDR::convertEtaPhi(l1t::TkMuonRef& tkmu) const {
101  float muRefEta = 0.;
102  float muRefPhi = 0.;
103 
104  if (tkmu->muonDetector() != emtfRegion_) {
105  if (tkmu->muRef().isNull())
106  return std::make_pair(muRefEta, muRefPhi);
107 
108  muRefEta = tkmu->muRef()->hwEta() * etaScale_;
109  muRefPhi = static_cast<float>(l1t::MicroGMTConfiguration::calcGlobalPhi(
110  tkmu->muRef()->hwPhi(), tkmu->muRef()->trackFinderType(), tkmu->muRef()->processor()));
111  muRefPhi = muRefPhi * phiScale_;
112  } else {
113  if (tkmu->emtfTrk().isNull())
114  return std::make_pair(muRefEta, muRefPhi);
115 
116  muRefEta = tkmu->emtfTrk()->Eta();
117  muRefPhi = angle_units::operators::convertDegToRad(tkmu->emtfTrk()->Phi_glob());
118  }
119  muRefPhi = reco::reduceRange(muRefPhi);
120 
121  return std::make_pair(muRefEta, muRefPhi);
122 }
123 
125  if (minDR_ < 0.)
126  return true;
127 
128  auto [eta1, phi1] = convertEtaPhi(r1);
129  auto [eta2, phi2] = convertEtaPhi(r2);
130  return (reco::deltaR2(eta1, phi1, eta2, phi2) > minDR_ * minDR_);
131 }
132 
133 // ------------ method called to produce the data ------------
135  const edm::EventSetup& iSetup,
136  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
137  // All HLT filters must create and fill an HLT filter object,
138  // recording any reconstructed physics objects satisfying (or not)
139  // this HLT filter, and place it in the Event.
140  bool accept(false);
141 
142  std::vector<l1t::TkMuonRef> coll1;
143  std::vector<l1t::TkMuonRef> coll2;
144 
145  if (getCollections(iEvent, coll1, coll2, filterproduct)) {
146  int n(0);
149 
150  for (unsigned int i1 = 0; i1 != coll1.size(); i1++) {
151  r1 = coll1[i1];
152  unsigned int I(0);
153  if (same_) {
154  I = i1 + 1;
155  }
156  for (unsigned int i2 = I; i2 != coll2.size(); i2++) {
157  r2 = coll2[i2];
158 
159  if (!computeDR(iEvent, r1, r2))
160  continue;
161 
162  n++;
165  }
166  }
167 
168  accept = accept || (n >= min_N_);
169  }
170 
171  return accept;
172 }
173 
constexpr double convertDegToRad(NumType degrees)
Definition: angle_units.h:27
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
static PFTauRenderPlugin instance
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken1_
This has all to be decided for Phase-2. Here is Thiago&#39;s proposal.
std::pair< float, float > convertEtaPhi(l1t::TkMuonRef &tkmu) const
bool computeDR(edm::Event &iEvent, l1t::TkMuonRef &c1, l1t::TkMuonRef &c2) const
const std::vector< edm::InputTag > originTag1_
std::string const & processName() const
Definition: Provenance.h:57
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
uint16_t size_type
std::string encode() const
Definition: InputTag.cc:159
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref&lt;C&gt;)
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 &)
const std::complex< double > I
Definition: I.h:8
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isNull() const
Checks for null.
Definition: Ref.h:235
static constexpr unsigned int emtfRegion_
static int calcGlobalPhi(int locPhi, tftype t, int proc)
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)
std::string const & moduleLabel() const
Definition: Provenance.h:55
bool saveTags() const
Definition: HLTFilter.h:46
bool getCollections(edm::Event &iEvent, std::vector< l1t::TkMuonRef > &coll1, std::vector< l1t::TkMuonRef > &coll2, trigger::TriggerFilterObjectWithRefs &filterproduct) const
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_
std::string const & productInstanceName() const
Definition: Provenance.h:58
tuple process
Definition: LaserDQM_cfg.py:3
Provenance const & getProvenance(BranchID const &theID) const
Definition: Event.cc:118