CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
lhef::JetInput Class Reference

#include <JetInput.h>

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.

Constructor & Destructor Documentation

lhef::JetInput::JetInput ( )

Definition at line 17 of file JetInput.cc.

17  :
18  partonicFinalState(false),
19  onlyHardProcess(false),
20  excludeResonances(false),
21  tausAsJets(false),
22  ptMin(0.0)
23 {
24 }
bool partonicFinalState
Definition: JetInput.h:88
bool excludeResonances
Definition: JetInput.h:90
double ptMin
Definition: JetInput.h:92
bool onlyHardProcess
Definition: JetInput.h:89
bool tausAsJets
Definition: JetInput.h:91
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().

26  :
27  partonicFinalState(params.getParameter<bool>("partonicFinalState")),
28  onlyHardProcess(params.getParameter<bool>("onlyHardProcess")),
29  excludeResonances(false),
30  tausAsJets(params.getParameter<bool>("tausAsJets")),
31  ptMin(0.0)
32 {
33  if (params.exists("ignoreParticleIDs"))
35  params.getParameter<std::vector<unsigned int> >(
36  "ignoreParticleIDs"));
37  if (params.exists("excludedResonances"))
39  params.getParameter<std::vector<unsigned int> >(
40  "excludedResonances"));
41  if (params.exists("excludedFromResonances"))
43  params.getParameter<std::vector<unsigned int> >(
44  "excludedFromResonances"));
45 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool partonicFinalState
Definition: JetInput.h:88
bool excludeResonances
Definition: JetInput.h:90
double ptMin
Definition: JetInput.h:92
void setExcludedFromResonances(const std::vector< unsigned int > &particleIds)
Definition: JetInput.cc:66
bool onlyHardProcess
Definition: JetInput.h:89
bool tausAsJets
Definition: JetInput.h:91
void setIgnoredParticles(const std::vector< unsigned int > &particleIDs)
Definition: JetInput.cc:51
void setExcludedResonances(const std::vector< unsigned int > &particleIds)
Definition: JetInput.cc:58
lhef::JetInput::~JetInput ( )

Definition at line 47 of file JetInput.cc.

48 {
49 }

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 customizeTrackingMonitorSeedNumber::idx, isHardProcess(), getDQMSummary::iter, lhef::partIdx(), and findQualityFiles::v.

Referenced by operator()().

205 {
206  unsigned int idx = partIdx(p, particle);
207 
208  if (invalid[idx])
209  return false;
210 
211  const HepMC::GenVertex *v = particle->production_vertex();
212  if (!v)
213  return false;
214  else if (isHardProcess(v, hardProcess))
215  return true;
216 
217  for(HepMC::GenVertex::particles_in_const_iterator iter =
218  v->particles_in_const_begin();
219  iter != v->particles_in_const_end(); ++iter)
220  if (fromHardProcess(invalid, p, *iter, hardProcess))
221  return true;
222 
223  return false;
224 }
static unsigned int partIdx(const JetInput::ParticleVector &p, const HepMC::GenParticle *particle)
Definition: JetInput.cc:130
bool isHardProcess(const HepMC::GenVertex *vertex, const VertexVector &hardProcess) const
Definition: JetInput.cc:121
bool fromHardProcess(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const VertexVector &hardProcess) const
Definition: JetInput.cc:201
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
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, customizeTrackingMonitorSeedNumber::idx, isExcludedFromResonances(), isHardProcess(), isIgnored(), isParton(), isResonance(), getDQMSummary::iter, kDirect, kIndirect, kNo, lhef::partIdx(), query::result, and findQualityFiles::v.

Referenced by operator()().

