CMS 3D CMS Logo

L1THLTTauMatching.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
6 
7 //
8 // class decleration
9 //
10 using namespace reco ;
11 using namespace std ;
12 using namespace edm ;
13 using namespace trigger;
14 
16  jetSrc ( consumes<PFTauCollection> (iConfig.getParameter<InputTag>("JetSrc" ) ) ),
17  tauTrigger( consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<InputTag>("L1TauTrigger") ) ),
18  mEt_Min ( iConfig.getParameter<double> ("EtMin" ) )
19 {
20  produces<PFTauCollection>();
21 }
23 
25 {
26 
27  unique_ptr<PFTauCollection> tauL2jets(new PFTauCollection);
28 
29  double deltaR = 1.0;
30  double matchingR = 0.5;
31 
32  // Getting HLT jets to be matched
34  iEvent.getByToken( jetSrc, tauJets );
35 
37  iEvent.getByToken(tauTrigger,l1TriggeredTaus);
38 
39  l1t::TauVectorRef tauCandRefVec;
40  l1TriggeredTaus->getObjects( trigger::TriggerL1Tau,tauCandRefVec);
41 
42  math::XYZPoint a(0.,0.,0.);
43 
44  for(unsigned int iL1Tau = 0; iL1Tau < tauCandRefVec.size(); iL1Tau++){
45  for(unsigned int iJet = 0; iJet < tauJets->size(); iJet++){
46  // Find the relative L2TauJets, to see if it has been reconstructed
47  const PFTau & myJet = (*tauJets)[iJet];
48  deltaR = ROOT::Math::VectorUtil::DeltaR(myJet.p4().Vect(), (tauCandRefVec[iL1Tau]->p4()).Vect());
49  if(deltaR < matchingR ) {
50  if(myJet.leadPFChargedHadrCand().isNonnull()){
51  a = myJet.leadPFChargedHadrCand()->vertex();
52  }
53  PFTau myPFTau(std::numeric_limits<int>::quiet_NaN(), myJet.p4(), a);
54  if(myJet.pt() > mEt_Min) {
55  tauL2jets->push_back(myPFTau);
56  }
57  break;
58  }
59  }
60  }
61 
62  iEvent.put(std::move(tauL2jets));
63 }
64 
66 {
68  desc.add<edm::InputTag>("L1TauTrigger", edm::InputTag("hltL1sDoubleIsoTau40er" ))->setComment("Name of trigger filter" );
69  desc.add<edm::InputTag>("JetSrc" , edm::InputTag("hltSelectedPFTausTrackPt1MediumIsolationReg"))->setComment("Input collection of PFTaus");
70  desc.add<double> ("EtMin",0.0)->setComment("Minimal pT of PFTau to match");
71  descriptions.setComment("This module produces collection of PFTaus matched to L1 Taus / Jets passing a HLT filter (Only p4 and vertex of returned PFTaus are set).");
72  descriptions.add ("L1THLTTauMatching",desc);
73 }
const edm::EDGetTokenT< reco::PFTauCollection > jetSrc
L1THLTTauMatching(const edm::ParameterSet &)
virtual double pt() const final
transverse momentum
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const PFCandidatePtr & leadPFChargedHadrCand() const
Definition: PFTau.cc:67
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
int iEvent
Definition: GenABIO.cc:230
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tauTrigger
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:652
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:169
const double mEt_Min
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void setComment(std::string const &value)
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
std::vector< TauRef > TauVectorRef
Definition: Tau.h:14
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
def move(src, dest)
Definition: eostools.py:510