fastjet::TrackJetPlugin Class Reference

#include <TrackJetPlugin.hh>

Inheritance diagram for fastjet::TrackJetPlugin:
Inheritance graph
[legend]
Collaboration diagram for fastjet::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(.

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

fastjet::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:
_radius the distance at which point a particle is no longer recombied into the jet
jet_recombination_scheme the recombination scheme used to sum the 4-vecors inside the jet
track_recombination_scheme the 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.

00060                                                                            {
00061     _radius  = radius;
00062     _radius2 = radius*radius;
00063     _jet_recombiner = JetDefinition::DefaultRecombiner(jet_recombination_scheme);
00064     _track_recombiner = JetDefinition::DefaultRecombiner(track_recombination_scheme);
00065   }

fastjet::TrackJetPlugin::TrackJetPlugin ( const TrackJetPlugin plugin  )  [inline]

copy constructor

Definition at line 68 of file TrackJetPlugin.hh.

00068                                                  {
00069     *this = plugin;
00070   }


Member Function Documentation

string fastjet::TrackJetPlugin::description (  )  const [virtual]

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

Implements JetDefinition::Plugin.

Definition at line 59 of file TrackJetPlugin.cc.

References R().

00059                                           {
00060   ostringstream desc;
00061   desc << "TrackJet algorithm with R = " << R();
00062   return desc.str();
00063 }

virtual double fastjet::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.

Referenced by description().

00078 {return _radius;}

void fastjet::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 _jet_recombiner, _radius2, _track_recombiner, ClusterSequence::jets(), PseudoJet::plain_distance(), ClusterSequence::plugin_record_iB_recombination(), ClusterSequence::plugin_record_ij_recombination(), JetDefinition::DefaultRecombiner::preprocess(), and JetDefinition::DefaultRecombiner::recombine().

00065                                                                      {
00066   // we first need to sort the particles in pt
00067   vector<TrackJetParticlePtr> particle_list;
00068 
00069   const vector<PseudoJet> & jets = clust_seq.jets();  
00070   int index=0;
00071   for (vector<PseudoJet>::const_iterator mom_it = jets.begin(); mom_it != jets.end(); mom_it++){
00072     particle_list.push_back(TrackJetParticlePtr(index, mom_it->perp2()));
00073     index++;
00074   }
00075 
00076   // sort the particles into decreasing pt
00077   sort(particle_list.begin(), particle_list.end());
00078 
00079 
00080   // if we're using a recombination scheme different from the E scheme,
00081   // we first need to update the particles' energy so that they
00082   // are massless (and rapidity = pseudorapidity)
00083   vector<PseudoJet> tuned_particles = clust_seq.jets();
00084   vector<PseudoJet> tuned_tracks = clust_seq.jets();
00085   for (vector<PseudoJet>::iterator pit = tuned_particles.begin();
00086        pit != tuned_particles.end(); pit++)
00087     _jet_recombiner.preprocess(*pit);
00088   for (vector<PseudoJet>::iterator pit = tuned_tracks.begin();
00089        pit != tuned_tracks.end(); pit++)
00090     _track_recombiner.preprocess(*pit);
00091 
00092 
00093   // we'll just need the particle indices for what follows
00094   list<int> sorted_pt_index;
00095   for (vector<TrackJetParticlePtr>::iterator mom_it = particle_list.begin();
00096        mom_it != particle_list.end(); mom_it++)
00097     sorted_pt_index.push_back(mom_it->index);
00098   
00099   // now start building the jets
00100   while (sorted_pt_index.size()){
00101     // note that here 'track' refers to the direction we're using to test if a particle belongs to the jet
00102     // 'jet' refers to the momentum of the jet
00103     // the difference between the two is in the recombination scheme used to compute the sum of 4-vectors
00104     int current_jet_index = sorted_pt_index.front();
00105     PseudoJet current_jet   = tuned_particles[current_jet_index];
00106     PseudoJet current_track = tuned_tracks[current_jet_index];
00107 
00108     // remove the first particle from the available ones    
00109     list<int>::iterator index_it = sorted_pt_index.begin();
00110     sorted_pt_index.erase(index_it);
00111 
00112     // start browsing the remaining ones
00113     index_it = sorted_pt_index.begin();
00114     while (index_it != sorted_pt_index.end()){
00115       const PseudoJet & current_particle = tuned_particles[*index_it];
00116       const PseudoJet & current_particle_track = tuned_tracks[*index_it];
00117 
00118       // check if the particle is within a distance R of the jet
00119       double distance2 = current_track.plain_distance(current_particle_track);
00120       if (distance2 <= _radius2){
00121         // add the particle to the jet
00122         PseudoJet new_track;
00123         PseudoJet new_jet;
00124         _jet_recombiner.recombine(current_jet, current_particle, new_jet);
00125         _track_recombiner.recombine(current_track, current_particle_track, new_track);
00126 
00127         int new_jet_index;
00128         clust_seq.plugin_record_ij_recombination(current_jet_index, *index_it, distance2, new_jet, new_jet_index);
00129 
00130         current_jet = new_jet;
00131         current_track = new_track;
00132         current_jet_index = new_jet_index;
00133 
00134         // particle has been clustered so remove it from the list
00135         sorted_pt_index.erase(index_it);
00136 
00137         // and don't forget to start again from the beginning
00138         //  as the jet axis may have changed
00139         index_it = sorted_pt_index.begin();
00140       } else {
00141         index_it++;
00142       }
00143     }
00144 
00145     // now we have a final jet, so cluster it with the beam
00146     clust_seq.plugin_record_iB_recombination(current_jet_index, _radius2);
00147   }
00148     
00149 }


Member Data Documentation

Definition at line 83 of file TrackJetPlugin.hh.

Referenced by run_clustering().

Definition at line 81 of file TrackJetPlugin.hh.

Definition at line 81 of file TrackJetPlugin.hh.

Referenced by run_clustering().

Definition at line 84 of file TrackJetPlugin.hh.

Referenced by run_clustering().


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

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1