257 {
258  unsigned int idx = partIdx(p, particle);
259  int id = particle->pdg_id();
260 
261  if (invalid[idx])
262  return kIndirect;
263 
264  if (particle->status() == 3 && isResonance(id))
265  return kDirect;
266 
267  if (!isIgnored(id) && excludedFromResonances.empty() && isParton(id))
268  return kNo;
269 
270  const HepMC::GenVertex *v = particle->production_vertex();
271  if (!v)
272  return kNo;
273  else if (v == signalVertex && excludedResonances.empty())
274  return kIndirect;
275  else if (v == signalVertex || isHardProcess(v, hardProcess))
276  return kNo;
277 
278  for(HepMC::GenVertex::particles_in_const_iterator iter =
279  v->particles_in_const_begin();
280  iter != v->particles_in_const_end(); ++iter) {
283  signalVertex, hardProcess);
284  switch(result) {
285  case kNo:
286  break;
287  case kDirect:
288  if ((*iter)->pdg_id() == id)
289  return kDirect;
290  if (!isExcludedFromResonances(id))
291  break;
292  case kIndirect:
293  return kIndirect;
294  }
295  }
296 
297  return kNo;
298 }
static unsigned int partIdx(const JetInput::ParticleVector &p, const HepMC::GenParticle *particle)
Definition: JetInput.cc:130
bool isHardProcess(const HepMC::GenVertex *vertex, const VertexVector &hardProcess) const
Definition: JetInput.cc:121
std::vector< unsigned int > excludedFromResonances
Definition: JetInput.h:87
bool isExcludedFromResonances(int pdgId) const
Definition: JetInput.cc:113
tuple result
Definition: query.py:137
bool isParton(int pdgId) const
Definition: JetInput.cc:73
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
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
bool isIgnored(int pdgId) const
Definition: JetInput.cc:108
std::vector< unsigned int > excludedResonances
Definition: JetInput.h:86
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 customizeTrackingMonitorSeedNumber::idx, getDQMSummary::iter, lhef::partIdx(), and findQualityFiles::v.

230 {
231  unsigned int idx = partIdx(p, particle);
232 
233  if (invalid[idx])
234  return false;
235 
236  const HepMC::GenVertex *v = particle->production_vertex();
237  if (!v)
238  return false;
239  else if (v == signalVertex)
240  return true;
241 
242  for(HepMC::GenVertex::particles_in_const_iterator iter =
243  v->particles_in_const_begin();
244  iter != v->particles_in_const_end(); ++iter)
245  if (fromSignalVertex(invalid, p, *iter, signalVertex))
246  return true;
247 
248  return false;
249 }
static unsigned int partIdx(const JetInput::ParticleVector &p, const HepMC::GenParticle *particle)
Definition: JetInput.cc:130
bool fromSignalVertex(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const HepMC::GenVertex *signalVertex) const
Definition: JetInput.cc:226
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
const std::vector<unsigned int>& lhef::JetInput::getExcludedFromResonances ( ) const
inline

Definition at line 35 of file JetInput.h.

References excludedFromResonances.

36  { return excludedFromResonances; }
std::vector< unsigned int > excludedFromResonances
Definition: JetInput.h:87
const std::vector<unsigned int>& lhef::JetInput::getExcludedResonances ( ) const
inline

Definition at line 33 of file JetInput.h.

References excludedResonances.

34  { return excludedResonances; }
std::vector< unsigned int > excludedResonances
Definition: JetInput.h:86
bool lhef::JetInput::getExcludeResonances ( ) const
inline

Definition at line 29 of file JetInput.h.

References excludeResonances.

29 { return excludeResonances; }
bool excludeResonances
Definition: JetInput.h:90
bool lhef::JetInput::getHardProcessOnly ( ) const
inline

Definition at line 27 of file JetInput.h.

References onlyHardProcess.

27 { return onlyHardProcess; }
bool onlyHardProcess
Definition: JetInput.h:89
const std::vector<unsigned int>& lhef::JetInput::getIgnoredParticles ( ) const
inline

Definition at line 31 of file JetInput.h.

References ignoreParticleIDs.

32  { return ignoreParticleIDs; }
std::vector< unsigned int > ignoreParticleIDs
Definition: JetInput.h:85
bool lhef::JetInput::getPartonicFinalState ( ) const
inline

Definition at line 26 of file JetInput.h.

References partonicFinalState.

