FastJet 3.0beta1
|
00001 //STARTHEADER 00002 // $Id: ClusterSequenceActiveAreaExplicitGhosts.hh 2269 2011-06-21 12:37:08Z salam $ 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_CLUSTERSEQUENCEACTIVEAREAEXPLICITGHOSTS_HH_ 00032 #define __FASTJET_CLUSTERSEQUENCEACTIVEAREAEXPLICITGHOSTS_HH_ 00033 00034 #include "fastjet/PseudoJet.hh" 00035 #include "fastjet/ClusterSequenceAreaBase.hh" 00036 #include "fastjet/GhostedAreaSpec.hh" 00037 #include "fastjet/internal/LimitedWarning.hh" 00038 #include<iostream> 00039 #include<vector> 00040 #include <cstdio> 00041 00042 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh 00043 00044 //====================================================================== 00045 /// @ingroup sec_area_classes 00046 /// \class ClusterSequenceActiveAreaExplicitGhosts 00047 /// Like ClusterSequence with computation of the active jet area with the 00048 /// addition of explicit ghosts 00049 /// 00050 /// Class that behaves essentially like ClusterSequence except 00051 /// that it also provides access to the area of a jet (which 00052 /// will be a random quantity... Figure out what to do about seeds 00053 /// later...) 00054 /// 00055 /// This class should not be used directly. Rather use 00056 /// ClusterSequenceArea with the appropriate AreaDefinition 00057 class ClusterSequenceActiveAreaExplicitGhosts : 00058 public ClusterSequenceAreaBase { 00059 public: 00060 /// constructor using a GhostedAreaSpec to specify how the area is 00061 /// to be measured 00062 template<class L> ClusterSequenceActiveAreaExplicitGhosts 00063 (const std::vector<L> & pseudojets, 00064 const JetDefinition & jet_def, 00065 const GhostedAreaSpec & ghost_spec, 00066 const bool & writeout_combinations = false) 00067 : ClusterSequenceAreaBase() { 00068 std::vector<L> * ghosts = NULL; 00069 _initialise(pseudojets,jet_def,&ghost_spec,ghosts,0.0, 00070 writeout_combinations); } 00071 00072 template<class L> ClusterSequenceActiveAreaExplicitGhosts 00073 (const std::vector<L> & pseudojets, 00074 const JetDefinition & jet_def, 00075 const std::vector<L> & ghosts, 00076 double ghost_area, 00077 const bool & writeout_combinations = false) 00078 : ClusterSequenceAreaBase() { 00079 const GhostedAreaSpec * ghost_spec = NULL; 00080 _initialise(pseudojets,jet_def,ghost_spec,&ghosts,ghost_area, 00081 writeout_combinations); } 00082 00083 00084 /// does the actual work of initialisation 00085 template<class L> void _initialise 00086 (const std::vector<L> & pseudojets, 00087 const JetDefinition & jet_def, 00088 const GhostedAreaSpec * ghost_spec, 00089 const std::vector<L> * ghosts, 00090 double ghost_area, 00091 const bool & writeout_combinations); 00092 00093 //vector<PseudoJet> constituents (const PseudoJet & jet) const; 00094 00095 /// returns the number of hard particles (i.e. those supplied by the user). 00096 unsigned int n_hard_particles() const; 00097 00098 /// returns the area of a jet 00099 virtual double area (const PseudoJet & jet) const; 00100 00101 /// returns a four vector corresponding to the sum (E-scheme) of the 00102 /// ghost four-vectors composing the jet area, normalised such that 00103 /// for a small contiguous area the p_t of the extended_area jet is 00104 /// equal to area of the jet. 00105 virtual PseudoJet area_4vector (const PseudoJet & jet) const; 00106 00107 /// true if a jet is made exclusively of ghosts 00108 virtual bool is_pure_ghost(const PseudoJet & jet) const; 00109 00110 /// true if the entry in the history index corresponds to a 00111 /// ghost; if hist_ix does not correspond to an actual particle 00112 /// (i.e. hist_ix < 0), then the result is false. 00113 bool is_pure_ghost(int history_index) const; 00114 00115 /// this class does have explicit ghosts 00116 virtual bool has_explicit_ghosts() const {return true;} 00117 00118 /// return the total area, corresponding to a given Selector, that 00119 /// consists of unclustered ghosts 00120 /// 00121 /// The selector needs to apply jet by jet 00122 virtual double empty_area(const Selector & selector) const; 00123 00124 /// returns the total area under study 00125 double total_area () const; 00126 00127 /// returns the largest squared transverse momentum among 00128 /// all ghosts 00129 double max_ghost_perp2() const {return _max_ghost_perp2;} 00130 00131 /// returns true if there are any particles whose transverse momentum 00132 /// if so low that there's a risk of the ghosts having modified the 00133 /// clustering sequence 00134 bool has_dangerous_particles() const {return _has_dangerous_particles;} 00135 00136 /// get the area of the ghosts 00137 //double ghost_area() const{return _ghost_area;} 00138 00139 private: 00140 00141 int _n_ghosts; 00142 double _ghost_area; 00143 std::vector<bool> _is_pure_ghost; 00144 std::vector<double> _areas; 00145 std::vector<PseudoJet> _area_4vectors; 00146 00147 // things related to checks for dangerous particles 00148 double _max_ghost_perp2; 00149 bool _has_dangerous_particles; 00150 static LimitedWarning _warnings; 00151 00152 //static int _n_warn_dangerous_particles; 00153 //static const int _max_warn_dangerous_particles = 5; 00154 00155 00156 unsigned int _initial_hard_n; 00157 00158 /// adds the "ghost" momenta, which will be used to estimate 00159 /// the jet area 00160 void _add_ghosts(const GhostedAreaSpec & ghost_spec); 00161 00162 /// another way of adding ghosts 00163 template<class L> void _add_ghosts ( 00164 const std::vector<L> & ghosts, 00165 double ghost_area); 00166 00167 /// routine to be called after the processing is done so as to 00168 /// establish summary information on all the jets (areas, whether 00169 /// pure ghost, etc.) 00170 void _post_process(); 00171 00172 }; 00173 00174 00175 //---------------------------------------------------------------------- 00176 // initialise from some generic type... Has to be made available 00177 // here in order for the template aspect of it to work... 00178 template<class L> void ClusterSequenceActiveAreaExplicitGhosts::_initialise 00179 (const std::vector<L> & pseudojets, 00180 const JetDefinition & jet_def, 00181 const GhostedAreaSpec * ghost_spec, 00182 const std::vector<L> * ghosts, 00183 double ghost_area, 00184 const bool & writeout_combinations) { 00185 // don't reserve space yet -- will be done below 00186 00187 // insert initial jets this way so that any type L that can be 00188 // converted to a pseudojet will work fine (basically PseudoJet 00189 // and any type that has [] subscript access to the momentum 00190 // components, such as CLHEP HepLorentzVector). 00191 for (unsigned int i = 0; i < pseudojets.size(); i++) { 00192 PseudoJet mom(pseudojets[i]); 00193 //mom.set_user_index(0); // for user's particles (user index now lost...) 00194 _jets.push_back(mom); 00195 _is_pure_ghost.push_back(false); 00196 } 00197 00198 _initial_hard_n = _jets.size(); 00199 00200 if (ghost_spec != NULL) { 00201 //std::cout << "about to reserve " << (_jets.size()+ghost_spec->n_ghosts())*2 << std::endl; 00202 _jets.reserve((_jets.size()+ghost_spec->n_ghosts())); 00203 _add_ghosts(*ghost_spec); 00204 } else { 00205 _jets.reserve(_jets.size()+ghosts->size()); 00206 _add_ghosts(*ghosts, ghost_area); 00207 } 00208 00209 if (writeout_combinations) { 00210 std::cout << "# Printing particles including ghosts\n"; 00211 for (unsigned j = 0; j < _jets.size(); j++) { 00212 printf("%5u %20.13f %20.13f %20.13e\n", 00213 j,_jets[j].rap(),_jets[j].phi_02pi(),_jets[j].kt2()); 00214 } 00215 std::cout << "# Finished printing particles including ghosts\n"; 00216 } 00217 00218 // this will ensure that we can still point to jets without 00219 // difficulties arising! 00220 //std::cout << _jets.size() << " " << _jets.size()*2 << " " << _jets.max_size() << std::endl; 00221 _jets.reserve(_jets.size()*2); //GPS tmp removed 00222 00223 // run the clustering 00224 _initialise_and_run(jet_def,writeout_combinations); 00225 00226 // set up all other information 00227 _post_process(); 00228 } 00229 00230 00231 inline unsigned int ClusterSequenceActiveAreaExplicitGhosts::n_hard_particles() const {return _initial_hard_n;} 00232 00233 00234 //---------------------------------------------------------------------- 00235 /// add an explicitly specified bunch of ghosts 00236 template<class L> void ClusterSequenceActiveAreaExplicitGhosts::_add_ghosts ( 00237 const std::vector<L> & ghosts, 00238 double ghost_area) { 00239 00240 00241 for (unsigned i = 0; i < ghosts.size(); i++) { 00242 _is_pure_ghost.push_back(true); 00243 _jets.push_back(ghosts[i]); 00244 } 00245 // and record some info about ghosts 00246 _ghost_area = ghost_area; 00247 _n_ghosts = ghosts.size(); 00248 } 00249 00250 00251 FASTJET_END_NAMESPACE 00252 00253 #endif // __FASTJET_CLUSTERSEQUENCEACTIVEAREAEXPLICITGHOSTS_HH_