CMS 3D CMS Logo

Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes

lhef::JetInput Class Reference

#include <JetInput.h>

List of all members.

Public Types

typedef std::vector< bool > ParticleBitmap
typedef std::vector< const
HepMC::GenParticle * > 
ParticleVector
enum  ResonanceState { kNo = 0, kDirect, kIndirect }
typedef std::vector< const
HepMC::GenVertex * > 
VertexVector

Public Member Functions

bool fromHardProcess (ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const VertexVector &hardProcess) const
ResonanceState fromResonance (ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const HepMC::GenVertex *signalVertex, const VertexVector &hardProcess) const
bool fromSignalVertex (ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const HepMC::GenVertex *signalVertex) const
const std::vector< unsigned int > & getExcludedFromResonances () const
const std::vector< unsigned int > & getExcludedResonances () const
bool getExcludeResonances () const
bool getHardProcessOnly () const
const std::vector< unsigned int > & getIgnoredParticles () const
bool getPartonicFinalState () const
double getPtMin () const
bool getTausAndJets () const
bool hasPartonChildren (ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle) const
bool isExcludedFromResonances (int pdgId) const
bool isHardProcess (const HepMC::GenVertex *vertex, const VertexVector &hardProcess) const
bool isIgnored (int pdgId) const
bool isParton (int pdgId) const
bool isResonance (int pdgId) const
 JetInput ()
 JetInput (const edm::ParameterSet &params)
ParticleVector operator() (const HepMC::GenEvent *event) const
void setExcludedFromResonances (const std::vector< unsigned int > &particleIds)
void setExcludedResonances (const std::vector< unsigned int > &particleIds)
void setExcludeResonances (bool flag=true)
void setHardProcessOnly (bool flag=true)
void setIgnoredParticles (const std::vector< unsigned int > &particleIDs)
void setPartonicFinalState (bool flag=true)
void setPtMin (double ptMin)
void setTausAsJets (bool flag=true)
 ~JetInput ()

Static Public Member Functions

static bool isHadron (int pdgId)

Private Member Functions

int testPartonChildren (ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenVertex *v) const

Private Attributes

std::vector< unsigned int > excludedFromResonances
std::vector< unsigned int > excludedResonances
bool excludeResonances
std::vector< unsigned int > ignoreParticleIDs
bool onlyHardProcess
bool partonicFinalState
double ptMin
bool tausAsJets

Detailed Description

Definition at line 14 of file JetInput.h.


Member Typedef Documentation

typedef std::vector<bool> lhef::JetInput::ParticleBitmap

Definition at line 16 of file JetInput.h.

typedef std::vector<const HepMC::GenParticle*> lhef::JetInput::ParticleVector

Definition at line 17 of file JetInput.h.

typedef std::vector<const HepMC::GenVertex*> lhef::JetInput::VertexVector

Definition at line 18 of file JetInput.h.


Member Enumeration Documentation

Enumerator:
kNo 
kDirect 
kIndirect 

Definition at line 57 of file JetInput.h.

                            {
                kNo = 0,
                kDirect,
                kIndirect
        };

Constructor & Destructor Documentation

lhef::JetInput::JetInput ( )

Definition at line 17 of file JetInput.cc.

                   :
        partonicFinalState(false),
        onlyHardProcess(false),
        excludeResonances(false),
        tausAsJets(false),
        ptMin(0.0)
{
}
lhef::JetInput::JetInput ( const edm::ParameterSet params)

