CMS 3D CMS Logo

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

HLTJetCollectionsForElePlusJets< T > Class Template Reference

#include <HLTJetCollectionsForElePlusJets.h>

Inheritance diagram for HLTJetCollectionsForElePlusJets< T >:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

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

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 hltElectronTag
double minDeltaR_
edm::InputTag sourceJetTag

Detailed Description

template<typename T>
class HLTJetCollectionsForElePlusJets< T >

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.

Date:
2012/02/12 10:24:17
Revision:
1.2
Author:
Lukasz Kreczko

Definition at line 40 of file HLTJetCollectionsForElePlusJets.h.


Constructor & Destructor Documentation

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

Definition at line 24 of file HLTJetCollectionsForElePlusJets.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"))
{
  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
  produces<TCollectionVector> ();
}

Definition at line 41 of file HLTJetCollectionsForElePlusJets.cc.

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

}

Member Function Documentation

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

Reimplemented from edm::EDProducer.

Definition at line 49 of file HLTJetCollectionsForElePlusJets.cc.

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

                                                                                                    {
    edm::ParameterSetDescription desc;
    desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
    desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
   // 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(std::string("hlt")+std::string(typeid(HLTJetCollectionsForElePlusJets<T>).name()),desc);
}
template<typename T >
void HLTJetCollectionsForElePlusJets< T >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 71 of file HLTJetCollectionsForElePlusJets.cc.

References edm::Event::getByLabel(), i, j, position, edm::Event::put(), trigger::TriggerCluster, trigger::TriggerElectron, x, detailsBasic3DVector::y, and z.

{
  using namespace edm;
  using namespace std;

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

  typedef edm::RefVector<TCollection> TRefVector;

  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
  
  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<TCollection> theJetCollectionHandle;
  iEvent.getByLabel(sourceJetTag, theJetCollectionHandle);
  //const TCollection* theJetCollection = theJetCollectionHandle.product();
  
  const TCollection & theJetCollection = *theJetCollectionHandle;
  
  //std::auto_ptr< TCollection >  theFilteredJetCollection(new TCollection);
  
  std::auto_ptr < TCollectionVector > allSelections(new TCollectionVector());
  
 //bool foundSolution(false);

    for (unsigned int i = 0; i < ElePs.size(); i++) {

       // bool VBFJetPair = false;
        //std::vector<int> store_jet;
        TRefVector refVector;

        for (unsigned int j = 0; j < theJetCollection.size(); j++) {
            TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
            double DR = ElePs[i].DeltaR(JetP);

            if (DR > minDeltaR_)
        refVector.push_back(TRef(theJetCollectionHandle, j));
        }
    allSelections->push_back(refVector);
    }

    //iEvent.put(theFilteredJetCollection);
    iEvent.put(allSelections);
  
  return;
  
}

Member Data Documentation

template<typename T >
edm::InputTag HLTJetCollectionsForElePlusJets< T >::hltElectronTag [private]

Definition at line 49 of file HLTJetCollectionsForElePlusJets.h.

template<typename T >
double HLTJetCollectionsForElePlusJets< T >::minDeltaR_ [private]

Definition at line 56 of file HLTJetCollectionsForElePlusJets.h.

template<typename T >
edm::InputTag HLTJetCollectionsForElePlusJets< T >::sourceJetTag [private]

Definition at line 50 of file HLTJetCollectionsForElePlusJets.h.