FastJet
3.1.3
|
Like ClusterSequence with computation of the active jet area. More...
#include <fastjet/ClusterSequenceActiveArea.hh>
Public Types | |
enum | mean_pt_strategies { median =0, non_ghost_median, pttot_over_areatot, pttot_over_areatot_cut, mean_ratio_cut, play, median_4vector } |
enum providing a variety of tentative strategies for estimating the background (e.g. More... | |
Public Types inherited from fastjet::ClusterSequence | |
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 | |
ClusterSequenceActiveArea () | |
default constructor | |
template<class L > | |
ClusterSequenceActiveArea (const std::vector< L > &pseudojets, const JetDefinition &jet_def_in, const GhostedAreaSpec &ghost_spec, const bool &writeout_combinations=false) | |
constructor based on JetDefinition and GhostedAreaSpec | |
virtual double | area (const PseudoJet &jet) const |
return the area associated with the given jet; this base class returns 0. More... | |
virtual double | area_error (const PseudoJet &jet) const |
return the error (uncertainty) associated with the determination of the area of this jet; this base class returns 0. More... | |
virtual PseudoJet | area_4vector (const PseudoJet &jet) const |
return a PseudoJet whose 4-vector is defined by the following integral More... | |
double | pt_per_unit_area (mean_pt_strategies strat=median, double range=2.0) const |
return the transverse momentum per unit area according to one of the above strategies; for some strategies (those with "cut" in their name) the parameter "range" allows one to exclude a subset of the jets for the background estimation, those that have pt/area > median(pt/area)*range. More... | |
virtual double | empty_area (const Selector &selector) const |
rewrite the empty area from the parent class, so as to use all info at our disposal return the total area, corresponding to a given Selector, that consists of ghost jets or unclustered ghosts More... | |
virtual double | n_empty_jets (const Selector &selector) const |
return the true number of empty jets (replaces ClusterSequenceAreaBase::n_empty_jets(...)) | |
Public Member Functions inherited from fastjet::ClusterSequenceAreaBase | |
template<class L > | |
ClusterSequenceAreaBase (const std::vector< L > &pseudojets, const JetDefinition &jet_def_in, const bool &writeout_combinations=false) | |
a constructor which just carries out the construction of the parent class | |
ClusterSequenceAreaBase () | |
default constructor | |
virtual | ~ClusterSequenceAreaBase () |
destructor | |
virtual bool | is_pure_ghost (const PseudoJet &) const |
true if a jet is made exclusively of ghosts More... | |
virtual bool | has_explicit_ghosts () const |
returns true if ghosts are explicitly included within jets for this ClusterSequence; More... | |
double | empty_area_from_jets (const std::vector< PseudoJet > &all_jets, const Selector &selector) const |
return the total area, corresponding to the given Selector, that is free of jets, based on the supplied all_jets More... | |
double | median_pt_per_unit_area (const Selector &selector) const |
the median of (pt/area) for jets contained within the selector range, making use also of the info on n_empty_jets More... | |
double | median_pt_per_unit_area_4vector (const Selector &selector) const |
the median of (pt/area_4vector) for jets contained within the selector range, making use also of the info on n_empty_jets More... | |
double | median_pt_per_unit_something (const Selector &selector, bool use_area_4vector) const |
the function that does the work for median_pt_per_unit_area and median_pt_per_unit_area_4vector: More... | |
virtual void | get_median_rho_and_sigma (const Selector &selector, bool use_area_4vector, double &median, double &sigma, double &mean_area) const |
using jets withing the selector range (and with 4-vector areas if use_area_4vector), calculate the median pt/area, as well as an "error" (uncertainty), which is defined as the 1-sigma half-width of the distribution of pt/A, obtained by looking for the point below which we have (1-0.6827)/2 of the jets (including empty jets). More... | |
virtual void | get_median_rho_and_sigma (const std::vector< PseudoJet > &all_jets, const Selector &selector, bool use_area_4vector, double &median, double &sigma, double &mean_area, bool all_are_inclusive=false) const |
a more advanced version of get_median_rho_and_sigma, which allows one to use any "view" of the event containing all jets (so that, e.g. More... | |
virtual void | get_median_rho_and_sigma (const Selector &selector, bool use_area_4vector, double &median, double &sigma) const |
same as the full version of get_median_rho_and_error, but without access to the mean_area More... | |
virtual void | parabolic_pt_per_unit_area (double &a, double &b, const Selector &selector, double exclude_above=-1.0, bool use_area_4vector=false) const |
fits a form pt_per_unit_area(y) = a + b*y^2 in the selector range. More... | |
std::vector< PseudoJet > | subtracted_jets (const double rho, const double ptmin=0.0) const |
return a vector of all subtracted jets, using area_4vector, given rho. More... | |
std::vector< PseudoJet > | subtracted_jets (const Selector &selector, const double ptmin=0.0) const |
return a vector of subtracted jets, using area_4vector. More... | |
PseudoJet | subtracted_jet (const PseudoJet &jet, const double rho) const |
return a subtracted jet, using area_4vector, given rho | |
PseudoJet | subtracted_jet (const PseudoJet &jet, const Selector &selector) const |
return a subtracted jet, using area_4vector; note that this is potentially inefficient if repeatedly used for many different jets, because rho will be recalculated each time around. More... | |
double | subtracted_pt (const PseudoJet &jet, const double rho, bool use_area_4vector=false) const |
return the subtracted pt, given rho | |
double | subtracted_pt (const PseudoJet &jet, const Selector &selector, bool use_area_4vector=false) const |
return the subtracted pt; note that this is potentially inefficient if repeatedly used for many different jets, because rho will be recalculated each time around. More... | |
Public Member Functions inherited from fastjet::ClusterSequence | |
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 | |
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 |
Protected Member Functions | |
void | _resize_and_zero_AA () |
void | _initialise_AA (const JetDefinition &jet_def, const GhostedAreaSpec &ghost_spec, const bool &writeout_combinations, bool &continue_running) |
void | _run_AA (const GhostedAreaSpec &ghost_spec) |
void | _postprocess_AA (const GhostedAreaSpec &ghost_spec) |
run the postprocessing for the active area (and derived classes) | |
void | _initialise_and_run_AA (const JetDefinition &jet_def, const GhostedAreaSpec &ghost_spec, const bool &writeout_combinations=false) |
does the initialisation and running specific to the active areas class More... | |
void | _transfer_ghost_free_history (const ClusterSequenceActiveAreaExplicitGhosts &clust_seq) |
transfer the history (and jet-momenta) from clust_seq to our own internal structure while removing ghosts | |
void | _transfer_areas (const std::vector< int > &unique_hist_order, const ClusterSequenceActiveAreaExplicitGhosts &) |
transfer areas from the ClusterSequenceActiveAreaExplicitGhosts object into our internal area bookkeeping... More... | |
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 | |
Protected Member Functions inherited from fastjet::ClusterSequenceAreaBase | |
void | _check_selector_good_for_median (const Selector &selector) const |
check the selector is suited for the computations i.e. applies jet by jet and has a finite area | |
Protected Member Functions inherited from fastjet::ClusterSequence | |
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 | |
std::valarray< double > | _average_area |
child classes benefit from having these at their disposal | |
std::valarray< double > | _average_area2 |
std::valarray< PseudoJet > | _average_area_4vector |
Protected Attributes inherited from fastjet::ClusterSequence | |
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... | |
Additional Inherited Members | |
Static Public Member Functions inherited from fastjet::ClusterSequence | |
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... | |
Like ClusterSequence with computation of the active jet area.
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...)
This class should not be used directly. Rather use ClusterSequenceArea with the appropriate AreaDefinition
Definition at line 63 of file ClusterSequenceActiveArea.hh.
enum providing a variety of tentative strategies for estimating the background (e.g.
non-jet) activity in a highly populated event; the one that has been most extensively tested is median.
These strategies are OBSOLETE and deprecated (see comment for pt_per_unit_area).
Definition at line 90 of file ClusterSequenceActiveArea.hh.
|
inlinevirtual |
return the area associated with the given jet; this base class returns 0.
Reimplemented from fastjet::ClusterSequenceAreaBase.
Definition at line 76 of file ClusterSequenceActiveArea.hh.
|
inlinevirtual |
return the error (uncertainty) associated with the determination of the area of this jet; this base class returns 0.
Reimplemented from fastjet::ClusterSequenceAreaBase.
Definition at line 78 of file ClusterSequenceActiveArea.hh.
|
inlinevirtual |
return a PseudoJet whose 4-vector is defined by the following integral
drap d PseudoJet("rap,phi,pt=one") *
where PseudoJet("rap,phi,pt=one") is a 4-vector with the given rapidity (rap), azimuth (phi) and pt=1, while Theta("rap,phi inside jet boundary") is a function that is 1 when rap,phi define a direction inside the jet boundary and 0 otherwise.
This base class returns a null 4-vector.
Reimplemented from fastjet::ClusterSequenceAreaBase.
Definition at line 81 of file ClusterSequenceActiveArea.hh.
double fastjet::ClusterSequenceActiveArea::pt_per_unit_area | ( | mean_pt_strategies | strat = median , |
double | range = 2.0 |
||
) | const |
return the transverse momentum per unit area according to one of the above strategies; for some strategies (those with "cut" in their name) the parameter "range" allows one to exclude a subset of the jets for the background estimation, those that have pt/area > median(pt/area)*range.
NB: This call is OBSOLETE and deprecated; use a JetMedianBackgroundEstimator or GridMedianBackgroundEstimator instead.
Definition at line 276 of file ClusterSequenceActiveArea.cc.
|
virtual |
rewrite the empty area from the parent class, so as to use all info at our disposal return the total area, corresponding to a given Selector, that consists of ghost jets or unclustered ghosts
The selector passed as an argument needs to apply jet by jet.
Reimplemented from fastjet::ClusterSequenceAreaBase.
Reimplemented in fastjet::ClusterSequencePassiveArea.
Definition at line 445 of file ClusterSequenceActiveArea.cc.
|
protected |
does the initialisation and running specific to the active areas class
global routine for running active area
Definition at line 53 of file ClusterSequenceActiveArea.cc.
|
protected |
transfer areas from the ClusterSequenceActiveAreaExplicitGhosts object into our internal area bookkeeping...
Definition at line 566 of file ClusterSequenceActiveArea.cc.