#include <PhysicsTools/CandUtils/interface/Thrust.h>
Public Types | |
typedef math::XYZVector | Vector |
spatial vector | |
Public Member Functions | |
const Vector & | axis () const |
thrust axis (with magnitude = 1) | |
double | thrust () const |
thrust value (in the range [0.5, 1.0]) | |
template<typename const_iterator> | |
Thrust (const_iterator begin, const_iterator end) | |
constructor from first and last iterators | |
Private Member Functions | |
Vector | axis (const ThetaPhi &tp) const |
Vector | axis (double theta, double phi) 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< Vector > | p_ |
double | pSum_ |
double | thrust_ |
Classes | |
struct | ThetaPhi |
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
Definition at line 40 of file Thrust.h.
typedef math::XYZVector Thrust::Vector |
Thrust::Thrust | ( | const_iterator | begin, | |
const_iterator | end | |||
) | [inline] |
constructor from first and last iterators
Definition at line 46 of file Thrust.h.
00046 : 00047 thrust_(0), axis_(0, 0, 0), pSum_(0), 00048 n_(end - begin), p_(n_) { 00049 if (n_ == 0) return; 00050 std::vector<const reco::Candidate*> cands; 00051 for(const_iterator i = begin; i != end; ++i) { 00052 cands.push_back(&*i); 00053 } 00054 init(cands); 00055 }
Definition at line 76 of file Thrust.h.
References axis(), Thrust::ThetaPhi::phi, and Thrust::ThetaPhi::theta.
Thrust::Vector Thrust::axis | ( | double | theta, | |
double | phi | |||
) | const [private] |
Thrust::ThetaPhi Thrust::finalAxis | ( | ThetaPhi | best | ) | const [private] |
Definition at line 64 of file Thrust.cc.
References a, axis(), b, c, geometryDiff::epsilon, parabola(), Thrust::ThetaPhi::phi, pi, pi2, pi_2, pi_4, Thrust::ThetaPhi::theta, and thrust().
Referenced by init().
00064 { 00065 static const double epsilon = 0.0001; 00066 double maxChange1, maxChange2, a, b, c, thr; 00067 int mandCt = 3, maxCt = 1000; 00068 bool done; 00069 do { 00070 parabola(a, b, c, 00071 axis(best), 00072 axis(best.theta + epsilon, best.phi), 00073 axis(best.theta - epsilon, best.phi)); 00074 maxChange1 = 10*(b < 0 ? -1 : 1); 00075 if (a != 0) maxChange1 = - b/(2*a); 00076 while (fabs(maxChange1 * epsilon) > pi_4) maxChange1 /= 2; 00077 if (maxChange1 == 0 && (best.theta == 0 || best.theta == pi)) { 00078 best.phi += pi_2; 00079 if (best.phi > pi2) best.phi -= pi2; 00080 parabola(a, b, c, 00081 axis(best), 00082 axis(best.theta + epsilon, best.phi), 00083 axis(best.theta - epsilon, best.phi)); 00084 maxChange1 = 10 * (b < 0 ? -1 : 1); 00085 if (a != 0) maxChange1 = - b / (2 * a); 00086 } 00087 do { 00088 // L.L.: fixed odd behavoir adding epsilon (???) 00089 thr = thrust(axis(best.theta + maxChange1 * epsilon, best.phi)) + epsilon; 00090 if (thr < c) maxChange1 /= 2; 00091 } while (thr < c); 00092 00093 best.theta += maxChange1 * epsilon; 00094 if (best.theta > pi) { 00095 best.theta = pi - (best.theta - pi); 00096 best.phi += pi; 00097 if (best.phi > pi2) best.phi -= pi2; 00098 } 00099 if (best.theta < 0) { 00100 best.theta *= -1; 00101 best.phi += pi; 00102 if (best.phi > pi2) best.phi -= pi2; 00103 } 00104 parabola(a, b, c, 00105 axis(best), 00106 axis(best.theta, best.phi + epsilon), 00107 axis(best.theta, best.phi - epsilon)); 00108 maxChange2 = 10 * (b < 0 ? -1 : 1); 00109 if (a != 0) maxChange2 = - b / (2 * a); 00110 while (fabs(maxChange2 * epsilon) > pi_4) { maxChange2 /= 2; } 00111 do { 00112 // L.L.: fixed odd behavoir adding epsilon 00113 thr = thrust(axis(best.theta, best.phi + maxChange2 * epsilon)) + epsilon; 00114 if (thr < c) maxChange2 /= 2; 00115 } while (thr < c); 00116 best.phi += maxChange2 * epsilon; 00117 if (best.phi > pi2) best.phi -= pi2; 00118 if (best.phi < 0) best.phi += pi2; 00119 if (mandCt > 0) mandCt --; 00120 maxCt --; 00121 done = (fabs(maxChange1) > 1 || fabs(maxChange2) > 1 || mandCt) && (maxCt > 0); 00122 } while (done); 00123 00124 return best; 00125 }
void Thrust::init | ( | const std::vector< const reco::Candidate * > & | cands | ) | [private] |
Definition at line 7 of file Thrust.cc.
References axis(), axis_, finalAxis(), i, initialAxis(), p_, pSum_, r, t, thrust(), and thrust_.
Referenced by Thrust().
00007 { 00008 int i = 0; 00009 for(std::vector<const Candidate*>::const_iterator t = cands.begin(); 00010 t != cands.end(); ++t, ++i) 00011 pSum_ += (p_[i] = (*t)->momentum()).r(); 00012 axis_ = axis(finalAxis(initialAxis())); 00013 if (axis_.z() < 0) axis_ *= -1; 00014 thrust_ = thrust(axis_); 00015 }
Thrust::ThetaPhi Thrust::initialAxis | ( | ) | const [private] |
Definition at line 17 of file Thrust.cc.
References a, b, c, funct::cos(), i, index, j, max, phi, pi, pi2, r, funct::sin(), funct::sqrt(), thrust(), and z.
Referenced by init().
00017 { 00018 static const int nSegsTheta = 10, nSegsPhi = 10, nSegs = nSegsTheta * nSegsPhi; 00019 int i, j; 00020 double thr[nSegs], max = 0; 00021 int indI = 0, indJ = 0, index = -1; 00022 for (i = 0; i < nSegsTheta ; ++i) { 00023 double z = cos(pi * i / (nSegsTheta - 1)); 00024 double r = sqrt(1 - z*z); 00025 for (j = 0; j < nSegsPhi ; ++j) { 00026 double phi = pi2 * j / nSegsPhi; 00027 thr[i * nSegsPhi + j] = thrust(Vector(r*cos(phi), r*sin(phi), z)); 00028 if (thr[i*nSegsPhi + j] > max) { 00029 index = i*nSegsPhi + j; 00030 indI = i; indJ = j; 00031 max = thr[index]; 00032 } 00033 } 00034 } 00035 00036 // take max and one point on either size, fitting to a parabola and 00037 // extrapolating to the real max. Do this separately for each dimension. 00038 // y = a x^2 + b x + c. At the max, x = 0, on either side, x = +/-1. 00039 // do phi first 00040 double a, b, c = max; 00041 int ind1 = indJ + 1; 00042 if (ind1 >= nSegsPhi) ind1 -= nSegsPhi; 00043 int ind2 = indJ - 1; 00044 if (ind2 < 0) ind2 += nSegsPhi; 00045 a = (thr[ind1] + thr[ind2] - 2*c) / 2; 00046 b = thr[ind1] - a - c; 00047 double maxPhiInd = 0; 00048 if (a != 0) maxPhiInd = -b/(2*a); 00049 double maxThetaInd; 00050 if (indI == 0 || indI == (nSegsTheta - 1)) 00051 maxThetaInd = indI; 00052 else { 00053 ind1 = indI + 1; 00054 ind2 = indI - 1; 00055 a = (thr[ind1] + thr[ind2] - 2*c) / 2; 00056 b = thr[ind1] - a - c; 00057 maxThetaInd = 0; 00058 if (a != 0) maxThetaInd = - b/(2*a); 00059 } 00060 return ThetaPhi(pi*(maxThetaInd + indI) / (nSegsTheta - 1), 00061 pi2*(maxPhiInd + indJ) / nSegsPhi); 00062 }
double Thrust::thrust | ( | const Vector & | theAxis | ) | const [private] |
Definition at line 140 of file Thrust.cc.
References i, n_, p_, pSum_, HLT_VtxMuL3::result, and sum().
00140 { 00141 double result = 0; 00142 double sum = 0; 00143 for (unsigned int i = 0; i < n_; ++i) 00144 sum += fabs(axis.Dot(p_[i])); 00145 if (pSum_ > 0) result = sum / pSum_; 00146 return result; 00147 }
double Thrust::thrust | ( | ) | const [inline] |
thrust value (in the range [0.5, 1.0])
Definition at line 57 of file Thrust.h.
References thrust_.
Referenced by finalAxis(), init(), initialAxis(), and parabola().
00057 { return thrust_; }
Vector Thrust::axis_ [private] |
const unsigned int Thrust::n_ [private] |
std::vector<Vector> Thrust::p_ [private] |
double Thrust::pSum_ [private] |
double Thrust::thrust_ [private] |