Definition at line 26 of file JetInput.cc.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), setExcludedFromResonances(), setExcludedResonances(), and setIgnoredParticles().

                                                :
        partonicFinalState(params.getParameter<bool>("partonicFinalState")),
        onlyHardProcess(params.getParameter<bool>("onlyHardProcess")),
        excludeResonances(false),
        tausAsJets(params.getParameter<bool>("tausAsJets")),
        ptMin(0.0)
{
        if (params.exists("ignoreParticleIDs"))
                setIgnoredParticles(
                        params.getParameter<std::vector<unsigned int> >(
                                                        "ignoreParticleIDs"));
        if (params.exists("excludedResonances"))
                setExcludedResonances(
                        params.getParameter<std::vector<unsigned int> >(
                                                        "excludedResonances"));
        if (params.exists("excludedFromResonances"))
                setExcludedFromResonances(
                        params.getParameter<std::vector<unsigned int> >(
                                                "excludedFromResonances"));
}
lhef::JetInput::~JetInput ( )

Definition at line 47 of file JetInput.cc.

{
}

Member Function Documentation

bool lhef::JetInput::fromHardProcess ( ParticleBitmap invalid,
const ParticleVector p,
const HepMC::GenParticle *  particle,
const VertexVector hardProcess 
) const

Definition at line 201 of file JetInput.cc.

References UserOptions_cff::idx, isHardProcess(), lhef::partIdx(), and v.

Referenced by operator()().

{
        unsigned int idx = partIdx(p, particle);

        if (invalid[idx])
                return false;

        const HepMC::GenVertex *v = particle->production_vertex();
        if (!v)
                return false;
        else if (isHardProcess(v, hardProcess))
                return true;

        for(HepMC::GenVertex::particles_in_const_iterator iter =
                                        v->particles_in_const_begin();
            iter != v->particles_in_const_end(); ++iter)
                if (fromHardProcess(invalid, p, *iter, hardProcess))
                        return true;

        return false;
}
JetInput::ResonanceState lhef::JetInput::fromResonance ( ParticleBitmap invalid,
const ParticleVector p,
const HepMC::GenParticle *  particle,
const HepMC::GenVertex *  signalVertex,
const VertexVector hardProcess 
) const

Definition at line 252 of file JetInput.cc.

References excludedFromResonances, excludedResonances, UserOptions_cff::idx, isExcludedFromResonances(), isHardProcess(), isIgnored(), isParton(), isResonance(), kDirect, kIndirect, kNo, lhef::partIdx(), query::result, and v.

Referenced by operator()().

{
        unsigned int idx = partIdx(p, particle);
        int id = particle->pdg_id();

        if (invalid[idx])
                return kIndirect;

        if (particle->status() == 3 && isResonance(id))
                return kDirect;

        if (!isIgnored(id) && excludedFromResonances.empty() && isParton(id))
                return kNo;

        const HepMC::GenVertex *v = particle->production_vertex();
        if (!v)
                return kNo;
        else if (v == signalVertex && excludedResonances.empty())
                return kIndirect;
        else if (v == signalVertex || isHardProcess(v, hardProcess))
                return kNo;

        for(HepMC::GenVertex::particles_in_const_iterator iter =
                                        v->particles_in_const_begin();
            iter != v->particles_in_const_end(); ++iter) {
                ResonanceState result =
                        fromResonance(invalid, p, *iter,
                                      signalVertex, hardProcess);
                switch(result) {
                    case kNo:
                        break;
                    case kDirect:
                        if ((*iter)->pdg_id() == id)
                                return kDirect;
                        if (!isExcludedFromResonances(id))
                                break;
                    case kIndirect:
                        return kIndirect;
                }
        }

        return kNo;
}
bool lhef::JetInput::fromSignalVertex ( ParticleBitmap invalid,
const ParticleVector p,
const HepMC::GenParticle *  particle,
const HepMC::GenVertex *  signalVertex 
) const

Definition at line 226 of file JetInput.cc.

References UserOptions_cff::idx, lhef::partIdx(), and v.

{
        unsigned int idx = partIdx(p, particle);

        if (invalid[idx])
                return false;

        const HepMC::GenVertex *v = particle->production_vertex();
        if (!v)
                return false;
        else if (v == signalVertex)
                return true;

        for(HepMC::GenVertex::particles_in_const_iterator iter =
                                        v->particles_in_const_begin();
            iter != v->particles_in_const_end(); ++iter)
                if (fromSignalVertex(invalid, p, *iter, signalVertex))
                        return true;

        return false;
}
const std::vector<unsigned int>& lhef::JetInput::getExcludedFromResonances ( ) const [inline]

