#include <SISConePlugin.hh>
Inheritance diagram for fastjet::SISConePlugin:
Public Types | |
enum | SplitMergeScale { SM_pt, SM_Et, SM_mt, SM_pttilde } |
enum for the different split-merge scale choices; Note that order _must_ be the same as in siscone More... | |
Public Member Functions | |
SISConePlugin (double cone_radius, double overlap_threshold=0.5, int n_pass_max=0, double protojet_ptmin=0.0, bool caching=false, SplitMergeScale split_merge_scale=SM_pttilde, double split_merge_stopping_scale=0.0) | |
Constructor for the SISCone Plugin class. | |
SISConePlugin (double cone_radius, double overlap_threshold, int n_pass_max, double protojet_ptmin, bool caching, bool split_merge_on_transverse_mass) | |
Backwards compatible constructor for the SISCone Plugin class. | |
SISConePlugin (double cone_radius, double overlap_threshold, int n_pass_max, bool caching) | |
backwards compatible constructor for the SISCone Plugin class (avoid using this in future). | |
SISConePlugin (const SISConePlugin &plugin) | |
copy constructor | |
double | cone_radius () const |
the cone radius | |
double | overlap_threshold () const |
Fraction of overlap energy in a jet above which jets are merged and below which jets are split. | |
int | n_pass_max () const |
the maximum number of passes of stable-cone searching (<=0 is same as infinity). | |
double | protojet_ptmin () const |
minimum pt for a protojet to be considered in the split-merge step of the algorithm | |
double | protojet_or_ghost_ptmin () const |
return the scale to be passed to SISCone as the protojet_ptmin -- if we have a ghost separation scale that is above the protojet_ptmin, then the ghost_separation_scale becomes the relevant one to use here | |
SplitMergeScale | split_merge_scale () const |
indicates scale used in split-merge | |
void | set_split_merge_scale (SplitMergeScale sms) |
sets scale used in split-merge | |
bool | split_merge_on_transverse_mass () const |
indicates whether the split-merge orders on transverse mass or not. | |
void | set_split_merge_on_transverse_mass (bool val) |
void | set_split_merge_stopping_scale (double scale) |
set the "split_merge_stopping_scale": if the scale variable for all protojets is below this, then stop the split-merge procedure and keep only those jets found so far. | |
double | split_merge_stopping_scale () |
return the value of the split_merge_stopping_scale (see set_split_merge_stopping_scale(. | |
bool | caching () const |
indicates whether caching is turned on or not. | |
virtual std::string | description () const |
return a textual description of the jet-definition implemented in this plugin | |
virtual void | run_clustering (ClusterSequence &) const |
given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:
| |
virtual double | R () const |
the plugin mechanism's standard way of accessing the jet radius | |
virtual bool | supports_ghosted_passive_areas () const |
return true since there is specific support for the measurement of passive areas, in the sense that areas determined from all particles below the ghost separation scale will be a passive area. | |
virtual void | set_ghost_separation_scale (double scale) const |
set the ghost separation scale for passive area determinations _just_ in the next run (strictly speaking that makes the routine a non const, so related internal info must be stored as a mutable) | |
virtual double | ghost_separation_scale () const |
Private Attributes | |
double | _cone_radius |
double | _overlap_threshold |
int | _n_pass_max |
double | _protojet_ptmin |
bool | _caching |
SplitMergeScale | _split_merge_scale |
double | _split_merge_stopping_scale |
double | _ghost_sep_scale |
Static Private Attributes | |
static std::auto_ptr< SISConePlugin > | stored_plugin |
static std::auto_ptr< std::vector< PseudoJet > > | stored_particles |
static std::auto_ptr< siscone::Csiscone > | stored_siscone |
As of 2006-12-26, this plugin is beta, as is the SISCone code itself.
SISCone uses geometrical techniques to exhaustively consider all possible distinct cones. It then finds out which ones are stable and sends the result to the Tevatron Run-II type split-merge procedure for overlapping cones.
Four parameters govern the "physics" of the algorithm:
One parameter governs some internal algorithmic shortcuts:
The final jets can be accessed by requestion the inclusive_jets(...) from the ClusterSequence object. Note that these PseudoJets have their user_index() set to the index of the pass in which they were found (first pass = 0). NB: This does not currently work for jets that consist of a single particle.
For further information on the details of the algorithm see the SISCone paper; for documentation about the implementation, see the siscone/doc/html/index.html file.
Definition at line 76 of file SISConePlugin.hh.
|
enum for the different split-merge scale choices; Note that order _must_ be the same as in siscone
Definition at line 81 of file SISConePlugin.hh. 00081 {SM_pt, 00082 SM_Et, 00083 00084 SM_mt, 00085 00086 SM_pttilde 00087 00088 };
|
|
Constructor for the SISCone Plugin class. Note: though the default value here for the overlap_threshold is 0.5 (for backwards compatibility), there is a strong recommendation to use a higher value, e.g. 0.75, especially in environments with a substantial amount of underlying event or pileup. Definition at line 98 of file SISConePlugin.hh. Referenced by run_clustering(). 00104 : 00105 _cone_radius (cone_radius ), 00106 _overlap_threshold (overlap_threshold ), 00107 _n_pass_max (n_pass_max ), 00108 _protojet_ptmin (protojet_ptmin), 00109 _caching (caching), 00110 _split_merge_scale (split_merge_scale), 00111 _split_merge_stopping_scale (split_merge_stopping_scale), 00112 _ghost_sep_scale (0.0) {}
|
|
Backwards compatible constructor for the SISCone Plugin class.
Definition at line 115 of file SISConePlugin.hh. 00120 : 00121 _cone_radius (cone_radius ), 00122 _overlap_threshold (overlap_threshold ), 00123 _n_pass_max (n_pass_max ), 00124 _protojet_ptmin (protojet_ptmin), 00125 _caching (caching), 00126 _split_merge_scale (split_merge_on_transverse_mass ? SM_mt : SM_pttilde), 00127 _ghost_sep_scale (0.0) {}
|
|
backwards compatible constructor for the SISCone Plugin class (avoid using this in future).
Definition at line 131 of file SISConePlugin.hh. 00134 : 00135 _cone_radius (cone_radius ), 00136 _overlap_threshold (overlap_threshold ), 00137 _n_pass_max (n_pass_max ), 00138 _protojet_ptmin (0.0), 00139 _caching (caching), 00140 _split_merge_scale (SM_mt), 00141 _ghost_sep_scale (0.0) {}
|
|
copy constructor
Definition at line 144 of file SISConePlugin.hh. 00144 {
00145 *this = plugin;
00146 }
|
|
indicates whether caching is turned on or not.
Definition at line 195 of file SISConePlugin.hh. Referenced by description(), and run_clustering(). 00195 {return _caching ;}
|
|
the cone radius
Definition at line 149 of file SISConePlugin.hh. Referenced by description(), and run_clustering(). 00149 {return _cone_radius ;}
|
|
return a textual description of the jet-definition implemented in this plugin
Implements fastjet::JetDefinition::Plugin. Definition at line 23 of file SISConePlugin.cc. References _split_merge_stopping_scale, caching(), cone_radius(), n_pass_max(), overlap_threshold(), protojet_ptmin(), and split_merge_scale(). 00023 { 00024 ostringstream desc; 00025 00026 const string on = "on"; 00027 const string off = "off"; 00028 00029 string sm_scale_string = "split-merge uses " + 00030 split_merge_scale_name(Esplit_merge_scale(split_merge_scale())); 00031 00032 desc << "SISCone jet algorithm with " ; 00033 desc << "cone_radius = " << cone_radius () << ", "; 00034 desc << "overlap_threshold = " << overlap_threshold () << ", "; 00035 desc << "n_pass_max = " << n_pass_max () << ", "; 00036 desc << "protojet_ptmin = " << protojet_ptmin() << ", "; 00037 desc << sm_scale_string << ", "; 00038 desc << "caching turned " << (caching() ? on : off); 00039 desc << ", SM stop scale = " << _split_merge_stopping_scale; 00040 00041 00042 // create a fake siscone object so that we can find out more about it 00043 Csiscone siscone; 00044 if (siscone.merge_identical_protocones) { 00045 desc << ", and (IR unsafe) merge_indentical_protocones=true" ; 00046 } 00047 00048 desc << ", SISCone code v" << siscone_version(); 00049 00050 return desc.str(); 00051 }
|
|
Reimplemented from fastjet::JetDefinition::Plugin. Definition at line 216 of file SISConePlugin.hh. Referenced by run_clustering(). 00216 {return _ghost_sep_scale;}
|
|
the maximum number of passes of stable-cone searching (<=0 is same as infinity).
Definition at line 157 of file SISConePlugin.hh. Referenced by description(), and run_clustering(). 00157 {return _n_pass_max ;}
|
|
Fraction of overlap energy in a jet above which jets are merged and below which jets are split.
Definition at line 153 of file SISConePlugin.hh. Referenced by description(), and run_clustering(). 00153 {return _overlap_threshold ;}
|
|
return the scale to be passed to SISCone as the protojet_ptmin -- if we have a ghost separation scale that is above the protojet_ptmin, then the ghost_separation_scale becomes the relevant one to use here
Definition at line 167 of file SISConePlugin.hh. Referenced by run_clustering(). 00167 {return std::max(_protojet_ptmin, 00168 _ghost_sep_scale);}
|
|
minimum pt for a protojet to be considered in the split-merge step of the algorithm
Definition at line 161 of file SISConePlugin.hh. Referenced by description(). 00161 {return _protojet_ptmin ;}
|
|
the plugin mechanism's standard way of accessing the jet radius
Implements fastjet::JetDefinition::Plugin. Definition at line 201 of file SISConePlugin.hh. 00201 {return cone_radius();}
|
|
given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:
..)
Implements fastjet::JetDefinition::Plugin. Definition at line 61 of file SISConePlugin.cc. References fastjet::SISConeExtras::_jet_def_plugin, fastjet::SISConeExtras::_most_ambiguous_split, fastjet::SISConeExtras::_protocones, _split_merge_stopping_scale, caching(), cone_radius(), fastjet::PseudoJet::E(), ghost_separation_scale(), fastjet::have_same_momentum(), fastjet::ClusterSequence::jets(), n_pass_max(), overlap_threshold(), fastjet::ClusterSequence::plugin_associate_extras(), fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), protojet_or_ghost_ptmin(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), fastjet::PseudoJet::pz(), fastjet::PseudoJet::set_user_index(), SISConePlugin(), split_merge_scale(), stored_particles, stored_plugin, and stored_siscone. 00061 { 00062 00063 00064 Csiscone * siscone; 00065 Csiscone local_siscone; 00066 00067 unsigned n = clust_seq.jets().size(); 00068 00069 bool new_siscone = true; // by default we'll be running it 00070 00071 if (caching()) { 00072 00073 // Establish if we have a cached run with the same R, npass and 00074 // particles. If not then do any tidying up / reallocation that's 00075 // necessary for the next round of caching, otherwise just set 00076 // relevant pointers so that we can reuse and old run. 00077 if (stored_siscone.get() != 0) { 00078 new_siscone = !(stored_plugin->cone_radius() == cone_radius() 00079 && stored_plugin->n_pass_max() == n_pass_max() 00080 && stored_particles->size() == n); 00081 if (!new_siscone) { 00082 for(unsigned i = 0; i < n; i++) { 00083 // only check momentum because indices will be correctly dealt 00084 // with anyway when extracting the clustering order. 00085 new_siscone |= !have_same_momentum(clust_seq.jets()[i], 00086 (*stored_particles)[i]); 00087 } 00088 } 00089 } 00090 00091 // allocate the new siscone, etc., if need be 00092 if (new_siscone) { 00093 stored_siscone .reset( new Csiscone ); 00094 stored_particles.reset( new vector<PseudoJet>(clust_seq.jets())); 00095 stored_plugin .reset( new SISConePlugin(*this) ); 00096 } 00097 00098 siscone = stored_siscone.get(); 00099 } else { 00100 siscone = &local_siscone; 00101 } 00102 00103 // make sure stopping scale is set in siscone 00104 siscone->SM_var2_hardest_cut_off = _split_merge_stopping_scale*_split_merge_stopping_scale; 00105 // when running with ghosts for passive areas, do not put the 00106 // ghosts into the stable-cone search (not relevant) 00107 siscone->stable_cone_soft_pt2_cutoff = ghost_separation_scale() 00108 * ghost_separation_scale(); 00109 00110 if (new_siscone) { 00111 // transfer fastjet initial particles into the siscone type 00112 vector<Cmomentum> siscone_momenta(n); 00113 for(unsigned i = 0; i < n; i++) { 00114 const PseudoJet & p = clust_seq.jets()[i]; // shorthand 00115 siscone_momenta[i] = Cmomentum(p.px(), p.py(), p.pz(), p.E()); 00116 } 00117 00118 // run the jet finding 00119 //cout << "plg sms: " << split_merge_scale() << endl; 00120 siscone->compute_jets(siscone_momenta, cone_radius(), overlap_threshold(), 00121 n_pass_max(), protojet_or_ghost_ptmin(), 00122 Esplit_merge_scale(split_merge_scale())); 00123 } else { 00124 // just run the overlap part of the jets. 00125 //cout << "plg rcmp sms: " << split_merge_scale() << endl; 00126 siscone->recompute_jets(overlap_threshold(), protojet_or_ghost_ptmin(), 00127 Esplit_merge_scale(split_merge_scale())); 00128 } 00129 00130 // extract the jets [in reverse order -- to get nice ordering in pt at end] 00131 int njet = siscone->jets.size(); 00132 00133 for (int ijet = njet-1; ijet >= 0; ijet--) { 00134 const Cjet & jet = siscone->jets[ijet]; // shorthand 00135 00136 // Successively merge the particles that make up the cone jet 00137 // until we have all particles in it. Start off with the zeroth 00138 // particle. 00139 int jet_k = jet.contents[0]; 00140 for (unsigned ipart = 1; ipart < jet.contents.size(); ipart++) { 00141 // take the last result of the merge 00142 int jet_i = jet_k; 00143 // and the next element of the jet 00144 int jet_j = jet.contents[ipart]; 00145 // and merge them (with a fake dij) 00146 double dij = 0.0; 00147 00148 // create the new jet by hand so that we can adjust its user index 00149 PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j]; 00150 00151 // set the user index to be the pass in which the jet was discovered 00152 newjet.set_user_index(jet.pass); 00153 00154 clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k); 00155 } 00156 // we have merged all the jet's particles into a single object, so now 00157 // "declare" it to be a beam (inclusive) jet. 00158 // [NB: put a sensible looking d_iB just to be nice...] 00159 double d_iB = clust_seq.jets()[jet_k].perp2(); 00160 clust_seq.plugin_record_iB_recombination(jet_k, d_iB); 00161 } 00162 00163 // now copy the list of protocones into an "extras" objects 00164 SISConeExtras * extras = new SISConeExtras; 00165 for (unsigned ipass = 0; ipass < siscone->protocones_list.size(); ipass++) { 00166 for (unsigned ipc = 0; ipc < siscone->protocones_list[ipass].size(); ipc++) { 00167 double rap = siscone->protocones_list[ipass][ipc].eta; 00168 double phi = siscone->protocones_list[ipass][ipc].phi; 00169 PseudoJet protocone(cos(phi),sin(phi),sinh(rap),cosh(rap)); 00170 //PseudoJet protocone(siscone->protocones_list[ipass][ipc]); 00171 protocone.set_user_index(ipass); 00172 extras->_protocones.push_back(protocone); 00173 } 00174 } 00175 extras->_most_ambiguous_split = siscone->most_ambiguous_split; 00176 00177 // tell it what the jet definition was 00178 extras->_jet_def_plugin = this; 00179 00180 // give the extras object to the cluster sequence. 00181 clust_seq.plugin_associate_extras(auto_ptr<ClusterSequence::Extras>(extras)); 00182 }
|
|
set the ghost separation scale for passive area determinations _just_ in the next run (strictly speaking that makes the routine a non const, so related internal info must be stored as a mutable)
Reimplemented from fastjet::JetDefinition::Plugin. Definition at line 212 of file SISConePlugin.hh. 00212 { 00213 _ghost_sep_scale = scale; 00214 }
|
|
Definition at line 178 of file SISConePlugin.hh. 00178 { 00179 _split_merge_scale = val ? SM_mt : SM_pt;}
|
|
sets scale used in split-merge
Definition at line 173 of file SISConePlugin.hh. 00173 {_split_merge_scale = sms;}
|
|
set the "split_merge_stopping_scale": if the scale variable for all protojets is below this, then stop the split-merge procedure and keep only those jets found so far. This is useful in determination of areas of hard jets because it can be used to avoid running the split-merging on the pure ghost-part of the event. Definition at line 187 of file SISConePlugin.hh. 00187 { 00188 _split_merge_stopping_scale = scale;}
|
|
indicates whether the split-merge orders on transverse mass or not. retained for backwards compatibility with 2.1.0b3 Definition at line 177 of file SISConePlugin.hh. 00177 {return _split_merge_scale == SM_mt ;}
|
|
indicates scale used in split-merge
Definition at line 171 of file SISConePlugin.hh. Referenced by description(), and run_clustering(). 00171 {return _split_merge_scale;}
|
|
return the value of the split_merge_stopping_scale (see set_split_merge_stopping_scale(. ..) for description) Definition at line 192 of file SISConePlugin.hh. 00192 {return _split_merge_stopping_scale;}
|
|
return true since there is specific support for the measurement of passive areas, in the sense that areas determined from all particles below the ghost separation scale will be a passive area.
Reimplemented from fastjet::JetDefinition::Plugin. Definition at line 207 of file SISConePlugin.hh. 00207 {return true;}
|
|
Definition at line 222 of file SISConePlugin.hh. |
|
Definition at line 219 of file SISConePlugin.hh. |
|
Definition at line 226 of file SISConePlugin.hh. |
|
Definition at line 220 of file SISConePlugin.hh. |
|
Definition at line 219 of file SISConePlugin.hh. |
|
Definition at line 221 of file SISConePlugin.hh. |
|
Definition at line 223 of file SISConePlugin.hh. |
|
Definition at line 224 of file SISConePlugin.hh. Referenced by description(), and run_clustering(). |
|
Definition at line 20 of file SISConePlugin.cc. Referenced by run_clustering(). |
|
Definition at line 19 of file SISConePlugin.cc. Referenced by run_clustering(). |
|
Definition at line 21 of file SISConePlugin.cc. Referenced by run_clustering(). |