CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
Thrust Class Reference

#include <Thrust.h>

Classes

struct  ThetaPhi
 

Public Types

typedef math::XYZVector Vector
 spatial vector More...
 

Public Member Functions

const Vectoraxis () const
 thrust axis (with magnitude = 1) More...
 
template<typename const_iterator >
 Thrust (const_iterator begin, const_iterator end)
 constructor from first and last iterators More...
 
double thrust () const
 thrust value (in the range [0.5, 1.0]) More...
 

Private Member Functions

Vector axis (double theta, double phi) const
 
Vector axis (const ThetaPhi &tp) const
 
ThetaPhi finalAxis (ThetaPhi) const
 
void init (const std::vector< const reco::Candidate * > &)
 
ThetaPhi initialAxis () const
 
void parabola (double &a, double &b, double &c, const Vector &, const Vector &, const Vector &) const
 
double thrust (const Vector &theAxis) const
 

Private Attributes

Vector axis_
 
const unsigned int n_
 
std::vector< Vectorp_
 
double pSum_
 
double thrust_
 

Detailed Description

Utility to compute thrust value and axis from a collection of candidates.

Ported from BaBar implementation.

The thrust axis is the vector T which maximises the following expression:

  sum_i=1...N ( | T . p_i | )

t = ------------------------------— sum_i=1...N ( (p_i . _p_i)^(1/2) )

where p_i, i=1...N are the particle momentum vectors. The thrust value is the maximum value of t.

The thrust axis has a two-fold ambiguity due to taking the absolute value of the dot product. This computation returns by convention a thurs axis with a positive component along the z-direction is positive.

The thrust measure the alignment of a collection of particles along a common axis. The lower the thrust, the more spherical the event is. The higher the thrust, the more jet-like the

Author
Luca Lista, INFN

Definition at line 38 of file Thrust.h.

Member Typedef Documentation

spatial vector

Definition at line 41 of file Thrust.h.

Constructor & Destructor Documentation

template<typename const_iterator >
Thrust::Thrust ( const_iterator  begin,
const_iterator  end 
)
inline

constructor from first and last iterators

Definition at line 44 of file Thrust.h.

References end, mps_fire::i, init(), and n_.

44  :
45  thrust_(0), axis_(0, 0, 0), pSum_(0),
46  n_(end - begin), p_(n_) {
47  if (n_ == 0) return;
48  std::vector<const reco::Candidate*> cands;
49  for(const_iterator i = begin; i != end; ++i) {
50  cands.push_back(&*i);
51  }
52  init(cands);
53  }
double pSum_
Definition: Thrust.h:62
double thrust_
Definition: Thrust.h:60
const unsigned int n_
Definition: Thrust.h:63
std::vector< Vector > p_
Definition: Thrust.h:64
#define end
Definition: vmac.h:37
Vector axis_
Definition: Thrust.h:61
#define begin
Definition: vmac.h:30
void init(const std::vector< const reco::Candidate * > &)
Definition: Thrust.cc:6

Member Function Documentation

const Vector& Thrust::axis ( ) const
inline

thrust axis (with magnitude = 1)

Definition at line 57 of file Thrust.h.

References axis_.

Referenced by axis().

57 { return axis_; }
Vector axis_
Definition: Thrust.h:61
Thrust::Vector Thrust::axis ( double  theta,
double  phi 
) const
private

Definition at line 134 of file Thrust.cc.

References funct::cos(), and funct::sin().

