CMS 3D CMS Logo

Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes

HLTJetCollectionsForLeptonPlusJets< jetType > Class Template Reference

#include <HLTJetCollectionsForLeptonPlusJets.h>

Inheritance diagram for HLTJetCollectionsForLeptonPlusJets< jetType >:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 HLTJetCollectionsForLeptonPlusJets (const edm::ParameterSet &)
 ~HLTJetCollectionsForLeptonPlusJets ()

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)

Private Member Functions

virtual void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

edm::InputTag hltLeptonTag
double minDeltaR_
edm::InputTag sourceJetTag

Detailed Description

template<typename jetType>
class HLTJetCollectionsForLeptonPlusJets< jetType >

This class is an EDProducer implementing an HLT trigger for lepton and jet objects, cutting on variables relating to the jet 4-momentum representation. The producer checks for overlaps between leptons and jets and if a combination of one lepton + jets cleaned against this leptons satisfy the cuts. These jets are then added to a cleaned jet collection which is put into the event.

Date:
2012/02/06 15:09:21
Revision:
1.4
Author:
Lukasz Kreczko

Definition at line 39 of file HLTJetCollectionsForLeptonPlusJets.h.


Constructor & Destructor Documentation

template<typename jetType >
HLTJetCollectionsForLeptonPlusJets< jetType >::HLTJetCollectionsForLeptonPlusJets ( const edm::ParameterSet iConfig) [explicit]

Definition at line 25 of file HLTJetCollectionsForLeptonPlusJets.cc.

                                                                                                             :
  hltLeptonTag(iConfig.getParameter< edm::InputTag > ("HltLeptonTag")),
  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
{
  using namespace edm;
  using namespace std;
  typedef vector<RefVector<vector<jetType>,jetType,refhelper::FindUsingAdvance<vector<jetType>,jetType> > > JetCollectionVector;
  produces<JetCollectionVector> ();
}
template<typename jetType >
HLTJetCollectionsForLeptonPlusJets< jetType >::~HLTJetCollectionsForLeptonPlusJets ( )

Definition at line 37 of file HLTJetCollectionsForLeptonPlusJets.cc.

{
   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

template<typename jetType >
void HLTJetCollectionsForLeptonPlusJets< jetType >::fillDescriptions ( edm::ConfigurationDescriptions descriptions) [static]

Reimplemented from edm::EDProducer.

Definition at line 46 of file HLTJetCollectionsForLeptonPlusJets.cc.

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

                                                                                                        {
    edm::ParameterSetDescription desc;
    desc.add<edm::InputTag> ("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
    desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("caloJetCollection"));
    desc.add<double> ("minDeltaR", 0.5);
    descriptions.add(std::string("hlt")+std::string(typeid(HLTJetCollectionsForLeptonPlusJets<jetType>).name()),desc);
}
template<typename jetType >
void HLTJetCollectionsForLeptonPlusJets< jetType >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 63 of file HLTJetCollectionsForLeptonPlusJets.cc.

References deltaR(), edm::Event::getByLabel(), HLTJetCollectionsForLeptonPlusJets< jetType >::hltLeptonTag, j, HLTJetCollectionsForLeptonPlusJets< jetType >::minDeltaR_, p4, position, edm::Event::put(), HLTJetCollectionsForLeptonPlusJets< jetType >::sourceJetTag, trigger::TriggerCluster, trigger::TriggerElectron, and trigger::TriggerMuon.

{
  using namespace edm;
  using namespace std;
  
  typedef vector<RefVector<vector<jetType>,jetType,refhelper::FindUsingAdvance<vector<jetType>,jetType> > > JetCollectionVector;
  typedef vector<jetType> JetCollection;
  typedef edm::RefVector<JetCollection> JetRefVector;
  typedef edm::Ref<JetCollection> JetRef;

  Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
  iEvent.getByLabel(hltLeptonTag,PrevFilterOutput);
 
  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
  vector<Ref<reco::RecoEcalCandidateCollection> > clusCands;
  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);

  vector<Ref<reco::ElectronCollection> > eleCands;
  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
  
  vector<reco::RecoChargedCandidateRef> muonCands;
  PrevFilterOutput->getObjects(trigger::TriggerMuon,muonCands);

  Handle<JetCollection> theJetCollectionHandle;
  iEvent.getByLabel(sourceJetTag, theJetCollectionHandle);
  
  const JetCollection & theJetCollection = *theJetCollectionHandle;
  
  auto_ptr < JetCollectionVector > allSelections(new JetCollectionVector());
  
 if(!clusCands.empty()){ //try trigger cluster
    for(size_t candNr=0;candNr<clusCands.size();candNr++){  
        JetRefVector refVector;
        for (unsigned int j = 0; j < theJetCollection.size(); j++) {
          if (deltaR(clusCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
        }
    allSelections->push_back(refVector);
    }
 }

 if(!eleCands.empty()){ //try trigger cluster
    for(size_t candNr=0;candNr<eleCands.size();candNr++){  
        JetRefVector refVector;
        for (unsigned int j = 0; j < theJetCollection.size(); j++) {
          if (deltaR(eleCands[candNr]->superCluster()->position(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
        }
    allSelections->push_back(refVector);
    }
 }

 if(!muonCands.empty()){ //try trigger cluster
    for(size_t candNr=0;candNr<muonCands.size();candNr++){  
        JetRefVector refVector;
        for (unsigned int j = 0; j < theJetCollection.size(); j++) {
          if (deltaR(muonCands[candNr]->p4(),theJetCollection[j]) > minDeltaR_) refVector.push_back(JetRef(theJetCollectionHandle, j));
        }
    allSelections->push_back(refVector);
    }
 }




 iEvent.put(allSelections);
  
  return;
  
}

Member Data Documentation

template<typename jetType >
edm::InputTag HLTJetCollectionsForLeptonPlusJets< jetType >::hltLeptonTag [private]
template<typename jetType >
double HLTJetCollectionsForLeptonPlusJets< jetType >::minDeltaR_ [private]
template<typename jetType >
edm::InputTag HLTJetCollectionsForLeptonPlusJets< jetType >::sourceJetTag [private]