FastJet 3.0beta1
CASubJetTagger.cc
00001 //STARTHEADER
00002 // $Id: CASubJetTagger.cc 2474 2011-07-26 13:13:54Z cacciari $
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/CASubJetTagger.hh>
00032 #include <fastjet/ClusterSequence.hh>
00033 
00034 #include <algorithm>
00035 #include <cmath>
00036 #include <sstream>
00037 
00038 using namespace std;
00039 
00040 FASTJET_BEGIN_NAMESPACE
00041 
00042 
00043 LimitedWarning CASubJetTagger::_non_ca_warnings;
00044 
00045 // the tagger's description
00046 //----------------------------------------------------------------------
00047 string CASubJetTagger::description() const{
00048   ostringstream oss;
00049   oss << "CASubJetTagger with z_threshold=" << _z_threshold ;
00050   if (_absolute_z_cut) oss << " (defined wrt original jet)";
00051   oss << " and scale choice ";
00052   switch (_scale_choice) {
00053   case kt2_distance:         oss << "kt2_distance";         break;
00054   case jade_distance:        oss << "jade_distance";        break;
00055   case jade2_distance:       oss << "jade2_distance";       break;
00056   case plain_distance:       oss << "plain_distance";       break;
00057   case mass_drop_distance:   oss << "mass_drop_distance";   break;
00058   case dot_product_distance: oss << "dot_product_distance"; break;
00059   default:
00060     throw fastjet::Error("unrecognized scale choice");
00061   }
00062 
00063   return oss.str();
00064 }
00065 
00066 // run the tagger on the given cs/jet
00067 // returns the tagged PseudoJet if successful, 0 otherwise
00068 //----------------------------------------------------------------------
00069 PseudoJet CASubJetTagger::result(const fastjet::PseudoJet & jet) const{
00070   // make sure that the jet results from a Cambridge/Aachen clustering
00071   if (jet.validated_cs()->jet_def().jet_algorithm() != cambridge_algorithm)
00072     _non_ca_warnings.warn("CASubJetTagger should only be applied on jets from a Cambridge/Aachen clustering; use it with other algorithms at your own risk");
00073 
00074   // recurse in the jet to find the max distance
00075   JetAux aux;
00076   aux.jet          = PseudoJet();
00077   aux.aux_distance = -numeric_limits<double>::max();
00078   aux.delta_r      = 0.0;
00079   aux.z            = 1.0;
00080   _recurse_through_jet(jet, aux, jet); // last arg remains original jet
00081 
00082   // create the result and its associated structure
00083   PseudoJet result = aux.jet;
00084 
00085   // the tagger is considered to have failed if aux has never been set
00086   // (in which case it will not have parents).
00087   if (result == PseudoJet()) return result;
00088 
00089   // otherwise sort out the structure
00090   CASubJetTaggerStructure * s = new CASubJetTaggerStructure(result);
00091 //  s->_original_jet = jet;
00092   s->_scale_choice = _scale_choice;
00093   s->_distance     = aux.aux_distance;
00094   s->_absolute_z   = _absolute_z_cut;
00095   s->_z            = aux.z;
00096 
00097   result.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(s));
00098 
00099   return result;
00100 }
00101 
00102 
00103 ///----------------------------------------------------------------------
00104 /// work through the jet, establishing a distance at each branching
00105 inline void CASubJetTagger::_recurse_through_jet(const fastjet::PseudoJet & jet, JetAux &aux, const PseudoJet & original_jet) const {
00106 
00107   fastjet::PseudoJet parent1, parent2;
00108   if (! jet.has_parents(parent1, parent2)) return;
00109 
00110   /// make sure the objects are not _too_ close together
00111   if (parent1.squared_distance(parent2) < _dr2_min) return;
00112 
00113   // distance
00114   double dist=0.0;
00115   switch (_scale_choice) {
00116   case kt2_distance:
00117     // a standard (LI) kt distance
00118     dist = parent1.kt_distance(parent2);
00119     break;
00120   case jade_distance:
00121     // something a bit like a mass: pti ptj Delta R_ij^2
00122     dist = parent1.perp()*parent2.perp()*parent1.squared_distance(parent2);
00123     break;
00124   case jade2_distance:
00125     // something a bit like a mass*deltaR^2: pti ptj Delta R_ij^4
00126     dist = parent1.perp()*parent2.perp()*pow(parent1.squared_distance(parent2),2);
00127     break;
00128   case plain_distance:
00129     // Delta R_ij^2
00130     dist = parent1.squared_distance(parent2);
00131     break;
00132   case mass_drop_distance:
00133     // Delta R_ij^2
00134     dist = jet.m() - std::max(parent1.m(),parent2.m());
00135     break;
00136   case dot_product_distance:
00137     // parent1 . parent2 
00138     // ( = jet.m2() - parent1.m2() - parent2.m() in a 
00139     // 4-vector recombination scheme) 
00140     dist = dot_product(parent1, parent2);
00141     break;
00142   default:
00143     throw fastjet::Error("unrecognized scale choice");
00144   }
00145 
00146   // check the z cut
00147   bool zcut1 = true;
00148   bool zcut2 = true;
00149   double z2 = 0.0;
00150 
00151   // not very efficient -- sort out later
00152   if (parent1.perp2() < parent2.perp2()) std::swap(parent1,parent2);
00153 
00154   if (_absolute_z_cut) {
00155     z2    = parent2.perp() / original_jet.perp();
00156     zcut1 = parent1.perp() / original_jet.perp() >= _z_threshold;
00157   } else {
00158     z2    = parent2.perp()/(parent1.perp()+parent2.perp());
00159   }
00160   zcut2 = z2 >= _z_threshold;
00161 
00162   if (zcut1 && zcut2){
00163     if (dist > aux.aux_distance){
00164       aux.jet          = jet;
00165       aux.aux_distance = dist;
00166       aux.delta_r      = sqrt(parent1.squared_distance(parent2));
00167       aux.z            = z2; // the softest
00168     }
00169   }    
00170 
00171   if (zcut1) _recurse_through_jet(parent1, aux, original_jet);
00172   if (zcut2) _recurse_through_jet(parent2, aux, original_jet);
00173 }
00174 
00175 FASTJET_END_NAMESPACE
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends