rivet is hosted by Hepforge, IPPP Durham

#include <JetShape.hh>

List of all members.

Public Member Functions

void clear ()
 Reset projection between events.
void calc (const Jets &jets)
 Do the calculation directly on a supplied collection of Jet objects.
size_t numBins () const
 Number of equidistant radius bins.
size_t numJets () const
 Number of jets which passed cuts.
double rMin () const
 $ r_\text{min} $ value.
double rMax () const
 $ r_\text{max} $ value.
double ptMin () const
 $ p_\perp^\text{min} $ value.
double ptMax () const
 $ p_\perp^\text{max} $ value.
double rBinMin (size_t rbin) const
 Central $ r $ value for bin rbin.
double rBinMax (size_t rbin) const
 Central $ r $ value for bin rbin.
double rBinMid (size_t rbin) const
 Central $ r $ value for bin rbin.
double diffJetShape (size_t ijet, size_t rbin) const
 Return value of differential jet shape profile histo bin.
double intJetShape (size_t ijet, size_t rbin) const
 Return value of integrated jet shape profile histo bin.
bool before (const Projection &p) const
virtual const std::set< PdgIdPairbeamPairs () const
virtual std::string name () const
 Get the name of the projection.
ProjectionaddPdgIdPair (PdgId beam1, PdgId beam2)
 Add a colliding beam pair.
LoggetLog () const
 Get a Log object based on the getName() property of the calling projection object.
void setName (const std::string &name)
 Used by derived classes to set their name.
void markAsOwned () const
Constructors etc.
 JetShape (const JetAlg &jetalg, double rmin, double rmax, size_t nbins, double ptmin=0, double ptmax=MAXDOUBLE, double absrapmin=-MAXDOUBLE, double absrapmax=-MAXDOUBLE, RapScheme rapscheme=RAPIDITY)
 Constructor from histo range and number of bins.
 JetShape (const JetAlg &jetalg, vector< double > binedges, double ptmin=0, double ptmax=MAXDOUBLE, double absrapmin=-MAXDOUBLE, double absrapmax=-MAXDOUBLE, RapScheme rapscheme=RAPIDITY)
 Constructor from vector of bin edges.
 DEFAULT_RIVET_PROJ_CLONE (JetShape)
 Clone on the heap.
Standard constructors and destructors.
virtual unique_ptr< Projectionclone () const =0
 Clone on the heap.
Projection "getting" functions
std::set< ConstProjectionPtrgetProjections () const
 Get the contained projections, including recursion.
template<typename PROJ >
const PROJ & getProjection (const std::string &name) const
const ProjectiongetProjection (const std::string &name) const
template<typename PROJ >
const PROJ & get (const std::string &name) const
Projection applying functions
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const Projection &proj) const
 Apply the supplied projection on event evt.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const PROJ &proj) const
 Apply the supplied projection on event evt.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const std::string &name) const
template<typename PROJ >
const PROJ & apply (const Event &evt, const Projection &proj) const
template<typename PROJ >
const PROJ & apply (const Event &evt, const PROJ &proj) const
template<typename PROJ >
const PROJ & apply (const Event &evt, const std::string &name) const

Protected Member Functions

void project (const Event &e)
 Apply the projection to the event.
int compare (const Projection &p) const
 Compare projections.
Cmp< ProjectionmkNamedPCmp (const Projection &otherparent, const std::string &pname) const
Cmp< ProjectionmkPCmp (const Projection &otherparent, const std::string &pname) const
ProjectionHandlergetProjHandler () const
 Get a reference to the ProjectionHandler for this thread.
const Projection_applyProjection (const Event &evt, const std::string &name) const
const Projection_applyProjection (const Event &evt, const Projection &proj) const
Projection registration functions
template<typename PROJ >
const PROJ & declareProjection (const PROJ &proj, const std::string &name)
 Register a contained projection.
template<typename PROJ >
const PROJ & declare (const PROJ &proj, const std::string &name)
 Register a contained projection (user-facing version)
template<typename PROJ >
const PROJ & addProjection (const PROJ &proj, const std::string &name)
 Register a contained projection (user-facing version)
const Projection_declareProjection (const Projection &proj, const std::string &name)
 Untemplated function to do the work...

Protected Attributes

bool _allowProjReg
 Flag to forbid projection registration in analyses until the init phase.

Private Attributes

Jet shape parameters
vector< double > _binedges
 Vector of radius bin edges.
