CMS 3D CMS Logo

Public Member Functions | Static Public Member Functions | Private Attributes

HLTJetL1MatchProducer< T > Class Template Reference

#include <HLTJetL1MatchProducer.h>

Inheritance diagram for HLTJetL1MatchProducer< T >:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual void beginJob ()
 HLTJetL1MatchProducer (const edm::ParameterSet &)
virtual void produce (edm::Event &, const edm::EventSetup &)
 ~HLTJetL1MatchProducer ()

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)

Private Attributes

double DeltaR_
edm::InputTag jetsInput_
edm::InputTag L1CenJets_
edm::InputTag L1ForJets_
edm::InputTag L1TauJets_

Detailed Description

template<typename T>
class HLTJetL1MatchProducer< T >

Definition at line 13 of file HLTJetL1MatchProducer.h.


Constructor & Destructor Documentation

template<typename T >
HLTJetL1MatchProducer< T >::HLTJetL1MatchProducer ( const edm::ParameterSet iConfig) [explicit]

Definition at line 21 of file HLTJetL1MatchProducer.cc.

{
  jetsInput_ = iConfig.template getParameter<edm::InputTag>("jetsInput");
  L1TauJets_ = iConfig.template getParameter<edm::InputTag>("L1TauJets");
  L1CenJets_ = iConfig.template getParameter<edm::InputTag>("L1CenJets");
  L1ForJets_ = iConfig.template getParameter<edm::InputTag>("L1ForJets");
  DeltaR_ = iConfig.template getParameter<double>("DeltaR");

  typedef std::vector<T> TCollection;
  produces<TCollection> ();

}
template<typename T >
HLTJetL1MatchProducer< T >::~HLTJetL1MatchProducer ( )

Definition at line 41 of file HLTJetL1MatchProducer.cc.

{

}

Member Function Documentation

template<typename T >
void HLTJetL1MatchProducer< T >::beginJob ( void  ) [virtual]

Reimplemented from edm::EDProducer.

Definition at line 35 of file HLTJetL1MatchProducer.cc.

{

}
template<typename T >
void HLTJetL1MatchProducer< T >::fillDescriptions ( edm::ConfigurationDescriptions descriptions) [static]

Reimplemented from edm::EDProducer.

Definition at line 47 of file HLTJetL1MatchProducer.cc.

References edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::add(), and mergeVDriftHistosByStation::name.

                                                                                          {
  edm::ParameterSetDescription desc;
  desc.add<edm::InputTag>("jetsInput",edm::InputTag("hltAntiKT5PFJets"));
  desc.add<edm::InputTag>("L1TauJets",edm::InputTag("hltL1extraParticles","Tau"));
  desc.add<edm::InputTag>("L1CenJets",edm::InputTag("hltL1extraParticles","Central"));
  desc.add<edm::InputTag>("L1ForJets",edm::InputTag("hltL1extraParticles","Forward"));
  desc.add<double>("DeltaR",0.5);
  descriptions.add(std::string("hlt")+std::string(typeid(HLTJetL1MatchProducer<T>).name()),desc);
}
template<typename T >
void HLTJetL1MatchProducer< T >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 58 of file HLTJetL1MatchProducer.cc.

References SiPixelRawToDigiRegional_cfi::deltaPhi, edm::Event::getByLabel(), trackerHitRTTI::isMatched(), analyzePatCleaning_cfg::jets, edm::Event::put(), query::result, and mathSSE::sqrt().

{

  typedef std::vector<T> TCollection;

  edm::Handle<TCollection> jets;
  iEvent.getByLabel(jetsInput_, jets);

  std::auto_ptr<TCollection> result (new TCollection);


  edm::Handle<l1extra::L1JetParticleCollection> l1TauJets;
  iEvent.getByLabel(L1TauJets_,l1TauJets);

  edm::Handle<l1extra::L1JetParticleCollection> l1CenJets;
  iEvent.getByLabel(L1CenJets_,l1CenJets);

  edm::Handle<l1extra::L1JetParticleCollection> l1ForJets;
  iEvent.getByLabel(L1ForJets_,l1ForJets);

  typename TCollection::const_iterator jet_iter;
  for (jet_iter = jets->begin(); jet_iter != jets->end(); ++jet_iter) {

    bool isMatched=false;

    //std::cout << "FL: l1TauJets.size  = " << l1TauJets->size() << std::endl;
    for (unsigned int jetc=0;jetc<l1TauJets->size();++jetc)
    {
      const double deltaeta=jet_iter->eta()-(*l1TauJets)[jetc].eta();
      const double deltaphi=deltaPhi(jet_iter->phi(),(*l1TauJets)[jetc].phi());
      //std::cout << "FL: sqrt(2) = " << sqrt(2) << std::endl;
      if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
    }

    for (unsigned int jetc=0;jetc<l1CenJets->size();++jetc)
    {
      const double deltaeta=jet_iter->eta()-(*l1CenJets)[jetc].eta();
      const double deltaphi=deltaPhi(jet_iter->phi(),(*l1CenJets)[jetc].phi());
      if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
    }

    for (unsigned int jetc=0;jetc<l1ForJets->size();++jetc)
    {
      const double deltaeta=jet_iter->eta()-(*l1ForJets)[jetc].eta();
      const double deltaphi=deltaPhi(jet_iter->phi(),(*l1ForJets)[jetc].phi());
      if (sqrt(deltaeta*deltaeta+deltaphi*deltaphi) < DeltaR_) isMatched=true;
    }


    if (isMatched==true) result->push_back(*jet_iter);

  } // jet_iter

  iEvent.put( result);

}

Member Data Documentation

template<typename T >
double HLTJetL1MatchProducer< T >::DeltaR_ [private]

Definition at line 26 of file HLTJetL1MatchProducer.h.

template<typename T >
edm::InputTag HLTJetL1MatchProducer< T >::jetsInput_ [private]

Definition at line 21 of file HLTJetL1MatchProducer.h.

template<typename T >
edm::InputTag HLTJetL1MatchProducer< T >::L1CenJets_ [private]

Definition at line 23 of file HLTJetL1MatchProducer.h.

template<typename T >
edm::InputTag HLTJetL1MatchProducer< T >::L1ForJets_ [private]

Definition at line 24 of file HLTJetL1MatchProducer.h.

template<typename T >
edm::InputTag HLTJetL1MatchProducer< T >::L1TauJets_ [private]

Definition at line 22 of file HLTJetL1MatchProducer.h.