fastjet::ClusterSequenceActiveAreaExplicitGhosts Class Reference

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>

Inheritance diagram for fastjet::ClusterSequenceActiveAreaExplicitGhosts:

Inheritance graph
fastjet::ClusterSequenceAreaBasefastjet::ClusterSequence
[legend]
Collaboration diagram for fastjet::ClusterSequenceActiveAreaExplicitGhosts:

Collaboration graph
fastjet::ClusterSequenceAreaBasefastjet::ClusterSequencefastjet::JetDefinitionfastjet::JetDefinition::DefaultRecombinerfastjet::JetDefinition::Recombinerfastjet::JetDefinition::PluginLimitedWarning
[legend]
List of all members.

Public Member Functions

template<class L>
 ClusterSequenceActiveAreaExplicitGhosts (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const GhostedAreaSpec &area_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 *area_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 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 &area_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)
 add an explicitly specified bunch of 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
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
 allow for warnings

Detailed Description

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 48 of file ClusterSequenceActiveAreaExplicitGhosts.hh.


Constructor & Destructor Documentation

template<class L>
fastjet::ClusterSequenceActiveAreaExplicitGhosts::ClusterSequenceActiveAreaExplicitGhosts ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const GhostedAreaSpec area_spec,
const bool &  writeout_combinations = false 
) [inline]

constructor using a GhostedAreaSpec to specify how the area is to be measured

Definition at line 54 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

00058            : ClusterSequenceAreaBase() {
00059            std::vector<L> * ghosts = NULL;
00060            _initialise(pseudojets,jet_def,&area_spec,ghosts,0.0,
00061                        writeout_combinations); }

template<class L>
fastjet::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 64 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

00069            : ClusterSequenceAreaBase() {
00070            const GhostedAreaSpec * area_spec = NULL;
00071            _initialise(pseudojets,jet_def,area_spec,&ghosts,ghost_area,
00072                        writeout_combinations); }


Member Function Documentation

template<class L>
void fastjet::ClusterSequenceActiveAreaExplicitGhosts::_initialise ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const GhostedAreaSpec area_spec,
const std::vector< L > *  ghosts,
double  ghost_area,
const bool &  writeout_combinations 
) [inline]

does the actual work of initialisation

Definition at line 162 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

00167                                               {
00168   // don't reserve space yet -- will be done below
00169 
00170   // insert initial jets this way so that any type L that can be
00171   // converted to a pseudojet will work fine (basically PseudoJet
00172   // and any type that has [] subscript access to the momentum
00173   // components, such as CLHEP HepLorentzVector).
00174   for (unsigned int i = 0; i < pseudojets.size(); i++) {
00175     PseudoJet mom(pseudojets[i]);
00176     //mom.set_user_index(0); // for user's particles (user index now lost...)
00177     _jets.push_back(mom);
00178     _is_pure_ghost.push_back(false);
00179   }
00180 
00181   _initial_hard_n = _jets.size();
00182 
00183   if (area_spec != NULL) {
00184     _add_ghosts(*area_spec);
00185   } else {
00186     _add_ghosts(*ghosts, ghost_area);
00187   }
00188 
00189   if (writeout_combinations) {
00190     std::cout << "# Printing particles including ghosts\n";
00191     for (unsigned j = 0; j < _jets.size(); j++) {
00192       printf("%5u %20.13f %20.13f %20.13e\n",
00193                j,_jets[j].rap(),_jets[j].phi_02pi(),_jets[j].kt2());
00194     }
00195     std::cout << "# Finished printing particles including ghosts\n";
00196   }
00197 
00198   // this will ensure that we can still point to jets without
00199   // difficulties arising!
00200   _jets.reserve(_jets.size()*2);
00201 
00202   // run the clustering
00203   _initialise_and_run(jet_def,writeout_combinations);
00204 
00205   // set up all other information
00206   _post_process();
00207 }

unsigned int fastjet::ClusterSequenceActiveAreaExplicitGhosts::n_hard_particles (  )  const [inline]

returns the number of hard particles (i.e. those supplied by the user).

Definition at line 210 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

References _initial_hard_n.

00210 {return _initial_hard_n;}

double fastjet::ClustSeqActAreaEG::area ( const PseudoJet jet  )  const [virtual]

returns the area of a jet

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 66 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _areas, and fastjet::PseudoJet::cluster_hist_index().

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), and empty_area().

00066                                                            {
00067   return _areas[jet.cluster_hist_index()];
00068 }

PseudoJet fastjet::ClustSeqActAreaEG::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 fastjet::ClusterSequenceAreaBase.

Definition at line 80 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _area_4vectors, and fastjet::PseudoJet::cluster_hist_index().

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas().

00080                                                                       {
00081   return _area_4vectors[jet.cluster_hist_index()];
00082 }

