ZFinder Class Reference

#include <ZFinder.hh>

Inheritance diagram for ZFinder:

Inheritance graph
[legend]

Collaboration diagram for ZFinder:

Collaboration graph
[legend]

List of all members.


Detailed Description

Chain together different projections as convenience for finding Z's from two leptons in the final state

Definition at line 17 of file ZFinder.hh.


Public Types

typedef Particle entity_type
typedef ParticleVector collection_type

Public Member Functions

const FinalStateremainingFinalState () const
const FinalStateconstituentsFinalState () const
virtual const ParticleVectorparticles () const
 Get the final-state particles.
template<typename F>
const ParticleVectorparticles (F sorter) const
 Get the final-state particles, ordered by supplied sorting function object.
const ParticleVectorparticlesByPt () const
 Get the final-state particles, ordered by $ p_T $.
const ParticleVectorparticlesByE () const
 Get the final-state particles, ordered by $ E $.
const ParticleVectorparticlesByEt () const
 Get the final-state particles, ordered by $ E_T $.
virtual size_t size () const
 Access the projected final-state particles.
virtual bool empty () const
 Is this final state empty?
virtual bool isEmpty () const
const collection_typeentities () const
 Template-usable interface common to JetAlg.
bool before (const Projection &p) const
virtual const std::set< BeamPairbeamPairs () const
virtual std::string name () const
 Get the name of the projection.
ProjectionaddBeamPair (const ParticleName &beam1, const ParticleName &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
 ZFinder (const FinalState &fs, PdgId pid, double m2_min, double m2_max, double dRmax)
 ZFinder (double etaMin, double etaMax, double pTmin, PdgId pid, double m2_min, double m2_max, double dRmax)
 ZFinder (const std::vector< std::pair< double, double > > &etaRanges, double pTmin, PdgId pid, double m2_min, const double m2_max, double dRmax)
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)
 Apply the projection on the supplied event.
int compare (const Projection &p) const
 Compare projections.
bool accept (const Particle &p) const
 Decide if a particle is to be accepted or not.
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

vector< pair< double, double > > _etaRanges
 The ranges allowed for pseudorapidity.
double _ptmin
 The minimum allowed transverse momentum.
ParticleVector _theParticles
 The final-state particles.
bool _allowProjReg
 Flag to forbid projection registration in analyses until the init phase.

Private Member Functions

void _init (const std::vector< std::pair< double, double > > &etaRanges, double pTmin, PdgId pid, double m2_min, double m2_max, double dRmax)
 Common implementation of constructor operation, taking FS params.
void _init (const FinalState &fs, PdgId pid, double m2_min, double m2_max, double dRmax)
 Common implementation of constructor operation, taking FS.

Friends

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

Member Typedef Documentation

typedef Particle entity_type [inherited]

Definition at line 75 of file FinalState.hh.

typedef ParticleVector collection_type [inherited]

Definition at line 76 of file FinalState.hh.


Constructor & Destructor Documentation

ZFinder ( const FinalState fs,
PdgId  pid,
double  m2_min,
double  m2_max,
double  dRmax 
)

Constructor taking a FinalState and type of the leptons, mass window, and maximum dR of photons around leptons to take into account for Z reconstruction.

Definition at line 13 of file ZFinder.cc.

References ZFinder::_init().

Referenced by ZFinder::clone().

00016                                  {
00017     _init(fs, pid, m2_min, m2_max, dRmax);
00018   }

ZFinder ( double  etaMin,
double  etaMax,
double  pTmin,
PdgId  pid,
double  m2_min,
double  m2_max,
double  dRmax 
)

Constructor taking single eta/pT bounds and type of the leptons, mass window, and maximum dR of photons around leptons to take into account for Z reconstruction.

Definition at line 21 of file ZFinder.cc.

References ZFinder::_init().

00025                                  {
00026     vector<pair<double, double> > etaRanges;
00027     etaRanges += std::make_pair(etaMin, etaMax);
00028     _init(etaRanges, pTmin, pid, m2_min, m2_max, dRmax);
00029   }

ZFinder ( const std::vector< std::pair< double, double > > &  etaRanges,
double  pTmin,
PdgId  pid,
double  m2_min,
const double  m2_max,
double  dRmax 
)

Constructor taking multiple eta/pT bounds and type of the leptons, mass window, and maximum dR of photons around leptons to take into account for Z reconstruction.

Definition at line 32 of file ZFinder.cc.

References ZFinder::_init().

