CMS 3D CMS Logo

Public Member Functions | Private Attributes

L1HLTTauMatching Class Reference

#include <L1HLTTauMatching.h>

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

List of all members.

Public Member Functions

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

Private Attributes

std::vector
< l1extra::L1JetParticleRef
jetCandRefVec
edm::InputTag jetSrc
double mEt_Min
std::vector
< l1extra::L1JetParticleRef
objL1CandRefVec
l1extra::L1JetParticleRef tauCandRef
std::vector
< l1extra::L1JetParticleRef
tauCandRefVec
edm::InputTag tauTrigger

Detailed Description

Definition at line 21 of file L1HLTTauMatching.h.


Constructor & Destructor Documentation

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

Definition at line 17 of file L1HLTTauMatching.cc.

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

{
  jetSrc = iConfig.getParameter<InputTag>("JetSrc");
  tauTrigger = iConfig.getParameter<InputTag>("L1TauTrigger");
  mEt_Min = iConfig.getParameter<double>("EtMin");
  
  produces<PFTauCollection>();
}
L1HLTTauMatching::~L1HLTTauMatching ( )

Definition at line 25 of file L1HLTTauMatching.cc.

{ }

Member Function Documentation

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

Implements edm::EDProducer.

Definition at line 27 of file L1HLTTauMatching.cc.

References a, reco::deltaR(), f, edm::Event::getByLabel(), edm::Ref< C, T, F >::isNonnull(), PatBasicFWLiteJetAnalyzer_Selector_cfg::jetSrc, reco::PFTau::leadPFChargedHadrCand(), reco::LeafCandidate::p4(), reco::LeafCandidate::pt(), edm::Event::put(), dt_dqm_sourceclient_common_cff::reco, trigger::TriggerL1CenJet, and trigger::TriggerL1TauJet.

{
        
        using namespace edm;
        using namespace std;
        using namespace reco;
        using namespace trigger;
        using namespace l1extra;
        
          auto_ptr<PFTauCollection> tauL2jets(new PFTauCollection);
        
        double deltaR = 1.0;
        double matchingR = 0.5;
        //Getting HLT jets to be matched
        edm::Handle<PFTauCollection > tauJets;
        iEvent.getByLabel( jetSrc, tauJets );

//              std::cout <<"Size of input jet collection "<<tauJets->size()<<std::endl;
                
        edm::Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus;
    iEvent.getByLabel(tauTrigger,l1TriggeredTaus);
                
                
                tauCandRefVec.clear();
                jetCandRefVec.clear();
                
                l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec);
                l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec);
                math::XYZPoint a(0.,0.,0.);
                CaloJet::Specific f;
                
                for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++)
                {  
                        for(unsigned int iJet=0;iJet<tauJets->size();iJet++)
                        {
                                //Find the relative L2TauJets, to see if it has been reconstructed
                                const PFTau &  myJet = (*tauJets)[iJet];
                                deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
                                if(deltaR < matchingR ) {
                                        //               LeafCandidate myLC(myJet);
                    if(myJet.leadPFChargedHadrCand().isNonnull()){
                        a =  myJet.leadPFChargedHadrCand()->vertex();  
                    }
                                  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
                                        if(myJet.pt() > mEt_Min) {
                                                //                tauL2LC->push_back(myLC);
                                                tauL2jets->push_back(myPFTau);
                                        }
                                        break;
                                }
                        }
                }  
                
                for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++)
                {  
                        for(unsigned int iJet=0;iJet<tauJets->size();iJet++)
                        {
                                const PFTau &  myJet = (*tauJets)[iJet];
                                //Find the relative L2TauJets, to see if it has been reconstructed
                                deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
                                if(deltaR < matchingR ) {
                                        //               LeafCandidate myLC(myJet);
                    if(myJet.leadPFChargedHadrCand().isNonnull()){
                        a =  myJet.leadPFChargedHadrCand()->vertex();  
                    }
                    
                                  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(),a);
                                        if(myJet.pt() > mEt_Min) {
                                                //tauL2LC->push_back(myLC);
                                                tauL2jets->push_back(myPFTau);
                                        }
                                        break;
                                }
                        }
                }
        

//std::cout <<"Size of L1HLT matched jets "<<tauL2jets->size()<<std::endl;

iEvent.put(tauL2jets);
// iEvent.put(tauL2LC);
}

Member Data Documentation

Definition at line 29 of file L1HLTTauMatching.h.

Definition at line 33 of file L1HLTTauMatching.h.

double L1HLTTauMatching::mEt_Min [private]

Definition at line 35 of file L1HLTTauMatching.h.

Definition at line 30 of file L1HLTTauMatching.h.

Definition at line 31 of file L1HLTTauMatching.h.

Definition at line 28 of file L1HLTTauMatching.h.

Definition at line 34 of file L1HLTTauMatching.h.