Main Page | Namespace List | Class Hierarchy | Class List | Directories | File List | Namespace Members | Class Members | File Members

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=0.5, 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 = 0.5,
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 fastjet::JetDefinition::Plugin.

Definition at line 45 of file CDFJetCluPlugin.cc.

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

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

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 fastjet::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 fastjet::JetDefinition::Plugin.

Definition at line 60 of file CDFJetCluPlugin.cc.

References adjacency_cut(), cone_radius(), fastjet::PseudoJet::E(), iratch(), fastjet::ClusterSequence::jets(), max_iterations(), overlap_threshold(), fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), fastjet::PseudoJet::px(), fastjet::PseudoJet::py(), fastjet::PseudoJet::pz(), and seed_threshold().

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

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

int fastjet::CDFJetCluPlugin::_adjacency_cut [private]
 

Definition at line 93 of file CDFJetCluPlugin.hh.

double fastjet::CDFJetCluPlugin::_cone_radius [private]
 

Definition at line 92 of file CDFJetCluPlugin.hh.

int fastjet::CDFJetCluPlugin::_iratch [private]
 

Definition at line 95 of file CDFJetCluPlugin.hh.

int fastjet::CDFJetCluPlugin::_max_iterations [private]
 

Definition at line 94 of file CDFJetCluPlugin.hh.

double fastjet::CDFJetCluPlugin::_overlap_threshold [private]
 

Definition at line 96 of file CDFJetCluPlugin.hh.

double fastjet::CDFJetCluPlugin::_seed_threshold [private]
 

Definition at line 91 of file CDFJetCluPlugin.hh.


The documentation for this class was generated from the following files:
Generated on Fri Aug 15 13:45:47 2008 for fastjet by  doxygen 1.4.2