CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
InvariantMassAlgorithm Class Reference

#include <InvariantMassAlgorithm.h>

Public Member Functions

float getMinimumClusterDR (edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const math::XYZVector &cluster_3vec)
 
 InvariantMassAlgorithm (const edm::ParameterSet &parameters)
 
 InvariantMassAlgorithm ()
 
std::pair< double,
reco::TauMassTagInfo
tag (edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const edm::Handle< reco::BasicClusterCollection > &clus_handle)
 
 ~InvariantMassAlgorithm ()
 

Private Attributes

double cluster_jet_matching_cone
 
double cluster_track_matching_cone
 
double inv_mass_cut
 
double leading_trk_pt
 
double matching_cone
 
double signal_cone
 
TrackDetectorAssociatortrackAssociator_
 
TrackAssociatorParameters trackAssociatorParameters_
 

Detailed Description

Definition at line 26 of file InvariantMassAlgorithm.h.

Constructor & Destructor Documentation

InvariantMassAlgorithm::InvariantMassAlgorithm ( const edm::ParameterSet parameters)

Definition at line 23 of file InvariantMassAlgorithm.cc.

References edm::ParameterSet::getParameter().

24 {
25  matching_cone = parameters.getParameter<double>("MatchingCone");
26  leading_trk_pt = parameters.getParameter<double>("LeadingTrackPt");
27  signal_cone = parameters.getParameter<double>("SignalCone");
28  cluster_jet_matching_cone = parameters.getParameter<double>("ClusterSelectionCone");
29  cluster_track_matching_cone = parameters.getParameter<double>("ClusterTrackMatchingCone");
30  inv_mass_cut = parameters.getParameter<double>("InvariantMassCutoff");
31 
32 
33 // TrackAssociator parameters
34  edm::ParameterSet tk_ass_pset = parameters.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
36 
39 }
TrackDetectorAssociator * trackAssociator_
T getParameter(std::string const &) const
void useDefaultPropagator()
use the default propagator
void loadParameters(const edm::ParameterSet &)
TrackAssociatorParameters trackAssociatorParameters_
InvariantMassAlgorithm::InvariantMassAlgorithm ( )

Definition at line 11 of file InvariantMassAlgorithm.cc.

12 {
13  matching_cone = 0.1; // check with simone
14  leading_trk_pt = 1.0; // check with simone
15  signal_cone = 0.07;// check with simone
18  inv_mass_cut = 2.5;
19 }
InvariantMassAlgorithm::~InvariantMassAlgorithm ( )

Definition at line 43 of file InvariantMassAlgorithm.cc.

43  {
45 }
TrackDetectorAssociator * trackAssociator_

Member Function Documentation

float InvariantMassAlgorithm::getMinimumClusterDR ( edm::Event theEvent,
const edm::EventSetup theEventSetup,
const reco::IsolatedTauTagInfoRef tauRef,
const math::XYZVector cluster_3vec 
)

Definition at line 86 of file InvariantMassAlgorithm.cc.

References edm::RefVector< C, T, F >::begin(), Geom::deltaR2(), edm::RefVector< C, T, F >::end(), info, metsig::jet, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), testEve_cfg::tracks, and TrackDetMatchInfo::trkGlobPosAtEcal.

