00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033 #include "InputGenJetsParticleSelector.h"
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035
00036 #include <memory>
00037 #include "PhysicsTools/CandUtils/interface/pdgIdUtils.h"
00038
00039 using namespace std;
00040
00041 InputGenJetsParticleSelector::InputGenJetsParticleSelector(const edm::ParameterSet ¶ms ):
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 }
00056
00057 InputGenJetsParticleSelector::~InputGenJetsParticleSelector(){}
00058
00059 void InputGenJetsParticleSelector::setIgnoredParticles(const std::vector<unsigned int> &particleIDs)
00060 {
00061 ignoreParticleIDs = particleIDs;
00062 std::sort(ignoreParticleIDs.begin(), ignoreParticleIDs.end());
00063 }
00064
00065 void InputGenJetsParticleSelector::setExcludeFromResonancePids(const std::vector<unsigned int> &particleIDs)
00066 {
00067 excludeFromResonancePids = particleIDs;
00068 std::sort( excludeFromResonancePids.begin(), excludeFromResonancePids.end());
00069 }
00070
00071 bool InputGenJetsParticleSelector::isParton(int pdgId) const{
00072 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
00073 return (pdgId > 0 && pdgId < 6) || pdgId == 7 ||
00074 pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
00075
00076
00077 }
00078
00079 bool InputGenJetsParticleSelector::isHadron(int pdgId)
00080 {
00081 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
00082 return (pdgId > 100 && pdgId < 900) ||
00083 (pdgId > 1000 && pdgId < 9000);
00084 }
00085
00086 bool InputGenJetsParticleSelector::isResonance(int pdgId)
00087 {
00088
00089 pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
00090 return (pdgId > 21 && pdgId <= 42) || pdgId == 6 || pdgId == 8 ;
00091 }
00092
00093 bool InputGenJetsParticleSelector::isIgnored(int pdgId) const
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 }
00102
00103 bool InputGenJetsParticleSelector::isExcludedFromResonance(int pdgId) const
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 }
00113
00114 static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p,
00115
00116 const reco::GenParticle *particle)
00117 {
00118 InputGenJetsParticleSelector::ParticleVector::const_iterator pos =
00119 std::lower_bound(p.begin(), p.end(), particle);
00120 if (pos == p.end() || *pos != particle)
00121 throw cms::Exception("CorruptedData")
00122 << "reco::GenEvent corrupted: Unlisted particles"
00123 " in decay tree." << std::endl;
00124
00125 return pos - p.begin();
00126 }
00127
00128 static void invalidateTree(InputGenJetsParticleSelector::ParticleBitmap &invalid,
00129 const InputGenJetsParticleSelector::ParticleVector &p,
00130 const reco::GenParticle *particle)
00131 {
00132 unsigned int npart=particle->numberOfDaughters();
00133 if (!npart) return;
00134
00135 for (unsigned int i=0;i<npart;++i){
00136 unsigned int idx=partIdx(p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
00137 if (invalid[idx])
00138 continue;
00139 invalid[idx] = true;
00140
00141 invalidateTree(invalid, p, dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
00142 }
00143 }
00144
00145
00146 int InputGenJetsParticleSelector::testPartonChildren
00147 (InputGenJetsParticleSelector::ParticleBitmap &invalid,
00148 const InputGenJetsParticleSelector::ParticleVector &p,
00149 const reco::GenParticle *particle) const
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 }
00169
00170 InputGenJetsParticleSelector::ResonanceState
00171 InputGenJetsParticleSelector::fromResonance(ParticleBitmap &invalid,
00172 const ParticleVector &p,
00173 const reco::GenParticle *particle) const
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 }
00211
00212
00213 bool InputGenJetsParticleSelector::hasPartonChildren(ParticleBitmap &invalid,
00214 const ParticleVector &p,
00215 const reco::GenParticle *particle) const {
00216 return testPartonChildren(invalid, p, particle) > 0;
00217 }
00218
00219
00220
00221 void InputGenJetsParticleSelector::produce (edm::Event &evt, const edm::EventSetup &evtSetup){
00222
00223
00224 std::auto_ptr<reco::GenParticleRefVector> selected_ (new reco::GenParticleRefVector);
00225
00226 edm::Handle<reco::GenParticleCollection> genParticles;
00227
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
00253 invalid[i] = true;
00254 }
00255 else if (!hasPartonChildren(invalid, particles,
00256 particle)) {
00257 selected[i] = true;
00258 invalidateTree(invalid, particles,particle);
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
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
00286 count++;
00287 }
00288 }
00289 evt.put(selected_);
00290 }
00291
00292