![]() |
![]() |
#include <GeneratorInterface/LHEInterface/interface/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 (const edm::ParameterSet ¶ms) | |
JetInput () | |
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 |
Definition at line 14 of file JetInput.h.
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.
lhef::JetInput::JetInput | ( | ) |
Definition at line 17 of file JetInput.cc.
00017 : 00018 partonicFinalState(false), 00019 onlyHardProcess(false), 00020 excludeResonances(false), 00021 tausAsJets(false), 00022 ptMin(0.0) 00023 { 00024 }
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().
00026 : 00027 partonicFinalState(params.getParameter<bool>("partonicFinalState")), 00028 onlyHardProcess(params.getParameter<bool>("onlyHardProcess")), 00029 excludeResonances(false), 00030 tausAsJets(params.getParameter<bool>("tausAsJets")), 00031 ptMin(0.0) 00032 { 00033 if (params.exists("ignoreParticleIDs")) 00034 setIgnoredParticles( 00035 params.getParameter<std::vector<unsigned int> >( 00036 "ignoreParticleIDs")); 00037 if (params.exists("excludedResonances")) 00038 setExcludedResonances( 00039 params.getParameter<std::vector<unsigned int> >( 00040 "excludedResonances")); 00041 if (params.exists("excludedFromResonances")) 00042 setExcludedFromResonances( 00043 params.getParameter<std::vector<unsigned int> >( 00044 "excludedFromResonances")); 00045 }
lhef::JetInput::~JetInput | ( | ) |
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 isHardProcess(), iter, lhef::partIdx(), and v.
Referenced by operator()().
00205 { 00206 unsigned int idx = partIdx(p, particle); 00207 00208 if (invalid[idx]) 00209 return false; 00210 00211 const HepMC::GenVertex *v = particle->production_vertex(); 00212 if (!v) 00213 return false; 00214 else if (isHardProcess(v, hardProcess)) 00215 return true; 00216 00217 for(HepMC::GenVertex::particles_in_const_iterator iter = 00218 v->particles_in_const_begin(); 00219 iter != v->particles_in_const_end(); ++iter) 00220 if (fromHardProcess(invalid, p, *iter, hardProcess)) 00221 return true; 00222 00223 return false; 00224 }
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, isExcludedFromResonances(), isHardProcess(), isIgnored(), isParton(), isResonance(), iter, kDirect, kIndirect, kNo, lhef::partIdx(), HLT_VtxMuL3::result, and v.
Referenced by operator()().
00257 { 00258 unsigned int idx = partIdx(p, particle); 00259 int id = particle->pdg_id(); 00260 00261 if (invalid[idx]) 00262 return kIndirect; 00263 00264 if (particle->status() == 3 && isResonance(id)) 00265 return kDirect; 00266 00267 if (!isIgnored(id) && excludedFromResonances.empty() && isParton(id)) 00268 return kNo; 00269 00270 const HepMC::GenVertex *v = particle->production_vertex(); 00271 if (!v) 00272 return kNo; 00273 else if (v == signalVertex && excludedResonances.empty()) 00274 return kIndirect; 00275 else if (v == signalVertex || isHardProcess(v, hardProcess)) 00276 return kNo; 00277 00278 for(HepMC::GenVertex::particles_in_const_iterator iter = 00279 v->particles_in_const_begin(); 00280 iter != v->particles_in_const_end(); ++iter) { 00281 ResonanceState result = 00282 fromResonance(invalid, p, *iter, 00283 signalVertex, hardProcess); 00284 switch(result) { 00285 case kNo: 00286 break; 00287 case kDirect: 00288 if ((*iter)->pdg_id() == id) 00289 return kDirect; 00290 if (!isExcludedFromResonances(id)) 00291 break; 00292 case kIndirect: 00293 return kIndirect; 00294 } 00295 } 00296 00297 return kNo; 00298 }
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 iter, lhef::partIdx(), and v.
00230 { 00231 unsigned int idx = partIdx(p, particle); 00232 00233 if (invalid[idx]) 00234 return false; 00235 00236 const HepMC::GenVertex *v = particle->production_vertex(); 00237 if (!v) 00238 return false; 00239 else if (v == signalVertex) 00240 return true; 00241 00242 for(HepMC::GenVertex::particles_in_const_iterator iter = 00243 v->particles_in_const_begin(); 00244 iter != v->particles_in_const_end(); ++iter) 00245 if (fromSignalVertex(invalid, p, *iter, signalVertex)) 00246 return true; 00247 00248 return false; 00249 }
const std::vector<unsigned int>& lhef::JetInput::getExcludedFromResonances | ( | ) | const [inline] |
Definition at line 35 of file JetInput.h.
References excludedFromResonances.
00036 { return excludedFromResonances; }
const std::vector<unsigned int>& lhef::JetInput::getExcludedResonances | ( | ) | const [inline] |
Definition at line 33 of file JetInput.h.
References excludedResonances.
00034 { return excludedResonances; }
bool lhef::JetInput::getExcludeResonances | ( | ) | const [inline] |
Definition at line 29 of file JetInput.h.
References excludeResonances.
00029 { return excludeResonances; }
bool lhef::JetInput::getHardProcessOnly | ( | ) | const [inline] |
Definition at line 27 of file JetInput.h.
References onlyHardProcess.
00027 { return onlyHardProcess; }
const std::vector<unsigned int>& lhef::JetInput::getIgnoredParticles | ( | ) | const [inline] |
Definition at line 31 of file JetInput.h.
References ignoreParticleIDs.
00032 { return ignoreParticleIDs; }
bool lhef::JetInput::getPartonicFinalState | ( | ) | const [inline] |
Definition at line 26 of file JetInput.h.
References partonicFinalState.
00026 { return partonicFinalState; }
double lhef::JetInput::getPtMin | ( | ) | const [inline] |
bool lhef::JetInput::getTausAndJets | ( | ) | const [inline] |
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()().
00197 { 00198 return testPartonChildren(invalid, p, particle->end_vertex()) > 0; 00199 }
Definition at line 113 of file JetInput.cc.
References excludedFromResonances, and lhef::isContained().
Referenced by fromResonance().
00114 { 00115 if (excludedFromResonances.empty()) 00116 return true; 00117 00118 return isContained(excludedFromResonances, pdgId); 00119 }
Definition at line 82 of file JetInput.cc.
Referenced by testPartonChildren().
00083 { 00084 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000; 00085 return (pdgId > 100 && pdgId < 900) || 00086 (pdgId > 1000 && pdgId < 9000); 00087 }
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()().
00123 { 00124 std::vector<const HepMC::GenVertex*>::const_iterator pos = 00125 std::lower_bound(hardProcess.begin(), 00126 hardProcess.end(), vertex); 00127 return pos != hardProcess.end() && *pos == vertex; 00128 }
Definition at line 108 of file JetInput.cc.
References ignoreParticleIDs, and lhef::isContained().
Referenced by fromResonance(), and operator()().
00109 { 00110 return isContained(ignoreParticleIDs, pdgId); 00111 }
Definition at line 73 of file JetInput.cc.
References tausAsJets.
Referenced by fromResonance(), operator()(), and testPartonChildren().
00074 { 00075 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000; 00076 return (pdgId > 0 && pdgId < 6) || pdgId == 7 || 00077 pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21; 00078 // tops are not considered "regular" partons 00079 // but taus eventually are (since they may hadronize later) 00080 }
Definition at line 97 of file JetInput.cc.
References excludedResonances, and lhef::isContained().
Referenced by fromResonance().
00098 { 00099 if (excludedResonances.empty()) { 00100 // gauge bosons and tops 00101 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000; 00102 return (pdgId > 21 && pdgId <= 39) || pdgId == 6 || pdgId == 8; 00103 } 00104 00105 return isContained(excludedResonances, pdgId); 00106 }
JetInput::ParticleVector lhef::JetInput::operator() | ( | const HepMC::GenEvent * | event | ) | const |
Definition at line 300 of file JetInput.cc.
References lat::endl(), excludeResonances, fromHardProcess(), fromResonance(), hasPartonChildren(), i, reco::method::invalid, lhef::invalidateTree(), isHardProcess(), isIgnored(), isParton(), iter, onlyHardProcess, partonicFinalState, ptMin, pv, HLT_VtxMuL3::result, size, python::multivaluedict::sort(), v, and dimuonsSequences_cff::veto.
00302 { 00303 if (!event->signal_process_vertex()) 00304 throw cms::Exception("InvalidHepMCEvent") 00305 << "HepMC event is lacking signal vertex." 00306 << std::endl; 00307 00308 VertexVector hardProcess, toLookAt; 00309 std::pair<HepMC::GenParticle*,HepMC::GenParticle*> beamParticles = 00310 event->beam_particles(); 00311 toLookAt.push_back(event->signal_process_vertex()); 00312 while(!toLookAt.empty()) { 00313 std::vector<const HepMC::GenVertex*> next; 00314 for(std::vector<const HepMC::GenVertex*>::const_iterator v = 00315 toLookAt.begin(); v != toLookAt.end(); ++v) { 00316 if (!*v || isHardProcess(*v, hardProcess)) 00317 continue; 00318 00319 bool veto = false; 00320 for(HepMC::GenVertex::particles_in_const_iterator iter = 00321 (*v)->particles_in_const_begin(); 00322 iter != (*v)->particles_in_const_end(); ++iter) { 00323 if (*iter == beamParticles.first || 00324 *iter == beamParticles.second) { 00325 veto = true; 00326 break; 00327 } 00328 } 00329 if (veto) 00330 continue; 00331 00332 hardProcess.push_back(*v); 00333 std::sort(hardProcess.begin(), hardProcess.end()); 00334 00335 for(HepMC::GenVertex::particles_in_const_iterator iter = 00336 (*v)->particles_in_const_begin(); 00337 iter != (*v)->particles_in_const_end(); ++iter) { 00338 const HepMC::GenVertex *pv = 00339 (*iter)->production_vertex(); 00340 if (pv) 00341 next.push_back(pv); 00342 } 00343 } 00344 00345 toLookAt = next; 00346 std::sort(toLookAt.begin(), toLookAt.end()); 00347 } 00348 00349 ParticleVector particles; 00350 for(HepMC::GenEvent::particle_const_iterator iter = event->particles_begin(); 00351 iter != event->particles_end(); ++iter) 00352 particles.push_back(*iter); 00353 00354 std::sort(particles.begin(), particles.end()); 00355 unsigned int size = particles.size(); 00356 00357 ParticleBitmap selected(size, false); 00358 ParticleBitmap invalid(size, false); 00359 00360 for(unsigned int i = 0; i < size; i++) { 00361 const HepMC::GenParticle *particle = particles[i]; 00362 if (invalid[i]) 00363 continue; 00364 if (particle->status() == 1) 00365 selected[i] = true; 00366 00367 if (partonicFinalState && isParton(particle->pdg_id())) { 00368 if (!particle->end_vertex() && 00369 particle->status() != 1) 00370 // some brokenness in event... 00371 invalid[i] = true; 00372 else if (!hasPartonChildren(invalid, particles, 00373 particle)) { 00374 selected[i] = true; 00375 invalidateTree(invalid, particles, 00376 particle->end_vertex()); 00377 } 00378 } 00379 00380 if (onlyHardProcess && 00381 !fromHardProcess(invalid, particles, 00382 particle, hardProcess) && 00383 !isHardProcess(particle->end_vertex(), hardProcess)) 00384 invalid[i] = true; 00385 } 00386 00387 ParticleVector result; 00388 for(unsigned int i = 0; i < size; i++) { 00389 const HepMC::GenParticle *particle = particles[i]; 00390 if (!selected[i] || invalid[i]) 00391 continue; 00392 00393 if (excludeResonances && 00394 fromResonance(invalid, particles, particle, 00395 event->signal_process_vertex(), 00396 hardProcess)) { 00397 invalid[i] = true; 00398 continue; 00399 } 00400 00401 if (isIgnored(particle->pdg_id())) 00402 continue; 00403 00404 if (particle->momentum().perp() >= ptMin) 00405 result.push_back(particle); 00406 } 00407 00408 return result; 00409 }
Definition at line 66 of file JetInput.cc.
References excludedFromResonances, and python::multivaluedict::sort().
Referenced by JetInput().
00068 { 00069 excludedFromResonances = particleIDs; 00070 std::sort(excludedFromResonances.begin(), excludedFromResonances.end()); 00071 }
Definition at line 58 of file JetInput.cc.
References excludedResonances, setExcludeResonances(), and python::multivaluedict::sort().
Referenced by JetInput().
00060 { 00061 setExcludeResonances(true); 00062 excludedResonances = particleIDs; 00063 std::sort(excludedResonances.begin(), excludedResonances.end()); 00064 }
Definition at line 43 of file JetInput.h.
References excludeResonances.
Referenced by setExcludedResonances().
00043 { excludeResonances = flag; }
Definition at line 40 of file JetInput.h.
References onlyHardProcess.
00041 { onlyHardProcess = flag; }
Definition at line 51 of file JetInput.cc.
References ignoreParticleIDs, and python::multivaluedict::sort().
Referenced by JetInput().
00053 { 00054 ignoreParticleIDs = particleIDs; 00055 std::sort(ignoreParticleIDs.begin(), ignoreParticleIDs.end()); 00056 }
Definition at line 38 of file JetInput.h.
References partonicFinalState.
00039 { partonicFinalState = flag; }
void lhef::JetInput::setPtMin | ( | double | ptMin | ) | [inline] |
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 isHadron(), isParton(), iter, lhef::partIdx(), and HLT_VtxMuL3::result.
Referenced by hasPartonChildren().
00168 { 00169 if (!v) 00170 return 0; 00171 00172 for(HepMC::GenVertex::particles_out_const_iterator iter = 00173 v->particles_out_const_begin(); 00174 iter != v->particles_out_const_end(); ++iter) { 00175 unsigned int idx = partIdx(p, *iter); 00176 00177 if (invalid[idx]) 00178 continue; 00179 00180 if (isParton((*iter)->pdg_id())) 00181 return 1; 00182 if (isHadron((*iter)->pdg_id())) 00183 return -1; 00184 00185 const HepMC::GenVertex *v = (*iter)->end_vertex(); 00186 int result = testPartonChildren(invalid, p, v); 00187 if (result) 00188 return result; 00189 } 00190 00191 return 0; 00192 }
std::vector<unsigned int> lhef::JetInput::excludedFromResonances [private] |
Definition at line 87 of file JetInput.h.
Referenced by fromResonance(), getExcludedFromResonances(), isExcludedFromResonances(), and setExcludedFromResonances().
std::vector<unsigned int> lhef::JetInput::excludedResonances [private] |
Definition at line 86 of file JetInput.h.
Referenced by fromResonance(), getExcludedResonances(), isResonance(), and setExcludedResonances().
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] |
bool lhef::JetInput::tausAsJets [private] |
Definition at line 91 of file JetInput.h.
Referenced by getTausAndJets(), isParton(), and setTausAsJets().