#include <RecoBTag/InvariantMass/src/InvariantMass.cc>
Public Member Functions | |
InvariantMass (const edm::ParameterSet &) | |
double | operator() (const T1 &t1, const T2 &t2) const |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
~InvariantMass () | |
Private Member Functions | |
reco::BasicClusterRefProd | getSelectedClusters (reco::BasicClusterRefProd bclus, reco::BasicClusterRefProd eclus, reco::Jet &jet) |
Private Attributes | |
std::string | jetTrackSrc |
InvariantMassAlgorithm * | m_algo |
std::string | m_ecalBClSrc |
EDProducer of the tagged TauJet with the InvariantMassAlgorithm. It returns two collections: base collection is the JetTag, and extended Collection which is the IsolatedTauTagInfo. The method implemented in the IsolatedTauTagInfo class are used to compute the discriminator variable. A trick is used to link the IsolatedTauTagInfo to a smart reference to the JetTag.
Description: <one line="" class="" summary>="">
Implementation: <Notes on="" implementation>="">
Definition at line 6 of file InvariantMass.h.
InvariantMass< T1, T2 >::InvariantMass | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 50 of file InvariantMass.cc.
References edm::ParameterSet::getParameter().
{ jetTrackSrc = iConfig.getParameter<string>("JetTrackSrc"); m_ecalBClSrc = iConfig.getParameter<string>("ecalbcl"); m_algo = new InvariantMassAlgorithm(iConfig); produces<reco::JetTagCollection>(); produces<reco::TauMassTagInfoCollection>(); }
InvariantMass< T1, T2 >::~InvariantMass | ( | ) |
Definition at line 62 of file InvariantMass.cc.
{ delete m_algo; }
reco::BasicClusterRefProd InvariantMass< T1, T2 >::getSelectedClusters | ( | reco::BasicClusterRefProd | bclus, |
reco::BasicClusterRefProd | eclus, | ||
reco::Jet & | jet | ||
) | [private] |
double InvariantMass< T1, T2 >::operator() | ( | const T1 & | t1, |
const T2 & | t2 | ||
) | const [inline] |
Definition at line 7 of file InvariantMass.h.
{
return ( t1.momentum() + t2.momentum() ).mass();
}
void InvariantMass< T1, T2 >::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 72 of file InvariantMass.cc.
References edm::Event::getByLabel(), i, metsig::jet, parseEventContent::prod, edm::Event::put(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), and dt_dqm_sourceclient_common_cff::reco.
{ using namespace edm; using namespace reco; Handle<IsolatedTauTagInfoCollection> isolatedTaus; iEvent.getByLabel(jetTrackSrc, isolatedTaus); std::auto_ptr<JetTagCollection> tagCollection; std::auto_ptr<TauMassTagInfoCollection> extCollection( new TauMassTagInfoCollection() ); // Island basic cluster collection Handle<BasicClusterCollection> barrelBasicClusterHandle; iEvent.getByLabel(m_ecalBClSrc, "islandBarrelBasicClusters", barrelBasicClusterHandle); Handle<BasicClusterCollection> endcapBasicClusterHandle; iEvent.getByLabel(m_ecalBClSrc, "islandEndcapBasicClusters", endcapBasicClusterHandle); if (isolatedTaus->empty()) { tagCollection.reset( new JetTagCollection() ); } else { RefToBaseProd<reco::Jet> prod( isolatedTaus->begin()->jet() ); tagCollection.reset( new JetTagCollection(RefToBaseProd<reco::Jet>(prod)) ); for (unsigned int i = 0; i < isolatedTaus->size(); ++i) { IsolatedTauTagInfoRef tauRef(isolatedTaus, i); const Jet & jet = *(tauRef->jet()); math::XYZVector jetDir(jet.px(),jet.py(),jet.pz()); pair<double,TauMassTagInfo> jetTauPair; if (jetDir.eta() < 1.2) // barrel jetTauPair = m_algo->tag(iEvent, iSetup, tauRef, barrelBasicClusterHandle); else // endcap jetTauPair = m_algo->tag(iEvent, iSetup, tauRef, endcapBasicClusterHandle); tagCollection->setValue( i, jetTauPair.first ); extCollection->push_back( jetTauPair.second ); } } iEvent.put( tagCollection ); iEvent.put( extCollection ); }
std::string InvariantMass< T1, T2 >::jetTrackSrc [private] |
Definition at line 37 of file InvariantMass.h.
InvariantMassAlgorithm* InvariantMass< T1, T2 >::m_algo [private] |
Definition at line 36 of file InvariantMass.h.
std::string InvariantMass< T1, T2 >::m_ecalBClSrc [private] |
Definition at line 38 of file InvariantMass.h.