CMS 3D CMS Logo

Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes

InputGenJetsParticleSelector Class Reference

#include <InputGenJetsParticleSelector.h>

Inheritance diagram for InputGenJetsParticleSelector:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Types

typedef std::vector< bool > ParticleBitmap
typedef std::vector< const
reco::GenParticle * > 
ParticleVector
enum  ResonanceState { kNo = 0, kDirect, kIndirect }

Public Member Functions

ResonanceState fromResonance (ParticleBitmap &invalid, const ParticleVector &p, const reco::GenParticle *particle) const
bool getExcludeResonances () 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 reco::GenParticle *particle) const
 InputGenJetsParticleSelector (const edm::ParameterSet &)
bool isIgnored (int pdgId) const
bool isParton (int pdgId) const
virtual void produce (edm::Event &evt, const edm::EventSetup &evtSetup)
void setExcludeResonances (bool flag=true)
void setPartonicFinalState (bool flag=true)
void setPtMin (double ptMin)
void setTausAsJets (bool flag=true)
 ~InputGenJetsParticleSelector ()

Static Public Member Functions

static bool isHadron (int pdgId)
static bool isResonance (int pdgId)

Private Member Functions

 InputGenJetsParticleSelector ()
bool isExcludedFromResonance (int pdgId) const
void setExcludeFromResonancePids (const std::vector< unsigned int > &particleIDs)
void setIgnoredParticles (const std::vector< unsigned int > &particleIDs)
int testPartonChildren (ParticleBitmap &invalid, const ParticleVector &p, const reco::GenParticle *particle) const

Private Attributes

std::vector< unsigned int > excludeFromResonancePids
bool excludeResonances
std::vector< unsigned int > ignoreParticleIDs
edm::InputTag inTag
bool partonicFinalState
double ptMin
bool tausAsJets

Detailed Description

Definition at line 28 of file InputGenJetsParticleSelector.h.


Member Typedef Documentation

Definition at line 31 of file InputGenJetsParticleSelector.h.

Definition at line 32 of file InputGenJetsParticleSelector.h.


Member Enumeration Documentation

Enumerator:
kNo 
kDirect 
kIndirect 

Definition at line 62 of file InputGenJetsParticleSelector.h.

                      {
    kNo = 0,
    kDirect,
    kIndirect
  };

Constructor & Destructor Documentation

InputGenJetsParticleSelector::InputGenJetsParticleSelector ( const edm::ParameterSet params)

Definition at line 42 of file InputGenJetsParticleSelector.cc.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), setExcludeFromResonancePids(), and setIgnoredParticles().

                                                                                        :
  inTag(params.getParameter<edm::InputTag>("src")),
  partonicFinalState(params.getParameter<bool>("partonicFinalState")),
  excludeResonances(params.getParameter<bool>("excludeResonances")),
  tausAsJets(params.getParameter<bool>("tausAsJets")),
  ptMin(0.0){
  if (params.exists("ignoreParticleIDs"))
    setIgnoredParticles(params.getParameter<std::vector<unsigned int> >
                        ("ignoreParticleIDs"));
  setExcludeFromResonancePids(params.getParameter<std::vector<unsigned int> >
                        ("excludeFromResonancePids"));

  produces <reco::GenParticleRefVector> ();
      
}
InputGenJetsParticleSelector::~InputGenJetsParticleSelector ( )

Definition at line 58 of file InputGenJetsParticleSelector.cc.

{}
InputGenJetsParticleSelector::InputGenJetsParticleSelector ( ) [inline, private]

Definition at line 76 of file InputGenJetsParticleSelector.h.

{} //should not be used!

Member Function Documentation

InputGenJetsParticleSelector::ResonanceState InputGenJetsParticleSelector::fromResonance ( ParticleBitmap invalid,
const ParticleVector p,
const reco::GenParticle particle 
) const

Definition at line 172 of file InputGenJetsParticleSelector.cc.

References i, UserOptions_cff::idx, isExcludedFromResonance(), isIgnored(), isParton(), isResonance(), kDirect, kIndirect, kNo, reco::CompositeRefCandidateT< D >::mother(), reco::CompositeRefCandidateT< D >::numberOfMothers(), partIdx(), reco::LeafCandidate::pdgId(), benchmark_cfg::pdgId, query::result, and reco::LeafCandidate::status().

Referenced by produce().

{
  unsigned int idx = partIdx(p, particle);
  int id = particle->pdgId();

  if (invalid[idx]) return kIndirect;
      
  if (isResonance(id) && particle->status() == 3){
    return kDirect;
  }

  
  if (!isIgnored(id) && (isParton(id)))
    return kNo;
  
  
  
  unsigned int nMo=particle->numberOfMothers();
  if (!nMo)
    return kNo;
  

  for(unsigned int i=0;i<nMo;++i){
    ResonanceState result = fromResonance(invalid,p,dynamic_cast<const reco::GenParticle*>(particle->mother(i)));
    switch(result) {
    case kNo:
      break;
    case kDirect:
      if (dynamic_cast<const reco::GenParticle*>(particle->mother(i))->pdgId()==id || isResonance(id))
        return kDirect;
      if(!isExcludedFromResonance(id))
        break;
    case kIndirect:
      return kIndirect;
    }
  }
return kNo;
}
bool InputGenJetsParticleSelector::getExcludeResonances ( ) const [inline]

