CMS 3D CMS Logo

Public Member Functions | Private Types | Private Attributes

L2TauJetsProvider Class Reference

#include <L2TauJetsProvider.h>

Inheritance diagram for L2TauJetsProvider:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

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

Private Types

typedef std::vector
< edm::InputTag
vtag

Private Attributes

std::vector
< l1extra::L1JetParticleRef
jetCandRefVec
vtag jetSrc
edm::InputTag l1ParticlesJet
edm::InputTag l1ParticlesTau
double mEt_Min
std::map< int, const
reco::CaloJet
myL2L1JetsMap
std::vector
< l1extra::L1JetParticleRef
objL1CandRefVec
l1extra::L1JetParticleRef tauCandRef
std::vector
< l1extra::L1JetParticleRef
tauCandRefVec
edm::InputTag tauTrigger

Detailed Description

Definition at line 23 of file L2TauJetsProvider.h.


Member Typedef Documentation

typedef std::vector<edm::InputTag> L2TauJetsProvider::vtag [private]

Definition at line 35 of file L2TauJetsProvider.h.


Constructor & Destructor Documentation

L2TauJetsProvider::L2TauJetsProvider ( const edm::ParameterSet iConfig) [explicit]

Definition at line 15 of file L2TauJetsProvider.cc.

References edm::ParameterSet::getParameter(), and PatBasicFWLiteJetAnalyzer_Selector_cfg::jetSrc.

{
  jetSrc = iConfig.getParameter<vtag>("JetSrc");
  l1ParticlesTau = iConfig.getParameter<InputTag>("L1ParticlesTau");
  l1ParticlesJet = iConfig.getParameter<InputTag>("L1ParticlesJet");
  tauTrigger = iConfig.getParameter<InputTag>("L1TauTrigger");
  mEt_Min = iConfig.getParameter<double>("EtMin");
  
  produces<CaloJetCollection>();
}
L2TauJetsProvider::~L2TauJetsProvider ( )

Definition at line 26 of file L2TauJetsProvider.cc.

{ }

Member Function Documentation

void L2TauJetsProvider::produce ( edm::Event iEvent,
const edm::EventSetup iES 
) [virtual]

Implements edm::EDProducer.

Definition at line 28 of file L2TauJetsProvider.cc.

References edm::HandleBase::clear(), reco::deltaR(), edm::Event::getByLabel(), i, j, PatBasicFWLiteJetAnalyzer_Selector_cfg::jetSrc, p4, reco::LeafCandidate::pt(), edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, alignCSCRings::s, trigger::TriggerL1CenJet, and trigger::TriggerL1TauJet.

{

 using namespace edm;
 using namespace std;
 using namespace reco;
 using namespace trigger;
 using namespace l1extra;


 //Getting all the L1Seeds

 
 //Getting the Collections of L2ReconstructedJets from L1Seeds
 //and removing the collinear jets
 myL2L1JetsMap.clear();
 int iL1Jet = 0;
 for( vtag::const_iterator s = jetSrc.begin(); s != jetSrc.end(); ++ s ) {
   edm::Handle<CaloJetCollection> tauJets;
   iEvent.getByLabel( * s, tauJets );
   CaloJetCollection::const_iterator iTau = tauJets->begin();
   if(iTau != tauJets->end()){
     //Create a Map to associate to every Jet its L1SeedId, i.e. 0,1,2 or 3
     myL2L1JetsMap.insert(std::pair<int, const CaloJet>(iL1Jet, *(iTau)));
   }
   iL1Jet++;
 }
 std::auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection);

 //Loop over the jetSrc to select the proper jets
  
  double deltaR = 0.1;
  double matchingR = 0.01;
  //Loop over the Map to find which jets has fired the trigger
  //myL1Tau is the Collection of L1TauCandidates (from 0 to max  4 elements)
  //get the list of trigger candidates from the HLTL1SeedGT filter
  edm::Handle< L1JetParticleCollection > tauColl ; 
  edm::Handle< L1JetParticleCollection > jetColl ; 
  
  iEvent.getByLabel( l1ParticlesTau, tauColl );
  iEvent.getByLabel(l1ParticlesJet, jetColl);

  const L1JetParticleCollection  myL1Tau = *(tauColl.product());

  const L1JetParticleCollection  myL1Jet = *(jetColl.product());  

    L1JetParticleCollection myL1Obj;
    myL1Obj.reserve(8);
    
    for(unsigned int i=0;i<myL1Tau.size();i++)
      {
        myL1Obj.push_back(myL1Tau[i]);
      }
    for(unsigned int j=0;j<myL1Jet.size();j++)
      {
        myL1Obj.push_back(myL1Jet[j]);
      }


    edm::Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus;
    if(iEvent.getByLabel(tauTrigger,l1TriggeredTaus)){
    
      
      tauCandRefVec.clear();
      jetCandRefVec.clear();
      
      l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec);
      l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec);
      
      for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++)
        {  
          for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++)
            {
              //Find the relative L2TauJets, to see if it has been reconstructed
            std::map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet);
            if(myL2itr!=myL2L1JetsMap.end()){
              //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger
              if(&tauCandRefVec[iL1Tau]) 
                deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
              if(deltaR < matchingR ) {
              //              Getting back from the map the L2TauJet
                const CaloJet myL2TauJet = myL2itr->second;
                if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet);
                myL2L1JetsMap.erase(myL2itr->first);
                break;
                
              }
            }
            
          }
      }
      
    for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++)
      {  
        for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++)
          {
            //Find the relative L2TauJets, to see if it has been reconstructed
            std::map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet);
            if(myL2itr!=myL2L1JetsMap.end()){
              //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger
              if(&jetCandRefVec[iL1Tau])
                deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
              if(deltaR < matchingR) {
                // Getting back from the map the L2TauJet
                const CaloJet myL2TauJet = myL2itr->second;
                
                if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet);
                myL2L1JetsMap.erase(myL2itr->first);
                break;
                
              }
            }
            
          }
      }
    
  }
  //  std::cout <<"Size of L2 jets "<<tauL2jets->size()<<std::endl;

  iEvent.put(tauL2jets);

}

Member Data Documentation

Definition at line 31 of file L2TauJetsProvider.h.

Definition at line 36 of file L2TauJetsProvider.h.

Definition at line 38 of file L2TauJetsProvider.h.

Definition at line 37 of file L2TauJetsProvider.h.

double L2TauJetsProvider::mEt_Min [private]

Definition at line 40 of file L2TauJetsProvider.h.

std::map<int, const reco::CaloJet> L2TauJetsProvider::myL2L1JetsMap [private]

Definition at line 41 of file L2TauJetsProvider.h.

Definition at line 32 of file L2TauJetsProvider.h.

Definition at line 33 of file L2TauJetsProvider.h.

Definition at line 30 of file L2TauJetsProvider.h.

Definition at line 39 of file L2TauJetsProvider.h.