CMS 3D CMS Logo

TtFullHadHypothesis.cc
Go to the documentation of this file.
3 
6  : jetsToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))) {
7  getMatch_ = false;
8  if (cfg.exists("match")) {
9  getMatch_ = true;
10  matchToken_ = consumes<std::vector<std::vector<int> > >(cfg.getParameter<edm::InputTag>("match"));
11  }
12  if (cfg.exists("jetCorrectionLevel")) {
13  jetCorrectionLevel_ = cfg.getParameter<std::string>("jetCorrectionLevel");
14  }
15  produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
16  produces<int>("Key");
17 }
18 
23 
24  std::vector<std::vector<int> > matchVec;
25  if (getMatch_) {
27  evt.getByToken(matchToken_, matchHandle);
28  matchVec = *matchHandle;
29  } else {
30  std::vector<int> dummyMatch;
31  for (unsigned int i = 0; i < 4; ++i)
32  dummyMatch.push_back(-1);
33  matchVec.push_back(dummyMatch);
34  }
35 
36  // declare unique_ptr for products
37  std::unique_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > > pOut(
38  new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > >);
39  std::unique_ptr<int> pKey(new int);
40 
41  // go through given vector of jet combinations
42  unsigned int idMatch = 0;
43  typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
44  for (MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
45  // reset pointers
47  // build hypothesis
48  buildHypo(evt, jets, *match, idMatch++);
49  pOut->push_back(std::make_pair(hypo(), *match));
50  }
51  // feed out hyps and matches
52  evt.put(std::move(pOut));
53 
54  // build and feed out key
55  buildKey();
56  *pKey = key();
57  evt.put(std::move(pKey), "Key");
58 }
59 
62  lightQ_.reset();
63  lightQBar_.reset();
64  b_.reset();
65  bBar_.reset();
66  lightP_.reset();
67  lightPBar_.reset();
68 }
69 
72  // check for sanity of the hypothesis
73  if (!lightQ_ || !lightQBar_ || !b_ || !bBar_ || !lightP_ || !lightPBar_)
74  return reco::CompositeCandidate();
75 
76  // setup transient references
77  reco::CompositeCandidate hyp, top, w, topBar, wBar;
78 
79  AddFourMomenta addFourMomenta;
80  // build up the top bar branch
83  addFourMomenta.set(wBar);
86  addFourMomenta.set(topBar);
87 
88  // build up the top branch that decays hadronically
91  addFourMomenta.set(w);
94  addFourMomenta.set(top);
95 
96  // build ttbar hypotheses
99  addFourMomenta.set(hyp);
100 
101  return hyp;
102 }
103 
106  // jetCorrectionLevel was not configured
107  if (jetCorrectionLevel_.empty())
108  throw cms::Exception("Configuration")
109  << "Unconfigured jetCorrectionLevel. Please use an appropriate, non-empty string.\n";
110 
111  // quarkType is unknown
112  if (!(quarkType == "wQuarkMix" || quarkType == "udsQuark" || quarkType == "cQuark" || quarkType == "bQuark"))
113  throw cms::Exception("Configuration") << quarkType << " is unknown as a quarkType for the jetCorrectionLevel.\n";
114 
115  // combine correction level; start with a ':' even if
116  // there is no flavor tag to be added, as it is needed
117  // by makeCandidate to disentangle the correction tag
118  // from a potential flavor tag, which can be empty
120  if (level == "L5Flavor:" || level == "L6UE:" || level == "L7Parton:") {
121  if (quarkType == "wQuarkMix") {
122  level += "wMix";
123  }
124  if (quarkType == "udsQuark") {
125  level += "uds";
126  }
127  if (quarkType == "cQuark") {
128  level += "charm";
129  }
130  if (quarkType == "bQuark") {
131  level += "bottom";
132  }
133  } else {
134  level += "none";
135  }
136  return level;
137 }
138 
140 std::unique_ptr<reco::ShallowClonePtrCandidate> TtFullHadHypothesis::makeCandidate(
141  const edm::Handle<std::vector<pat::Jet> >& handle, const int& idx, const std::string& correctionLevel) {
143  // disentangle the correction from the potential flavor tag
144  // by the separating ':'; the flavor tag can be empty though
145  std::string step = correctionLevel.substr(0, correctionLevel.find(':'));
146  std::string flavor = correctionLevel.substr(1 + correctionLevel.find(':'));
147  float corrFactor = 1.;
148  if (flavor == "wMix")
149  corrFactor = 0.75 * ptr->jecFactor(step, "uds") + 0.25 * ptr->jecFactor(step, "charm");
150  else
151  corrFactor = ptr->jecFactor(step, flavor);
152  return std::make_unique<reco::ShallowClonePtrCandidate>(ptr, ptr->charge(), ptr->p4() * corrFactor, ptr->vertex());
153 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
float jecFactor(const std::string &level, const std::string &flavor="none", const std::string &set="") const
void set(reco::Candidate &c) const
set up a candidate
std::unique_ptr< reco::ShallowClonePtrCandidate > makeCandidate(const edm::Handle< C > &handle, const int &idx)
use one object in a collection to set a ShallowClonePtrCandidate
T w() const
static const std::string LightPBar
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
static const std::string LightQ
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:526
std::unique_ptr< reco::ShallowClonePtrCandidate > b_
Definition: HeavyIon.h:7
static const std::string Top
Definition: Jet.py:1
std::unique_ptr< reco::ShallowClonePtrCandidate > lightP_
static const std::string LightP
static const std::string WMinus
static const std::string B
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
std::unique_ptr< reco::ShallowClonePtrCandidate > bBar_
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQ_
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
void produce(edm::Event &, const edm::EventSetup &) override
produce the event hypothesis as CompositeCandidate and Key
edm::EDGetTokenT< std::vector< std::vector< int > > > matchToken_
std::string jetCorrectionLevel_
static const std::string WPlus
std::unique_ptr< reco::ShallowClonePtrCandidate > lightPBar_
static const std::string BBar
static const std::string TopBar
virtual void buildKey()=0
build the event hypothesis key
void resetCandidates()
reset candidate pointers before hypo build process
HLT enums.
int key() const
return key
TtFullHadHypothesis(const edm::ParameterSet &cfg)
default constructor
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
static const std::string LightQBar
step
Definition: StallMonitor.cc:83
reco::CompositeCandidate hypo()
return event hypothesis
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
input label for all necessary collections
def move(src, dest)
Definition: eostools.py:511
std::unique_ptr< reco::ShallowClonePtrCandidate > lightQBar_
int charge() const final
electric charge