CMS 3D CMS Logo

Public Member Functions | Private Attributes

PartonHadronDecayGenEvtSelector Class Reference

#include <PartonHadronDecayGenEvtSelector.h>

Inheritance diagram for PartonHadronDecayGenEvtSelector:
BaseHiGenEvtSelector

List of all members.

Public Member Functions

bool filter (HepMC::GenEvent *)
 PartonHadronDecayGenEvtSelector (const edm::ParameterSet &pset)
bool selectParticle (HepMC::GenParticle *par, int status, int pdg, double ptMin, double etaMax)
bool selectParticle (HepMC::GenParticle *par, int status, int pdg, double etaMax, double etaMin, double pMin, double ptMax, double ptMin)
virtual ~PartonHadronDecayGenEvtSelector ()

Private Attributes

double decayEtaMax_
double decayEtaMin_
int decayId_
int decayNtrig_
double decayPMin_
double decayPtMax_
double decayPtMin_
int decayStatus_
std::vector< double > hadronEtaMax_
std::vector< double > hadronEtaMin_
std::vector< int > hadronId_
std::vector< double > hadronPMin_
std::vector< double > hadronPtMax_
std::vector< double > hadronPtMin_
std::vector< int > hadronStatus_
std::vector< double > partonEtaMax_
std::vector< int > partonId_
std::vector< double > partonPtMin_
std::vector< int > partonStatus_

Detailed Description

Definition at line 7 of file PartonHadronDecayGenEvtSelector.h.


Constructor & Destructor Documentation

PartonHadronDecayGenEvtSelector::PartonHadronDecayGenEvtSelector ( const edm::ParameterSet pset)

Definition at line 6 of file PartonHadronDecayGenEvtSelector.cc.

References decayEtaMax_, decayEtaMin_, decayId_, decayNtrig_, decayPMin_, decayPtMax_, decayPtMin_, decayStatus_, Exception, edm::ParameterSet::getParameter(), hadronEtaMax_, hadronEtaMin_, hadronId_, hadronPMin_, hadronPtMax_, hadronPtMin_, hadronStatus_, edm::errors::LogicError, partonEtaMax_, partonId_, partonPtMin_, partonStatus_, and ptmin.

                                                                                            : BaseHiGenEvtSelector(pset)
{

  hadronId_       = pset.getParameter<vector<int> >("hadrons");
  hadronStatus_   = pset.getParameter<vector<int> >("hadronStatus");
  hadronEtaMax_   = pset.getParameter<vector<double> >("hadronEtaMax");
  hadronEtaMin_   = pset.getParameter<vector<double> >("hadronEtaMin");
  hadronPMin_     = pset.getParameter<vector<double> >("hadronPMin");
  hadronPtMax_    = pset.getParameter<vector<double> >("hadronPtMax");
  hadronPtMin_    = pset.getParameter<vector<double> >("hadronPtMin");
  
  decayId_        = pset.getParameter<int>("decays");
  decayStatus_    = pset.getParameter<int>("decayStatus");
  decayEtaMax_    = pset.getParameter<double>("decayEtaMax");
  decayEtaMin_    = pset.getParameter<double>("decayEtaMin");
  decayPMin_      = pset.getParameter<double>("decayPMin");
  decayPtMax_     = pset.getParameter<double>("decayPtMax");
  decayPtMin_     = pset.getParameter<double>("decayPtMin");
  decayNtrig_     = pset.getParameter<int>("decayNtrig");

  partonId_       = pset.getParameter<vector<int> >("partons");
  partonStatus_   = pset.getParameter<vector<int> >("partonStatus");
  partonEtaMax_   = pset.getParameter<vector<double> >("partonEtaMax");
  partonPtMin_    = pset.getParameter<vector<double> >("partonPtMin");
     
  int id     = hadronId_.size();
  int st     = hadronStatus_.size();
  int etamax = hadronEtaMax_.size();
  int etamin = hadronEtaMin_.size();
  int pmin   = hadronPMin_.size();
  int ptmax  = hadronPtMax_.size();
  int ptmin  = hadronPtMin_.size();
  
  if( id!=st || id!=etamax || id!=etamin || id!=ptmax || id!=ptmin || id!=pmin)
    {
      throw edm::Exception(edm::errors::LogicError)<<"Hadron selection parameters: "<<id<<st<<etamax<<etamin<<pmin<<ptmax<<ptmin<<endl;
    }
  

  id     = partonId_.size();
  st     = partonStatus_.size();
  etamax = partonEtaMax_.size();
  ptmin  = partonPtMin_.size();
  
  if( id!=st || id!=etamax || id!=ptmin )
    {
      throw edm::Exception(edm::errors::LogicError)<<"Parton selection parameters: "<<id<<st<<etamax<<ptmin<<endl;
    }
  
}
virtual PartonHadronDecayGenEvtSelector::~PartonHadronDecayGenEvtSelector ( ) [inline, virtual]

