CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Protected Attributes

TtFullHadHypothesis Class Reference

#include <TtFullHadHypothesis.h>

Inheritance diagram for TtFullHadHypothesis:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper TtFullHadHypGenMatch TtFullHadHypKinFit

List of all members.

Public Member Functions

 TtFullHadHypothesis (const edm::ParameterSet &cfg)
 default constructor
 ~TtFullHadHypothesis ()
 default destructor

Protected Member Functions

virtual void buildHypo (edm::Event &event, const edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &jetPartonAssociation, const unsigned int iComb)=0
 build event hypothesis from the reco objects of a full-hadronic event
virtual void buildKey ()=0
 build the event hypothesis key
reco::CompositeCandidate hypo ()
 return event hypothesis
bool isValid (const int &idx, const edm::Handle< std::vector< pat::Jet > > &jets)
 check if index is in valid range of selected jets
std::string jetCorrectionLevel (const std::string &quarkType)
 helper function to construct the proper correction level string for corresponding quarkType
int key () const
 return key
virtual void produce (edm::Event &, const edm::EventSetup &)
 produce the event hypothesis as CompositeCandidate and Key
void resetCandidates ()
 reset candidate pointers before hypo build process
void setCandidate (const edm::Handle< std::vector< pat::Jet > > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone, const std::string &correctionLevel)
 use one object in a jet collection to set a ShallowClonePtrCandidate with proper jet corrections
