fastjet 2.4.3
Public Member Functions | Private Attributes

TrackJetPlugin Class Reference

#include <TrackJetPlugin.hh>

Inheritance diagram for TrackJetPlugin:
Inheritance graph
[legend]
Collaboration diagram for TrackJetPlugin:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 TrackJetPlugin (double radius, RecombinationScheme jet_recombination_scheme=pt_scheme, RecombinationScheme track_recombination_scheme=pt_scheme)
 Main constructor for the TrackJet Plugin class.
 TrackJetPlugin (const TrackJetPlugin &plugin)
 copy constructor
virtual std::string description () const
 return a textual description of the jet-definition implemented in this plugin
virtual void run_clustering (ClusterSequence &) const
 given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:

  • plugin_do_ij_recombination(...)
  • plugin_do_iB_recombination(...)

virtual double R () const
 the plugin mechanism's standard way of accessing the jet radius here we return the R of the last alg in the list

Private Attributes

double _radius
double _radius2
JetDefinition::DefaultRecombiner _jet_recombiner
JetDefinition::DefaultRecombiner _track_recombiner

Detailed Description

Definition at line 44 of file TrackJetPlugin.hh.


Constructor & Destructor Documentation

TrackJetPlugin::TrackJetPlugin ( double  radius,
RecombinationScheme  jet_recombination_scheme = pt_scheme,
RecombinationScheme  track_recombination_scheme = pt_scheme 
) [inline]

Main constructor for the TrackJet Plugin class.

The argument is an initialised list of jet algorithms

Parameters:
_radiusthe distance at which point a particle is no longer recombied into the jet
jet_recombination_schemethe recombination scheme used to sum the 4-vecors inside the jet
track_recombination_schemethe recombination scheme used to sum the 4-vecors when accumulating track into a the jet Both recombiners are defaulted to pt_scheme recomb as for the Rivet implementation.

Definition at line 58 of file TrackJetPlugin.hh.

References _jet_recombiner, _radius, _radius2, and _track_recombiner.

                                                                           {
    _radius  = radius;
    _radius2 = radius*radius;
    _jet_recombiner = JetDefinition::DefaultRecombiner(jet_recombination_scheme);
    _track_recombiner = JetDefinition::DefaultRecombiner(track_recombination_scheme);
  }
TrackJetPlugin::TrackJetPlugin ( const TrackJetPlugin plugin) [inline]

copy constructor

Definition at line 68 of file TrackJetPlugin.hh.

                                                 {
    *this = plugin;
  }

Member Function Documentation

string TrackJetPlugin::description ( ) const [virtual]

return a textual description of the jet-definition implemented in this plugin

Implements JetDefinition::Plugin.

Definition at line 84 of file GhostedAreaSpec.cc.

                                          {

  ostringstream ostr;
  ostr << "ghosts of area " << actual_ghost_area() 
       << " (had requested " << ghost_area() << ")"
       << ", placed up to y = " << ghost_maxrap() 
       << ", scattered wrt to perfect grid by (rel) " << grid_scatter() 
       << ", mean_ghost_kt = " << mean_ghost_kt()
       << ", rel kt_scatter =  " << kt_scatter()
       << ", n repetitions of ghost distributions =  " << repeat();
  return ostr.str();
}
virtual double TrackJetPlugin::R ( ) const [inline, virtual]

the plugin mechanism's standard way of accessing the jet radius here we return the R of the last alg in the list

Implements JetDefinition::Plugin.

Definition at line 78 of file TrackJetPlugin.hh.

References _radius.

{return _radius;}
void TrackJetPlugin::run_clustering ( ClusterSequence ) const [virtual]

given a ClusterSequence that has been filled up with initial particles, the following function should fill up the rest of the ClusterSequence, using the following member functions of ClusterSequence:

  • plugin_do_ij_recombination(...)
  • plugin_do_iB_recombination(...)

Implements JetDefinition::Plugin.

Definition at line 65 of file TrackJetPlugin.cc.

