FastJet 3.0beta1
JHTopTagger.cc
00001 //STARTHEADER
00002 // $Id: JHTopTagger.cc 2493 2011-08-03 15:48:16Z salam $
00003 //
00004 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin Salam and Gregory Soyez
00005 //
00006 //----------------------------------------------------------------------
00007 // This file is part of FastJet.
00008 //
00009 //  FastJet is free software; you can redistribute it and/or modify
00010 //  it under the terms of the GNU General Public License as published by
00011 //  the Free Software Foundation; either version 2 of the License, or
00012 //  (at your option) any later version.
00013 //
00014 //  The algorithms that underlie FastJet have required considerable
00015 //  development and are described in hep-ph/0512210. If you use
00016 //  FastJet as part of work towards a scientific publication, please
00017 //  include a citation to the FastJet paper.
00018 //
00019 //  FastJet is distributed in the hope that it will be useful,
00020 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00021 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00022 //  GNU General Public License for more details.
00023 //
00024 //  You should have received a copy of the GNU General Public License
00025 //  along with FastJet; if not, write to the Free Software
00026 //  Foundation, Inc.:
00027 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00028 //----------------------------------------------------------------------
00029 //ENDHEADER
00030 
00031 #include <fastjet/tools/JHTopTagger.hh>
00032 #include <fastjet/Error.hh>
00033 #include <fastjet/JetDefinition.hh>
00034 #include <fastjet/ClusterSequence.hh>
00035 #include <sstream>
00036 #include <limits>
00037 
00038 FASTJET_BEGIN_NAMESPACE
00039 
00040 using namespace std;
00041 
00042 //----------------------------------------------------------------------
00043 // JHTopTagger class implementation
00044 //----------------------------------------------------------------------
00045 
00046 LimitedWarning JHTopTagger::_warnings_nonca;
00047 
00048 //------------------------------------------------------------------------
00049 // description of the tagger
00050 string JHTopTagger::description() const{ 
00051   ostringstream oss;
00052   oss << "JHTopTagger with delta_p=" << _delta_p << ", delta_r=" << _delta_r
00053       << ", cos_theta_W_max=" << _cos_theta_W_max
00054       << " and mW = " << _mW;
00055   oss << description_of_selectors();
00056   return oss.str();
00057 }
00058 
00059 //------------------------------------------------------------------------
00060 // returns the tagged PseudoJet if successful, 0 otherwise
00061 //  - jet   the PseudoJet to tag
00062 PseudoJet JHTopTagger::result(const PseudoJet & jet) const{
00063   // make sure that there is a "regular" cluster sequence associated
00064   // with the jet. Note that we also check it is valid (to avoid a
00065   // more criptic error later on)
00066   if (!jet.has_valid_cluster_sequence()){
00067     throw Error("JHTopTagger can only be applied on jets having an associated (and valid) ClusterSequence");
00068   }
00069 
00070   // warn if the jet has not been clustered with a Cambridge/Aachen
00071   // algorithm
00072   if (! jet.validated_cs()->jet_def().jet_algorithm() == cambridge_algorithm)
00073     _warnings_nonca.warn("JHTopTagger should only be applied on jets from a Cambridge/Aachen clustering; use it with other algorithms at your own risk.");
00074 
00075 
00076   // do the first splitting
00077   vector<PseudoJet> split0 = _split_once(jet, jet);
00078   if (split0.size() == 0) return PseudoJet();
00079 
00080   // now try a second splitting on each of the resulting objects
00081   vector<PseudoJet> subjets;
00082   for (unsigned i = 0; i < 2; i++) {
00083     vector<PseudoJet> split1 = _split_once(split0[i], jet);
00084     if (split1.size() > 0) {
00085       subjets.push_back(split1[0]);
00086       subjets.push_back(split1[1]);
00087     } else {
00088       subjets.push_back(split0[i]);
00089     }
00090   }
00091 
00092   // make sure things make sense
00093   if (subjets.size() < 3) return PseudoJet();
00094 
00095   // now find the pair of objects closest in mass to the W
00096   double dmW_min = numeric_limits<double>::max();
00097   int ii=-1, jj=-1;
00098   for (unsigned i = 0 ; i < subjets.size()-1; i++) {
00099     for (unsigned j = i+1 ; j < subjets.size(); j++) {
00100       double dmW = abs(_mW - (subjets[i]+subjets[j]).m());
00101       if (dmW < dmW_min) {
00102         dmW_min = dmW; ii = i; jj = j;
00103       }
00104     }
00105   }
00106 
00107   // order the subjets in the following order:
00108   //  - hardest of the W subjets
00109   //  - softest of the W subjets
00110   //  - hardest of the remaining subjets
00111   //  - softest of the remaining subjets (if any)
00112   if (ii>0) std::swap(subjets[ii], subjets[0]);
00113   if (jj>1) std::swap(subjets[jj], subjets[1]);
00114   if (subjets[0].perp2() < subjets[1].perp2()) std::swap(subjets[0], subjets[1]);
00115   if ((subjets.size()>3) && (subjets[2].perp2() < subjets[3].perp2())) 
00116     std::swap(subjets[2], subjets[3]);
00117   
00118   // create the result and its structure
00119   const JetDefinition::Recombiner *rec
00120     = jet.associated_cluster_sequence()->jet_def().recombiner();
00121 
00122   PseudoJet W = join(subjets[0], subjets[1], *rec);
00123   PseudoJet non_W;
00124   if (subjets.size()>3) {
00125     non_W = join(subjets[2], subjets[3], *rec);
00126   } else {
00127     non_W = join(subjets[2], *rec);
00128   }
00129   PseudoJet result = join<JHTopTaggerStructure>(W, non_W, *rec);
00130   JHTopTaggerStructure *s = (JHTopTaggerStructure*) result.structure_non_const_ptr();
00131   s->_cos_theta_w = _cos_theta_W(result);
00132 
00133   // if the polarisation angle does not pass the cut, consider that
00134   // the tagging has failed
00135   //
00136   // Note that we could perhaps ensure this cut before constructing
00137   // the result structure but this has the advantage that the top
00138   // 4-vector is already available and does not have to de re-computed
00139   if (s->_cos_theta_w >= _cos_theta_W_max ||
00140       ! _top_selector.pass(result) || ! _W_selector.pass(W)
00141       ) {
00142     result *= 0.0;
00143   }
00144 
00145   return result;
00146 
00147   // // old version
00148   // PseudoJet result = join<JHTopTaggerStructure>(subjets, *rec);
00149   // JHTopTaggerStructure *s = (JHTopTaggerStructure*) result.structure_non_const_ptr();
00150   // //  s->_original_jet = jet;
00151   // s->_W = join(subjets[0], subjets[1], *rec);
00152   // if (subjets.size()>3)
00153   //   s->_non_W = join(subjets[2], subjets[3], *rec);
00154   // else
00155   //   s->_non_W = join(subjets[2], *rec);
00156   // s->_cos_theta_w = _cos_theta_W(result);
00157   // 
00158   // // if the polarisation angle does not pass the cut, consider that
00159   // // the tagging has failed
00160   // //
00161   // // Note that we could perhaps ensure this cut before constructing
00162   // // the result structure but this has the advantage that the top
00163   // // 4-vector is already available and does not have to de re-computed
00164   // if (s->_cos_theta_w >= _cos_theta_W_max)
00165   //   return PseudoJet();
00166   // 
00167   // return result;
00168 }
00169 
00170 // runs the Johns Hopkins decomposition procedure
00171 vector<PseudoJet> JHTopTagger::_split_once(const PseudoJet & jet_to_split,
00172                                            const PseudoJet & reference_jet) const{
00173   PseudoJet this_jet = jet_to_split;
00174   PseudoJet p1, p2;
00175   vector<PseudoJet> result;
00176   while (this_jet.has_parents(p1, p2)) {
00177     if (p2.perp2() > p1.perp2()) std::swap(p1,p2); // order with hardness
00178     if (p1.perp() < _delta_p * reference_jet.perp()) break; // harder is too soft wrt original jet
00179     if ( (abs(p2.rap()-p1.rap()) + abs(p2.delta_phi_to(p1))) < _delta_r) break; // distance is too small
00180     if (p2.perp() < _delta_p * reference_jet.perp()) {
00181       this_jet = p1; // softer is too soft wrt original, so ignore it
00182       continue; 
00183     }
00184     //result.push_back(this_jet);
00185     result.push_back(p1);
00186     result.push_back(p2);
00187     break;
00188   }
00189   return result;
00190 }
00191 
00192 
00193 // compute the W helicity angle
00194 //
00195 // The helicity angle is a standard observable in top decays, used to
00196 // determine the Lorentz structure of the top- W coupling [13]. It is
00197 // defined as the angle, measured in the rest frame of the
00198 // reconstructed W, between the reconstructed top's flight direction
00199 // and one of the W decay products. Normally, it is studied in
00200 // semi-leptonic top decays, where the charge of the lepton uniquely
00201 // identifies these decay products. In hadronic top decays there is an
00202 // ambiguity which we resolve by choosing the lower pT subjet, as
00203 // measured in the lab frame.
00204 double JHTopTagger::_cos_theta_W(const PseudoJet & result) const{
00205   // the two jets of interest: top and lower-pt prong of W
00206   const PseudoJet & W  = result.structure_of<JHTopTagger>().W();
00207   vector<PseudoJet> W_pieces = W.pieces();
00208   assert(W_pieces[0].perp2() >= W_pieces[1].perp2());
00209   PseudoJet W2  = W_pieces[1];
00210   PseudoJet top = result;
00211   
00212   // transform these jets into jets in the rest frame of the W
00213   W2.unboost(W);
00214   top.unboost(W);
00215 
00216   return (W2.px()*top.px() + W2.py()*top.py() + W2.pz()*top.pz())/
00217     sqrt(W2.modp2() * top.modp2());
00218 }
00219 
00220 
00221 FASTJET_END_NAMESPACE
00222 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends