fastjet 2.4.3
|
Class that behaves essentially like ClusterSequence except that it also provides access to the area of a jet (which will be a random quantity... More...
#include <ClusterSequenceActiveAreaExplicitGhosts.hh>
Public Member Functions | |
template<class L > | |
ClusterSequenceActiveAreaExplicitGhosts (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const GhostedAreaSpec &ghost_spec, const bool &writeout_combinations=false) | |
constructor using a GhostedAreaSpec to specify how the area is to be measured | |
template<class L > | |
ClusterSequenceActiveAreaExplicitGhosts (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const std::vector< L > &ghosts, double ghost_area, const bool &writeout_combinations=false) | |
template<class L > | |
void | _initialise (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const GhostedAreaSpec *ghost_spec, const std::vector< L > *ghosts, double ghost_area, const bool &writeout_combinations) |
does the actual work of initialisation | |
unsigned int | n_hard_particles () const |
returns the number of hard particles (i.e. those supplied by the user). | |
virtual double | area (const PseudoJet &jet) const |
returns the area of a jet | |
virtual PseudoJet | area_4vector (const PseudoJet &jet) const |
returns a four vector corresponding to the sum (E-scheme) of the ghost four-vectors composing the jet area, normalised such that for a small contiguous area the p_t of the extended_area jet is equal to area of the jet. | |
virtual bool | is_pure_ghost (const PseudoJet &jet) const |
true if a jet is made exclusively of ghosts | |
bool | is_pure_ghost (int history_index) const |
true if the entry in the history index corresponds to a ghost; if hist_ix does not correspond to an actual particle (i.e. | |
virtual bool | has_explicit_ghosts () const |
this class does have explicit ghosts | |
virtual double | empty_area (const RangeDefinition &range) const |
return the total area, up to |y|<maxrap, that consists of unclustered ghosts | |
double | total_area () const |
returns the total area under study | |
double | max_ghost_perp2 () const |
returns the largest squared transverse momentum among all ghosts | |
bool | has_dangerous_particles () const |
returns true if there are any particles whose transverse momentum if so low that there's a risk of the ghosts having modified the clustering sequence | |
Private Member Functions | |
void | _add_ghosts (const GhostedAreaSpec &ghost_spec) |
adds the "ghost" momenta, which will be used to estimate the jet area | |
template<class L > | |
void | _add_ghosts (const std::vector< L > &ghosts, double ghost_area) |
another way of adding ghosts | |
void | _post_process () |
routine to be called after the processing is done so as to establish summary information on all the jets (areas, whether pure ghost, etc.) | |
Private Attributes | |
int | _n_ghosts |
get the area of the ghosts | |
double | _ghost_area |
std::vector< bool > | _is_pure_ghost |
std::vector< double > | _areas |
std::vector< PseudoJet > | _area_4vectors |
double | _max_ghost_perp2 |
bool | _has_dangerous_particles |
unsigned int | _initial_hard_n |
Static Private Attributes | |
static LimitedWarning | _warnings |
handle warning messages |
Class that behaves essentially like ClusterSequence except that it also provides access to the area of a jet (which will be a random quantity...
Figure out what to do about seeds later...)
Definition at line 49 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
ClusterSequenceActiveAreaExplicitGhosts::ClusterSequenceActiveAreaExplicitGhosts | ( | const std::vector< L > & | pseudojets, |
const JetDefinition & | jet_def, | ||
const GhostedAreaSpec & | ghost_spec, | ||
const bool & | writeout_combinations = false |
||
) | [inline] |
constructor using a GhostedAreaSpec to specify how the area is to be measured
Definition at line 55 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
References _initialise().
: ClusterSequenceAreaBase() { std::vector<L> * ghosts = NULL; _initialise(pseudojets,jet_def,&ghost_spec,ghosts,0.0, writeout_combinations); }
ClusterSequenceActiveAreaExplicitGhosts::ClusterSequenceActiveAreaExplicitGhosts | ( | const std::vector< L > & | pseudojets, |
const JetDefinition & | jet_def, | ||
const std::vector< L > & | ghosts, | ||
double | ghost_area, | ||
const bool & | writeout_combinations = false |
||
) | [inline] |
Definition at line 65 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
References _initialise().
: ClusterSequenceAreaBase() { const GhostedAreaSpec * ghost_spec = NULL; _initialise(pseudojets,jet_def,ghost_spec,&ghosts,ghost_area, writeout_combinations); }
void ClusterSequenceActiveAreaExplicitGhosts::_add_ghosts | ( | const GhostedAreaSpec & | ghost_spec | ) | [private] |
adds the "ghost" momenta, which will be used to estimate the jet area
void ClusterSequenceActiveAreaExplicitGhosts::_add_ghosts | ( | const std::vector< L > & | ghosts, |
double | ghost_area | ||
) | [private] |
another way of adding ghosts
add an explicitly specified bunch of ghosts
Definition at line 222 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
References _ghost_area, _is_pure_ghost, ClusterSequence::_jets, and _n_ghosts.
{ for (unsigned i = 0; i < ghosts.size(); i++) { _is_pure_ghost.push_back(true); _jets.push_back(ghosts[i]); } // and record some info about ghosts _ghost_area = ghost_area; _n_ghosts = ghosts.size(); }
void ClusterSequenceActiveAreaExplicitGhosts::_initialise | ( | const std::vector< L > & | pseudojets, |
const JetDefinition & | jet_def, | ||
const GhostedAreaSpec * | ghost_spec, | ||
const std::vector< L > * | ghosts, | ||
double | ghost_area, | ||
const bool & | writeout_combinations | ||
) |
does the actual work of initialisation
Definition at line 169 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by ClusterSequenceActiveAreaExplicitGhosts().
{ // don't reserve space yet -- will be done below // insert initial jets this way so that any type L that can be // converted to a pseudojet will work fine (basically PseudoJet // and any type that has [] subscript access to the momentum // components, such as CLHEP HepLorentzVector). for (unsigned int i = 0; i < pseudojets.size(); i++) { PseudoJet mom(pseudojets[i]); //mom.set_user_index(0); // for user's particles (user index now lost...) _jets.push_back(mom); _is_pure_ghost.push_back(false); } _initial_hard_n = _jets.size(); if (ghost_spec != NULL) { _add_ghosts(*ghost_spec); } else { _add_ghosts(*ghosts, ghost_area); } if (writeout_combinations) { std::cout << "# Printing particles including ghosts\n"; for (unsigned j = 0; j < _jets.size(); j++) { printf("%5u %20.13f %20.13f %20.13e\n", j,_jets[j].rap(),_jets[j].phi_02pi(),_jets[j].kt2()); } std::cout << "# Finished printing particles including ghosts\n"; } // this will ensure that we can still point to jets without // difficulties arising! _jets.reserve(_jets.size()*2); // run the clustering _initialise_and_run(jet_def,writeout_combinations); // set up all other information _post_process(); }
void ClusterSequenceActiveAreaExplicitGhosts::_post_process | ( | ) | [private] |
routine to be called after the processing is done so as to establish summary information on all the jets (areas, whether pure ghost, etc.)
virtual double ClusterSequenceActiveAreaExplicitGhosts::area | ( | const PseudoJet & | jet | ) | const [virtual] |
returns the area of a jet
Reimplemented from ClusterSequenceAreaBase.
Referenced by ClusterSequenceActiveArea::_transfer_areas().
virtual PseudoJet ClusterSequenceActiveAreaExplicitGhosts::area_4vector | ( | const PseudoJet & | jet | ) | const [virtual] |
returns a four vector corresponding to the sum (E-scheme) of the ghost four-vectors composing the jet area, normalised such that for a small contiguous area the p_t of the extended_area jet is equal to area of the jet.
Reimplemented from ClusterSequenceAreaBase.
Referenced by ClusterSequenceActiveArea::_transfer_areas().
virtual double ClusterSequenceActiveAreaExplicitGhosts::empty_area | ( | const RangeDefinition & | range | ) | const [virtual] |
return the total area, up to |y|<maxrap, that consists of unclustered ghosts
Reimplemented from ClusterSequenceAreaBase.
bool ClusterSequenceActiveAreaExplicitGhosts::has_dangerous_particles | ( | ) | const [inline] |
returns true if there are any particles whose transverse momentum if so low that there's a risk of the ghosts having modified the clustering sequence
Definition at line 124 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
References _has_dangerous_particles.
Referenced by ClusterSequenceActiveArea::_run_AA(), and ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E().
{return _has_dangerous_particles;}
virtual bool ClusterSequenceActiveAreaExplicitGhosts::has_explicit_ghosts | ( | ) | const [inline, virtual] |
this class does have explicit ghosts
Reimplemented from ClusterSequenceAreaBase.
Definition at line 108 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
{return true;}
bool ClusterSequenceActiveAreaExplicitGhosts::is_pure_ghost | ( | int | history_index | ) | const |
true if the entry in the history index corresponds to a ghost; if hist_ix does not correspond to an actual particle (i.e.
hist_ix < 0), then the result is false.
virtual bool ClusterSequenceActiveAreaExplicitGhosts::is_pure_ghost | ( | const PseudoJet & | jet | ) | const [virtual] |
true if a jet is made exclusively of ghosts
Reimplemented from ClusterSequenceAreaBase.
Referenced by ClusterSequenceActiveArea::_transfer_areas(), and ClusterSequenceActiveArea::_transfer_ghost_free_history().
double ClusterSequenceActiveAreaExplicitGhosts::max_ghost_perp2 | ( | ) | const [inline] |
returns the largest squared transverse momentum among all ghosts
Definition at line 119 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
References _max_ghost_perp2.
{return _max_ghost_perp2;}
unsigned int ClusterSequenceActiveAreaExplicitGhosts::n_hard_particles | ( | ) | const [inline] |
returns the number of hard particles (i.e. those supplied by the user).
Definition at line 217 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
References _initial_hard_n.
{return _initial_hard_n;}
double ClusterSequenceActiveAreaExplicitGhosts::total_area | ( | ) | const |
returns the total area under study
std::vector<PseudoJet> ClusterSequenceActiveAreaExplicitGhosts::_area_4vectors [private] |
Definition at line 135 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
std::vector<double> ClusterSequenceActiveAreaExplicitGhosts::_areas [private] |
Definition at line 134 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
double ClusterSequenceActiveAreaExplicitGhosts::_ghost_area [private] |
Definition at line 132 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by _add_ghosts().
Definition at line 139 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by has_dangerous_particles().
unsigned int ClusterSequenceActiveAreaExplicitGhosts::_initial_hard_n [private] |
Definition at line 146 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by n_hard_particles().
std::vector<bool> ClusterSequenceActiveAreaExplicitGhosts::_is_pure_ghost [private] |
Definition at line 133 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by _add_ghosts().
double ClusterSequenceActiveAreaExplicitGhosts::_max_ghost_perp2 [private] |
Definition at line 138 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by max_ghost_perp2().
int ClusterSequenceActiveAreaExplicitGhosts::_n_ghosts [private] |
get the area of the ghosts
Definition at line 131 of file ClusterSequenceActiveAreaExplicitGhosts.hh.
Referenced by _add_ghosts().
LimitedWarning ClusterSequenceActiveAreaExplicitGhosts::_warnings [static, private] |
handle warning messages
allow for warnings
Reimplemented from ClusterSequenceAreaBase.
Definition at line 140 of file ClusterSequenceActiveAreaExplicitGhosts.hh.