CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VZeroFinder.cc
Go to the documentation of this file.
2 
4 
6 
8 
9 #include <utility>
10 /*****************************************************************************/
12  (const edm::EventSetup& es,
13  const edm::ParameterSet& pset)
14 {
15  // Get track-pair level cuts
16  maxDca = pset.getParameter<double>("maxDca");
17 
18  // Get mother cuts
19  minCrossingRadius = pset.getParameter<double>("minCrossingRadius");
20  maxCrossingRadius = pset.getParameter<double>("maxCrossingRadius");
21  maxImpactMother = pset.getParameter<double>("maxImpactMother");
22 
23  // Get magnetic field
25  es.get<IdealMagneticFieldRecord>().get(magField);
26  theMagField = magField.product();
27 }
28 
29 /*****************************************************************************/
31 {
32 }
33 
34 /*****************************************************************************/
36  (const reco::Track& track)
37 {
38  GlobalPoint position(track.vertex().x(),
39  track.vertex().y(),
40  track.vertex().z());
41 
42  GlobalVector momentum(track.momentum().x(),
43  track.momentum().y(),
44  track.momentum().z());
45 
47  track.charge(),theMagField);
48 
49  return gtp;
50 }
51 
52 /*****************************************************************************/
54 {
55  double pt = p.perp();
56  return GlobalVector( -pt*cos(a), -pt*sin(a), p.z());
57 }
58 
59 /*****************************************************************************/
61  const reco::Track& negTrack,
64 {
65 /*
66  LogTrace("MinBiasTracking") << " [VZeroFinder] tracks"
67  << " +" << posTrack.algoName() << " " << posTrack.d0()
68  << " -" << negTrack.algoName() << " " << negTrack.d0();
69 */
70 
71  // Get trajectories
74 
75  // Two track minimum distance
77 
78  // Closest approach
79  theMinimum.calculate(posGtp,negGtp);
80 
81  // Closest points
82  // std::pair<GlobalPoint, GlobalPoint> points = theMinimum.points();
83 
84  // Momenta at closest point
85  std::pair<GlobalVector,GlobalVector> momenta;
86  momenta.first = rotate(posGtp.momentum(), theMinimum.firstAngle() );
87  momenta.second = rotate(negGtp.momentum(), theMinimum.secondAngle());
88 
89  // dcaR
90  float dca = theMinimum.distance();
91  GlobalPoint crossing = theMinimum.crossingPoint();
92 
93 /*
94  LogTrace("MinBiasTracking") << fixed << setprecision(2)
95  << " [VZeroFinder] dca = "<<dca<<" (<"<<maxDca<<")";
96  LogTrace("MinBiasTracking") << fixed << setprecision(2)
97  << " [VZeroFinder] crossR = "<<crossing.perp()<<" ("
98  <<minCrossingRadius<<"< <"<<maxCrossingRadius<<")";
99 */
100 
101  if(dca < maxDca &&
102  crossing.perp() > minCrossingRadius &&
103  crossing.perp() < maxCrossingRadius)
104  {
105  // Momentum of the mother
106  GlobalVector momentum = momenta.first + momenta.second;
107  float impact = -1.;
108 
109  if(vertices->size() > 0)
110  {
111  // Impact parameter of the mother wrt vertices, choose smallest
112  for(reco::VertexCollection::const_iterator
113  vertex = vertices->begin(); vertex!= vertices->end(); vertex++)
114  {
115  GlobalVector r(crossing.x(),
116  crossing.y(),
117  crossing.z() - vertex->position().z());
118  GlobalVector p(momentum.x(),momentum.y(),momentum.z());
119 
120  GlobalVector b = r - (r*p)*p / p.mag2();
121 
122  float im = b.mag();
123  if(im < impact || vertex == vertices->begin())
124  impact = im;
125  }
126  }
127  else
128  {
129  // Impact parameter of the mother in the plane
130  GlobalVector r_(crossing.x(),crossing.y(),0);
131  GlobalVector p_(momentum.x(),momentum.y(),0);
132 
133  GlobalVector b_ = r_ - (r_*p_)*p_ / p_.mag2();
134  impact = b_.mag();
135  }
136 
137 /*
138  LogTrace("MinBiasTracking") << fixed << setprecision(2)
139  << " [VZeroFinder] impact = "<<impact<<" (<"<<maxImpactMother<<")";
140 */
141 
142  if(impact < maxImpactMother)
143  {
144  // Armenteros variables
145  float pt =
146  (momenta.first.cross(momenta.second)).mag()/momentum.mag();
147  float alpha =
148  (momenta.first.mag2() - momenta.second.mag2())/momentum.mag2();
149 
150  // Fill data
151  data.dca = dca;
152  data.crossingPoint = crossing;
153  data.impactMother = impact;
154  data.armenterosPt = pt;
155  data.armenterosAlpha = alpha;
156  data.momenta = momenta;
157 
158 /*
159  LogTrace("MinBiasTracking") << fixed << setprecision(2)
160  << " [VZeroFinder] success {alpha = "<<alpha<<", pt = "<<pt<<"}";
161 */
162 
163  return true;
164  }
165  }
166 
167  return false;
168 }
169 
virtual float distance() const
T getParameter(std::string const &) const
math::GlobalPoint crossingPoint
Definition: VZero.h:18
virtual bool calculate(const TrajectoryStateOnSurface &sta, const TrajectoryStateOnSurface &stb)
float alpha
Definition: AMPTWrapper.h:95
T mag2() const
Definition: PV3DBase.h:66
T perp() const
Definition: PV3DBase.h:72
virtual GlobalPoint crossingPoint() const
GlobalTrajectoryParameters getGlobalTrajectoryParameters(const reco::Track &track)
Definition: VZeroFinder.cc:36
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:63
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:628
float armenterosPt
Definition: VZero.h:17
std::pair< GlobalVector, GlobalVector > momenta
Definition: VZero.h:19
float maxCrossingRadius
Definition: VZeroFinder.h:38
float maxImpactMother
Definition: VZeroFinder.h:38
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:640
T mag() const
Definition: PV3DBase.h:67
float dca
Definition: VZero.h:17
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
float armenterosAlpha
Definition: VZero.h:17
float minCrossingRadius
Definition: VZeroFinder.h:38
GlobalVector rotate(const GlobalVector &p, double a)
Definition: VZeroFinder.cc:53
float maxDca
Definition: VZeroFinder.h:38
float impactMother
Definition: VZero.h:17
const T & get() const
Definition: EventSetup.h:55
double b
Definition: hdecay.h:120
#define begin
Definition: vmac.h:30
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
VZeroFinder(const edm::EventSetup &es, const edm::ParameterSet &pset)
Definition: VZeroFinder.cc:12
int charge() const
track electric charge
Definition: TrackBase.h:520
bool checkTrackPair(const reco::Track &posTrack, const reco::Track &negTrack, const reco::VertexCollection *vertices, reco::VZeroData &data)
Definition: VZeroFinder.cc:60
T x() const
Definition: PV3DBase.h:62
Global3DVector GlobalVector
Definition: GlobalVector.h:10