template<typename C >
void setCandidate (const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
 use one object in a collection to set a ShallowClonePtrCandidate

Protected Attributes

reco::ShallowClonePtrCandidateb_
reco::ShallowClonePtrCandidatebBar_
bool getMatch_
std::string jetCorrectionLevel_
edm::InputTag jets_
 input label for all necessary collections
int key_
 hypothesis key (to be set by the buildKey function)
reco::ShallowClonePtrCandidatelightP_
reco::ShallowClonePtrCandidatelightPBar_
reco::ShallowClonePtrCandidatelightQ_
reco::ShallowClonePtrCandidatelightQBar_
edm::InputTag match_

Detailed Description

Definition at line 29 of file TtFullHadHypothesis.h.


Constructor & Destructor Documentation

TtFullHadHypothesis::TtFullHadHypothesis ( const edm::ParameterSet cfg) [explicit]

default constructor

Definition at line 5 of file TtFullHadHypothesis.cc.

References edm::ParameterSet::exists(), getMatch_, edm::ParameterSet::getParameter(), jetCorrectionLevel_, and match_.

                                                                  :
  jets_(cfg.getParameter<edm::InputTag>("jets")),
  lightQ_(0), lightQBar_(0), b_(0), bBar_(0), lightP_(0), lightPBar_(0)
{
  getMatch_ = false;
  if( cfg.exists("match") ) {
    getMatch_ = true;
    match_ = cfg.getParameter<edm::InputTag>("match");
  }
  if( cfg.exists("jetCorrectionLevel") ) {
    jetCorrectionLevel_ = cfg.getParameter<std::string>("jetCorrectionLevel");
  }
  produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
  produces<int>("Key");
}
TtFullHadHypothesis::~TtFullHadHypothesis ( )

default destructor

Definition at line 22 of file TtFullHadHypothesis.cc.

References b_, bBar_, lightP_, lightPBar_, lightQ_, and lightQBar_.

{
  if( lightQ_    ) delete lightQ_;
  if( lightQBar_ ) delete lightQBar_;
  if( b_         ) delete b_;
  if( bBar_      ) delete bBar_;
  if( lightP_    ) delete lightP_;
  if( lightPBar_ ) delete lightPBar_;
}

Member Function Documentation

virtual void TtFullHadHypothesis::buildHypo ( edm::Event event,
const edm::Handle< std::vector< pat::Jet > > &  jets,
std::vector< int > &  jetPartonAssociation,
const unsigned int  iComb 
) [protected, pure virtual]

build event hypothesis from the reco objects of a full-hadronic event

Implemented in TtFullHadHypGenMatch, and TtFullHadHypKinFit.

Referenced by produce().

virtual void TtFullHadHypothesis::buildKey ( ) [protected, pure virtual]

build the event hypothesis key

Implemented in TtFullHadHypGenMatch, and TtFullHadHypKinFit.

Referenced by produce().

reco::CompositeCandidate TtFullHadHypothesis::hypo ( ) [protected]

return event hypothesis

Definition at line 90 of file TtFullHadHypothesis.cc.

References reco::CompositeCandidate::addDaughter(), TtFullHadDaughter::B, b_, TtFullHadDaughter::BBar, bBar_, TtFullHadDaughter::LightP, lightP_, TtFullHadDaughter::LightPBar, lightPBar_, TtFullHadDaughter::LightQ, lightQ_, TtFullHadDaughter::LightQBar, lightQBar_, AddFourMomenta::set(), TtFullHadDaughter::Top, TtFullHadDaughter::TopBar, TtFullHadDaughter::WMinus, and TtFullHadDaughter::WPlus.

Referenced by produce().

{
  // check for sanity of the hypothesis
  if( !lightQ_ || !lightQBar_ || !b_ || 
      !bBar_ || !lightP_ || !lightPBar_ )
    return reco::CompositeCandidate();
  
  // setup transient references
  reco::CompositeCandidate hyp, top, w, topBar, wBar;

  AddFourMomenta addFourMomenta;  
  // build up the top bar branch
  wBar  .addDaughter(*lightP_,    TtFullHadDaughter::LightP    );
  wBar  .addDaughter(*lightPBar_, TtFullHadDaughter::LightPBar );
  addFourMomenta.set( wBar );
  topBar.addDaughter( wBar,  TtFullHadDaughter::WMinus );
  topBar.addDaughter(*bBar_, TtFullHadDaughter::BBar   );
  addFourMomenta.set( topBar );
  
  // build up the top branch that decays hadronically
  w  .addDaughter(*lightQ_,    TtFullHadDaughter::LightQ    );
  w  .addDaughter(*lightQBar_, TtFullHadDaughter::LightQBar );
  addFourMomenta.set( w );
  top.addDaughter( w,  TtFullHadDaughter::WPlus );
  top.addDaughter(*b_, TtFullHadDaughter::B     );
  addFourMomenta.set( top );

  // build ttbar hypotheses
  hyp.addDaughter( topBar, TtFullHadDaughter::TopBar );
  hyp.addDaughter( top,    TtFullHadDaughter::Top    );
  addFourMomenta.set( hyp );

  return hyp;
}
bool TtFullHadHypothesis::isValid ( const int &  idx,
const edm::Handle< std::vector< pat::Jet > > &  jets 
) [inline, protected]

check if index is in valid range of selected jets

Definition at line 55 of file TtFullHadHypothesis.h.

References analyzePatCleaning_cfg::jets.

Referenced by TtFullHadHypGenMatch::buildHypo().

{ return (0<=idx && idx<(int)jets->size()); };
std::string TtFullHadHypothesis::jetCorrectionLevel ( const std::string &  quarkType) [protected]

helper function to construct the proper correction level string for corresponding quarkType

helper function to construct the proper correction level string for corresponding quarkType, for unknown quarkTypes an empty string is returned

Definition at line 127 of file TtFullHadHypothesis.cc.

References jetCorrectionLevel_, and testEve_cfg::level.

Referenced by TtFullHadHypGenMatch::buildHypo().

{
  // jetCorrectionLevel was not configured
  if(jetCorrectionLevel_.empty())
    throw cms::Exception("Configuration")
      << "Unconfigured jetCorrectionLevel. Please use an appropriate, non-empty string.\n";

  // quarkType is unknown
  if( !(quarkType=="wQuarkMix" ||
        quarkType=="udsQuark" ||
        quarkType=="cQuark" ||
        quarkType=="bQuark") )
    throw cms::Exception("Configuration")
      << quarkType << " is unknown as a quarkType for the jetCorrectionLevel.\n";

  // combine correction level; start with a ':' even if 
  // there is no flavor tag to be added, as it is needed
  // by setCandidate to disentangle the correction tag 
  // from a potential flavor tag, which can be empty
  std::string level=jetCorrectionLevel_+":";
  if( level=="had:" || level=="ue:" || level=="part:" ){
    if(quarkType=="wQuarkMix"){level+="wMix";}
    if(quarkType=="udsQuark" ){level+="uds";}
    if(quarkType=="cQuark"   ){level+="charm";}
    if(quarkType=="bQuark"   ){level+="bottom";}
  }
  else{
    level+="none";
  }
  return level;
}
int TtFullHadHypothesis::key ( ) const [inline, protected]

return key

Definition at line 51 of file TtFullHadHypothesis.h.

References key_.

Referenced by produce().

{ return key_; };
void TtFullHadHypothesis::produce ( edm::Event evt,
const edm::EventSetup setup 
) [protected, virtual]

produce the event hypothesis as CompositeCandidate and Key

Implements edm::EDProducer.

Definition at line 34 of file TtFullHadHypothesis.cc.

References buildHypo(), buildKey(), edm::Event::getByLabel(), getMatch_, hypo(), i, analyzePatCleaning_cfg::jets, jets_, key(), match(), match_, edm::Event::put(), and resetCandidates().

{
  edm::Handle<std::vector<pat::Jet> > jets;
  evt.getByLabel(jets_, jets);
  
  std::vector<std::vector<int> > matchVec;
  if( getMatch_ ) {
    edm::Handle<std::vector<std::vector<int> > > matchHandle;
    evt.getByLabel(match_, matchHandle);
    matchVec = *matchHandle;
  }
  else {
    std::vector<int> dummyMatch;
    for(unsigned int i = 0; i < 4; ++i) 
      dummyMatch.push_back( -1 );
    matchVec.push_back( dummyMatch );
  }

  // declare auto_ptr for products
  std::auto_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >
    pOut( new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > );
  std::auto_ptr<int> pKey(new int);

  // go through given vector of jet combinations
  unsigned int idMatch = 0;
  typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
  for(MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
    // reset pointers
    resetCandidates();
    // build hypothesis
    buildHypo(evt, jets, *match, idMatch++);
    pOut->push_back( std::make_pair(hypo(), *match) );
  }
  // feed out hyps and matches
  evt.put(pOut);

  // build and feed out key
  buildKey();
  *pKey=key();
  evt.put(pKey, "Key");
}
void TtFullHadHypothesis::resetCandidates ( ) [protected]

reset candidate pointers before hypo build process

Definition at line 78 of file TtFullHadHypothesis.cc.

References b_, bBar_, lightP_, lightPBar_, lightQ_, and lightQBar_.

Referenced by produce().

{
  lightQ_    = 0;
  lightQBar_ = 0;
  b_         = 0;
  bBar_      = 0;
  lightP_    = 0;
  lightPBar_ = 0;
}
void TtFullHadHypothesis::setCandidate ( const edm::Handle< std::vector< pat::Jet > > &  handle,
const int &  idx,
reco::ShallowClonePtrCandidate *&  clone,
const std::string &  correctionLevel 
) [protected]

use one object in a jet collection to set a ShallowClonePtrCandidate with proper jet corrections

Definition at line 161 of file TtFullHadHypothesis.cc.

References patZpeak::handle, and launcher::step.

{
  edm::Ptr<pat::Jet> ptr = edm::Ptr<pat::Jet>(handle, idx);
  // disentangle the correction from the potential flavor tag 
  // by the separating ':'; the flavor tag can be empty though
  std::string step   = correctionLevel.substr(0,correctionLevel.find(":"));
  std::string flavor = correctionLevel.substr(1+correctionLevel.find(":"));
  float corrFactor = 1.;
  if(flavor=="wMix")
    corrFactor = 0.75*ptr->jecFactor(step, "uds") + 0.25*ptr->jecFactor(step, "charm");
  else
    corrFactor = ptr->jecFactor(step, flavor);
  clone = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), ptr->p4()*corrFactor, ptr->vertex() );
}
template<typename C >
void TtFullHadHypothesis::setCandidate ( const edm::Handle< C > &  handle,
const int &  idx,
reco::ShallowClonePtrCandidate *&  clone 
) [protected]

