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 
5 
6 /* class CaloRecoTauDiscriminationByDeltaE
7  * created : September 23 2010,
8  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
9  * based on H+ tau ID by Lauri Wendland
10  */
11 
13 #include "TLorentzVector.h"
14 
15 using namespace reco;
16 using namespace std;
17 using namespace edm;
18 
20  public:
22  deltaEmin = iConfig.getParameter<double>("deltaEmin");
23  deltaEmax = iConfig.getParameter<double>("deltaEmax");
24  chargedPionMass = 0.139;
25  booleanOutput = iConfig.getParameter<bool>("BooleanOutput");
26  }
27 
29 
30  void beginEvent(const edm::Event&, const edm::EventSetup&);
31  double discriminate(const reco::CaloTauRef&);
32 
33  private:
34  double DeltaE(const CaloTauRef&);
35 
37 
38  double deltaEmin,deltaEmax;
40 };
41 
43 }
44 
46 
47  double dE = DeltaE(tau);
48  if(booleanOutput) return ( dE > deltaEmin && dE < deltaEmax ? 1. : 0. );
49  return dE;
50 }
51 
53  double tracksE = 0;
54  reco::TrackRefVector signalTracks = tau->signalTracks();
55  for(size_t i = 0; i < signalTracks.size(); ++i){
56  TLorentzVector p4;
57  p4.SetXYZM(signalTracks[i]->px(),
58  signalTracks[i]->py(),
59  signalTracks[i]->pz(),
60  chargedPionMass);
61  tracksE += p4.E();
62  }
63  if(tau->leadTrackHCAL3x3hitsEtSum() == 0) return -1; // electron
64  return tracksE/tau->leadTrackHCAL3x3hitsEtSum() - 1.0;
65 }
66 
68 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
int iEvent
Definition: GenABIO.cc:243
double p4[4]
Definition: TauolaWrapper.h:92
CaloRecoTauDiscriminationByDeltaE(const ParameterSet &iConfig)
DEFINE_FWK_MODULE(CosmicTrackingParticleSelector)
void beginEvent(const edm::Event &, const edm::EventSetup &)
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89