CMS 3D CMS Logo

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

reco::modules::JetFlavourIdentifier Class Reference

Inheritance diagram for reco::modules::JetFlavourIdentifier:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Types

enum  DEFINITION_T {
  PHYSICS = 0, ALGO, NEAREST_STATUS2, NEAREST_STATUS3,
  HEAVIEST, N_DEFINITIONS, NULL_DEF
}

Public Member Functions

 JetFlavourIdentifier (const edm::ParameterSet &)
 ~JetFlavourIdentifier ()

Private Member Functions

void fillLeptons (const std::vector< const reco::Candidate * > &, JetFlavour::Leptons &, int, int)
std::vector< const
reco::Candidate * > 
findCandidates (const reco::Candidate *, int)
JetFlavour::Leptons findLeptons (const GenParticleRef &)
virtual void produce (edm::Event &, const edm::EventSetup &)

Static Private Member Functions

static int heaviestFlavour (int)

Private Attributes

DEFINITION_T definition
bool leptonInfo_
bool physDefinition
InputTag sourceByRefer_
math::XYZTLorentzVector thePartonLV
Handle
< JetMatchedPartonsCollection
theTagByRef

Detailed Description

Definition at line 85 of file JetFlavourIdentifier.cc.


Member Enumeration Documentation

Enumerator:
PHYSICS 
ALGO 
NEAREST_STATUS2 
NEAREST_STATUS3 
HEAVIEST 
N_DEFINITIONS 
NULL_DEF 

Definition at line 88 of file JetFlavourIdentifier.cc.


Constructor & Destructor Documentation

JetFlavourIdentifier::JetFlavourIdentifier ( const edm::ParameterSet iConfig)

Definition at line 116 of file JetFlavourIdentifier.cc.

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

{
    produces<JetFlavourMatchingCollection>();
    sourceByRefer_ = iConfig.getParameter<InputTag>("srcByReference");
    physDefinition = iConfig.getParameter<bool>("physicsDefinition");
    leptonInfo_ = iConfig.exists("leptonInfo") ? iConfig.getParameter<bool>("leptonInfo") : true;
    // If we have a definition of which parton to identify, use it,
    // otherwise we default to the "old" behavior of either "physics" or "algorithmic".
    // Furthermore, if the specified definition is not sensible for the given jet,
    // then the "physDefinition" switch is used to identify the flavour of the jet. 
    if ( iConfig.exists("definition") ) {
      definition = static_cast<DEFINITION_T>( iConfig.getParameter<int>("definition") );
    } else {
      definition = NULL_DEF;
    }
}
JetFlavourIdentifier::~JetFlavourIdentifier ( )

Definition at line 135 of file JetFlavourIdentifier.cc.

{
}

Member Function Documentation

void JetFlavourIdentifier::fillLeptons ( const std::vector< const reco::Candidate * > &  cands,
JetFlavour::Leptons leptons,
int  rank,
int  flavour 
) [private]

test for neutrinos because of conversions and dalitz pions

Definition at line 342 of file JetFlavourIdentifier.cc.

References abs, reco::JetFlavour::Leptons::electron, i, j, max(), reco::JetFlavour::Leptons::muon, benchmark_cfg::pdgId, and reco::JetFlavour::Leptons::tau.

{
  for(unsigned int j = 0; j < cands.size(); j++) {
    for(unsigned int i = 0; i < cands[j]->numberOfDaughters(); i++) {
      int pdgId = std::abs(cands[j]->daughter(i)->pdgId());

//      for(int z = 1; z <= rank; z *= 10) std::cout << " ------ ";
//      std::cout << pdgId << std::endl;

      if (pdgId == 12)
        leptons.electron += rank;
      else if (pdgId == 14)
        leptons.muon += rank;
      else if (pdgId == 16)
        leptons.tau += rank;
      else {
        int heaviest = heaviestFlavour(pdgId);
        int heaviest_ = heaviest < 10 ? heaviest : 0;
        if (!heaviest || (flavour < 4 ? (heaviest_ < 4) : (heaviest >= 4))) {
          std::vector<const reco::Candidate*> newcands = findCandidates(cands[j]->daughter(i), heaviest);
          if (pdgId <= 110) newcands.push_back(cands[j]->daughter(i));
          fillLeptons(newcands, leptons, rank * 10, std::max(heaviest_, flavour));
        }
      }
    }
  }
}
std::vector< const reco::Candidate * > JetFlavourIdentifier::findCandidates ( const reco::Candidate cand,
int  partonFlavour 
) [private]

Definition at line 306 of file JetFlavourIdentifier.cc.

References abs, filterCSVwithJSON::copy, reco::Candidate::daughter(), reco::flavour(), i, reco::Candidate::numberOfDaughters(), reco::Candidate::p4(), reco::Candidate::pdgId(), and benchmark_cfg::pdgId.

