rivet is hosted by Hepforge, IPPP Durham
CMS_2014_I1305624.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/FinalState.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 
00006 namespace Rivet {
00007 
00008 
00009   namespace {
00010 
00011     /// Number of event shape variables
00012     /// @todo Move into the EventShape class
00013     const int NEVTVAR = 5;
00014 
00015     /// Number of leading jet pT thresholds
00016     /// @todo Move into the analysis class
00017     const int NJETPTMN = 5;
00018     /// Leading jet pT thresholds
00019     /// @todo Move into the analysis class
00020     const double LEADINGPTTHRESHOLD[NJETPTMN] = { 110.0, 170.0, 250.0, 320.0, 390.0 };
00021 
00022 
00023     // Helpers for event shape calculations in hidden namespace; implementation at bottom of file
00024     /// @todo Why a class? Improve/remove this junk
00025     class EventShape {
00026     public:
00027 
00028       /// Constructor from vectors of four-vectors as input objects in the event to calculate the event shapes
00029       EventShape(const vector<double>& px_vector, const vector<double>& py_vector, const vector<double>& pz_vector,
00030                  const vector<double>& e_vector, double eta_central, int irap, int nmn)
00031         : _object_px(px_vector), _object_py(py_vector), _object_pz(pz_vector),
00032           _object_e(e_vector), _eta_c(eta_central), _irap(irap), _nmnjet(nmn)
00033       {   }
00034 
00035       /// @brief Returns the values of the five event shapes
00036       ///
00037       /// Event shape indices:
00038       /// 0. central transverse thrust
00039       /// 1. central total jet broadening
00040       /// 2. central total jet mass
00041       /// 3. central total transverse jet mass
00042       /// 4. central three-jet resolution threshold
00043       vector<double> getEventShapes() {
00044         _calculate(); ///< @todo There should be some test for success/failure!!
00045         return _event_shapes;
00046       }
00047 
00048       /// Returns the global thrust axis Nx, Ny, Nz=0
00049       vector<double> getThrustAxis() {
00050         _calculate(); ///< @todo There should be some test for success/failure!!
00051         return _thrust_axis;
00052       }
00053 
00054       /// Returns the central thrust axis Nx, Ny, Nz=0
00055       vector<double> getThrustAxisC() {
00056         _calculate(); ///< @todo There should be some test for success/failure!!
00057         return _thrust_axis_c;
00058       }
00059 
00060       // /// @brief Choice of the central region
00061       // void setEtaC(double eta_central) { _eta_c = eta_central; }
00062 
00063       // // Whether to use the rapidity y (rap==1)  or the pseudorapidity eta (rap==0)
00064       // void setRapType(int irap) { _irap = irap; }
00065 
00066 
00067     private:
00068 
00069       /// Calculate everything
00070       int _calculate();
00071 
00072       /// Returns the difference in phi between two vectors
00073       double _delta_phi(double, double);
00074 
00075       /// The Lorentz scalar product
00076       double _lorentz_sp(const vector<double>&, const vector<double>&);
00077 
00078       // Calculates the three-jet resolutions
00079       double _three_jet_res(const vector<double>&, const vector<double>&, const vector<double>&, const vector<double>&, int);
00080 
00081       // Calculates the thrust axis and the tau values
00082       vector<double> _thrust(const vector<double>&, const vector<double>&);
00083 
00084 
00085       vector<double> _object_px, _object_py, _object_pz, _object_p;
00086       vector<double> _object_pt, _object_e, _object_phi, _object_eta;
00087       vector<double> _event_shapes;
00088       vector<double> _thrust_axis, _thrust_axis_c;
00089 
00090       double _eta_c;
00091       int _irap;
00092       size_t _nmnjet;
00093 
00094     };
00095 
00096   }
00097 
00098 
00099 
00100 
00101   class CMS_2014_I1305624 : public Analysis {
00102   public:
00103 
00104     /// Constructor
00105     CMS_2014_I1305624()
00106       : Analysis("CMS_2014_I1305624")
00107     {    }
00108 
00109 
00110     /// @name Analysis methods
00111 
00112     /// Book histograms and initialise projections before the run
00113     void init() {
00114       const FastJets jets(FinalState(Cuts::abseta < 2.6), FastJets::ANTIKT, 0.5);
00115       declare(jets, "Jets");
00116 
00117       for (int ij=0; ij < NJETPTMN; ij++) {
00118         _h_thrustc[ij] = bookHisto1D(1, 1, ij+1);
00119         _h_broadt[ij] = bookHisto1D(1, 2, ij+1);
00120         _h_tot3dmass[ij] = bookHisto1D(1, 3, ij+1);
00121         _h_tottrnsmass[ij] = bookHisto1D(1, 4, ij+1);
00122         _h_y23c[ij] = bookHisto1D(1, 5, ij+1);
00123         //
00124         _alow1[ij] = _h_thrustc[ij]->xMin();
00125         _alow2[ij] = _h_broadt[ij]->xMin();
00126         _alow3[ij] = _h_tot3dmass[ij]->xMin();
00127         _alow4[ij] = _h_tottrnsmass[ij]->xMin();
00128         _alow5[ij] = _h_y23c[ij]->xMin();
00129         //
00130         _ahgh1[ij] = _h_thrustc[ij]->xMax();
00131         _ahgh2[ij] = _h_broadt[ij]->xMax();
00132         _ahgh3[ij] = _h_tot3dmass[ij]->xMax();
00133         _ahgh4[ij] = _h_tottrnsmass[ij]->xMax();
00134         _ahgh5[ij] = _h_y23c[ij]->xMax();
00135       }
00136     }
00137 
00138 
00139     /// Perform the per-event analysis
00140     void analyze(const Event& event) {
00141 
00142       const Jets& jets = apply<FastJets>(event, "Jets").jetsByPt(30.0*GeV);
00143       if (jets.size() < 2) vetoEvent;
00144       if (jets[0].abseta() > 2.4 || jets[1].abseta() > 2.4) vetoEvent;
00145 
00146       const double leadingpt = jets[0].pT();
00147       if (leadingpt < 110*GeV) vetoEvent;
00148 
00149       vector<double> jtpx, jtpy, jtpz, jten;
00150       foreach (const Jet& j, jets) {
00151         if (j.abseta() < 2.4) {
00152           jtpx.push_back(j.px());
00153           jtpy.push_back(j.py());
00154           jtpz.push_back(j.pz());
00155           jten.push_back(j.E());
00156         }
00157       }
00158 
00159       EventShape eventshape(jtpx, jtpy, jtpz, jten, 2.4, 0, 2);
00160       const vector<double> eventvar = eventshape.getEventShapes();
00161       if (eventvar[NEVTVAR] < 0) vetoEvent; // Jets are not only one hemisphere
00162 
00163       const double weight = event.weight();
00164       for (int ij = NJETPTMN-1; ij >= 0; --ij) {
00165         if (leadingpt/GeV > LEADINGPTTHRESHOLD[ij]) {
00166           if (inRange(eventvar[0], _alow1[ij], _ahgh1[ij])) _h_thrustc[ij]->fill(eventvar[0], weight);
00167           if (inRange(eventvar[2], _alow3[ij], _ahgh3[ij])) _h_tot3dmass[ij]->fill(eventvar[2], weight);
00168           if (inRange(eventvar[3], _alow4[ij], _ahgh4[ij])) _h_tottrnsmass[ij]->fill(eventvar[3], weight);
00169           if (eventvar[NEVTVAR] >= 3) {
00170             if (inRange(eventvar[1], _alow2[ij], _ahgh2[ij])) _h_broadt[ij]->fill(eventvar[1], weight);
00171             if (inRange(eventvar[4], _alow5[ij], _ahgh5[ij])) _h_y23c[ij]->fill(eventvar[4], weight);
00172           }
00173           break;
00174         }
00175       }
00176 
00177     }
00178 
00179 
00180     /// Normalise histograms etc., after the run
00181     void finalize() {
00182       for (int ij = 0; ij < NJETPTMN; ij++) {
00183         normalize(_h_thrustc[ij]);
00184         normalize(_h_broadt[ij]);
00185         normalize(_h_tot3dmass[ij]);
00186         normalize(_h_tottrnsmass[ij]);
00187         normalize(_h_y23c[ij]);
00188       }
00189     }
00190 
00191 
00192   private:
00193 
00194     /// @name Histograms
00195     //@{
00196     Histo1DPtr _h_thrustc[NJETPTMN];
00197     Histo1DPtr _h_broadt[NJETPTMN];
00198     Histo1DPtr _h_tot3dmass[NJETPTMN];
00199     Histo1DPtr _h_tottrnsmass[NJETPTMN];
00200     Histo1DPtr _h_y23c[NJETPTMN];
00201     //@}
00202 
00203     // Data members
00204     double _alow1[NJETPTMN], _alow2[NJETPTMN], _alow3[NJETPTMN], _alow4[NJETPTMN], _alow5[NJETPTMN];
00205     double _ahgh1[NJETPTMN], _ahgh2[NJETPTMN], _ahgh3[NJETPTMN], _ahgh4[NJETPTMN], _ahgh5[NJETPTMN];
00206 
00207   };
00208 
00209 
00210   DECLARE_RIVET_PLUGIN(CMS_2014_I1305624);
00211 
00212 
00213 
00214   /////////////////////
00215 
00216 
00217   namespace {
00218 
00219     // EventShape helper class method implementations:
00220 
00221     int EventShape::_calculate() {
00222       if (!_event_shapes.empty() && !_thrust_axis.empty() && !_thrust_axis_c.empty())
00223         return 1; //< return success if this appears to already have been run
00224 
00225       const size_t length = (size_t) _object_px.size();
00226 
00227       if (((size_t) _object_py.size() != length) ||
00228           ((size_t) _object_pz.size() != length) ||
00229           ((size_t) _object_e.size() != length)) {
00230         /// @todo Change to exception or assert
00231         // cout << "ERROR!!!! Input vectors differ in size! Change that please!" << endl;
00232         // cout<<"py_size: "<<_object_py.size()<<" ,pz_size: "<<_object_pz.size()
00233         //     <<" ,px_size: "<<_object_px.size()<<" ,E_size: "<<_object_e.size()<<endl;
00234         return 0;
00235       }
00236 
00237       if (!_object_p.empty()) {
00238         _object_p.clear();
00239         _object_pt.clear();
00240         _object_eta.clear();
00241         _object_phi.clear();
00242         _event_shapes.clear();
00243         _thrust_axis.clear();
00244         _thrust_axis_c.clear();
00245       }
00246 
00247       for (size_t j = 0; j < length; j++) {
00248         _object_p.push_back(0.);
00249         _object_pt.push_back(0.);
00250         _object_eta.push_back(0.);
00251         _object_phi.push_back(0.);
00252       }
00253 
00254       for (int j = 0; j < NEVTVAR; j++) {
00255         _event_shapes.push_back(-50.);
00256       }
00257 
00258       _event_shapes.push_back(double(_object_px.size())); //< WTF?
00259 
00260       for (size_t j = 0; j < 3; j++) {
00261         _thrust_axis.push_back(0.);
00262         _thrust_axis_c.push_back(0.);
00263       }
00264 
00265       double theta = 0;
00266 
00267       for (size_t k = 0; k < length; k++) {
00268         _object_p[k] = sqrt(pow(_object_px[k],2) + pow(_object_py[k],2) + pow(_object_pz[k],2));
00269         _object_pt[k] = sqrt(pow(_object_px[k],2) + pow(_object_py[k],2));
00270         if (_object_p[k] > _object_e[k] + 1e-4) {
00271           /// @todo Change to exception or assert
00272           // cout << "ERROR!!! object " << k <<" has P = " << _object_p[k]
00273           //      << " which is bigger than E = " << _object_e[k] <<" "
00274           //      << _object_px[k] <<" "<< _object_py[k] <<" "
00275           //      << _object_pz[k] <<" of total length "<< length
00276           //      << endl;
00277           return 0;
00278         }
00279 
00280         //to prevent a division by zero
00281         if (_irap == 0) {
00282           if (fabs(_object_pz[k]) > 1e-5) {
00283             theta = atan(_object_pt[k]/(_object_pz[k]));
00284           } else {
00285             theta = M_PI/2;
00286           }
00287           if (theta < 0.) theta = theta + M_PI;
00288           _object_eta[k] = -log(tan(0.5*theta));
00289         }
00290         if (_irap == 1) {
00291           if (_object_pz[k] == _object_e[k]) {
00292             /// @todo Change to exception or assert
00293             // cout << "ERROR!!! object "<<k<<" has Pz "<< _object_pz[k] <<" which is equal to E = "<< _object_e[k] <<endl;
00294             return 0;
00295           }
00296           _object_eta[k]=0.5*log((_object_e[k]+_object_pz[k])/(_object_e[k]-_object_pz[k]));
00297         }
00298         if (_irap != 0 && _irap != 1) {
00299           /// @todo Change to exception or assert
00300           // cout << "ERROR!!!, The choice to use the rapidity y or the pseudorapidity eta is not set correctly! Change that please!" << endl;
00301           return 0;
00302         }
00303         _object_phi[k] = atan2(_object_py[k], _object_px[k]);
00304       }
00305 
00306       vector<double> object_px_in, object_py_in, object_pz_in, object_pt_in, object_e_in, object_et_in, object_eta_in;
00307       vector<double> object_px_out, object_py_out, object_pz_out, object_e_out, object_pt_out, object_eta_out;
00308       if (!object_px_in.empty()) { //< FFS, this is impossible: it's only just been created!
00309         object_px_in.clear();
00310         object_py_in.clear();
00311         object_pz_in.clear();
00312         object_pt_in.clear();
00313         object_e_in.clear();
00314         object_et_in.clear();
00315         object_eta_in.clear();
00316         object_px_out.clear();
00317         object_py_out.clear();
00318         object_pz_out.clear();
00319         object_pt_out.clear();
00320         object_e_out.clear();
00321         object_eta_out.clear();
00322       }
00323 
00324       size_t nin = 0;
00325 
00326       for (size_t j = 0; j < length; j++) {
00327         if (fabs(_object_eta[j]) < _eta_c) {
00328           object_px_in.push_back(_object_px[j]);
00329           object_py_in.push_back(_object_py[j]);
00330           object_pz_in.push_back(_object_pz[j]);
00331           object_e_in.push_back(_object_e[j]);
00332           object_pt_in.push_back(sqrt(pow(_object_px[j],2)+pow(_object_py[j],2)));
00333           object_et_in.push_back(sqrt((pow(_object_e[j],2)*pow(_object_pt[j],2))/(pow(_object_pt[j],2)+pow(_object_pz[j],2))));
00334           object_eta_in.push_back(_object_eta[j]);
00335           nin += 1;
00336       } else {
00337           object_px_out.push_back(_object_px[j]);
00338           object_py_out.push_back(_object_py[j]);
00339           object_pz_out.push_back(_object_pz[j]);
00340           object_e_out.push_back(_object_e[j]);
00341           object_pt_out.push_back(sqrt(pow(_object_px[j],2)+pow(_object_py[j],2)));
00342           object_eta_out.push_back(_object_eta[j]);
00343         }
00344       }
00345 
00346       if (object_px_in.size() != nin) {
00347         /// @todo Change to exception or assert
00348         cout<<"ERROR!!! wrong dimension of 'in' momenta"<<endl;
00349         //return 0; ///< @todo Why not do this?
00350       }
00351       const size_t nout = length - nin;
00352 
00353       if (nin < _nmnjet) {
00354         for (int i = 0; i < NEVTVAR; i++) {
00355           _event_shapes[i] = -50.0;
00356         }
00357       }
00358 
00359       _event_shapes[NEVTVAR] = nin;
00360 
00361       if (nin >= _nmnjet) {
00362         double p_sum_c = 0; //GMA
00363         double pt_sum_c = 0;
00364         double eta_cw=0;
00365         double px_sum_in = 0;
00366         double py_sum_in = 0;
00367         for (size_t j = 0; j < nin; j++) {
00368           pt_sum_c += object_pt_in[j];
00369           p_sum_c += sqrt(pow(object_pt_in[j],2.) + pow(object_pz_in[j], 2.0)); //GMA
00370           eta_cw += object_pt_in[j]*object_eta_in[j];
00371           px_sum_in += object_px_in[j];
00372           py_sum_in += object_py_in[j];
00373         }
00374         eta_cw /= pt_sum_c;
00375 
00376         double expTerm = 0;
00377         for (size_t j = 0; j < nout; j++) {
00378           expTerm += object_pt_out[j] * exp(-fabs(object_eta_out[j]-eta_cw));
00379         }
00380         expTerm /= pt_sum_c;
00381 
00382         //the central global transverse thrust centrthr is calculated
00383         double centrthr = 0;
00384         vector<double> thrust_central = _thrust(object_px_in, object_py_in);
00385 
00386         for (size_t l=0; l<3; l++) _thrust_axis_c[l] = thrust_central[l];
00387         //the variable which gets resummed is not thrust
00388         //but tau = 1 - thrust - see calculation
00389         centrthr = thrust_central[3];
00390         _event_shapes[0] = centrthr;
00391 
00392         double alpha_c = atan2(_thrust_axis_c[1], _thrust_axis_c[0]);
00393         //central jet masses
00394         //define two jet masses in region U and D
00395         double cenjm_up = 0;
00396         double cenjm_down= 0;
00397         double dot_product = 0;
00398 
00399         vector<double> up_sum;
00400         vector<double> down_sum;
00401         for (size_t j=0; j<4;j++) {
00402           up_sum.push_back(0.);
00403           down_sum.push_back(0.);
00404         }
00405         for (size_t i=0;i<nin;i++) {
00406           dot_product = object_px_in[i] * _thrust_axis_c[0] + object_py_in[i] * _thrust_axis_c[1];
00407           if (dot_product >= 0) {
00408             up_sum[0]+=object_px_in[i];
00409             up_sum[1]+=object_py_in[i];
00410             up_sum[2]+=object_pz_in[i];
00411             up_sum[3]+=object_e_in[i];
00412           } else {
00413             down_sum[0]+=object_px_in[i];
00414             down_sum[1]+=object_py_in[i];
00415             down_sum[2]+=object_pz_in[i];
00416             down_sum[3]+=object_e_in[i];
00417           }
00418         }
00419         cenjm_up = _lorentz_sp(up_sum, up_sum) / pow(p_sum_c, 2.); //GMA pow(pt_sum_c,2);
00420         cenjm_down = _lorentz_sp(down_sum, down_sum) / pow(p_sum_c, 2.); //GMA pow(pt_sum_c,2);
00421 
00422         //central total jet mass centotjm
00423         double centotjm=0;
00424         centotjm = cenjm_up + cenjm_down;
00425 
00426         _event_shapes[2]=centotjm;
00427 
00428         double centrjm_up=0, centrjm_down=0;
00429         vector<double> upsum;
00430         vector<double> downsum;
00431         for (size_t j = 0; j < 3; j++) {
00432           upsum.push_back(0.);
00433           downsum.push_back(0.);
00434         }
00435         for (size_t i = 0; i < nin; i++) {
00436           dot_product = object_px_in[i]*_thrust_axis_c[0]+object_py_in[i]*_thrust_axis_c[1];
00437           if (dot_product >= 0) {
00438             upsum[0] += object_px_in[i];
00439             upsum[1] += object_py_in[i];
00440             upsum[2] += object_et_in[i];
00441           } else {
00442             downsum[0] += object_px_in[i];
00443             downsum[1] += object_py_in[i];
00444             downsum[2] += object_et_in[i];
00445           }
00446         }
00447         centrjm_up = _lorentz_sp(upsum, upsum) / pow(pt_sum_c, 2);
00448         centrjm_down = _lorentz_sp(downsum, downsum) / pow(pt_sum_c, 2);
00449         double centottrjm = centrjm_up + centrjm_down;
00450 
00451         _event_shapes[3] = centottrjm;
00452 
00453         //central three-jet resolution threshold
00454         double ceny3=0;
00455         if (nin < 3) {
00456           ceny3 = -1.0;
00457         } else {
00458           ceny3 = _three_jet_res(object_px_in, object_py_in, object_pz_in, object_e_in, _irap);
00459         }
00460 
00461         _event_shapes[4] = ceny3;
00462 
00463         //the central jet broadenings in the up and down region
00464         double cenbroad_up=0;
00465         double cenbroad_down=0;
00466 
00467         double eta_up=0;
00468         size_t num_up=0;
00469         double eta_down =0;
00470         size_t num_down =0;
00471         double phi_temp =0;
00472         double phi_up_aver =0;
00473         double phi_down_aver =0;
00474         double pt_sum_up =0;
00475         double pt_sum_down =0;
00476         double dot_product_b =0;
00477         vector<double> phi_up;
00478         vector<double> phi_down;
00479         double py_rot =0;
00480         double px_rot =0;
00481 
00482         for (size_t j = 0; j < 4; j++) {
00483           up_sum.push_back(0.);
00484           down_sum.push_back(0.);
00485         }
00486 
00487         for (size_t i=0;i<nin;i++) {
00488           dot_product_b =sqrt(object_px_in[i]*_thrust_axis_c[0] + object_py_in[i]*_thrust_axis_c[1]);
00489           if (dot_product_b>=0){
00490             pt_sum_up += object_pt_in[i];
00491             //rotate the coordinate system so that
00492             //the central thrust axis is e_x
00493             px_rot = cos(alpha_c)*object_px_in[i]+sin(alpha_c)*object_py_in[i];
00494             py_rot = - sin(alpha_c)*object_px_in[i]+cos(alpha_c)*object_py_in[i];
00495             //calculate the eta and phi in the rotated system
00496             eta_up += object_pt_in[i]*object_eta_in[i];
00497             phi_temp = atan2(py_rot,px_rot);
00498 
00499             if(phi_temp > M_PI/2){
00500               phi_temp = phi_temp - M_PI/2;
00501             }
00502             if (phi_temp < -M_PI/2){
00503               phi_temp = phi_temp + M_PI/2;
00504             }
00505             phi_up.push_back(phi_temp);
00506             phi_up_aver += object_pt_in[i]*phi_temp;
00507             num_up += 1;
00508           } else {
00509             eta_down += object_pt_in[i]*object_eta_in[i];
00510             pt_sum_down += object_pt_in[i];
00511             px_rot = cos(alpha_c)*object_px_in[i]+sin(alpha_c)*object_py_in[i];
00512             py_rot = - sin(alpha_c)*object_px_in[i]+cos(alpha_c)*object_py_in[i];
00513             phi_temp = atan2(py_rot,px_rot);
00514             if (phi_temp > M_PI/2) {
00515               //if phi is bigger than pi/2 in the new system calculate
00516               //the difference to the thrust axis
00517               phi_temp = M_PI -phi_temp;
00518             }
00519             if (phi_temp<-M_PI/2) {
00520               //if phi is smaller than
00521               phi_temp = -M_PI-phi_temp;
00522             }
00523             phi_down.push_back(phi_temp);
00524             //calculate the pt-weighted phi
00525             phi_down_aver += object_pt_in[i]*phi_temp;
00526             num_down += 1;
00527           }
00528         }
00529         if (num_up!=0){
00530           eta_up = eta_up/pt_sum_up;
00531           phi_up_aver = phi_up_aver/pt_sum_up;
00532         }
00533         if (num_down!=0) {
00534           eta_down = eta_down/pt_sum_down;
00535           phi_down_aver = phi_down_aver/pt_sum_down;
00536         }
00537 
00538         size_t index_up=0, index_down=0;
00539         for (size_t i = 0; i < nin; i++) {
00540           dot_product_b = object_px_in[i]*_thrust_axis_c[0] + object_py_in[i]*_thrust_axis_c[1];
00541           if (dot_product_b >= 0) {
00542             //calculate the broadenings of the regions with the rotated system
00543             //and the pt-weighted average of phi in the rotated system
00544             cenbroad_up += object_pt_in[i]*sqrt(pow(object_eta_in[i]-eta_up, 2) +
00545                                                 pow(_delta_phi(phi_up[index_up], phi_up_aver), 2));
00546             index_up += 1;
00547           } else {
00548             cenbroad_down += object_pt_in[i]*sqrt(pow(object_eta_in[i]-eta_down, 2)+
00549                                                   pow(_delta_phi(phi_down[index_down], phi_down_aver), 2));
00550             index_down += 1;
00551           }
00552         }
00553 
00554         if (index_up == 0 || index_down ==0) _event_shapes[NEVTVAR] *= -1;
00555 
00556         cenbroad_up=cenbroad_up/(2*pt_sum_c);
00557         cenbroad_down=cenbroad_down/(2*pt_sum_c);
00558 
00559         //central total jet broadening
00560         double centotbroad = 0;
00561         centotbroad = cenbroad_up + cenbroad_down;
00562 
00563         _event_shapes[1] = centotbroad;
00564 
00565         for (int ij = 0; ij < 5; ij++) {
00566           if (_event_shapes[ij] < 1.e-20) _event_shapes[ij] = 1.e-20;
00567           _event_shapes[ij] = log(_event_shapes[ij]);
00568         }
00569       }
00570 
00571       return 1;
00572     }
00573 
00574 
00575     double EventShape::_three_jet_res(const vector<double>& in_object_px, const vector<double>& in_object_py, const vector<double>& in_object_pz, const vector<double>& in_object_e, int irap) {
00576 
00577       size_t y3_length = (size_t)in_object_px.size();
00578       if (((size_t) in_object_py.size()!=y3_length) ||
00579           ((size_t) in_object_pz.size()!=y3_length) ||
00580           (in_object_e.size()!=y3_length)) {
00581         // cout << "ERROR!!!! Input vectors differ in size! Change that please!" << endl;
00582         // cout<<"py_size: "<<in_object_py.size()<<" ,pz_size: "<<in_object_pz.size()
00583         //     <<" ,px_size: "<<in_object_px.size()<<" , E_size: "<<in_object_e.size() <<endl;
00584         return 0.0;
00585       }
00586 
00587       vector<double> in_object_p, in_object_pt, in_object_eta, in_object_phi;
00588       if (!in_object_p.empty()) {
00589         in_object_p.clear();
00590         in_object_pt.clear();
00591         in_object_eta.clear();
00592         in_object_phi.clear();
00593       }
00594       for (size_t j = 0; j < y3_length; j++) {
00595         in_object_p.push_back(0.);
00596         in_object_pt.push_back(0.);
00597         in_object_eta.push_back(0.);
00598         in_object_phi.push_back(0.);
00599       }
00600       double theta_y3_1st = 0;
00601       for (size_t k =0; k<y3_length; k++) {
00602         in_object_p[k] = sqrt(pow(in_object_px[k],2) + pow(in_object_py[k],2) + pow(in_object_pz[k],2));
00603         in_object_pt[k] = sqrt(pow(in_object_px[k],2) + pow(in_object_py[k],2));
00604 
00605         //calculates the pseudorapidity to prevent a division by zero
00606         if (irap == 0) {
00607           if (fabs(in_object_pz[k]) > 1E-5) {
00608             theta_y3_1st = atan(in_object_pt[k]/(in_object_pz[k]));
00609           } else {
00610             theta_y3_1st = M_PI/2;
00611           }
00612           if (theta_y3_1st<0.) theta_y3_1st = theta_y3_1st + M_PI;
00613           in_object_eta[k] = - log(tan(0.5*theta_y3_1st));
00614         }
00615         //calculates the real rapidity
00616         if (irap == 1) {
00617           in_object_eta[k]=0.5*log((in_object_e[k]+in_object_pz[k])/(in_object_e[k]-in_object_pz[k]));
00618         }
00619         in_object_phi[k] = atan2(in_object_py[k], in_object_px[k]);
00620       }
00621 
00622       //the three-jet resolution
00623       //threshold y3
00624       double y3 = 0;
00625 
00626       //vector which will be filled with the
00627       //minimum of the distances
00628       double max_dmin_temp=0;
00629 
00630       double max_dmin = 0;
00631 
00632       //distance input object k, beam
00633       double distance_jB = 0;
00634       double distance_jB_min = 0;
00635       //distance of input object k to l
00636       double distance_jk = 0;
00637       double distance_jk_min = 0;
00638       //as we search the minimum of the distances
00639       //give them values which are for sure higher
00640       //than those we evaluate first in the for-loups
00641 
00642       size_t index_jB = 0;
00643       size_t index_j_jk = 0;
00644       size_t index_k_jk = 0;
00645 
00646       //to decide later if the minmum is a jB or jk
00647       int decide_jB = -1;
00648 
00649       vector<double> input_pt, input_px, input_py, input_pz;
00650       vector<double> input_p, input_e, input_phi, input_eta;
00651 
00652       if (!input_pt.empty()) {
00653         input_pt.clear();
00654         input_px.clear();
00655         input_px.clear();
00656         input_pz.clear();
00657         input_p.clear();
00658         input_e.clear();
00659         input_phi.clear();
00660         input_eta.clear();
00661       }
00662 
00663       for (size_t j = 0; j < y3_length; j++){
00664         input_pt.push_back(in_object_pt[j]);
00665         input_px.push_back(in_object_px[j]);
00666         input_py.push_back(in_object_py[j]);
00667         input_pz.push_back(in_object_pz[j]);
00668         input_p.push_back(in_object_p[j]);
00669         input_e.push_back(in_object_e[j]);
00670         input_phi.push_back(in_object_phi[j]);
00671         input_eta.push_back(in_object_eta[j]);
00672       }
00673       if (y3_length<3) {
00674         return -1;
00675       } else {
00676         size_t rest = y3_length;
00677         for (size_t i = 0; i<y3_length; i++) {
00678           //make the minima at the initialization step
00679           //of each looping bigger than the first values
00680           distance_jB_min = 0.36*pow(input_pt[0],2) + 10;
00681           //DELTA PHIs wanted not the pure difference
00682           distance_jk_min = min(pow(input_pt[1], 2), pow(input_pt[0], 2)) *
00683             (pow(input_eta[1]-input_eta[0], 2) +
00684              pow(_delta_phi(input_phi[1], input_phi[0]), 2)) + 10;
00685           //do the procedure only until we have only 2 objects left anymore
00686           if (rest > 2) {
00687             for (size_t j=0; j<rest;j++) {
00688               //calculate the distance between object j and the beam
00689               distance_jB = 0.36*pow(input_pt[j], 2);
00690               if(distance_jB < distance_jB_min){
00691                 distance_jB_min = distance_jB;
00692                 index_jB = j;
00693               }
00694               if (j > 0) {
00695                 for(size_t k=0; k<j;k++){
00696                   //calculate the distance in delta eta and delta phi between object i and object j
00697                   distance_jk = min(pow(input_pt[j], 2),pow(input_pt[k], 2))*
00698                     (pow(input_eta[j]-input_eta[k], 2)+
00699                      pow(_delta_phi(input_phi[j],input_phi[k]), 2));
00700                   if (distance_jk<distance_jk_min) {
00701                     distance_jk_min = distance_jk;
00702                     index_j_jk = j;
00703                     index_k_jk =k;
00704                   }
00705                 }
00706               }
00707             }
00708             //decide if the minimum is from a jB or jk combination
00709             if (distance_jk_min<distance_jB_min) {
00710               max_dmin_temp = max(distance_jk_min,max_dmin_temp);
00711               decide_jB = 0;
00712             } else {
00713               max_dmin_temp = max(distance_jB_min,max_dmin_temp);
00714               decide_jB=1;
00715             }
00716             //if we have only three jets left calculate
00717             //the maxima of the dmin's
00718             //if the minimum is a jB eliminate the input object
00719             if (decide_jB == 1) {
00720               //if index_jB is the last one nothing is to do
00721               if (index_jB != rest-1) {
00722                 for (size_t i=index_jB; i<rest-1;i++) {
00723                   input_pt[i]=input_pt[i+1];
00724                   input_phi[i]=input_phi[i+1];
00725                   input_eta[i]=input_eta[i+1];
00726                   input_px[i]=input_px[i+1];
00727                   input_py[i]=input_py[i+1];
00728                   input_pz[i]=input_pz[i+1];
00729                   input_e[i]=input_e[i+1];
00730                 }
00731               }
00732             }
00733             //if the minimum is a jk combine both input objects
00734             if(decide_jB==0) {
00735               input_px[index_k_jk] = input_px[index_k_jk]+input_px[index_j_jk];
00736               input_py[index_k_jk] = input_py[index_k_jk]+input_py[index_j_jk];
00737               input_pz[index_k_jk] = input_pz[index_k_jk]+input_pz[index_j_jk];
00738               input_e[index_k_jk] = input_e[index_k_jk]+input_e[index_j_jk];
00739               input_p[index_k_jk] = sqrt(pow(input_px[index_k_jk], 2)+
00740                                          pow(input_py[index_k_jk], 2)+
00741                                          pow(input_pz[index_k_jk], 2));
00742               //calculate the pt, eta and phi of the new combined momenta k_jk
00743               input_pt[index_k_jk] = sqrt(pow(input_px[index_k_jk], 2)+
00744                                           pow(input_py[index_k_jk], 2));
00745               //in the case of pseudorapidity
00746               if (irap == 0) {
00747                 double theta_new =0;
00748                 if (fabs(input_pz[index_k_jk]) > 1E-5){
00749                   theta_new = atan(input_pt[index_k_jk]/(input_pz[index_k_jk]));
00750                 } else {
00751                   theta_new = M_PI/2;
00752                 }
00753                 if (theta_new < 0) {
00754                   theta_new = theta_new + M_PI;
00755                 }
00756                 input_eta[index_k_jk] = - log(tan(0.5*theta_new));
00757               }
00758               //in the real rapidity y is wanted
00759               if (irap == 1) {
00760                 input_eta[index_k_jk] = 0.5 * log((input_e[index_k_jk]+
00761                                                    input_pz[index_k_jk]) /
00762                                                   (input_e[index_k_jk] -
00763                                                    input_pz[index_k_jk]));
00764               }
00765               input_phi[index_k_jk] = atan2(input_py[index_k_jk], input_px[index_k_jk]);
00766               if (index_j_jk != rest-1) {
00767                 for (size_t i = index_j_jk; i<rest-1;i++) {
00768                   input_pt[i] = input_pt[i+1];
00769                   input_phi[i] = input_phi[i+1];
00770                   input_eta[i] = input_eta[i+1];
00771                   input_px[i] = input_px[i+1];
00772                   input_py[i] = input_py[i+1];
00773                   input_pz[i] = input_pz[i+1];
00774                   input_e[i] = input_e[i+1];
00775                 }
00776               }
00777             }
00778           }
00779           if (rest == 3) max_dmin = max_dmin_temp;
00780           rest = rest-1;
00781         }
00782       }
00783 
00784       double et2 = 0;
00785       et2 = input_pt[0] + input_pt[1];
00786       y3 = max_dmin/pow(et2,2);
00787 
00788       return y3;
00789     }
00790 
00791 
00792     vector<double> EventShape::_thrust(const vector<double>& input_px, const vector<double>& input_py) {
00793 
00794       double thrustmax_calc = 0;
00795       double temp_calc = 0;
00796       size_t length_thrust_calc = 0;
00797       vector<double> thrust_values, thrust_axis_calc;
00798       vector<double> p_thrust_max_calc, p_dec_1_calc,  p_dec_2_calc, p_pt_beam_calc;
00799 
00800       if (!thrust_values.empty()){
00801         thrust_values.clear();
00802         thrust_axis_calc.clear();
00803         p_thrust_max_calc.clear();
00804         p_dec_1_calc.clear();
00805         p_dec_2_calc.clear();
00806         p_pt_beam_calc.clear();
00807       }
00808 
00809       for (size_t j = 0; j < 3; j++){
00810         p_pt_beam_calc.push_back(0.);
00811         p_dec_1_calc.push_back(0.);
00812         p_dec_2_calc.push_back(0.);
00813         p_thrust_max_calc.push_back(0.);
00814         thrust_axis_calc.push_back(0.);
00815       }
00816 
00817       for (size_t j = 0; j < 4; j++) {
00818         thrust_values.push_back(0.);
00819       }
00820 
00821       length_thrust_calc = input_px.size();
00822       if (input_py.size() != length_thrust_calc) {
00823         /// @todo Change to exception or assert
00824         cout<<"ERROR in thrust calculation!!! Size of input vectors differs. Change that please!"<<endl;
00825         return thrust_values;
00826       }
00827 
00828       double pt_sum_calc =0;
00829       for(size_t k=0;k<length_thrust_calc;k++){
00830         pt_sum_calc+=sqrt(pow(input_px[k],2)+pow(input_py[k],2));
00831         for(size_t j = 0; j < 3; j++){
00832           p_thrust_max_calc[j]=0;
00833         }
00834         //get a vector perpendicular to the beam axis and
00835         //perpendicular to the momentum of particle k
00836         //per default beam axis b = (0,0,1)
00837         p_pt_beam_calc[0] = input_py[k]*1;
00838         p_pt_beam_calc[1] = - input_px[k]*1;
00839         p_pt_beam_calc[2] = 0.; // GMA p_pt_beam_calc[3] = 0.;
00840         for(size_t i=0;i<length_thrust_calc;i++){
00841           if(i!=k){
00842             if((input_px[i]*p_pt_beam_calc[0]+input_py[i]*p_pt_beam_calc[1])>=0){
00843               p_thrust_max_calc[0]= p_thrust_max_calc[0] + input_px[i];
00844               p_thrust_max_calc[1]= p_thrust_max_calc[1] + input_py[i];
00845             }else{
00846               p_thrust_max_calc[0]= p_thrust_max_calc[0] - input_px[i];
00847               p_thrust_max_calc[1]= p_thrust_max_calc[1] - input_py[i];
00848             }
00849           }
00850         }
00851         p_dec_1_calc[0] = p_thrust_max_calc[0] + input_px[k];
00852         p_dec_1_calc[1] = p_thrust_max_calc[1] + input_py[k];
00853         p_dec_1_calc[2] = 0;
00854         p_dec_2_calc[0] = p_thrust_max_calc[0] - input_px[k];
00855         p_dec_2_calc[1] = p_thrust_max_calc[1] - input_py[k];
00856         p_dec_2_calc[2] = 0;
00857         temp_calc = pow(p_dec_1_calc[0], 2) + pow(p_dec_1_calc[1], 2);
00858 
00859         if (temp_calc>thrustmax_calc) {
00860           thrustmax_calc =temp_calc;
00861           for (size_t i=0; i<3; i++) {
00862             thrust_axis_calc[i] = p_dec_1_calc[i]/sqrt(thrustmax_calc);
00863           }
00864         }
00865         temp_calc = pow(p_dec_2_calc[0], 2)+pow(p_dec_2_calc[1], 2);
00866         if (temp_calc > thrustmax_calc) {
00867           thrustmax_calc =temp_calc;
00868           for (size_t i=0; i<3; i++) {
00869             thrust_axis_calc[i] = p_dec_2_calc[i]/sqrt(thrustmax_calc);
00870           }
00871         }
00872       }
00873       for (size_t j = 0; j < 3; j++) thrust_values[j] = thrust_axis_calc[j];
00874       const double thrust_calc = sqrt(thrustmax_calc)/pt_sum_calc;
00875 
00876       // the variable which gets returned is not the thrust but tau=1-thrust
00877       thrust_values[3] = 1 - thrust_calc;
00878       if (thrust_values[3] < 1e-20) thrust_values[3] = 1e-20;
00879 
00880       return thrust_values;
00881     }
00882 
00883 
00884     double EventShape::_delta_phi(double phi1, double phi2) {
00885       double dphi = fabs(phi2 - phi1);
00886       if (dphi > M_PI) dphi = 2*M_PI - dphi;
00887       return dphi;
00888     }
00889 
00890 
00891     // Returns the scalar product between two 4 momenta
00892     double EventShape::_lorentz_sp(const vector<double>& a, const vector<double>& b) {
00893       size_t dim = (size_t) a.size();
00894       if (a.size()!=b.size()) {
00895         cout<<"ERROR!!! Dimension of input vectors are different! Change that please!"<<endl;
00896         return 0;
00897       } else {
00898         double l_dot_product=a[dim-1]*b[dim-1];
00899         for(size_t i=0; i<dim-1;i++){
00900           l_dot_product-=a[i]*b[i];
00901         }
00902         return l_dot_product;
00903       }
00904     }
00905 
00906 
00907   }
00908 
00909 
00910 }