Definition at line 41 of file InputGenJetsParticleSelector.h.

References excludeResonances.

{ return excludeResonances; }
const std::vector<unsigned int>& InputGenJetsParticleSelector::getIgnoredParticles ( ) const [inline]

Definition at line 44 of file InputGenJetsParticleSelector.h.

References ignoreParticleIDs.

    { return ignoreParticleIDs; }
bool InputGenJetsParticleSelector::getPartonicFinalState ( ) const [inline]

Definition at line 40 of file InputGenJetsParticleSelector.h.

References partonicFinalState.

{ return partonicFinalState; }
double InputGenJetsParticleSelector::getPtMin ( ) const [inline]

Definition at line 43 of file InputGenJetsParticleSelector.h.

References ptMin.

{ return ptMin; }
bool InputGenJetsParticleSelector::getTausAndJets ( ) const [inline]

Definition at line 42 of file InputGenJetsParticleSelector.h.

References tausAsJets.

{ return tausAsJets; }
bool InputGenJetsParticleSelector::hasPartonChildren ( ParticleBitmap invalid,
const ParticleVector p,
const reco::GenParticle particle 
) const

Definition at line 214 of file InputGenJetsParticleSelector.cc.

References testPartonChildren().

Referenced by produce().

                                                                                            {
  return testPartonChildren(invalid, p, particle) > 0;
}
bool InputGenJetsParticleSelector::isExcludedFromResonance ( int  pdgId) const [private]

Definition at line 104 of file InputGenJetsParticleSelector.cc.

References excludeFromResonancePids, benchmark_cfg::pdgId, and pos.

Referenced by fromResonance().

{
  pdgId = pdgId > 0 ? pdgId : -pdgId;
  std::vector<unsigned int>::const_iterator pos =
    std::lower_bound(excludeFromResonancePids.begin(),
                     excludeFromResonancePids.end(),
                     (unsigned int)pdgId);
  return pos != excludeFromResonancePids.end() && *pos == (unsigned int)pdgId;
 
}
bool InputGenJetsParticleSelector::isHadron ( int  pdgId) [static]

Definition at line 80 of file InputGenJetsParticleSelector.cc.

References benchmark_cfg::pdgId.

{
  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
  return (pdgId > 100 && pdgId < 900) ||
    (pdgId > 1000 && pdgId < 9000);
}
bool InputGenJetsParticleSelector::isIgnored ( int  pdgId) const

Definition at line 94 of file InputGenJetsParticleSelector.cc.

References ignoreParticleIDs, benchmark_cfg::pdgId, and pos.

Referenced by fromResonance(), and produce().

{
  pdgId = pdgId > 0 ? pdgId : -pdgId;
  std::vector<unsigned int>::const_iterator pos =
    std::lower_bound(ignoreParticleIDs.begin(),
                     ignoreParticleIDs.end(),
                     (unsigned int)pdgId);
  return pos != ignoreParticleIDs.end() && *pos == (unsigned int)pdgId;
}
bool InputGenJetsParticleSelector::isParton ( int  pdgId) const

Definition at line 72 of file InputGenJetsParticleSelector.cc.

References benchmark_cfg::pdgId, and tausAsJets.

Referenced by fromResonance(), and produce().

                                                          {
  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
  return (pdgId > 0 && pdgId < 6) || pdgId == 7 ||
    pdgId == 9 || (tausAsJets && pdgId == 15) || pdgId == 21;
  // tops are not considered "regular" partons
  // but taus eventually are (since they may hadronize later)
}
bool InputGenJetsParticleSelector::isResonance ( int  pdgId) [static]

Definition at line 87 of file InputGenJetsParticleSelector.cc.

References benchmark_cfg::pdgId.

Referenced by fromResonance().