Definition at line 35 of file JetInput.h.

References excludedFromResonances.

const std::vector<unsigned int>& lhef::JetInput::getExcludedResonances ( ) const [inline]

Definition at line 33 of file JetInput.h.

References excludedResonances.

        { return excludedResonances; }
bool lhef::JetInput::getExcludeResonances ( ) const [inline]

Definition at line 29 of file JetInput.h.

References excludeResonances.

{ return excludeResonances; }
bool lhef::JetInput::getHardProcessOnly ( ) const [inline]

Definition at line 27 of file JetInput.h.

References onlyHardProcess.

{ return onlyHardProcess; }
const std::vector<unsigned int>& lhef::JetInput::getIgnoredParticles ( ) const [inline]

Definition at line 31 of file JetInput.h.

References ignoreParticleIDs.

        { return ignoreParticleIDs; }
bool lhef::JetInput::getPartonicFinalState ( ) const [inline]

Definition at line 26 of file JetInput.h.

References partonicFinalState.

{ return partonicFinalState; }
double lhef::JetInput::getPtMin ( ) const [inline]

Definition at line 30 of file JetInput.h.

References ptMin.

{ return ptMin; }
bool lhef::JetInput::getTausAndJets ( ) const [inline]

Definition at line 28 of file JetInput.h.

References tausAsJets.

{ return tausAsJets; }
bool lhef::JetInput::hasPartonChildren ( ParticleBitmap invalid,
const ParticleVector p,
const HepMC::GenParticle *  particle 
) const

Definition at line 194 of file JetInput.cc.

References testPartonChildren().

Referenced by operator()().

{
        return testPartonChildren(invalid, p, particle->end_vertex()) > 0;
}
bool lhef::JetInput::isExcludedFromResonances ( int  pdgId) const

Definition at line 113 of file JetInput.cc.

References excludedFromResonances, and lhef::isContained().

Referenced by fromResonance().

{
        if (excludedFromResonances.empty())
                return true;

        return isContained(excludedFromResonances, pdgId);
}
bool lhef::JetInput::isHadron ( int  pdgId) [static]

Definition at line 82 of file JetInput.cc.

References benchmark_cfg::pdgId.

Referenced by testPartonChildren().

{
        pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
        return (pdgId > 100 && pdgId < 900) ||
               (pdgId > 1000 && pdgId < 9000);
}
bool lhef::JetInput::isHardProcess ( const HepMC::GenVertex *  vertex,
const VertexVector hardProcess 
) const

Definition at line 121 of file JetInput.cc.

References pos.

Referenced by fromHardProcess(), fromResonance(), and operator()().

{
        std::vector<const HepMC::GenVertex*>::const_iterator pos =
                        std::lower_bound(hardProcess.begin(),
                                         hardProcess.end(), vertex);
        return pos != hardProcess.end() && *pos == vertex;
}
bool lhef::JetInput::isIgnored ( int  pdgId) const

Definition at line 108 of file JetInput.cc.

References ignoreParticleIDs, and lhef::isContained().

Referenced by fromResonance(), and operator()().

bool lhef::JetInput::isParton ( int  pdgId) const

Definition at line 73 of file JetInput.cc.

References benchmark_cfg::pdgId, and tausAsJets.

Referenced by fromResonance(), operator()(), and testPartonChildren().

{
        pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
        return (pdgId > 0 && pdgId < 6) || pdgId == 7 ||
               pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
        // tops are not considered "regular" partons
        // but taus eventually are (since they may hadronize later)
}
bool lhef::JetInput::isResonance ( int  pdgId) const

Definition at line 97 of file JetInput.cc.

References excludedResonances, lhef::isContained(), and benchmark_cfg::pdgId.

