18 template <
typename T1,
typename T2,
typename T3>
32 std::vector<T1Ref>& coll1,
33 std::vector<T2Ref>& coll2,
34 std::vector<T3Ref>& coll3,
57 template <
typename T1,
typename T2,
typename T3>
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")) {
78 <<
") and \"MaxInvMass\" (" <<
max_InvMass_.size() <<
") differ";
81 throw cms::Exception(
"Configuration") <<
"invalid value for parameter \"MaxDR\" (must be >= 0): " <<
max_DR_;
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")});
95 desc.add<
int>(
"triggerType1", 0);
96 desc.add<
int>(
"triggerType2", 0);
97 desc.add<
int>(
"triggerType3", 0);
99 desc.add<vector<double>>(
"MinInvMass", {0});
100 desc.add<vector<double>>(
"MaxInvMass", {1e12});
102 desc.add<
double>(
"MaxDR", 1
e4);
103 desc.add<
int>(
"MinN", 0);
105 desc.add<
bool>(
"is1and2Same",
false);
106 desc.add<
bool>(
"is2and3Same",
false);
110 template <
typename T1,
typename T2,
typename T3>
112 std::vector<T1Ref>& coll1,
113 std::vector<T2Ref>& coll2,
114 std::vector<T3Ref>& coll3,
117 if (
iEvent.getByToken(inputToken1_, handle1) and
iEvent.getByToken(inputToken2_, handle2) and
118 iEvent.getByToken(inputToken3_, handle3)) {
129 for (
unsigned int i = 0;
i < originTag1_.size(); ++
i) {
135 const auto& prov =
iEvent.getStableProvenance(pid);
140 if (tagOld.
encode() != tagNew.encode()) {
145 for (
unsigned int i = 0;
i < originTag2_.size(); ++
i) {
151 const auto& prov =
iEvent.getStableProvenance(pid);
156 if (tagOld.
encode() != tagNew.encode()) {
161 for (
unsigned int i = 0;
i < originTag3_.size(); ++
i) {
167 const auto& prov =
iEvent.getStableProvenance(pid);
172 if (tagOld.
encode() != tagNew.encode()) {
185 template <
typename T1,
typename T2,
typename T3>
193 std::vector<T1Ref> coll1;
194 std::vector<T2Ref> coll2;
195 std::vector<T3Ref> coll3;
198 if (getCollections(
iEvent, coll1, coll2, coll3, filterproduct)) {
204 for (
unsigned int i1 = 0;
i1 != coll1.size();
i1++) {
207 unsigned int i2 = is1and2Same_ ?
i1 + 1 : 0;
208 for (;
i2 != coll2.size();
i2++) {
211 dauAB_p4 = dauA_p4 + dauB_p4;
213 unsigned int i3 = is2and3Same_ ?
i2 + 1 : 0;
214 for (;
i3 != coll3.size();
i3++) {
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;
232 filterproduct.
addObject(triggerType3_, r3);
239 return (
n >= min_N_);
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 &)
edm::Ref< std::vector< T2 > > T2Ref
static PFTauRenderPlugin instance
const std::vector< edm::InputTag > originTag2_
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 std::vector< double > min_InvMass_
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
edm::Ref< std::vector< T1 > > T1Ref
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
const std::vector< double > max_InvMass_
const std::vector< edm::InputTag > originTag3_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > inputToken2_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
HLTTripletMass< reco::RecoChargedCandidate, reco::RecoChargedCandidate, reco::RecoEcalCandidate > HLT3MuonMuonPhotonMass