{
  // gauge bosons and tops
  pdgId = (pdgId > 0 ? pdgId : -pdgId) % 10000;
  return (pdgId > 21 && pdgId <= 42) || pdgId == 6 || pdgId == 8 ;  //BUG! war 21. 22=gamma..
}
void InputGenJetsParticleSelector::produce ( edm::Event evt,
const edm::EventSetup evtSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 222 of file InputGenJetsParticleSelector.cc.

References prof2calltree::count, excludeResonances, fromResonance(), genParticleCandidates2GenParticles_cfi::genParticles, edm::Event::getByLabel(), hasPartonChildren(), i, UserOptions_cff::idx, inTag, align::invalid, invalidateTree(), isIgnored(), isParton(), reco::CompositeRefCandidateT< D >::numberOfDaughters(), partonicFinalState, reco::LeafCandidate::pdgId(), reco::LeafCandidate::pt(), ptMin, edm::Event::put(), findQualityFiles::size, python::multivaluedict::sort(), and reco::LeafCandidate::status().

                                                                                       {

 
  std::auto_ptr<reco::GenParticleRefVector> selected_ (new reco::GenParticleRefVector);
    
  edm::Handle<reco::GenParticleCollection> genParticles;
  //  evt.getByLabel("genParticles", genParticles );
  evt.getByLabel(inTag, genParticles );
    
    
  ParticleVector particles;
  for (reco::GenParticleCollection::const_iterator iter=genParticles->begin();iter!=genParticles->end();++iter){
    particles.push_back(&*iter); 
  }
  
  std::sort(particles.begin(), particles.end());
  unsigned int size = particles.size();
  
  ParticleBitmap selected(size, false);
  ParticleBitmap invalid(size, false);

  for(unsigned int i = 0; i < size; i++) {
    const reco::GenParticle *particle = particles[i];
    if (invalid[i])
      continue;
    if (particle->status() == 1)
      selected[i] = true;
    if (partonicFinalState && isParton(particle->pdgId())) {
          
      if (particle->numberOfDaughters()==0 &&
          particle->status() != 1) {
        // some brokenness in event...
        invalid[i] = true;
      }
      else if (!hasPartonChildren(invalid, particles,
                                  particle)) {
        selected[i] = true;
        invalidateTree(invalid, particles,particle); //this?!?
      }
    }
        
  }
 unsigned int count=0;
  for(size_t idx=0;idx<genParticles->size();++idx){ 
    const reco::GenParticle *particle = particles[idx];
    if (!selected[idx] || invalid[idx]){
      continue;
    }
        
    if (excludeResonances &&
        fromResonance(invalid, particles, particle)) {
      invalid[idx] = true;
      //cout<<"[INPUTSELECTOR] Invalidates FROM RESONANCE!: ["<<setw(4)<<idx<<"] "<<particle->pdgId()<<" "<<particle->pt()<<endl;
      continue;
    }
        
    if (isIgnored(particle->pdgId())){
      continue;
    }

   
    if (particle->pt() >= ptMin){
      edm::Ref<reco::GenParticleCollection> particleRef(genParticles,idx);
      selected_->push_back(particleRef);
      //cout<<"Finally we have: ["<<setw(4)<<idx<<"] "<<setw(4)<<particle->pdgId()<<" "<<particle->pt()<<endl;
      count++;
    }
  }
  evt.put(selected_);
}
void InputGenJetsParticleSelector::setExcludeFromResonancePids ( const std::vector< unsigned int > &  particleIDs) [private]
void InputGenJetsParticleSelector::setExcludeResonances ( bool  flag = true) [inline]

Definition at line 49 of file InputGenJetsParticleSelector.h.

References excludeResonances.

void InputGenJetsParticleSelector::setIgnoredParticles ( const std::vector< unsigned int > &  particleIDs) [private]
void InputGenJetsParticleSelector::setPartonicFinalState ( bool  flag = true) [inline]

Definition at line 47 of file InputGenJetsParticleSelector.h.

References partonicFinalState.

void InputGenJetsParticleSelector::setPtMin ( double  ptMin) [inline]

Definition at line 52 of file InputGenJetsParticleSelector.h.

References ptMin.

{ this->ptMin = ptMin; }
void InputGenJetsParticleSelector::setTausAsJets ( bool  flag = true) [inline]

Definition at line 51 of file InputGenJetsParticleSelector.h.

References tausAsJets.

int InputGenJetsParticleSelector::testPartonChildren ( InputGenJetsParticleSelector::ParticleBitmap invalid,
const ParticleVector p,
const reco::GenParticle particle 
) const [private]

Definition at line 148 of file InputGenJetsParticleSelector.cc.

References reco::CompositeRefCandidateT< D >::daughter(), i, UserOptions_cff::idx, npart, reco::CompositeRefCandidateT< D >::numberOfDaughters(), partIdx(), and query::result.

Referenced by hasPartonChildren().

{
  unsigned int npart=particle->numberOfDaughters();
  if (!npart) {return 0;}

  for (unsigned int i=0;i<npart;++i){
    unsigned int idx = partIdx(p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
    if (invalid[idx])
      continue;
    if (isParton((particle->daughter(i)->pdgId()))){
      return 1;
    }
    if (isHadron((particle->daughter(i)->pdgId()))){
      return -1;
    }
    int result = testPartonChildren(invalid,p,dynamic_cast<const reco::GenParticle*>(particle->daughter(i)));
    if (result) return result;
  }
  return 0;
}

Member Data Documentation

std::vector<unsigned int> InputGenJetsParticleSelector::excludeFromResonancePids [private]
std::vector<unsigned int> InputGenJetsParticleSelector::ignoreParticleIDs [private]

Definition at line 78 of file InputGenJetsParticleSelector.h.

Referenced by produce().

Definition at line 92 of file InputGenJetsParticleSelector.h.

Referenced by getPtMin(), produce(), and setPtMin().

Definition at line 91 of file InputGenJetsParticleSelector.h.

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