CMS 3D CMS Logo

HLTPFTauPairLeadTrackDzMatchFilter.cc
Go to the documentation of this file.
1 
3 
11 
12 
14 
15  tauSrc_ = conf.getParameter<edm::InputTag>("tauSrc");
16  tauSrcToken_ = consumes<reco::PFTauCollection>(tauSrc_);
17  tauMinPt_ = conf.getParameter<double>("tauMinPt");
18  tauMaxEta_ = conf.getParameter<double>("tauMaxEta");
19  tauMinDR_ = conf.getParameter<double>("tauMinDR");
20  tauLeadTrackMaxDZ_ = conf.getParameter<double>("tauLeadTrackMaxDZ");
21  triggerType_ = conf.getParameter<int>("triggerType");
22 
23  // set the minimum DR between taus, so that one never has a chance
24  // to create a pair out of the same Tau replicated with two
25  // different vertices
26  if (tauMinDR_ < 0.1) tauMinDR_ = 0.1;
27 
28 }
29 
31 
33 
36 
37  desc.add<edm::InputTag>("tauSrc",edm::InputTag("hltPFTaus") );
38  desc.add<double>("tauMinPt",5.0);
39  desc.add<double>("tauMaxEta",2.5);
40  desc.add<double>("tauMinDR",0.1);
41  desc.add<double>("tauLeadTrackMaxDZ",0.2);
42  desc.add<int>("triggerType",trigger::TriggerTau);
43  descriptions.add("hltPFTauPairLeadTrackDzMatchFilter",desc);
44 }
45 
47 
48  using namespace std;
49  using namespace reco;
50 
51  // The resuilting filter object to store in the Event
52  if(saveTags() ) filterproduct.addCollectionTag(tauSrc_);
53 
54  // Ref to Candidate object to be recorded in the filter object
55  PFTauRef ref;
56 
57  // Pick up taus
59  ev.getByToken( tauSrcToken_, tausHandle );
60  const PFTauCollection & taus = *tausHandle;
61  const size_t n_taus = taus.size();
62 
63  // Combine taus into pairs and check the dz matching
64  size_t npairs = 0, nfail_dz = 0;
65  if(n_taus > 1) for(size_t t1 = 0; t1 < n_taus; ++t1) {
66  if( taus[t1].leadPFChargedHadrCand().isNull() ||
67  taus[t1].leadPFChargedHadrCand()->trackRef().isNull() ||
68  taus[t1].pt() < tauMinPt_ ||
69  std::abs(taus[t1].eta() ) > tauMaxEta_ )
70  continue;
71 
72  float mindz = 99.f;
73  for(size_t t2 = t1+1; t2 < n_taus; ++t2){
74  if( taus[t2].leadPFChargedHadrCand().isNull() ||
75  taus[t2].leadPFChargedHadrCand()->trackRef().isNull() ||
76  taus[t2].pt() < tauMinPt_ ||
77  std::abs(taus[t2].eta() ) > tauMaxEta_ )
78  continue;
79 
80  float dr2 = reco::deltaR2(taus[t1].eta(), taus[t1].phi(),
81  taus[t2].eta(), taus[t2].phi() );
82  float dz = ( taus[t1].leadPFChargedHadrCand()->trackRef()->vz() -
83  taus[t2].leadPFChargedHadrCand()->trackRef()->vz() );
84 
85  // skip pairs of taus that are close
86  if ( dr2 < tauMinDR_*tauMinDR_ ) {
87  continue;
88  }
89 
90  if (std::abs(dz) < std::abs(mindz)) mindz = dz;
91 
92  // do not form a pair if dz is too large
93  if ( std::abs(dz) > tauLeadTrackMaxDZ_ ) {
94  ++nfail_dz;
95  continue;
96  }
97 
98  // add references to both jets
99  ref = PFTauRef(tausHandle, t1);
100  filterproduct.addObject(triggerType_, ref);
101 
102  ref = PFTauRef(tausHandle, t2);
103  filterproduct.addObject(triggerType_, ref);
104 
105  ++npairs;
106 
107  }
108 
109  }
110 
111  // return truth if at least one good pair found
112  return npairs>0;
113 
114 }
T getParameter(std::string const &) const
std::vector< PFTau > PFTauCollection
collection of PFTau objects
Definition: PFTauFwd.h:9
bool hltFilter(edm::Event &ev, const edm::EventSetup &es, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HLTPFTauPairLeadTrackDzMatchFilter(const edm::ParameterSet &conf)
bool ev
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Ref< PFTauCollection > PFTauRef
presistent reference to a PFTau
Definition: PFTauFwd.h:13
ParameterDescriptionBase * add(U const &iLabel, T const &value)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::PFTauCollection > tauSrcToken_
fixed size matrix
bool saveTags() const
Definition: HLTFilter.h:45