#include <RecoJets/JetProducers/src/InputGenJetsParticleSelector.h>
Definition at line 29 of file InputGenJetsParticleSelector.h.
typedef std::vector<bool> InputGenJetsParticleSelector::ParticleBitmap |
Definition at line 32 of file InputGenJetsParticleSelector.h.
typedef std::vector<const reco::GenParticle*> InputGenJetsParticleSelector::ParticleVector |
Definition at line 33 of file InputGenJetsParticleSelector.h.
InputGenJetsParticleSelector::InputGenJetsParticleSelector | ( | const edm::ParameterSet & | params | ) |
Definition at line 41 of file InputGenJetsParticleSelector.cc.
References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), setExcludeFromResonancePids(), and setIgnoredParticles().
00041 : 00042 inTag(params.getParameter<edm::InputTag>("src")), 00043 partonicFinalState(params.getParameter<bool>("partonicFinalState")), 00044 excludeResonances(params.getParameter<bool>("excludeResonances")), 00045 tausAsJets(params.getParameter<bool>("tausAsJets")), 00046 ptMin(0.0){ 00047 if (params.exists("ignoreParticleIDs")) 00048 setIgnoredParticles(params.getParameter<std::vector<unsigned int> > 00049 ("ignoreParticleIDs")); 00050 setExcludeFromResonancePids(params.getParameter<std::vector<unsigned int> > 00051 ("excludeFromResonancePids")); 00052 00053 produces <reco::GenParticleRefVector> (); 00054 00055 }
InputGenJetsParticleSelector::~InputGenJetsParticleSelector | ( | ) |
InputGenJetsParticleSelector::InputGenJetsParticleSelector | ( | ) | [inline, private] |
InputGenJetsParticleSelector::ResonanceState InputGenJetsParticleSelector::fromResonance | ( | ParticleBitmap & | invalid, | |
const ParticleVector & | p, | |||
const reco::GenParticle * | particle | |||
) | const |
Definition at line 171 of file InputGenJetsParticleSelector.cc.
References i, isExcludedFromResonance(), isIgnored(), isParton(), isResonance(), kDirect, kIndirect, kNo, reco::CompositeRefCandidateT< D >::mother(), reco::CompositeRefCandidateT< D >::numberOfMothers(), partIdx(), reco::Particle::pdgId(), HLT_VtxMuL3::result, and reco::Particle::status().
Referenced by produce().
00174 { 00175 unsigned int idx = partIdx(p, particle); 00176 int id = particle->pdgId(); 00177 00178 if (invalid[idx]) return kIndirect; 00179 00180 if (isResonance(id) && particle->status() == 3){ 00181 return kDirect; 00182 } 00183 00184 00185 if (!isIgnored(id) && (isParton(id))) 00186 return kNo; 00187 00188 00189 00190 unsigned int nMo=particle->numberOfMothers(); 00191 if (!nMo) 00192 return kNo; 00193 00194 00195 for(unsigned int i=0;i<nMo;++i){ 00196 ResonanceState result = fromResonance(invalid,p,dynamic_cast<const reco::GenParticle*>(particle->mother(i))); 00197 switch(result) { 00198 case kNo: 00199 break; 00200 case kDirect: 00201 if (dynamic_cast<const reco::GenParticle*>(particle->mother(i))->pdgId()==id || isResonance(id)) 00202 return kDirect; 00203 if(!isExcludedFromResonance(id)) 00204 break; 00205 case kIndirect: 00206 return kIndirect; 00207 } 00208 } 00209 return kNo; 00210 }
bool InputGenJetsParticleSelector::getExcludeResonances | ( | ) | const [inline] |
Definition at line 42 of file InputGenJetsParticleSelector.h.
References excludeResonances.
00042 { return excludeResonances; }
const std::vector<unsigned int>& InputGenJetsParticleSelector::getIgnoredParticles | ( | ) | const [inline] |
Definition at line 45 of file InputGenJetsParticleSelector.h.
References ignoreParticleIDs.
00046 { return ignoreParticleIDs; }
bool InputGenJetsParticleSelector::getPartonicFinalState | ( | ) | const [inline] |
Definition at line 41 of file InputGenJetsParticleSelector.h.
References partonicFinalState.
00041 { return partonicFinalState; }
double InputGenJetsParticleSelector::getPtMin | ( | ) | const [inline] |
Definition at line 44 of file InputGenJetsParticleSelector.h.
References ptMin.
00044 { return ptMin; }
bool InputGenJetsParticleSelector::getTausAndJets | ( | ) | const [inline] |
Definition at line 43 of file InputGenJetsParticleSelector.h.
References tausAsJets.
00043 { return tausAsJets; }
bool InputGenJetsParticleSelector::hasPartonChildren | ( | ParticleBitmap & | invalid, | |
const ParticleVector & | p, | |||
const reco::GenParticle * | particle | |||
) | const |
Definition at line 213 of file InputGenJetsParticleSelector.cc.
References testPartonChildren().
Referenced by produce().
00215 { 00216 return testPartonChildren(invalid, p, particle) > 0; 00217 }
Definition at line 103 of file InputGenJetsParticleSelector.cc.
References excludeFromResonancePids, and int.
Referenced by fromResonance().
00104 { 00105 pdgId = pdgId > 0 ? pdgId : -pdgId; 00106 std::vector<unsigned int>::const_iterator pos = 00107 std::lower_bound(excludeFromResonancePids.begin(), 00108 excludeFromResonancePids.end(), 00109 (unsigned int)pdgId); 00110 return pos != excludeFromResonancePids.end() && *pos == (unsigned int)pdgId; 00111 00112 }
Definition at line 79 of file InputGenJetsParticleSelector.cc.
Referenced by testPartonChildren().
00080 { 00081 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000; 00082 return (pdgId > 100 && pdgId < 900) || 00083 (pdgId > 1000 && pdgId < 9000); 00084 }
Definition at line 93 of file InputGenJetsParticleSelector.cc.
References ignoreParticleIDs, and int.
Referenced by fromResonance(), and produce().
00094 { 00095 pdgId = pdgId > 0 ? pdgId : -pdgId; 00096 std::vector<unsigned int>::const_iterator pos = 00097 std::lower_bound(ignoreParticleIDs.begin(), 00098 ignoreParticleIDs.end(), 00099 (unsigned int)pdgId); 00100 return pos != ignoreParticleIDs.end() && *pos == (unsigned int)pdgId; 00101 }
Definition at line 71 of file InputGenJetsParticleSelector.cc.
References tausAsJets.
Referenced by fromResonance(), produce(), and testPartonChildren().
00071 { 00072 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000; 00073 return (pdgId > 0 && pdgId < 6) || pdgId == 7 || 00074 pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21; 00075 // tops are not considered "regular" partons 00076 // but taus eventually are (since they may hadronize later) 00077 }
Definition at line 86 of file InputGenJetsParticleSelector.cc.
Referenced by fromResonance().
00087 { 00088 // gauge bosons and tops 00089 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000; 00090 return (pdgId > 21 && pdgId <= 42) || pdgId == 6 || pdgId == 8 ; //BUG! war 21. 22=gamma.. 00091 }
void InputGenJetsParticleSelector::produce | ( | edm::Event & | evt, | |
const edm::EventSetup & | evtSetup | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 221 of file InputGenJetsParticleSelector.cc.
References count, excludeResonances, fromResonance(), genParticles_cfi::genParticles, edm::Event::getByLabel(), hasPartonChildren(), i, inTag, reco::method::invalid, invalidateTree(), isIgnored(), isParton(), iter, reco::CompositeRefCandidateT< D >::numberOfDaughters(), partonicFinalState, reco::Particle::pdgId(), reco::Particle::pt(), ptMin, edm::Event::put(), size, python::multivaluedict::sort(), and reco::Particle::status().
00221 { 00222 00223 00224 std::auto_ptr<reco::GenParticleRefVector> selected_ (new reco::GenParticleRefVector); 00225 00226 edm::Handle<reco::GenParticleCollection> genParticles; 00227 // evt.getByLabel("genParticles", genParticles ); 00228 evt.getByLabel(inTag, genParticles ); 00229 00230 00231 ParticleVector particles; 00232 for (reco::GenParticleCollection::const_iterator iter=genParticles->begin();iter!=genParticles->end();++iter){ 00233 particles.push_back(&*iter); 00234 } 00235 00236 std::sort(particles.begin(), particles.end()); 00237 unsigned int size = particles.size(); 00238 00239 ParticleBitmap selected(size, false); 00240 ParticleBitmap invalid(size, false); 00241 00242 for(unsigned int i = 0; i < size; i++) { 00243 const reco::GenParticle *particle = particles[i]; 00244 if (invalid[i]) 00245 continue; 00246 if (particle->status() == 1) 00247 selected[i] = true; 00248 if (partonicFinalState && isParton(particle->pdgId())) { 00249 00250 if (particle->numberOfDaughters()==0 && 00251 particle->status() != 1) { 00252 // some brokenness in event... 00253 invalid[i] = true; 00254 } 00255 else if (!hasPartonChildren(invalid, particles, 00256 particle)) { 00257 selected[i] = true; 00258 invalidateTree(invalid, particles,particle); //this?!? 00259 } 00260 } 00261 00262 } 00263 unsigned int count=0; 00264 for(size_t idx=0;idx<genParticles->size();++idx){ 00265 const reco::GenParticle *particle = particles[idx]; 00266 if (!selected[idx] || invalid[idx]){ 00267 continue; 00268 } 00269 00270 if (excludeResonances && 00271 fromResonance(invalid, particles, particle)) { 00272 invalid[idx] = true; 00273 //cout<<"[INPUTSELECTOR] Invalidates FROM RESONANCE!: ["<<setw(4)<<idx<<"] "<<particle->pdgId()<<" "<<particle->pt()<<endl; 00274 continue; 00275 } 00276 00277 if (isIgnored(particle->pdgId())){ 00278 continue; 00279 } 00280 00281 00282 if (particle->pt() >= ptMin){ 00283 edm::Ref<reco::GenParticleCollection> particleRef(genParticles,idx); 00284 selected_->push_back(particleRef); 00285 //cout<<"Finally we have: ["<<setw(4)<<idx<<"] "<<setw(4)<<particle->pdgId()<<" "<<particle->pt()<<endl; 00286 count++; 00287 } 00288 } 00289 evt.put(selected_); 00290 }
void InputGenJetsParticleSelector::setExcludeFromResonancePids | ( | const std::vector< unsigned int > & | particleIDs | ) | [private] |
Definition at line 65 of file InputGenJetsParticleSelector.cc.
References excludeFromResonancePids, and python::multivaluedict::sort().
Referenced by InputGenJetsParticleSelector().
00066 { 00067 excludeFromResonancePids = particleIDs; 00068 std::sort( excludeFromResonancePids.begin(), excludeFromResonancePids.end()); 00069 }
Definition at line 50 of file InputGenJetsParticleSelector.h.
References excludeResonances.
00051 { excludeResonances = flag; }
void InputGenJetsParticleSelector::setIgnoredParticles | ( | const std::vector< unsigned int > & | particleIDs | ) | [private] |
Definition at line 59 of file InputGenJetsParticleSelector.cc.
References ignoreParticleIDs, and python::multivaluedict::sort().
Referenced by InputGenJetsParticleSelector().
00060 { 00061 ignoreParticleIDs = particleIDs; 00062 std::sort(ignoreParticleIDs.begin(), ignoreParticleIDs.end()); 00063 }
Definition at line 48 of file InputGenJetsParticleSelector.h.
References partonicFinalState.
00049 { partonicFinalState = flag; }
void InputGenJetsParticleSelector::setPtMin | ( | double | ptMin | ) | [inline] |
Definition at line 52 of file InputGenJetsParticleSelector.h.
References tausAsJets.
00052 { tausAsJets = flag; }
int InputGenJetsParticleSelector::testPartonChildren | ( | InputGenJetsParticleSelector::ParticleBitmap & | invalid, | |
const ParticleVector & | p, | |||
const reco::GenParticle * | particle | |||
) | const [private] |
Definition at line 147 of file InputGenJetsParticleSelector.cc.
References reco::CompositeRefCandidateT< D >::daughter(), i, isHadron(), isParton(), npart, reco::CompositeRefCandidateT< D >::numberOfDaughters(), partIdx(), and HLT_VtxMuL3::result.
Referenced by hasPartonChildren().
00150 { 00151 unsigned int npart=particle->numberOfDaughters(); 00152 if (!npart) {return 0;} 00153 00154 for (unsigned int i=0;i<npart;++i){ 00155 unsigned int idx = partIdx(p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i))); 00156 if (invalid[idx]) 00157 continue; 00158 if (isParton((particle->daughter(i)->pdgId()))){ 00159 return 1; 00160 } 00161 if (isHadron((particle->daughter(i)->pdgId()))){ 00162 return -1; 00163 } 00164 int result = testPartonChildren(invalid,p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i))); 00165 if (result) return result; 00166 } 00167 return 0; 00168 }
std::vector<unsigned int> InputGenJetsParticleSelector::excludeFromResonancePids [private] |
Definition at line 85 of file InputGenJetsParticleSelector.h.
Referenced by isExcludedFromResonance(), and setExcludeFromResonancePids().
Definition at line 91 of file InputGenJetsParticleSelector.h.
Referenced by getExcludeResonances(), produce(), and setExcludeResonances().
std::vector<unsigned int> InputGenJetsParticleSelector::ignoreParticleIDs [private] |
Definition at line 84 of file InputGenJetsParticleSelector.h.
Referenced by getIgnoredParticles(), isIgnored(), and setIgnoredParticles().
Definition at line 90 of file InputGenJetsParticleSelector.h.
Referenced by getPartonicFinalState(), produce(), and setPartonicFinalState().
double InputGenJetsParticleSelector::ptMin [private] |
Definition at line 93 of file InputGenJetsParticleSelector.h.
Referenced by getPtMin(), and produce().
bool InputGenJetsParticleSelector::tausAsJets [private] |
Definition at line 92 of file InputGenJetsParticleSelector.h.
Referenced by getTausAndJets(), isParton(), and setTausAsJets().