fastjet::CDFJetCluPlugin Class Reference

a plugin for fastjet-v2.1 that provides an interface to the CDF jetclu algorithm More...

#include <CDFJetCluPlugin.hh>

Inheritance diagram for fastjet::CDFJetCluPlugin:
Inheritance graph
[legend]
Collaboration diagram for fastjet::CDFJetCluPlugin:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 CDFJetCluPlugin (double cone_radius, double overlap_threshold, double seed_threshold=1.0, int iratch=1)
 a compact constructor
 CDFJetCluPlugin (double seed_threshold, double cone_radius, int adjacency_cut, int max_iterations, int iratch, double overlap_threshold)
 a constructor that looks like the one provided by CDF
double seed_threshold () const
double cone_radius () const
int adjacency_cut () const
int max_iterations () const
int iratch () const
double overlap_threshold () const
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

Private Member Functions

void _insert_unique (PseudoJet &jet, std::map< double, int > &jetmap) const
 given a jet try inserting its energy into the map -- if that energy entry already exists, modify the jet infinitesimally so as ensure that the jet energy is unique

Private Attributes

double _seed_threshold
double _cone_radius
int _adjacency_cut
int _max_iterations
int _iratch
double _overlap_threshold

Detailed Description

a plugin for fastjet-v2.1 that provides an interface to the CDF jetclu algorithm

Definition at line 44 of file CDFJetCluPlugin.hh.


Constructor & Destructor Documentation

fastjet::CDFJetCluPlugin::CDFJetCluPlugin ( double  cone_radius,
double  overlap_threshold,
double  seed_threshold = 1.0,
int  iratch = 1 
) [inline]

a compact constructor

Definition at line 47 of file CDFJetCluPlugin.hh.

00050                                         : 
00051     _seed_threshold    ( seed_threshold    ),    
00052     _cone_radius       ( cone_radius       ),
00053     _adjacency_cut     (   2               ),
00054     _max_iterations    ( 100               ),
00055     _iratch            ( iratch            ),
00056     _overlap_threshold ( overlap_threshold )  {}

fastjet::CDFJetCluPlugin::CDFJetCluPlugin ( double  seed_threshold,
double  cone_radius,
int  adjacency_cut,
int  max_iterations,
int  iratch,
double  overlap_threshold 
) [inline]

a constructor that looks like the one provided by CDF

Definition at line 59 of file CDFJetCluPlugin.hh.


Member Function Documentation

void fastjet::CDFJetCluPlugin::_insert_unique ( PseudoJet jet,
std::map< double, int > &  jetmap 
) const [private]

given a jet try inserting its energy into the map -- if that energy entry already exists, modify the jet infinitesimally so as ensure that the jet energy is unique

int fastjet::CDFJetCluPlugin::adjacency_cut (  )  const [inline]

Definition at line 76 of file CDFJetCluPlugin.hh.

Referenced by description(), and run_clustering().

00076 {return _adjacency_cut     ;}

double fastjet::CDFJetCluPlugin::cone_radius (  )  const [inline]

Definition at line 75 of file CDFJetCluPlugin.hh.

Referenced by description(), and run_clustering().

00075 {return _cone_radius       ;}

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

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

Implements JetDefinition::Plugin.

Definition at line 46 of file CDFJetCluPlugin.cc.

References adjacency_cut(), cone_radius(), iratch(), max_iterations(), overlap_threshold(), and seed_threshold().

00046                                            {
00047   ostringstream desc;
00048   
00049   desc << "CDF JetClu jet algorithm with " 
00050        << "seed_threshold = "     << seed_threshold    () << ", "
00051        << "cone_radius = "        << cone_radius       () << ", "
00052        << "adjacency_cut = "      << adjacency_cut     () << ", " 
00053        << "max_iterations = "     << max_iterations    () << ", "
00054        << "iratch = "             << iratch            () << ", "
00055        << "overlap_threshold = "  << overlap_threshold () ;
00056 
00057   return desc.str();
00058 }

int fastjet::CDFJetCluPlugin::iratch (  )  const [inline]

Definition at line 78 of file CDFJetCluPlugin.hh.

Referenced by description(), and run_clustering().

00078 {return _iratch            ;}

int fastjet::CDFJetCluPlugin::max_iterations (  )  const [inline]

Definition at line 77 of file CDFJetCluPlugin.hh.

Referenced by description(), and run_clustering().

00077 {return _max_iterations    ;}

double fastjet::CDFJetCluPlugin::overlap_threshold (  )  const [inline]

Definition at line 79 of file CDFJetCluPlugin.hh.

Referenced by description(), and run_clustering().

00079 {return _overlap_threshold ;}

virtual double fastjet::CDFJetCluPlugin::R (  )  const [inline, virtual]

the plugin mechanism's standard way of accessing the jet radius

Implements JetDefinition::Plugin.

Definition at line 86 of file CDFJetCluPlugin.hh.

00086 {return cone_radius();}

void fastjet::CDFJetCluPlugin::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 61 of file CDFJetCluPlugin.cc.

References adjacency_cut(), cone_radius(), iratch(), ClusterSequence::jets(), max_iterations(), overlap_threshold(), ClusterSequence::plugin_record_iB_recombination(), ClusterSequence::plugin_record_ij_recombination(), seed_threshold(), and sort_indices().

