CMS 3D CMS Logo

Classes | Typedefs | Functions

/data/refman/pasoursint/CMSSW_5_3_6/src/RecoTauTag/RecoTau/plugins/TauDiscriminationAgainstCaloMuon.cc File Reference

#include "RecoTauTag/RecoTau/interface/TauDiscriminationProducerBase.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/MuonReco/interface/Muon.h"
#include "DataFormats/MuonReco/interface/MuonFwd.h"
#include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
#include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
#include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
#include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/CaloTowers/interface/CaloTower.h"
#include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
#include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "TrackingTools/Records/interface/TransientTrackRecord.h"
#include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToPoint.h"
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
#include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
#include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
#include "Geometry/Records/interface/CaloGeometryRecord.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "DataFormats/Math/interface/deltaR.h"
#include <TVector3.h>
#include <TMath.h>
#include <string>
#include "FWCore/Framework/interface/MakerMacros.h"

Go to the source code of this file.

Classes

class  TauDiscriminationAgainstCaloMuon< TauType, TauDiscriminator >
class  TauLeadTrackExtractor< T >
class  TauLeadTrackExtractor< reco::PFTau >

Typedefs

typedef
TauDiscriminationAgainstCaloMuon
< CaloTau,
CaloTauDiscriminator
CaloRecoTauDiscriminationAgainstCaloMuon
typedef
TauDiscriminationAgainstCaloMuon
< PFTau, PFTauDiscriminator
PFRecoTauDiscriminationAgainstCaloMuon

Functions

double compEcalEnergySum (const EcalRecHitCollection &ecalRecHits, const CaloSubdetectorGeometry *detGeometry, const reco::TransientTrack &transientTrack, double dR, const GlobalPoint &eventVertexPosition)
double compHcalEnergySum (const HBHERecHitCollection &hcalRecHits, const CaloSubdetectorGeometry *hbGeometry, const CaloSubdetectorGeometry *heGeometry, const reco::TransientTrack &transientTrack, double dR, const GlobalPoint &eventVertexPosition)
 DEFINE_FWK_MODULE (PFRecoTauDiscriminationAgainstCaloMuon)

Typedef Documentation

Definition at line 315 of file TauDiscriminationAgainstCaloMuon.cc.

Definition at line 314 of file TauDiscriminationAgainstCaloMuon.cc.


Function Documentation

double compEcalEnergySum ( const EcalRecHitCollection ecalRecHits,
const CaloSubdetectorGeometry detGeometry,
const reco::TransientTrack transientTrack,
double  dR,
const GlobalPoint eventVertexPosition 
)

Definition at line 185 of file TauDiscriminationAgainstCaloMuon.cc.