26 { return partonicFinalState; }
bool partonicFinalState
Definition: JetInput.h:88
double lhef::JetInput::getPtMin ( ) const
inline

Definition at line 30 of file JetInput.h.

References ptMin.

30 { return ptMin; }
double ptMin
Definition: JetInput.h:92
bool lhef::JetInput::getTausAndJets ( ) const
inline

Definition at line 28 of file JetInput.h.

References tausAsJets.

28 { return tausAsJets; }
bool tausAsJets
Definition: JetInput.h:91
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()().

197 {
198  return testPartonChildren(invalid, p, particle->end_vertex()) > 0;
199 }
int testPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenVertex *v) const
Definition: JetInput.cc:165
bool lhef::JetInput::isExcludedFromResonances ( int  pdgId) const

Definition at line 113 of file JetInput.cc.

References excludedFromResonances, and lhef::isContained().

Referenced by fromResonance().

114 {
115  if (excludedFromResonances.empty())
116  return true;
117 
119 }
static bool isContained(const std::vector< unsigned int > &list, int id)
Definition: JetInput.cc:89
std::vector< unsigned int > excludedFromResonances
Definition: JetInput.h:87
bool lhef::JetInput::isHadron ( int  pdgId)
static

Definition at line 82 of file JetInput.cc.

References benchmark_cfg::pdgId.

Referenced by testPartonChildren().

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

Definition at line 121 of file JetInput.cc.

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

123 {
124  std::vector<const HepMC::GenVertex*>::const_iterator pos =
125  std::lower_bound(hardProcess.begin(),
126  hardProcess.end(), vertex);
127  return pos != hardProcess.end() && *pos == vertex;
128 }
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()().

109 {
111 }
static bool isContained(const std::vector< unsigned int > &list, int id)
Definition: JetInput.cc:89
std::vector< unsigned int > ignoreParticleIDs
Definition: JetInput.h:85
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().

74 {
75  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
76  return (pdgId > 0 && pdgId < 6) || pdgId == 7 ||
77  pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
78  // tops are not considered "regular" partons
79  // but taus eventually are (since they may hadronize later)
80 }
bool tausAsJets
Definition: JetInput.h:91
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().

98 {
99  if (excludedResonances.empty()) {
100  // gauge bosons and tops
101  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
102  return (pdgId > 21 && pdgId <= 39) || pdgId == 6 || pdgId == 8;
103  }
104 
106 }
static bool isContained(const std::vector< unsigned int > &list, int id)
Definition: JetInput.cc:89
std::vector< unsigned int > excludedResonances
Definition: JetInput.h:86
JetInput::ParticleVector lhef::JetInput::operator() ( const HepMC::GenEvent *  event) const

Definition at line 300 of file JetInput.cc.

References excludeResonances, fromHardProcess(), fromResonance(), GenParticle::GenParticle, hasPartonChildren(), i, align::invalid, lhef::invalidateTree(), isHardProcess(), isIgnored(), isParton(), getDQMSummary::iter, GetRecoTauVFromDQM_MC_cff::next, onlyHardProcess, partonicFinalState, ptMin, MetAnalyzer::pv(), query::result, findQualityFiles::size, python.multivaluedict::sort(), findQualityFiles::v, and TriggerAnalyzer::veto.