{
  std::vector<const reco::Candidate*> cands;
  
  for(unsigned int i = 0; i < cand->numberOfDaughters(); i++) {
/*
    std::cout << "DeltaR - " << partonFlavour << " ";
    if (DeltaR(thePartonLV, cand->daughter(i)->p4()) > 0.7) std::cout << "(";
    std::cout << cand->daughter(i)->pdgId() << ": " << DeltaR(thePartonLV, cand->daughter(i)->p4());
    if (DeltaR(thePartonLV, cand->daughter(i)->p4()) > 0.7) std::cout << ")";
    std::cout << std::endl;
*/

    if (DeltaR(thePartonLV, cand->daughter(i)->p4()) < 0.7) {
      int pdgId = std::abs(cand->daughter(i)->pdgId());
      int flavour = heaviestFlavour(pdgId);
      if (flavour == partonFlavour ||
          (flavour >= 10 && partonFlavour >= 10)) {
//        std::cout << "<------- " << std::endl;
        std::vector<const reco::Candidate*> newcands = findCandidates(cand->daughter(i), partonFlavour);
//        std::cout << " ------->" << std::endl;
        std::copy(newcands.begin(), newcands.end(), std::back_inserter(cands));
      }
      if (partonFlavour >= 10)
        cands.push_back(cand->daughter(i));
    }
  }

  if (cands.empty() && std::abs(cand->pdgId()) > 110 &&
      !(partonFlavour >= 4 && partonFlavour < 10 &&
        heaviestFlavour(cand->pdgId()) < 4))
    cands.push_back(cand);

  return cands;
}
JetFlavour::Leptons JetFlavourIdentifier::findLeptons ( const GenParticleRef parton) [private]

first daughter of the parton should be an MC particle (pdgId==92,93)

lookup particles with parton flavour and weak decay

count leptons of candidates

Definition at line 284 of file JetFlavourIdentifier.cc.

References abs.

{
  JetFlavour::Leptons theLeptons;

  thePartonLV = parton->p4();

  const reco::Candidate *mcstring = parton->daughter(0);
  int partonFlavour = std::abs(parton->pdgId());
//  std::cout << "parton DeltaR: " << DeltaR(thePartonLV, parton->p4()) << std::endl;

  std::vector<const reco::Candidate*> candidates = findCandidates(mcstring, partonFlavour);
//  std::cout << "Candidates are:" << std::endl;
//  for(unsigned int j = 0; j < candidates.size(); j++) std::cout << "   --> " << candidates[j]->pdgId() << std::endl;

  fillLeptons(candidates, theLeptons, 1, partonFlavour);

  return theLeptons;
}
int JetFlavourIdentifier::heaviestFlavour ( int  pdgId) [static, private]

Definition at line 371 of file JetFlavourIdentifier.cc.

References abs, reco::flavour(), and benchmark_cfg::pdgId.

