Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Backport of HLTJetTagWithMatching in 73X (#6985) #7087

Merged
merged 2 commits into from
Jan 13, 2015
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
150 changes: 150 additions & 0 deletions HLTrigger/btau/src/HLTJetTagWithMatching.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,150 @@
/** \class HLTJetTagWithMatching
*
* This class is an HLTFilter (a spcialized EDFilter) implementing
* tagged multi-jet trigger for b and tau.
* It should be run after the normal multi-jet trigger.
* It is like HLTJetTag, but it loop on jets collection instead of jetTag collection.
*
*/

#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/Common/interface/RefToBase.h"
#include "DataFormats/JetReco/interface/CaloJetCollection.h"
#include "DataFormats/BTauReco/interface/JetTag.h"
#include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
#include "HLTrigger/HLTcore/interface/HLTFilter.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "DataFormats/Math/interface/deltaR.h"

#include "HLTJetTagWithMatching.h"

#include<vector>
#include<string>
#include<typeinfo>

//
// constructors and destructor
//

template<typename T>
HLTJetTagWithMatching<T>::HLTJetTagWithMatching(const edm::ParameterSet & config) : HLTFilter(config),
m_Jets (config.getParameter<edm::InputTag>("Jets") ),
m_JetTags(config.getParameter<edm::InputTag>("JetTags") ),
m_MinTag (config.getParameter<double> ("MinTag") ),
m_MaxTag (config.getParameter<double> ("MaxTag") ),
m_MinJets(config.getParameter<int> ("MinJets") ),
m_TriggerType(config.getParameter<int> ("TriggerType") ),
m_deltaR (config.getParameter<double> ("deltR") )
{
m_JetsToken = consumes<std::vector<T> >(m_Jets),
m_JetTagsToken = consumes<reco::JetTagCollection>(m_JetTags),

edm::LogInfo("") << " (HLTJetTagWithMatching) trigger cuts: " << std::endl
<< "\ttype of jets used: " << m_Jets.encode() << std::endl
<< "\ttype of tagged jets used: " << m_JetTags.encode() << std::endl
<< "\tmin/max tag value: [" << m_MinTag << ".." << m_MaxTag << "]" << std::endl
<< "\tmin no. tagged jets: " << m_MinJets
<< "\tTriggerType: " << m_TriggerType << std::endl;
}

template<typename T>
HLTJetTagWithMatching<T>::~HLTJetTagWithMatching()
{
}

template<typename T> float HLTJetTagWithMatching<T>::findCSV(const typename std::vector<T>::const_iterator & jet, const reco::JetTagCollection & jetTags, float minDr=0.1 ){
float tmpCSV = -20 ;
for (reco::JetTagCollection::const_iterator jetb = jetTags.begin(); (jetb!=jetTags.end()); ++jetb) {
float tmpDr = reco::deltaR(*jet,*(jetb->first));

if (tmpDr < minDr) {
minDr = tmpDr ;
tmpCSV= jetb->second;
}

}
return tmpCSV;

}


template<typename T>
void
HLTJetTagWithMatching<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
makeHLTFilterDescription(desc);
desc.add<edm::InputTag>("Jets",edm::InputTag("hltJetCollection"));
desc.add<edm::InputTag>("JetTags",edm::InputTag("HLTJetTagWithMatchingCollection"));
desc.add<double>("MinTag",2.0);
desc.add<double>("MaxTag",999999.0);
desc.add<int>("MinJets",1);
desc.add<int>("TriggerType",0);
desc.add<double>("deltaR",0.1);
descriptions.add(std::string("hlt")+std::string(typeid(HLTJetTagWithMatching<T>).name()),desc);
}

//
// member functions
//


