fastjet 2.4.5
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
fastjet::GhostedAreaSpec Class Reference

Class that defines the parameters that go into the measurement of active jet areas. More...

#include <GhostedAreaSpec.hh>

Collaboration diagram for fastjet::GhostedAreaSpec:
Collaboration graph
[legend]

List of all members.

Public Member Functions

 GhostedAreaSpec ()
 default constructor
 GhostedAreaSpec (double ghost_maxrap, int repeat=gas::def_repeat, double ghost_area=gas::def_ghost_area, double grid_scatter=gas::def_grid_scatter, double kt_scatter=gas::def_kt_scatter, double mean_ghost_kt=gas::def_mean_ghost_kt)
 explicit constructor
 GhostedAreaSpec (double ghost_minrap, double ghost_maxrap, int repeat=gas::def_repeat, double ghost_area=gas::def_ghost_area, double grid_scatter=gas::def_grid_scatter, double kt_scatter=gas::def_kt_scatter, double mean_ghost_kt=gas::def_mean_ghost_kt)
 explicit constructor
void _initialize ()
 does the initialization of actual ghost parameters
double ghost_etamax () const
double ghost_maxrap () const
double ghost_area () const
double grid_scatter () const
double kt_scatter () const
double mean_ghost_kt () const
int repeat () const
double actual_ghost_area () const
int n_ghosts () const
void set_ghost_area (double val)
void set_ghost_etamax (double val)
void set_ghost_maxrap (double val)
void set_grid_scatter (double val)
void set_kt_scatter (double val)
void set_mean_ghost_kt (double val)
void set_repeat (int val)
int nphi () const
 return nphi (ghosts layed out (-nrap, 0..nphi-1), (-nrap+1,0..nphi-1), ...
int nrap () const
void get_random_status (std::vector< int > &__iseed) const
 get all relevant information about the status of the random number generator, so that it can be reset subsequently with set_random_status.
void set_random_status (const std::vector< int > &__iseed)
 set the status of the random number generator, as obtained previously with get_random_status.
void checkpoint_random ()
void restore_checkpoint_random ()
std::string description () const
 for a summary
void add_ghosts (std::vector< PseudoJet > &) const
 push the ghost 4-momenta onto the back of the vector of PseudoJets
double random_at_own_risk () const
 very deprecated public access to a random number from the internal generator
BasicRandom< double > & generator_at_own_risk () const
 very deprecated public access to the generator itself

Private Member Functions

double _our_rand () const

Private Attributes

double _ghost_maxrap
double _ghost_rap_offset
int _repeat
double _ghost_area
double _grid_scatter
double _kt_scatter
double _mean_ghost_kt
double _actual_ghost_area
double _dphi
double _drap
int _n_ghosts
int _nphi
int _nrap
std::vector< int > _random_checkpoint

Static Private Attributes

static BasicRandom< double > _random_generator

Detailed Description

Class that defines the parameters that go into the measurement of active jet areas.

Definition at line 58 of file GhostedAreaSpec.hh.


Constructor & Destructor Documentation

fastjet::GhostedAreaSpec::GhostedAreaSpec ( ) [inline]
fastjet::GhostedAreaSpec::GhostedAreaSpec ( double  ghost_maxrap,
int  repeat = gas::def_repeat,
double  ghost_area = gas::def_ghost_area,
double  grid_scatter = gas::def_grid_scatter,
double  kt_scatter = gas::def_kt_scatter,
double  mean_ghost_kt = gas::def_mean_ghost_kt 
) [inline, explicit]
fastjet::GhostedAreaSpec::GhostedAreaSpec ( double  ghost_minrap,
double  ghost_maxrap,
int  repeat = gas::def_repeat,
double  ghost_area = gas::def_ghost_area,
double  grid_scatter = gas::def_grid_scatter,
double  kt_scatter = gas::def_kt_scatter,
double  mean_ghost_kt = gas::def_mean_ghost_kt 
) [inline, explicit]

Member Function Documentation

void fastjet::GhostedAreaSpec::_initialize ( )

does the initialization of actual ghost parameters

sets the detailed parameters for the ghosts (which may not be quite the same as those requested -- this is in order for things to fit in nicely into 2pi etc...

Definition at line 45 of file GhostedAreaSpec.cc.

References twopi.

                                  {
  // add on area-measuring dummy particles
  _drap = sqrt(_ghost_area);
  _dphi = _drap;
  _nphi = int(ceil(twopi/_dphi)); _dphi = twopi/_nphi;
  _nrap = int(ceil(_ghost_maxrap/_drap)); _drap = _ghost_maxrap / _nrap;
  _actual_ghost_area = _dphi * _drap;
  _n_ghosts   = (2*_nrap+1)*_nphi;

  // checkpoint the status of the random number generator.
  checkpoint_random();
  //_random_generator.info(cerr);
}
double fastjet::GhostedAreaSpec::_our_rand ( ) const [inline, private]

Definition at line 186 of file GhostedAreaSpec.hh.

{return _random_generator();};
double fastjet::GhostedAreaSpec::actual_ghost_area ( ) const [inline]
void fastjet::GhostedAreaSpec::add_ghosts ( std::vector< PseudoJet > &  ) const

push the ghost 4-momenta onto the back of the vector of PseudoJets

adds the ghost 4-momenta to the vector of PseudoJet's

Definition at line 61 of file GhostedAreaSpec.cc.

References fastjet::d0::inline_maths::phi().

Referenced by fastjet::ClusterSequenceActiveAreaExplicitGhosts::_add_ghosts(), fastjet::ClusterSequence1GhostPassiveArea::_run_1GPA(), and main().

                                                                {
  // add momenta for ghosts
  for (int irap = -_nrap; irap <= _nrap; irap++) {
    for (int iphi = 0; iphi < _nphi; iphi++) {
     
      // include random offsets for all quantities
      double phi = (iphi+0.5) * _dphi + _dphi*(_our_rand()-0.5)*_grid_scatter;
      double rap = irap * _drap + _drap*(_our_rand()-0.5)*_grid_scatter
                                                         + _ghost_rap_offset ;
      double kt = _mean_ghost_kt*(1+(_our_rand()-0.5)*_kt_scatter);

      double pminus = kt*exp(-rap);
      double pplus  = kt*exp(+rap);
      double px = kt*sin(phi);
      double py = kt*cos(phi);
      //cout << kt<<" "<<rap<<" "<<phi<<"\n";
      //if (phi>=twopi || phi < 0.0) cout << "Hey: "<< phi-twopi<<"\n";
      PseudoJet mom(px,py,0.5*(pplus-pminus),0.5*(pplus+pminus));
      event.push_back(mom);
    }
  }
}
void fastjet::GhostedAreaSpec::checkpoint_random ( ) [inline]
string fastjet::GhostedAreaSpec::description ( ) const

for a summary

Definition at line 84 of file GhostedAreaSpec.cc.

                                          {

  ostringstream ostr;
  ostr << "ghosts of area " << actual_ghost_area() 
       << " (had requested " << ghost_area() << ")"
       << ", placed up to y = " << ghost_maxrap() 
       << ", scattered wrt to perfect grid by (rel) " << grid_scatter() 
       << ", mean_ghost_kt = " << mean_ghost_kt()
       << ", rel kt_scatter =  " << kt_scatter()
       << ", n repetitions of ghost distributions =  " << repeat();
  return ostr.str();
}
BasicRandom<double>& fastjet::GhostedAreaSpec::generator_at_own_risk ( ) const [inline]

very deprecated public access to the generator itself

Definition at line 163 of file GhostedAreaSpec.hh.

                                                             {
    return _random_generator;}
void fastjet::GhostedAreaSpec::get_random_status ( std::vector< int > &  __iseed) const [inline]

get all relevant information about the status of the random number generator, so that it can be reset subsequently with set_random_status.

Definition at line 139 of file GhostedAreaSpec.hh.

                                                                {
    _random_generator.get_status(__iseed);}
double fastjet::GhostedAreaSpec::ghost_area ( ) const [inline]

Definition at line 112 of file GhostedAreaSpec.hh.

{return _ghost_area   ;};
double fastjet::GhostedAreaSpec::ghost_etamax ( ) const [inline]

Definition at line 110 of file GhostedAreaSpec.hh.

{return _ghost_maxrap;};
double fastjet::GhostedAreaSpec::ghost_maxrap ( ) const [inline]

Definition at line 111 of file GhostedAreaSpec.hh.

Referenced by fastjet::ClusterSequenceActiveArea::_initialise_AA().

{return _ghost_maxrap;};
double fastjet::GhostedAreaSpec::grid_scatter ( ) const [inline]

Definition at line 113 of file GhostedAreaSpec.hh.

{return _grid_scatter;};
double fastjet::GhostedAreaSpec::kt_scatter ( ) const [inline]

Definition at line 114 of file GhostedAreaSpec.hh.

{return _kt_scatter  ;};
double fastjet::GhostedAreaSpec::mean_ghost_kt ( ) const [inline]
int fastjet::GhostedAreaSpec::n_ghosts ( ) const [inline]
int fastjet::GhostedAreaSpec::nphi ( ) const [inline]

return nphi (ghosts layed out (-nrap, 0..nphi-1), (-nrap+1,0..nphi-1), ...

(nrap,0..nphi-1)

Definition at line 133 of file GhostedAreaSpec.hh.

{return _nphi;}
int fastjet::GhostedAreaSpec::nrap ( ) const [inline]

Definition at line 134 of file GhostedAreaSpec.hh.

{return _nrap;}
double fastjet::GhostedAreaSpec::random_at_own_risk ( ) const [inline]

very deprecated public access to a random number from the internal generator

Definition at line 161 of file GhostedAreaSpec.hh.

{return _our_rand();};
int fastjet::GhostedAreaSpec::repeat ( ) const [inline]
void fastjet::GhostedAreaSpec::restore_checkpoint_random ( ) [inline]
void fastjet::GhostedAreaSpec::set_ghost_area ( double  val) [inline]

Definition at line 123 of file GhostedAreaSpec.hh.

{_ghost_area    = val; _initialize();};
void fastjet::GhostedAreaSpec::set_ghost_etamax ( double  val) [inline]

Definition at line 124 of file GhostedAreaSpec.hh.

void fastjet::GhostedAreaSpec::set_ghost_maxrap ( double  val) [inline]

Definition at line 125 of file GhostedAreaSpec.hh.

void fastjet::GhostedAreaSpec::set_grid_scatter ( double  val) [inline]

Definition at line 126 of file GhostedAreaSpec.hh.

Referenced by main().

{_grid_scatter   = val; };
void fastjet::GhostedAreaSpec::set_kt_scatter ( double  val) [inline]

Definition at line 127 of file GhostedAreaSpec.hh.

{_kt_scatter     = val; };
void fastjet::GhostedAreaSpec::set_mean_ghost_kt ( double  val) [inline]

Definition at line 128 of file GhostedAreaSpec.hh.

{_mean_ghost_kt  = val; };
void fastjet::GhostedAreaSpec::set_random_status ( const std::vector< int > &  __iseed) [inline]

set the status of the random number generator, as obtained previously with get_random_status.

Note that the random generator is a static member of the class, i.e. common to all instances of the class --- so if you modify the random for this instance, you modify it for all instances.

Definition at line 147 of file GhostedAreaSpec.hh.

                                                                {
    _random_generator.set_status(__iseed);}
void fastjet::GhostedAreaSpec::set_repeat ( int  val) [inline]

Definition at line 129 of file GhostedAreaSpec.hh.

{_repeat         = val; };

Member Data Documentation

Definition at line 178 of file GhostedAreaSpec.hh.

Definition at line 178 of file GhostedAreaSpec.hh.

Definition at line 178 of file GhostedAreaSpec.hh.

Definition at line 172 of file GhostedAreaSpec.hh.

Definition at line 169 of file GhostedAreaSpec.hh.

Definition at line 170 of file GhostedAreaSpec.hh.

Definition at line 173 of file GhostedAreaSpec.hh.

Definition at line 174 of file GhostedAreaSpec.hh.

Definition at line 175 of file GhostedAreaSpec.hh.

Definition at line 179 of file GhostedAreaSpec.hh.

Definition at line 179 of file GhostedAreaSpec.hh.

Definition at line 179 of file GhostedAreaSpec.hh.

Definition at line 182 of file GhostedAreaSpec.hh.

Definition at line 183 of file GhostedAreaSpec.hh.

Definition at line 171 of file GhostedAreaSpec.hh.


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