References ClusterSequence::jets(), ClusterSequence::plugin_record_iB_recombination(), and ClusterSequence::plugin_record_ij_recombination().

                                                                     {
  // we first need to sort the particles in pt
  vector<TrackJetParticlePtr> particle_list;

  const vector<PseudoJet> & jets = clust_seq.jets();  
  int index=0;
  for (vector<PseudoJet>::const_iterator mom_it = jets.begin(); mom_it != jets.end(); mom_it++){
    particle_list.push_back(TrackJetParticlePtr(index, mom_it->perp2()));
    index++;
  }

  // sort the particles into decreasing pt
  sort(particle_list.begin(), particle_list.end());


  // if we're using a recombination scheme different from the E scheme,
  // we first need to update the particles' energy so that they
  // are massless (and rapidity = pseudorapidity)
  vector<PseudoJet> tuned_particles = clust_seq.jets();
  vector<PseudoJet> tuned_tracks = clust_seq.jets();
  for (vector<PseudoJet>::iterator pit = tuned_particles.begin();
       pit != tuned_particles.end(); pit++)
    _jet_recombiner.preprocess(*pit);
  for (vector<PseudoJet>::iterator pit = tuned_tracks.begin();
       pit != tuned_tracks.end(); pit++)
    _track_recombiner.preprocess(*pit);


  // we'll just need the particle indices for what follows
  list<int> sorted_pt_index;
  for (vector<TrackJetParticlePtr>::iterator mom_it = particle_list.begin();
       mom_it != particle_list.end(); mom_it++)
    sorted_pt_index.push_back(mom_it->index);
  
  // now start building the jets
  while (sorted_pt_index.size()){
    // note that here 'track' refers to the direction we're using to test if a particle belongs to the jet
    // 'jet' refers to the momentum of the jet
    // the difference between the two is in the recombination scheme used to compute the sum of 4-vectors
    int current_jet_index = sorted_pt_index.front();
    PseudoJet current_jet   = tuned_particles[current_jet_index];
    PseudoJet current_track = tuned_tracks[current_jet_index];

    // remove the first particle from the available ones    
    list<int>::iterator index_it = sorted_pt_index.begin();
    sorted_pt_index.erase(index_it);

    // start browsing the remaining ones
    index_it = sorted_pt_index.begin();
    while (index_it != sorted_pt_index.end()){
      const PseudoJet & current_particle = tuned_particles[*index_it];
      const PseudoJet & current_particle_track = tuned_tracks[*index_it];

      // check if the particle is within a distance R of the jet
      double distance2 = current_track.plain_distance(current_particle_track);
      if (distance2 <= _radius2){
        // add the particle to the jet
        PseudoJet new_track;
        PseudoJet new_jet;
        _jet_recombiner.recombine(current_jet, current_particle, new_jet);
        _track_recombiner.recombine(current_track, current_particle_track, new_track);

        int new_jet_index;
        clust_seq.plugin_record_ij_recombination(current_jet_index, *index_it, distance2, new_jet, new_jet_index);

        current_jet = new_jet;
        current_track = new_track;
        current_jet_index = new_jet_index;

        // particle has been clustered so remove it from the list
        sorted_pt_index.erase(index_it);

        // and don't forget to start again from the beginning
        //  as the jet axis may have changed
        index_it = sorted_pt_index.begin();
      } else {
        index_it++;
      }
    }

    // now we have a final jet, so cluster it with the beam
    clust_seq.plugin_record_iB_recombination(current_jet_index, _radius2);
  }
    
}

Member Data Documentation

Definition at line 83 of file TrackJetPlugin.hh.

Referenced by TrackJetPlugin().

double TrackJetPlugin::_radius [private]

Definition at line 81 of file TrackJetPlugin.hh.

Referenced by R(), and TrackJetPlugin().

double TrackJetPlugin::_radius2 [private]

Definition at line 81 of file TrackJetPlugin.hh.

Referenced by TrackJetPlugin().

Definition at line 84 of file TrackJetPlugin.hh.

Referenced by TrackJetPlugin().


The documentation for this class was generated from the following files: