FastJet 3.0.2
CASubJetTagger.hh
00001 //STARTHEADER
00002 // $Id: CASubJetTagger.hh 2616 2011-09-30 18:03:40Z salam $
00003 //
00004 // Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. 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, see <http://www.gnu.org/licenses/>.
00026 //----------------------------------------------------------------------
00027 //ENDHEADER
00028 
00029 #ifndef __CASUBJET_TAGGER_HH__
00030 #define __CASUBJET_TAGGER_HH__
00031 
00032 #include <fastjet/PseudoJet.hh>
00033 #include <fastjet/WrappedStructure.hh>
00034 #include <fastjet/tools/Transformer.hh>
00035 #include "fastjet/LimitedWarning.hh"
00036 
00037 FASTJET_BEGIN_NAMESPACE
00038 
00039 class CASubJetTagger;
00040 class CASubJetTaggerStructure;
00041 
00042 //----------------------------------------------------------------------
00043 /// @ingroup tools_taggers
00044 /// \class CASubJetTagger
00045 /// clean (almost parameter-free) tagger searching for the element in
00046 /// the clustering history that maximises a chosen distance
00047 ///
00048 /// class to help us get a clean (almost parameter-free) handle on
00049 /// substructure inside a C/A jet. It follows the logic described in
00050 /// arXiv:0906.0728 (and is inspired by the original Cambridge
00051 /// algorithm paper in its use of separate angular and dimensionful
00052 /// distances), but provides some extra flexibility.
00053 ///
00054 /// It searches for all splittings that pass a symmetry cut (zcut) and
00055 /// then selects the one with the largest auxiliary scale choice
00056 /// (e.g. jade distance of the splitting, kt distance of the
00057 /// splitting, etc.)
00058 ///
00059 /// By default, the zcut is calculated from the fraction of the child
00060 /// pt carried by the parent jet. If one calls set_absolute_z_cut the
00061 /// fraction of transverse momentum will be computed wrt the original
00062 /// jet.
00063 ///
00064 /// original code copyright (C) 2009 by Gavin Salam, released under
00065 /// the GPL.
00066 ///
00067 /// \section desc Options
00068 /// 
00069 ///  - the distance choice: options are
00070 ///        kt2_distance        : usual min(kti^2,ktj^2)DeltaR_{ij}^2
00071 ///        jade_distance       : kti . ktj DeltaR_{ij}^2 (LI version of jade)
00072 ///        jade2_distance      : kti . ktj DeltaR_{ij}^4 (LI version of jade * DR^2)
00073 ///        plain_distance      :           DeltaR_{ij}^2
00074 ///        mass_drop_distance  : m_jet - max(m_parent1,m_parent2)
00075 ///        dot_product_distance: parent1.parent2
00076 ///    (kt2_distance by default)
00077 ///
00078 ///  - the z cut (0 by default)
00079 ///
00080 ///  - by calling set_absolute_z_cut(), one can ask that the pt
00081 ///    fraction if calculated wrt the original jet
00082 ///
00083 ///  - by calling set_dr_min(drmin), one can ask that only the
00084 ///    recombinations where the 2 objects are (geometrically) distant
00085 ///    by at least drmin are kept in the maximisation.
00086 ///
00087 /// \section input Input conditions
00088 /// 
00089 ///  - the jet must have been obtained from a Cambridge/Aachen cluster
00090 ///    sequence
00091 ///
00092 /// \section output Output/structure
00093 /// 
00094 ///  - the element of the cluster sequence maximising the requested
00095 ///    distance (and satisfying the zcut) is returned.
00096 ///
00097 ///  - if the original jet has no parents, it will be returned
00098 ///
00099 ///  - the value of the "z" and distance corresponding to that history
00100 ///    element are stored and accessible through
00101 ///      result.structure_of<CASubJetTagger>().z();
00102 ///      result.structure_of<CASubJetTagger>().distance();
00103 ///
00104 class CASubJetTagger : public Transformer {
00105 public:
00106 
00107   /// the available choices of auxiliary scale with respect to which
00108   /// to order the splittings
00109   enum ScaleChoice {
00110     kt2_distance,        //< usual min(kti^2,ktj^2)DeltaR_{ij}^2
00111     jade_distance,       //< kti . ktj DeltaR_{ij}^2 (LI version of jade)
00112     jade2_distance,      //< kti . ktj DeltaR_{ij}^4 (LI version of jade * DR^2)
00113     plain_distance,      //<           DeltaR_{ij}^2
00114     mass_drop_distance,  //<    m_jet - max(m_parent1,m_parent2)
00115     dot_product_distance //<  parent1.parent2
00116   };
00117 
00118   /// just constructs
00119   CASubJetTagger(ScaleChoice scale_choice = jade_distance,
00120                  double      z_threshold  = 0.1)
00121     : _scale_choice(scale_choice), _z_threshold(z_threshold),
00122       _dr2_min(0.0), _absolute_z_cut(false){};
00123 
00124   /// sets a minimum delta R below which spliting will be ignored
00125   /// (only relevant if set prior to calling run())
00126   void set_dr_min(double drmin) {_dr2_min = drmin*drmin;}
00127 
00128   /// returns a textual description of the tagger
00129   virtual std::string description() const;
00130 
00131   /// If (abs_z_cut) is set to false (the default) then for a
00132   /// splitting to be considered, each subjet must satisfy 
00133   /// 
00134   ///        p_{t,sub} > z_threshold * p_{t,parent}
00135   ///
00136   /// whereas if it is set to true, then each subject must satisfy
00137   ///
00138   ///        p_{t,sub} > z_threshold * p_{t,original-jet}
00139   ///
00140   /// where parent is the immediate parent of the splitting, and
00141   /// original jet is the one supplied to the run() function.
00142   ///
00143   /// Only relevant is called prior to run().
00144   void set_absolute_z_cut(bool abs_z_cut=true) {_absolute_z_cut = abs_z_cut;}
00145   
00146   /// runs the tagger on the given jet and
00147   /// returns the tagged PseudoJet if successful, or a PseudoJet==0 otherwise
00148   /// (standard access is through operator()).
00149   virtual PseudoJet result(const PseudoJet & jet) const;
00150 
00151   /// the type of Structure returned
00152   typedef CASubJetTaggerStructure StructureType;
00153 
00154 protected:
00155   /// class that contains the result internally
00156   class JetAux {
00157   public:
00158     PseudoJet jet;          //< the subjet (immediate parent of splitting)
00159     double    aux_distance; //< the auxiliary distance between its two subjets
00160     double    delta_r;      //< the angular distance between its two subjets
00161     double    z;            //< the transverse momentum fraction
00162   };
00163 
00164 
00165   void _recurse_through_jet(const PseudoJet & current_jet, 
00166                             JetAux &aux_max,
00167                             const PseudoJet & original_jet) const;
00168 
00169   ScaleChoice _scale_choice;
00170   double      _z_threshold;
00171   double      _dr2_min;
00172   bool        _absolute_z_cut;
00173 
00174   static  LimitedWarning _non_ca_warnings;
00175 };
00176 
00177 
00178 //------------------------------------------------------------------------
00179 /// @ingroup tools_taggers
00180 /// the structure returned by a CASubJetTagger
00181 ///
00182 /// Since this is directly an element of the ClusterSequence, we keep
00183 /// basically the original ClusterSequenceStructure (wrapped for
00184 /// memory-management reasons) and add information about the pt
00185 /// fraction and distance of the subjet structure
00186 class CASubJetTaggerStructure : public WrappedStructure{
00187 public:
00188   /// default ctor
00189   ///  \param result_jet   the jet for which we have to keep the
00190   ///                      structure
00191   CASubJetTaggerStructure(const PseudoJet & result_jet)
00192     : WrappedStructure(result_jet.structure_shared_ptr()){}
00193 
00194   /// returns the scale choice asked for the maximisation
00195   CASubJetTagger::ScaleChoice scale_choice() const {return _scale_choice;}
00196 
00197   /// returns the value of the distance measure (corresponding to
00198   /// ScaleChoice) for this jet's splitting
00199   double distance() const {return _distance;}
00200 
00201   /// returns the pt fraction contained by the softer of the two component
00202   /// pieces of this jet (normalised relative to this jet)
00203   double z() const {return _z;}
00204 
00205   /// returns the pt fraction contained by the softer of the two component
00206   /// pieces of this jet (normalised relative to the original jet)
00207   bool absolute_z() const {return _absolute_z;}
00208 
00209 //  /// returns the original jet (before tagging)
00210 //  const PseudoJet & original() const {return _original_jet;}
00211 
00212 protected:
00213   CASubJetTagger::ScaleChoice _scale_choice; ///< the user scale choice 
00214   double _distance;  ///< the maximal distance associated with the result
00215   bool _absolute_z;  ///< whether z is computed wrt to the original jet or not
00216   double _z;         ///< the transverse momentum fraction
00217 //  PseudoJet _original_jet;  ///< the original jet (before tagging)
00218 
00219   friend class CASubJetTagger; ///< to allow setting the internal information
00220 };
00221 
00222 FASTJET_END_NAMESPACE
00223 
00224 #endif // __CASUBJET_HH__
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends