CMS 3D CMS Logo

PFRecoTauDiscriminationByNProngs.cc
Go to the documentation of this file.
6 
7 /* class PFRecoTauDiscriminationByNProngs
8  * created : August 30 2010,
9  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
10  * based on H+ tau ID by Lauri Wendland
11  * Modified April 16 2014 by S.Lehti
12  */
13 
14 using namespace reco;
15 using namespace std;
16 using namespace edm;
17 
19  public:
22 
23  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
24  double discriminate(const reco::PFTauRef&) const override;
25 
26  private:
27  std::unique_ptr<tau::RecoTauQualityCuts> qcuts_;
28  std::unique_ptr<tau::RecoTauVertexAssociator> vertexAssociator_;
29 
30  uint32_t minN,maxN;
33 };
34 
37  qualityCuts(iConfig.getParameterSet("qualityCuts"))
38 {
39  minN = iConfig.getParameter<uint32_t>("MinN");
40  maxN = iConfig.getParameter<uint32_t>("MaxN");
41  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
42 
43  qcuts_.reset(new tau::RecoTauQualityCuts(qualityCuts.getParameterSet("signalQualityCuts")));
44  vertexAssociator_.reset(new tau::RecoTauVertexAssociator(qualityCuts,consumesCollector()));
45 }
46 
48  vertexAssociator_->setEvent(iEvent);
49 }
50 
52 
53  reco::VertexRef pv = vertexAssociator_->associatedVertex(*tau);
54  const PFCandidatePtr leadingTrack = tau->leadPFChargedHadrCand();
55 
56  uint np = 0;
57  if(leadingTrack.isNonnull() && pv.isNonnull()){
58  qcuts_->setPV(pv);
59  qcuts_->setLeadTrack(tau->leadPFChargedHadrCand());
60 
61  for(auto const& cand : tau->signalPFChargedHadrCands() ) {
62  if ( qcuts_->filterCandRef(cand) ) np++;
63  }
64  }
65 
66  bool accepted = false;
67  if(maxN == 0){
68  if(np == 1 || np == 3) accepted = true;
69  }else{
70  if(np >= minN && np <= maxN) accepted = true;
71  }
72 
73  if(!accepted) np = 0;
74  if(booleanOutput) return accepted;
75  return np;
76 }
77 
T getParameter(std::string const &) const
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:251
double discriminate(const reco::PFTauRef &) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::unique_ptr< tau::RecoTauVertexAssociator > vertexAssociator_
int iEvent
Definition: GenABIO.cc:230
int np
Definition: AMPTWrapper.h:33
def pv(vc)
Definition: MetAnalyzer.py:7
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
void beginEvent(const edm::Event &, const edm::EventSetup &) override
ParameterSet const & getParameterSet(std::string const &) const
def uint(string)
fixed size matrix
HLT enums.
std::unique_ptr< tau::RecoTauQualityCuts > qcuts_