bool fastjet::ClustSeqActAreaEG::is_pure_ghost ( const PseudoJet jet  )  const [virtual]

true if a jet is made exclusively of ghosts

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 85 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _is_pure_ghost, and fastjet::PseudoJet::cluster_hist_index().

Referenced by fastjet::ClusterSequenceActiveArea::_transfer_areas(), fastjet::ClusterSequenceActiveArea::_transfer_ghost_free_history(), and empty_area().

00086 {
00087   return _is_pure_ghost[jet.cluster_hist_index()];
00088 }

bool fastjet::ClustSeqActAreaEG::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.

Definition at line 91 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _is_pure_ghost.

00092 {
00093   return hist_ix >= 0 ? _is_pure_ghost[hist_ix] : false;
00094 }

double fastjet::ClustSeqActAreaEG::empty_area ( const RangeDefinition range  )  const [virtual]

return the total area, up to |y|<maxrap, that consists of unclustered ghosts

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 97 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _ghost_area, area(), fastjet::RangeDefinition::is_in_range(), is_pure_ghost(), and fastjet::ClusterSequence::unclustered_particles().

00097                                                                         {
00098   vector<PseudoJet> unclust = unclustered_particles();
00099   double area = 0.0;
00100   for (unsigned iu = 0; iu < unclust.size();  iu++) {
00101     if (is_pure_ghost(unclust[iu]) && range.is_in_range(unclust[iu])) {
00102       area += _ghost_area;
00103     }
00104   }
00105   return area;
00106 }

double fastjet::ClustSeqActAreaEG::total_area (  )  const

returns the total area under study

Definition at line 73 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _ghost_area, and _n_ghosts.

00073                                             {
00074   return _n_ghosts * _ghost_area;
00075 }

double fastjet::ClusterSequenceActiveAreaExplicitGhosts::max_ghost_perp2 (  )  const [inline]

returns the largest squared transverse momentum among all ghosts

Definition at line 115 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

00115 {return _max_ghost_perp2;}

bool fastjet::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 120 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E().

00120 {return _has_dangerous_particles;}

void fastjet::ClustSeqActAreaEG::_add_ghosts ( const GhostedAreaSpec area_spec  )  [private]

adds the "ghost" momenta, which will be used to estimate the jet area

Definition at line 46 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _ghost_area, _initial_hard_n, _is_pure_ghost, fastjet::ClusterSequence::_jets, _n_ghosts, fastjet::GhostedAreaSpec::actual_ghost_area(), fastjet::GhostedAreaSpec::add_ghosts(), and fastjet::GhostedAreaSpec::n_ghosts().

00047                                                             {
00048 
00049   // add the ghosts to the jets
00050   area_spec.add_ghosts(_jets);
00051 
00052   // now add labelling...
00053   for (unsigned i = _initial_hard_n; i < _jets.size(); i++) {
00054     //_jets[i].set_user_index(1);
00055     _is_pure_ghost.push_back(true);
00056   }
00057 
00058   // and record some info from the area_spec
00059   _ghost_area = area_spec.actual_ghost_area();
00060   _n_ghosts   = area_spec.n_ghosts();
00061 }

template<class L>
void fastjet::ClusterSequenceActiveAreaExplicitGhosts::_add_ghosts ( const std::vector< L > &  ghosts,
double  ghost_area 
) [inline, private]

add an explicitly specified bunch of ghosts

Definition at line 215 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

References _ghost_area, _is_pure_ghost, fastjet::ClusterSequence::_jets, and _n_ghosts.

00217                                              {
00218 
00219   
00220   for (unsigned i = 0; i < ghosts.size(); i++) {
00221     _is_pure_ghost.push_back(true);
00222     _jets.push_back(ghosts[i]);
00223   }
00224   // and record some info about ghosts
00225   _ghost_area = ghost_area;
00226   _n_ghosts   = ghosts.size();
00227 }

void fastjet::ClustSeqActAreaEG::_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.

)

Definition at line 110 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References _area_4vectors, _areas, _ghost_area, _has_dangerous_particles, fastjet::ClusterSequence::_history, fastjet::ClusterSequence::_initial_n, _is_pure_ghost, fastjet::ClusterSequence::_jet_def, fastjet::ClusterSequence::_jets, _max_ghost_perp2, _warnings, fastjet::ClusterSequence::BeamJet, fastjet::JetDefinition::recombiner(), and LimitedWarning::warn().

