CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/PhysicsTools/CandUtils/src/Thrust.cc

Go to the documentation of this file.
00001 // $Id: Thrust.cc,v 1.15 2011/10/12 15:30:11 gowdy Exp $
00002 #include "PhysicsTools/CandUtils/interface/Thrust.h"
00003 #include <cmath>
00004 using namespace reco;
00005 const double pi = M_PI, pi2 = 2 * pi, pi_2 = pi / 2, pi_4 = pi / 4;
00006 
00007 void Thrust::init(const std::vector<const Candidate*>& cands) {
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 }
00016 
00017 Thrust::ThetaPhi Thrust::initialAxis() const {
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 }
00063 
00064 Thrust::ThetaPhi Thrust::finalAxis(ThetaPhi best) const {
00065   static const double epsilon = 0.0001;
00066   double maxChange1=0.0, maxChange2=0.0, a=0.0, b=0.0, c=0.0, thr=0.0;
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 }
00126 
00127 void Thrust::parabola(double & a, double & b, double & c, 
00128                        const Vector & a1, const Vector & a2, const Vector & a3) const {
00129   double t1 = thrust(a1), t2 = thrust(a2), t3 = thrust(a3);
00130   a = (t2 - 2 * c + t3) / 2;
00131   b = t2 - a - c;
00132   c = t1;
00133 }
00134 
00135 Thrust::Vector Thrust::axis(double theta, double phi) const {
00136   double theSin = sin(theta);
00137   return Vector(theSin * cos(phi), theSin * sin(phi), cos(theta));
00138 }
00139 
00140 double Thrust::thrust(const Vector & axis) const {
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 }