class that is intended to hold a full definition of the jet clusterer More...
#include <JetDefinition.hh>
Classes | |
class | DefaultRecombiner |
A class that will provide the recombination scheme facilities and/or allow a user to extend these facilities. More... | |
class | Plugin |
a class that allows a user to introduce their own "plugin" jet finder More... | |
class | Recombiner |
An abstract base class that will provide the recombination scheme facilities and/or allow a user to extend these facilities. More... | |
Public Member Functions | |
JetDefinition (JetAlgorithm jet_algorithm, double R, RecombinationScheme recomb_scheme=E_scheme, Strategy strategy=Best) | |
constructor with alternative ordering or arguments -- note that we have not provided a default jet finder, to avoid ambiguous JetDefinition() constructor. | |
JetDefinition (JetAlgorithm jet_algorithm, RecombinationScheme recomb_scheme=E_scheme, Strategy strategy=Best) | |
constructor for algorithms that have no free parameters (e.g. | |
JetDefinition (JetAlgorithm jet_algorithm, double R, double xtra_param, RecombinationScheme recomb_scheme=E_scheme, Strategy strategy=Best) | |
constructor for algorithms that require R + one extra parameter to be set (the gen-kt series for example) | |
JetDefinition (JetAlgorithm jet_algorithm, double R, const Recombiner *recombiner, Strategy strategy=Best) | |
constructor in a form that allows the user to provide a pointer to an external recombiner class (which must remain valid for the life of the JetDefinition object). | |
JetDefinition (JetAlgorithm jet_algorithm, const Recombiner *recombiner, Strategy strategy=Best) | |
constructor for case with 0 parameters (ee_kt_algorithm) and and external recombiner | |
JetDefinition (JetAlgorithm jet_algorithm, double R, double xtra_param, const Recombiner *recombiner, Strategy strategy=Best) | |
constructor allowing the extra parameter to be set and a pointer to a recombiner | |
JetDefinition () | |
a default constructor | |
JetDefinition (const Plugin *plugin) | |
constructor based on a pointer to a user's plugin; the object pointed to must remain valid for the whole duration of existence of the JetDefinition and any related ClusterSequences | |
JetDefinition (JetAlgorithm jet_algorithm, double R, Strategy strategy, RecombinationScheme recomb_scheme=E_scheme, int nparameters=1) | |
constructor to fully specify a jet-definition (together with information about how algorithically to run it). | |
void | set_recombination_scheme (RecombinationScheme) |
set the recombination scheme to the one provided | |
void | set_recombiner (const Recombiner *recomb) |
set the recombiner class to the one provided | |
const Plugin * | plugin () const |
return a pointer to the plugin | |
JetAlgorithm | jet_algorithm () const |
return information about the definition... | |
JetAlgorithm | jet_finder () const |
same as above for backward compatibility | |
double | R () const |
double | extra_param () const |
Strategy | strategy () const |
RecombinationScheme | recombination_scheme () const |
void | set_jet_algorithm (JetAlgorithm njf) |
(re)set the jet finder | |
void | set_jet_finder (JetAlgorithm njf) |
same as above for backward compatibility | |
void | set_extra_param (double xtra_param) |
(re)set the general purpose extra parameter | |
const Recombiner * | recombiner () const |
return a pointer to the currently defined recombiner. | |
std::string | description () const |
return a textual description of the current jet definition |
class that is intended to hold a full definition of the jet clusterer
Definition at line 172 of file JetDefinition.hh.
fastjet::JetDefinition::JetDefinition | ( | JetAlgorithm | jet_algorithm, | |
double | R, | |||
RecombinationScheme | recomb_scheme = E_scheme , |
|||
Strategy | strategy = Best | |||
) | [inline] |
constructor with alternative ordering or arguments -- note that we have not provided a default jet finder, to avoid ambiguous JetDefinition() constructor.
Definition at line 190 of file JetDefinition.hh.
{ *this = JetDefinition(jet_algorithm, R, strategy, recomb_scheme, 1); }
fastjet::JetDefinition::JetDefinition | ( | JetAlgorithm | jet_algorithm, | |
RecombinationScheme | recomb_scheme = E_scheme , |
|||
Strategy | strategy = Best | |||
) | [inline] |
constructor for algorithms that have no free parameters (e.g.
ee_kt_algorithm)
Definition at line 199 of file JetDefinition.hh.
{ double dummyR = 0.0; *this = JetDefinition(jet_algorithm, dummyR, strategy, recomb_scheme, 0); }
fastjet::JetDefinition::JetDefinition | ( | JetAlgorithm | jet_algorithm, | |
double | R, | |||
const Recombiner * | recombiner, | |||
Strategy | strategy = Best | |||
) | [inline] |
constructor in a form that allows the user to provide a pointer to an external recombiner class (which must remain valid for the life of the JetDefinition object).
Definition at line 221 of file JetDefinition.hh.
References fastjet::external_scheme.
{ *this = JetDefinition(jet_algorithm, R, external_scheme, strategy); _recombiner = recombiner; }
fastjet::JetDefinition::JetDefinition | ( | JetAlgorithm | jet_algorithm, | |
double | R, | |||
Strategy | strategy, | |||
RecombinationScheme | recomb_scheme = E_scheme , |
|||
int | nparameters = 1 | |||
) |
constructor to fully specify a jet-definition (together with information about how algorithically to run it).
the ordering of arguments here is old and deprecated (except as the common constructor for internal use)
Definition at line 42 of file JetDefinition.cc.
References fastjet::ee_genkt_algorithm, fastjet::ee_kt_algorithm, fastjet::genkt_algorithm, fastjet::plugin_strategy, set_extra_param(), and set_recombination_scheme().
: _jet_algorithm(jet_algorithm), _Rparam(R), _strategy(strategy) { // set R parameter or ensure its sensibleness, as appropriate if (jet_algorithm == ee_kt_algorithm) { _Rparam = 4.0; // introduce a fictional R that ensures that // our clustering sequence will not produce // "beam" jets except when only a single particle remains. // Any value > 2 would have done here } else { // We maintain some limit on R because particles with pt=0, m=0 // can have rapidities O(10000) and one doesn't want the // clustering to start including them as if their rapidities were // physical. assert ( R < 1000.0); } // cross-check the number of parameters that were declared in setting up the // algorithm (passed internally from the public constructors) ostringstream oss; switch (jet_algorithm) { case ee_kt_algorithm: if (nparameters != 0) oss << "ee_kt_algorithm should be constructed with 0 parameters but was called with " << nparameters << " parameter(s)\n"; break; case genkt_algorithm: case ee_genkt_algorithm: if (nparameters != 2) oss << "(ee_)genkt_algorithm should be constructed with 2 parameters but was called with " << nparameters << " parameter(s)\n"; break; default: if (nparameters != 1) oss << "The jet algorithm you requested (" << jet_algorithm << ") should be constructed with 1 parameter but was called with " << nparameters << " parameter(s)\n"; } if (oss.str() != "") throw Error(oss.str()); // make sure the strategy requested is sensible assert (_strategy != plugin_strategy); _plugin = NULL; set_recombination_scheme(recomb_scheme); set_extra_param(0.0); // make sure it's defined }
const Recombiner* fastjet::JetDefinition::recombiner | ( | ) | const [inline] |
return a pointer to the currently defined recombiner.
Warning: the pointer may be to an internal recombiner (for default recombination schemes), in which case if the JetDefinition becomes invalid (e.g. is deleted), the pointer will then point to an object that no longer exists.
Note also that if you copy a JetDefinition with a default recombination scheme, then the two copies will have distinct recombiners, and return different recombiner() pointers.
Definition at line 332 of file JetDefinition.hh.
Referenced by description().
{
return _recombiner == 0 ? & _default_recombiner : _recombiner;}