00036                                  {
00037     _init(etaRanges, pTmin, pid, m2_min, m2_max, dRmax);
00038   }


Member Function Documentation

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

Clone on the heap.

Reimplemented from FinalState.

Definition at line 54 of file ZFinder.hh.

References ZFinder::ZFinder().

00054                                             {
00055       return new ZFinder(*this);
00056     }

const FinalState & remainingFinalState (  )  const

Access to the remaining particles, after the Z and clustered photons have been removed from the full final state (e.g. for running a jet finder on it)

Definition at line 75 of file ZFinder.cc.

Referenced by MC_ZJETS::init(), D0_2009_S8349509::init(), D0_2009_S8202443::init(), and D0_2008_S7863608::init().

00076   {
00077     return getProjection<FinalState>("RFS");
00078   }

const FinalState & constituentsFinalState (  )  const

Access to the Z constituent leptons final state (e.g. for more fine-grained cuts on the leptons)

Definition at line 81 of file ZFinder.cc.

Referenced by MC_ZJETS::analyze(), D0_2009_S8349509::analyze(), D0_2008_S7863608::analyze(), and D0_2007_S7075677::analyze().

00082   {
00083     return getProjection<FinalState>("IMFS");
00084   }

void project ( const Event e  )  [protected, virtual]

Apply the projection on the supplied event.

Reimplemented from FinalState.

Definition at line 97 of file ZFinder.cc.

References FinalState::_theParticles, Log::DEBUG, Projection::getLog(), Particle::momentum(), Projection::name(), FinalState::particles(), Particle::pdgId(), Particle::setMomentum(), Rivet::PID::threeCharge(), and Rivet::PID::Z().

00097                                       {
00098     _theParticles.clear();
00099 
00100     const FinalState& imfs=applyProjection<FinalState>(e, "IMFS");
00101     if (imfs.particles().size() != 2) return;
00102     FourMomentum pZ = imfs.particles()[0].momentum() + imfs.particles()[1].momentum();
00103     const int z3charge = PID::threeCharge(imfs.particles()[0].pdgId()) + PID::threeCharge(imfs.particles()[1].pdgId());
00104     assert(z3charge == 0);
00105 
00106     stringstream msg;
00107     msg << "Z reconstructed from: " << endl
00108         << "   " << imfs.particles()[0].momentum() << " " << imfs.particles()[0].pdgId() << endl
00109         << " + " << imfs.particles()[1].momentum() << " " << imfs.particles()[1].pdgId() << endl;
00110 
00111     // Add in clustered photons
00112     const FinalState& photons = applyProjection<FinalState>(e, "CPhotons");
00113     foreach (const Particle& photon, photons.particles()) {
00114       msg << " + " << photon.momentum() << " " << photon.pdgId() << endl;
00115       pZ += photon.momentum();
00116     }
00117     msg << " = " << pZ;
00118     getLog() << Log::DEBUG << msg.str() << endl;
00119 
00120     Particle Z;
00121     Z.setMomentum(pZ);
00122     _theParticles.push_back(Z);
00123     getLog() << Log::DEBUG << name() << " found " << _theParticles.size()
00124              << " particles." << endl;
00125   }

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

Compare projections.

Reimplemented from FinalState.

Definition at line 86 of file ZFinder.cc.

References Rivet::cmp(), Rivet::EQUIVALENT, and Projection::mkNamedPCmp().

00086                                                 {
00087     PCmp cmp = mkNamedPCmp(p, "IMFS");
00088     if (cmp != EQUIVALENT) return cmp;
00089 
00090     cmp = mkNamedPCmp(p, "CPhotons");
00091     if (cmp != EQUIVALENT) return cmp;
00092 
00093     return EQUIVALENT;
00094   }

void _init ( const std::vector< std::pair< double, double > > &  etaRanges,
double  pTmin,
PdgId  pid,
double  m2_min,
double  m2_max,
double  dRmax 
) [private]

Common implementation of constructor operation, taking FS params.

Definition at line 41 of file ZFinder.cc.

Referenced by ZFinder::ZFinder().

00044                                     {
00045     FinalState fs(etaRanges, pTmin);
00046     _init(fs, pid, m2_min, m2_max, dRmax);
00047   }

void _init ( const FinalState fs,
PdgId  pid,
double  m2_min,
double  m2_max,
double  dRmax 
) [private]

Common implementation of constructor operation, taking FS.

Definition at line 50 of file ZFinder.cc.

References ProjectionApplier::addProjection(), VetoedFinalState::addVetoOnThisFinalState(), FinalState::FinalState(), and Projection::setName().

00054   {
00055     setName("ZFinder");
00056 
00057     addProjection(fs, "FS");
00058 
00059     InvMassFinalState imfs(fs, std::make_pair(pid, -pid), m2_min, m2_max);
00060     addProjection(imfs, "IMFS");
00061  
00062     ClusteredPhotons cphotons(FinalState(), imfs, dRmax);
00063     addProjection(cphotons, "CPhotons");
00064 
00065     VetoedFinalState remainingFS;
00066     remainingFS.addVetoOnThisFinalState(imfs);
00067     remainingFS.addVetoOnThisFinalState(cphotons);
00068     addProjection(remainingFS, "RFS");
00069   }

virtual const ParticleVector& particles (  )  const [inline, virtual, inherited]

Get the final-state particles.

Definition at line 39 of file FinalState.hh.

References FinalState::_theParticles.

Referenced by UA5_1986_S1583476::analyze(), UA5_1982_S875503::analyze(), UA1_1990_S2044935::analyze(), STAR_2009_UE_HELEN::analyze(), STAR_2008_S7993412::analyze(), STAR_2006_S6500200::analyze(), SFM_1984_S1178091::analyze(), PDG_HADRON_MULTIPLICITIES_RATIOS::analyze(), PDG_HADRON_MULTIPLICITIES::analyze(), OPAL_1998_S3780481::analyze(), MC_ZJETS::analyze(), MC_WJETS::analyze(), MC_TTBAR::analyze(), MC_SUSY::analyze(), MC_PHOTONJETUE::analyze(), MC_PHOTONJETS::analyze(), MC_LEADINGJETS::analyze(), MC_DIJET::analyze(), H1_2000_S4129130::analyze(), H1_1995_S3167097::analyze(), H1_1994_S2919893::analyze(), E735_1998_S3905616::analyze(), DELPHI_2002_069_CONF_603::analyze(), DELPHI_1995_S3137023::analyze(), D0_2009_S8349509::analyze(), D0_2009_S8202443::analyze(), D0_2008_S7863608::analyze(), D0_2008_S7837160::analyze(), D0_2008_S7719523::analyze(), D0_2008_S7554427::analyze(), D0_2008_S6879055::analyze(), D0_2007_S7075677::analyze(), D0_2006_S6438750::analyze(), D0_2001_S4674421::analyze(), D0_1998_S3711838::analyze(), CDF_2009_S8383952::analyze(), CDF_2009_S8233977::analyze(), CDF_2008_S8095620::analyze(), CDF_2008_S7541902::analyze(), CDF_2008_S7540469::analyze(), CDF_2008_NOTE_9351::analyze(), CDF_2008_LEADINGJETS::analyze(), CDF_2006_S6653332::analyze(), CDF_2002_S4796047::analyze(), CDF_2000_S4155203::analyze(), CDF_1991_S2313472::analyze(), CDF_1990_S2089246::analyze(), CDF_1988_S1865951::analyze(), ATLAS_2010_S8591806::analyze(), ALEPH_1996_S3486095::analyze(), Thrust::calc(), Sphericity::calc(), FinalState::entities(), FinalState::particlesByE(), FinalState::particlesByEt(), FinalState::particlesByPt(), ZFinder::project(), WFinder::project(), VetoedFinalState::project(), TriggerUA5::project(), TriggerCDFRun0Run1::project(), TotalVisibleMomentum::project(), SVertex::project(), NeutralFinalState::project(), Multiplicity::project(), MergedFinalState::project(), LossyFinalState::project(), LeadingParticlesFinalState::project(), KtJets::project(), JetShape::project(), InvMassFinalState::project(), IdentifiedFinalState::project(), Hemispheres::project(), HadronicFinalState::project(), FoxWolframMoments::project(), FinalStateHCM::project(), FinalState::project(), FastJets::project(), DISLepton::project(), ClusteredPhotons::project(), ChargedLeptons::project(), ChargedFinalState::project(), and CentralEtHCM::project().

00039 { return _theParticles; }

const ParticleVector& particles ( sorter  )  const [inline, inherited]

Get the final-state particles, ordered by supplied sorting function object.

Definition at line 43 of file FinalState.hh.

References FinalState::_theParticles.

00043                                                     {
00044       std::sort(_theParticles.begin(), _theParticles.end(), sorter);
00045       return _theParticles;
00046     }

const ParticleVector& particlesByPt (  )  const [inline, inherited]

