![]() |
![]() |
#include <RecoTauTag/HLTProducers/interface/L2TauJetsProvider.h>
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 |
Definition at line 23 of file L2TauJetsProvider.h.
typedef std::vector<edm::InputTag> L2TauJetsProvider::vtag [private] |
Definition at line 35 of file L2TauJetsProvider.h.
L2TauJetsProvider::L2TauJetsProvider | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 15 of file L2TauJetsProvider.cc.
References edm::ParameterSet::getParameter(), jetSrc, l1ParticlesJet, l1ParticlesTau, mEt_Min, and tauTrigger.
00016 { 00017 jetSrc = iConfig.getParameter<vtag>("JetSrc"); 00018 l1ParticlesTau = iConfig.getParameter<InputTag>("L1ParticlesTau"); 00019 l1ParticlesJet = iConfig.getParameter<InputTag>("L1ParticlesJet"); 00020 tauTrigger = iConfig.getParameter<InputTag>("L1TauTrigger"); 00021 mEt_Min = iConfig.getParameter<double>("EtMin"); 00022 00023 produces<CaloJetCollection>(); 00024 }
L2TauJetsProvider::~L2TauJetsProvider | ( | ) |
void L2TauJetsProvider::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iES | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 28 of file L2TauJetsProvider.cc.
References deltaR(), edm::Event::getByLabel(), i, j, jetCandRefVec, jetColl, jetSrc, l1ParticlesJet, l1ParticlesTau, mEt_Min, myL2L1JetsMap, p4, edm::Handle< T >::product(), reco::Particle::pt(), edm::Event::put(), HcalSimpleRecAlgoImpl::reco(), s, std, tauCandRefVec, tauTrigger, trigger::TriggerL1CenJet, and trigger::TriggerL1TauJet.
00029 { 00030 00031 using namespace edm; 00032 using namespace std; 00033 using namespace reco; 00034 using namespace trigger; 00035 using namespace l1extra; 00036 00037 00038 //Getting all the L1Seeds 00039 00040 00041 //Getting the Collections of L2ReconstructedJets from L1Seeds 00042 //and removing the collinear jets 00043 myL2L1JetsMap.clear(); 00044 int iL1Jet = 0; 00045 for( vtag::const_iterator s = jetSrc.begin(); s != jetSrc.end(); ++ s ) { 00046 edm::Handle<CaloJetCollection> tauJets; 00047 iEvent.getByLabel( * s, tauJets ); 00048 CaloJetCollection::const_iterator iTau = tauJets->begin(); 00049 if(iTau != tauJets->end()){ 00050 //Create a Map to associate to every Jet its L1SeedId, i.e. 0,1,2 or 3 00051 myL2L1JetsMap.insert(pair<int, const CaloJet>(iL1Jet, *(iTau))); 00052 } 00053 iL1Jet++; 00054 } 00055 auto_ptr<CaloJetCollection> tauL2jets(new CaloJetCollection); 00056 00057 //Loop over the jetSrc to select the proper jets 00058 00059 double deltaR = 0.1; 00060 double matchingR = 0.01; 00061 //Loop over the Map to find which jets has fired the trigger 00062 //myL1Tau is the Collection of L1TauCandidates (from 0 to max 4 elements) 00063 //get the list of trigger candidates from the HLTL1SeedGT filter 00064 Handle< L1JetParticleCollection > tauColl ; 00065 Handle< L1JetParticleCollection > jetColl ; 00066 00067 iEvent.getByLabel( l1ParticlesTau, tauColl ); 00068 iEvent.getByLabel(l1ParticlesJet, jetColl); 00069 00070 const L1JetParticleCollection myL1Tau = *(tauColl.product()); 00071 00072 const L1JetParticleCollection myL1Jet = *(jetColl.product()); 00073 00074 L1JetParticleCollection myL1Obj; 00075 myL1Obj.reserve(8); 00076 00077 for(unsigned int i=0;i<myL1Tau.size();i++) 00078 { 00079 myL1Obj.push_back(myL1Tau[i]); 00080 } 00081 for(unsigned int j=0;j<myL1Jet.size();j++) 00082 { 00083 myL1Obj.push_back(myL1Jet[j]); 00084 } 00085 00086 00087 Handle<trigger::TriggerFilterObjectWithRefs> l1TriggeredTaus; 00088 if(iEvent.getByLabel(tauTrigger,l1TriggeredTaus)){ 00089 00090 00091 tauCandRefVec.clear(); 00092 jetCandRefVec.clear(); 00093 00094 l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec); 00095 l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec); 00096 00097 for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++) 00098 { 00099 for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++) 00100 { 00101 //Find the relative L2TauJets, to see if it has been reconstructed 00102 map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet); 00103 if(myL2itr!=myL2L1JetsMap.end()){ 00104 //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger 00105 if(&tauCandRefVec[iL1Tau]) 00106 deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect()); 00107 if(deltaR < matchingR ) { 00108 // Getting back from the map the L2TauJet 00109 const CaloJet myL2TauJet = myL2itr->second; 00110 if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet); 00111 myL2L1JetsMap.erase(myL2itr->first); 00112 break; 00113 00114 } 00115 } 00116 00117 } 00118 } 00119 00120 for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++) 00121 { 00122 for(unsigned int iJet=0;iJet<myL1Obj.size();iJet++) 00123 { 00124 //Find the relative L2TauJets, to see if it has been reconstructed 00125 map<int, const reco::CaloJet>::const_iterator myL2itr = myL2L1JetsMap.find(iJet); 00126 if(myL2itr!=myL2L1JetsMap.end()){ 00127 //Calculate the DeltaR between L1TauCandidate and L1Tau which fired the trigger 00128 if(&jetCandRefVec[iL1Tau]) 00129 deltaR = ROOT::Math::VectorUtil::DeltaR(myL1Obj[iJet].p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect()); 00130 if(deltaR < matchingR) { 00131 // Getting back from the map the L2TauJet 00132 const CaloJet myL2TauJet = myL2itr->second; 00133 00134 if(myL2TauJet.pt() > mEt_Min) tauL2jets->push_back(myL2TauJet); 00135 myL2L1JetsMap.erase(myL2itr->first); 00136 break; 00137 00138 } 00139 } 00140 00141 } 00142 } 00143 00144 } 00145 // cout <<"Size of L2 jets "<<tauL2jets->size()<<endl; 00146 00147 iEvent.put(tauL2jets); 00148 00149 }
std::vector<l1extra::L1JetParticleRef> L2TauJetsProvider::jetCandRefVec [private] |
vtag L2TauJetsProvider::jetSrc [private] |
Definition at line 36 of file L2TauJetsProvider.h.
Referenced by L2TauJetsProvider(), and produce().
Definition at line 38 of file L2TauJetsProvider.h.
Referenced by L2TauJetsProvider(), and produce().
Definition at line 37 of file L2TauJetsProvider.h.
Referenced by L2TauJetsProvider(), and produce().
double L2TauJetsProvider::mEt_Min [private] |
Definition at line 40 of file L2TauJetsProvider.h.
Referenced by L2TauJetsProvider(), and produce().
std::map<int, const reco::CaloJet> L2TauJetsProvider::myL2L1JetsMap [private] |
std::vector<l1extra::L1JetParticleRef> L2TauJetsProvider::objL1CandRefVec [private] |
Definition at line 32 of file L2TauJetsProvider.h.
Definition at line 33 of file L2TauJetsProvider.h.
std::vector<l1extra::L1JetParticleRef> L2TauJetsProvider::tauCandRefVec [private] |
edm::InputTag L2TauJetsProvider::tauTrigger [private] |
Definition at line 39 of file L2TauJetsProvider.h.
Referenced by L2TauJetsProvider(), and produce().