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

fastjet::PxConePlugin Class Reference

PxConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the fortran pxcone iterative cone algorithm with midpoint seeds. More...

#include <PxConePlugin.hh>

Inheritance diagram for fastjet::PxConePlugin:

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

Collaboration graph
[legend]
List of all members.

Public Member Functions

 PxConePlugin (double cone_radius, double min_jet_energy=5.0, double overlap_threshold=0.5, bool E_scheme_jets=false)
 constructor for the PxConePlugin, whose arguments have the following meaning:
double cone_radius () const
 the cone radius
double min_jet_energy () const
 minimum jet energy (protojets below this are thrown own before merging/splitting) -- called epslon in pxcone
double overlap_threshold () const
 Maximum fraction of overlap energy in a jet -- called ovlim in pxcone.
bool E_scheme_jets () const
 if true then the final jets are returned as the E-scheme recombination of the particle momenta (by default, pxcone returns massless jets with a mean phi,eta type of recombination); regardless of what is returned, the internal pxcone jet-finding procedure is unaffected.
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(.


Private Attributes

double _cone_radius
double _min_jet_energy
double _overlap_threshold
bool _E_scheme_jets

Detailed Description

PxConePlugin is a plugin for fastjet (v2.1 upwards) that provides an interface to the fortran pxcone iterative cone algorithm with midpoint seeds.

Pxcone was written by Luis del Pozo and Michael H. Seymour. It is not a "supported" program, so if you encounter problems, you are on your own...

Note that pxcone sometimes encounters non-stable iterations; in such cases it returns an error -- the plugin propagates this by throwing a fastjet::Error exception; if the user wishes to have robust code, they should catch this exception.

Pxcone has a hard-coded limit (by default 4000) on the maximum number of particles and protojets; if the number of particles or protojets exceeds this, again a fastjet::Error exception will be thrown.

The functionality of pxcone is described at http://www.hep.man.ac.uk/u/wplano/ConeJet.ps

Definition at line 65 of file PxConePlugin.hh.


Constructor & Destructor Documentation

fastjet::PxConePlugin::PxConePlugin double  cone_radius,
double  min_jet_energy = 5.0,
double  overlap_threshold = 0.5,
bool  E_scheme_jets = false
[inline]
 

constructor for the PxConePlugin, whose arguments have the following meaning:

  • the cone_radius is as usual in cone algorithms

  • stables cones (protojets) below min_jet_energy are discarded before calling the splitting procedure to resolve overlaps (called epslon in pxcone).

  • when two protojets overlap, if (overlapping_Et)/(Et_of_softer_protojet) < overlap_threshold the overlapping energy is split between the two protojets; otherwise the less energetic protojet is discarded. Called ovlim in pxcone.

  • pxcone carries out p-scheme recombination, and the resulting jets are massless; setting E_scheme_jets = true (default false) doesn't change the jet composition, but the final momentum sum for the jets is carried out by direct four-vector addition instead of p-scheme recombination.

Definition at line 89 of file PxConePlugin.hh.


Member Function Documentation

double fastjet::PxConePlugin::cone_radius  )  const [inline]
 

the cone radius

Definition at line 102 of file PxConePlugin.hh.

Referenced by description(), and run_clustering().

00102 {return _cone_radius        ;}

string fastjet::PxConePlugin::description  )  const [virtual]
 

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

Implements fastjet::JetDefinition::Plugin.

Definition at line 44 of file PxConePlugin.cc.

References cone_radius(), E_scheme_jets(), min_jet_energy(), and overlap_threshold().

00044                                         {
00045   ostringstream desc;
00046   
00047   desc << "PxCone jet algorithm with " 
00048        << "cone_radius = "        << cone_radius        () << ", "
00049        << "min_jet_energy = "     << min_jet_energy     () << ", "
00050        << "overlap_threshold  = " << overlap_threshold  () << ", "
00051        << "E_scheme_jets  = "     << E_scheme_jets      () 
00052        << " (NB: non-standard version of PxCone, containing small bug fixes by Gavin Salam)";
00053 
00054   return desc.str();
00055 }

bool fastjet::PxConePlugin::E_scheme_jets  )  const [inline]
 

if true then the final jets are returned as the E-scheme recombination of the particle momenta (by default, pxcone returns massless jets with a mean phi,eta type of recombination); regardless of what is returned, the internal pxcone jet-finding procedure is unaffected.

Definition at line 116 of file PxConePlugin.hh.

Referenced by description(), and run_clustering().

00116 {return _E_scheme_jets      ;}

double fastjet::PxConePlugin::min_jet_energy  )  const [inline]
 

minimum jet energy (protojets below this are thrown own before merging/splitting) -- called epslon in pxcone

Definition at line 106 of file PxConePlugin.hh.

Referenced by description(), and run_clustering().

00106 {return _min_jet_energy     ;}

double fastjet::PxConePlugin::overlap_threshold  )  const [inline]
 

Maximum fraction of overlap energy in a jet -- called ovlim in pxcone.

Definition at line 109 of file PxConePlugin.hh.

Referenced by description(), and run_clustering().

00109 {return _overlap_threshold  ;}

void fastjet::PxConePlugin::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 58 of file PxConePlugin.cc.

References cone_radius(), E_scheme_jets(), fastjet::ClusterSequence::jets(), min_jet_energy(), overlap_threshold(), fastjet::ClusterSequence::plugin_record_iB_recombination(), fastjet::ClusterSequence::plugin_record_ij_recombination(), and pxcone.