References edm::SortedCollection< T, SORT >::begin(), deltaR(), dir, CaloRecHits_cff::ecalRecHit, edm::SortedCollection< T, SORT >::end(), PV3DBase< T, VectorTag, FrameTag >::eta(), reco::TrackBase::eta(), CaloSubdetectorGeometry::getGeometry(), CaloCellGeometry::getPosition(), PV3DBase< T, VectorTag, FrameTag >::phi(), reco::TrackBase::phi(), TrajectoryStateClosestToPoint::position(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), reco::TransientTrack::track(), reco::TransientTrack::trajectoryStateClosestToPoint(), PV3DBase< T, VectorTag, FrameTag >::x(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, VectorTag, FrameTag >::y(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and PV3DBase< T, VectorTag, FrameTag >::z().

Referenced by TauDiscriminationAgainstCaloMuon< TauType, TauDiscriminator >::discriminate().

{
  double ecalEnergySum = 0.;
  for ( EcalRecHitCollection::const_iterator ecalRecHit = ecalRecHits.begin();
        ecalRecHit != ecalRecHits.end(); ++ecalRecHit ) {
    const CaloCellGeometry* cellGeometry = detGeometry->getGeometry(ecalRecHit->detid());
    
    if ( !cellGeometry ) {
      edm::LogError ("compEcalEnergySum") 
        << " Failed to access ECAL geometry for detId = " << ecalRecHit->detid().rawId()
        << " --> skipping !!";
      continue;
    }

    const GlobalPoint& cellPosition = cellGeometry->getPosition();

//--- CV: speed up computation by requiring eta-phi distance
//        between cell position and track direction to be dR < 0.5
    Vector3DBase<float, GlobalTag> cellPositionRelVertex = (cellPosition) - eventVertexPosition;
    if ( deltaR(cellPositionRelVertex.eta(), cellPositionRelVertex.phi(), 
                transientTrack.track().eta(), transientTrack.track().phi()) > 0.5 ) continue;

    TrajectoryStateClosestToPoint dcaPosition = transientTrack.trajectoryStateClosestToPoint(cellPosition);
    
    Vector3DBase<float, GlobalTag> d = (cellPosition - dcaPosition.position());

    TVector3 d3(d.x(), d.y(), d.z());
    TVector3 dir(transientTrack.track().px(), transientTrack.track().py(), transientTrack.track().pz());

    double dPerp = d3.Cross(dir.Unit()).Mag();
    double dParl = TVector3(cellPosition.x(), cellPosition.y(), cellPosition.z()).Dot(dir.Unit());
    
    if ( dPerp < dR && dParl > 100. ) {
      ecalEnergySum += ecalRecHit->energy();
    }
  }
  
  return ecalEnergySum;
}
double compHcalEnergySum ( const HBHERecHitCollection hcalRecHits,
const CaloSubdetectorGeometry hbGeometry,
const CaloSubdetectorGeometry heGeometry,
const reco::TransientTrack transientTrack,
double  dR,
const GlobalPoint eventVertexPosition 
)

Definition at line 228 of file TauDiscriminationAgainstCaloMuon.cc.

References edm::SortedCollection< T, SORT >::begin(), deltaR(), dir, edm::SortedCollection< T, SORT >::end(), PV3DBase< T, VectorTag, FrameTag >::eta(), reco::TrackBase::eta(), CaloSubdetectorGeometry::getGeometry(), CaloCellGeometry::getPosition(), PV3DBase< T, VectorTag, FrameTag >::phi(), reco::TrackBase::phi(), TrajectoryStateClosestToPoint::position(), reco::TrackBase::px(), reco::TrackBase::py(), reco::TrackBase::pz(), reco::TransientTrack::track(), reco::TransientTrack::trajectoryStateClosestToPoint(), PV3DBase< T, VectorTag, FrameTag >::x(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, VectorTag, FrameTag >::y(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and PV3DBase< T, VectorTag, FrameTag >::z().

Referenced by TauDiscriminationAgainstCaloMuon< TauType, TauDiscriminator >::discriminate().

{
  double hcalEnergySum = 0.;
  for ( HBHERecHitCollection::const_iterator hcalRecHit = hcalRecHits.begin();
        hcalRecHit != hcalRecHits.end(); ++hcalRecHit ) {
    const CaloCellGeometry* hbCellGeometry = hbGeometry->getGeometry(hcalRecHit->detid());
    const CaloCellGeometry* heCellGeometry = heGeometry->getGeometry(hcalRecHit->detid());

    const GlobalPoint* cellPosition = 0;
    if ( hbCellGeometry ) cellPosition = &(hbCellGeometry->getPosition());
    if ( heCellGeometry ) cellPosition = &(heCellGeometry->getPosition());

    if ( !cellPosition ) {
      edm::LogError ("compHcalEnergySum") 
        << " Failed to access HCAL geometry for detId = " << hcalRecHit->detid().rawId()
        << " --> skipping !!";
      continue;
    }

//--- CV: speed up computation by requiring eta-phi distance
//        between cell position and track direction to be dR < 0.5
    Vector3DBase<float, GlobalTag> cellPositionRelVertex = (*cellPosition) - eventVertexPosition;
    if ( deltaR(cellPositionRelVertex.eta(), cellPositionRelVertex.phi(), 
                transientTrack.track().eta(), transientTrack.track().phi()) > 0.5 ) continue;

    TrajectoryStateClosestToPoint dcaPosition = transientTrack.trajectoryStateClosestToPoint(*cellPosition);
    
    Vector3DBase<float, GlobalTag> d = ((*cellPosition) - dcaPosition.position());

    TVector3 d3(d.x(), d.y(), d.z());
    TVector3 dir(transientTrack.track().px(), transientTrack.track().py(), transientTrack.track().pz());

    double dPerp = d3.Cross(dir.Unit()).Mag();
    double dParl = TVector3(cellPosition->x(), cellPosition->y(), cellPosition->z()).Dot(dir.Unit());

    if ( dPerp < dR && dParl > 100. ) {
      hcalEnergySum += hcalRecHit->energy();
    }
  }
  
  return hcalEnergySum;
}
DEFINE_FWK_MODULE ( PFRecoTauDiscriminationAgainstCaloMuon  )