00061                                                                       {
00062  
00063   // create the physics towers needed by the CDF code
00064   vector<PhysicsTower> towers;
00065   towers.reserve(clust_seq.jets().size());
00066 
00067   // create a map to identify jets (actually just the input particles)...
00068   //map<double,int> jetmap;
00069 
00070   for (unsigned i = 0; i < clust_seq.jets().size(); i++) {
00071     PseudoJet particle(clust_seq.jets()[i]);
00072     //_insert_unique(particle, jetmap);
00073     LorentzVector fourvect(particle.px(), particle.py(),
00074                            particle.pz(), particle.E());
00075     PhysicsTower tower(fourvect);
00076     // add tracking information for later
00077     tower.fjindex = i;
00078     towers.push_back(tower);
00079   }
00080 
00081   // prepare the CDF algorithm
00082   JetCluAlgorithm j(seed_threshold(), cone_radius(), adjacency_cut(),
00083                     max_iterations(), iratch(), overlap_threshold());
00084     
00085   // run the CDF algorithm
00086   std::vector<Cluster> jets;
00087   j.run(towers,jets);
00088 
00089 
00090   // now transfer the jets back into our own structure -- we will
00091   // mimic the cone code with a sequential recombination sequence in
00092   // which the jets are built up by adding one particle at a time
00093 
00094   // NB: with g++-4.0, the reverse iterator code gave problems, so switch
00095   //     to indices instead
00096   //for(vector<Cluster>::const_reverse_iterator jetIter = jets.rbegin(); 
00097   //                                    jetIter != jets.rend(); jetIter++) {
00098   //  const vector<PhysicsTower> & tower_list = jetIter->towerList;
00099   //  int jet_k = jetmap[tower_list[0].fourVector.E];
00100   //
00101   //  int ntow = int(jetIter->towerList.size());
00102 
00103   for(int iCDFjets = jets.size()-1; iCDFjets >= 0; iCDFjets--) {
00104 
00105     const vector<PhysicsTower> & tower_list = jets[iCDFjets].towerList;
00106     int ntow = int(tower_list.size());
00107     
00108     // 2008-09-04: sort the towerList (according to fjindex) so as
00109     //             to have a consistent order for particles in jet
00110     //             (necessary because addition of ultra-soft particles
00111     //             sometimes often modifies the order, while maintaining
00112     //             the same overall set)
00113     vector<int>    jc_indices(ntow);
00114     vector<double> fj_indices(ntow); // use double: benefit from existing routine
00115     for (int itow = 0; itow < ntow; itow++) {
00116       jc_indices[itow] = itow;
00117       fj_indices[itow] = tower_list[itow].fjindex;
00118     }
00119     sort_indices(jc_indices, fj_indices);
00120 
00121     int jet_k = tower_list[jc_indices[0]].fjindex;
00122   
00123     for (int itow = 1; itow < ntow; itow++) {
00124       if (tower_list[jc_indices[itow]].Et() > 1e-50) {
00125       }
00126       int jet_i = jet_k;
00127       // retrieve our index for the jet
00128       int jet_j;
00129       jet_j = tower_list[jc_indices[itow]].fjindex;
00130 
00131       // safety check
00132       assert (jet_j >= 0 && jet_j < int(towers.size()));
00133 
00134       // do a fake recombination step with dij=0
00135       double dij = 0.0;
00136 
00137       // JetClu does E-scheme recombination so we can stick with the
00138       // simple option
00139       clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
00140 
00141     }
00142   
00143     // NB: put a sensible looking d_iB just to be nice...
00144     double d_iB = clust_seq.jets()[jet_k].perp2();
00145     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00146   }
00147 
00148 
00149   // following code is for testing only
00150   //cout << endl;
00151   //for(vector<Cluster>::const_iterator jetIter = jets.begin(); 
00152   //                                    jetIter != jets.end(); jetIter++) {
00153   //  cout << jetIter->fourVector.pt() << " " << jetIter->fourVector.y() << endl;
00154   //}
00155   //cout << "-----------------------------------------------------\n";
00156   //vector<PseudoJet> ourjets(clust_seq.inclusive_jets());
00157   //for (vector<PseudoJet>::const_iterator ourjet = ourjets.begin();
00158   //     ourjet != ourjets.end(); ourjet++) {
00159   //  cout << ourjet->perp() << " " << ourjet->rap() << endl;
00160   //}
00161   //cout << endl;
00162 }

double fastjet::CDFJetCluPlugin::seed_threshold (  )  const [inline]

Definition at line 74 of file CDFJetCluPlugin.hh.

Referenced by description(), and run_clustering().

00074 {return _seed_threshold    ;}


Member Data Documentation

Definition at line 93 of file CDFJetCluPlugin.hh.

Definition at line 92 of file CDFJetCluPlugin.hh.

Definition at line 95 of file CDFJetCluPlugin.hh.

Definition at line 94 of file CDFJetCluPlugin.hh.

Definition at line 96 of file CDFJetCluPlugin.hh.

Definition at line 91 of file CDFJetCluPlugin.hh.


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

Generated on 26 Feb 2010 for fastjet by  doxygen 1.6.1