CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecoTauDiscriminationByLeadingObjectPtCut.cc
Go to the documentation of this file.
2 
3 /*
4  * class PFRecoTauDiscriminationByLeadingObjectPtCut
5  * created : October 08 2008,
6  * revised : Wed Aug 19 17:13:04 PDT 2009
7  * Authors : Simone Gennai (SNS), Evan Friis (UC Davis)
8  */
9 
10 using namespace reco;
11 
13  public:
15  chargedOnly_ = iConfig.getParameter<bool>("UseOnlyChargedHadrons");
16  minPtLeadObject_ = iConfig.getParameter<double>("MinPtLeadingObject");
17  }
19  double discriminate(const PFTauRef& pfTau) override;
20  private:
23 };
24 
26 {
27  double leadObjectPt = -1.;
28  if( chargedOnly_ )
29  {
30  // consider only charged hadrons. note that the leadPFChargedHadrCand is the highest pt
31  // charged signal cone object above the quality cut level (typically 0.5 GeV).
32  if( thePFTauRef->leadPFChargedHadrCand().isNonnull() )
33  {
34  leadObjectPt = thePFTauRef->leadPFChargedHadrCand()->pt();
35  }
36  }
37  else
38  {
39  // If using the 'leading pion' option, require that:
40  // 1) at least one charged hadron exists above threshold (thePFTauRef->leadPFChargedHadrCand().isNonnull())
41  // 2) the lead PFCand exists. In the case that the highest pt charged hadron is above the PFRecoTauProducer threshold
42  // (typically 5 GeV), the leadPFCand and the leadPFChargedHadrCand are the same object. If the leadPFChargedHadrCand
43  // is below 5GeV, but there exists a neutral PF particle > 5 GeV, it is set to be the leadPFCand
44  if( thePFTauRef->leadPFCand().isNonnull() && thePFTauRef->leadPFChargedHadrCand().isNonnull() )
45  {
46  leadObjectPt = thePFTauRef->leadPFCand()->pt();
47  }
48  }
49 
50  return ( leadObjectPt > minPtLeadObject_ ? 1. : 0. );
51 }
52 
54 
55 /*
56 void PFRecoTauDiscriminationByLeadingPionPtCut::produce(edm::Event& iEvent,const edm::EventSetup& iEventSetup){
57  edm::Handle<PFTauCollection> thePFTauCollection;
58  iEvent.getByLabel(PFTauProducer_,thePFTauCollection);
59 
60 
61  auto_ptr<PFTauDiscriminator> thePFTauDiscriminatorByLeadingPionPtCut(new PFTauDiscriminator(PFTauRefProd(thePFTauCollection)));
62 
63  //loop over the PFTau candidates
64  for(size_t iPFTau=0;iPFTau<thePFTauCollection->size();++iPFTau) {
65  PFTauRef thePFTauRef(thePFTauCollection,iPFTau);
66  PFTau thePFTau=*thePFTauRef;
67  double theleadTrackPtCutDiscriminator = 0.;
68  // fill the AssociationVector object
69  if (!thePFTau.leadPFCand() || !thePFTau.leadPFChargedHadrCand())
70  {
71  theleadTrackPtCutDiscriminator=0.;
72  }
73  else if(thePFTau.leadPFCand()->pt() > minPtLeadTrack_) theleadTrackPtCutDiscriminator=1.;
74 
75  thePFTauDiscriminatorByLeadingPionPtCut->setValue(iPFTau,theleadTrackPtCutDiscriminator);
76  }
77 
78  iEvent.put(thePFTauDiscriminatorByLeadingPionPtCut);
79 
80 }
81 
82 */
T getParameter(std::string const &) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250