JetShape Class Reference

Calculate the jet shape. More...

#include <JetShape.hh>

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

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.
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.
virtual const Projectionclone () const
 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
 Get the named projection, specifying return type via a template argument.
const ProjectiongetProjection (const std::string &name) const
Projection applying functions

template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const PROJ &proj) const
 Apply the supplied projection on event.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const Projection &proj) const
 Apply the supplied projection on event.
template<typename PROJ >
const PROJ & applyProjection (const Event &evt, const std::string &name) const
 Apply the named projection on event.

Protected Member Functions

void project (const Event &e)
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.
Projection registration functions

template<typename PROJ >
const PROJ & addProjection (const PROJ &proj, const std::string &name)
const Projection_addProjection (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

Calculate the jet shape.

Calculate the differential and integral jet shapes in $P_{\perp}$ for a given set of jets. This particular jet shape projection calculates jet shapes relative to jet centroids, using only the particles associated to each jet, for the hardest $ n $ jets.

The rapidity scheme ($ \eta $ or $ y $) has to be specified when invoking the constructor.

The differential jet shape around a given jet axis at distance interval $ r \pm \delta{r}/2 $ is defined as

\[ \rho(r) = \frac{1}{\delta r} \frac{1}{N_\mathrm{jets}} \sum_\mathrm{jets} \frac{P_\perp(r - \delta r/2, r+\delta r/2)}{p_\perp(0, R)} \]

with $ 0 \le r \le R $ and $ P_\perp(r_1, r_2) = \sum_{\in [r_1, r_2)} p_\perp $.

The integral jet shape around a given jet axes until distance $ r $ is defined as

\[ \Psi(r) = \frac{1}{N_\mathrm{jets}} \sum_\mathrm{jets} \frac{P_\perp(0, r)}{p_\perp(0, R)} \]

with $ 0 \le r \le R $ and $ P_\perp(r_1, r_2) = \sum_{\in [r_1, r_2)} p_\perp $.

The constructor expects also the binning in radius $ r $ to be supplied.

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.

References JetShape::_binedges, JetShape::_ptcuts, JetShape::_rapcuts, ProjectionApplier::addProjection(), Rivet::linspace(), and Projection::setName().

Referenced by JetShape::clone().

00014     : _rapscheme(rapscheme)
00015   {
00016     setName("JetShape");
00017     _binedges = linspace(rmin, rmax, nbins);
00018     _ptcuts = make_pair(ptmin, ptmax);
00019     _rapcuts = make_pair(absrapmin, absrapmax);
00020     addProjection(jetalg, "Jets");
00021   }

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.

References JetShape::_ptcuts, JetShape::_rapcuts, ProjectionApplier::addProjection(), and Projection::setName().

00030     : _binedges(binedges), _rapscheme(rapscheme)
00031   {
00032     setName("JetShape");
00033     _ptcuts = make_pair(ptmin, ptmax);
00034     _rapcuts = make_pair(absrapmin, absrapmax);
00035     addProjection(jetalg, "Jets");
00036   }


Member Function Documentation

const Projection & _addProjection ( const Projection proj,
const std::string &  name 
) [protected, inherited]

Untemplated function to do the work...

Definition at line 33 of file ProjectionApplier.cc.

References ProjectionApplier::_allowProjReg, ProjectionApplier::getProjHandler(), ProjectionApplier::name(), Projection::name(), and ProjectionHandler::registerProjection().

Referenced by ProjectionApplier::addProjection().

00034                                                                              {
00035     if (!_allowProjReg) {
00036       cerr << "Trying to register projection '"
00037            << proj.name() << "' before init phase in '" << this->name() << "'." << endl;
00038       exit(2);
00039     }
00040     const Projection& reg = getProjHandler().registerProjection(*this, proj, name);
00041     return reg;
00042   }

Projection& addPdgIdPair ( PdgId  beam1,
PdgId  beam2 
) [inline, inherited]

Add a colliding beam pair.

Definition at line 107 of file Projection.hh.

References Projection::_beamPairs.

Referenced by Projection::Projection().

00107                                                        {
00108       _beamPairs.insert(PdgIdPair(beam1, beam2));
00109       return *this;
00110     }

const PROJ& addProjection ( 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.

Definition at line 113 of file ProjectionApplier.hh.

References ProjectionApplier::_addProjection().

Referenced by ZFinder::_init(), WFinder::_init(), VetoedFinalState::addVetoOnThisFinalState(), BeamThrust::BeamThrust(), CDF_2009_S8057893::CDF_2009_S8057893::init(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), ClusteredPhotons::ClusteredPhotons(), DISFinalState::DISFinalState(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), FinalState::FinalState(), FoxWolframMoments::FoxWolframMoments(), FParameter::FParameter(), HadronicFinalState::HadronicFinalState(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), ZEUS_2001_S4815815::init(), UA5_1989_S1926373::init(), UA5_1988_S1867512::init(), UA5_1987_S1640666::init(), UA5_1986_S1583476::init(), UA5_1982_S875503::init(), UA1_1990_S2044935::init(), TASSO_1990_S2148048::init(), STAR_2009_UE_HELEN::init(), STAR_2008_S7993412::init(), STAR_2008_S7869363::init(), STAR_2006_S6870392::init(), STAR_2006_S6860818::init(), STAR_2006_S6500200::init(), SFM_1984_S1178091::init(), PDG_HADRON_MULTIPLICITIES_RATIOS::init(), PDG_HADRON_MULTIPLICITIES::init(), OPAL_2004_S6132243::init(), OPAL_2001_S4553896::init(), OPAL_1998_S3780481::init(), OPAL_1993_S2692198::init(), MC_ZZJETS::init(), MC_ZJETS::init(), MC_WWJETS::init(), MC_WPOL::init(), MC_WJETS::init(), MC_VH2BB::init(), MC_TTBAR::init(), MC_SUSY::init(), MC_PHOTONJETUE::init(), MC_PHOTONJETS::init(), MC_LEADINGJETS::init(), MC_JETS::init(), MC_HJETS::init(), MC_GENERIC::init(), MC_DIPHOTON::init(), MC_DIJET::init(), LHCB_2010_S8758301::init(), JADE_OPAL_2000_S4300807::init(), JADE_1998_S3612880::init(), H1_2000_S4129130::init(), H1_1995_S3167097::init(), H1_1994_S2919893::init(), ExampleAnalysis::init(), E735_1998_S3905616::init(), DELPHI_2003_WUD_03_11::init(), DELPHI_2002_069_CONF_603::init(), DELPHI_1996_S3430090::init(), DELPHI_1995_S3137023::init(), D0_2010_S8821313::init(), D0_2010_S8671338::init(), D0_2010_S8570965::init(), D0_2010_S8566488::init(), D0_2009_S8349509::init(), D0_2009_S8320160::init(), D0_2009_S8202443::init(), D0_2008_S7863608::init(), D0_2008_S7837160::init(), D0_2008_S7719523::init(), D0_2008_S7662670::init(), D0_2008_S7554427::init(), D0_2008_S6879055::init(), D0_2007_S7075677::init(), D0_2006_S6438750::init(), D0_2004_S5992206::init(), D0_2001_S4674421::init(), D0_2000_S4480767::init(), D0_1996_S3324664::init(), D0_1996_S3214044::init(), CMS_2011_S9120041::init(), CMS_2011_S9088458::init(), CMS_2011_S9086218::init(), CMS_2011_S8978280::init(), CMS_2011_S8968497::init(), CMS_2011_S8957746::init(), CMS_2011_S8950903::init(), CMS_2011_S8884919::init(), CMS_2010_S8656010::init(), CMS_2010_S8547297::init(), CDF_2010_S8591881_QCD::init(), CDF_2010_S8591881_DY::init(), CDF_2009_S8436959::init(), CDF_2009_S8383952::init(), CDF_2009_S8233977::init(), CDF_2009_NOTE_9936::init(), CDF_2008_S8095620::init(), CDF_2008_S8093652::init(), CDF_2008_S7828950::init(), CDF_2008_S7782535::init(), CDF_2008_S7541902::init(), CDF_2008_S7540469::init(), CDF_2008_NOTE_9351::init(), CDF_2008_LEADINGJETS::init(), CDF_2007_S7057202::init(), CDF_2006_S6653332::init(), CDF_2006_S6450792::init(), CDF_2005_S6217184::init(), CDF_2005_S6080774::init(), CDF_2004_S5839831::init(), CDF_2001_S4751469::init(), CDF_2001_S4563131::init(), CDF_2001_S4517016::init(), CDF_2000_S4266730::init(), CDF_2000_S4155203::init(), CDF_1998_S3618439::init(), CDF_1997_S3541940::init(), CDF_1996_S3418421::init(), CDF_1996_S3349578::init(), CDF_1996_S3108457::init(), CDF_1994_S2952106::init(), CDF_1993_S2742446::init(), CDF_1990_S2089246::init(), CDF_1988_S1865951::init(), BELLE_2006_S6265367::init(), ATLAS_2011_S9225137::init(), ATLAS_2011_S9212183::init(), ATLAS_2011_S9131140::init(), ATLAS_2011_S9128077::init(), ATLAS_2011_S9126244::init(), ATLAS_2011_S9120807::init(), ATLAS_2011_S9108483::init(), ATLAS_2011_S9041966::init(), ATLAS_2011_S9019561::init(), ATLAS_2011_S9002537::init(), ATLAS_2011_S8994773::init(), ATLAS_2011_S8983313::init(), ATLAS_2011_S8971293::init(), ATLAS_2011_S8924791::init(), ATLAS_2011_I925932::init(), ATLAS_2011_I919017::init(), ATLAS_2011_CONF_2011_098::init(), ATLAS_2011_CONF_2011_090::init(), ATLAS_2010_S8919674::init(), ATLAS_2010_S8918562::init(), ATLAS_2010_S8914702::init(), ATLAS_2010_S8894728::init(), ATLAS_2010_S8817804::init(), ATLAS_2010_S8591806::init(), ATLAS_2010_CONF_2010_049::init(), ALICE_2011_S8945144::init(), ALICE_2011_S8909580::init(), ALICE_2010_S8706239::init(), ALICE_2010_S8625980::init(), ALICE_2010_S8624100::init(), ALEPH_2004_S5765862::init(), ALEPH_1996_S3486095::init(), ALEPH_1996_S3196992::init(), ALEPH_1991_S2435284::init(), IsolationProjection< PROJ1, PROJ2, EST >::IsolationProjection(), JetAlg::JetAlg(), JetShape::JetShape(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LeptonClusters::LeptonClusters(), LossyFinalState< ConstRandomFilter >::LossyFinalState(), MergedFinalState::MergedFinalState(), MissingMomentum::MissingMomentum(), Multiplicity::Multiplicity(), NeutralFinalState::NeutralFinalState(), NonHadronicFinalState::NonHadronicFinalState(), ParisiTensor::ParisiTensor(), Sphericity::Sphericity(), Spherocity::Spherocity(), SVertex::SVertex(), Thrust::Thrust(), TotalVisibleMomentum::TotalVisibleMomentum(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerCDFRun2::TriggerCDFRun2(), TriggerUA5::TriggerUA5(), VetoedFinalState::VetoedFinalState(), and VisibleFinalState::VisibleFinalState().

00113                                                                        {
00114       const Projection& reg = _addProjection(proj, name);
00115       const PROJ& rtn = dynamic_cast<const PROJ&>(reg);
00116       return rtn;
00117     }

const PROJ& applyProjection ( const Event evt,
const std::string &  name 
) const [inline, inherited]

Apply the named projection on event.

Definition at line 81 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

00081                                                                                {
00082       return pcast<PROJ>(_applyProjection(evt, name));
00083     }

const PROJ& applyProjection ( const Event evt,
const Projection proj 
) const [inline, inherited]

Apply the supplied projection on event.

Definition at line 74 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

00074                                                                                 {
00075       return pcast<PROJ>(_applyProjection(evt, proj));
00076     }

const PROJ& applyProjection ( const Event evt,
const PROJ &  proj 
) const [inline, inherited]

Apply the supplied projection on event.

Definition at line 67 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

Referenced by DISFinalState::project().

00067                                                                           {
00068       return pcast<PROJ>(_applyProjection(evt, proj));
00069     }

const set< PdgIdPair > beamPairs (  )  const [virtual, inherited]

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

Definition at line 39 of file Projection.cc.

References Projection::_beamPairs, Projection::beamPairs(), Projection::getLog(), ProjectionApplier::getProjections(), Rivet::intersection(), and Log::TRACE.

Referenced by Projection::beamPairs().

00039                                                    {
00040     set<PdgIdPair> ret = _beamPairs;
00041     set<ConstProjectionPtr> projs = getProjections();
00042     for (set<ConstProjectionPtr>::const_iterator ip = projs.begin(); ip != projs.end(); ++ip) {
00043       ConstProjectionPtr p = *ip;
00044       getLog() << Log::TRACE << "Proj addr = " << p << endl;
00045       if (p) ret = intersection(ret, p->beamPairs());
00046     }
00047     return ret;
00048   }

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 28 of file Projection.cc.

References Projection::compare().

Referenced by less< const Rivet::Projection * >::operator()().

00028                                                    {
00029     const std::type_info& thisid = typeid(*this);
00030     const std::type_info& otherid = typeid(p);
00031     if (thisid == otherid) {
00032       return compare(p) < 0;
00033     } else {
00034       return thisid.before(otherid);
00035     }
00036   }

void calc ( const Jets jets  ) 

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

Todo:
Introduce a better (i.e. more safe and general) eta/y selection mechanism: MomentumFilter

Definition at line 62 of file JetShape.cc.

References JetShape::_binedges, JetShape::_diffjetshapes, JetShape::_ptcuts, JetShape::_rapcuts, JetShape::_rapscheme, JetShape::clear(), Rivet::deltaR(), FourVector::eta(), Rivet::index_between(), Rivet::inRange(), Rivet::integral(), Particle::momentum(), Jet::momentum(), MSG_DEBUG, JetShape::numBins(), Jet::particles(), Rivet::PSEUDORAPIDITY, FourMomentum::pT(), FourMomentum::rapidity(), and Rivet::RAPIDITY.

Referenced by CDF_2008_S7782535::analyze(), and JetShape::project().

00062                                       {
00063     clear();
00064 
00065     foreach (const Jet& j, jets) {
00066       // Apply jet cuts
00067       FourMomentum pj = j.momentum();
00068       if (!inRange(pj.pT(), _ptcuts)) continue;
00069       /// @todo Introduce a better (i.e. more safe and general) eta/y selection mechanism: MomentumFilter
00070       if (_rapscheme == PSEUDORAPIDITY && !inRange(fabs(pj.eta()), _rapcuts)) continue;
00071       if (_rapscheme == RAPIDITY && !inRange(fabs(pj.rapidity()), _rapcuts)) continue;
00072 
00073       // Fill bins
00074       vector<double> bins(numBins(), 0.0);
00075       foreach (const Particle& p, j.particles()) {
00076         const double dR = deltaR(pj, p.momentum(), _rapscheme);
00077         const int dRindex = index_between(dR, _binedges);
00078         if (dRindex == -1) continue; //< Out of histo range
00079         bins[dRindex] += p.momentum().pT();
00080       }
00081 
00082       // Add bin vector for this jet to the diffjetshapes container
00083       _diffjetshapes += bins;
00084     }
00085 
00086     // Normalize to total pT
00087     foreach (vector<double>& binsref, _diffjetshapes) {
00088       double integral = 0.0;
00089       for (size_t i = 0; i < numBins(); ++i) {
00090         integral += binsref[i];
00091       }
00092       if (integral > 0) {
00093         for (size_t i = 0; i < numBins(); ++i) {
00094           binsref[i] /= integral;
00095         }
00096       } else {
00097         // It's just-about conceivable that a jet would have no particles in the given Delta(r) range...
00098         MSG_DEBUG("No pT contributions in jet Delta(r) range: weird!");
00099       }
00100     }
00101 
00102   }

void clear (  ) 

Reset projection between events.

Definition at line 57 of file JetShape.cc.

References JetShape::_diffjetshapes.

Referenced by JetShape::calc().

00057                        {
00058     _diffjetshapes.clear();
00059   }

virtual const Projection* clone (  )  const [inline, virtual]

Clone on the heap.

Implements Projection.

Definition at line 65 of file JetShape.hh.

References JetShape::JetShape().

00065                                             {
00066       return new JetShape(*this);
00067     }

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

Compare projections.

Implements Projection.

Definition at line 39 of file JetShape.cc.

References JetShape::_binedges, JetShape::_rapcuts, Rivet::cmp(), Rivet::EQUIVALENT, Projection::mkNamedPCmp(), JetShape::numBins(), JetShape::ptMax(), and JetShape::ptMin().

00039                                                  {
00040     const int jcmp = mkNamedPCmp(p, "Jets");
00041     if (jcmp != EQUIVALENT) return jcmp;
00042     const JetShape& other = pcast<JetShape>(p);
00043     const int ptcmp = cmp(ptMin(), other.ptMin()) || cmp(ptMax(), other.ptMax());
00044     if (ptcmp != EQUIVALENT) return ptcmp;
00045     const int rapcmp = cmp(_rapcuts.first, other._rapcuts.first) || cmp(_rapcuts.second, other._rapcuts.second);
00046     if (rapcmp != EQUIVALENT) return rapcmp;
00047     int bincmp = cmp(numBins(), other.numBins());
00048     if (bincmp != EQUIVALENT) return bincmp;
00049     for (size_t i = 0; i < _binedges.size(); ++i) {
00050       bincmp = cmp(_binedges[i], other._binedges[i]);
00051       if (bincmp != EQUIVALENT) return bincmp;
00052     }
00053     return EQUIVALENT;
00054   }

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

Return value of differential jet shape profile histo bin.

Definition at line 133 of file JetShape.hh.

References JetShape::_diffjetshapes, Rivet::inRange(), JetShape::numBins(), and JetShape::numJets().

Referenced by CDF_2005_S6217184::analyze(), and ATLAS_2011_S8924791::analyze().

00133                                                         {
00134       assert(inRange(ijet, 0, numJets()));
00135       assert(inRange(rbin, 0, numBins()));
00136       return _diffjetshapes[ijet][rbin];
00137     }

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 114 of file Projection.hh.

References Projection::name().

Referenced by Projection::beamPairs(), InvMassFinalState::calc(), VetoedFinalState::project(), UnstableFinalState::project(), LossyFinalState< ConstRandomFilter >::project(), IsolationProjection< PROJ1, PROJ2, EST >::project(), InitialQuarks::project(), ChargedFinalState::project(), and TotalVisibleMomentum::TotalVisibleMomentum().

00114                         {
00115       string logname = "Rivet.Projection." + name();
00116       return Log::getLog(logname);
00117     }

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 57 of file ProjectionApplier.hh.

References ProjectionHandler::getProjection(), and ProjectionApplier::getProjHandler().

00057                                                                  {
00058       return getProjHandler().getProjection(*this, name);
00059     }

const PROJ& getProjection ( const std::string &  name  )  const [inline, inherited]

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

Definition at line 49 of file ProjectionApplier.hh.

References ProjectionHandler::getProjection(), and ProjectionApplier::getProjHandler().

Referenced by ProjectionApplier::_applyProjection(), Rivet::pcmp(), and Hemispheres::project().

00049                                                            {
00050       const Projection& p = getProjHandler().getProjection(*this, name);
00051       return pcast<PROJ>(p);
00052     }

std::set<ConstProjectionPtr> getProjections (  )  const [inline, inherited]

Get the contained projections, including recursion.

Definition at line 42 of file ProjectionApplier.hh.

References ProjectionHandler::DEEP, ProjectionHandler::getChildProjections(), and ProjectionApplier::getProjHandler().

Referenced by Projection::beamPairs().

00042                                                       {
00043       return getProjHandler().getChildProjections(*this, ProjectionHandler::DEEP);
00044     }

ProjectionHandler& getProjHandler (  )  const [inline, protected, inherited]
double intJetShape ( size_t  ijet,
size_t  rbin 
) const [inline]

Return value of integrated jet shape profile histo bin.

Definition at line 140 of file JetShape.hh.

References JetShape::_diffjetshapes, Rivet::inRange(), JetShape::numBins(), and JetShape::numJets().

Referenced by CDF_2008_S7782535::analyze(), CDF_2005_S6217184::analyze(), and ATLAS_2011_S8924791::analyze().

00140                                                        {
00141       assert(inRange(ijet, 0, numJets()));
00142       assert(inRange(rbin, 0, numBins()));
00143       double rtn  = 0;
00144       for (size_t i = 0; i <= rbin; ++i) {
00145         rtn += _diffjetshapes[ijet][i];
00146       }
00147       return rtn;
00148     }

Cmp< Projection > mkNamedPCmp ( const Projection otherparent,
const std::string &  pname 
) const [protected, inherited]
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.

Definition at line 57 of file Projection.cc.

References Rivet::pcmp().

00058                                                                 {
00059     return pcmp(*this, otherparent, pname);
00060   }

virtual std::string name (  )  const [inline, virtual, inherited]
size_t numBins (  )  const [inline]
size_t numJets (  )  const [inline]

Number of jets which passed cuts.

Definition at line 89 of file JetShape.hh.

References JetShape::_diffjetshapes.

Referenced by CDF_2008_S7782535::analyze(), CDF_2005_S6217184::analyze(), ATLAS_2011_S8924791::analyze(), JetShape::diffJetShape(), and JetShape::intJetShape().

00089                            {
00090       return _diffjetshapes.size();
00091     }

void project ( const Event e  )  [protected, virtual]
Todo:
Provide int and diff jet shapes with some sort of area normalisation?

Apply the projection to the event.

Implements Projection.

Definition at line 105 of file JetShape.cc.

References JetShape::_ptcuts, JetShape::_rapcuts, JetShape::_rapscheme, and JetShape::calc().

00105                                        {
00106     const Jets jets = applyProjection<JetAlg>(e, "Jets").jets(_ptcuts.first, _ptcuts.second,
00107                                                               -_rapcuts.second, _rapcuts.second, _rapscheme);
00108     calc(jets);
00109   }

double ptMax (  )  const [inline]

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

Definition at line 109 of file JetShape.hh.

References JetShape::_ptcuts.

Referenced by JetShape::compare().

00109                          {
00110       return _ptcuts.second;
00111     }

double ptMin (  )  const [inline]

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

Definition at line 104 of file JetShape.hh.

References JetShape::_ptcuts.

Referenced by JetShape::compare().

00104                          {
00105       return _ptcuts.first;
00106     }

double rBinMax ( size_t  rbin  )  const [inline]

Central $ r $ value for bin rbin.

Definition at line 120 of file JetShape.hh.

References JetShape::_binedges, Rivet::inRange(), and JetShape::numBins().

Referenced by CDF_2008_S7782535::analyze(), and CDF_2005_S6217184::analyze().

00120                                       {
00121       assert(inRange(rbin, 0, numBins()));
00122       return _binedges[rbin+1];
00123     }

double rBinMid ( size_t  rbin  )  const [inline]

Central $ r $ value for bin rbin.

Definition at line 126 of file JetShape.hh.

References JetShape::_binedges, Rivet::inRange(), and JetShape::numBins().

Referenced by CDF_2005_S6217184::analyze(), and ATLAS_2011_S8924791::analyze().

00126                                       {
00127       assert(inRange(rbin, 0, numBins()));
00128       //cout << _binedges << endl;
00129       return (_binedges[rbin] + _binedges[rbin+1])/2.0;
00130     }

double rBinMin ( size_t  rbin  )  const [inline]

Central $ r $ value for bin rbin.

Definition at line 114 of file JetShape.hh.

References JetShape::_binedges, Rivet::inRange(), and JetShape::numBins().

00114                                       {
00115       assert(inRange(rbin, 0, numBins()));
00116       return _binedges[rbin];
00117     }

double rMax (  )  const [inline]

$ r_\text{max} $ value.

Definition at line 99 of file JetShape.hh.

References JetShape::_binedges.

00099                         {
00100       return _binedges.back();
00101     }

double rMin (  )  const [inline]

$ r_\text{min} $ value.

Definition at line 94 of file JetShape.hh.

References JetShape::_binedges.

00094                         {
00095       return _binedges.front();
00096     }

void setName ( const std::string &  name  )  [inline, inherited]

Used by derived classes to set their name.

Definition at line 120 of file Projection.hh.

References Projection::_name.

Referenced by ZFinder::_init(), WFinder::_init(), FastJets::_init1(), FastJets::_init2(), FastJets::_init3(), Beam::Beam(), BeamThrust::BeamThrust(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), ClusteredPhotons::ClusteredPhotons(), ConstLossyFinalState::ConstLossyFinalState(), DISFinalState::DISFinalState(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), FinalState::FinalState(), FoxWolframMoments::FoxWolframMoments(), FParameter::FParameter(), HadronicFinalState::HadronicFinalState(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), InitialQuarks::InitialQuarks(), IsolationProjection< PROJ1, PROJ2, EST >::IsolationProjection(), JetAlg::JetAlg(), JetShape::JetShape(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LeptonClusters::LeptonClusters(), LossyFinalState< ConstRandomFilter >::LossyFinalState(), MergedFinalState::MergedFinalState(), MissingMomentum::MissingMomentum(), Multiplicity::Multiplicity(), NeutralFinalState::NeutralFinalState(), NonHadronicFinalState::NonHadronicFinalState(), ParisiTensor::ParisiTensor(), PVertex::PVertex(), Sphericity::Sphericity(), Spherocity::Spherocity(), SVertex::SVertex(), Thrust::Thrust(), TotalVisibleMomentum::TotalVisibleMomentum(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerCDFRun2::TriggerCDFRun2(), TriggerUA5::TriggerUA5(), UnstableFinalState::UnstableFinalState(), VetoedFinalState::VetoedFinalState(), and VisibleFinalState::VisibleFinalState().

00120                                         {
00121       _name = name;
00122     }


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 140 of file ProjectionApplier.hh.

Referenced by ProjectionApplier::_addProjection(), and Analysis::Analysis().

vector<double> _binedges [private]
vector< vector<double> > _diffjetshapes [private]

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

Definition at line 192 of file JetShape.hh.

Referenced by JetShape::calc(), JetShape::clear(), JetShape::diffJetShape(), JetShape::intJetShape(), and JetShape::numJets().

pair<double, double> _ptcuts [private]

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

Definition at line 177 of file JetShape.hh.

Referenced by JetShape::calc(), JetShape::JetShape(), JetShape::project(), JetShape::ptMax(), and JetShape::ptMin().

pair<double, double> _rapcuts [private]

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

Definition at line 180 of file JetShape.hh.

Referenced by JetShape::calc(), JetShape::compare(), JetShape::JetShape(), and JetShape::project().

RapScheme _rapscheme [private]

Rapidity scheme.

Definition at line 183 of file JetShape.hh.

Referenced by JetShape::calc(), and JetShape::project().


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