CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Thrust.cc
Go to the documentation of this file.
1 // $Id: Thrust.cc,v 1.15 2011/10/12 15:30:11 gowdy Exp $
3 #include <cmath>
4 using namespace reco;
5 const double pi = M_PI, pi2 = 2 * pi, pi_2 = pi / 2, pi_4 = pi / 4;
6 
7 void Thrust::init(const std::vector<const Candidate*>& cands) {
8  int i = 0;
9  for(std::vector<const Candidate*>::const_iterator t = cands.begin();
10  t != cands.end(); ++t, ++i)
11  pSum_ += (p_[i] = (*t)->momentum()).r();
12  axis_ = axis(finalAxis(initialAxis()));
13  if (axis_.z() < 0) axis_ *= -1;
14  thrust_ = thrust(axis_);
15 }
16 
18  static const int nSegsTheta = 10, nSegsPhi = 10, nSegs = nSegsTheta * nSegsPhi;
19  int i, j;
20  double thr[nSegs], max = 0;
21  int indI = 0, indJ = 0, index = -1;
22  for (i = 0; i < nSegsTheta ; ++i) {
23  double z = cos(pi * i / (nSegsTheta - 1));
24  double r = sqrt(1 - z*z);
25  for (j = 0; j < nSegsPhi ; ++j) {
26  double phi = pi2 * j / nSegsPhi;
27  thr[i * nSegsPhi + j] = thrust(Vector(r*cos(phi), r*sin(phi), z));
28  if (thr[i*nSegsPhi + j] > max) {
29  index = i*nSegsPhi + j;
30  indI = i; indJ = j;
31  max = thr[index];
32  }
33  }
34  }
35 
36  // take max and one point on either size, fitting to a parabola and
37  // extrapolating to the real max. Do this separately for each dimension.
38  // y = a x^2 + b x + c. At the max, x = 0, on either side, x = +/-1.
39  // do phi first
40  double a, b, c = max;
41  int ind1 = indJ + 1;
42  if (ind1 >= nSegsPhi) ind1 -= nSegsPhi;
43  int ind2 = indJ - 1;
44  if (ind2 < 0) ind2 += nSegsPhi;
45  a = (thr[ind1] + thr[ind2] - 2*c) / 2;
46  b = thr[ind1] - a - c;
47  double maxPhiInd = 0;
48  if (a != 0) maxPhiInd = -b/(2*a);
49  double maxThetaInd;
50  if (indI == 0 || indI == (nSegsTheta - 1))
51  maxThetaInd = indI;
52  else {
53  ind1 = indI + 1;
54  ind2 = indI - 1;
55  a = (thr[ind1] + thr[ind2] - 2*c) / 2;
56  b = thr[ind1] - a - c;
57  maxThetaInd = 0;
58  if (a != 0) maxThetaInd = - b/(2*a);
59  }
60  return ThetaPhi(pi*(maxThetaInd + indI) / (nSegsTheta - 1),
61  pi2*(maxPhiInd + indJ) / nSegsPhi);
62 }
63 
65  static const double epsilon = 0.0001;
66  double maxChange1=0.0, maxChange2=0.0, a=0.0, b=0.0, c=0.0, thr=0.0;
67  int mandCt = 3, maxCt = 1000;
68  bool done;
69  do {
70  parabola(a, b, c,
71  axis(best),
72  axis(best.theta + epsilon, best.phi),
73  axis(best.theta - epsilon, best.phi));
74  maxChange1 = 10*(b < 0 ? -1 : 1);
75  if (a != 0) maxChange1 = - b/(2*a);
76  while (fabs(maxChange1 * epsilon) > pi_4) maxChange1 /= 2;
77  if (maxChange1 == 0 && (best.theta == 0 || best.theta == pi)) {
78  best.phi += pi_2;
79  if (best.phi > pi2) best.phi -= pi2;
80  parabola(a, b, c,
81  axis(best),
82  axis(best.theta + epsilon, best.phi),
83  axis(best.theta - epsilon, best.phi));
84  maxChange1 = 10 * (b < 0 ? -1 : 1);
85  if (a != 0) maxChange1 = - b / (2 * a);
86  }
87  do {
88  // L.L.: fixed odd behavoir adding epsilon (???)
89  thr = thrust(axis(best.theta + maxChange1 * epsilon, best.phi)) + epsilon;
90  if (thr < c) maxChange1 /= 2;
91  } while (thr < c);
92 
93  best.theta += maxChange1 * epsilon;
94  if (best.theta > pi) {
95  best.theta = pi - (best.theta - pi);
96  best.phi += pi;
97  if (best.phi > pi2) best.phi -= pi2;
98  }
99  if (best.theta < 0) {
100  best.theta *= -1;
101  best.phi += pi;
102  if (best.phi > pi2) best.phi -= pi2;
103  }
104  parabola(a, b, c,
105  axis(best),
106  axis(best.theta, best.phi + epsilon),
107  axis(best.theta, best.phi - epsilon));
108  maxChange2 = 10 * (b < 0 ? -1 : 1);
109  if (a != 0) maxChange2 = - b / (2 * a);
110  while (fabs(maxChange2 * epsilon) > pi_4) { maxChange2 /= 2; }
111  do {
112  // L.L.: fixed odd behavoir adding epsilon
113  thr = thrust(axis(best.theta, best.phi + maxChange2 * epsilon)) + epsilon;
114  if (thr < c) maxChange2 /= 2;
115  } while (thr < c);
116  best.phi += maxChange2 * epsilon;
117  if (best.phi > pi2) best.phi -= pi2;
118  if (best.phi < 0) best.phi += pi2;
119  if (mandCt > 0) mandCt --;
120  maxCt --;
121  done = (fabs(maxChange1) > 1 || fabs(maxChange2) > 1 || mandCt) && (maxCt > 0);
122  } while (done);
123 
124  return best;
125 }
126 
127 void Thrust::parabola(double & a, double & b, double & c,
128  const Vector & a1, const Vector & a2, const Vector & a3) const {
129  double t1 = thrust(a1), t2 = thrust(a2), t3 = thrust(a3);
130  a = (t2 - 2 * c + t3) / 2;
131  b = t2 - a - c;
132  c = t1;
133 }
134 
135 Thrust::Vector Thrust::axis(double theta, double phi) const {
136  double theSin = sin(theta);
137  return Vector(theSin * cos(phi), theSin * sin(phi), cos(theta));
138 }
139 
140 double Thrust::thrust(const Vector & axis) const {
141  double result = 0;
142  double sum = 0;
143  for (unsigned int i = 0; i < n_; ++i)
144  sum += fabs(axis.Dot(p_[i]));
145  if (pSum_ > 0) result = sum / pSum_;
146  return result;
147 }
int i
Definition: DBlmapReader.cc:9
const double pi_2
Definition: Thrust.cc:5
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:28
const Vector & axis() const
thrust axis (with magnitude = 1)
Definition: Thrust.h:59
double phi
Definition: Thrust.h:70
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
ThetaPhi initialAxis() const
Definition: Thrust.cc:17
double theta
Definition: Thrust.h:70
Geom::Theta< T > theta() const
const double pi2
Definition: Thrust.cc:5
double double double z
math::XYZVector Vector
spatial vector
Definition: Thrust.h:43
const T & max(const T &a, const T &b)
T sqrt(T t)
Definition: SSEVec.h:46
const double pi_4
Definition: Thrust.cc:5
tuple result
Definition: query.py:137
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const reco::Candidate::LorentzVector & axis_
double thrust() const
thrust value (in the range [0.5, 1.0])
Definition: Thrust.h:57
int j
Definition: DBlmapReader.cc:9
#define M_PI
Definition: BFit3D.cc:3
ThetaPhi finalAxis(ThetaPhi) const
Definition: Thrust.cc:64
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:127
void init(const std::vector< const reco::Candidate * > &)
Definition: Thrust.cc:7
const double epsilon
double pi
Definition: DDAxes.h:10