CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoTauTag/RecoTau/plugins/CaloRecoTauDiscriminationByDeltaE.cc

Go to the documentation of this file.
00001 #include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
00002 
00003 #include "RecoTauTag/TauTagTools/interface/PFTauQualityCutWrapper.h"
00004 #include "FWCore/Utilities/interface/InputTag.h"
00005 
00006 /* class CaloRecoTauDiscriminationByDeltaE
00007  * created : September 23 2010,
00008  * contributors : Sami Lehti (sami.lehti@cern.ch ; HIP, Helsinki)
00009  * based on H+ tau ID by Lauri Wendland
00010  */
00011 
00012 #include "DataFormats/TrackReco/interface/Track.h"
00013 #include "TLorentzVector.h"
00014 
00015 using namespace reco;
00016 using namespace std;
00017 using namespace edm;
00018 
00019 class CaloRecoTauDiscriminationByDeltaE : public CaloTauDiscriminationProducerBase  {
00020     public:
00021         explicit CaloRecoTauDiscriminationByDeltaE(const ParameterSet& iConfig):CaloTauDiscriminationProducerBase(iConfig){
00022                 deltaEmin               = iConfig.getParameter<double>("deltaEmin");
00023                 deltaEmax               = iConfig.getParameter<double>("deltaEmax");
00024                 chargedPionMass         = 0.139;
00025                 booleanOutput           = iConfig.getParameter<bool>("BooleanOutput");
00026         }
00027 
00028         ~CaloRecoTauDiscriminationByDeltaE(){}
00029 
00030         void beginEvent(const edm::Event&, const edm::EventSetup&);
00031         double discriminate(const reco::CaloTauRef&);
00032 
00033     private:
00034         double DeltaE(const CaloTauRef&);
00035 
00036         double chargedPionMass;
00037 
00038         double deltaEmin,deltaEmax;
00039         bool booleanOutput;
00040 };
00041 
00042 void CaloRecoTauDiscriminationByDeltaE::beginEvent(const Event& iEvent, const EventSetup& iSetup){
00043 }
00044 
00045 double CaloRecoTauDiscriminationByDeltaE::discriminate(const CaloTauRef& tau){
00046 
00047         double dE = DeltaE(tau);
00048         if(booleanOutput) return ( dE > deltaEmin && dE < deltaEmax ? 1. : 0. );
00049         return dE;
00050 }
00051 
00052 double CaloRecoTauDiscriminationByDeltaE::DeltaE(const CaloTauRef& tau){
00053         double tracksE = 0;
00054         reco::TrackRefVector signalTracks = tau->signalTracks();
00055         for(size_t i = 0; i < signalTracks.size(); ++i){
00056                 TLorentzVector p4;
00057                 p4.SetXYZM(signalTracks[i]->px(),
00058                signalTracks[i]->py(),
00059                signalTracks[i]->pz(),
00060                chargedPionMass);
00061                 tracksE += p4.E();
00062         }
00063         if(tau->leadTrackHCAL3x3hitsEtSum() == 0) return -1; // electron
00064         return tracksE/tau->leadTrackHCAL3x3hitsEtSum() - 1.0;
00065 }
00066 
00067 DEFINE_FWK_MODULE(CaloRecoTauDiscriminationByDeltaE);
00068