FastJet 3.0beta1
|
00001 //STARTHEADER 00002 // $Id: ClusterSequenceActiveArea.hh 2366 2011-07-06 09:57:52Z soyez $ 00003 // 00004 // Copyright (c) 2005-2006, Matteo Cacciari and Gavin Salam 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 #ifndef __FASTJET_CLUSTERSEQUENCEACTIVEAREA_HH__ 00032 #define __FASTJET_CLUSTERSEQUENCEACTIVEAREA_HH__ 00033 00034 00035 #include "fastjet/PseudoJet.hh" 00036 #include "fastjet/ClusterSequenceAreaBase.hh" 00037 #include "fastjet/ClusterSequenceActiveAreaExplicitGhosts.hh" 00038 #include<iostream> 00039 #include<vector> 00040 00041 //------------ backwards compatibility with version 2.1 ------------- 00042 // for backwards compatibility make ActiveAreaSpec name available 00043 //#include "fastjet/ActiveAreaSpec.hh" 00044 //#include "fastjet/ClusterSequenceWithArea.hh" 00045 //-------------------------------------------------------------------- 00046 00047 00048 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 00049 00050 //using namespace std; 00051 00052 /// @ingroup sec_area_classes 00053 /// \class ClusterSequenceActiveArea 00054 /// Like ClusterSequence with computation of the active jet area 00055 /// 00056 /// Class that behaves essentially like ClusterSequence except 00057 /// that it also provides access to the area of a jet (which 00058 /// will be a random quantity... Figure out what to do about seeds 00059 /// later...) 00060 /// 00061 /// This class should not be used directly. Rather use 00062 /// ClusterSequenceArea with the appropriate AreaDefinition 00063 class ClusterSequenceActiveArea : public ClusterSequenceAreaBase { 00064 public: 00065 00066 /// default constructor 00067 ClusterSequenceActiveArea() {} 00068 00069 /// constructor based on JetDefinition and GhostedAreaSpec 00070 template<class L> ClusterSequenceActiveArea 00071 (const std::vector<L> & pseudojets, 00072 const JetDefinition & jet_def, 00073 const GhostedAreaSpec & ghost_spec, 00074 const bool & writeout_combinations = false) ; 00075 00076 virtual double area (const PseudoJet & jet) const { 00077 return _average_area[jet.cluster_hist_index()];}; 00078 virtual double area_error (const PseudoJet & jet) const { 00079 return _average_area2[jet.cluster_hist_index()];}; 00080 00081 virtual PseudoJet area_4vector (const PseudoJet & jet) const { 00082 return _average_area_4vector[jet.cluster_hist_index()];}; 00083 00084 /// enum providing a variety of tentative strategies for estimating 00085 /// the background (e.g. non-jet) activity in a highly populated event; the 00086 /// one that has been most extensively tested is median. 00087 enum mean_pt_strategies{median=0, non_ghost_median, pttot_over_areatot, 00088 pttot_over_areatot_cut, mean_ratio_cut, play, 00089 median_4vector}; 00090 00091 /// return the transverse momentum per unit area according to one 00092 /// of the above strategies; for some strategies (those with "cut" 00093 /// in their name) the parameter "range" allows one to exclude a 00094 /// subset of the jets for the background estimation, those that 00095 /// have pt/area > median(pt/area)*range. 00096 /// 00097 /// NB: This call is OBSOLETE; use media_pt_per_unit_area from the 00098 // ClusterSequenceAreaBase class instead 00099 double pt_per_unit_area(mean_pt_strategies strat=median, 00100 double range=2.0 ) const; 00101 00102 // following code removed -- now dealt with by AreaBase class (and 00103 // this definition here conflicts with it). 00104 // /// fits a form pt_per_unit_area(y) = a + b*y^2 in the range 00105 // /// abs(y)<raprange (for negative raprange, it defaults to 00106 // /// _safe_rap_for_area). 00107 // void parabolic_pt_per_unit_area(double & a,double & b, double raprange=-1.0, 00108 // double exclude_above=-1.0, 00109 // bool use_area_4vector=false ) const; 00110 // 00111 /// rewrite the empty area from the parent class, so as to use 00112 /// all info at our disposal 00113 /// return the total area, corresponding to a given Selector, that 00114 /// consists of ghost jets or unclustered ghosts 00115 /// 00116 /// The selector passed as an argument needs to apply jet by jet. 00117 virtual double empty_area(const Selector & selector) const; 00118 00119 /// return the true number of empty jets (replaces 00120 /// ClusterSequenceAreaBase::n_empty_jets(...)) 00121 virtual double n_empty_jets(const Selector & selector) const; 00122 00123 protected: 00124 void _resize_and_zero_AA (); 00125 void _initialise_AA(const JetDefinition & jet_def, 00126 const GhostedAreaSpec & ghost_spec, 00127 const bool & writeout_combinations, 00128 bool & continue_running); 00129 00130 void _run_AA(const GhostedAreaSpec & ghost_spec); 00131 00132 void _postprocess_AA(const GhostedAreaSpec & ghost_spec); 00133 00134 /// does the initialisation and running specific to the active 00135 /// areas class 00136 void _initialise_and_run_AA (const JetDefinition & jet_def, 00137 const GhostedAreaSpec & ghost_spec, 00138 const bool & writeout_combinations = false); 00139 00140 /// transfer the history (and jet-momenta) from clust_seq to our 00141 /// own internal structure while removing ghosts 00142 void _transfer_ghost_free_history( 00143 const ClusterSequenceActiveAreaExplicitGhosts & clust_seq); 00144 00145 00146 /// transfer areas from the ClusterSequenceActiveAreaExplicitGhosts 00147 /// object into our internal area bookkeeping... 00148 void _transfer_areas(const std::vector<int> & unique_hist_order, 00149 const ClusterSequenceActiveAreaExplicitGhosts & ); 00150 00151 /// child classes benefit from having these at their disposal 00152 std::valarray<double> _average_area, _average_area2; 00153 std::valarray<PseudoJet> _average_area_4vector; 00154 00155 /// returns true if there are any particles whose transverse momentum 00156 /// if so low that there's a risk of the ghosts having modified the 00157 /// clustering sequence 00158 bool has_dangerous_particles() const {return _has_dangerous_particles;} 00159 00160 private: 00161 00162 00163 double _non_jet_area, _non_jet_area2, _non_jet_number; 00164 00165 double _maxrap_for_area; // max rap where we put ghosts 00166 double _safe_rap_for_area; // max rap where we trust jet areas 00167 00168 bool _has_dangerous_particles; 00169 00170 00171 /// routine for extracting the tree in an order that will be independent 00172 /// of any degeneracies in the recombination sequence that don't 00173 /// affect the composition of the final jets 00174 void _extract_tree(std::vector<int> &) const; 00175 /// do the part of the extraction associated with pos, working 00176 /// through its children and their parents 00177 void _extract_tree_children(int pos, std::valarray<bool> &, const std::valarray<int> &, std::vector<int> &) const; 00178 /// do the part of the extraction associated with the parents of pos. 00179 void _extract_tree_parents (int pos, std::valarray<bool> &, const std::valarray<int> &, std::vector<int> &) const; 00180 00181 /// check if two jets have the same momentum to within the 00182 /// tolerance (and if pt's are not the same we're forgiving and 00183 /// look to see if the energy is the same) 00184 void _throw_unless_jets_have_same_perp_or_E(const PseudoJet & jet, 00185 const PseudoJet & refjet, 00186 double tolerance, 00187 const ClusterSequenceActiveAreaExplicitGhosts & jets_ghosted_seq 00188 ) const; 00189 00190 /// since we are playing nasty games with seeds, we should warn 00191 /// the user a few times 00192 //static int _n_seed_warnings; 00193 //const static int _max_seed_warnings = 10; 00194 00195 // record the number of repeats 00196 int _ghost_spec_repeat; 00197 00198 /// a class for our internal storage of ghost jets 00199 class GhostJet : public PseudoJet { 00200 public: 00201 GhostJet(const PseudoJet & j, double a) : PseudoJet(j), area(a){} 00202 double area; 00203 }; 00204 00205 std::vector<GhostJet> _ghost_jets; 00206 std::vector<GhostJet> _unclustered_ghosts; 00207 }; 00208 00209 00210 00211 00212 template<class L> ClusterSequenceActiveArea::ClusterSequenceActiveArea 00213 (const std::vector<L> & pseudojets, 00214 const JetDefinition & jet_def, 00215 const GhostedAreaSpec & ghost_spec, 00216 const bool & writeout_combinations) { 00217 00218 // transfer the initial jets (type L) into our own array 00219 _transfer_input_jets(pseudojets); 00220 00221 // run the clustering for active areas 00222 _initialise_and_run_AA(jet_def, ghost_spec, writeout_combinations); 00223 00224 } 00225 00226 00227 00228 FASTJET_END_NAMESPACE 00229 00230 #endif // __FASTJET_CLUSTERSEQUENCEACTIVEAREA_HH__