![]() |
![]() |
#include <HLTJetCollForElePlusJets.h>
Public Member Functions | |
HLTJetCollForElePlusJets (const edm::ParameterSet &) | |
~HLTJetCollForElePlusJets () | |
Static Public Member Functions | |
static void | fillDescriptions (edm::ConfigurationDescriptions &descriptions) |
Private Member Functions | |
virtual void | beginJob () |
virtual void | endJob () |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
edm::InputTag | hltElectronTag |
double | maxAbsJetEta_ |
double | minDeltaEta_ |
double | minDeltaR_ |
double | minJetPt_ |
unsigned int | minNJets_ |
double | minSoftJetPt_ |
edm::InputTag | sourceJetTag |
This class is an EDProducer implementing an HLT trigger for electron and jet objects, cutting on variables relating to the jet 4-momentum representation. The producer checks for overlaps between electrons and jets and if a combination of one electron + jets cleaned against this electrons satisfy the cuts. These jets are then added to a cleaned jet collection which is put into the event.
Definition at line 39 of file HLTJetCollForElePlusJets.h.
HLTJetCollForElePlusJets::HLTJetCollForElePlusJets | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 18 of file HLTJetCollForElePlusJets.cc.
: hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")), sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")), minJetPt_(iConfig.getParameter<double> ("MinJetPt")), maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")), minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")), minDeltaR_(iConfig.getParameter< double > ("minDeltaR")), //Only for VBF minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")), minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta")) { produces<reco::CaloJetCollection>(); }
HLTJetCollForElePlusJets::~HLTJetCollForElePlusJets | ( | ) |
Definition at line 33 of file HLTJetCollForElePlusJets.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void HLTJetCollForElePlusJets::beginJob | ( | void | ) | [private, virtual] |
void HLTJetCollForElePlusJets::endJob | ( | void | ) | [private, virtual] |
void HLTJetCollForElePlusJets::fillDescriptions | ( | edm::ConfigurationDescriptions & | descriptions | ) | [static] |
Reimplemented from edm::EDProducer.
Definition at line 40 of file HLTJetCollForElePlusJets.cc.
References edm::ParameterSetDescription::add(), and edm::ConfigurationDescriptions::add().
{ edm::ParameterSetDescription desc; desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs")); desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("caloJetCollection")); desc.add<double> ("MinJetPt", 30.); desc.add<double> ("MaxAbsJetEta", 2.6); desc.add<unsigned int> ("MinNJets", 1); desc.add<double> ("minDeltaR", 0.5); //Only for VBF desc.add<double> ("MinSoftJetPt", 25.); desc.add<double> ("MinDeltaEta", -1.); descriptions.add("hltJetCollForElePlusJets", desc); }
void HLTJetCollForElePlusJets::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 62 of file HLTJetCollForElePlusJets.cc.
References abs, edm::Event::getByLabel(), hltElectronTag, i, j, gen::k, maxAbsJetEta_, minDeltaEta_, minDeltaR_, minJetPt_, minNJets_, minSoftJetPt_, position, edm::RefVector< C, T, F >::push_back(), edm::Event::put(), python::multivaluedict::sort(), sourceJetTag, trigger::TriggerCluster, trigger::TriggerElectron, x, detailsBasic3DVector::y, and z.
{ using namespace edm; edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput; iEvent.getByLabel(hltElectronTag,PrevFilterOutput); //its easier on the if statement flow if I try everything at once, shouldnt add to timing std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands; PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands); std::vector<edm::Ref<reco::ElectronCollection> > eleCands; PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands); //prepare the collection of 3-D vector for electron momenta std::vector<TVector3> ElePs; if(!clusCands.empty()){ //try trigger cluster for(size_t candNr=0;candNr<clusCands.size();candNr++){ TVector3 positionVector( clusCands[candNr]->superCluster()->position().x(), clusCands[candNr]->superCluster()->position().y(), clusCands[candNr]->superCluster()->position().z()); ElePs.push_back(positionVector); } }else if(!eleCands.empty()){ // try trigger electrons for(size_t candNr=0;candNr<eleCands.size();candNr++){ TVector3 positionVector( eleCands[candNr]->superCluster()->position().x(), eleCands[candNr]->superCluster()->position().y(), eleCands[candNr]->superCluster()->position().z()); ElePs.push_back(positionVector); } } edm::Handle<reco::CaloJetCollection> theCaloJetCollectionHandle; iEvent.getByLabel(sourceJetTag, theCaloJetCollectionHandle); //const reco::CaloJetCollection* theCaloJetCollection = theCaloJetCollectionHandle.product(); const reco::CaloJetCollection & theCaloJetCollection = *theCaloJetCollectionHandle; std::auto_ptr< reco::CaloJetCollection > theFilteredCaloJetCollection(new reco::CaloJetCollection); std::auto_ptr < std::vector<reco::CaloJetRefVector> > allSelections(new std::vector<reco::CaloJetRefVector>()); bool foundSolution(false); for (unsigned int i = 0; i < ElePs.size(); i++) { bool VBFJetPair = false; std::vector<int> store_jet; reco::CaloJetRefVector refVector; for (unsigned int j = 0; j < theCaloJetCollection.size(); j++) { TVector3 JetP(theCaloJetCollection[j].px(), theCaloJetCollection[j].py(), theCaloJetCollection[j].pz()); double DR = ElePs[i].DeltaR(JetP); if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR > minDeltaR_) { store_jet.push_back(j); // The VBF part of the filter if ( minDeltaEta_ > 0 ) { for ( unsigned int k = j+1; k < theCaloJetCollection.size(); k++ ) { TVector3 SoftJetP(theCaloJetCollection[k].px(), theCaloJetCollection[k].py(), theCaloJetCollection[k].pz()); double softDR = ElePs[i].DeltaR(SoftJetP); if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR > minDeltaR_) if ( std::abs(SoftJetP.Eta() - JetP.Eta()) > minDeltaEta_ ) { store_jet.push_back(k); VBFJetPair = true; } } } } } // Now remove duplicates from the jet collection to store std::sort( store_jet.begin(), store_jet.end() ); store_jet.erase( unique( store_jet.begin(), store_jet.end() ), store_jet.end() ); // Now save the cleaned jets for ( unsigned int ijet = 0; ijet < store_jet.size(); ijet++ ) { //store all selections refVector.push_back(reco::CaloJetRef(theCaloJetCollectionHandle, store_jet.at(ijet))); //store first selection which matches the criteria if(!foundSolution) theFilteredCaloJetCollection->push_back(theCaloJetCollection[store_jet.at(ijet)]); } //store all selections allSelections->push_back(refVector); if (theFilteredCaloJetCollection->size() >= minNJets_ && minDeltaEta_ < 0) foundSolution = true; else if (VBFJetPair && minDeltaEta_ > 0) foundSolution = true; else if (!foundSolution) theFilteredCaloJetCollection->clear(); } iEvent.put(theFilteredCaloJetCollection); return; }
Definition at line 50 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
double HLTJetCollForElePlusJets::maxAbsJetEta_ [private] |
Definition at line 54 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
double HLTJetCollForElePlusJets::minDeltaEta_ [private] |
Definition at line 60 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
double HLTJetCollForElePlusJets::minDeltaR_ [private] |
Definition at line 57 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
double HLTJetCollForElePlusJets::minJetPt_ [private] |
Definition at line 53 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
unsigned int HLTJetCollForElePlusJets::minNJets_ [private] |
Definition at line 55 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
double HLTJetCollForElePlusJets::minSoftJetPt_ [private] |
Definition at line 59 of file HLTJetCollForElePlusJets.h.
Referenced by produce().
Definition at line 51 of file HLTJetCollForElePlusJets.h.
Referenced by produce().