fastjet 2.4.5
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
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
[legend]
Collaboration diagram for fastjet::ClusterSequenceActiveAreaExplicitGhosts:
Collaboration graph
[legend]

List of all members.

Public Member Functions

template<class L >
 ClusterSequenceActiveAreaExplicitGhosts (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const GhostedAreaSpec &ghost_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 *ghost_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 bool has_explicit_ghosts () const
 this class does have explicit ghosts
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 &ghost_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)
 another way of adding 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
 get the area of the 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
 handle warning messages

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


Constructor & Destructor Documentation

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

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

Definition at line 55 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

           : ClusterSequenceAreaBase() {
           std::vector<L> * ghosts = NULL;
           _initialise(pseudojets,jet_def,&ghost_spec,ghosts,0.0,
                       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 65 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

           : ClusterSequenceAreaBase() {
           const GhostedAreaSpec * ghost_spec = NULL;
           _initialise(pseudojets,jet_def,ghost_spec,&ghosts,ghost_area,
                       writeout_combinations); }

Member Function Documentation

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

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

Definition at line 46 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

References fastjet::GhostedAreaSpec::actual_ghost_area(), fastjet::GhostedAreaSpec::add_ghosts(), and fastjet::GhostedAreaSpec::n_ghosts().

                                                             {

  // add the ghosts to the jets
  ghost_spec.add_ghosts(_jets);

  // now add labelling...
  for (unsigned i = _initial_hard_n; i < _jets.size(); i++) {
    //_jets[i].set_user_index(1);
    _is_pure_ghost.push_back(true);
  }

  // and record some info from the ghost_spec
  _ghost_area = ghost_spec.actual_ghost_area();
  _n_ghosts   = ghost_spec.n_ghosts();
}
template<class L >
void fastjet::ClusterSequenceActiveAreaExplicitGhosts::_add_ghosts ( const std::vector< L > &  ghosts,
double  ghost_area 
) [private]

another way of adding ghosts

add an explicitly specified bunch of ghosts

Definition at line 222 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

                                             {

  
  for (unsigned i = 0; i < ghosts.size(); i++) {
    _is_pure_ghost.push_back(true);
    _jets.push_back(ghosts[i]);
  }
  // and record some info about ghosts
  _ghost_area = ghost_area;
  _n_ghosts   = ghosts.size();
}
template<class L >
void fastjet::ClusterSequenceActiveAreaExplicitGhosts::_initialise ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const GhostedAreaSpec ghost_spec,
const std::vector< L > *  ghosts,
double  ghost_area,
const bool &  writeout_combinations 
)

does the actual work of initialisation

Definition at line 169 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

                                              {
  // don't reserve space yet -- will be done below

  // insert initial jets this way so that any type L that can be
  // converted to a pseudojet will work fine (basically PseudoJet
  // and any type that has [] subscript access to the momentum
  // components, such as CLHEP HepLorentzVector).
  for (unsigned int i = 0; i < pseudojets.size(); i++) {
    PseudoJet mom(pseudojets[i]);
    //mom.set_user_index(0); // for user's particles (user index now lost...)
    _jets.push_back(mom);
    _is_pure_ghost.push_back(false);
  }

  _initial_hard_n = _jets.size();

  if (ghost_spec != NULL) {
    _add_ghosts(*ghost_spec);
  } else {
    _add_ghosts(*ghosts, ghost_area);
  }

  if (writeout_combinations) {
    std::cout << "# Printing particles including ghosts\n";
    for (unsigned j = 0; j < _jets.size(); j++) {
      printf("%5u %20.13f %20.13f %20.13e\n",
               j,_jets[j].rap(),_jets[j].phi_02pi(),_jets[j].kt2());
    }
    std::cout << "# Finished printing particles including ghosts\n";
  }

  // this will ensure that we can still point to jets without
  // difficulties arising!
  _jets.reserve(_jets.size()*2);

  // run the clustering
  _initialise_and_run(jet_def,writeout_combinations);

  // set up all other information
  _post_process();
}
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.

                                      {

  // first check for danger signals.
  // Establish largest ghost transverse momentum
  _max_ghost_perp2 = 0.0;
  for (int i = 0; i < _initial_n; i++) {
    if (_is_pure_ghost[i] && _jets[i].perp2() > _max_ghost_perp2) 
      _max_ghost_perp2 = _jets[i].perp2();
  }

  // now find out if any of the particles are close to danger
  double danger_ratio = numeric_limits<double>::epsilon();
  danger_ratio = danger_ratio * danger_ratio;
  _has_dangerous_particles = false;
  for (int i = 0; i < _initial_n; i++) {
    if (!_is_pure_ghost[i] && 
        danger_ratio * _jets[i].perp2() <=  _max_ghost_perp2) {
      _has_dangerous_particles = true;
      break;
    }
  }

  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");

  // sort out sizes
  _areas.resize(_history.size());
  _area_4vectors.resize(_history.size());
  _is_pure_ghost.resize(_history.size());
  
  // First set up areas for the initial particles (ghost=_ghost_area,
  // real particles = 0); recall that _initial_n here is the number of
  // particles including ghosts
  for (int i = 0; i < _initial_n; i++) {
    if (_is_pure_ghost[i]) {
      _areas[i] = _ghost_area;
      // normalise pt to be _ghost_area (NB we make use of fact that
      // for initial particles, jet and clust_hist index are the same).
      _area_4vectors[i] = (_ghost_area/_jets[i].perp()) * _jets[i];
    } else {
      _areas[i] = 0;
      _area_4vectors[i] = PseudoJet(0.0,0.0,0.0,0.0);
    }
  }
  
  // next follow the branching through and set up the areas 
  // and ghost-nature at each step of the clustering (rather than
  // each jet).
  for (unsigned i = _initial_n; i < _history.size(); i++) {
    if (_history[i].parent2 == BeamJet) {
      _is_pure_ghost[i]  = _is_pure_ghost[_history[i].parent1];
      _areas[i]          = _areas[_history[i].parent1];
      _area_4vectors[i] = _area_4vectors[_history[i].parent1];
    } else {
      _is_pure_ghost[i]  = _is_pure_ghost[_history[i].parent1] && 
                           _is_pure_ghost[_history[i].parent2]   ;
      _areas[i]          = _areas[_history[i].parent1] + 
                           _areas[_history[i].parent2]  ;                          
      _jet_def.recombiner()->recombine(_area_4vectors[_history[i].parent1], 
                           _area_4vectors[_history[i].parent2],
                           _area_4vectors[i]);  
//      _area_4vectors[i] = _area_4vectors[_history[i].parent1] + 
//                          _area_4vectors[_history[i].parent2]  ;
    }

  }
  
}
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 fastjet::PseudoJet::cluster_hist_index().

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

                                                           {
  return _areas[jet.cluster_hist_index()];
}
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 fastjet::PseudoJet::cluster_hist_index().

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

                                                                      {
  return _area_4vectors[jet.cluster_hist_index()];
}
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 fastjet::RangeDefinition::is_in_range().

                                                                        {
  vector<PseudoJet> unclust = unclustered_particles();
  double area = 0.0;
  for (unsigned iu = 0; iu < unclust.size();  iu++) {
    if (is_pure_ghost(unclust[iu]) && range.is_in_range(unclust[iu])) {
      area += _ghost_area;
    }
  }
  return area;
}
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 124 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Referenced by fastjet::ClusterSequenceActiveArea::_run_AA(), and fastjet::ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E().

virtual bool fastjet::ClusterSequenceActiveAreaExplicitGhosts::has_explicit_ghosts ( ) const [inline, virtual]

this class does have explicit ghosts

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 108 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

{return true;}
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.

{
  return hist_ix >= 0 ? _is_pure_ghost[hist_ix] : false;
}
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 fastjet::PseudoJet::cluster_hist_index().

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

{
  return _is_pure_ghost[jet.cluster_hist_index()];
}
double fastjet::ClusterSequenceActiveAreaExplicitGhosts::max_ghost_perp2 ( ) const [inline]

returns the largest squared transverse momentum among all ghosts

Definition at line 119 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

{return _max_ghost_perp2;}
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 217 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

{return _initial_hard_n;}
double fastjet::ClustSeqActAreaEG::total_area ( ) const

returns the total area under study

Definition at line 73 of file ClusterSequenceActiveAreaExplicitGhosts.cc.

                                            {
  return _n_ghosts * _ghost_area;
}

Member Data Documentation

Definition at line 135 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Definition at line 134 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Definition at line 132 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Definition at line 139 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Definition at line 146 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Definition at line 133 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

Definition at line 138 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

get the area of the ghosts

Definition at line 131 of file ClusterSequenceActiveAreaExplicitGhosts.hh.

handle warning messages

allow for warnings

Reimplemented from fastjet::ClusterSequenceAreaBase.

Definition at line 140 of file ClusterSequenceActiveAreaExplicitGhosts.hh.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines