CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Types | Private Member Functions | Private Attributes
TtSemiLepHitFitProducer< LeptonCollection > Class Template Reference

#include <TtSemiLepHitFitProducer.h>

Inheritance diagram for TtSemiLepHitFitProducer< LeptonCollection >:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Classes

struct  FitResult
 

Public Member Functions

 TtSemiLepHitFitProducer (const edm::ParameterSet &)
 
 ~TtSemiLepHitFitProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

typedef hitfit::RunHitFit
< pat::Electron, pat::Muon,
pat::Jet, pat::MET
PatHitFit
 

Private Member Functions

virtual void produce (edm::Event &, const edm::EventSetup &)
 

Private Attributes

std::string bTagAlgo_
 input tag for b-tagging algorithm More...
 
hitfit::LeptonTranslatorBase
< pat::Electron
electronTranslator_
 
PatHitFitHitFit
 
edm::FileInPath hitfitBJetResolution_
 
edm::FileInPath hitfitDefault_
 
edm::FileInPath hitfitElectronResolution_
 
edm::FileInPath hitfitMETResolution_
 
edm::FileInPath hitfitMuonResolution_
 
edm::FileInPath hitfitUdscJetResolution_
 
double jes_
 jet energy scale More...
 
double jesB_
 
std::string jetCorrectionLevel_
 jet correction level More...
 
edm::EDGetTokenT< std::vector
< pat::Jet > > 
jetsToken_
 
hitfit::JetTranslatorBase
< pat::Jet
jetTranslator_
 
edm::EDGetTokenT
< LeptonCollection > 
lepsToken_
 
double maxBTagValueNonBJet_
 max value of bTag for a non-b-jet More...
 
double maxEtaEle_
 maximum eta value for electrons, needed to limited range in which resolutions are provided More...
 
double maxEtaJet_
 maximum eta value for jets, needed to limited range in which resolutions are provided More...
 
double maxEtaMu_
 maximum eta value for muons, needed to limited range in which resolutions are provided More...
 
int maxNComb_
 maximal number of combinations to be written to the event More...
 
int maxNJets_
 maximal number of jets (-1 possible to indicate 'all') More...
 
edm::EDGetTokenT< std::vector
< pat::MET > > 
metsToken_
 
hitfit::METTranslatorBase
< pat::MET
metTranslator_
 
double minBTagValueBJet_
 min value of bTag for a b-jet More...
 
double mTop_
 
hitfit::LeptonTranslatorBase
< pat::Muon
muonTranslator_
 
double mW_
 constraints More...
 
bool useBTag_
 switch to tell whether to use b-tagging or not More...
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

template<typename LeptonCollection>
class TtSemiLepHitFitProducer< LeptonCollection >

Definition at line 20 of file TtSemiLepHitFitProducer.h.

Member Typedef Documentation

template<typename LeptonCollection >
typedef hitfit::RunHitFit<pat::Electron,pat::Muon,pat::Jet,pat::MET> TtSemiLepHitFitProducer< LeptonCollection >::PatHitFit
private

Definition at line 83 of file TtSemiLepHitFitProducer.h.

Constructor & Destructor Documentation

template<typename LeptonCollection >
TtSemiLepHitFitProducer< LeptonCollection >::TtSemiLepHitFitProducer ( const edm::ParameterSet cfg)
explicit

Definition at line 101 of file TtSemiLepHitFitProducer.h.

References TtSemiLepHitFitProducer< LeptonCollection >::electronTranslator_, edm::FileInPath::fullPath(), TtSemiLepHitFitProducer< LeptonCollection >::HitFit, TtSemiLepHitFitProducer< LeptonCollection >::hitfitDefault_, TtSemiLepHitFitProducer< LeptonCollection >::jetTranslator_, TtSemiLepHitFitProducer< LeptonCollection >::maxEtaEle_, TtSemiLepHitFitProducer< LeptonCollection >::maxEtaJet_, TtSemiLepHitFitProducer< LeptonCollection >::maxEtaMu_, TtSemiLepHitFitProducer< LeptonCollection >::metTranslator_, TtSemiLepHitFitProducer< LeptonCollection >::mTop_, TtSemiLepHitFitProducer< LeptonCollection >::muonTranslator_, and TtSemiLepHitFitProducer< LeptonCollection >::mW_.