Get the final-state particles, ordered by $ p_T $.

Definition at line 49 of file FinalState.hh.

References Rivet::cmpParticleByPt(), and FinalState::particles().

Referenced by MC_PHOTONJETUE::analyze().

00049                                                 {
00050       return particles(cmpParticleByPt);
00051     }

const ParticleVector& particlesByE (  )  const [inline, inherited]

Get the final-state particles, ordered by $ E $.

Definition at line 54 of file FinalState.hh.

References Rivet::cmpParticleByE(), and FinalState::particles().

00054                                                {
00055       return particles(cmpParticleByE);
00056     }

const ParticleVector& particlesByEt (  )  const [inline, inherited]

Get the final-state particles, ordered by $ E_T $.

Definition at line 59 of file FinalState.hh.

References Rivet::cmpParticleByEt(), and FinalState::particles().

00059                                                 {
00060       return particles(cmpParticleByEt);
00061     }

virtual size_t size (  )  const [inline, virtual, inherited]

virtual bool empty (  )  const [inline, virtual, inherited]

virtual bool isEmpty (  )  const [inline, virtual, inherited]

Deprecated:
Is this final state empty?

Definition at line 70 of file FinalState.hh.

References FinalState::_theParticles.

00070 { return _theParticles.empty(); }

const collection_type& entities (  )  const [inline, inherited]

Template-usable interface common to JetAlg.

Definition at line 79 of file FinalState.hh.

References FinalState::particles().

00079                                             {
00080       return particles();
00081     }

bool accept ( const Particle p  )  const [protected, inherited]

Decide if a particle is to be accepted or not.

Definition at line 99 of file FinalState.cc.

References FinalState::_etaRanges, FinalState::_ptmin, FourVector::eta(), Rivet::eta(), Particle::genParticle(), Particle::momentum(), FourMomentum::pT(), and Rivet::pT().

Referenced by LeadingParticlesFinalState::project(), InvMassFinalState::project(), and FinalState::project().

00099                                                  {
00100     // Not having s.c. == 1 should never happen!
00101     assert(p.genParticle().status() == 1);
00102 
00103     // Check pT cut
00104     if (_ptmin > 0.0) {
00105       const double pT = p.momentum().pT();
00106       if (pT < _ptmin) return false;
00107     }
00108 
00109     // Check eta cuts
00110     if (!_etaRanges.empty()) {
00111       bool eta_pass = false;
00112       const double eta = p.momentum().eta();
00113       typedef pair<double,double> EtaPair;
00114       foreach (const EtaPair& etacuts, _etaRanges) {
00115         if (eta > etacuts.first && eta < etacuts.second) {
00116           eta_pass = true;
00117           break;
00118         }
00119       }
00120       if (!eta_pass) return false;
00121     }
00122  
00123     return true;
00124   }

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   }

const set< BeamPair > 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<BeamPair> 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   }

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

Projection& addBeamPair ( const ParticleName beam1,
const ParticleName beam2 
) [inline, inherited]

Add a colliding beam pair.

Definition at line 105 of file Projection.hh.

References Projection::_beamPairs.

Referenced by Projection::Projection().

00105                                                                                   {
00106       _beamPairs.insert(BeamPair(beam1, beam2));
00107       return *this;
00108     }

Log& getLog (  )  const [inline, inherited]

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

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.

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

Get the contained projections, including recursion.

Definition at line 43 of file ProjectionApplier.hh.

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

Referenced by Projection::beamPairs().

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

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

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

Definition at line 50 of file ProjectionApplier.hh.

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

Referenced by VetoedFinalState::compare(), Rivet::pcmp(), and Hemispheres::project().

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

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

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

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

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

Apply the supplied projection on event.

Definition at line 68 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

Referenced by HadronicFinalState::project(), and FinalStateHCM::project().

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

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

Apply the supplied projection on event.

Definition at line 75 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

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

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

Apply the named projection on event.

Definition at line 82 of file ProjectionApplier.hh.

References ProjectionApplier::_applyProjection().

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

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

Get a reference to the ProjectionHandler for this thread.

Definition at line 95 of file ProjectionApplier.hh.

References ProjectionApplier::_projhandler.

Referenced by ProjectionApplier::_addProjection(), ProjectionApplier::getProjection(), ProjectionApplier::getProjections(), and ProjectionApplier::~ProjectionApplier().

00095                                               {
00096       assert(_projhandler);
00097       return *_projhandler;
00098     }

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

References ProjectionApplier::_addProjection().