// ------------ method called to produce the data ------------
template<typename T>
bool
HLTJetTagWithMatching<T>::hltFilter(edm::Event& event, const edm::EventSetup& setup, trigger::TriggerFilterObjectWithRefs & filterproduct) const
{
using namespace std;
using namespace edm;
using namespace reco;

typedef vector<T> TCollection;
typedef Ref<TCollection> TRef;

edm::Handle<TCollection> h_Jets;
event.getByToken(m_JetsToken, h_Jets);
if (saveTags()) filterproduct.addCollectionTag(m_Jets);

edm::Handle<JetTagCollection> h_JetTags;
event.getByToken(m_JetTagsToken, h_JetTags);

// check if the product this one depends on is available
auto const & handle = h_JetTags;
auto const & dependent = handle->keyProduct();
if (not dependent.isNull() and not dependent.hasCache()) {
// only an empty AssociationVector can have a invalid dependent collection
edm::Provenance const & dependent_provenance = event.getProvenance(dependent.id());
if (dependent_provenance.constBranchDescription().dropped())
// FIXME the error message should be made prettier
throw edm::Exception(edm::errors::ProductNotFound) << "Product " << handle.provenance()->branchName() << " requires product " << dependent_provenance.branchName() << ", which has been dropped";
}

TRef jetRef;

// Look at all jets in decreasing order of Et.
int nJet = 0;
int nTag = 0;
for (typename TCollection::const_iterator jet = h_Jets->begin(); jet != h_Jets->end(); ++jet) {
jetRef = TRef(h_Jets,nJet);
LogTrace("") << "Jet " << nJet
<< " : Et = " << jet->et()
<< " , tag value = " << findCSV(jet,*h_JetTags,m_deltaR);
++nJet;
// Check if jet is tagged.
if ( (m_MinTag <= findCSV(jet,*h_JetTags,m_deltaR)) and (findCSV(jet,*h_JetTags,m_deltaR) <= m_MaxTag) ) {
++nTag;

// Store a reference to the jets which passed tagging cuts
filterproduct.addObject(m_TriggerType,jetRef);
}
}

// filter decision
bool accept = (nTag >= m_MinJets);

edm::LogInfo("") << " trigger accept ? = " << accept
<< " nTag/nJet = " << nTag << "/" << nJet << std::endl;

return accept;
}
51 changes: 51 additions & 0 deletions HLTrigger/btau/src/HLTJetTagWithMatching.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,51 @@
#ifndef HLTrigger_btau_HLTJetTagWithMatching_h
#define HLTrigger_btau_HLTJetTagWithMatching_h

/** \class HLTJetTagWithMatching
*
* This class is an HLTFilter (a spcialized EDFilter) implementing
* tagged multi-jet trigger for b and tau.
* It should be run after the normal multi-jet trigger.
* It is like HLTJetTag, but it loop on jets collection instead of jetTag collection.
*
*/

#include <string>

#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "HLTrigger/HLTcore/interface/HLTFilter.h"
#include "DataFormats/BTauReco/interface/JetTag.h"

namespace edm {
class ConfigurationDescriptions;
}

//
// class declaration
//


template<typename T>
class HLTJetTagWithMatching : public HLTFilter {

public:
explicit HLTJetTagWithMatching(const edm::ParameterSet & config);
~HLTJetTagWithMatching();
static float findCSV(const typename std::vector<T>::const_iterator & jet, const reco::JetTagCollection & jetTags, float minDr);
static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
virtual bool hltFilter(edm::Event & event, const edm::EventSetup & setup, trigger::TriggerFilterObjectWithRefs & filterproduct) const override;

private:
edm::InputTag m_Jets; // module label of input JetCollection
edm::EDGetTokenT<std::vector<T> > m_JetsToken;
edm::InputTag m_JetTags; // module label of input JetTagCollection
edm::EDGetTokenT<reco::JetTagCollection> m_JetTagsToken;
double m_MinTag, m_MaxTag; // tag descriminator cuts applied to each jet
int m_MinJets; // min. number of jets required to be tagged
int m_TriggerType;
double m_deltaR; // deltaR used to match jet with jetTags

};

#endif // HLTrigger_btau_HLTJetTagWithMatching_h
7 changes: 7 additions & 0 deletions HLTrigger/btau/src/modules.cc
Original file line number Diff line number Diff line change
Expand Up @@ -59,6 +59,13 @@ typedef HLTJetTag<reco:: PFJet> HLTPFJetTag;
DEFINE_FWK_MODULE(HLTCaloJetTag);
DEFINE_FWK_MODULE(HLTPFJetTag);

#include "HLTJetTagWithMatching.h"
#include "HLTJetTagWithMatching.cc"
typedef HLTJetTagWithMatching<reco::CaloJet> HLTCaloJetTagWithMatching;
typedef HLTJetTagWithMatching<reco:: PFJet> HLTPFJetTagWithMatching;
DEFINE_FWK_MODULE(HLTCaloJetTagWithMatching);
DEFINE_FWK_MODULE(HLTPFJetTagWithMatching);

#include "HLTCollectionProducer.h"
typedef HLTCollectionProducer<reco::CaloJet> HLTCaloJetCollectionProducer;
typedef HLTCollectionProducer<reco::PFJet> HLTPFJetCollectionProducer;
Expand Down