Definition at line 11 of file PartonHadronDecayGenEvtSelector.h.

{;}

Member Function Documentation

bool PartonHadronDecayGenEvtSelector::filter ( HepMC::GenEvent *  evt) [virtual]

Reimplemented from BaseHiGenEvtSelector.

Definition at line 59 of file PartonHadronDecayGenEvtSelector.cc.

References begin, decayEtaMax_, decayEtaMin_, decayId_, decayNtrig_, decayPMin_, decayPtMax_, decayPtMin_, decayStatus_, end, hadronEtaMax_, hadronEtaMin_, hadronId_, hadronPMin_, hadronPtMax_, hadronPtMin_, hadronStatus_, i, partonEtaMax_, partonId_, partonPtMin_, partonStatus_, and selectParticle().

{
  // loop over HepMC event, and search for  products of interest

  HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
  HepMC::GenEvent::particle_const_iterator end   = evt->particles_end();
  
  bool foundHadron   = false;
  bool foundDecay    = false;
  bool foundParton   = false;
  
  HepMC::GenEvent::particle_const_iterator it = begin;
  while( !foundParton && it != end )
    {
      for(unsigned i = 0; i < partonId_.size(); ++i)
        {
          if( selectParticle(*it,
                             partonStatus_[i], partonId_[i],  
                             partonPtMin_[i],partonEtaMax_[i]) ) foundParton = true;
        }
      ++it;
    }

  int foundtrig = 0;
  HepMC::GenEvent::particle_const_iterator it2 = begin;
 
  if(foundParton)
    {
      while( (!foundHadron || !foundDecay) && it2 != end )
        {
          
          for(unsigned i = 0; i < hadronId_.size(); ++i)
            {
              if( selectParticle(*it2, 
                                 hadronStatus_[i], hadronId_[i], 
                                 hadronEtaMax_[i],hadronEtaMin_[i], 
                                 hadronPMin_[i],
                                 hadronPtMax_[i],hadronPtMin_[i]) ) foundHadron = true;
            }
        
          if( selectParticle(*it2, 
                             decayStatus_, decayId_, 
                             decayEtaMax_,decayEtaMin_, 
                             decayPMin_,
                             decayPtMax_,decayPtMin_) ) foundtrig++;
          if(decayNtrig_ == foundtrig) foundDecay = true;
           
          ++it2;
        }
    }
  
  return (foundHadron && foundDecay && foundParton);
}
bool PartonHadronDecayGenEvtSelector::selectParticle ( HepMC::GenParticle *  par,
int  status,
int  pdg,
double  etaMax,
double  etaMin,
double  pMin,
double  ptMax,
double  ptMin 
) [inline]

Definition at line 17 of file PartonHadronDecayGenEvtSelector.h.

References abs, and PtMinSelector_cfg::ptMin.

                                                                                                                                          {
    return (par->status() == status &&
            abs(par->pdg_id()) == pdg &&
            par->momentum().eta() < etaMax && par->momentum().eta() > etaMin &&
            par->momentum().rho() > pMin &&
            par->momentum().perp() < ptMax && par->momentum().perp() > ptMin);
  }
bool PartonHadronDecayGenEvtSelector::selectParticle ( HepMC::GenParticle *  par,
int  status,
int  pdg,
double  ptMin,
double  etaMax 
) [inline]

Reimplemented from BaseHiGenEvtSelector.

Definition at line 14 of file PartonHadronDecayGenEvtSelector.h.

References abs.

Referenced by filter().

                                                                                                {
    return (par->status() == status && abs(par->pdg_id()) == pdg && par->momentum().perp() > ptMin && fabs(par->momentum().eta()) < etaMax);
  }

Member Data Documentation

Definition at line 42 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 43 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 40 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 47 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 44 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 45 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 46 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 41 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::hadronEtaMax_ [private]

Definition at line 34 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::hadronEtaMin_ [private]

Definition at line 35 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<int> PartonHadronDecayGenEvtSelector::hadronId_ [private]

Definition at line 32 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::hadronPMin_ [private]

Definition at line 36 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::hadronPtMax_ [private]

Definition at line 37 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::hadronPtMin_ [private]

Definition at line 38 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 33 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::partonEtaMax_ [private]

Definition at line 29 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<int> PartonHadronDecayGenEvtSelector::partonId_ [private]

Definition at line 27 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

std::vector<double> PartonHadronDecayGenEvtSelector::partonPtMin_ [private]

Definition at line 30 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().

Definition at line 28 of file PartonHadronDecayGenEvtSelector.h.

Referenced by filter(), and PartonHadronDecayGenEvtSelector().