pair< double, double > _ptcuts
 Lower and upper cuts on contributing jet $ p_\perp $.
pair< double, double > _rapcuts
 Lower and upper cuts on contributing jet (pseudo)rapidity.
RapScheme _rapscheme
 Rapidity scheme.
The projected jet shapes
vector< vector< double > > _diffjetshapes
 Jet shape histo -- first index is jet number, second is r bin.

Friends

class Event
 Event is a friend.
class Cmp< Projection >
 The Cmp specialization for Projection is a friend.

Detailed Description

Definition at line 45 of file JetShape.hh.


Constructor & Destructor Documentation

JetShape ( const JetAlg jetalg,
double  rmin,
double  rmax,
size_t  nbins,
double  ptmin = 0,
double  ptmax = MAXDOUBLE,
double  absrapmin = -MAXDOUBLE,
double  absrapmax = -MAXDOUBLE,
RapScheme  rapscheme = RAPIDITY 
)

Constructor from histo range and number of bins.

Definition at line 9 of file JetShape.cc.

    : _rapscheme(rapscheme)
  {
    setName("JetShape");
    _binedges = linspace(nbins, rmin, rmax);
    _ptcuts = make_pair(ptmin, ptmax);
    _rapcuts = make_pair(absrapmin, absrapmax);
    addProjection(jetalg, "Jets");
  }
JetShape ( const JetAlg jetalg,
vector< double >  binedges,
double  ptmin = 0,
double  ptmax = MAXDOUBLE,
double  absrapmin = -MAXDOUBLE,
double  absrapmax = -MAXDOUBLE,
RapScheme  rapscheme = RAPIDITY 
)

Constructor from vector of bin edges.

Definition at line 25 of file JetShape.cc.

    : _binedges(binedges), _rapscheme(rapscheme)
  {
    setName("JetShape");
    _ptcuts = make_pair(ptmin, ptmax);
    _rapcuts = make_pair(absrapmin, absrapmax);
    addProjection(jetalg, "Jets");
  }

Member Function Documentation

const Projection & _applyProjection ( const Event evt,
const std::string &  name 
) const [protected, inherited]

Non-templated version of string-based applyProjection, to work around header dependency issue.

Definition at line 22 of file ProjectionApplier.cc.

                                                                                  {
    return evt.applyProjection(getProjection(name));
  }
const Projection & _applyProjection ( const Event evt,
const Projection proj 
) const [protected, inherited]

Non-templated version of proj-based applyProjection, to work around header dependency issue.

Definition at line 28 of file ProjectionApplier.cc.

                                                                                      {
    return evt.applyProjection(proj);
  }
const Projection & _declareProjection ( const Projection proj,
const std::string &  name 
) [protected, inherited]

Untemplated function to do the work...

Definition at line 34 of file ProjectionApplier.cc.

                                                                             {
    if (!_allowProjReg) {
      cerr << "Trying to register projection '"
           << proj.name() << "' before init phase in '" << this->name() << "'." << endl;
      exit(2);
    }
    const Projection& reg = getProjHandler().registerProjection(*this, proj, name);
    return reg;
  }
Projection& addPdgIdPair ( PdgId  beam1,
PdgId  beam2 
) [inline, inherited]

Add a colliding beam pair.

Definition at line 108 of file Projection.hh.

                                                       {
      _beamPairs.insert(PdgIdPair(beam1, beam2));
      return *this;
    }
const PROJ& addProjection ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection (user-facing version)

Deprecated:
Use declareProjection() or declare()
Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 157 of file ProjectionApplier.hh.

{ return declareProjection(proj, name); }
const PROJ& apply ( const Event evt,
const Projection proj 
) const [inline, inherited]