use one object in a collection to set a ShallowClonePtrCandidate

Definition at line 96 of file TtFullHadHypothesis.h.

References patZpeak::handle.

Referenced by TtFullHadHypGenMatch::buildHypo(), and TtFullHadHypKinFit::buildHypo().

                                                                                                                  {
  typedef typename C::value_type O;
  edm::Ptr<O> ptr = edm::Ptr<O>(handle, idx);
  clone = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), ptr->p4(), ptr->vertex() );
}

Member Data Documentation

internal check whether the match information exists or not, if false a blind dummy match vector will be used internally

Definition at line 73 of file TtFullHadHypothesis.h.

Referenced by produce(), and TtFullHadHypothesis().

specify the desired jet correction level (the default should be L3Absolute-'abs')

Definition at line 79 of file TtFullHadHypothesis.h.

Referenced by jetCorrectionLevel(), and TtFullHadHypothesis().

input label for all necessary collections

Definition at line 75 of file TtFullHadHypothesis.h.

Referenced by produce().

int TtFullHadHypothesis::key_ [protected]

hypothesis key (to be set by the buildKey function)

Definition at line 81 of file TtFullHadHypothesis.h.

Referenced by TtFullHadHypGenMatch::buildKey(), TtFullHadHypKinFit::buildKey(), and key().

candidates for internal use for the creation of the hypothesis candidate

Definition at line 84 of file TtFullHadHypothesis.h.

Referenced by TtFullHadHypGenMatch::buildHypo(), TtFullHadHypKinFit::buildHypo(), hypo(), resetCandidates(), and ~TtFullHadHypothesis().

Definition at line 76 of file TtFullHadHypothesis.h.

Referenced by produce(), and TtFullHadHypothesis().