31 #include "fastjet/JetDefinition.hh" 
   32 #include "fastjet/Error.hh" 
   33 #include "fastjet/CompositeJetStructure.hh" 
   36 FASTJET_BEGIN_NAMESPACE      
 
   40 const 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.reset(other_jet_def._shared_recombiner);
 
  224   if (other_jd.recombination_scheme() != scheme) 
return false;
 
  235   if (_recombiner == 0){
 
  236     throw Error(
"tried to call JetDefinition::delete_recombiner_when_unused() for a JetDefinition without a user-defined recombination scheme");
 
  237   } 
else if (_shared_recombiner.get()) {
 
  238     throw Error(
"Error in JetDefinition::delete_recombiner_when_unused: the recombiner is already scheduled for deletion when unused (or was already set as shared)");
 
  241   _shared_recombiner.reset(_recombiner);
 
  248     throw Error(
"tried to call JetDefinition::delete_plugin_when_unused() for a JetDefinition without a plugin");
 
  251   _plugin_shared.reset(_plugin);
 
  257   switch(_recomb_scheme) {
 
  259     return "E scheme recombination";
 
  261     return "pt scheme recombination";
 
  263     return "pt2 scheme recombination";
 
  265     return "Et scheme recombination";
 
  267     return "Et2 scheme recombination";
 
  269     return "boost-invariant pt scheme recombination";
 
  271     return "boost-invariant pt2 scheme recombination";
 
  273     return "pt-ordered Winner-Takes-All recombination";
 
  280     return "|3-momentum|-ordered Winner-Takes-All recombination";
 
  283     err << 
"DefaultRecombiner: unrecognized recombination scheme "  
  285     throw Error(err.str());
 
  294   double weighta, weightb;
 
  296   switch(_recomb_scheme) {
 
  301     pab.
reset(pa.px()+pb.px(),
 
  317     weighta = pa.
perp2(); 
 
  318     weightb = pb.
perp2();
 
  324                       phard.
rap(), phard.
phi(), phard.
m());
 
  350     const PseudoJet & phard = a_hardest ? pa : pb;
 
  351     const PseudoJet & psoft = a_hardest ? pb : pa;
 
  356     double modp_hard = phard.
modp();
 
  357     double modp_ab = modp_hard + psoft.
modp();
 
  358     if (phard.
modp2()==0.0){
 
  359       pab.
reset(0.0, 0.0, 0.0, phard.
m());
 
  361       double scale = modp_ab/modp_hard;
 
  362       pab.
reset(phard.px()*scale, phard.py()*scale, phard.pz()*scale,
 
  363                 sqrt(modp_ab*modp_ab + phard.
m2()));
 
  368     err << 
"DefaultRecombiner: unrecognized recombination scheme "  
  370     throw Error(err.str());
 
  373   double perp_ab = pa.
perp() + pb.
perp();
 
  374   if (perp_ab != 0.0) { 
 
  375     double y_ab    = (weighta * pa.
rap() + weightb * pb.
rap())/(weighta+weightb);
 
  378     double phi_a = pa.
phi(), phi_b = pb.
phi();
 
  379     if (phi_a - phi_b > pi)  phi_b += twopi;
 
  380     if (phi_a - phi_b < -pi) phi_b -= twopi;
 
  381     double phi_ab = (weighta * phi_a + weightb * phi_b)/(weighta+weightb);
 
  391     pab.
reset(0.0, 0.0, 0.0, 0.0);
 
  397   switch(_recomb_scheme) {
 
  410       double newE = sqrt(p.
perp2()+p.pz()*p.pz());
 
  423       double rescale = p.E()/sqrt(p.
perp2()+p.pz()*p.pz());
 
  424       p.
reset_momentum(rescale*p.px(), rescale*p.py(), rescale*p.pz(), p.E());
 
  433     err << 
"DefaultRecombiner: unrecognized recombination scheme "  
  435     throw Error(err.str());
 
  440   throw Error(
"set_ghost_separation_scale not supported");
 
  459   if (pieces.size()>0){
 
  461     for (
unsigned int i=1; i<pieces.size(); i++)
 
  477   return join(vector<PseudoJet>(1,j1), recombiner);
 
  483   vector<PseudoJet> pieces;
 
  484   pieces.push_back(j1);
 
  485   pieces.push_back(j2);
 
  486   return join(pieces, recombiner);
 
  492   vector<PseudoJet> pieces;
 
  493   pieces.push_back(j1);
 
  494   pieces.push_back(j2);
 
  495   pieces.push_back(j3);
 
  496   return join(pieces, recombiner);
 
  502   vector<PseudoJet> pieces;
 
  503   pieces.push_back(j1);
 
  504   pieces.push_back(j2);
 
  505   pieces.push_back(j3);
 
  506   pieces.push_back(j4);
 
  507   return join(pieces, recombiner);
 
  513 FASTJET_END_NAMESPACE
 
pt weighted recombination of y,phi (and summing of pt's), with no preprocessing 
 
void set_recombination_scheme(RecombinationScheme)
set the recombination scheme to the one provided 
 
double rap() const 
returns the rapidity or some large value when the rapidity is infinite 
 
a version of cambridge with a special distance measure for particles whose pt is < extra_param(); thi...
 
like the k_t but with distance measures dij = min(kti^{2p},ktj^{2p}) Delta R_{ij}^2 / R^2 diB = 1/kti...
 
void delete_recombiner_when_unused()
calling this tells the JetDefinition to handle the deletion of the recombiner when it is no longer us...
 
the value for the jet algorithm in a JetDefinition for which no algorithm has yet been defined ...
 
pt^2 weighted recombination of y,phi (and summing of pt's) no preprocessing 
 
JetAlgorithm jet_algorithm() const 
return information about the definition... 
 
void delete_plugin_when_unused()
calling this causes the JetDefinition to handle the deletion of the plugin when it is no longer used ...
 
pt^2 weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless...
 
the e+e- genkt algorithm (R > 2 and p=1 gives ee_kt) 
 
double modp() const 
return the 3-vector modulus = sqrt(px^2+py^2+pz^2) 
 
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...
 
the longitudinally invariant kt algorithm 
 
static std::string algorithm_description(const JetAlgorithm jet_alg)
a short textual description of the algorithm jet_alg 
 
pt-based Winner-Takes-All (WTA) recombination: the result of the recombination has the rapidity...
 
The structure for a jet made of pieces. 
 
const Plugin * plugin() const 
return a pointer to the plugin 
 
void set_extra_param(double xtra_param)
(re)set the general purpose extra parameter 
 
double m2() const 
returns the squared invariant mass // like CLHEP 
 
any plugin algorithm supplied by the user 
 
mod-p-based Winner-Takes-All (WTA) recombination: the result of the recombination gets the 3-vector d...
 
virtual bool is_spherical() const 
returns true if the plugin implements an algorithm intended for use on a spherical geometry (e...
 
the plugin has been used... 
 
A class that will provide the recombination scheme facilities and/or allow a user to extend these fac...
 
virtual std::string description() const 
return a textual description of the recombination scheme implemented here 
 
std::string description() const 
return a textual description of the current jet definition 
 
pt^2 weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless...
 
virtual void preprocess(PseudoJet &p) const 
routine called to preprocess each input jet (to make all input jets compatible with the scheme requir...
 
virtual std::string description() const =0
return a textual description of the jet-definition implemented in this plugin 
 
base class corresponding to errors that can be thrown by FastJet 
 
std::string description_no_recombiner() const 
returns a description not including the recombiner information 
 
RecombinationScheme
The various recombination schemes. 
 
virtual void set_ghost_separation_scale(double scale) const 
set the ghost separation scale for passive area determinations in future runs (strictly speaking that...
 
void plus_equal(PseudoJet &pa, const PseudoJet &pb) const 
pa += pb in the given recombination scheme. 
 
static const double max_allowable_R
R values larger than max_allowable_R are not allowed. 
 
an implementation of C++0x shared pointers (or boost's) 
 
pt weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless b...
 
virtual std::string description() const =0
return a textual description of the recombination scheme implemented here 
 
pt weighted recombination of y,phi (and summing of pt's) with preprocessing to make things massless b...
 
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 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 phi() const 
returns phi (in the range 0..2pi) 
 
virtual void recombine(const PseudoJet &pa, const PseudoJet &pb, PseudoJet &pab) const 
recombine pa and pb and put result into pab 
 
void set_recombiner(const Recombiner *recomb)
set the recombiner class to the one provided 
 
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...
 
bool is_spherical() const 
returns true if the jet definition involves an algorithm intended for use on a spherical geometry (e...
 
void set_structure_shared_ptr(const SharedPtr< PseudoJetStructureBase > &structure)
set the associated structure 
 
double perp() const 
returns the scalar transverse momentum 
 
the longitudinally invariant variant of the cambridge algorithm (aka Aachen algoithm). 
 
double modp2() const 
return the squared 3-vector modulus = px^2+py^2+pz^2 
 
double m() const 
returns the invariant mass (If m2() is negative then -sqrt(-m2()) is returned, as in CLHEP) ...
 
An abstract base class that will provide the recombination scheme facilities and/or allow a user to e...
 
JetAlgorithm
the various families of jet-clustering algorithm 
 
double pt() const 
returns the scalar transverse momentum 
 
Strategy
the various options for the algorithmic strategy to adopt in clustering events with kt and cambridge ...
 
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
 
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 pt2() const 
returns the squared transverse momentum 
 
for the user's external scheme 
 
double perp2() const 
returns the squared transverse momentum 
 
class that is intended to hold a full definition of the jet clusterer 
 
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...