CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CaloRecoTauDiscriminationByDeltaE.cc
Go to the documentation of this file.
2 
4 
5 /* class CaloRecoTauDiscriminationByDeltaE
6  * created : September 23 2010,
7  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
8  * based on H+ tau ID by Lauri Wendland
9  */
10 
12 #include "TLorentzVector.h"
13 
14 using namespace reco;
15 using namespace std;
16 using namespace edm;
17 
19  public:
21  deltaEmin = iConfig.getParameter<double>("deltaEmin");
22  deltaEmax = iConfig.getParameter<double>("deltaEmax");
23  chargedPionMass = 0.139;
24  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
25  }
26 
28 
29  void beginEvent(const edm::Event&, const edm::EventSetup&) override;
30  double discriminate(const reco::CaloTauRef&) const override;
31 
32  private:
33  double DeltaE(const CaloTauRef&) const ;
34 
36 
37  double deltaEmin,deltaEmax;
39 };
40 
42 }
43 
45 
46  double dE = DeltaE(tau);
47  if(booleanOutput) return ( dE > deltaEmin && dE < deltaEmax ? 1. : 0. );
48  return dE;
49 }
50 
52  double tracksE = 0;
53  reco::TrackRefVector signalTracks = tau->signalTracks();
54  for(size_t i = 0; i < signalTracks.size(); ++i){
55  TLorentzVector p4;
56  p4.SetXYZM(signalTracks[i]->px(),
57  signalTracks[i]->py(),
58  signalTracks[i]->pz(),
59  chargedPionMass);
60  tracksE += p4.E();
61  }
62  if(tau->leadTrackHCAL3x3hitsEtSum() == 0) return -1; // electron
63  return tracksE/tau->leadTrackHCAL3x3hitsEtSum() - 1.0;
64 }
65 
67 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void beginEvent(const edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int iEvent
Definition: GenABIO.cc:230
double p4[4]
Definition: TauolaWrapper.h:92
CaloRecoTauDiscriminationByDeltaE(const ParameterSet &iConfig)
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
double discriminate(const reco::CaloTauRef &) const override