CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
InvariantMassAlgorithm.cc
Go to the documentation of this file.
2 #include <boost/regex.hpp>
3 
4 using namespace std;
5 using namespace edm;
6 using namespace reco;
7 
8 //
9 // -- Constructor
10 //
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
16  cluster_jet_matching_cone = 0.4;
17  cluster_track_matching_cone = 0.08;
18  inv_mass_cut = 2.5;
19 }
20 //
21 // -- Constructor
22 //
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");
35  trackAssociatorParameters_.loadParameters( tk_ass_pset );
36 
37  trackAssociator_ = new TrackDetectorAssociator();
38  trackAssociator_->useDefaultPropagator();
39 }
40 //
41 // -- Destructor
42 //
44  if(trackAssociator_) delete trackAssociator_;
45 }
46 //
47 // -- Tag
48 //
49 pair<double, reco::TauMassTagInfo> InvariantMassAlgorithm::tag(edm::Event& theEvent,
50  const edm::EventSetup& theEventSetup,
51  const reco::IsolatedTauTagInfoRef& tauRef,
52  const Handle<BasicClusterCollection>& clus_handle)
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,
76  leading_trk_pt,
77  signal_cone,
78  cluster_track_matching_cone,
79  inv_mass_cut);
80  }
81  return make_pair(discriminator, resultExtended);
82 }
83 //
84 // get Cluster Map
85 //
87  const edm::EventSetup& theEventSetup,
88  const reco::IsolatedTauTagInfoRef& tauRef,
89  const math::XYZVector& cluster_3vec)
90 {
91 
92  trackAssociatorParameters_.useEcal = true;
93  trackAssociatorParameters_.useHcal = false;
94  trackAssociatorParameters_.useHO = false;
95  trackAssociatorParameters_.useMuon = false;
96  trackAssociatorParameters_.dREcal = 0.03;
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))),
111  trackAssociatorParameters_);
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 }
T getParameter(std::string const &) const
dictionary parameters
Definition: Parameters.py:2
Base class for all types of Jets.
Definition: Jet.h:21
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
double double double z
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
std::pair< double, reco::TauMassTagInfo > tag(edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const edm::Handle< reco::BasicClusterCollection > &clus_handle)
float getMinimumClusterDR(edm::Event &theEvent, const edm::EventSetup &theEventSetup, const reco::IsolatedTauTagInfoRef &tauRef, const math::XYZVector &cluster_3vec)
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
virtual double px() const
x coordinate of momentum vector
float discriminator() const
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
void storeClusterTrackCollection(reco::BasicClusterRef clusterRef, float dr)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
virtual double py() const
y coordinate of momentum vector
void setIsolatedTauTag(const IsolatedTauTagInfoRef)