00110                                       {
00111 
00112   // first check for danger signals.
00113   // Establish largest ghost transverse momentum
00114   _max_ghost_perp2 = 0.0;
00115   for (int i = 0; i < _initial_n; i++) {
00116     if (_is_pure_ghost[i] && _jets[i].perp2() > _max_ghost_perp2) 
00117       _max_ghost_perp2 = _jets[i].perp2();
00118   }
00119 
00120   // now find out if any of the particles are close to danger
00121   double danger_ratio = numeric_limits<double>::epsilon();
00122   danger_ratio = danger_ratio * danger_ratio;
00123   _has_dangerous_particles = false;
00124   for (int i = 0; i < _initial_n; i++) {
00125     if (!_is_pure_ghost[i] && 
00126         danger_ratio * _jets[i].perp2() <=  _max_ghost_perp2) {
00127       _has_dangerous_particles = true;
00128       break;
00129     }
00130   }
00131 
00132   if (_has_dangerous_particles) _warnings.warn("ClusterSequenceActiveAreaExplicitGhosts: \n  ghosts not sufficiently soft wrt some of the input particles\n  a common cause is (unphysical?) input particles with pt=0 but finite rapidity");
00133 
00134   // sort out sizes
00135   _areas.resize(_history.size());
00136   _area_4vectors.resize(_history.size());
00137   _is_pure_ghost.resize(_history.size());
00138   
00139   // First set up areas for the initial particles (ghost=_ghost_area,
00140   // real particles = 0); recall that _initial_n here is the number of
00141   // particles including ghosts
00142   for (int i = 0; i < _initial_n; i++) {
00143     if (_is_pure_ghost[i]) {
00144       _areas[i] = _ghost_area;
00145       // normalise pt to be _ghost_area (NB we make use of fact that
00146       // for initial particles, jet and clust_hist index are the same).
00147       _area_4vectors[i] = (_ghost_area/_jets[i].perp()) * _jets[i];
00148     } else {
00149       _areas[i] = 0;
00150       _area_4vectors[i] = PseudoJet(0.0,0.0,0.0,0.0);
00151     }
00152   }
00153   
00154   // next follow the branching through and set up the areas 
00155   // and ghost-nature at each step of the clustering (rather than
00156   // each jet).
00157   for (unsigned i = _initial_n; i < _history.size(); i++) {
00158     if (_history[i].parent2 == BeamJet) {
00159       _is_pure_ghost[i]  = _is_pure_ghost[_history[i].parent1];
00160       _areas[i]          = _areas[_history[i].parent1];
00161       _area_4vectors[i] = _area_4vectors[_history[i].parent1];
00162     } else {
00163       _is_pure_ghost[i]  = _is_pure_ghost[_history[i].parent1] && 
00164                            _is_pure_ghost[_history[i].parent2]   ;
00165       _areas[i]          = _areas[_history[i].parent1] + 
00166                            _areas[_history[i].parent2]  ;                          
00167       _jet_def.recombiner()->recombine(_area_4vectors[_history[i].parent1], 
00168                            _area_4vectors[_history[i].parent2],
00169                            _area_4vectors[i]);  
00170 //      _area_4vectors[i] = _area_4vectors[_history[i].parent1] + 
00171 //                          _area_4vectors[_history[i].parent2]  ;
00172     }
00173 
00174   }
00175   
00176 }


Member Data Documentation

int fastjet::ClusterSequenceActiveAreaExplicitGhosts::_n_ghosts [private]

Definition at line 124 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _add_ghosts(), and total_area().

double fastjet::ClusterSequenceActiveAreaExplicitGhosts::_ghost_area [private]

Definition at line 125 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _add_ghosts(), _post_process(), empty_area(), and total_area().

std::vector<bool> fastjet::ClusterSequenceActiveAreaExplicitGhosts::_is_pure_ghost [private]

Definition at line 126 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _add_ghosts(), _post_process(), and is_pure_ghost().

std::vector<double> fastjet::ClusterSequenceActiveAreaExplicitGhosts::_areas [private]

Definition at line 127 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _post_process(), and area().

std::vector<PseudoJet> fastjet::ClusterSequenceActiveAreaExplicitGhosts::_area_4vectors [private]

Definition at line 128 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _post_process(), and area_4vector().

double fastjet::ClusterSequenceActiveAreaExplicitGhosts::_max_ghost_perp2 [private]

Definition at line 131 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _post_process().

bool fastjet::ClusterSequenceActiveAreaExplicitGhosts::_has_dangerous_particles [private]

Definition at line 132 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _post_process().

LimitedWarning fastjet::ClustSeqActAreaEG::_warnings [static, private]

allow for warnings

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 133 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _post_process().

unsigned int fastjet::ClusterSequenceActiveAreaExplicitGhosts::_initial_hard_n [private]

Definition at line 139 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by _add_ghosts(), and n_hard_particles().


The documentation for this class was generated from the following files:
Generated on Thu Jan 3 19:05:16 2008 for fastjet by  doxygen 1.5.2