FastJet
3.3.0
|
deals with clustering More...
#include <fastjet/ClusterSequence.hh>
Classes | |
class | Extras |
base class to store extra information that plugins may provide More... | |
struct | history_element |
a single element in the clustering history More... | |
Public Types | |
enum | JetType { Invalid =-3, InexistentParent = -2, BeamJet = -1 } |
typedef ClusterSequenceStructure | StructureType |
the structure type associated with a jet belonging to a ClusterSequence | |
Public Member Functions | |
ClusterSequence () | |
default constructor | |
template<class L > | |
ClusterSequence (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const bool &writeout_combinations=false) | |
create a ClusterSequence, starting from the supplied set of PseudoJets and clustering them with jet definition specified by jet_def (which also specifies the clustering strategy) More... | |
ClusterSequence (const ClusterSequence &cs) | |
copy constructor for a ClusterSequence | |
ClusterSequence & | operator= (const ClusterSequence &cs) |
explicit assignment operator for a ClusterSequence | |
std::vector< PseudoJet > | inclusive_jets (const double ptmin=0.0) const |
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin. More... | |
int | n_exclusive_jets (const double dcut) const |
return the number of jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut. More... | |
std::vector< PseudoJet > | exclusive_jets (const double dcut) const |
return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut. More... | |
std::vector< PseudoJet > | exclusive_jets (const int njets) const |
return a vector of all jets when the event is clustered (in the exclusive sense) to exactly njets. More... | |
std::vector< PseudoJet > | exclusive_jets_up_to (const int njets) const |
return a vector of all jets when the event is clustered (in the exclusive sense) to exactly njets. More... | |
double | exclusive_dmerge (const int njets) const |
return the dmin corresponding to the recombination that went from n+1 to n jets (sometimes known as d_{n n+1}). More... | |
double | exclusive_dmerge_max (const int njets) const |
return the maximum of the dmin encountered during all recombinations up to the one that led to an n-jet final state; identical to exclusive_dmerge, except in cases where the dmin do not increase monotonically. More... | |
double | exclusive_ymerge (int njets) const |
return the ymin corresponding to the recombination that went from n+1 to n jets (sometimes known as y_{n n+1}). More... | |
double | exclusive_ymerge_max (int njets) const |
same as exclusive_dmerge_max, but normalised to squared total energy | |
int | n_exclusive_jets_ycut (double ycut) const |
the number of exclusive jets at the given ycut | |
std::vector< PseudoJet > | exclusive_jets_ycut (double ycut) const |
the exclusive jets obtained at the given ycut | |
std::vector< PseudoJet > | exclusive_subjets (const PseudoJet &jet, const double dcut) const |
return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut. More... | |
int | n_exclusive_subjets (const PseudoJet &jet, const double dcut) const |
return the size of exclusive_subjets(...); still n ln n with same coefficient, but marginally more efficient than manually taking exclusive_subjets.size() | |
std::vector< PseudoJet > | exclusive_subjets (const PseudoJet &jet, int nsub) const |
return the list of subjets obtained by unclustering the supplied jet down to nsub subjets. More... | |
std::vector< PseudoJet > | exclusive_subjets_up_to (const PseudoJet &jet, int nsub) const |
return the list of subjets obtained by unclustering the supplied jet down to nsub subjets (or all constituents if there are fewer than nsub). More... | |
double | exclusive_subdmerge (const PseudoJet &jet, int nsub) const |
returns the dij that was present in the merging nsub+1 -> nsub subjets inside this jet. More... | |
double | exclusive_subdmerge_max (const PseudoJet &jet, int nsub) const |
returns the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet. More... | |
double | Q () const |
returns the sum of all energies in the event (relevant mainly for e+e-) | |
double | Q2 () const |
return Q()^2 | |
bool | object_in_jet (const PseudoJet &object, const PseudoJet &jet) const |
returns true iff the object is included in the jet. More... | |
bool | has_parents (const PseudoJet &jet, PseudoJet &parent1, PseudoJet &parent2) const |
if the jet has parents in the clustering, it returns true and sets parent1 and parent2 equal to them. More... | |
bool | has_child (const PseudoJet &jet, PseudoJet &child) const |
if the jet has a child then return true and give the child jet otherwise return false and set the child to zero | |
bool | has_child (const PseudoJet &jet, const PseudoJet *&childp) const |
Version of has_child that sets a pointer to the child if the child exists;. | |
bool | has_partner (const PseudoJet &jet, PseudoJet &partner) const |
if this jet has a child (and so a partner) return true and give the partner, otherwise return false and set the partner to zero | |
std::vector< PseudoJet > | constituents (const PseudoJet &jet) const |
return a vector of the particles that make up jet | |
void | print_jets_for_root (const std::vector< PseudoJet > &jets, std::ostream &ostr=std::cout) const |
output the supplied vector of jets in a format that can be read by an appropriate root script; the format is: jet-n jet-px jet-py jet-pz jet-E particle-n particle-rap particle-phi particle-pt particle-n particle-rap particle-phi particle-pt ... More... | |
void | print_jets_for_root (const std::vector< PseudoJet > &jets, const std::string &filename, const std::string &comment="") const |
print jets for root to the file labelled filename, with an optional comment at the beginning | |
void | add_constituents (const PseudoJet &jet, std::vector< PseudoJet > &subjet_vector) const |
add on to subjet_vector the constituents of jet (for internal use mainly) | |
Strategy | strategy_used () const |
return the enum value of the strategy used to cluster the event | |
std::string | strategy_string () const |
return the name of the strategy used to cluster the event | |
std::string | strategy_string (Strategy strategy_in) const |
return the name of the strategy associated with the enum strategy_in | |
const JetDefinition & | jet_def () const |
return a reference to the jet definition | |
void | delete_self_when_unused () |
by calling this routine you tell the ClusterSequence to delete itself when all the Pseudojets associated with it have gone out of scope. More... | |
bool | will_delete_self_when_unused () const |
return true if the object has been told to delete itself when unused | |
void | signal_imminent_self_deletion () const |
tell the ClusterSequence it's about to be self deleted (internal use only) | |
double | jet_scale_for_algorithm (const PseudoJet &jet) const |
returns the scale associated with a jet as required for this clustering algorithm (kt^2 for the kt-algorithm, 1 for the Cambridge algorithm). More... | |
void | plugin_record_ij_recombination (int jet_i, int jet_j, double dij, int &newjet_k) |
record the fact that there has been a recombination between jets()[jet_i] and jets()[jet_k], with the specified dij, and return the index (newjet_k) allocated to the new jet, whose momentum is assumed to be the 4-vector sum of that of jet_i and jet_j | |
void | plugin_record_ij_recombination (int jet_i, int jet_j, double dij, const PseudoJet &newjet, int &newjet_k) |
as for the simpler variant of plugin_record_ij_recombination, except that the new jet is attributed the momentum and user_index of newjet | |
void | plugin_record_iB_recombination (int jet_i, double diB) |
record the fact that there has been a recombination between jets()[jet_i] and the beam, with the specified diB; when looking for inclusive jets, any iB recombination will returned to the user as a jet. More... | |
void | plugin_associate_extras (Extras *extras_in) |
the plugin can associate some extra information with the ClusterSequence object by calling this function. More... | |
void | plugin_associate_extras (std::auto_ptr< Extras > extras_in) |
the plugin can associate some extra information with the ClusterSequence object by calling this function More... | |
bool | plugin_activated () const |
returns true when the plugin is allowed to run the show. | |
const Extras * | extras () const |
returns a pointer to the extras object (may be null) | |
template<class GBJ > | |
void | plugin_simple_N2_cluster () |
allows a plugin to run a templated clustering (nearest-neighbour heuristic) More... | |
const std::vector< PseudoJet > & | jets () const |
allow the user to access the internally stored _jets() array, which contains both the initial particles and the various intermediate and final stages of recombination. More... | |
const std::vector< history_element > & | history () const |
allow the user to access the raw internal history. More... | |
unsigned int | n_particles () const |
returns the number of particles that were provided to the clustering algorithm (helps the user find their way around the history and jets objects if they weren't paying attention beforehand). More... | |
std::vector< int > | particle_jet_indices (const std::vector< PseudoJet > &) const |
returns a vector of size n_particles() which indicates, for each of the initial particles (in the order in which they were supplied), which of the supplied jets it belongs to; if it does not belong to any of the supplied jets, the index is set to -1; | |
std::vector< int > | unique_history_order () const |
routine that returns an order in which to read the history such that clusterings that lead to identical jet compositions but different histories (because of degeneracies in the clustering order) will have matching constituents for each matching entry in the unique_history_order. More... | |
std::vector< PseudoJet > | unclustered_particles () const |
return the set of particles that have not been clustered. More... | |
std::vector< PseudoJet > | childless_pseudojets () const |
Return the list of pseudojets in the ClusterSequence that do not have children (and are not among the inclusive jets). More... | |
bool | contains (const PseudoJet &object) const |
returns true if the object (jet or particle) is contained by (ie belongs to) this cluster sequence. More... | |
void | transfer_from_sequence (const ClusterSequence &from_seq, const FunctionOfPseudoJet< PseudoJet > *action_on_jets=0) |
transfer the sequence contained in other_seq into our own; any plugin "extras" contained in the from_seq will be lost from there. More... | |
const SharedPtr< PseudoJetStructureBase > & | structure_shared_ptr () const |
retrieve a shared pointer to the wrapper to this ClusterSequence More... | |
template<> | |
void | _bj_set_jetinfo (EEBriefJet *const jetA, const int _jets_index) const |
template<> | |
double | _bj_dist (const EEBriefJet *const jeta, const EEBriefJet *const jetb) const |
Static Public Member Functions | |
static void | print_banner () |
This is the function that is automatically called during clustering to print the FastJet banner. More... | |
static std::ostream * | fastjet_banner_stream () |
returns a pointer to the stream to be used to print banners (cout by default). More... | |
Protected Member Functions | |
template<class L > | |
void | _transfer_input_jets (const std::vector< L > &pseudojets) |
transfer the vector<L> of input jets into our own vector<PseudoJet> _jets (with some reserved space for future growth). More... | |
void | _initialise_and_run (const JetDefinition &jet_def, const bool &writeout_combinations) |
This is what is called to do all the initialisation and then run the clustering (may be called by various constructors). More... | |
void | _initialise_and_run_no_decant () |
void | _decant_options (const JetDefinition &jet_def, const bool &writeout_combinations) |
fills in the various member variables with "decanted" options from the jet_definition and writeout_combinations variables | |
void | _decant_options_partial () |
assuming that the jet definition, writeout_combinations and _structure_shared_ptr have been set (e.g. More... | |
void | _fill_initial_history () |
fill out the history (and jet cross refs) related to the initial set of jets (assumed already to have been "transferred"), without any clustering | |
void | _do_ij_recombination_step (const int jet_i, const int jet_j, const double dij, int &newjet_k) |
carry out the recombination between the jets numbered jet_i and jet_j, at distance scale dij; return the index newjet_k of the result of the recombination of i and j. More... | |
void | _do_iB_recombination_step (const int jet_i, const double diB) |
carry out an recombination step in which _jets[jet_i] merges with the beam, More... | |
void | _set_structure_shared_ptr (PseudoJet &j) |
every time a jet is added internally during clustering, this should be called to set the jet's structure shared ptr to point to the CS (and the count of internally associated objects is also updated). More... | |
void | _update_structure_use_count () |
make sure that the CS's internal tally of the use count matches that of the _structure_shared_ptr | |
Strategy | _best_strategy () const |
returns a suggestion for the best strategy to use on event multiplicity, algorithm, R, etc. More... | |
void | get_subhist_set (std::set< const history_element *> &subhist, const PseudoJet &jet, double dcut, int maxjet) const |
set subhist to be a set pointers to history entries corresponding to the subjets of this jet; one stops going working down through the subjets either when More... | |
Protected Attributes | |
JetDefinition | _jet_def |
std::vector< PseudoJet > | _jets |
This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in the _history by looking at _jets[i].cluster_hist_index(). More... | |
std::vector< history_element > | _history |
this vector will contain the branching history; for each stage, _history[i].jetp_index indicates where to look in the _jets vector to get the physical PseudoJet. More... | |
bool | _writeout_combinations |
int | _initial_n |
double | _Rparam |
double | _R2 |
double | _invR2 |
double | _Qtot |
Strategy | _strategy |
JetAlgorithm | _jet_algorithm |
SharedPtr< PseudoJetStructureBase > | _structure_shared_ptr |
int | _structure_use_count_after_construction |
bool | _deletes_self_when_unused |
if true then the CS will delete itself when the last external object referring to it disappears. More... | |
deals with clustering
Definition at line 63 of file ClusterSequence.hh.
fastjet::ClusterSequence::ClusterSequence | ( | const std::vector< L > & | pseudojets, |
const JetDefinition & | jet_def, | ||
const bool & | writeout_combinations = false |
||
) |
create a ClusterSequence, starting from the supplied set of PseudoJets and clustering them with jet definition specified by jet_def (which also specifies the clustering strategy)
constructor of a jet-clustering sequence from a vector of four-momenta, with the jet definition specified by jet_def
Definition at line 972 of file ClusterSequence.hh.
vector< PseudoJet > fastjet::ClusterSequence::inclusive_jets | ( | const double | ptmin = 0.0 | ) | const |
return a vector of all jets (in the sense of the inclusive algorithm) with pt >= ptmin.
Time taken should be of the order of the number of jets returned.
Definition at line 876 of file ClusterSequence.cc.
int fastjet::ClusterSequence::n_exclusive_jets | ( | const double | dcut | ) | const |
return the number of jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.
Definition at line 928 of file ClusterSequence.cc.
vector< PseudoJet > fastjet::ClusterSequence::exclusive_jets | ( | const double | dcut | ) | const |
return a vector of all jets (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.
Definition at line 947 of file ClusterSequence.cc.
vector< PseudoJet > fastjet::ClusterSequence::exclusive_jets | ( | const int | njets | ) | const |
return a vector of all jets when the event is clustered (in the exclusive sense) to exactly njets.
If there are fewer than njets particles in the ClusterSequence an error is thrown
Definition at line 956 of file ClusterSequence.cc.
vector< PseudoJet > fastjet::ClusterSequence::exclusive_jets_up_to | ( | const int | njets | ) | const |
return a vector of all jets when the event is clustered (in the exclusive sense) to exactly njets.
If there are fewer than njets particles in the ClusterSequence the function just returns however many particles there were.
Definition at line 973 of file ClusterSequence.cc.
double fastjet::ClusterSequence::exclusive_dmerge | ( | const int | njets | ) | const |
return the dmin corresponding to the recombination that went from n+1 to n jets (sometimes known as d_{n n+1}).
return the dmin corresponding to the recombination that went from n+1 to n jets
If the number of particles in the event is <= njets, the function returns 0.
Definition at line 1041 of file ClusterSequence.cc.
double fastjet::ClusterSequence::exclusive_dmerge_max | ( | const int | njets | ) | const |
return the maximum of the dmin encountered during all recombinations up to the one that led to an n-jet final state; identical to exclusive_dmerge, except in cases where the dmin do not increase monotonically.
Definition at line 1053 of file ClusterSequence.cc.
|
inline |
return the ymin corresponding to the recombination that went from n+1 to n jets (sometimes known as y_{n n+1}).
Definition at line 137 of file ClusterSequence.hh.
std::vector< PseudoJet > fastjet::ClusterSequence::exclusive_subjets | ( | const PseudoJet & | jet, |
const double | dcut | ||
) | const |
return a vector of all subjets of the current jet (in the sense of the exclusive algorithm) that would be obtained when running the algorithm with the given dcut.
Time taken is O(m ln m), where m is the number of subjets that are found. If m gets to be of order of the total number of constituents in the jet, this could be substantially slower than just getting that list of constituents.
Definition at line 1065 of file ClusterSequence.cc.
std::vector< PseudoJet > fastjet::ClusterSequence::exclusive_subjets | ( | const PseudoJet & | jet, |
int | nsub | ||
) | const |
return the list of subjets obtained by unclustering the supplied jet down to nsub subjets.
Throws an error if there are fewer than nsub particles in the jet.
This requires nsub ln nsub time
Throws an error if there are fewer than nsub particles in the jet.
Definition at line 1101 of file ClusterSequence.cc.
std::vector< PseudoJet > fastjet::ClusterSequence::exclusive_subjets_up_to | ( | const PseudoJet & | jet, |
int | nsub | ||
) | const |
return the list of subjets obtained by unclustering the supplied jet down to nsub subjets (or all constituents if there are fewer than nsub).
This requires nsub ln nsub time
Definition at line 1118 of file ClusterSequence.cc.
double fastjet::ClusterSequence::exclusive_subdmerge | ( | const PseudoJet & | jet, |
int | nsub | ||
) | const |
returns the dij that was present in the merging nsub+1 -> nsub subjets inside this jet.
return the dij that was present in the merging nsub+1 -> nsub subjets inside this jet.
Returns 0 if there were nsub or fewer constituents in the jet.
If the jet has nsub or fewer constituents, it will return 0.
will be zero if nconst <= nsub, since highest will be an original particle have zero dij
Definition at line 1146 of file ClusterSequence.cc.
double fastjet::ClusterSequence::exclusive_subdmerge_max | ( | const PseudoJet & | jet, |
int | nsub | ||
) | const |
returns the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet.
return the maximum dij that occurred in the whole event at the stage that the nsub+1 -> nsub merge of subjets occurred inside this jet.
Returns 0 if there were nsub or fewer constituents in the jet.
If the jet has nsub or fewer constituents, it will return 0.
will be zero if nconst <= nsub, since highest will be an original particle have zero dij
Definition at line 1167 of file ClusterSequence.cc.
bool fastjet::ClusterSequence::object_in_jet | ( | const PseudoJet & | object, |
const PseudoJet & | jet | ||
) | const |
returns true iff the object is included in the jet.
NB: this is only sensible if the object is already registered within the cluster sequence, so you cannot use it with an input particle to the CS (since the particle won't have the history index set properly).
For nice clustering structures it should run in O(ln(N)) time but in worst cases (certain cone plugins) it can take O(n) time, where n is the number of particles in the jet.
Definition at line 1225 of file ClusterSequence.cc.
bool fastjet::ClusterSequence::has_parents | ( | const PseudoJet & | jet, |
PseudoJet & | parent1, | ||
PseudoJet & | parent2 | ||
) | const |
if the jet has parents in the clustering, it returns true and sets parent1 and parent2 equal to them.
if it has no parents it returns false and sets parent1 and parent2 to zero
Definition at line 1251 of file ClusterSequence.cc.
void fastjet::ClusterSequence::print_jets_for_root | ( | const std::vector< PseudoJet > & | jets, |
std::ostream & | ostr = std::cout |
||
) | const |
output the supplied vector of jets in a format that can be read by an appropriate root script; the format is: jet-n jet-px jet-py jet-pz jet-E particle-n particle-rap particle-phi particle-pt particle-n particle-rap particle-phi particle-pt ...
#END ... [i.e. above repeated]
void fastjet::ClusterSequence::delete_self_when_unused | ( | ) |
by calling this routine you tell the ClusterSequence to delete itself when all the Pseudojets associated with it have gone out of scope.
At the time you call this, there must be at least one jet or other object outside the CS that is associated with the CS (e.g. the result of inclusive_jets()).
NB: after having made this call, the user is still allowed to delete the CS. Jets associated with it will then simply not be able to access their substructure after that point.
Definition at line 1738 of file ClusterSequence.cc.
double fastjet::ClusterSequence::jet_scale_for_algorithm | ( | const PseudoJet & | jet | ) | const |
returns the scale associated with a jet as required for this clustering algorithm (kt^2 for the kt-algorithm, 1 for the Cambridge algorithm).
Intended mainly for internal use and not valid for plugin algorithms.
Definition at line 567 of file ClusterSequence.cc.
|
inline |
record the fact that there has been a recombination between jets()[jet_i] and the beam, with the specified diB; when looking for inclusive jets, any iB recombination will returned to the user as a jet.
Definition at line 338 of file ClusterSequence.hh.
|
inline |
the plugin can associate some extra information with the ClusterSequence object by calling this function.
The ClusterSequence takes ownership of the pointer (and responsibility for deleting it when the CS gets deleted).
Definition at line 360 of file ClusterSequence.hh.
|
inline |
the plugin can associate some extra information with the ClusterSequence object by calling this function
As of FJ v3.1, this is deprecated, in line with the deprecation of auto_ptr in C++11
Definition at line 372 of file ClusterSequence.hh.
|
inline |
allows a plugin to run a templated clustering (nearest-neighbour heuristic)
This has N^2 behaviour on "good" distance, but a worst case behaviour of N^3 (and many algs trigger the worst case behaviour)
For more details on how this works, see GenBriefJet below
Definition at line 391 of file ClusterSequence.hh.
|
inline |
allow the user to access the internally stored _jets() array, which contains both the initial particles and the various intermediate and final stages of recombination.
The first n_particles() entries are the original particles, in the order in which they were supplied to the ClusterSequence constructor. It can be useful to access them for example when examining whether a given input object is part of a specific jet, via the objects_in_jet(...) member function (which only takes PseudoJets that are registered in the ClusterSequence).
One of the other (internal uses) is related to the fact because we don't seem to be able to access protected elements of the class for an object that is not "this" (at least in case where "this" is of a slightly different kind from the object, both derived from ClusterSequence).
Definition at line 991 of file ClusterSequence.hh.
|
inline |
allow the user to access the raw internal history.
This is present (as for jets()) in part so that protected derived classes can access this information about other ClusterSequences.
A user who wishes to follow the details of the ClusterSequence can also make use of this information (and should consult the history_element documentation for more information), but should be aware that these internal structures may evolve in future FastJet versions.
Definition at line 995 of file ClusterSequence.hh.
|
inline |
returns the number of particles that were provided to the clustering algorithm (helps the user find their way around the history and jets objects if they weren't paying attention beforehand).
Definition at line 999 of file ClusterSequence.hh.
vector< int > fastjet::ClusterSequence::unique_history_order | ( | ) | const |
routine that returns an order in which to read the history such that clusterings that lead to identical jet compositions but different histories (because of degeneracies in the clustering order) will have matching constituents for each matching entry in the unique_history_order.
The order has the property that an entry's parents will always appear prior to that entry itself.
Roughly speaking the order is such that we first provide all steps that lead to the final jet containing particle 1; then we have the steps that lead to reconstruction of the jet containing the next-lowest-numbered unclustered particle, etc... [see GPS CCN28-12 for more info – of course a full explanation here would be better...]
Definition at line 1531 of file ClusterSequence.cc.
vector< PseudoJet > fastjet::ClusterSequence::unclustered_particles | ( | ) | const |
return the set of particles that have not been clustered.
For kt and cam/aachen algorithms this should always be null, but for cone type algorithms it can be non-null;
Definition at line 1586 of file ClusterSequence.cc.
vector< PseudoJet > fastjet::ClusterSequence::childless_pseudojets | ( | ) | const |
Return the list of pseudojets in the ClusterSequence that do not have children (and are not among the inclusive jets).
They may result from a clustering step or may be one of the pseudojets returned by unclustered_particles().
Definition at line 1600 of file ClusterSequence.cc.
bool fastjet::ClusterSequence::contains | ( | const PseudoJet & | object | ) | const |
returns true if the object (jet or particle) is contained by (ie belongs to) this cluster sequence.
Tests performed: if thejet's interface is this cluster sequence and its cluster history index is in a consistent range.
Definition at line 1615 of file ClusterSequence.cc.
void fastjet::ClusterSequence::transfer_from_sequence | ( | const ClusterSequence & | from_seq, |
const FunctionOfPseudoJet< PseudoJet > * | action_on_jets = 0 |
||
) |
transfer the sequence contained in other_seq into our own; any plugin "extras" contained in the from_seq will be lost from there.
It also sets the ClusterSequence pointers of the PseudoJets in the history to point to this ClusterSequence
When specified, the second argument is an action that will be applied on every jets in the resulting ClusterSequence
Definition at line 799 of file ClusterSequence.cc.
|
inline |
retrieve a shared pointer to the wrapper to this ClusterSequence
this may turn useful if you want to track when this ClusterSequence goes out of scope
Definition at line 537 of file ClusterSequence.hh.
|
static |
This is the function that is automatically called during clustering to print the FastJet banner.
Only the first call to this function will result in the printout of the banner. Users may wish to call this function themselves, during the initialization phase of their program, in order to ensure that the banner appears before other output. This call will not affect 3rd-party banners, e.g. those from plugins.
Definition at line 423 of file ClusterSequence.cc.
|
inlinestatic |
returns a pointer to the stream to be used to print banners (cout by default).
This function is used by plugins to determine where to direct their banners. Plugins should properly handle the case where the pointer is null.
Definition at line 574 of file ClusterSequence.hh.
|
protected |
transfer the vector<L> of input jets into our own vector<PseudoJet> _jets (with some reserved space for future growth).
Definition at line 936 of file ClusterSequence.hh.
|
protected |
This is what is called to do all the initialisation and then run the clustering (may be called by various constructors).
It assumes _jets contains the momenta to be clustered.
Definition at line 207 of file ClusterSequence.cc.
|
protected |
assuming that the jet definition, writeout_combinations and _structure_shared_ptr have been set (e.g.
in an initialiser list in the constructor), it handles the remaining decanting of options.
Definition at line 469 of file ClusterSequence.cc.
|
protected |
carry out the recombination between the jets numbered jet_i and jet_j, at distance scale dij; return the index newjet_k of the result of the recombination of i and j.
carries out the bookkeeping associated with the step of recombining jet_i and jet_j (assuming a distance dij) and returns the index of the recombined jet, newjet_k.
Definition at line 1658 of file ClusterSequence.cc.
|
protected |
carry out an recombination step in which _jets[jet_i] merges with the beam,
carries out the bookkeeping associated with the step of recombining jet_i with the beam
Definition at line 1699 of file ClusterSequence.cc.
|
protected |
every time a jet is added internally during clustering, this should be called to set the jet's structure shared ptr to point to the CS (and the count of internally associated objects is also updated).
This should not be called outside construction of a CS object.
Definition at line 1719 of file ClusterSequence.cc.
|
protected |
returns a suggestion for the best strategy to use on event multiplicity, algorithm, R, etc.
Definition at line 603 of file ClusterSequence.cc.
|
protected |
set subhist to be a set pointers to history entries corresponding to the subjets of this jet; one stops going working down through the subjets either when
return a set of pointers to history entries corresponding to the subjets of this jet; one stops going working down through the subjets either when
Definition at line 1191 of file ClusterSequence.cc.
|
protected |
This contains the physical PseudoJets; for each PseudoJet one can find the corresponding position in the _history by looking at _jets[i].cluster_hist_index().
Definition at line 682 of file ClusterSequence.hh.
|
protected |
this vector will contain the branching history; for each stage, _history[i].jetp_index indicates where to look in the _jets vector to get the physical PseudoJet.
Definition at line 688 of file ClusterSequence.hh.
|
mutableprotected |
if true then the CS will delete itself when the last external object referring to it disappears.
It is mutable so as to ensure that signal_imminent_self_deletion() [const] can make relevant changes.
Definition at line 714 of file ClusterSequence.hh.