90 {
91 
97 
98  const TrackRefVector tracks = tauRef->allTracks();
99  const Jet & jet = *(tauRef->jet());
100  math::XYZVector jet3Vec(jet.px(),jet.py(),jet.pz());
101  float min_dR = 999.9;
102  float deltaR1 = ROOT::Math::VectorUtil::DeltaR(cluster_3vec, jet3Vec);
103  if (deltaR1 > cluster_jet_matching_cone) return -1.0;
104 
106  it != tracks.end(); it++)
107  {
109  info = trackAssociator_->associate(theEvent, theEventSetup,
110  trackAssociator_->getFreeTrajectoryState(theEventSetup, (*(*it))),
112  math::XYZVector track3Vec(info.trkGlobPosAtEcal.x(),
113  info.trkGlobPosAtEcal.y(),
114  info.trkGlobPosAtEcal.z());
115 
116  float deltaR2 = ROOT::Math::VectorUtil::DeltaR(cluster_3vec, track3Vec);
117  if (deltaR2 < min_dR) min_dR = deltaR2;
118  }
119  return min_dR;
120 }
TrackDetectorAssociator * trackAssociator_
static FreeTrajectoryState getFreeTrajectoryState(const edm::EventSetup &, const reco::Track &)
get FreeTrajectoryState from different track representations
Base class for all types of Jets.
Definition: Jet.h:21
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
virtual double px() const
x coordinate of momentum vector
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
virtual double pz() const
z coordinate of momentum vector
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
virtual double py() const
y coordinate of momentum vector
TrackAssociatorParameters trackAssociatorParameters_
pair< double, reco::TauMassTagInfo > InvariantMassAlgorithm::tag ( edm::Event theEvent,
const edm::EventSetup theEventSetup,
const reco::IsolatedTauTagInfoRef tauRef,
const edm::Handle< reco::BasicClusterCollection > &  clus_handle 
)

Definition at line 49 of file InvariantMassAlgorithm.cc.

References gather_cfg::cout, reco::TauMassTagInfo::discriminator(), reco::TauMassTagInfo::setIsolatedTauTag(), funct::sin(), reco::TauMassTagInfo::storeClusterTrackCollection(), vdt::x, detailsBasic3DVector::y, and detailsBasic3DVector::z.

53 {
54 
55  TauMassTagInfo resultExtended;
56  resultExtended.setIsolatedTauTag(tauRef);
57 
58  double discriminator = tauRef->discriminator();
59  if (discriminator > 0.0) {
60 
61  int nSel = 0;
62  for (size_t ic = 0; ic < clus_handle->size(); ic++) {
63  math::XYZVector cluster3Vec((*clus_handle)[ic].x(),(*clus_handle)[ic].y(),(*clus_handle)[ic].z());
64  BasicClusterRef cluster(clus_handle, ic);
65  float et = (*clus_handle)[ic].energy() * sin(cluster3Vec.theta());
66  if ( (et < 0.04) || (et > 300.0) ) continue;
67  float delR = getMinimumClusterDR(theEvent, theEventSetup, tauRef, cluster3Vec);
68  if (delR != -1.0) {
69  nSel++;
70  resultExtended.storeClusterTrackCollection(cluster, delR);
71  }
72  }
73  if (0) cout << " Total Number of Clusters " << clus_handle->size() << " " << nSel << endl;
74 
75  discriminator = resultExtended.discriminator(matching_cone,
79  inv_mass_cut);
80  }
81  return make_pair(discriminator, resultExtended);
82 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double double double z
float getMinimumClusterDR(edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const math::XYZVector &cluster_3vec)
float discriminator() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
void storeClusterTrackCollection(reco::BasicClusterRef clusterRef, float dr)
tuple cout
Definition: gather_cfg.py:121
x
Definition: VDTMath.h:216
void setIsolatedTauTag(const IsolatedTauTagInfoRef)

Member Data Documentation

double InvariantMassAlgorithm::cluster_jet_matching_cone
private

Definition at line 50 of file InvariantMassAlgorithm.h.

double InvariantMassAlgorithm::cluster_track_matching_cone
private

Definition at line 51 of file InvariantMassAlgorithm.h.

double InvariantMassAlgorithm::inv_mass_cut
private

Definition at line 52 of file InvariantMassAlgorithm.h.

double InvariantMassAlgorithm::leading_trk_pt
private

Definition at line 48 of file InvariantMassAlgorithm.h.

double InvariantMassAlgorithm::matching_cone
private

Definition at line 47 of file InvariantMassAlgorithm.h.

double InvariantMassAlgorithm::signal_cone
private

Definition at line 49 of file InvariantMassAlgorithm.h.

TrackDetectorAssociator* InvariantMassAlgorithm::trackAssociator_
private

Definition at line 54 of file InvariantMassAlgorithm.h.

TrackAssociatorParameters InvariantMassAlgorithm::trackAssociatorParameters_
private

Definition at line 55 of file InvariantMassAlgorithm.h.