Referenced by fromResonance().

{
        if (excludedResonances.empty()) {
                // gauge bosons and tops
                pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
                return (pdgId > 21 && pdgId <= 39) || pdgId == 6 || pdgId == 8;
        }

        return isContained(excludedResonances, pdgId);
}
JetInput::ParticleVector lhef::JetInput::operator() ( const HepMC::GenEvent *  event) const

Definition at line 300 of file JetInput.cc.

References excludeResonances, fromHardProcess(), fromResonance(), configurableAnalysis::GenParticle, hasPartonChildren(), i, align::invalid, lhef::invalidateTree(), isHardProcess(), isIgnored(), isParton(), onlyHardProcess, partonicFinalState, ptMin, query::result, findQualityFiles::size, python::multivaluedict::sort(), and v.

{
        if (!event->signal_process_vertex())
                throw cms::Exception("InvalidHepMCEvent")
                        << "HepMC event is lacking signal vertex."
                        << std::endl;

        VertexVector hardProcess, toLookAt;
        std::pair<HepMC::GenParticle*,HepMC::GenParticle*> beamParticles =
                                                event->beam_particles();
        toLookAt.push_back(event->signal_process_vertex());
        while(!toLookAt.empty()) {
                std::vector<const HepMC::GenVertex*> next;
                for(std::vector<const HepMC::GenVertex*>::const_iterator v =
                        toLookAt.begin(); v != toLookAt.end(); ++v) {
                        if (!*v || isHardProcess(*v, hardProcess))
                                continue;

                        bool veto = false;
                        for(HepMC::GenVertex::particles_in_const_iterator iter =
                                        (*v)->particles_in_const_begin();
                            iter != (*v)->particles_in_const_end(); ++iter) {
                                if (*iter == beamParticles.first ||
                                    *iter == beamParticles.second) {
                                        veto = true;
                                        break;
                                }
                        }
                        if (veto)
                                continue;

                        hardProcess.push_back(*v);
                        std::sort(hardProcess.begin(), hardProcess.end());

                        for(HepMC::GenVertex::particles_in_const_iterator iter =
                                        (*v)->particles_in_const_begin();
                            iter != (*v)->particles_in_const_end(); ++iter) {
                                const HepMC::GenVertex *pv =
                                                (*iter)->production_vertex();
                                if (pv)
                                        next.push_back(pv);
                        }
                }

                toLookAt = next;
                std::sort(toLookAt.begin(), toLookAt.end());
        }

        ParticleVector particles;
        for(HepMC::GenEvent::particle_const_iterator iter = event->particles_begin();
            iter != event->particles_end(); ++iter)
                particles.push_back(*iter);

        std::sort(particles.begin(), particles.end());
        unsigned int size = particles.size();

        ParticleBitmap selected(size, false);
        ParticleBitmap invalid(size, false);

        for(unsigned int i = 0; i < size; i++) {
                const HepMC::GenParticle *particle = particles[i];
                if (invalid[i])
                        continue;
                if (particle->status() == 1)
                        selected[i] = true;

                if (partonicFinalState && isParton(particle->pdg_id())) {
                        if (!particle->end_vertex() &&
                            particle->status() != 1)
                                // some brokenness in event...
                                invalid[i] = true;
                        else if (!hasPartonChildren(invalid, particles,
                                                    particle)) {
                                selected[i] = true;
                                invalidateTree(invalid, particles,
                                               particle->end_vertex());
                        }
                }

                if (onlyHardProcess &&
                    !fromHardProcess(invalid, particles,
                                     particle, hardProcess) &&
                    !isHardProcess(particle->end_vertex(), hardProcess))
                        invalid[i] = true;
        }

        ParticleVector result;
        for(unsigned int i = 0; i < size; i++) {
                const HepMC::GenParticle *particle = particles[i];
                if (!selected[i] || invalid[i])
                        continue;

                if (excludeResonances &&
                    fromResonance(invalid, particles, particle,
                                  event->signal_process_vertex(),
                                  hardProcess)) {
                        invalid[i] = true;
                        continue;
                }

                if (isIgnored(particle->pdg_id()))
                        continue;

                if (particle->momentum().perp() >= ptMin)
                        result.push_back(particle);
        }

        return result;
}
void lhef::JetInput::setExcludedFromResonances ( const std::vector< unsigned int > &  particleIds)

