fastjet 2.4.3
Public Member Functions | Private Attributes

CMSIterativeConePlugin Class Reference

CMSIterativeConePlugin is a plugin for fastjet (v2.4 upwards) More...

#include <CMSIterativeConePlugin.hh>

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

List of all members.

Public Member Functions

 CMSIterativeConePlugin (double ConeRadius, double SeedThreshold=1.0)
 Main constructor for the CMSIterativeCone Plugin class.
 CMSIterativeConePlugin (const CMSIterativeConePlugin &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
virtual double seed_threshold () const
 get the seed threshold

Private Attributes

double theConeRadius
 cone radius
double theSeedThreshold
 seed threshold

Detailed Description

CMSIterativeConePlugin is a plugin for fastjet (v2.4 upwards)

Definition at line 47 of file CMSIterativeConePlugin.hh.


Constructor & Destructor Documentation

CMSIterativeConePlugin::CMSIterativeConePlugin ( double  ConeRadius,
double  SeedThreshold = 1.0 
) [inline]

Main constructor for the CMSIterativeCone Plugin class.

The arguments are ConeRadius the radius of the cone SeedThreshold a threshold for the seeds to iterate from

NOTE: to be more coherent with all other fastjet plugins, we've put the radius before the seed threshold. CMS does the opposite. In this way, we also put a default value of 0 for the seed threshold.

Definition at line 60 of file CMSIterativeConePlugin.hh.

                                                                       :
    theConeRadius(ConeRadius), theSeedThreshold(SeedThreshold){}
CMSIterativeConePlugin::CMSIterativeConePlugin ( const CMSIterativeConePlugin plugin) [inline]

copy constructor

Definition at line 64 of file CMSIterativeConePlugin.hh.

                                                                 {
    *this = plugin;
  }

Member Function Documentation

string CMSIterativeConePlugin::description ( ) const [virtual]

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

Implements JetDefinition::Plugin.

Definition at line 68 of file CMSIterativeConePlugin.cc.

                                                  {
  ostringstream desc;
  desc << "CMSIterativeCone plugin with R = " << theConeRadius << " and seed threshold = " << theSeedThreshold;
  return desc.str();
}
virtual double CMSIterativeConePlugin::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 74 of file CMSIterativeConePlugin.hh.

References theConeRadius.

{return theConeRadius;}
void CMSIterativeConePlugin::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 74 of file CMSIterativeConePlugin.cc.

References deltaR2(), ClusterSequence::jets(), M_PI, d0::inline_maths::phi(), ClusterSequence::plugin_record_iB_recombination(), and ClusterSequence::plugin_record_ij_recombination().

                                                                             {

  // This code is adapted from CMSIterativeConeAlgorithms.cc from the
  // CMSSW software. 
  // The adaptation is just meant to use 
  //   - the FastJet 4-vectors instead of the CMS ones
  //   - the FastJet clustering history structures instead of the 
  //     ProtoJet one used by CMS.

  //make a list of input objects ordered by ET
  //cout << "copying the list of particles" << endl;
  list<PseudoJet> input;
  for (unsigned int i=0 ; i<clust_seq.jets().size() ; i++) {
    input.push_back(clust_seq.jets()[i]);
  }
  NumericSafeGreaterByEt<PseudoJet> compCandidate;
  //cout << "sorting" << endl;
  input.sort(compCandidate);

  //find jets
  //cout << "launching the main loop" << endl;
  while( !input.empty() && (input.front().Et() > theSeedThreshold )) {
    //cone centre 
    double eta0=input.front().eta();
    double phi0=input.front().phi();
    //protojet properties
    double eta=0;
    double phi=0;
    double et=0;
    //list of towers in cone
    list< list<PseudoJet>::iterator> cone;
    for(int iteration=0;iteration<100;iteration++){
      //cout << "iterating" << endl;
      cone.clear();
      eta=0;
      phi=0;
      et=0;
      for(list<PseudoJet>::iterator inp=input.begin();
          inp!=input.end();inp++){
        const PseudoJet tower = *inp;   
        if( deltaR2(eta0,phi0,tower.eta(),tower.phi()) < 
            theConeRadius*theConeRadius) {
          double tower_et = tower.Et();   
          cone.push_back(inp);
          eta+= tower_et*tower.eta();
          double dphi=tower.phi()-phi0;
          if(dphi>M_PI) dphi-=2*M_PI;
          else if(dphi<=-M_PI) dphi+=2*M_PI;
          phi+=tower_et*dphi;
          et +=tower_et;
        }
      }
      eta=eta/et;
      phi=phi0+phi/et;
      if(phi>M_PI)phi-=2*M_PI;
      else if(phi<=-M_PI)phi+=2*M_PI;
      
      if(fabs(eta-eta0)<.001 && fabs(phi-phi0)<.001) break;//stable cone found
      eta0=eta;
      phi0=phi;
    }

    //cout << "make the jet final" << endl;

    //make a final jet and remove the jet constituents from the input list 
    //  InputCollection jetConstituents;     
    //  list< list<InputItem>::iterator>::const_iterator inp;
    //  for(inp=cone.begin();inp!=cone.end();inp++)  {
    //    jetConstituents.push_back(**inp);
    //    input.erase(*inp);
    //  }
    //  fOutput->push_back (ProtoJet (jetConstituents));
    //
    // IMPORTANT NOTE:
    //   while the stability of the stable cone is tested using the Et
    //   scheme recombination, the final jet uses E-scheme
    //   recombination.
    //
    // The technique used here is the same as what we already used for
    // SISCone except that we act on the 'cone' list.
    // We successively merge the particles that make up the cone jet
    // until we have all particles in it.  We start off with the zeroth
    // particle.
    list< list<PseudoJet>::iterator>::const_iterator inp;
    inp = cone.begin();
    int jet_k = (*inp)->cluster_hist_index();
    // gps tmp
    //float px=(*inp)->px(), py=(*inp)->py(), pz=(*inp)->pz(), E = (*inp)->E();

    // remove the particle from the list and jump to the next one
    input.erase(*inp);
    inp++;

    // now loop over the remaining particles
    while (inp != cone.end()){
      // take the last result of the merge
      int jet_i = jet_k;
      // and the next element of the jet
      int jet_j = (*inp)->cluster_hist_index();
      // and merge them (with a fake dij)
      double dij = 0.0;

      // create the new jet by hand so that we can adjust its user index
      // Note again the use of the E-scheme recombination here!
      PseudoJet newjet = clust_seq.jets()[jet_i] + clust_seq.jets()[jet_j];

      // gps tmp to try out floating issues
      //px+=(*inp)->px(), py+=(*inp)->py(), pz+=(*inp)->pz(), E += (*inp)->E();
      //PseudoJet newjet(px,py,pz,E);
        
      clust_seq.plugin_record_ij_recombination(jet_i, jet_j, dij, newjet, jet_k);

      // remove the particle from the list and jump to the next one
      input.erase(*inp);
      inp++;
    }

    // we have merged all the jet's particles into a single object, so now
    // "declare" it to be a beam (inclusive) jet.
    // [NB: put a sensible looking d_iB just to be nice...]
    double d_iB = clust_seq.jets()[jet_k].perp2();
    clust_seq.plugin_record_iB_recombination(jet_k, d_iB);


  } //loop over seeds ended

    
}
virtual double CMSIterativeConePlugin::seed_threshold ( ) const [inline, virtual]

get the seed threshold

Definition at line 77 of file CMSIterativeConePlugin.hh.

References theSeedThreshold.

{return theSeedThreshold;}

Member Data Documentation

cone radius

Definition at line 80 of file CMSIterativeConePlugin.hh.

Referenced by R().

seed threshold

Definition at line 81 of file CMSIterativeConePlugin.hh.

Referenced by seed_threshold().


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