{
  int flavour = 0;

  pdgId = std::abs(pdgId) % 100000;
  if (pdgId > 110) {
    while(pdgId % 10 > 0 && pdgId % 10 < 6) {
      pdgId /= 10;
      if (pdgId % 10 > flavour)
        flavour = pdgId % 10;
    }
  } else
    flavour = pdgId;

  return flavour;
}
void JetFlavourIdentifier::produce ( edm::Event iEvent,
const edm::EventSetup iEs 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 141 of file JetFlavourIdentifier.cc.

References lumi::ALGO, reco::MatchedPartons::algoDefinitionParton(), edm::Ref< C, T, F >::get(), edm::Event::getByLabel(), reco::MatchedPartons::heaviest(), edm::Ref< C, T, F >::isNonnull(), j, findQualityFiles::jj, reco::MatchedPartons::nearest_status2(), reco::MatchedPartons::nearest_status3(), sistrip::PHYSICS, reco::MatchedPartons::physicsDefinitionParton(), and edm::Event::put().

{
  // Get the JetMatchedPartons
  iEvent.getByLabel (sourceByRefer_, theTagByRef);

  // Create a JetFlavourMatchingCollection
  JetFlavourMatchingCollection *jfmc;
  if (!theTagByRef->empty()) {
    RefToBase<Jet> jj = theTagByRef->begin()->first;
    jfmc = new JetFlavourMatchingCollection(RefToBaseProd<Jet>(jj));
  } else {
    jfmc = new JetFlavourMatchingCollection();
  }
  auto_ptr<reco::JetFlavourMatchingCollection> jetFlavMatching(jfmc);

  // Loop over the matched partons and see which match. 
  for ( JetMatchedPartonsCollection::const_iterator j  = theTagByRef->begin();
                                                    j != theTagByRef->end();
                                                    j ++ ) {


    // Consider this match. 
    const MatchedPartons aMatch = (*j).second;

    // This will hold the 4-vector, vertex, flavour and the leptonian decays (0: no lepton, xyz: x leptons in layer 2, y in layer 1 and z in layer 0) of the requested object.
    math::XYZTLorentzVector thePartonLorentzVector(0,0,0,0);
    math::XYZPoint          thePartonVertex(0,0,0);
    int                     thePartonFlavour = 0;
    JetFlavour::Leptons     theLeptons;

    // get the partons based on which definition to use. 
    switch (definition) {
      case PHYSICS: {
        const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
        if (aPartPhy.isNonnull()) {
                thePartonLorentzVector = aPartPhy.get()->p4();
                thePartonVertex        = aPartPhy.get()->vertex();
                thePartonFlavour       = aPartPhy.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
        }
        break;
      }
      case ALGO: {
        const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
        if (aPartAlg.isNonnull()) {
                thePartonLorentzVector = aPartAlg.get()->p4();
                thePartonVertex        = aPartAlg.get()->vertex();
                thePartonFlavour       = aPartAlg.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
        }
        break;
      }
      case NEAREST_STATUS2 : {
        const GenParticleRef aPartN2 = aMatch.nearest_status2();
        if (aPartN2.isNonnull()) {
                thePartonLorentzVector = aPartN2.get()->p4();
                thePartonVertex        = aPartN2.get()->vertex();
                thePartonFlavour       = aPartN2.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartN2);
        }
        break;
      }
      case NEAREST_STATUS3: {
        const GenParticleRef aPartN3 = aMatch.nearest_status3();
        if (aPartN3.isNonnull()) {
                thePartonLorentzVector = aPartN3.get()->p4();
                thePartonVertex        = aPartN3.get()->vertex();
                thePartonFlavour       = aPartN3.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartN3);
        }
        break;
      }
      case HEAVIEST: {
        const GenParticleRef aPartHeaviest = aMatch.heaviest();
        if (aPartHeaviest.isNonnull()) {
                thePartonLorentzVector = aPartHeaviest.get()->p4();
                thePartonVertex        = aPartHeaviest.get()->vertex();
                thePartonFlavour       = aPartHeaviest.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartHeaviest);
        }
        break;
      }
      // Default case is backwards-compatible
      default:{
        if (physDefinition) {
          const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
          if (aPartPhy.isNonnull()) {
            thePartonLorentzVector = aPartPhy.get()->p4();
            thePartonVertex        = aPartPhy.get()->vertex();
            thePartonFlavour       = aPartPhy.get()->pdgId();
            if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
          }
        } else {
          const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
          if (aPartAlg.isNonnull()) {
            thePartonLorentzVector = aPartAlg.get()->p4();
            thePartonVertex        = aPartAlg.get()->vertex();
            thePartonFlavour       = aPartAlg.get()->pdgId();
            if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
          }
        }
      } break;
    } // end switch on definition

    // Now make sure we have a match. If the user specified "heaviest", for instance,
    // and there is no b- or c-quarks, then fall back to the "physDefinition" switch.
    
    if (thePartonFlavour == 0) {
      if (physDefinition) {
        const GenParticleRef aPartPhy = aMatch.physicsDefinitionParton();
        if (aPartPhy.isNonnull()) {
          thePartonLorentzVector = aPartPhy.get()->p4();
          thePartonVertex        = aPartPhy.get()->vertex();
          thePartonFlavour       = aPartPhy.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartPhy);
        }
      } else {
        const GenParticleRef aPartAlg = aMatch.algoDefinitionParton();
        if (aPartAlg.isNonnull()) {
          thePartonLorentzVector = aPartAlg.get()->p4();
          thePartonVertex        = aPartAlg.get()->vertex();
          thePartonFlavour       = aPartAlg.get()->pdgId();
          if (leptonInfo_) theLeptons = findLeptons(aPartAlg);
        }
      }
    }

/*
     std::cout << "Leptons of " <<thePartonFlavour << " Jet: " << std::endl;
     std::cout << "  electrons: " <<theLeptons.electron << std::endl;
     std::cout << "  muons    : " <<theLeptons.muon << std::endl;
     std::cout << "  tau      : " <<theLeptons.tau << std::endl;
*/
    // Add the jet->flavour match to the map. 
    (*jetFlavMatching)[(*j).first] = JetFlavour(thePartonLorentzVector, thePartonVertex, thePartonFlavour, theLeptons);
  }// end loop over jets


  // Put the object into the event. 
  iEvent.put(  jetFlavMatching );

}

Member Data Documentation

Definition at line 107 of file JetFlavourIdentifier.cc.

Definition at line 106 of file JetFlavourIdentifier.cc.

Definition at line 105 of file JetFlavourIdentifier.cc.

Definition at line 104 of file JetFlavourIdentifier.cc.

Definition at line 108 of file JetFlavourIdentifier.cc.

Definition at line 103 of file JetFlavourIdentifier.cc.