Apply the supplied projection on event evt (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 80 of file ProjectionApplier.hh.

{ return applyProjection<PROJ>(evt, proj); }
const PROJ& apply ( const Event evt,
const PROJ &  proj 
) const [inline, inherited]

Apply the supplied projection on event evt (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 92 of file ProjectionApplier.hh.

{ return applyProjection<PROJ>(evt, proj); }
const PROJ& apply ( const Event evt,
const std::string &  name 
) const [inline, inherited]

Apply the supplied projection on event evt (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 104 of file ProjectionApplier.hh.

{ return applyProjection<PROJ>(evt, name); }
const PROJ& applyProjection ( const Event evt,
const Projection proj 
) const [inline, inherited]

Apply the supplied projection on event evt.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 74 of file ProjectionApplier.hh.

                                                                                {
      return pcast<PROJ>(_applyProjection(evt, proj));
    }
const PROJ& applyProjection ( const Event evt,
const PROJ &  proj 
) const [inline, inherited]

Apply the supplied projection on event evt.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 86 of file ProjectionApplier.hh.

                                                                          {
      return pcast<PROJ>(_applyProjection(evt, proj));
    }
const PROJ& applyProjection ( const Event evt,
const std::string &  name 
) const [inline, inherited]

Apply the named projection on event evt.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 98 of file ProjectionApplier.hh.

                                                                               {
      return pcast<PROJ>(_applyProjection(evt, name));
    }
const set< PdgIdPair > beamPairs ( ) const [virtual, inherited]

Return the allowed beam pairs on which this projection can operate, not including recursion. Derived classes should ensure that all contained projections are registered in the _projections set for the beam constraint chaining to work.

Todo:
Remove the beam constraints system from projections.

Definition at line 35 of file Projection.cc.

                                                   {
    set<PdgIdPair> ret = _beamPairs;
    set<ConstProjectionPtr> projs = getProjections();
    for (set<ConstProjectionPtr>::const_iterator ip = projs.begin(); ip != projs.end(); ++ip) {
      ConstProjectionPtr p = *ip;
      getLog() << Log::TRACE << "Proj addr = " << p << endl;
      if (p) ret = intersection(ret, p->beamPairs());
    }
    return ret;
  }
bool before ( const Projection p) const [inherited]

Determine whether this object should be ordered before the object p given as argument. If p is of a different class than this, the before() function of the corresponding type_info objects is used. Otherwise, if the objects are of the same class, the virtual compare(const Projection &) will be returned.

Definition at line 24 of file Projection.cc.

                                                   {
    const std::type_info& thisid = typeid(*this);
    const std::type_info& otherid = typeid(p);
    if (thisid == otherid) {
      return compare(p) < 0;
    } else {
      return thisid.before(otherid);
    }
  }
void calc ( const Jets jets)

Do the calculation directly on a supplied collection of Jet objects.

Todo:
Use Cut for better eta/y selection

< Out of histo range

Definition at line 62 of file JetShape.cc.

                                      {
    clear();

    foreach (const Jet& j, jets) {
      // Apply jet cuts
      const FourMomentum& pj = j.momentum();
      if (!inRange(pj.pT(), _ptcuts)) continue;
      /// @todo Use Cut for better eta/y selection
      if (_rapscheme == PSEUDORAPIDITY && !inRange(fabs(pj.eta()), _rapcuts)) continue;
      if (_rapscheme == RAPIDITY && !inRange(fabs(pj.rapidity()), _rapcuts)) continue;

      // Fill bins
      vector<double> bins(numBins(), 0.0);
      foreach (const Particle& p, j.particles()) {
        const double dR = deltaR(pj, p.momentum(), _rapscheme);
        const int dRindex = binIndex(dR, _binedges);
        if (dRindex == -1) continue; ///< Out of histo range
        bins[dRindex] += p.pT();
      }

      // Add bin vector for this jet to the diffjetshapes container
      _diffjetshapes += bins;
    }

    // Normalize to total pT
    foreach (vector<double>& binsref, _diffjetshapes) {
      double integral = 0.0;
      for (size_t i = 0; i < numBins(); ++i) {
        integral += binsref[i];
      }
      if (integral > 0) {
        for (size_t i = 0; i < numBins(); ++i) {
          binsref[i] /= integral;
        }
      } else {
        // It's just-about conceivable that a jet would have no particles in the given Delta(r) range...
        MSG_DEBUG("No pT contributions in jet Delta(r) range: weird!");
      }
    }

  }
void clear ( )

Reset projection between events.

Definition at line 57 of file JetShape.cc.

                       {
    _diffjetshapes.clear();
  }
virtual unique_ptr<Projection> clone ( ) const [pure virtual, inherited]

Clone on the heap.

Implemented in JetAlg, AxesDefinition, and ParticleFinder.

int compare ( const Projection p) const [protected, virtual]

Compare projections.

Implements Projection.

Definition at line 39 of file JetShape.cc.

                                                 {
    const int jcmp = mkNamedPCmp(p, "Jets");
    if (jcmp != EQUIVALENT) return jcmp;
    const JetShape& other = pcast<JetShape>(p);
    const int ptcmp = cmp(ptMin(), other.ptMin()) || cmp(ptMax(), other.ptMax());
    if (ptcmp != EQUIVALENT) return ptcmp;
    const int rapcmp = cmp(_rapcuts.first, other._rapcuts.first) || cmp(_rapcuts.second, other._rapcuts.second);
    if (rapcmp != EQUIVALENT) return rapcmp;
    int bincmp = cmp(numBins(), other.numBins());
    if (bincmp != EQUIVALENT) return bincmp;
    for (size_t i = 0; i < _binedges.size(); ++i) {
      bincmp = cmp(_binedges[i], other._binedges[i]);
      if (bincmp != EQUIVALENT) return bincmp;
    }
    return EQUIVALENT;
  }
const PROJ& declare ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection (user-facing version)

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 151 of file ProjectionApplier.hh.

{ return declareProjection(proj, name); }
const PROJ& declareProjection ( const PROJ &  proj,
const std::string &  name 
) [inline, protected, inherited]

Register a contained projection.

The type of the argument is used to instantiate a new projection internally: this new object is applied to events rather than the argument object. Hence you are advised to only use locally-scoped Projection objects in your Projection and Analysis constructors, and to avoid polymorphism (e.g. handling ConcreteProjection via a pointer or reference to type Projection) since this will screw up the internal type management.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 142 of file ProjectionApplier.hh.

                                                                           {
      const Projection& reg = _declareProjection(proj, name);
      const PROJ& rtn = dynamic_cast<const PROJ&>(reg);
      return rtn;
    }

Clone on the heap.

double diffJetShape ( size_t  ijet,
size_t  rbin 
) const [inline]

Return value of differential jet shape profile histo bin.

Definition at line 131 of file JetShape.hh.

                                                        {
      assert(inRange(ijet, 0u, numJets()));
      assert(inRange(rbin, 0u, numBins()));
      return _diffjetshapes[ijet][rbin];
    }
const PROJ& get ( const std::string &  name) const [inline, inherited]

Get the named projection, specifying return type via a template argument (user-facing alias).

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 57 of file ProjectionApplier.hh.

{ return getProjection<PROJ>(name); }
Log& getLog ( ) const [inline, inherited]

Get a Log object based on the getName() property of the calling projection object.

Reimplemented from ProjectionApplier.

Definition at line 115 of file Projection.hh.

                        {
      string logname = "Rivet.Projection." + name();
      return Log::getLog(logname);
    }
const PROJ& getProjection ( const std::string &  name) const [inline, inherited]

Get the named projection, specifying return type via a template argument.

Todo:
Add SFINAE to require that PROJ inherit from Projection

Definition at line 50 of file ProjectionApplier.hh.

                                                           {
      const Projection& p = getProjHandler().getProjection(*this, name);
      return pcast<PROJ>(p);
    }
const Projection& getProjection ( const std::string &  name) const [inline, inherited]

Get the named projection (non-templated, so returns as a reference to a Projection base class).

Definition at line 61 of file ProjectionApplier.hh.

                                                                 {
      return getProjHandler().getProjection(*this, name);
    }
std::set<ConstProjectionPtr> getProjections ( ) const [inline, inherited]

Get the contained projections, including recursion.

Definition at line 43 of file ProjectionApplier.hh.

ProjectionHandler& getProjHandler ( ) const [inline, protected, inherited]

Get a reference to the ProjectionHandler for this thread.

Definition at line 122 of file ProjectionApplier.hh.

                                              {
      return _projhandler;
    }
double intJetShape ( size_t  ijet,
size_t  rbin 
) const [inline]

Return value of integrated jet shape profile histo bin.

Definition at line 138 of file JetShape.hh.

                                                       {
      assert(inRange(ijet, 0u, numJets()));
      assert(inRange(rbin, 0u, numBins()));
      double rtn  = 0;
      for (size_t i = 0; i <= rbin; ++i) {
        rtn += _diffjetshapes[ijet][i];
      }
      return rtn;
    }
void markAsOwned ( ) const [inline, inherited]

Mark object as owned by the _projhandler

Todo:
Huh? What's this for?

Definition at line 111 of file ProjectionApplier.hh.

{ _owned = true; }
Cmp< Projection > mkNamedPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]

Shortcut to make a named Cmp<Projection> comparison with the *this object automatically passed as one of the parent projections.

Definition at line 47 of file Projection.cc.

                                                                                                  {
    return pcmp(*this, otherparent, pname);
  }
Cmp< Projection > mkPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]

Shortcut to make a named Cmp<Projection> comparison with the *this object automatically passed as one of the parent projections.

Note:
Alias for mkNamedPCmp

Definition at line 51 of file Projection.cc.

                                                                                             {
    return pcmp(*this, otherparent, pname);
  }
virtual std::string name ( ) const [inline, virtual, inherited]

Get the name of the projection.

Implements ProjectionApplier.

Definition at line 102 of file Projection.hh.

                                   {
      return _name;
    }
size_t numBins ( ) const [inline]

Number of equidistant radius bins.

Definition at line 82 of file JetShape.hh.

                           {
      return _binedges.size() - 1;
    }
size_t numJets ( ) const [inline]

Number of jets which passed cuts.

Definition at line 87 of file JetShape.hh.

                           {
      return _diffjetshapes.size();
    }
void project ( const Event e) [protected, virtual]

Apply the projection to the event.

Todo:
Provide int and diff jet shapes with some sort of area normalisation?

Implements Projection.

Definition at line 105 of file JetShape.cc.

                                       {
    const Jets jets = applyProjection<JetAlg>(e, "Jets").jets(Cuts::ptIn(_ptcuts.first, _ptcuts.second) &
                                                              ((_rapscheme == PSEUDORAPIDITY) ?
                                                               Cuts::etaIn(-_rapcuts.second, _rapcuts.second) :
                                                               Cuts::rapIn(-_rapcuts.second, _rapcuts.second)) );
    calc(jets);
  }
double ptMax ( ) const [inline]

$ p_\perp^\text{max} $ value.

Definition at line 107 of file JetShape.hh.

                         {
      return _ptcuts.second;
    }
double ptMin ( ) const [inline]

$ p_\perp^\text{min} $ value.

Definition at line 102 of file JetShape.hh.

                         {
      return _ptcuts.first;
    }
double rBinMax ( size_t  rbin) const [inline]

Central $ r $ value for bin rbin.

Definition at line 118 of file JetShape.hh.

                                      {
      assert(inRange(rbin, 0u, numBins()));
      return _binedges[rbin+1];
    }
double rBinMid ( size_t  rbin) const [inline]

Central $ r $ value for bin rbin.

Definition at line 124 of file JetShape.hh.

                                      {
      assert(inRange(rbin, 0u, numBins()));
      //cout << _binedges << endl;
      return (_binedges[rbin] + _binedges[rbin+1])/2.0;
    }
double rBinMin ( size_t  rbin) const [inline]

Central $ r $ value for bin rbin.

Definition at line 112 of file JetShape.hh.

                                      {
      assert(inRange(rbin, 0u, numBins()));
      return _binedges[rbin];
    }
double rMax ( ) const [inline]

$ r_\text{max} $ value.

Definition at line 97 of file JetShape.hh.

                        {
      return _binedges.back();
    }
double rMin ( ) const [inline]

$ r_\text{min} $ value.

Definition at line 92 of file JetShape.hh.

                        {
      return _binedges.front();
    }
void setName ( const std::string &  name) [inline, inherited]

Used by derived classes to set their name.

Definition at line 121 of file Projection.hh.

                                        {
      _name = name;
    }

Friends And Related Function Documentation

friend class Cmp< Projection > [friend, inherited]

The Cmp specialization for Projection is a friend.

Definition at line 36 of file Projection.hh.

friend class Event [friend, inherited]

Event is a friend.

Definition at line 33 of file Projection.hh.


Member Data Documentation

bool _allowProjReg [protected, inherited]

Flag to forbid projection registration in analyses until the init phase.

Definition at line 176 of file ProjectionApplier.hh.

vector<double> _binedges [private]

Vector of radius bin edges.

Definition at line 172 of file JetShape.hh.

vector< vector<double> > _diffjetshapes [private]

Jet shape histo -- first index is jet number, second is r bin.

Definition at line 190 of file JetShape.hh.

pair<double, double> _ptcuts [private]

Lower and upper cuts on contributing jet $ p_\perp $.

Definition at line 175 of file JetShape.hh.

pair<double, double> _rapcuts [private]

Lower and upper cuts on contributing jet (pseudo)rapidity.

Definition at line 178 of file JetShape.hh.

RapScheme _rapscheme [private]

Rapidity scheme.

Definition at line 181 of file JetShape.hh.


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