302 {
303  if (!event->signal_process_vertex())
304  throw cms::Exception("InvalidHepMCEvent")
305  << "HepMC event is lacking signal vertex."
306  << std::endl;
307 
308  VertexVector hardProcess, toLookAt;
309  std::pair<HepMC::GenParticle*,HepMC::GenParticle*> beamParticles =
310  event->beam_particles();
311  toLookAt.push_back(event->signal_process_vertex());
312  while(!toLookAt.empty()) {
313  std::vector<const HepMC::GenVertex*> next;
314  for(std::vector<const HepMC::GenVertex*>::const_iterator v =
315  toLookAt.begin(); v != toLookAt.end(); ++v) {
316  if (!*v || isHardProcess(*v, hardProcess))
317  continue;
318 
319  bool veto = false;
320  for(HepMC::GenVertex::particles_in_const_iterator iter =
321  (*v)->particles_in_const_begin();
322  iter != (*v)->particles_in_const_end(); ++iter) {
323  if (*iter == beamParticles.first ||
324  *iter == beamParticles.second) {
325  veto = true;
326  break;
327  }
328  }
329  if (veto)
330  continue;
331 
332  hardProcess.push_back(*v);
333  std::sort(hardProcess.begin(), hardProcess.end());
334 
335  for(HepMC::GenVertex::particles_in_const_iterator iter =
336  (*v)->particles_in_const_begin();
337  iter != (*v)->particles_in_const_end(); ++iter) {
338  const HepMC::GenVertex *pv =
339  (*iter)->production_vertex();
340  if (pv)
341  next.push_back(pv);
342  }
343  }
344 
345  toLookAt = next;
346  std::sort(toLookAt.begin(), toLookAt.end());
347  }
348 
349  ParticleVector particles;
350  for(HepMC::GenEvent::particle_const_iterator iter = event->particles_begin();
351  iter != event->particles_end(); ++iter)
352  particles.push_back(*iter);
353 
354  std::sort(particles.begin(), particles.end());
355  unsigned int size = particles.size();
356 
357  ParticleBitmap selected(size, false);
358  ParticleBitmap invalid(size, false);
359 
360  for(unsigned int i = 0; i < size; i++) {
361  const HepMC::GenParticle *particle = particles[i];
362  if (invalid[i])
363  continue;
364  if (particle->status() == 1)
365  selected[i] = true;
366 
367  if (partonicFinalState && isParton(particle->pdg_id())) {
368  if (!particle->end_vertex() &&
369  particle->status() != 1)
370  // some brokenness in event...
371  invalid[i] = true;
372  else if (!hasPartonChildren(invalid, particles,
373  particle)) {
374  selected[i] = true;
375  invalidateTree(invalid, particles,
376  particle->end_vertex());
377  }
378  }
379 
380  if (onlyHardProcess &&
381  !fromHardProcess(invalid, particles,
382  particle, hardProcess) &&
383  !isHardProcess(particle->end_vertex(), hardProcess))
384  invalid[i] = true;
385  }
386 
388  for(unsigned int i = 0; i < size; i++) {
389  const HepMC::GenParticle *particle = particles[i];
390  if (!selected[i] || invalid[i])
391  continue;
392 
393  if (excludeResonances &&
394  fromResonance(invalid, particles, particle,
395  event->signal_process_vertex(),
396  hardProcess)) {
397  invalid[i] = true;
398  continue;
399  }
400 
401  if (isIgnored(particle->pdg_id()))
402  continue;
403 
404  if (particle->momentum().perp() >= ptMin)
405  result.push_back(particle);
406  }
407 
408  return result;
409 }
int i
Definition: DBlmapReader.cc:9
std::vector< bool > ParticleBitmap
Definition: JetInput.h:16
bool isHardProcess(const HepMC::GenVertex *vertex, const VertexVector &hardProcess) const
Definition: JetInput.cc:121
bool partonicFinalState
Definition: JetInput.h:88
bool excludeResonances
Definition: JetInput.h:90
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
tuple result
Definition: query.py:137
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
static void invalidateTree(JetInput::ParticleBitmap &invalid, const JetInput::ParticleVector &p, const HepMC::GenVertex *v)
Definition: JetInput.cc:143
bool isParton(int pdgId) const
Definition: JetInput.cc:73
bool onlyHardProcess
Definition: JetInput.h:89
ResonanceState fromResonance(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle, const HepMC::GenVertex *signalVertex, const VertexVector &hardProcess) const
Definition: JetInput.cc:252
std::vector< const HepMC::GenParticle * > ParticleVector
Definition: JetInput.h:17
std::vector< const HepMC::GenVertex * > VertexVector
Definition: JetInput.h:18
bool hasPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenParticle *particle) const
Definition: JetInput.cc:194
bool isIgnored(int pdgId) const
Definition: JetInput.cc:108
tuple size
Write out results.
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().

