CMS 3D CMS Logo

HLTJetTagWithMatching.cc
Go to the documentation of this file.
1 
10 #include <vector>
11 #include <string>
12 
24 
25 #include "HLTJetTagWithMatching.h"
26 
27 
28 //
29 // constructors and destructor
30 //
31 
32 template<typename T>
34  m_Jets (config.getParameter<edm::InputTag>("Jets") ),
35  m_JetTags(config.getParameter<edm::InputTag>("JetTags") ),
36  m_MinTag (config.getParameter<double> ("MinTag") ),
37  m_MaxTag (config.getParameter<double> ("MaxTag") ),
38  m_MinJets(config.getParameter<int> ("MinJets") ),
39  m_TriggerType(config.getParameter<int> ("TriggerType") ),
40  m_deltaR (config.getParameter<double> ("deltaR") )
41 {
42  m_JetsToken = consumes<std::vector<T> >(m_Jets),
43  m_JetTagsToken = consumes<reco::JetTagCollection>(m_JetTags),
44 
45  edm::LogInfo("") << " (HLTJetTagWithMatching) trigger cuts: " << std::endl
46  << "\ttype of jets used: " << m_Jets.encode() << std::endl
47  << "\ttype of tagged jets used: " << m_JetTags.encode() << std::endl
48  << "\tmin/max tag value: [" << m_MinTag << ".." << m_MaxTag << "]" << std::endl
49  << "\tmin no. tagged jets: " << m_MinJets
50  << "\tTriggerType: " << m_TriggerType << std::endl;
51 }
52 
53 template<typename T>
55 
56 template<typename T> float HLTJetTagWithMatching<T>::findCSV(const typename std::vector<T>::const_iterator & jet, const reco::JetTagCollection & jetTags, float minDr){
57  float tmpCSV = -20 ;
58  for (auto jetb = jetTags.begin(); (jetb!=jetTags.end()); ++jetb) {
59  float tmpDr = reco::deltaR(*jet,*(jetb->first));
60 
61  if (tmpDr < minDr) {
62  minDr = tmpDr ;
63  tmpCSV= jetb->second;
64  }
65 
66  }
67  return tmpCSV;
68 
69  }
70 
71 
72 template<typename T>
73 void
77  desc.add<edm::InputTag>("Jets",edm::InputTag("hltJetCollection"));
78  desc.add<edm::InputTag>("JetTags",edm::InputTag("HLTJetTagWithMatchingCollection"));
79  desc.add<double>("MinTag",2.0);
80  desc.add<double>("MaxTag",999999.0);
81  desc.add<int>("MinJets",1);
82  desc.add<int>("TriggerType",0);
83  desc.add<double>("deltaR",0.1);
84  descriptions.add(defaultModuleLabel<HLTJetTagWithMatching<T>>(), desc);
85 }
86 
87 //
88 // member functions
89 //
90 
91 
92 // ------------ method called to produce the data ------------
93 template<typename T>
94 bool
96 {
97  using namespace std;
98  using namespace edm;
99  using namespace reco;
100 
101  typedef vector<T> TCollection;
102  typedef Ref<TCollection> TRef;
103 
105  event.getByToken(m_JetsToken, h_Jets);
106  if (saveTags()) filterproduct.addCollectionTag(m_Jets);
107 
109  event.getByToken(m_JetTagsToken, h_JetTags);
110 
111  // check if the product this one depends on is available
112  auto const & handle = h_JetTags;
113  auto const & dependent = handle->keyProduct();
114  if (not dependent.isNull() and not dependent.hasCache()) {
115  // only an empty AssociationVector can have a invalid dependent collection
116  edm::Provenance const & dependent_provenance = event.getProvenance(dependent.id());
117  if (dependent_provenance.branchDescription().dropped())
118  // FIXME the error message should be made prettier
119  throw edm::Exception(edm::errors::ProductNotFound) << "Product " << handle.provenance()->branchName() << " requires product " << dependent_provenance.branchName() << ", which has been dropped";
120  }
121 
122  TRef jetRef;
123 
124  // Look at all jets in decreasing order of Et.
125  int nJet = 0;
126  int nTag = 0;
127  for (typename TCollection::const_iterator jet = h_Jets->begin(); jet != h_Jets->end(); ++jet) {
128  jetRef = TRef(h_Jets,nJet);
129  LogTrace("") << "Jet " << nJet
130  << " : Et = " << jet->et()
131  << " , tag value = " << findCSV(jet,*h_JetTags,m_deltaR);
132  ++nJet;
133  // Check if jet is tagged.
134  if ( (m_MinTag <= findCSV(jet,*h_JetTags,m_deltaR)) and (findCSV(jet,*h_JetTags,m_deltaR) <= m_MaxTag) ) {
135  ++nTag;
136 
137  // Store a reference to the jets which passed tagging cuts
138  filterproduct.addObject(m_TriggerType,jetRef);
139  }
140  }
141 
142  // filter decision
143  bool accept = (nTag >= m_MinJets);
144 
145  edm::LogInfo("") << " trigger accept ? = " << accept
146  << " nTag/nJet = " << nTag << "/" << nJet << std::endl;
147 
148  return accept;
149 }
const_iterator end() const
~HLTJetTagWithMatching() override
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: config.py:1
std::string defaultModuleLabel()
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
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>)
edm::EDGetTokenT< std::vector< T > > m_JetsToken
static float findCSV(const typename std::vector< T >::const_iterator &jet, const reco::JetTagCollection &jetTags, float minDr=0.1)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
#define LogTrace(id)
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
BranchDescription const & branchDescription() const
Definition: Provenance.h:45
HLTJetTagWithMatching(const edm::ParameterSet &config)
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool hltFilter(edm::Event &event, const edm::EventSetup &setup, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45
HLT enums.
edm::EDGetTokenT< reco::JetTagCollection > m_JetTagsToken
std::string const & branchName() const
Definition: Provenance.h:51
const_iterator begin() const
Definition: event.py:1