CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Thrust.cc
Go to the documentation of this file.
2 #include <cmath>
3 using namespace reco;
4 const double pi = M_PI, pi2 = 2 * pi, pi_2 = pi / 2, pi_4 = pi / 4;
5 
6 void Thrust::init(const std::vector<const Candidate*>& cands) {
7  int i = 0;
8  for (std::vector<const Candidate*>::const_iterator t = cands.begin(); t != cands.end(); ++t, ++i)
9  pSum_ += (p_[i] = (*t)->momentum()).r();
10  axis_ = axis(finalAxis(initialAxis()));
11  if (axis_.z() < 0)
12  axis_ *= -1;
13  thrust_ = thrust(axis_);
14 }
15 
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;
30  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)
43  ind1 -= nSegsPhi;
44  int ind2 = indJ - 1;
45  if (ind2 < 0)
46  ind2 += nSegsPhi;
47  a = (thr[ind1] + thr[ind2] - 2 * c) / 2;
48  b = thr[ind1] - a - c;
49  double maxPhiInd = 0;
50  if (a != 0)
51  maxPhiInd = -b / (2 * a);
52  double maxThetaInd;
53  if (indI == 0 || indI == (nSegsTheta - 1))
54  maxThetaInd = indI;
55  else {
56  ind1 = indI + 1;
57  ind2 = indI - 1;
58  a = (thr[ind1] + thr[ind2] - 2 * c) / 2;
59  b = thr[ind1] - a - c;
60  maxThetaInd = 0;
61  if (a != 0)
62  maxThetaInd = -b / (2 * a);
63  }
64  return ThetaPhi(pi * (maxThetaInd + indI) / (nSegsTheta - 1), pi2 * (maxPhiInd + indJ) / nSegsPhi);
65 }
66 
68  static const double epsilon = 0.0001;
69  double maxChange1 = 0.0, maxChange2 = 0.0, a = 0.0, b = 0.0, c = 0.0, thr = 0.0;
70  int mandCt = 3, maxCt = 1000;
71  bool done;
72  do {
73  parabola(a, b, c, axis(best), axis(best.theta + epsilon, best.phi), axis(best.theta - epsilon, best.phi));
74  maxChange1 = 10 * (b < 0 ? -1 : 1);
75  if (a != 0)
76  maxChange1 = -b / (2 * a);
77  while (fabs(maxChange1 * epsilon) > pi_4)
78  maxChange1 /= 2;
79  if (maxChange1 == 0 && (best.theta == 0 || best.theta == pi)) {
80  best.phi += pi_2;
81  if (best.phi > pi2)
82  best.phi -= pi2;
83  parabola(a, b, c, axis(best), axis(best.theta + epsilon, best.phi), axis(best.theta - epsilon, best.phi));
84  maxChange1 = 10 * (b < 0 ? -1 : 1);
85  if (a != 0)
86  maxChange1 = -b / (2 * a);
87  }
88  do {
89  // L.L.: fixed odd behavoir adding epsilon (???)
90  thr = thrust(axis(best.theta + maxChange1 * epsilon, best.phi)) + epsilon;
91  if (thr < c)
92  maxChange1 /= 2;
93  } while (thr < c);
94 
95  best.theta += maxChange1 * epsilon;
96  if (best.theta > pi) {
97  best.theta = pi - (best.theta - pi);
98  best.phi += pi;
99  if (best.phi > pi2)
100  best.phi -= pi2;
101  }
102  if (best.theta < 0) {
103  best.theta *= -1;
104  best.phi += pi;
105  if (best.phi > pi2)
106  best.phi -= pi2;
107  }
108  parabola(a, b, c, axis(best), axis(best.theta, best.phi + epsilon), axis(best.theta, best.phi - epsilon));
109  maxChange2 = 10 * (b < 0 ? -1 : 1);
110  if (a != 0)
111  maxChange2 = -b / (2 * a);
112  while (fabs(maxChange2 * epsilon) > pi_4) {
113  maxChange2 /= 2;
114  }
115  do {
116  // L.L.: fixed odd behavoir adding epsilon
117  thr = thrust(axis(best.theta, best.phi + maxChange2 * epsilon)) + epsilon;
118  if (thr < c)
119  maxChange2 /= 2;
120  } while (thr < c);
121  best.phi += maxChange2 * epsilon;
122  if (best.phi > pi2)
123  best.phi -= pi2;
124  if (best.phi < 0)
125  best.phi += pi2;
126  if (mandCt > 0)
127  mandCt--;
128  maxCt--;
129  done = (fabs(maxChange1) > 1 || fabs(maxChange2) > 1 || mandCt) && (maxCt > 0);
130  } while (done);
131 
132  return best;
133 }
134 
135 void Thrust::parabola(double& a, double& b, double& c, const Vector& a1, const Vector& a2, const Vector& a3) const {
136  double t1 = thrust(a1), t2 = thrust(a2), t3 = thrust(a3);
137  a = (t2 - 2 * c + t3) / 2;
138  b = t2 - a - c;
139  c = t1;
140 }
141 
142 Thrust::Vector Thrust::axis(double theta, double phi) const {
143  double theSin = sin(theta);
144  return Vector(theSin * cos(phi), theSin * sin(phi), cos(theta));
145 }
146 
147 double Thrust::thrust(const Vector& axis) const {
148  double result = 0;
149  double sum = 0;
150  for (unsigned int i = 0; i < n_; ++i)
151  sum += fabs(axis.Dot(p_[i]));
152  if (pSum_ > 0)
153  result = sum / pSum_;
154  return result;
155 }
const edm::EventSetup & c
const double pi_2
Definition: Thrust.cc:4
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:29
const Vector & axis() const
thrust axis (with magnitude = 1)
Definition: Thrust.h:56
double phi
Definition: Thrust.h:67
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
ThetaPhi initialAxis() const
Definition: Thrust.cc:16
double theta
Definition: Thrust.h:67
Geom::Theta< T > theta() const
const double pi2
Definition: Thrust.cc:4
tuple result
Definition: mps_fire.py:311
const Double_t pi
math::XYZVector Vector
spatial vector
Definition: Thrust.h:41
T sqrt(T t)
Definition: SSEVec.h:19
const double pi_4
Definition: Thrust.cc:4
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:54
#define M_PI
ThetaPhi finalAxis(ThetaPhi) const
Definition: Thrust.cc:67
double b
Definition: hdecay.h:118
double a
Definition: hdecay.h:119
void parabola(double &a, double &b, double &c, const Vector &, const Vector &, const Vector &) const
Definition: Thrust.cc:135
void init(const std::vector< const reco::Candidate * > &)
Definition: Thrust.cc:6