101  :
102  jetsToken_ (consumes<std::vector<pat::Jet> >(cfg.getParameter<edm::InputTag>("jets"))),
103  lepsToken_ (consumes<LeptonCollection>(cfg.getParameter<edm::InputTag>("leps"))),
104  metsToken_ (consumes<std::vector<pat::MET> >(cfg.getParameter<edm::InputTag>("mets"))),
105  maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
106  maxNComb_ (cfg.getParameter<int> ("maxNComb" )),
107  bTagAlgo_ (cfg.getParameter<std::string> ("bTagAlgo" )),
108  minBTagValueBJet_ (cfg.getParameter<double> ("minBDiscBJets" )),
109  maxBTagValueNonBJet_ (cfg.getParameter<double> ("maxBDiscLightJets" )),
110  useBTag_ (cfg.getParameter<bool> ("useBTagging" )),
111  mW_ (cfg.getParameter<double> ("mW" )),
112  mTop_ (cfg.getParameter<double> ("mTop" )),
113  jetCorrectionLevel_ (cfg.getParameter<std::string> ("jetCorrectionLevel" )),
114  jes_ (cfg.getParameter<double> ("jes" )),
115  jesB_ (cfg.getParameter<double> ("jesB" )),
116 
117  // The following five initializers read the config parameters for the
118  // ASCII text files which contains the physics object resolutions.
120  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/setting/RunHitFitConfiguration.txt")))),
122  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafElectronResolution.txt")))),
124  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafMuonResolution.txt")))),
126  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafUdscJetResolution.txt")))),
128  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafBJetResolution.txt")))),
130  edm::FileInPath(std::string("TopQuarkAnalysis/TopHitFit/data/resolution/tqafKtResolution.txt")))),
131 
132  // The following four initializers instantiate the translator between PAT objects
133  // and HitFit objects using the ASCII text files which contains the resolutions.
138 
139 {
140  // Create an instance of RunHitFit and initialize it.
146  mW_,
147  mW_,
148  mTop_);
149 
150  maxEtaMu_ = 2.4;
151  maxEtaEle_ = 2.5;
152  maxEtaJet_ = 2.5;
153 
154  edm::LogVerbatim( "TopHitFit" )
155  << "\n"
156  << "+++++++++++ TtSemiLepHitFitProducer ++++++++++++ \n"
157  << " Due to the eta ranges for which resolutions \n"
158  << " are provided in \n"
159  << " TopQuarkAnalysis/TopHitFit/data/resolution/ \n"
160  << " so far, the following cuts are currently \n"
161  << " implemented in the TtSemiLepHitFitProducer: \n"
162  << " |eta(muons )| <= " << maxEtaMu_ << " \n"
163  << " |eta(electrons)| <= " << maxEtaEle_ << " \n"
164  << " |eta(jets )| <= " << maxEtaJet_ << " \n"
165  << "++++++++++++++++++++++++++++++++++++++++++++++++ \n";
166 
167  produces< std::vector<pat::Particle> >("PartonsHadP");
168  produces< std::vector<pat::Particle> >("PartonsHadQ");
169  produces< std::vector<pat::Particle> >("PartonsHadB");
170  produces< std::vector<pat::Particle> >("PartonsLepB");
171  produces< std::vector<pat::Particle> >("Leptons");
172  produces< std::vector<pat::Particle> >("Neutrinos");
173 
174  produces< std::vector<std::vector<int> > >();
175  produces< std::vector<double> >("Chi2");
176  produces< std::vector<double> >("Prob");
177  produces< std::vector<double> >("MT");
178  produces< std::vector<double> >("SigMT");
179  produces< std::vector<int> >("Status");
180  produces< int >("NumberOfConsideredJets");
181 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int maxNComb_
maximal number of combinations to be written to the event
double maxEtaEle_
maximum eta value for electrons, needed to limited range in which resolutions are provided ...
double maxBTagValueNonBJet_
max value of bTag for a non-b-jet
hitfit::JetTranslatorBase< pat::Jet > jetTranslator_
hitfit::LeptonTranslatorBase< pat::Electron > electronTranslator_
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
hitfit::LeptonTranslatorBase< pat::Muon > muonTranslator_
double maxEtaJet_
maximum eta value for jets, needed to limited range in which resolutions are provided ...
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::FileInPath hitfitUdscJetResolution_
bool useBTag_
switch to tell whether to use b-tagging or not
hitfit::RunHitFit< pat::Electron, pat::Muon, pat::Jet, pat::MET > PatHitFit
int maxNJets_
maximal number of jets (-1 possible to indicate &#39;all&#39;)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
hitfit::METTranslatorBase< pat::MET > metTranslator_
edm::EDGetTokenT< LeptonCollection > lepsToken_
edm::FileInPath hitfitElectronResolution_
std::string jetCorrectionLevel_
jet correction level
double maxEtaMu_
maximum eta value for muons, needed to limited range in which resolutions are provided ...
std::string fullPath() const
Definition: FileInPath.cc:165
double minBTagValueBJet_
min value of bTag for a b-jet
std::string bTagAlgo_
input tag for b-tagging algorithm
double jes_
jet energy scale
template<typename LeptonCollection >
TtSemiLepHitFitProducer< LeptonCollection >::~TtSemiLepHitFitProducer ( )

Definition at line 184 of file TtSemiLepHitFitProducer.h.

185 {
186  delete HitFit;
187 }

Member Function Documentation

template<typename LeptonCollection >
void TtSemiLepHitFitProducer< LeptonCollection >::produce ( edm::Event evt,
const edm::EventSetup setup 
)
privatevirtual

Implements edm::EDProducer.

Definition at line 190 of file TtSemiLepHitFitProducer.h.

References funct::abs(), TtSemiLepHitFitProducer< LeptonCollection >::FitResult::Chi2, eta, create_public_lumi_plots::exp, edm::Event::getByToken(), TtSemiLepEvtPartons::HadB, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::HadB, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::HadP, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::HadQ, i, metsig::jet, hitfit::Lepjets_Event::jet(), TtSemiLepHitFitProducer< LeptonCollection >::FitResult::JetCombi, fwrapper::jets, hitfit::Lepjets_Event::lep(), TtSemiLepEvtPartons::LepB, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::LepB, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::LepL, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::LepN, TtSemiLepEvtPartons::LightQ, TtSemiLepEvtPartons::LightQBar, hitfit::Lepjets_Event::met(), TtSemiLepHitFitProducer< LeptonCollection >::FitResult::MT, nPartons, hitfit::Lepjets_Event_Lep::p(), TtSemiLepHitFitProducer< LeptonCollection >::FitResult::Prob, edm::Event::put(), query::result, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::SigMT, TtSemiLepHitFitProducer< LeptonCollection >::FitResult::Status, and hitfit::Lepjets_Event_Lep::type().

191 {
192  std::auto_ptr< std::vector<pat::Particle> > pPartonsHadP( new std::vector<pat::Particle> );
193  std::auto_ptr< std::vector<pat::Particle> > pPartonsHadQ( new std::vector<pat::Particle> );
194  std::auto_ptr< std::vector<pat::Particle> > pPartonsHadB( new std::vector<pat::Particle> );
195  std::auto_ptr< std::vector<pat::Particle> > pPartonsLepB( new std::vector<pat::Particle> );
196  std::auto_ptr< std::vector<pat::Particle> > pLeptons ( new std::vector<pat::Particle> );
197  std::auto_ptr< std::vector<pat::Particle> > pNeutrinos ( new std::vector<pat::Particle> );
198 
199  std::auto_ptr< std::vector<std::vector<int> > > pCombi ( new std::vector<std::vector<int> > );
200  std::auto_ptr< std::vector<double> > pChi2 ( new std::vector<double> );
201  std::auto_ptr< std::vector<double> > pProb ( new std::vector<double> );
202  std::auto_ptr< std::vector<double> > pMT ( new std::vector<double> );
203  std::auto_ptr< std::vector<double> > pSigMT ( new std::vector<double> );
204  std::auto_ptr< std::vector<int> > pStatus( new std::vector<int> );
205  std::auto_ptr< int > pJetsConsidered( new int );
206 
208  evt.getByToken(jetsToken_, jets);
209 
211  evt.getByToken(metsToken_, mets);
212 
214  evt.getByToken(lepsToken_, leps);
215 
216  // -----------------------------------------------------
217  // skip events with no appropriate lepton candidate in
218  // or empty MET or less jets than partons
219  // -----------------------------------------------------
220 
221  const unsigned int nPartons = 4;
222 
223  // Clear the internal state
224  HitFit->clear();
225 
226  // Add lepton into HitFit
227  bool foundLepton = false;
228  if(!leps->empty()) {
229  double maxEtaLep = maxEtaMu_;
230  if( !dynamic_cast<const reco::Muon*>(&((*leps)[0])) ) // assume electron if it is not a muon
231  maxEtaLep = maxEtaEle_;
232  for(unsigned iLep=0; iLep<(*leps).size() && !foundLepton; ++iLep) {
233  if(std::abs((*leps)[iLep].eta()) <= maxEtaLep) {
234  HitFit->AddLepton((*leps)[iLep]);
235  foundLepton = true;
236  }
237  }
238  }
239 
240  // Add jets into HitFit
241  unsigned int nJetsFound = 0;
242  for(unsigned iJet=0; iJet<(*jets).size() && (int)nJetsFound!=maxNJets_; ++iJet) {
243  double jet_eta = (*jets)[iJet].eta();
244  if ( (*jets)[iJet].isCaloJet() ) {
245  jet_eta = ((reco::CaloJet*) ((*jets)[iJet]).originalObject())->detectorP4().eta();
246  }
247  if(std::abs(jet_eta) <= maxEtaJet_) {
248  HitFit->AddJet((*jets)[iJet]);
249  nJetsFound++;
250  }
251  }
252  *pJetsConsidered = nJetsFound;
253 
254  // Add missing transverse energy into HitFit
255  if(!mets->empty())
256  HitFit->SetMet((*mets)[0]);
257 
258  if( !foundLepton || mets->empty() || nJetsFound<nPartons ) {
259  // the kinFit getters return empty objects here
260  pPartonsHadP->push_back( pat::Particle() );
261  pPartonsHadQ->push_back( pat::Particle() );
262  pPartonsHadB->push_back( pat::Particle() );
263  pPartonsLepB->push_back( pat::Particle() );
264  pLeptons ->push_back( pat::Particle() );
265  pNeutrinos ->push_back( pat::Particle() );
266  // indices referring to the jet combination
267  std::vector<int> invalidCombi;
268  for(unsigned int i = 0; i < nPartons; ++i)
269  invalidCombi.push_back( -1 );
270  pCombi->push_back( invalidCombi );
271  // chi2
272  pChi2->push_back( -1. );
273  // chi2 probability
274  pProb->push_back( -1. );
275  // fitted top mass
276  pMT->push_back( -1. );
277  pSigMT->push_back( -1. );
278  // status of the fitter
279  pStatus->push_back( -1 );
280  // feed out all products
281  evt.put(pCombi);
282  evt.put(pPartonsHadP, "PartonsHadP");
283  evt.put(pPartonsHadQ, "PartonsHadQ");
284  evt.put(pPartonsHadB, "PartonsHadB");
285  evt.put(pPartonsLepB, "PartonsLepB");
286  evt.put(pLeptons , "Leptons" );
287  evt.put(pNeutrinos , "Neutrinos" );
288  evt.put(pChi2 , "Chi2" );
289  evt.put(pProb , "Prob" );
290  evt.put(pMT , "MT" );
291  evt.put(pSigMT , "SigMT" );
292  evt.put(pStatus , "Status" );
293  evt.put(pJetsConsidered, "NumberOfConsideredJets");
294  return;
295  }
296 
297  std::list<FitResult> FitResultList;
298 
299  //
300  // BEGIN DECLARATION OF VARIABLES FROM KINEMATIC FIT
301  //
302 
303  // In this part are variables from the
304  // kinematic fit procedure
305 
306  // Number of all permutations of the event
307  size_t nHitFit = 0 ;
308 
309  // Number of jets in the event
310  size_t nHitFitJet = 0 ;
311 
312  // Results of the fit for all jet permutations of the event
313  std::vector<hitfit::Fit_Result> hitFitResult;
314 
315  //
316  // R U N H I T F I T
317  //
318  // Run the kinematic fit and get how many permutations are possible
319  // in the fit
320 
321  nHitFit = HitFit->FitAllPermutation();
322 
323  //
324  // BEGIN PART WHICH EXTRACTS INFORMATION FROM HITFIT
325  //
326 
327  // Get the number of jets
328  nHitFitJet = HitFit->GetUnfittedEvent()[0].njets();
329 
330  // Get the fit results for all permutations
331  hitFitResult = HitFit->GetFitAllPermutation();
332 
333  // Loop over all permutations and extract the information
334  for (size_t fit = 0 ; fit != nHitFit ; ++fit) {
335 
336  // Get the event after the fit
337  hitfit::Lepjets_Event fittedEvent = hitFitResult[fit].ev();
338 
339  /*
340  Get jet permutation according to TQAF convention
341  11 : leptonic b
342  12 : hadronic b
343  13 : hadronic W
344  14 : hadronic W
345  */
346  std::vector<int> hitCombi(4);
347  for (size_t jet = 0 ; jet != nHitFitJet ; ++jet) {
348  int jet_type = fittedEvent.jet(jet).type();
349 
350  switch(jet_type) {
351  case 11: hitCombi[TtSemiLepEvtPartons::LepB ] = jet;
352  break;
353  case 12: hitCombi[TtSemiLepEvtPartons::HadB ] = jet;
354  break;
355  case 13: hitCombi[TtSemiLepEvtPartons::LightQ ] = jet;
356  break;
357  case 14: hitCombi[TtSemiLepEvtPartons::LightQBar] = jet;
358  break;
359  }
360  }
361 
362  // Store the kinematic quantities in the corresponding containers.
363 
364  hitfit::Lepjets_Event_Jet hadP_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LightQ ]);
365  hitfit::Lepjets_Event_Jet hadQ_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LightQBar]);
366  hitfit::Lepjets_Event_Jet hadB_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::HadB ]);
367  hitfit::Lepjets_Event_Jet lepB_ = fittedEvent.jet(hitCombi[TtSemiLepEvtPartons::LepB ]);
368  hitfit::Lepjets_Event_Lep lepL_ = fittedEvent.lep(0);
369 
370  /*
372  std::string bTagAlgo_;
374  double minBTagValueBJet_;
376  double maxBTagValueNonBJet_;
378  bool useBTag_;
379  */
380 
381  if ( hitFitResult[fit].chisq() > 0 // only take into account converged fits
382  && (!useBTag_ || ( useBTag_ // use btag information if chosen
383  && jets->at(hitCombi[TtSemiLepEvtPartons::LightQ ]).bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_
384  && jets->at(hitCombi[TtSemiLepEvtPartons::LightQBar]).bDiscriminator(bTagAlgo_) < maxBTagValueNonBJet_
385  && jets->at(hitCombi[TtSemiLepEvtPartons::HadB ]).bDiscriminator(bTagAlgo_) > minBTagValueBJet_
386  && jets->at(hitCombi[TtSemiLepEvtPartons::LepB ]).bDiscriminator(bTagAlgo_) > minBTagValueBJet_
387  )
388  )
389  ) {
390  FitResult result;
391  result.Status = 0;
392  result.Chi2 = hitFitResult[fit].chisq();
393  result.Prob = exp(-1.0*(hitFitResult[fit].chisq())/2.0);
394  result.MT = hitFitResult[fit].mt();
395  result.SigMT= hitFitResult[fit].sigmt();
396  result.HadB = pat::Particle(reco::LeafCandidate(0,
397  math::XYZTLorentzVector(hadB_.p().x(), hadB_.p().y(),
398  hadB_.p().z(), hadB_.p().t()), math::XYZPoint()));
399  result.HadP = pat::Particle(reco::LeafCandidate(0,
400  math::XYZTLorentzVector(hadP_.p().x(), hadP_.p().y(),
401  hadP_.p().z(), hadP_.p().t()), math::XYZPoint()));
402  result.HadQ = pat::Particle(reco::LeafCandidate(0,
403  math::XYZTLorentzVector(hadQ_.p().x(), hadQ_.p().y(),
404  hadQ_.p().z(), hadQ_.p().t()), math::XYZPoint()));
405  result.LepB = pat::Particle(reco::LeafCandidate(0,
406  math::XYZTLorentzVector(lepB_.p().x(), lepB_.p().y(),
407  lepB_.p().z(), lepB_.p().t()), math::XYZPoint()));
408  result.LepL = pat::Particle(reco::LeafCandidate(0,
409  math::XYZTLorentzVector(lepL_.p().x(), lepL_.p().y(),
410  lepL_.p().z(), lepL_.p().t()), math::XYZPoint()));
411  result.LepN = pat::Particle(reco::LeafCandidate(0,
412  math::XYZTLorentzVector(fittedEvent.met().x(), fittedEvent.met().y(),
413  fittedEvent.met().z(), fittedEvent.met().t()), math::XYZPoint()));
414  result.JetCombi = hitCombi;
415 
416  FitResultList.push_back(result);
417  }
418 
419  }
420 
421  // sort results w.r.t. chi2 values
422  FitResultList.sort();
423 
424  // -----------------------------------------------------
425  // feed out result
426  // starting with the JetComb having the smallest chi2
427  // -----------------------------------------------------
428 
429  if( ((unsigned)FitResultList.size())<1 ) { // in case no fit results were stored in the list (all fits aborted)
430  pPartonsHadP->push_back( pat::Particle() );
431  pPartonsHadQ->push_back( pat::Particle() );
432  pPartonsHadB->push_back( pat::Particle() );
433  pPartonsLepB->push_back( pat::Particle() );
434  pLeptons ->push_back( pat::Particle() );
435  pNeutrinos ->push_back( pat::Particle() );
436  // indices referring to the jet combination
437  std::vector<int> invalidCombi;
438  for(unsigned int i = 0; i < nPartons; ++i)
439  invalidCombi.push_back( -1 );
440  pCombi->push_back( invalidCombi );
441  // chi2
442  pChi2->push_back( -1. );
443  // chi2 probability
444  pProb->push_back( -1. );
445  // fitted top mass
446  pMT->push_back( -1. );
447  pSigMT->push_back( -1. );
448  // status of the fitter
449  pStatus->push_back( -1 );
450  }
451  else {
452  unsigned int iComb = 0;
453  for(typename std::list<FitResult>::const_iterator result = FitResultList.begin(); result != FitResultList.end(); ++result) {
454  if(maxNComb_ >= 1 && iComb == (unsigned int) maxNComb_) break;
455  iComb++;
456  // partons
457  pPartonsHadP->push_back( result->HadP );
458  pPartonsHadQ->push_back( result->HadQ );
459  pPartonsHadB->push_back( result->HadB );
460  pPartonsLepB->push_back( result->LepB );
461  // lepton
462  pLeptons->push_back( result->LepL );
463  // neutrino
464  pNeutrinos->push_back( result->LepN );
465  // indices referring to the jet combination
466  pCombi->push_back( result->JetCombi );
467  // chi2
468  pChi2->push_back( result->Chi2 );
469  // chi2 probability
470  pProb->push_back( result->Prob );
471  // fitted top mass
472  pMT->push_back( result->MT );
473  pSigMT->push_back( result->SigMT );
474  // status of the fitter
475  pStatus->push_back( result->Status );
476  }
477  }
478  evt.put(pCombi);
479  evt.put(pPartonsHadP, "PartonsHadP");
480  evt.put(pPartonsHadQ, "PartonsHadQ");
481  evt.put(pPartonsHadB, "PartonsHadB");
482  evt.put(pPartonsLepB, "PartonsLepB");
483  evt.put(pLeptons , "Leptons" );
484  evt.put(pNeutrinos , "Neutrinos" );
485  evt.put(pChi2 , "Chi2" );
486  evt.put(pProb , "Prob" );
487  evt.put(pMT , "MT" );
488  evt.put(pSigMT , "SigMT" );
489  evt.put(pStatus , "Status" );
490  evt.put(pJetsConsidered, "NumberOfConsideredJets");
491 }
int maxNComb_
maximal number of combinations to be written to the event
double maxEtaEle_
maximum eta value for electrons, needed to limited range in which resolutions are provided ...
int i
Definition: DBlmapReader.cc:9
std::vector< Lepjets_Event > GetUnfittedEvent()
Return the unfitted events for all permutations.
Definition: RunHitFit.h:521
double maxBTagValueNonBJet_
max value of bTag for a non-b-jet
static const unsigned int nPartons
Jets made from CaloTowers.
Definition: CaloJet.h:29
std::vector< Fit_Result >::size_type FitAllPermutation()
Fit all permutations of the internal event. Returns the number of permutations.
Definition: RunHitFit.h:418
Represent a lepton in an instance of Lepjets_Event class. This class hold the following information: ...
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
void SetMet(const AMet &met, bool useObjRes=false)
Set the missing transverse energy of the internal event.
Definition: RunHitFit.h:371
Lepjets_Event_Jet & jet(std::vector< Lepjets_Event_Jet >::size_type i)
Return a reference to jet at index position i.
double maxEtaJet_
maximum eta value for jets, needed to limited range in which resolutions are provided ...
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int & type()
Return a reference to the type code.
std::vector< Fit_Result > GetFitAllPermutation()
Return the results of fitting all permutations of the internal event.
Definition: RunHitFit.h:531
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
bool useBTag_
switch to tell whether to use b-tagging or not
vector< PseudoJet > jets
Represent a simple event consisting of lepton(s) and jet(s). An instance of this class holds a list o...
Definition: Lepjets_Event.h:66
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
A class to represent a jet in an instance of Lepjets_Event class. The class is derived from the Lepje...
int maxNJets_
maximal number of jets (-1 possible to indicate &#39;all&#39;)
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
Fourvec & met()
Return a reference to the missing transverse energy.
Lepjets_Event_Lep & lep(std::vector< Lepjets_Event_Lep >::size_type i)
Return a reference to lepton at index position i.
edm::EDGetTokenT< LeptonCollection > lepsToken_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
Analysis-level particle class.
Definition: Particle.h:32
Fourvec & p()
Return a reference to the four-momentum.
double maxEtaMu_
maximum eta value for muons, needed to limited range in which resolutions are provided ...
double minBTagValueBJet_
min value of bTag for a b-jet
void AddLepton(const AElectron &electron, bool useObjRes=false)
Add one electron into the internal event.
Definition: RunHitFit.h:306
void clear()
Clear the internal event, fit results, and jets.
Definition: RunHitFit.h:287
std::string bTagAlgo_
input tag for b-tagging algorithm
void AddJet(const AJet &jet, bool useObjRes=false)
Add one jet into the internal event. This function will do nothing if the internal event has already ...
Definition: RunHitFit.h:351

Member Data Documentation

template<typename LeptonCollection >
std::string TtSemiLepHitFitProducer< LeptonCollection >::bTagAlgo_
private

input tag for b-tagging algorithm

Definition at line 48 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
hitfit::LeptonTranslatorBase<pat::Electron> TtSemiLepHitFitProducer< LeptonCollection >::electronTranslator_
private
template<typename LeptonCollection >
PatHitFit* TtSemiLepHitFitProducer< LeptonCollection >::HitFit
private
template<typename LeptonCollection >
edm::FileInPath TtSemiLepHitFitProducer< LeptonCollection >::hitfitBJetResolution_
private

Definition at line 89 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
edm::FileInPath TtSemiLepHitFitProducer< LeptonCollection >::hitfitDefault_
private
template<typename LeptonCollection >
edm::FileInPath TtSemiLepHitFitProducer< LeptonCollection >::hitfitElectronResolution_
private

Definition at line 86 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
edm::FileInPath TtSemiLepHitFitProducer< LeptonCollection >::hitfitMETResolution_
private

Definition at line 90 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
edm::FileInPath TtSemiLepHitFitProducer< LeptonCollection >::hitfitMuonResolution_
private

Definition at line 87 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
edm::FileInPath TtSemiLepHitFitProducer< LeptonCollection >::hitfitUdscJetResolution_
private

Definition at line 88 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::jes_
private

jet energy scale

Definition at line 64 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::jesB_
private

Definition at line 65 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
std::string TtSemiLepHitFitProducer< LeptonCollection >::jetCorrectionLevel_
private

jet correction level

Definition at line 61 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
edm::EDGetTokenT<std::vector<pat::Jet> > TtSemiLepHitFitProducer< LeptonCollection >::jetsToken_
private

Definition at line 31 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
hitfit::JetTranslatorBase<pat::Jet> TtSemiLepHitFitProducer< LeptonCollection >::jetTranslator_
private
template<typename LeptonCollection >
edm::EDGetTokenT<LeptonCollection> TtSemiLepHitFitProducer< LeptonCollection >::lepsToken_
private

Definition at line 32 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::maxBTagValueNonBJet_
private

max value of bTag for a non-b-jet

Definition at line 52 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::maxEtaEle_
private

maximum eta value for electrons, needed to limited range in which resolutions are provided

Definition at line 43 of file TtSemiLepHitFitProducer.h.

Referenced by TtSemiLepHitFitProducer< LeptonCollection >::TtSemiLepHitFitProducer().

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::maxEtaJet_
private

maximum eta value for jets, needed to limited range in which resolutions are provided

Definition at line 45 of file TtSemiLepHitFitProducer.h.

Referenced by TtSemiLepHitFitProducer< LeptonCollection >::TtSemiLepHitFitProducer().

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::maxEtaMu_
private

maximum eta value for muons, needed to limited range in which resolutions are provided

Definition at line 41 of file TtSemiLepHitFitProducer.h.

Referenced by TtSemiLepHitFitProducer< LeptonCollection >::TtSemiLepHitFitProducer().

template<typename LeptonCollection >
int TtSemiLepHitFitProducer< LeptonCollection >::maxNComb_
private

maximal number of combinations to be written to the event

Definition at line 38 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
int TtSemiLepHitFitProducer< LeptonCollection >::maxNJets_
private

maximal number of jets (-1 possible to indicate 'all')

Definition at line 36 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
edm::EDGetTokenT<std::vector<pat::MET> > TtSemiLepHitFitProducer< LeptonCollection >::metsToken_
private

Definition at line 33 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
hitfit::METTranslatorBase<pat::MET> TtSemiLepHitFitProducer< LeptonCollection >::metTranslator_
private
template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::minBTagValueBJet_
private

min value of bTag for a b-jet

Definition at line 50 of file TtSemiLepHitFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::mTop_
private
template<typename LeptonCollection >
hitfit::LeptonTranslatorBase<pat::Muon> TtSemiLepHitFitProducer< LeptonCollection >::muonTranslator_
private
template<typename LeptonCollection >
double TtSemiLepHitFitProducer< LeptonCollection >::mW_
private
template<typename LeptonCollection >
bool TtSemiLepHitFitProducer< LeptonCollection >::useBTag_
private

switch to tell whether to use b-tagging or not

Definition at line 54 of file TtSemiLepHitFitProducer.h.