CMS 3D CMS Logo

L1HLTTauMatching.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
8 
9 //
10 // class decleration
11 //
12 using namespace reco;
13 using namespace std;
14 using namespace edm;
15 using namespace l1extra;
16 
18  jetSrc( consumes<PFTauCollection>(iConfig.getParameter<InputTag>("JetSrc") ) ),
19  tauTrigger( consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger") ) ),
20  mEt_Min( iConfig.getParameter<double>("EtMin") )
21 {
22  produces<PFTauCollection>();
23 }
25 
27 {
28 
29  using namespace edm;
30  using namespace std;
31  using namespace reco;
32  using namespace trigger;
33  using namespace l1extra;
34 
35  unique_ptr<PFTauCollection> tauL2jets(new PFTauCollection);
36 
37  double deltaR = 1.0;
38  double matchingR = 0.5;
39  //Getting HLT jets to be matched
41  iEvent.getByToken( jetSrc, tauJets );
42 
43  // std::cout <<"Size of input jet collection "<<tauJets->size()<<std::endl;
44 
46  iEvent.getByToken(tauTrigger,l1TriggeredTaus);
47 
48  vector<l1extra::L1JetParticleRef> tauCandRefVec;
49  vector<l1extra::L1JetParticleRef> jetCandRefVec;
50  l1TriggeredTaus->getObjects( trigger::TriggerL1TauJet,tauCandRefVec);
51  l1TriggeredTaus->getObjects( trigger::TriggerL1CenJet,jetCandRefVec);
52 
53  math::XYZPoint a(0.,0.,0.);
54 
55  for( unsigned int iL1Tau=0; iL1Tau <tauCandRefVec.size();iL1Tau++)
56  {
57  for(unsigned int iJet=0;iJet<tauJets->size();iJet++)
58  {
59  //Find the relative L2TauJets, to see if it has been reconstructed
60  const PFTau & myJet = (*tauJets)[iJet];
61  deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
62  if(deltaR < matchingR ) {
63  // LeafCandidate myLC(myJet);
64  if(myJet.leadPFChargedHadrCand().isNonnull()){
65  a = myJet.leadPFChargedHadrCand()->vertex();
66  }
67  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
68  if(myJet.pt() > mEt_Min) {
69  // tauL2LC->push_back(myLC);
70  tauL2jets->push_back(myPFTau);
71  }
72  break;
73  }
74  }
75  }
76 
77  for(unsigned int iL1Tau=0; iL1Tau <jetCandRefVec.size();iL1Tau++)
78  {
79  for(unsigned int iJet=0;iJet<tauJets->size();iJet++)
80  {
81  const PFTau & myJet = (*tauJets)[iJet];
82  //Find the relative L2TauJets, to see if it has been reconstructed
83  deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (jetCandRefVec[iL1Tau]->p4()).Vect());
84  if(deltaR < matchingR ) {
85  // LeafCandidate myLC(myJet);
86  if(myJet.leadPFChargedHadrCand().isNonnull()){
87  a = myJet.leadPFChargedHadrCand()->vertex();
88  }
89 
90  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(),a);
91  if(myJet.pt() > mEt_Min) {
92  // tauL2LC->push_back(myLC);
93  tauL2jets->push_back(myPFTau);
94  }
95  break;
96  }
97  }
98  }
99 
100 
101  //std::cout <<"Size of L1HLT matched jets "<<tauL2jets->size()<<std::endl;
102 
103  iEvent.put(std::move(tauL2jets));
104  // iEvent.put(std::move(tauL2LC));
105 }
106 
108 {
110  desc.add<edm::InputTag>("L1TauTrigger",edm::InputTag("hltL1sDoubleIsoTau40er"))->setComment("Name of trigger filter");
111  desc.add<edm::InputTag>("JetSrc",edm::InputTag("hltSelectedPFTausTrackPt1MediumIsolationReg"))->setComment("Input collection of PFTaus");
112  desc.add<double>("EtMin",0.0)->setComment("Minimal pT of PFTau to match");
113  descriptions.setComment("This module produces collection of PFTaus matched to L1ExtraTaus/Jets passing a HLT filter (Only p4 and vertex of returned PFTaus are set).");
114  descriptions.add("L1HLTJetsMatching",desc);
115 }
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
~L1HLTTauMatching() override
const PFCandidatePtr & leadPFChargedHadrCand() const
Definition: PFTau.cc:67
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
double pt() const final
transverse momentum
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
L1HLTTauMatching(const edm::ParameterSet &)
const Point & vertex() const override
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:656
int iEvent
Definition: GenABIO.cc:230
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tauTrigger
const edm::EDGetTokenT< reco::PFTauCollection > jetSrc
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void setComment(std::string const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
const double mEt_Min
def move(src, dest)
Definition: eostools.py:510