Go to the documentation of this file.00001 #ifndef GeneratorInterface_LHEInterface_JetInput_h
00002 #define GeneratorInterface_LHEInterface_JetInput_h
00003
00004 #include <vector>
00005
00006 #include <HepMC/GenEvent.h>
00007 #include <HepMC/GenParticle.h>
00008 #include <HepMC/GenVertex.h>
00009
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011
00012 namespace lhef {
00013
00014 class JetInput {
00015 public:
00016 typedef std::vector<bool> ParticleBitmap;
00017 typedef std::vector<const HepMC::GenParticle*> ParticleVector;
00018 typedef std::vector<const HepMC::GenVertex*> VertexVector;
00019
00020 JetInput();
00021 JetInput(const edm::ParameterSet ¶ms);
00022 ~JetInput();
00023
00024 ParticleVector operator () (const HepMC::GenEvent *event) const;
00025
00026 bool getPartonicFinalState() const { return partonicFinalState; }
00027 bool getHardProcessOnly() const { return onlyHardProcess; }
00028 bool getTausAndJets() const { return tausAsJets; }
00029 bool getExcludeResonances() const { return excludeResonances; }
00030 double getPtMin() const { return ptMin; }
00031 const std::vector<unsigned int> &getIgnoredParticles() const
00032 { return ignoreParticleIDs; }
00033 const std::vector<unsigned int> &getExcludedResonances() const
00034 { return excludedResonances; }
00035 const std::vector<unsigned int> &getExcludedFromResonances() const
00036 { return excludedFromResonances; }
00037
00038 void setPartonicFinalState(bool flag = true)
00039 { partonicFinalState = flag; }
00040 void setHardProcessOnly(bool flag = true)
00041 { onlyHardProcess = flag; }
00042 void setTausAsJets(bool flag = true) { tausAsJets = flag; }
00043 void setExcludeResonances(bool flag = true) { excludeResonances = flag; }
00044 void setPtMin(double ptMin) { this->ptMin = ptMin; }
00045 void setIgnoredParticles(const std::vector<unsigned int> &particleIDs);
00046 void setExcludedResonances(const std::vector<unsigned int> &particleIds);
00047 void setExcludedFromResonances(const std::vector<unsigned int> &particleIds);
00048
00049 bool isParton(int pdgId) const;
00050 static bool isHadron(int pdgId);
00051 bool isResonance(int pdgId) const;
00052 bool isExcludedFromResonances(int pdgId) const;
00053 bool isHardProcess(const HepMC::GenVertex *vertex,
00054 const VertexVector &hardProcess) const;
00055 bool isIgnored(int pdgId) const;
00056
00057 enum ResonanceState {
00058 kNo = 0,
00059 kDirect,
00060 kIndirect
00061 };
00062
00063 bool hasPartonChildren(ParticleBitmap &invalid,
00064 const ParticleVector &p,
00065 const HepMC::GenParticle *particle) const;
00066 bool fromHardProcess(ParticleBitmap &invalid,
00067 const ParticleVector &p,
00068 const HepMC::GenParticle *particle,
00069 const VertexVector &hardProcess) const;
00070 bool fromSignalVertex(ParticleBitmap &invalid,
00071 const ParticleVector &p,
00072 const HepMC::GenParticle *particle,
00073 const HepMC::GenVertex *signalVertex) const;
00074 ResonanceState fromResonance(ParticleBitmap &invalid,
00075 const ParticleVector &p,
00076 const HepMC::GenParticle *particle,
00077 const HepMC::GenVertex *signalVertex,
00078 const VertexVector &hardProcess) const;
00079
00080 private:
00081 int testPartonChildren(ParticleBitmap &invalid,
00082 const ParticleVector &p,
00083 const HepMC::GenVertex *v) const;
00084
00085 std::vector<unsigned int> ignoreParticleIDs;
00086 std::vector<unsigned int> excludedResonances;
00087 std::vector<unsigned int> excludedFromResonances;
00088 bool partonicFinalState;
00089 bool onlyHardProcess;
00090 bool excludeResonances;
00091 bool tausAsJets;
00092 double ptMin;
00093 };
00094
00095 }
00096
00097 #endif // GeneratorCommon_LHEInterface_JetInput_h