134  {
135  double theSin = sin(theta);
136  return Vector(theSin * cos(phi), theSin * sin(phi), cos(theta));
137 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
math::XYZVector Vector
spatial vector
Definition: Thrust.h:41
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Vector Thrust::axis ( const ThetaPhi tp) const
inlineprivate

Definition at line 74 of file Thrust.h.

References a, axis(), b, EnergyCorrector::c, init(), parabola(), Thrust::ThetaPhi::phi, and Thrust::ThetaPhi::theta.

74  {
75  return axis(tp.theta, tp.phi);
76  }
const Vector & axis() const
thrust axis (with magnitude = 1)
Definition: Thrust.h:57
Thrust::ThetaPhi Thrust::finalAxis ( ThetaPhi  best) const
private

Definition at line 63 of file Thrust.cc.

References a, b, EnergyCorrector::c, geometryDiff::epsilon, Thrust::ThetaPhi::phi, pi, pi2, pi_2, pi_4, and Thrust::ThetaPhi::theta.

63  {
64  static const double epsilon = 0.0001;
65  double maxChange1=0.0, maxChange2=0.0, a=0.0, b=0.0, c=0.0, thr=0.0;
66  int mandCt = 3, maxCt = 1000;
67  bool done;
68  do {
69  parabola(a, b, c,
70  axis(best),
71  axis(best.theta + epsilon, best.phi),
72  axis(best.theta - epsilon, best.phi));
73  maxChange1 = 10*(b < 0 ? -1 : 1);
74  if (a != 0) maxChange1 = - b/(2*a);
75  while (fabs(maxChange1 * epsilon) > pi_4) maxChange1 /= 2;
76  if (maxChange1 == 0 && (best.theta == 0 || best.theta == pi)) {
77  best.phi += pi_2;
78  if (best.phi > pi2) best.phi -= pi2;
79  parabola(a, b, c,
80  axis(best),
81  axis(best.theta + epsilon, best.phi),
82  axis(best.theta - epsilon, best.phi));
83  maxChange1 = 10 * (b < 0 ? -1 : 1);
84  if (a != 0) maxChange1 = - b / (2 * a);
85  }
86  do {
87  // L.L.: fixed odd behavoir adding epsilon (???)
88  thr = thrust(axis(best.theta + maxChange1 * epsilon, best.phi)) + epsilon;
89  if (thr < c) maxChange1 /= 2;
90  } while (thr < c);
91 
92  best.theta += maxChange1 * epsilon;
93  if (best.theta > pi) {
94  best.theta = pi - (best.theta - pi);
95  best.phi += pi;
96  if (best.phi > pi2) best.phi -= pi2;
97  }
98  if (best.theta < 0) {
99  best.theta *= -1;
100  best.phi += pi;
101  if (best.phi > pi2) best.phi -= pi2;
102  }
103  parabola(a, b, c,
104  axis(best),
105  axis(best.theta, best.phi + epsilon),
106  axis(best.theta, best.phi - epsilon));
107  maxChange2 = 10 * (b < 0 ? -1 : 1);
108  if (a != 0) maxChange2 = - b / (2 * a);
109  while (fabs(maxChange2 * epsilon) > pi_4) { maxChange2 /= 2; }
110  do {
111  // L.L.: fixed odd behavoir adding epsilon
112  thr = thrust(axis(best.theta, best.phi + maxChange2 * epsilon)) + epsilon;
113  if (thr < c) maxChange2 /= 2;
114  } while (thr < c);
115  best.phi += maxChange2 * epsilon;
116  if (best.phi > pi2) best.phi -= pi2;
117  if (best.phi < 0) best.phi += pi2;
118  if (mandCt > 0) mandCt --;
119  maxCt --;
120  done = (fabs(maxChange1) > 1 || fabs(maxChange2) > 1 || mandCt) && (maxCt > 0);
121  } while (done);
122 
123  return best;
124 }
const double pi_2
Definition: Thrust.cc:4
const Vector & axis() const
thrust axis (with magnitude = 1)
Definition: Thrust.h:57
const double pi2
Definition: Thrust.cc:4
const double pi_4
Definition: Thrust.cc:4
double thrust() const
thrust value (in the range [0.5, 1.0])
Definition: Thrust.h:55
const double pi
Definition: Thrust.cc:4
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
void parabola(double &a, double &b, double &c, const Vector &, const Vector &, const Vector &) const
Definition: Thrust.cc:126
void Thrust::init ( const std::vector< const reco::Candidate * > &  cands)
private

Definition at line 6 of file Thrust.cc.

References axis_, mps_fire::i, alignCSCRings::r, and lumiQTWidget::t.

Referenced by axis(), and Thrust().

6  {
7  int i = 0;
8  for(std::vector<const Candidate*>::const_iterator t = cands.begin();
9  t != cands.end(); ++t, ++i)
10  pSum_ += (p_[i] = (*t)->momentum()).r();
12  if (axis_.z() < 0) axis_ *= -1;
13  thrust_ = thrust(axis_);
14 }
const Vector & axis() const
thrust axis (with magnitude = 1)
Definition: Thrust.h:57
ThetaPhi initialAxis() const
Definition: Thrust.cc:16
double pSum_
Definition: Thrust.h:62
double thrust_
Definition: Thrust.h:60
std::vector< Vector > p_
Definition: Thrust.h:64
double thrust() const
thrust value (in the range [0.5, 1.0])
Definition: Thrust.h:55
Vector axis_
Definition: Thrust.h:61
ThetaPhi finalAxis(ThetaPhi) const
Definition: Thrust.cc:63
Thrust::ThetaPhi Thrust::initialAxis ( ) const
private

Definition at line 16 of file Thrust.cc.

References a, b, EnergyCorrector::c, funct::cos(), mps_fire::i, diffTreeTool::index, SiStripPI::max, pi, pi2, alignCSCRings::r, funct::sin(), and mathSSE::sqrt().

16  {
17  static const int nSegsTheta = 10, nSegsPhi = 10, nSegs = nSegsTheta * nSegsPhi;
18  int i, j;
19  double thr[nSegs], max = 0;
20  int indI = 0, indJ = 0, index = -1;
21  for (i = 0; i < nSegsTheta ; ++i) {
22  double z = cos(pi * i / (nSegsTheta - 1));
23  double r = sqrt(1 - z*z);
24  for (j = 0; j < nSegsPhi ; ++j) {
25  double phi = pi2 * j / nSegsPhi;
26  thr[i * nSegsPhi + j] = thrust(Vector(r*cos(phi), r*sin(phi), z));
27  if (thr[i*nSegsPhi + j] > max) {
28  index = i*nSegsPhi + j;
29  indI = i; indJ = j;
30  max = thr[index];
31  }
32  }
33  }
34 
35  // take max and one point on either size, fitting to a parabola and
36  // extrapolating to the real max. Do this separately for each dimension.
37  // y = a x^2 + b x + c. At the max, x = 0, on either side, x = +/-1.
38  // do phi first
39  double a, b, c = max;
40  int ind1 = indJ + 1;
41  if (ind1 >= nSegsPhi) ind1 -= nSegsPhi;
42  int ind2 = indJ - 1;
43  if (ind2 < 0) ind2 += nSegsPhi;
44  a = (thr[ind1] + thr[ind2] - 2*c) / 2;
45  b = thr[ind1] - a - c;
46  double maxPhiInd = 0;
47  if (a != 0) maxPhiInd = -b/(2*a);
48  double maxThetaInd;
49  if (indI == 0 || indI == (nSegsTheta - 1))
50  maxThetaInd = indI;
51  else {
52  ind1 = indI + 1;
53  ind2 = indI - 1;
54  a = (thr[ind1] + thr[ind2] - 2*c) / 2;
55  b = thr[ind1] - a - c;
56  maxThetaInd = 0;
57  if (a != 0) maxThetaInd = - b/(2*a);
58  }
59  return ThetaPhi(pi*(maxThetaInd + indI) / (nSegsTheta - 1),
60  pi2*(maxPhiInd + indJ) / nSegsPhi);
61 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const double pi2
Definition: Thrust.cc:4
math::XYZVector Vector
spatial vector
Definition: Thrust.h:41
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double thrust() const
thrust value (in the range [0.5, 1.0])
Definition: Thrust.h:55
const double pi
Definition: Thrust.cc:4
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
void Thrust::parabola ( double &  a,
double &  b,
double &  c,
const Vector a1,
const Vector a2,
const Vector a3 
) const
private

Definition at line 126 of file Thrust.cc.

References EnergyCorrector::c, and reco::t2.

Referenced by axis().

127  {
128  double t1 = thrust(a1), t2 = thrust(a2), t3 = thrust(a3);
129  a = (t2 - 2 * c + t3) / 2;
130  b = t2 - a - c;
131  c = t1;
132 }
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
double thrust() const
thrust value (in the range [0.5, 1.0])
Definition: Thrust.h:55
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
double Thrust::thrust ( ) const
inline

thrust value (in the range [0.5, 1.0])

Definition at line 55 of file Thrust.h.

References thrust_.

55 { return thrust_; }
double thrust_
Definition: Thrust.h:60
double Thrust::thrust ( const Vector theAxis) const
private

Definition at line 139 of file Thrust.cc.

References mps_fire::i, and mps_fire::result.

139  {
140  double result = 0;
141  double sum = 0;
142  for (unsigned int i = 0; i < n_; ++i)
143  sum += fabs(axis.Dot(p_[i]));
144  if (pSum_ > 0) result = sum / pSum_;
145  return result;
146 }
const Vector & axis() const
thrust axis (with magnitude = 1)
Definition: Thrust.h:57
double pSum_
Definition: Thrust.h:62
const unsigned int n_
Definition: Thrust.h:63
std::vector< Vector > p_
Definition: Thrust.h:64

Member Data Documentation

Vector Thrust::axis_
private

Definition at line 61 of file Thrust.h.

Referenced by axis().

const unsigned int Thrust::n_
private

Definition at line 63 of file Thrust.h.

Referenced by Thrust().

std::vector<Vector> Thrust::p_
private

Definition at line 64 of file Thrust.h.

double Thrust::pSum_
private

Definition at line 62 of file Thrust.h.

double Thrust::thrust_
private

Definition at line 60 of file Thrust.h.

Referenced by thrust().