68 {
69  excludedFromResonances = particleIDs;
71 }
std::vector< unsigned int > excludedFromResonances
Definition: JetInput.h:87
void lhef::JetInput::setExcludedResonances ( const std::vector< unsigned int > &  particleIds)

Definition at line 58 of file JetInput.cc.

References excludedResonances, setExcludeResonances(), and python.multivaluedict::sort().

Referenced by JetInput().

60 {
62  excludedResonances = particleIDs;
64 }
void setExcludeResonances(bool flag=true)
Definition: JetInput.h:43
std::vector< unsigned int > excludedResonances
Definition: JetInput.h:86
void lhef::JetInput::setExcludeResonances ( bool  flag = true)
inline

Definition at line 43 of file JetInput.h.

References excludeResonances, and archive::flag.

Referenced by setExcludedResonances().

bool excludeResonances
Definition: JetInput.h:90
void lhef::JetInput::setHardProcessOnly ( bool  flag = true)
inline

Definition at line 40 of file JetInput.h.

References archive::flag, and onlyHardProcess.

41  { onlyHardProcess = flag; }
bool onlyHardProcess
Definition: JetInput.h:89
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().

53 {
54  ignoreParticleIDs = particleIDs;
56 }
std::vector< unsigned int > ignoreParticleIDs
Definition: JetInput.h:85
void lhef::JetInput::setPartonicFinalState ( bool  flag = true)
inline

Definition at line 38 of file JetInput.h.

References archive::flag, and partonicFinalState.

bool partonicFinalState
Definition: JetInput.h:88
void lhef::JetInput::setPtMin ( double  ptMin)
inline

Definition at line 44 of file JetInput.h.

References ptMin.

44 { this->ptMin = ptMin; }
double ptMin
Definition: JetInput.h:92
void lhef::JetInput::setTausAsJets ( bool  flag = true)
inline

Definition at line 42 of file JetInput.h.

References archive::flag, and tausAsJets.

42 { tausAsJets = flag; }
bool tausAsJets
Definition: JetInput.h:91
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 customizeTrackingMonitorSeedNumber::idx, isHadron(), isParton(), getDQMSummary::iter, lhef::partIdx(), and query::result.

Referenced by hasPartonChildren().

168 {
169  if (!v)
170  return 0;
171 
172  for(HepMC::GenVertex::particles_out_const_iterator iter =
173  v->particles_out_const_begin();
174  iter != v->particles_out_const_end(); ++iter) {
175  unsigned int idx = partIdx(p, *iter);
176 
177  if (invalid[idx])
178  continue;
179 
180  if (isParton((*iter)->pdg_id()))
181  return 1;
182  if (isHadron((*iter)->pdg_id()))
183  return -1;
184 
185  const HepMC::GenVertex *v = (*iter)->end_vertex();
186  int result = testPartonChildren(invalid, p, v);
187  if (result)
188  return result;
189  }
190 
191  return 0;
192 }
static unsigned int partIdx(const JetInput::ParticleVector &p, const HepMC::GenParticle *particle)
Definition: JetInput.cc:130
tuple result
Definition: query.py:137
static bool isHadron(int pdgId)
Definition: JetInput.cc:82
bool isParton(int pdgId) const
Definition: JetInput.cc:73
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
int testPartonChildren(ParticleBitmap &invalid, const ParticleVector &p, const HepMC::GenVertex *v) const
Definition: JetInput.cc:165

Member Data Documentation

std::vector<unsigned int> lhef::JetInput::excludedFromResonances
private
std::vector<unsigned int> lhef::JetInput::excludedResonances
private
bool lhef::JetInput::excludeResonances
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().

bool lhef::JetInput::onlyHardProcess
private

Definition at line 89 of file JetInput.h.

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

bool lhef::JetInput::partonicFinalState
private

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().

bool lhef::JetInput::tausAsJets
private

Definition at line 91 of file JetInput.h.

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