CMS 3D CMS Logo

TtFullHadHypKinFit.cc
Go to the documentation of this file.
2 
3 
5  TtFullHadHypothesis( cfg ),
6  statusToken_ (consumes<std::vector<int> > (cfg.getParameter<edm::InputTag>("status" ))),
7  lightQToken_ (consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightQTag" ))),
8  lightQBarToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightQBarTag"))),
9  bToken_ (consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("bTag" ))),
10  bBarToken_ (consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("bBarTag" ))),
11  lightPToken_ (consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightPTag" ))),
12  lightPBarToken_(consumes<std::vector<pat::Particle> >(cfg.getParameter<edm::InputTag>("lightPBarTag")))
13 {
14 }
15 
17 
18 void
20  const edm::Handle<std::vector<pat::Jet> >& jets,
21  std::vector<int>& match, const unsigned int iComb)
22 {
24  evt.getByToken(statusToken_, status);
25  if( (*status)[iComb] != 0 ){
26  // create empty hypothesis if kinematic fit did not converge
27  return;
28  }
29 
36 
37  evt.getByToken(lightQToken_ , lightQ );
38  evt.getByToken(lightQBarToken_, lightQBar);
39  evt.getByToken(bToken_ , b );
40  evt.getByToken(bBarToken_ , bBar );
41  evt.getByToken(lightPToken_ , lightP );
42  evt.getByToken(lightPBarToken_, lightPBar);
43 
44  // -----------------------------------------------------
45  // add jets
46  // -----------------------------------------------------
47  if( !( lightQ->empty() || lightQBar->empty() || b->empty() || bBar->empty() ||
48  lightP->empty() || lightPBar->empty() ) ) {
49  setCandidate(lightQ , iComb, lightQ_ );
50  setCandidate(lightQBar, iComb, lightQBar_);
51  setCandidate(b , iComb, b_ );
52  setCandidate(bBar , iComb, bBar_ );
53  setCandidate(lightP , iComb, lightP_ );
54  setCandidate(lightPBar, iComb, lightPBar_);
55  }
56 }
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
edm::EDGetTokenT< std::vector< pat::Particle > > lightQToken_
edm::EDGetTokenT< std::vector< int > > statusToken_
reco::ShallowClonePtrCandidate * lightPBar_
Definition: HeavyIon.h:7
edm::EDGetTokenT< std::vector< pat::Particle > > bBarToken_
edm::EDGetTokenT< std::vector< pat::Particle > > lightQBarToken_
vector< PseudoJet > jets
reco::ShallowClonePtrCandidate * lightQBar_
reco::ShallowClonePtrCandidate * lightP_
edm::EDGetTokenT< std::vector< pat::Particle > > lightPToken_
TtFullHadHypKinFit(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< pat::Particle > > bToken_
virtual void buildHypo(edm::Event &, const edm::Handle< std::vector< pat::Jet > > &, std::vector< int > &, const unsigned int iComb)
build event hypothesis from the reco objects of a full-hadronic event
double b
Definition: hdecay.h:120
reco::ShallowClonePtrCandidate * bBar_
HLT enums.
reco::ShallowClonePtrCandidate * b_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
reco::ShallowClonePtrCandidate * lightQ_
edm::EDGetTokenT< std::vector< pat::Particle > > lightPBarToken_