CMS 3D CMS Logo

HLTDoublet.cc
Go to the documentation of this file.
1 
11 #include "HLTDoublet.h"
12 
14 
17 
19 
21 
23 #include <cmath>
24 
25 //
26 // constructors and destructor
27 //
28 template <typename T1, typename T2>
30  : HLTFilter(iConfig),
31  originTag1_(iConfig.template getParameter<std::vector<edm::InputTag>>("originTag1")),
32  originTag2_(iConfig.template getParameter<std::vector<edm::InputTag>>("originTag2")),
33  inputTag1_(iConfig.template getParameter<edm::InputTag>("inputTag1")),
34  inputTag2_(iConfig.template getParameter<edm::InputTag>("inputTag2")),
35  inputToken1_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag1_)),
36  inputToken2_(consumes<trigger::TriggerFilterObjectWithRefs>(inputTag2_)),
37  triggerType1_(iConfig.template getParameter<int>("triggerType1")),
38  triggerType2_(iConfig.template getParameter<int>("triggerType2")),
39  min_Dphi_(iConfig.template getParameter<double>("MinDphi")),
40  max_Dphi_(iConfig.template getParameter<double>("MaxDphi")),
41  min_Deta_(iConfig.template getParameter<double>("MinDeta")),
42  max_Deta_(iConfig.template getParameter<double>("MaxDeta")),
43  min_Minv_(iConfig.template getParameter<double>("MinMinv")),
44  max_Minv_(iConfig.template getParameter<double>("MaxMinv")),
45  min_DelR_(iConfig.template getParameter<double>("MinDelR")),
46  max_DelR_(iConfig.template getParameter<double>("MaxDelR")),
47  min_Pt_(iConfig.template getParameter<double>("MinPt")),
48  max_Pt_(iConfig.template getParameter<double>("MaxPt")),
49  min_N_(iConfig.template getParameter<int>("MinN")),
50  same_(inputTag1_.encode() == inputTag2_.encode()), // same collections to be compared?
51  cutdphi_(min_Dphi_ <= max_Dphi_), // cut active?
52  cutdeta_(min_Deta_ <= max_Deta_), // cut active?
53  cutminv_(min_Minv_ <= max_Minv_), // cut active?
54  cutdelr_(min_DelR_ <= max_DelR_), // cut active?
55  cutpt_(min_Pt_ <= max_Pt_) // cut active?
56 {
57  LogDebug("") << "InputTags and cuts : " << inputTag1_.encode() << " " << inputTag2_.encode() << triggerType1_ << " "
58  << triggerType2_ << " Dphi [" << min_Dphi_ << " " << max_Dphi_ << "]"
59  << " Deta [" << min_Deta_ << " " << max_Deta_ << "]"
60  << " Minv [" << min_Minv_ << " " << max_Minv_ << "]"
61  << " DelR [" << min_DelR_ << " " << max_DelR_ << "]"
62  << " Pt [" << min_Pt_ << " " << max_Pt_ << "]"
63  << " MinN =" << min_N_ << " same/dphi/deta/minv/delr/pt " << same_ << cutdphi_ << cutdeta_ << cutminv_
64  << cutdelr_ << cutpt_;
65 }
66 
67 template <typename T1, typename T2>
69 template <typename T1, typename T2>
73  std::vector<edm::InputTag> originTag1(1, edm::InputTag("hltOriginal1"));
74  std::vector<edm::InputTag> originTag2(1, edm::InputTag("hltOriginal2"));
75  desc.add<std::vector<edm::InputTag>>("originTag1", originTag1);
76  desc.add<std::vector<edm::InputTag>>("originTag2", originTag2);
77  desc.add<edm::InputTag>("inputTag1", edm::InputTag("hltFiltered1"));
78  desc.add<edm::InputTag>("inputTag2", edm::InputTag("hltFiltered22"));
79  desc.add<int>("triggerType1", 0);
80  desc.add<int>("triggerType2", 0);
81  desc.add<double>("MinDphi", +1.0);
82  desc.add<double>("MaxDphi", -1.0);
83  desc.add<double>("MinDeta", +1.0);
84  desc.add<double>("MaxDeta", -1.0);
85  desc.add<double>("MinMinv", +1.0);
86  desc.add<double>("MaxMinv", -1.0);
87  desc.add<double>("MinDelR", +1.0);
88  desc.add<double>("MaxDelR", -1.0);
89  desc.add<double>("MinPt", +1.0);
90  desc.add<double>("MaxPt", -1.0);
91  desc.add<int>("MinN", 1);
92  descriptions.add(defaultModuleLabel<HLTDoublet<T1, T2>>(), desc);
93 }
94 
95 //
96 // member functions
97 //
98 
99 // ------------ method called to produce the data ------------
100 template <typename T1, typename T2>
102  const edm::EventSetup& iSetup,
103  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
104  using namespace std;
105  using namespace edm;
106  using namespace reco;
107  using namespace trigger;
108 
109  // All HLT filters must create and fill an HLT filter object,
110  // recording any reconstructed physics objects satisfying (or not)
111  // this HLT filter, and place it in the Event.
112 
113  bool accept(false);
114 
115  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 0 " << std::endl;
116 
117  std::vector<T1Ref> coll1;
118  std::vector<T2Ref> coll2;
119 
120  // get hold of pre-filtered object collections
121  Handle<TriggerFilterObjectWithRefs> handle1, handle2;
122  if (iEvent.getByToken(inputToken1_, handle1) && iEvent.getByToken(inputToken2_, handle2)) {
123  handle1->getObjects(triggerType1_, coll1);
124  handle2->getObjects(triggerType2_, coll2);
125  const size_type n1(coll1.size());
126  const size_type n2(coll2.size());
127 
128  if (saveTags()) {
129  InputTag tagOld;
130  for (unsigned int i = 0; i < originTag1_.size(); ++i) {
131  filterproduct.addCollectionTag(originTag1_[i]);
132  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 1a/" << i << " " << originTag1_[i].encode()
133  << std::endl;
134  }
135  tagOld = InputTag();
136  for (size_type i1 = 0; i1 != n1; ++i1) {
137  const ProductID pid(coll1[i1].id());
138  const string& label(iEvent.getProvenance(pid).moduleLabel());
139  const string& instance(iEvent.getProvenance(pid).productInstanceName());
140  const string& process(iEvent.getProvenance(pid).processName());
142  if (tagOld.encode() != tagNew.encode()) {
143  filterproduct.addCollectionTag(tagNew);
144  tagOld = tagNew;
145  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 1b " << tagNew.encode() << std::endl;
146  }
147  }
148  for (unsigned int i = 0; i < originTag2_.size(); ++i) {
149  filterproduct.addCollectionTag(originTag2_[i]);
150  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 2a/" << i << " " << originTag2_[i].encode()
151  << std::endl;
152  }
153  tagOld = InputTag();
154  for (size_type i2 = 0; i2 != n2; ++i2) {
155  const ProductID pid(coll2[i2].id());
156  const string& label(iEvent.getProvenance(pid).moduleLabel());
157  const string& instance(iEvent.getProvenance(pid).productInstanceName());
158  const string& process(iEvent.getProvenance(pid).processName());
160  if (tagOld.encode() != tagNew.encode()) {
161  filterproduct.addCollectionTag(tagNew);
162  tagOld = tagNew;
163  LogVerbatim("HLTDoublet") << " XXX " << moduleLabel() << " 2b " << tagNew.encode() << std::endl;
164  }
165  }
166  }
167 
168  int n(0);
169  T1Ref r1;
170  T2Ref r2;
172  for (unsigned int i1 = 0; i1 != n1; i1++) {
173  r1 = coll1[i1];
174  p1 = r1->p4();
175  unsigned int I(0);
176  if (same_) {
177  I = i1 + 1;
178  }
179  for (unsigned int i2 = I; i2 != n2; i2++) {
180  r2 = coll2[i2];
181  p2 = r2->p4();
182 
183  double Dphi(std::abs(p1.phi() - p2.phi()));
184  if (Dphi > M_PI)
185  Dphi = 2.0 * M_PI - Dphi;
186 
187  double Deta(std::abs(p1.eta() - p2.eta()));
188 
189  p = p1 + p2;
190  double Minv(std::abs(p.mass()));
191  double Pt(p.pt());
192 
193  double DelR(sqrt(Dphi * Dphi + Deta * Deta));
194 
195  if (((!cutdphi_) || ((min_Dphi_ <= Dphi) && (Dphi <= max_Dphi_))) &&
196  ((!cutdeta_) || ((min_Deta_ <= Deta) && (Deta <= max_Deta_))) &&
197  ((!cutminv_) || ((min_Minv_ <= Minv) && (Minv <= max_Minv_))) &&
198  ((!cutdelr_) || ((min_DelR_ <= DelR) && (DelR <= max_DelR_))) &&
199  ((!cutpt_) || ((min_Pt_ <= Pt) && (Pt <= max_Pt_)))) {
200  n++;
201  filterproduct.addObject(triggerType1_, r1);
202  filterproduct.addObject(triggerType2_, r2);
203  }
204  }
205  }
206  // filter decision
207  accept = (n >= min_N_);
208  }
209 
210  return accept;
211 }
#define LogDebug(id)
const edm::InputTag inputTag2_
Definition: HLTDoublet.h:48
const double min_Pt_
Definition: HLTDoublet.h:57
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
const double min_Deta_
Definition: HLTDoublet.h:54
HLTDoublet(const edm::ParameterSet &)
Definition: HLTDoublet.cc:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HLTDoublet.cc:70
const bool same_
Definition: HLTDoublet.h:61
const bool cutpt_
Definition: HLTDoublet.h:62
static PFTauRenderPlugin instance
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const double min_Dphi_
Definition: HLTDoublet.h:53
const edm::InputTag inputTag1_
Definition: HLTDoublet.h:47
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken1_
Definition: HLTDoublet.h:49
const double min_Minv_
Definition: HLTDoublet.h:55
std::string defaultModuleLabel()
const double max_Deta_
Definition: HLTDoublet.h:54
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:30
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<C>)
char const * label
int iEvent
Definition: GenABIO.cc:224
const int triggerType2_
Definition: HLTDoublet.h:52
const bool cutdphi_
Definition: HLTDoublet.h:62
T sqrt(T t)
Definition: SSEVec.h:19
const double min_DelR_
Definition: HLTDoublet.h:56
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const bool cutdelr_
Definition: HLTDoublet.h:62
const std::complex< double > I
Definition: I.h:8
const int min_N_
Definition: HLTDoublet.h:58
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const double max_Minv_
Definition: HLTDoublet.h:55
const std::string * moduleLabel() const
Definition: HLTFilter.cc:61
double p2[4]
Definition: TauolaWrapper.h:90
#define M_PI
const double max_DelR_
Definition: HLTDoublet.h:56
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:25
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_
Definition: HLTDoublet.h:50
const double max_Pt_
Definition: HLTDoublet.h:57
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~HLTDoublet() override
std::string const & moduleLabel() const
Definition: Provenance.h:55
def encode(args, files)
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:46
HLT enums.
double p1[4]
Definition: TauolaWrapper.h:89
const std::vector< edm::InputTag > originTag2_
Definition: HLTDoublet.h:46
const std::vector< edm::InputTag > originTag1_
Definition: HLTDoublet.h:45
const int triggerType1_
Definition: HLTDoublet.h:51
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
Definition: HLTDoublet.cc:101
const bool cutminv_
Definition: HLTDoublet.h:62
const double max_Dphi_
Definition: HLTDoublet.h:53
std::string const & productInstanceName() const
Definition: Provenance.h:58
Provenance getProvenance(BranchID const &theID) const
Definition: Event.cc:114
const bool cutdeta_
Definition: HLTDoublet.h:62
math::PtEtaPhiELorentzVectorF LorentzVector