31#include "fastjet/JetDefinition.hh"
32#include "fastjet/Error.hh"
33#include "fastjet/CompositeJetStructure.hh"
36FASTJET_BEGIN_NAMESPACE
40const double JetDefinition::max_allowable_R = 1000.0;
50 _jet_algorithm(jet_algorithm_in), _Rparam(R_in), _strategy(strategy_in) {
65 oss <<
"Requested R = " << R_in <<
" for jet definition is larger than max_allowable_R = " <<
max_allowable_R;
66 throw Error(oss.str());
73 if (nparameters != (
int) nparameters_expected){
75 oss <<
"The jet algorithm you requested ("
76 << jet_algorithm_in <<
") should be constructed with " << nparameters_expected
77 <<
" parameter(s) but was called with " << nparameters <<
" parameter(s)\n";
78 throw Error(oss.str());
131 return "uninitialised JetDefinition (jet_algorithm=undefined_jet_algorithm)" ;
136 case 0: name <<
" (NB: no R)";
break;
137 case 1: name <<
" with R = " << R();
break;
140 name <<
" with R = " << R();
143 name <<
"and a special hack whereby particles with kt < "
144 << extra_param() <<
"are treated as passive ghosts";
146 name <<
", p = " << extra_param();
158 case kt_algorithm:
return "Longitudinally invariant kt algorithm";
161 case genkt_algorithm:
return "Longitudinally invariant generalised kt algorithm";
167 throw Error(
"JetDefinition::algorithm_description(): unrecognized jet_algorithm");
187 if (_shared_recombiner) _shared_recombiner.reset();
194 assert(other_jet_def._recombiner ||
198 if (other_jet_def._recombiner == 0){
204 _recombiner = other_jet_def._recombiner;
209 _shared_recombiner = other_jet_def._shared_recombiner;
225 if (other_jd.recombination_scheme() != scheme)
return false;
236 if (_recombiner == 0){
237 throw Error(
"tried to call JetDefinition::delete_recombiner_when_unused() for a JetDefinition without a user-defined recombination scheme");
238 }
else if (_shared_recombiner.get()) {
239 throw Error(
"Error in JetDefinition::delete_recombiner_when_unused: the recombiner is already scheduled for deletion when unused (or was already set as shared)");
242 _shared_recombiner.reset(_recombiner);
249 throw Error(
"tried to call JetDefinition::delete_plugin_when_unused() for a JetDefinition without a plugin");
252 _plugin_shared.reset(_plugin);
258 switch(_recomb_scheme) {
260 return "E scheme recombination";
262 return "pt scheme recombination";
264 return "pt2 scheme recombination";
266 return "Et scheme recombination";
268 return "Et2 scheme recombination";
270 return "boost-invariant pt scheme recombination";
272 return "boost-invariant pt2 scheme recombination";
274 return "pt-ordered Winner-Takes-All recombination";
281 return "|3-momentum|-ordered Winner-Takes-All recombination";
284 err <<
"DefaultRecombiner: unrecognized recombination scheme "
286 throw Error(err.str());
295 double weighta, weightb;
297 switch(_recomb_scheme) {
302 pab.
reset(pa.px()+pb.px(),
318 weighta = pa.
perp2();
319 weightb = pb.
perp2();
325 phard.
rap(), phard.
phi(), phard.
m());
351 const PseudoJet & phard = a_hardest ? pa : pb;
352 const PseudoJet & psoft = a_hardest ? pb : pa;
357 double modp_hard = phard.
modp();
358 double modp_ab = modp_hard + psoft.
modp();
359 if (phard.
modp2()==0.0){
360 pab.
reset(0.0, 0.0, 0.0, phard.
m());
362 double scale = modp_ab/modp_hard;
363 pab.
reset(phard.px()*scale, phard.py()*scale, phard.pz()*scale,
364 sqrt(modp_ab*modp_ab + phard.
m2()));
369 err <<
"DefaultRecombiner: unrecognized recombination scheme "
371 throw Error(err.str());
374 double perp_ab = pa.
perp() + pb.
perp();
375 if (perp_ab != 0.0) {
376 double y_ab = (weighta * pa.
rap() + weightb * pb.
rap())/(weighta+weightb);
379 double phi_a = pa.
phi(), phi_b = pb.
phi();
380 if (phi_a - phi_b > pi) phi_b += twopi;
381 if (phi_a - phi_b < -pi) phi_b -= twopi;
382 double phi_ab = (weighta * phi_a + weightb * phi_b)/(weighta+weightb);
392 pab.
reset(0.0, 0.0, 0.0, 0.0);
398 switch(_recomb_scheme) {
411 double newE = sqrt(p.
perp2()+p.pz()*p.pz());
424 double rescale = p.E()/sqrt(p.
perp2()+p.pz()*p.pz());
425 p.
reset_momentum(rescale*p.px(), rescale*p.py(), rescale*p.pz(), p.E());
434 err <<
"DefaultRecombiner: unrecognized recombination scheme "
436 throw Error(err.str());
441 throw Error(
"set_ghost_separation_scale not supported");
460 if (pieces.size()>0){
462 for (
unsigned int i=1; i<pieces.size(); i++)
478 return join(vector<PseudoJet>(1,j1),
recombiner);
484 vector<PseudoJet> pieces;
485 pieces.push_back(j1);
486 pieces.push_back(j2);
493 vector<PseudoJet> pieces;
494 pieces.push_back(j1);
495 pieces.push_back(j2);
496 pieces.push_back(j3);
503 vector<PseudoJet> pieces;
504 pieces.push_back(j1);
505 pieces.push_back(j2);
506 pieces.push_back(j3);
507 pieces.push_back(j4);
The structure for a jet made of pieces.
base class corresponding to errors that can be thrown by FastJet
A class that will provide the recombination scheme facilities and/or allow a user to extend these fac...
virtual void recombine(const PseudoJet &pa, const PseudoJet &pb, PseudoJet &pab) const override
recombine pa and pb and put result into pab
virtual std::string description() const override
return a textual description of the recombination scheme implemented here
virtual void preprocess(PseudoJet &p) const override
routine called to preprocess each input jet (to make all input jets compatible with the scheme requir...
virtual void set_ghost_separation_scale(double scale) const
set the ghost separation scale for passive area determinations in future runs (strictly speaking that...
virtual bool is_spherical() const
returns true if the plugin implements an algorithm intended for use on a spherical geometry (e....
virtual std::string description() const =0
return a textual description of the jet-definition implemented in this plugin
An abstract base class that will provide the recombination scheme facilities and/or allow a user to e...
virtual std::string description() const =0
return a textual description of the recombination scheme implemented here
void plus_equal(PseudoJet &pa, const PseudoJet &pb) const
pa += pb in the given recombination scheme.
class that is intended to hold a full definition of the jet clusterer
std::string description() const
return a textual description of the current jet definition
bool is_spherical() const
returns true if the jet definition involves an algorithm intended for use on a spherical geometry (e....
const Plugin * plugin() const
return a pointer to the plugin
JetAlgorithm jet_algorithm() const
return information about the definition...
static std::string algorithm_description(const JetAlgorithm jet_alg)
a short textual description of the algorithm jet_alg
const Recombiner * recombiner() const
returns a pointer to the currently defined recombiner.
static unsigned int n_parameters_for_algorithm(const JetAlgorithm jet_alg)
the number of parameters associated to a given jet algorithm
void delete_recombiner_when_unused()
calling this tells the JetDefinition to handle the deletion of the recombiner when it is no longer us...
void set_extra_param(double xtra_param)
(re)set the general purpose extra parameter
bool has_same_recombiner(const JetDefinition &other_jd) const
returns true if the current jet definitions shares the same recombiner as the one passed as an argume...
std::string description_no_recombiner() const
returns a description not including the recombiner information
static const double max_allowable_R
R values larger than max_allowable_R are not allowed.
void set_recombiner(const Recombiner *recomb)
set the recombiner class to the one provided
void set_recombination_scheme(RecombinationScheme)
set the recombination scheme to the one provided
void delete_plugin_when_unused()
calling this causes the JetDefinition to handle the deletion of the plugin when it is no longer used
Class to contain pseudojets, including minimal information of use to jet-clustering routines.
double rap() const
returns the rapidity or some large value when the rapidity is infinite
double modp() const
return the 3-vector modulus = sqrt(px^2+py^2+pz^2)
void reset_momentum(double px, double py, double pz, double E)
reset the 4-momentum according to the supplied components but leave all other information (indices,...
double phi() const
returns phi (in the range 0..2pi)
double perp() const
returns the scalar transverse momentum
double perp2() const
returns the squared transverse momentum
double m2() const
returns the squared invariant mass // like CLHEP
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure_in)
set the associated structure
void reset_PtYPhiM(double pt_in, double y_in, double phi_in, double m_in=0.0)
reset the PseudoJet according to the specified pt, rapidity, azimuth and mass (also resetting indices...
double pt() const
returns the scalar transverse momentum
double m() const
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP)
double modp2() const
return the squared 3-vector modulus = px^2+py^2+pz^2
void reset(double px, double py, double pz, double E)
reset the 4-momentum according to the supplied components and put the user and history indices back t...
double pt2() const
returns the squared transverse momentum
An implementation of shared pointers that is broadly similar to C++11 shared_ptr (https://en....
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
@ plugin_strategy
the plugin has been used...
RecombinationScheme
The various recombination schemes.
@ E_scheme
summing the 4-momenta
@ BIpt2_scheme
pt^2 weighted recombination of y,phi (and summing of pt's) no preprocessing
@ pt2_scheme
pt^2 weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless...
@ pt_scheme
pt weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless b...
@ external_scheme
for the user's external scheme
@ Et_scheme
pt weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless b...
@ WTA_pt_scheme
pt-based Winner-Takes-All (WTA) recombination: the result of the recombination has the rapidity,...
@ WTA_modp_scheme
mod-p-based Winner-Takes-All (WTA) recombination: the result of the recombination gets the 3-vector d...
@ BIpt_scheme
pt weighted recombination of y,phi (and summing of pt's), with no preprocessing
@ Et2_scheme
pt^2 weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless...
JetAlgorithm
the various families of jet-clustering algorithm
@ ee_genkt_algorithm
the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt)
@ genkt_algorithm
like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti...
@ undefined_jet_algorithm
the value for the jet algorithm in a JetDefinition for which no algorithm has yet been defined
@ plugin_algorithm
any plugin algorithm supplied by the user
@ ee_kt_algorithm
the e+e- kt algorithm
@ cambridge_algorithm
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm).
@ cambridge_for_passive_algorithm
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
@ antikt_algorithm
like the k_t but with distance measures dij = min(1/kti^2,1/ktj^2) Delta R_{ij}^2 / R^2 diB = 1/kti^2
@ kt_algorithm
the longitudinally invariant kt algorithm