fastjet 2.4.5
Public Member Functions | Protected Member Functions | Private Member Functions
fastjet::ClusterSequence1GhostPassiveArea Class Reference

Class that behaves essentially like ClusterSequence except that it also provides access to the area of a jet (which will be a random quantity... More...

#include <ClusterSequence1GhostPassiveArea.hh>

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

List of all members.

Public Member Functions

 ClusterSequence1GhostPassiveArea ()
template<class L >
 ClusterSequence1GhostPassiveArea (const std::vector< L > &pseudojets, const JetDefinition &jet_def, const GhostedAreaSpec &area_spec, const bool &writeout_combinations=false)
 constructor based on JetDefinition and 1GhostPassiveAreaSpec
virtual double n_empty_jets (const RangeDefinition &range) const
 return an estimate for the number of empty jets -- one uses the AreaBase one rather than the ActiveArea one (which for which we do not have the information).

Protected Member Functions

void _initialise_and_run_1GPA (const JetDefinition &jet_def, const GhostedAreaSpec &area_spec, const bool &writeout_combinations=false)
 does the initialisation and running specific to the passive areas class

Private Member Functions

void _run_1GPA (const GhostedAreaSpec &area_spec)
 routine for running a passive area one ghost at a time.

Detailed Description

Class that behaves essentially like ClusterSequence except that it also provides access to the area of a jet (which will be a random quantity...

Figure out what to do about seeds later...)

Definition at line 49 of file ClusterSequence1GhostPassiveArea.hh.


Constructor & Destructor Documentation

fastjet::ClusterSequence1GhostPassiveArea::ClusterSequence1GhostPassiveArea ( ) [inline]

Definition at line 52 of file ClusterSequence1GhostPassiveArea.hh.

{}
template<class L >
fastjet::ClusterSequence1GhostPassiveArea::ClusterSequence1GhostPassiveArea ( const std::vector< L > &  pseudojets,
const JetDefinition jet_def,
const GhostedAreaSpec area_spec,
const bool &  writeout_combinations = false 
)

constructor based on JetDefinition and 1GhostPassiveAreaSpec

Definition at line 84 of file ClusterSequence1GhostPassiveArea.hh.

                                     {

  // transfer the initial jets (type L) into our own array
  _transfer_input_jets(pseudojets);

  // run the clustering for passive areas
  _initialise_and_run_1GPA(jet_def, area_spec, writeout_combinations);

}

Member Function Documentation

void fastjet::ClusterSequence1GhostPassiveArea::_initialise_and_run_1GPA ( const JetDefinition jet_def,
const GhostedAreaSpec area_spec,
const bool &  writeout_combinations = false 
) [protected]

does the initialisation and running specific to the passive areas class

global routine for initialising and running a general passive area

Definition at line 40 of file ClusterSequence1GhostPassiveArea.cc.

                                                    {

  bool continue_running;
  _initialise_AA(jet_def,  area_spec, writeout_combinations, continue_running);
  if (continue_running) {
    _run_1GPA(area_spec);
    _postprocess_AA(area_spec);
  }
}
void fastjet::ClusterSequence1GhostPassiveArea::_run_1GPA ( const GhostedAreaSpec area_spec) [private]

routine for running a passive area one ghost at a time.

Definition at line 56 of file ClusterSequence1GhostPassiveArea.cc.

References fastjet::GhostedAreaSpec::actual_ghost_area(), fastjet::GhostedAreaSpec::add_ghosts(), and fastjet::GhostedAreaSpec::repeat().

                                                                                   {
    // record the input jets as they are currently
  vector<PseudoJet> input_jets(_jets);

  // code for testing the unique tree
  vector<int> unique_tree;

  // initialise our temporary average area^2
  valarray<double> lcl_average_area2(0.0, _average_area.size());
  valarray<double> last_average_area(0.0, _average_area.size());

  // run the clustering multiple times so as to get areas of all the jets
  for (int irepeat = 0; irepeat < area_spec.repeat(); irepeat++) {

    // first establish list of all ghosts
    vector<PseudoJet> all_ghosts;
    area_spec.add_ghosts(all_ghosts);

    // then run many subsets of ghosts (actually each subset contains just one ghost)
    for (unsigned ig = 0; ig < all_ghosts.size(); ig++) {
      vector<PseudoJet> some_ghosts;
      some_ghosts.push_back(all_ghosts[ig]);
      ClusterSequenceActiveAreaExplicitGhosts clust_seq(input_jets, jet_def(), 
                                                       some_ghosts, area_spec.actual_ghost_area());

      if (irepeat == 0 && ig == 0) {
        // take the non-ghost part of the history and put into our own
        // history.
        _transfer_ghost_free_history(clust_seq);
        // get the "unique" order that will be used for transferring all areas. 
        unique_tree = unique_history_order();
      }
      
      // transfer areas from clust_seq into our object
      _transfer_areas(unique_tree, clust_seq);
    }
    // keep track of fluctuations in area
    lcl_average_area2 += (_average_area - last_average_area)*
                         (_average_area - last_average_area);
    last_average_area = _average_area;
  }
  _average_area2 = lcl_average_area2;
}
virtual double fastjet::ClusterSequence1GhostPassiveArea::n_empty_jets ( const RangeDefinition range) const [inline, virtual]

return an estimate for the number of empty jets -- one uses the AreaBase one rather than the ActiveArea one (which for which we do not have the information).

Reimplemented from fastjet::ClusterSequenceActiveArea.

Definition at line 64 of file ClusterSequence1GhostPassiveArea.hh.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines