CMS 3D CMS Logo

HLTTripletMass.cc
Go to the documentation of this file.
1 //
2 // HLTFilter selecting events with a minimum number of tri-object candidates passing invariant-mass cuts
3 //
4 
15 #include <string>
16 #include <vector>
17 
18 template <typename T1, typename T2, typename T3>
19 class HLTTripletMass : public HLTFilter {
23 
24 public:
25  explicit HLTTripletMass(const edm::ParameterSet&);
26  ~HLTTripletMass() override = default;
27  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
28  bool hltFilter(edm::Event&,
29  const edm::EventSetup&,
30  trigger::TriggerFilterObjectWithRefs& filterproduct) const override;
32  std::vector<T1Ref>& coll1,
33  std::vector<T2Ref>& coll2,
34  std::vector<T3Ref>& coll3,
35  trigger::TriggerFilterObjectWithRefs& filterproduct) const;
36 
37 private:
38  // configuration
39  const std::vector<edm::InputTag> originTag1_; // input tag identifying originals 1st product
40  const std::vector<edm::InputTag> originTag2_; // input tag identifying originals 2nd product
41  const std::vector<edm::InputTag> originTag3_; // input tag identifying originals 3rd product
45  const int triggerType1_;
46  const int triggerType2_;
47  const int triggerType3_;
48  const std::vector<double> min_InvMass_; // minimum invariant mass of pair
49  const std::vector<double> max_InvMass_; // maximum invariant mass of pair
50  const double max_DR_; // maximum deltaR between the p4 of 3rd product and p4 of (1+2) product
51  const double max_DR2_; // maximum deltaR^2 between the p4 of 3rd product and p4 of (1+2) product
52  const int min_N_;
53  const bool is1and2Same_;
54  const bool is2and3Same_;
55 };
56 
57 template <typename T1, typename T2, typename T3>
59  : HLTFilter(iConfig),
60  originTag1_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag1")),
61  originTag2_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag2")),
62  originTag3_(iConfig.getParameter<std::vector<edm::InputTag>>("originTag3")),
63  inputToken1_(consumes(iConfig.getParameter<edm::InputTag>("inputTag1"))),
64  inputToken2_(consumes(iConfig.getParameter<edm::InputTag>("inputTag2"))),
65  inputToken3_(consumes(iConfig.getParameter<edm::InputTag>("inputTag3"))),
66  triggerType1_(iConfig.getParameter<int>("triggerType1")),
67  triggerType2_(iConfig.getParameter<int>("triggerType2")),
68  triggerType3_(iConfig.getParameter<int>("triggerType3")),
69  min_InvMass_(iConfig.getParameter<vector<double>>("MinInvMass")),
70  max_InvMass_(iConfig.getParameter<vector<double>>("MaxInvMass")),
71  max_DR_(iConfig.getParameter<double>("MaxDR")),
72  max_DR2_(max_DR_ * max_DR_),
73  min_N_(iConfig.getParameter<int>("MinN")),
74  is1and2Same_(iConfig.getParameter<bool>("is1and2Same")),
75  is2and3Same_(iConfig.getParameter<bool>("is2and3Same")) {
76  if (min_InvMass_.size() != max_InvMass_.size()) {
77  throw cms::Exception("Configuration") << "size of \"MinInvMass\" (" << min_InvMass_.size()
78  << ") and \"MaxInvMass\" (" << max_InvMass_.size() << ") differ";
79  }
80  if (max_DR_ < 0) {
81  throw cms::Exception("Configuration") << "invalid value for parameter \"MaxDR\" (must be >= 0): " << max_DR_;
82  }
83 }
84 
85 template <typename T1, typename T2, typename T3>
88  makeHLTFilterDescription(desc);
89  desc.add<std::vector<edm::InputTag>>("originTag1", {edm::InputTag("hltOriginal1")});
90  desc.add<std::vector<edm::InputTag>>("originTag2", {edm::InputTag("hltOriginal2")});
91  desc.add<std::vector<edm::InputTag>>("originTag3", {edm::InputTag("hltOriginal3")});
92  desc.add<edm::InputTag>("inputTag1", edm::InputTag("hltFiltered1"));
93  desc.add<edm::InputTag>("inputTag2", edm::InputTag("hltFiltered2"));
94  desc.add<edm::InputTag>("inputTag3", edm::InputTag("hltFiltered3"));
95  desc.add<int>("triggerType1", 0);
96  desc.add<int>("triggerType2", 0);
97  desc.add<int>("triggerType3", 0);
98 
99  desc.add<vector<double>>("MinInvMass", {0});
100  desc.add<vector<double>>("MaxInvMass", {1e12});
101 
102  desc.add<double>("MaxDR", 1e4);
103  desc.add<int>("MinN", 0);
104 
105  desc.add<bool>("is1and2Same", false);
106  desc.add<bool>("is2and3Same", false);
107  descriptions.addWithDefaultLabel(desc);
108 }
109 
110 template <typename T1, typename T2, typename T3>
112  std::vector<T1Ref>& coll1,
113  std::vector<T2Ref>& coll2,
114  std::vector<T3Ref>& coll3,
115  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
116  edm::Handle<trigger::TriggerFilterObjectWithRefs> handle1, handle2, handle3;
117  if (iEvent.getByToken(inputToken1_, handle1) and iEvent.getByToken(inputToken2_, handle2) and
118  iEvent.getByToken(inputToken3_, handle3)) {
119  // get hold of pre-filtered object collections
120  handle1->getObjects(triggerType1_, coll1);
121  handle2->getObjects(triggerType2_, coll2);
122  handle3->getObjects(triggerType3_, coll3);
123  const trigger::size_type n1(coll1.size());
124  const trigger::size_type n2(coll2.size());
125  const trigger::size_type n3(coll3.size());
126 
127  if (saveTags()) {
128  edm::InputTag tagOld;
129  for (unsigned int i = 0; i < originTag1_.size(); ++i) {
130  filterproduct.addCollectionTag(originTag1_[i]);
131  }
132  tagOld = edm::InputTag();
133  for (trigger::size_type i1 = 0; i1 != n1; ++i1) {
134  const edm::ProductID pid(coll1[i1].id());
135  const auto& prov = iEvent.getStableProvenance(pid);
136  const std::string& label(prov.moduleLabel());
137  const std::string& instance(prov.productInstanceName());
138  const std::string& process(prov.processName());
140  if (tagOld.encode() != tagNew.encode()) {
141  filterproduct.addCollectionTag(tagNew);
142  tagOld = tagNew;
143  }
144  }
145  for (unsigned int i = 0; i < originTag2_.size(); ++i) {
146  filterproduct.addCollectionTag(originTag2_[i]);
147  }
148  tagOld = edm::InputTag();
149  for (trigger::size_type i2 = 0; i2 != n2; ++i2) {
150  const edm::ProductID pid(coll2[i2].id());
151  const auto& prov = iEvent.getStableProvenance(pid);
152  const std::string& label(prov.moduleLabel());
153  const std::string& instance(prov.productInstanceName());
154  const std::string& process(prov.processName());
156  if (tagOld.encode() != tagNew.encode()) {
157  filterproduct.addCollectionTag(tagNew);
158  tagOld = tagNew;
159  }
160  }
161  for (unsigned int i = 0; i < originTag3_.size(); ++i) {
162  filterproduct.addCollectionTag(originTag3_[i]);
163  }
164  tagOld = edm::InputTag();
165  for (trigger::size_type i3 = 0; i3 != n3; ++i3) {
166  const edm::ProductID pid(coll3[i3].id());
167  const auto& prov = iEvent.getStableProvenance(pid);
168  const std::string& label(prov.moduleLabel());
169  const std::string& instance(prov.productInstanceName());
170  const std::string& process(prov.processName());
172  if (tagOld.encode() != tagNew.encode()) {
173  filterproduct.addCollectionTag(tagNew);
174  tagOld = tagNew;
175  }
176  }
177  }
178 
179  return true;
180  } else
181  return false;
182 }
183 
184 // ------------ method called to produce the data ------------
185 template <typename T1, typename T2, typename T3>
187  const edm::EventSetup& iSetup,
188  trigger::TriggerFilterObjectWithRefs& filterproduct) const {
189  // All HLT filters must create and fill an HLT filter object,
190  // recording any reconstructed physics objects satisfying (or not)
191  // this HLT filter, and place it in the Event.
192 
193  std::vector<T1Ref> coll1;
194  std::vector<T2Ref> coll2;
195  std::vector<T3Ref> coll3;
196 
197  int n(0);
198  if (getCollections(iEvent, coll1, coll2, coll3, filterproduct)) {
199  T1Ref r1;
200  T2Ref r2;
201  T3Ref r3;
202 
203  reco::Particle::LorentzVector dauA_p4, dauB_p4, dauAB_p4, dauC_p4;
204  for (unsigned int i1 = 0; i1 != coll1.size(); i1++) {
205  r1 = coll1[i1];
206  dauA_p4 = reco::Particle::LorentzVector(r1->px(), r1->py(), r1->pz(), r1->energy());
207  unsigned int i2 = is1and2Same_ ? i1 + 1 : 0;
208  for (; i2 != coll2.size(); i2++) {
209  r2 = coll2[i2];
210  dauB_p4 = reco::Particle::LorentzVector(r2->px(), r2->py(), r2->pz(), r2->energy());
211  dauAB_p4 = dauA_p4 + dauB_p4;
212 
213  unsigned int i3 = is2and3Same_ ? i2 + 1 : 0;
214  for (; i3 != coll3.size(); i3++) {
215  r3 = coll3[i3];
216  dauC_p4 = reco::Particle::LorentzVector(r3->px(), r3->py(), r3->pz(), r3->energy());
217  if (reco::deltaR2(dauAB_p4, dauC_p4) > max_DR2_) {
218  continue;
219  }
220  bool passesMassCut = false;
221  auto const mass_ABC = (dauC_p4 + dauAB_p4).mass();
222  for (unsigned int j = 0; j < max_InvMass_.size(); j++) {
223  if ((mass_ABC >= min_InvMass_[j]) and (mass_ABC < max_InvMass_[j])) {
224  passesMassCut = true;
225  break;
226  }
227  }
228  if (passesMassCut) {
229  n++;
230  filterproduct.addObject(triggerType1_, r1);
231  filterproduct.addObject(triggerType2_, r2);
232  filterproduct.addObject(triggerType3_, r3);
233  }
234  }
235  }
236  }
237  }
238 
239  return (n >= min_N_);
240 }
241 
244 
const int triggerType3_
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
const std::vector< edm::InputTag > originTag1_
HLTTripletMass(const edm::ParameterSet &)
const int triggerType2_
edm::Ref< std::vector< T2 > > T2Ref
static PFTauRenderPlugin instance
std::string encode() const
Definition: InputTag.cc:159
const std::vector< edm::InputTag > originTag2_
const double max_DR2_
bool getCollections(edm::Event &iEvent, std::vector< T1Ref > &coll1, std::vector< T2Ref > &coll2, std::vector< T3Ref > &coll3, trigger::TriggerFilterObjectWithRefs &filterproduct) const
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken3_
const double max_DR_
const std::vector< double > min_InvMass_
uint16_t size_type
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken1_
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
~HLTTripletMass() override=default
edm::Ref< std::vector< T3 > > T3Ref
char const * label
const int min_N_
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const int triggerType1_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::Ref< std::vector< T1 > > T1Ref
const bool is2and3Same_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
const std::vector< double > max_InvMass_
HLT enums.
const bool is1and2Same_
const std::vector< edm::InputTag > originTag3_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
HLTTripletMass< reco::RecoChargedCandidate, reco::RecoChargedCandidate, reco::RecoEcalCandidate > HLT3MuonMuonPhotonMass