CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TauMassTagInfo.cc
Go to the documentation of this file.
6 
7 using namespace edm;
8 using namespace reco;
9 using namespace std;
10 
11 float reco::TauMassTagInfo::discriminator(double matching_cone, double leading_trk_pt,
12  double signal_cone, double cluster_track_cone,
13  double m_cut) const {
14  float discriminator = 0.0;
15  double invariantMass = getInvariantMass(matching_cone,leading_trk_pt, signal_cone,cluster_track_cone);
16  if (invariantMass >= 0.0 && invariantMass < m_cut ) discriminator = 1.0;
17  return discriminator;
18 }
19 //
20 // -- Set IsolatedTauTag
21 //
23  isolatedTau = isolationRef;
24 }
25 //
26 // -- Get IsolatedTauTag
27 //
29  return isolatedTau;
30 }
31 //
32 // -- Set Cluster Collection
33 //
35 
36  clusterMap.insert(clusterRef, dr);
37 }
38 //
39 // -- Calculate 4 momentum vector from tracks only
40 //
41 bool reco::TauMassTagInfo::calculateTrkP4(double matching_cone,
42  double leading_trk_pt, double signal_cone,
43  math::XYZTLorentzVector& p4) const {
44 
45 
46  const TrackRef leadTk= isolatedTau->leadingSignalTrack(matching_cone,leading_trk_pt);
47  if (!leadTk) {
48  std::cout <<" TauMassTagInfo:: No Leading Track !! " << std::endl;
49  return false;
50  }
51  math::XYZVector momentum = (*leadTk).momentum();
52  const RefVector<TrackCollection> signalTracks =
53  isolatedTau->tracksInCone(momentum,signal_cone,1.0);
54 // if (signalTracks.size() == 0 || signalTracks.size()%2 == 0) return false;
55  if (signalTracks.size() == 0) return false;
56 
57  double px_inv = 0.0;
58  double py_inv = 0.0;
59  double pz_inv = 0.0;
60  double e_inv = 0.0;
61  for (RefVector<TrackCollection>::const_iterator itrack = signalTracks.begin();
62  itrack != signalTracks.end(); itrack++) {
63  double p = (*itrack)->p();
64  double energy = sqrt(p*p + 0.139*0.139); // use correct value!
65  px_inv += (*itrack)->px();
66  py_inv += (*itrack)->py();
67  pz_inv += (*itrack)->pz();
68  e_inv += energy;
69  }
70 
71  p4.SetPx(px_inv);
72  p4.SetPy(py_inv);
73  p4.SetPz(pz_inv);
74  p4.SetE(e_inv);
75 
76  return true;
77 }
78 //
79 // -- Get Invariant Mass
80 //
81 double reco::TauMassTagInfo::getInvariantMassTrk(double matching_cone,
82  double leading_trk_pt, double signal_cone) const
83 {
85  if (!calculateTrkP4(matching_cone,leading_trk_pt, signal_cone, totalP4)) return -1.0;
86  return totalP4.M();
87 }
88 //
89 // -- Get Invariant Mass
90 //
91 double reco::TauMassTagInfo::getInvariantMass(double matching_cone,
92  double leading_trk_pt, double signal_cone,
93  double track_cone) const
94 {
95 
97  if (!calculateTrkP4(matching_cone, leading_trk_pt, signal_cone, totalP4)) return -1.0;
98 
99  // Add Clusters away from tracks
100  for (ClusterTrackAssociationCollection::const_iterator mapIter = clusterMap.begin();
101  mapIter != clusterMap.end(); mapIter++) {
102  const reco::BasicClusterRef & iclus = mapIter->key;
103  float dr = mapIter->val;
104  if (dr > track_cone) {
105  math::XYZVector clus3Vec(iclus->x(), iclus->y(), iclus->z());
106  double e = iclus->energy();
107  double theta = clus3Vec.theta();
108  double phi = clus3Vec.phi();
109  double px = e * sin(theta) * cos(phi);
110  double py = e * sin(theta) * sin(phi);
111  double pz = e * cos(theta);
112  math::XYZTLorentzVector p4(px,py,pz,e);
113  totalP4 += p4;
114  }
115  }
116  return totalP4.M();
117 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
tuple discriminator
key_type key() const
Accessor for product key.
Definition: Ref.h:264
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:255
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:250
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const IsolatedTauTagInfoRef & getIsolatedTauTag() const
T sqrt(T t)
Definition: SSEVec.h:48
double p4[4]
Definition: TauolaWrapper.h:92
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double getInvariantMass(double matching_cone, double leading_trk_pt, double signal_cone, double cluster_track_cone) const
double getInvariantMassTrk(double matching_cone, double leading_trk_pt, double signal_cone) const
float discriminator() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
bool calculateTrkP4(double matching_cone, double leading_trk_pt, double signal_cone, math::XYZTLorentzVector &p4) const
Geom::Phi< T > phi() const
void storeClusterTrackCollection(reco::BasicClusterRef clusterRef, float dr)
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
tuple cout
Definition: gather_cfg.py:121
void setIsolatedTauTag(const IsolatedTauTagInfoRef)