CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TtFullHadHypothesis.cc
Go to the documentation of this file.
3 
6  jetsToken_(consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
7  lightQ_(0), lightQBar_(0), b_(0), bBar_(0), lightP_(0), lightPBar_(0)
8 {
9  getMatch_ = false;
10  if( cfg.exists("match") ) {
11  getMatch_ = true;
12  matchToken_ = consumes<std::vector<std::vector<int> > >(cfg.getParameter<edm::InputTag>("match"));
13  }
14  if( cfg.exists("jetCorrectionLevel") ) {
15  jetCorrectionLevel_ = cfg.getParameter<std::string>("jetCorrectionLevel");
16  }
17  produces<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >();
18  produces<int>("Key");
19 }
20 
23 {
24  if( lightQ_ ) delete lightQ_;
25  if( lightQBar_ ) delete lightQBar_;
26  if( b_ ) delete b_;
27  if( bBar_ ) delete bBar_;
28  if( lightP_ ) delete lightP_;
29  if( lightPBar_ ) delete lightPBar_;
30 }
31 
33 void
35 {
37  evt.getByToken(jetsToken_, jets);
38 
39  std::vector<std::vector<int> > matchVec;
40  if( getMatch_ ) {
42  evt.getByToken(matchToken_, matchHandle);
43  matchVec = *matchHandle;
44  }
45  else {
46  std::vector<int> dummyMatch;
47  for(unsigned int i = 0; i < 4; ++i)
48  dummyMatch.push_back( -1 );
49  matchVec.push_back( dummyMatch );
50  }
51 
52  // declare auto_ptr for products
53  std::auto_ptr<std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > >
54  pOut( new std::vector<std::pair<reco::CompositeCandidate, std::vector<int> > > );
55  std::auto_ptr<int> pKey(new int);
56 
57  // go through given vector of jet combinations
58  unsigned int idMatch = 0;
59  typedef std::vector<std::vector<int> >::iterator MatchVecIterator;
60  for(MatchVecIterator match = matchVec.begin(); match != matchVec.end(); ++match) {
61  // reset pointers
63  // build hypothesis
64  buildHypo(evt, jets, *match, idMatch++);
65  pOut->push_back( std::make_pair(hypo(), *match) );
66  }
67  // feed out hyps and matches
68  evt.put(pOut);
69 
70  // build and feed out key
71  buildKey();
72  *pKey=key();
73  evt.put(pKey, "Key");
74 }
75 
77 void
79 {
80  lightQ_ = 0;
81  lightQBar_ = 0;
82  b_ = 0;
83  bBar_ = 0;
84  lightP_ = 0;
85  lightPBar_ = 0;
86 }
87 
91 {
92  // check for sanity of the hypothesis
93  if( !lightQ_ || !lightQBar_ || !b_ ||
94  !bBar_ || !lightP_ || !lightPBar_ )
95  return reco::CompositeCandidate();
96 
97  // setup transient references
98  reco::CompositeCandidate hyp, top, w, topBar, wBar;
99 
100  AddFourMomenta addFourMomenta;
101  // build up the top bar branch
104  addFourMomenta.set( wBar );
105  topBar.addDaughter( wBar, TtFullHadDaughter::WMinus );
107  addFourMomenta.set( topBar );
108 
109  // build up the top branch that decays hadronically
112  addFourMomenta.set( w );
115  addFourMomenta.set( top );
116 
117  // build ttbar hypotheses
118  hyp.addDaughter( topBar, TtFullHadDaughter::TopBar );
120  addFourMomenta.set( hyp );
121 
122  return hyp;
123 }
124 
128 {
129  // jetCorrectionLevel was not configured
130  if(jetCorrectionLevel_.empty())
131  throw cms::Exception("Configuration")
132  << "Unconfigured jetCorrectionLevel. Please use an appropriate, non-empty string.\n";
133 
134  // quarkType is unknown
135  if( !(quarkType=="wQuarkMix" ||
136  quarkType=="udsQuark" ||
137  quarkType=="cQuark" ||
138  quarkType=="bQuark") )
139  throw cms::Exception("Configuration")
140  << quarkType << " is unknown as a quarkType for the jetCorrectionLevel.\n";
141 
142  // combine correction level; start with a ':' even if
143  // there is no flavor tag to be added, as it is needed
144  // by setCandidate to disentangle the correction tag
145  // from a potential flavor tag, which can be empty
147  if( level=="L5Flavor:" || level=="L6UE:" || level=="L7Parton:" ){
148  if(quarkType=="wQuarkMix"){level+="wMix";}
149  if(quarkType=="udsQuark" ){level+="uds";}
150  if(quarkType=="cQuark" ){level+="charm";}
151  if(quarkType=="bQuark" ){level+="bottom";}
152  }
153  else{
154  level+="none";
155  }
156  return level;
157 }
158 
160 void
161 TtFullHadHypothesis::setCandidate(const edm::Handle<std::vector<pat::Jet> >& handle, const int& idx,
162  reco::ShallowClonePtrCandidate*& clone, const std::string& correctionLevel)
163 {
165  // disentangle the correction from the potential flavor tag
166  // by the separating ':'; the flavor tag can be empty though
167  std::string step = correctionLevel.substr(0,correctionLevel.find(":"));
168  std::string flavor = correctionLevel.substr(1+correctionLevel.find(":"));
169  float corrFactor = 1.;
170  if(flavor=="wMix")
171  corrFactor = 0.75*ptr->jecFactor(step, "uds") + 0.25*ptr->jecFactor(step, "charm");
172  else
173  corrFactor = ptr->jecFactor(step, flavor);
174  clone = new reco::ShallowClonePtrCandidate( ptr, ptr->charge(), ptr->p4()*corrFactor, ptr->vertex() );
175 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
int key() const
return key
void setCandidate(const edm::Handle< C > &handle, const int &idx, reco::ShallowClonePtrCandidate *&clone)
use one object in a collection to set a ShallowClonePtrCandidate
const double w
Definition: UKUtility.cc:23
virtual void produce(edm::Event &, const edm::EventSetup &)
produce the event hypothesis as CompositeCandidate and Key
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
static const std::string LightPBar
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::string jetCorrectionLevel(const std::string &quarkType)
helper function to construct the proper correction level string for corresponding quarkType ...
static const std::string LightQ
reco::ShallowClonePtrCandidate * lightPBar_
static const std::string Top
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
vector< PseudoJet > jets
static const std::string LightP
reco::ShallowClonePtrCandidate * lightQBar_
tuple handle
Definition: patZpeak.py:22
static const std::string WMinus
reco::ShallowClonePtrCandidate * lightP_
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
void addDaughter(const Candidate &, const std::string &s="")
add a clone of the passed candidate as daughter
edm::EDGetTokenT< std::vector< std::vector< int > > > matchToken_
std::string jetCorrectionLevel_
static const std::string WPlus
static const std::string BBar
static const std::string TopBar
string top
Definition: fff_deleter.py:272
virtual void buildKey()=0
build the event hypothesis key
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
reco::ShallowClonePtrCandidate * bBar_
void resetCandidates()
reset candidate pointers before hypo build process
~TtFullHadHypothesis()
default destructor
reco::ShallowClonePtrCandidate * b_
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:6
void set(reco::Candidate &c) const
set up a candidate
static const std::string LightQBar
tuple level
Definition: testEve_cfg.py:34
reco::CompositeCandidate hypo()
return event hypothesis
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
input label for all necessary collections
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
reco::ShallowClonePtrCandidate * lightQ_