1 #ifndef HLTrigger_HLTfilters_plugins_HLTDoubletSinglet_h
2 #define HLTrigger_HLTfilters_plugins_HLTDoubletSinglet_h
45 template <
typename T1,
typename T2,
typename T3>
97 template <
typename T1,
typename T2,
typename T3>
100 originTag1_(iConfig.getParameter<std::
vector<edm::
InputTag>>(
"originTag1")),
101 originTag2_(iConfig.getParameter<std::
vector<edm::
InputTag>>(
"originTag2")),
102 originTag3_(iConfig.getParameter<std::
vector<edm::
InputTag>>(
"originTag3")),
103 inputTag1_(iConfig.getParameter<edm::
InputTag>(
"inputTag1")),
104 inputTag2_(iConfig.getParameter<edm::
InputTag>(
"inputTag2")),
105 inputTag3_(iConfig.getParameter<edm::
InputTag>(
"inputTag3")),
106 inputToken1_(consumes(inputTag1_)),
107 inputToken2_(consumes(inputTag2_)),
108 inputToken3_(consumes(inputTag3_)),
109 triggerType1_(iConfig.getParameter<int>(
"triggerType1")),
110 triggerType2_(iConfig.getParameter<int>(
"triggerType2")),
111 triggerType3_(iConfig.getParameter<int>(
"triggerType3")),
112 min_Dphi_(iConfig.getParameter<double>(
"MinDphi")),
113 max_Dphi_(iConfig.getParameter<double>(
"MaxDphi")),
114 min_Deta_(iConfig.getParameter<double>(
"MinDeta")),
115 max_Deta_(iConfig.getParameter<double>(
"MaxDeta")),
116 min_Minv_(iConfig.getParameter<double>(
"MinMinv")),
117 max_Minv_(iConfig.getParameter<double>(
"MaxMinv")),
118 min_DelR_(iConfig.getParameter<double>(
"MinDelR")),
119 max_DelR_(iConfig.getParameter<double>(
"MaxDelR")),
120 min_Pt_(iConfig.getParameter<double>(
"MinPt")),
121 max_Pt_(iConfig.getParameter<double>(
"MaxPt")),
122 min_N_(iConfig.getParameter<int>(
"MinN")),
126 min_DelR2_(min_DelR_ < 0 ? 0 : min_DelR_ * min_DelR_),
127 max_DelR2_(max_DelR_ < 0 ? 0 : max_DelR_ * max_DelR_),
128 cutdphi_(min_Dphi_ <= max_Dphi_),
129 cutdeta_(min_Deta_ <= max_Deta_),
130 cutminv_(min_Minv_ <= max_Minv_),
131 cutdelr_(min_DelR_ <= max_DelR_),
132 cutpt_(min_Pt_ <= max_Pt_)
146 <<
"Warning: The deltaR requirement is active, but its range is invalid: DelR [" <<
min_DelR_ <<
" "
150 template <
typename T1,
typename T2,
typename T3>
153 template <
typename T1,
typename T2,
typename T3>
156 makeHLTFilterDescription(desc);
157 desc.
add<std::vector<edm::InputTag>>(
"originTag1", {
edm::InputTag(
"hltOriginal1")});
158 desc.
add<std::vector<edm::InputTag>>(
"originTag2", {
edm::InputTag(
"hltOriginal2")});
159 desc.
add<std::vector<edm::InputTag>>(
"originTag3", {
edm::InputTag(
"hltOriginal3")});
163 desc.
add<
int>(
"triggerType1", 0);
164 desc.
add<
int>(
"triggerType2", 0);
165 desc.
add<
int>(
"triggerType3", 0);
166 desc.
add<
double>(
"MinDphi", +1.0);
167 desc.
add<
double>(
"MaxDphi", -1.0);
168 desc.
add<
double>(
"MinDeta", +1.0);
169 desc.
add<
double>(
"MaxDeta", -1.0);
170 desc.
add<
double>(
"MinMinv", +1.0);
171 desc.
add<
double>(
"MaxMinv", -1.0);
172 desc.
add<
double>(
"MinDelR", +1.0);
173 desc.
add<
double>(
"MaxDelR", -1.0);
174 desc.
add<
double>(
"MinPt", +1.0);
175 desc.
add<
double>(
"MaxPt", -1.0);
176 desc.
add<
int>(
"MinN", 1);
185 template <
typename T1,
typename T2,
typename T3>
191 using namespace reco;
192 using namespace trigger;
200 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 0 ";
203 std::vector<T1Ref> coll1;
204 auto const& objsWithRefs1 = iEvent.
get(inputToken1_);
205 objsWithRefs1.getObjects(triggerType1_, coll1);
206 std::vector<T2Ref> coll2;
207 auto const& objsWithRefs2 = iEvent.
get(inputToken2_);
208 objsWithRefs2.getObjects(triggerType2_, coll2);
209 std::vector<T3Ref> coll3;
210 auto const& objsWithRefs3 = iEvent.
get(inputToken3_);
211 objsWithRefs3.getObjects(triggerType3_, coll3);
219 for (
size_t i = 0;
i < originTag1_.size(); ++
i) {
221 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 1a/" << i <<
" "
222 << originTag1_[
i].encode();
228 const string&
label(prov.moduleLabel());
229 const string&
instance(prov.productInstanceName());
230 const string&
process(prov.processName());
232 if (tagOld.encode() != tagNew.encode()) {
235 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 1b " << tagNew.encode();
238 for (
size_t i = 0;
i < originTag2_.size(); ++
i) {
240 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 2a/" << i <<
" "
241 << originTag2_[
i].encode();
247 const string&
label(prov.moduleLabel());
248 const string&
instance(prov.productInstanceName());
249 const string&
process(prov.processName());
251 if (tagOld.encode() != tagNew.encode()) {
254 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 2b " << tagNew.encode();
257 for (
size_t i = 0;
i < originTag3_.size(); ++
i) {
259 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 3a/" << i <<
" "
260 << originTag3_[
i].encode();
266 const string&
label(prov.moduleLabel());
267 const string&
instance(prov.productInstanceName());
268 const string&
process(prov.processName());
270 if (tagOld.encode() != tagNew.encode()) {
273 LogVerbatim(
"HLTDoubletSinglet") <<
" moduleLabel: " << moduleLabel() <<
" 3b " << tagNew.encode();
281 for (
size_t i1 = 0; i1 != n1; i1++) {
284 auto const i2_min = (same12_ ? i1 + 1 : 0);
285 for (
size_t i2 = i2_min; i2 != n2; i2++) {
289 auto const i3_min = (same23_ ? i2_min + 1 : (same13_ ? i1 + 1 : 0));
290 for (
size_t i3 = i3_min; i3 != n3; i3++) {
296 if (cutdphi_ && (min_Dphi_ > dPhi13 || dPhi13 > max_Dphi_))
299 if (cutdphi_ && (min_Dphi_ > dPhi23 || dPhi23 > max_Dphi_))
303 auto const dEta13(
std::abs(p1.eta() - p3.eta()));
304 if (cutdeta_ && (min_Deta_ > dEta13 || dEta13 > max_Deta_))
306 auto const dEta23(
std::abs(p2.eta() - p3.eta()));
307 if (cutdeta_ && (min_Deta_ > dEta23 || dEta23 > max_Deta_))
311 auto const delR2_13(dPhi13 * dPhi13 + dEta13 * dEta13);
312 if (cutdelr_ && (min_DelR2_ > delR2_13 || delR2_13 > max_DelR2_))
314 auto const delR2_23(dPhi23 * dPhi23 + dEta23 * dEta23);
315 if (cutdelr_ && (min_DelR2_ > delR2_23 || delR2_23 > max_DelR2_))
320 auto const mInv13(
std::abs(p13.mass()));
321 if (cutminv_ && (min_Minv_ > mInv13 || mInv13 > max_Minv_))
323 auto const pt13(p13.pt());
324 if (cutpt_ && (min_Pt_ > pt13 || pt13 > max_Pt_))
328 auto const mInv23(
std::abs(p23.mass()));
329 if (cutminv_ && (min_Minv_ > mInv23 || mInv23 > max_Minv_))
331 auto const pt23(p23.pt());
332 if (cutpt_ && (min_Pt_ > pt23 || pt23 > max_Pt_))
336 filterproduct.
addObject(triggerType1_, r1);
337 filterproduct.
addObject(triggerType2_, r2);
338 filterproduct.
addObject(triggerType3_, r3);
345 return (n >= min_N_);
348 #endif //HLTrigger_HLTfilters_plugins_HLTDoubletSinglet_h
Log< level::Info, true > LogVerbatim
static PFTauRenderPlugin instance
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
std::vector< T1 > T1Collection
std::string defaultModuleLabel()
const std::vector< edm::InputTag > originTag1_
edm::Ref< T2Collection > T2Ref
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
~HLTDoubletSinglet() override
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken3_
const edm::InputTag inputTag1_
edm::Ref< T1Collection > T1Ref
const edm::InputTag inputTag3_
Abs< T >::type abs(const T &t)
edm::Ref< T3Collection > T3Ref
bool get(ProductID const &oid, Handle< PROD > &result) const
std::vector< T3 > T3Collection
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const std::string * moduleLabel() const
const std::vector< edm::InputTag > originTag3_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::InputTag inputTag2_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken1_
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
const std::vector< edm::InputTag > originTag2_
std::vector< T2 > T2Collection
Log< level::Warning, false > LogWarning
HLTDoubletSinglet(const edm::ParameterSet &)
StableProvenance const & getStableProvenance(BranchID const &theID) const
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_