CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetInput.h
Go to the documentation of this file.
1 #ifndef GeneratorInterface_LHEInterface_JetInput_h
2 #define GeneratorInterface_LHEInterface_JetInput_h
3 
4 #include <vector>
5 
6 #include <HepMC/GenEvent.h>
7 #include <HepMC/GenParticle.h>
8 #include <HepMC/GenVertex.h>
9 
11 
12 namespace lhef {
13 
14 class JetInput {
15  public:
16  typedef std::vector<bool> ParticleBitmap;
17  typedef std::vector<const HepMC::GenParticle*> ParticleVector;
18  typedef std::vector<const HepMC::GenVertex*> VertexVector;
19 
20  JetInput();
21  JetInput(const edm::ParameterSet &params);
22  ~JetInput();
23 
24  ParticleVector operator () (const HepMC::GenEvent *event) const;
25 
26  bool getPartonicFinalState() const { return partonicFinalState; }
27  bool getHardProcessOnly() const { return onlyHardProcess; }
28  bool getTausAndJets() const { return tausAsJets; }
29  bool getExcludeResonances() const { return excludeResonances; }
30  double getPtMin() const { return ptMin; }
31  const std::vector<unsigned int> &getIgnoredParticles() const
32  { return ignoreParticleIDs; }
33  const std::vector<unsigned int> &getExcludedResonances() const
34  { return excludedResonances; }
35  const std::vector<unsigned int> &getExcludedFromResonances() const
36  { return excludedFromResonances; }
37 
38  void setPartonicFinalState(bool flag = true)
40  void setHardProcessOnly(bool flag = true)
41  { onlyHardProcess = flag; }
42  void setTausAsJets(bool flag = true) { tausAsJets = flag; }
44  void setPtMin(double ptMin) { this->ptMin = ptMin; }
45  void setIgnoredParticles(const std::vector<unsigned int> &particleIDs);
46  void setExcludedResonances(const std::vector<unsigned int> &particleIds);
47  void setExcludedFromResonances(const std::vector<unsigned int> &particleIds);
48 
49  bool isParton(int pdgId) const;
50  static bool isHadron(int pdgId);
51  bool isResonance(int pdgId) const;
52  bool isExcludedFromResonances(int pdgId) const;
53  bool isHardProcess(const HepMC::GenVertex *vertex,
54  const VertexVector &hardProcess) const;
55  bool isIgnored(int pdgId) const;
56 
58  kNo = 0,
61  };
62 
64  const ParticleVector &p,
65  const HepMC::GenParticle *particle) const;
67  const ParticleVector &p,
68  const HepMC::GenParticle *particle,
69  const VertexVector &hardProcess) const;
71  const ParticleVector &p,
72  const HepMC::GenParticle *particle,
73  const HepMC::GenVertex *signalVertex) const;
75  const ParticleVector &p,
76  const HepMC::GenParticle *particle,
77  const HepMC::GenVertex *signalVertex,
78  const VertexVector &hardProcess) const;
79 
80  private:
82  const ParticleVector &p,
83  const HepMC::GenVertex *v) const;
84 
85  std::vector<unsigned int> ignoreParticleIDs;
86  std::vector<unsigned int> excludedResonances;
87  std::vector<unsigned int> excludedFromResonances;
91  bool tausAsJets;
92  double ptMin;
93 };
94 
95 } // namespace lhef
96 
97 #endif // GeneratorCommon_LHEInterface_JetInput_h
const std::vector< unsigned int > & getExcludedFromResonances() const
Definition: JetInput.h:35
long int flag
Definition: mlp_lapack.h:47
void setExcludeResonances(bool flag=true)
Definition: JetInput.h:43
std::vector< bool > ParticleBitmap
Definition: JetInput.h:16
bool isHardProcess(const HepMC::GenVertex *vertex, const VertexVector &hardProcess) const
Definition: JetInput.cc:121
void setTausAsJets(bool flag=true)
Definition: JetInput.h:42
std::vector< unsigned int > excludedFromResonances
Definition: JetInput.h:87
void setPtMin(double ptMin)
Definition: JetInput.h:44
bool fromSignalVertex(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const HepMC::GenVertex *signalVertex) const
Definition: JetInput.cc:226
bool partonicFinalState
Definition: JetInput.h:88
bool excludeResonances
Definition: JetInput.h:90
bool getHardProcessOnly() const
Definition: JetInput.h:27
bool fromHardProcess(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const VertexVector &hardProcess) const
Definition: JetInput.cc:201
double ptMin
Definition: JetInput.h:92
bool isExcludedFromResonances(int pdgId) const
Definition: JetInput.cc:113
static bool isHadron(int pdgId)
Definition: JetInput.cc:82
ParticleVector operator()(const HepMC::GenEvent *event) const
Definition: JetInput.cc:300
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
void setExcludedFromResonances(const std::vector< unsigned int > &particleIds)
Definition: JetInput.cc:66
bool isParton(int pdgId) const
Definition: JetInput.cc:73
std::vector< unsigned int > ignoreParticleIDs
Definition: JetInput.h:85
const std::vector< unsigned int > & getIgnoredParticles() const
Definition: JetInput.h:31
bool getExcludeResonances() const
Definition: JetInput.h:29
bool onlyHardProcess
Definition: JetInput.h:89
bool isResonance(int pdgId) const
Definition: JetInput.cc:97
ResonanceState fromResonance(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const HepMC::GenVertex *signalVertex, const VertexVector &hardProcess) const
Definition: JetInput.cc:252
void setPartonicFinalState(bool flag=true)
Definition: JetInput.h:38
std::vector< const HepMC::GenParticle * > ParticleVector
Definition: JetInput.h:17
int testPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenVertex *v) const
Definition: JetInput.cc:165
std::vector< const HepMC::GenVertex * > VertexVector
Definition: JetInput.h:18
bool tausAsJets
Definition: JetInput.h:91
void setIgnoredParticles(const std::vector< unsigned int > &particleIDs)
Definition: JetInput.cc:51
bool hasPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle) const
Definition: JetInput.cc:194
const std::vector< unsigned int > & getExcludedResonances() const
Definition: JetInput.h:33
void setHardProcessOnly(bool flag=true)
Definition: JetInput.h:40
bool getTausAndJets() const
Definition: JetInput.h:28
bool getPartonicFinalState() const
Definition: JetInput.h:26
void setExcludedResonances(const std::vector< unsigned int > &particleIds)
Definition: JetInput.cc:58
bool isIgnored(int pdgId) const
Definition: JetInput.cc:108
std::vector< unsigned int > excludedResonances
Definition: JetInput.h:86
double getPtMin() const
Definition: JetInput.h:30
mathSSE::Vec4< T > v