Definition at line 66 of file JetInput.cc.

References excludedFromResonances, and python::multivaluedict::sort().

Referenced by JetInput().

void lhef::JetInput::setExcludedResonances ( const std::vector< unsigned int > &  particleIds)
void lhef::JetInput::setExcludeResonances ( bool  flag = true) [inline]

Definition at line 43 of file JetInput.h.

References excludeResonances.

Referenced by setExcludedResonances().

void lhef::JetInput::setHardProcessOnly ( bool  flag = true) [inline]

Definition at line 40 of file JetInput.h.

References onlyHardProcess.

void lhef::JetInput::setIgnoredParticles ( const std::vector< unsigned int > &  particleIDs)

Definition at line 51 of file JetInput.cc.

References ignoreParticleIDs, and python::multivaluedict::sort().

Referenced by JetInput().

{
        ignoreParticleIDs = particleIDs;
        std::sort(ignoreParticleIDs.begin(), ignoreParticleIDs.end());
}
void lhef::JetInput::setPartonicFinalState ( bool  flag = true) [inline]

Definition at line 38 of file JetInput.h.

References partonicFinalState.

void lhef::JetInput::setPtMin ( double  ptMin) [inline]

Definition at line 44 of file JetInput.h.

References ptMin.

{ this->ptMin = ptMin; }
void lhef::JetInput::setTausAsJets ( bool  flag = true) [inline]

Definition at line 42 of file JetInput.h.

References tausAsJets.

int lhef::JetInput::testPartonChildren ( JetInput::ParticleBitmap invalid,
const ParticleVector p,
const HepMC::GenVertex *  v 
) const [private]

Definition at line 165 of file JetInput.cc.

References UserOptions_cff::idx, isHadron(), isParton(), lhef::partIdx(), and query::result.

Referenced by hasPartonChildren().

{
        if (!v)
                return 0;

        for(HepMC::GenVertex::particles_out_const_iterator iter =
                                        v->particles_out_const_begin();
            iter != v->particles_out_const_end(); ++iter) {
                unsigned int idx = partIdx(p, *iter);

                if (invalid[idx])
                        continue;

                if (isParton((*iter)->pdg_id()))
                        return 1;
                if (isHadron((*iter)->pdg_id()))
                        return -1;

                const HepMC::GenVertex *v = (*iter)->end_vertex();
                int result = testPartonChildren(invalid, p, v);
                if (result)
                        return result;
        }

        return 0;
}

Member Data Documentation

std::vector<unsigned int> lhef::JetInput::excludedFromResonances [private]
std::vector<unsigned int> lhef::JetInput::excludedResonances [private]

Definition at line 90 of file JetInput.h.

Referenced by getExcludeResonances(), operator()(), and setExcludeResonances().

std::vector<unsigned int> lhef::JetInput::ignoreParticleIDs [private]

Definition at line 85 of file JetInput.h.

Referenced by getIgnoredParticles(), isIgnored(), and setIgnoredParticles().

Definition at line 89 of file JetInput.h.

Referenced by getHardProcessOnly(), operator()(), and setHardProcessOnly().

Definition at line 88 of file JetInput.h.

Referenced by getPartonicFinalState(), operator()(), and setPartonicFinalState().

double lhef::JetInput::ptMin [private]

Definition at line 92 of file JetInput.h.

Referenced by getPtMin(), operator()(), and setPtMin().

Definition at line 91 of file JetInput.h.

Referenced by getTausAndJets(), isParton(), and setTausAsJets().