Referenced by ZFinder::_init(), WFinder::_init(), VetoedFinalState::addVetoOnThisFinalState(), CDF_2009_S8057893::CDF_2009_S8057893::init(), CentralEtHCM::CentralEtHCM(), ChargedFinalState::ChargedFinalState(), ChargedLeptons::ChargedLeptons(), ClusteredPhotons::ClusteredPhotons(), DISKinematics::DISKinematics(), DISLepton::DISLepton(), FinalState::FinalState(), FinalStateHCM::FinalStateHCM(), FoxWolframMoments::FoxWolframMoments(), HadronicFinalState::HadronicFinalState(), Hemispheres::Hemispheres(), IdentifiedFinalState::IdentifiedFinalState(), ZEUS_2001_S4815815::init(), UA5_1989_S1926373::init(), UA5_1988_S1867512::init(), UA5_1986_S1583476::init(), UA5_1982_S875503::init(), UA1_1990_S2044935::init(), STAR_2009_UE_HELEN::init(), STAR_2008_S7993412::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_1998_S3780481::init(), MC_ZJETS::init(), MC_WJETS::init(), MC_TTBAR::init(), MC_SUSY::init(), MC_PHOTONJETUE::init(), MC_PHOTONJETS::init(), MC_LEADINGJETS::init(), MC_JETS::init(), MC_DIPHOTON::init(), MC_DIJET::init(), JADE_OPAL_2000_S4300807::init(), H1_2000_S4129130::init(), H1_1995_S3167097::init(), H1_1994_S2919893::init(), ExampleAnalysis::init(), E735_1998_S3905616::init(), DELPHI_2002_069_CONF_603::init(), DELPHI_1995_S3137023::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_1998_S3711838::init(), D0_1996_S3324664::init(), D0_1996_S3214044::init(), CDF_2009_S8436959::init(), CDF_2009_S8383952::init(), CDF_2009_S8233977::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_2002_S4796047::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_1991_S2313472::init(), CDF_1990_S2089246::init(), CDF_1988_S1865951::init(), BELLE_2006_S6265367::init(), ATLAS_2010_S8591806::init(), ALEPH_2004_S5765862::init(), ALEPH_1996_S3486095::init(), ALEPH_1996_S3196992::init(), ALEPH_1991_S2435284::init(), IsolationProjection::IsolationProjection(), JetAlg::JetAlg(), JetShape::JetShape(), KtJets::KtJets(), LeadingParticlesFinalState::LeadingParticlesFinalState(), LossyFinalState::LossyFinalState(), MergedFinalState::MergedFinalState(), Multiplicity::Multiplicity(), NeutralFinalState::NeutralFinalState(), ParisiTensor::ParisiTensor(), Sphericity::Sphericity(), SVertex::SVertex(), Thrust::Thrust(), TotalVisibleMomentum::TotalVisibleMomentum(), TriggerCDFRun0Run1::TriggerCDFRun0Run1(), TriggerUA5::TriggerUA5(), and VetoedFinalState::VetoedFinalState().

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

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, Log::ERROR, ProjectionApplier::getLog(), ProjectionApplier::getProjHandler(), ProjectionApplier::name(), Projection::name(), and ProjectionHandler::registerProjection().

Referenced by ProjectionApplier::addProjection().

00034                                                                              {
00035     if (!_allowProjReg) {
00036       getLog() << Log::ERROR << "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   }


Friends And Related Function Documentation

friend class Event [friend, inherited]

Event is a friend.

Definition at line 31 of file Projection.hh.

friend class Cmp< Projection > [friend, inherited]

The Cmp specialization for Projection is a friend.

Definition at line 34 of file Projection.hh.

friend class Projectionhandler [friend, inherited]

Definition at line 23 of file ProjectionApplier.hh.


Member Data Documentation

vector<pair<double,double> > _etaRanges [protected, inherited]

The ranges allowed for pseudorapidity.

Definition at line 99 of file FinalState.hh.

Referenced by FinalState::accept(), FinalState::compare(), FinalState::FinalState(), and FinalState::project().

double _ptmin [protected, inherited]

The minimum allowed transverse momentum.

Definition at line 102 of file FinalState.hh.

Referenced by FinalState::accept(), FinalState::compare(), and FinalState::project().

ParticleVector _theParticles [mutable, protected, inherited]

bool _allowProjReg [protected, inherited]

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

Definition at line 141 of file ProjectionApplier.hh.

Referenced by ProjectionApplier::_addProjection(), and AnalysisHandler::init().


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