00058                                                                    {
00059  
00060   // only have hh mode
00061   int mode = 2;
00062 
00063   int    ntrak = clust_seq.jets().size(), itkdm = 4;
00064   double ptrak[ntrak][4];
00065   for (int i = 0; i < ntrak; i++) {
00066     ptrak[i][0] = clust_seq.jets()[i].px();
00067     ptrak[i][1] = clust_seq.jets()[i].py();
00068     ptrak[i][2] = clust_seq.jets()[i].pz();
00069     ptrak[i][3] = clust_seq.jets()[i].E();
00070   }  
00071 
00072   // max number of allowed jets
00073   int mxjet = ntrak;
00074   int njet;
00075   double pjet[mxjet][5];
00076   int    ipass[ntrak], ijmul[mxjet], ierr;
00077 
00078   // run pxcone
00079   pxcone(
00080     mode   ,    // 1=>e+e-, 2=>hadron-hadron
00081     ntrak  ,    // Number of particles
00082     itkdm  ,    // First dimension of PTRAK array: 
00083     &(ptrak[0][0])  ,    // Array of particle 4-momenta (Px,Py,Pz,E)
00084     cone_radius()  ,    // Cone size (half angle) in radians
00085     min_jet_energy() ,    // Minimum Jet energy (GeV)
00086     overlap_threshold()  ,    // Maximum fraction of overlap energy in a jet
00087     mxjet  ,    // Maximum possible number of jets
00088     njet   ,    // Number of jets found
00089     &(pjet[0][0]),       // 5-vectors of jets
00090     ipass,      // Particle k belongs to jet number IPASS(k)-1
00091                 // IPASS = -1 if not assosciated to a jet
00092     ijmul,      // Jet i contains IJMUL[i] particles
00093     ierr        // = 0 if all is OK ;   = -1 otherwise
00094     );
00095 
00096   if (ierr != 0) throw Error("An error occurred while running PXCONE");
00097 
00098   // now transfer information back 
00099   valarray<int> last_index_created(njet);
00100 
00101   vector<vector<int> > jet_particle_content(njet);
00102 
00103   // get a list of particles in each jet
00104   for (int itrak = 0; itrak < ntrak; itrak++) {
00105     int jet_i = ipass[itrak] - 1;
00106     if (jet_i >= 0) jet_particle_content[jet_i].push_back(itrak);
00107   }
00108 
00109   // now transfer the jets back into our own structure -- we will
00110   // mimic the cone code with a sequential recombination sequence in
00111   // which the jets are built up by adding one particle at a time
00112   for(int ipxjet = njet-1; ipxjet >= 0; ipxjet--) {
00113     const vector<int> & jet_trak_list = jet_particle_content[ipxjet];
00114     int jet_k = jet_trak_list[0];
00115   
00116     for (unsigned ilist = 1; ilist < jet_trak_list.size(); ilist++) {
00117       int jet_i = jet_k;
00118       // retrieve our misappropriated index for the jet
00119       int jet_j = jet_trak_list[ilist];
00120       // do a fake recombination step with dij=0
00121       double dij = 0.0;
00122       //clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
00123       if (ilist != jet_trak_list.size()-1 || E_scheme_jets()) {
00124         // our E-scheme recombination in cases where it doesn't matter
00125         clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, jet_k);
00126       } else {
00127         // put in pxcone's momentum for the last recombination so that the
00128         // final inclusive jet corresponds exactly to PXCONE's
00129         clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, 
00130                                 PseudoJet(pjet[ipxjet][0],pjet[ipxjet][1],
00131                                           pjet[ipxjet][2],pjet[ipxjet][3]),
00132                                                  jet_k);
00133       }
00134     }
00135   
00136     // NB: put a sensible looking d_iB just to be nice...
00137     double d_iB = clust_seq.jets()[jet_k].perp2();
00138     clust_seq.plugin_record_iB_recombination(jet_k, d_iB);
00139   }
00140 
00141 
00143   //cout << endl;
00144   //for (int ijet = 0; ijet < njet; ijet++) {
00145   //  PseudoJet jet(pjet[ijet][0],pjet[ijet][1],pjet[ijet][2],pjet[ijet][3]);
00146   //  cout << jet.perp() << " " << jet.rap() << endl;
00147   //}
00148   //cout << "-----------------------------------------------------\n";
00149   //vector<PseudoJet> ourjets(clust_seq.inclusive_jets());
00150   //for (vector<PseudoJet>::const_iterator ourjet = ourjets.begin();
00151   //     ourjet != ourjets.end(); ourjet++) {
00152   //  cout << ourjet->perp() << " " << ourjet->rap() << endl;
00153   //}
00155 }


Member Data Documentation

double fastjet::PxConePlugin::_cone_radius [private]
 

Definition at line 124 of file PxConePlugin.hh.

bool fastjet::PxConePlugin::_E_scheme_jets [private]
 

Definition at line 128 of file PxConePlugin.hh.

double fastjet::PxConePlugin::_min_jet_energy [private]
 

Definition at line 125 of file PxConePlugin.hh.

double fastjet::PxConePlugin::_overlap_threshold [private]
 

Definition at line 126 of file PxConePlugin.hh.


The documentation for this class was generated from the following files:
Generated on Mon Apr 2 20:58:24 2007 for fastjet by  doxygen 1.4.2