CMS 3D CMS Logo

Classes | Public Member Functions | Private Member Functions | Private Attributes

TtSemiLepKinFitProducer< LeptonCollection > Class Template Reference

#include <TtSemiLepKinFitProducer.h>

Inheritance diagram for TtSemiLepKinFitProducer< LeptonCollection >:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Classes

struct  KinFitResult

Public Member Functions

 TtSemiLepKinFitProducer (const edm::ParameterSet &)
 ~TtSemiLepKinFitProducer ()

Private Member Functions

TtSemiLepKinFitter::Constraint constraint (unsigned)
std::vector
< TtSemiLepKinFitter::Constraint
constraints (std::vector< unsigned > &)
bool doBTagging (bool &useBTag_, edm::Handle< std::vector< pat::Jet > > &jets, std::vector< int > &combi, std::string &bTagAlgo_, double &minBTagValueBJets_, double &maxBTagValueNonBJets_)
TtSemiLepKinFitter::Param param (unsigned)
virtual void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

std::vector< edm::ParameterSetbResolutions_
std::string bTagAlgo_
 input tag for b-tagging algorithm
std::vector< unsigned > constraints_
 constrains
TtSemiLepKinFitterfitter
std::vector< double > jetEnergyResolutionEtaBinning_
std::vector< double > jetEnergyResolutionScaleFactors_
 scale factors for jet energy resolution
unsigned int jetParam_
edm::InputTag jets_
unsigned int lepParam_
std::vector< edm::ParameterSetlepResolutions_
edm::InputTag leps_
edm::InputTag match_
double maxBTagValueNonBJet_
 max value of bTag for a non-b-jet
double maxDeltaS_
 maximal chi2 equivalent
double maxF_
 maximal deviation for contstraints
int maxNComb_
 maximal number of combinations to be written to the event
int maxNJets_
 maximal number of jets (-1 possible to indicate 'all')
unsigned int maxNrIter_
 maximal number of iterations to be performed for the fit
unsigned int metParam_
std::vector< edm::ParameterSetmetResolutions_
edm::InputTag mets_
double minBTagValueBJet_
 min value of bTag for a b-jet
double mTop_
double mW_
std::vector< edm::ParameterSetudscResolutions_
 config-file-based object resolutions
bool useBTag_
 switch to tell whether to use b-tagging or not
bool useOnlyMatch_
 switch to use only a combination given by another hypothesis

Detailed Description

template<typename LeptonCollection>
class TtSemiLepKinFitProducer< LeptonCollection >

Definition at line 13 of file TtSemiLepKinFitProducer.h.


Constructor & Destructor Documentation

template<typename LeptonCollection >
TtSemiLepKinFitProducer< LeptonCollection >::TtSemiLepKinFitProducer ( const edm::ParameterSet cfg) [explicit]

Definition at line 94 of file TtSemiLepKinFitProducer.h.

References TtSemiLepKinFitProducer< LeptonCollection >::bResolutions_, TtSemiLepKinFitProducer< LeptonCollection >::constraints(), TtSemiLepKinFitProducer< LeptonCollection >::constraints_, Exception, edm::ParameterSet::exists(), TtSemiLepKinFitProducer< LeptonCollection >::fitter, edm::ParameterSet::getParameter(), TtSemiLepKinFitProducer< LeptonCollection >::jetEnergyResolutionEtaBinning_, TtSemiLepKinFitProducer< LeptonCollection >::jetEnergyResolutionScaleFactors_, TtSemiLepKinFitProducer< LeptonCollection >::jetParam_, TtSemiLepKinFitProducer< LeptonCollection >::lepParam_, TtSemiLepKinFitProducer< LeptonCollection >::lepResolutions_, TtSemiLepKinFitProducer< LeptonCollection >::maxDeltaS_, TtSemiLepKinFitProducer< LeptonCollection >::maxF_, TtSemiLepKinFitProducer< LeptonCollection >::maxNrIter_, TtSemiLepKinFitProducer< LeptonCollection >::metParam_, TtSemiLepKinFitProducer< LeptonCollection >::metResolutions_, TtSemiLepKinFitProducer< LeptonCollection >::mTop_, TtSemiLepKinFitProducer< LeptonCollection >::mW_, TtSemiLepKinFitProducer< LeptonCollection >::param(), and TtSemiLepKinFitProducer< LeptonCollection >::udscResolutions_.

                                                                                            :
  jets_                    (cfg.getParameter<edm::InputTag>("jets")),
  leps_                    (cfg.getParameter<edm::InputTag>("leps")),
  mets_                    (cfg.getParameter<edm::InputTag>("mets")),
  match_                   (cfg.getParameter<edm::InputTag>("match")),
  useOnlyMatch_            (cfg.getParameter<bool>         ("useOnlyMatch"        )),
  bTagAlgo_                (cfg.getParameter<std::string>  ("bTagAlgo"            )),
  minBTagValueBJet_        (cfg.getParameter<double>       ("minBDiscBJets"       )),
  maxBTagValueNonBJet_     (cfg.getParameter<double>       ("maxBDiscLightJets"   )),
  useBTag_                 (cfg.getParameter<bool>         ("useBTagging"         )),
  maxNJets_                (cfg.getParameter<int>          ("maxNJets"            )),
  maxNComb_                (cfg.getParameter<int>          ("maxNComb"            )),
  maxNrIter_               (cfg.getParameter<unsigned>     ("maxNrIter"           )),
  maxDeltaS_               (cfg.getParameter<double>       ("maxDeltaS"           )),
  maxF_                    (cfg.getParameter<double>       ("maxF"                )),
  jetParam_                (cfg.getParameter<unsigned>     ("jetParametrisation"  )),
  lepParam_                (cfg.getParameter<unsigned>     ("lepParametrisation"  )),
  metParam_                (cfg.getParameter<unsigned>     ("metParametrisation"  )),
  constraints_             (cfg.getParameter<std::vector<unsigned> >("constraints")),
  mW_                      (cfg.getParameter<double>       ("mW"                  )),
  mTop_                    (cfg.getParameter<double>       ("mTop"                )),
  jetEnergyResolutionScaleFactors_(cfg.getParameter<std::vector<double> >("jetEnergyResolutionScaleFactors")),
  jetEnergyResolutionEtaBinning_  (cfg.getParameter<std::vector<double> >("jetEnergyResolutionEtaBinning")),
  udscResolutions_(0), bResolutions_(0), lepResolutions_(0), metResolutions_(0)
{
  if(cfg.exists("udscResolutions") && cfg.exists("bResolutions") && cfg.exists("lepResolutions") && cfg.exists("metResolutions")){
    udscResolutions_ = cfg.getParameter<std::vector<edm::ParameterSet> >("udscResolutions");
    bResolutions_    = cfg.getParameter<std::vector<edm::ParameterSet> >("bResolutions"   );
    lepResolutions_  = cfg.getParameter<std::vector<edm::ParameterSet> >("lepResolutions" );
    metResolutions_  = cfg.getParameter<std::vector<edm::ParameterSet> >("metResolutions" );
  }
  else if(cfg.exists("udscResolutions") || cfg.exists("bResolutions") || cfg.exists("lepResolutions") || cfg.exists("metResolutions") ){
    throw cms::Exception("Configuration") << "Parameters 'udscResolutions', 'bResolutions', 'lepResolutions', 'metResolutions' should be used together.\n";
  }

  fitter = new TtSemiLepKinFitter(param(jetParam_), param(lepParam_), param(metParam_), maxNrIter_, maxDeltaS_, maxF_,
                                  constraints(constraints_), mW_, mTop_, &udscResolutions_, &bResolutions_, &lepResolutions_, &metResolutions_,
                                  &jetEnergyResolutionScaleFactors_, &jetEnergyResolutionEtaBinning_);

  produces< std::vector<pat::Particle> >("PartonsHadP");
  produces< std::vector<pat::Particle> >("PartonsHadQ");
  produces< std::vector<pat::Particle> >("PartonsHadB");
  produces< std::vector<pat::Particle> >("PartonsLepB");
  produces< std::vector<pat::Particle> >("Leptons");
  produces< std::vector<pat::Particle> >("Neutrinos");

  produces< std::vector<std::vector<int> > >();
  produces< std::vector<double> >("Chi2");
  produces< std::vector<double> >("Prob");
  produces< std::vector<int> >("Status");

  produces<int>("NumberOfConsideredJets");
}
template<typename LeptonCollection >
TtSemiLepKinFitProducer< LeptonCollection >::~TtSemiLepKinFitProducer ( )

Definition at line 149 of file TtSemiLepKinFitProducer.h.

{
  delete fitter;
}

Member Function Documentation

template<typename LeptonCollection >
TtFullHadKinFitter::Constraint TtFullHadKinFitter::KinFit::constraint ( unsigned  val) [private]
template<typename LeptonCollection >
std::vector< TtFullHadKinFitter::Constraint > TtFullHadKinFitter::KinFit::constraints ( std::vector< unsigned > &  val) [private]

Definition at line 427 of file TtSemiLepKinFitProducer.h.

References i, and query::result.

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

{
  std::vector<TtSemiLepKinFitter::Constraint> result;
  for(unsigned i=0; i<val.size(); ++i){
    result.push_back(constraint(val[i]));
  }
  return result; 
}
template<typename LeptonCollection >
bool TtSemiLepKinFitProducer< LeptonCollection >::doBTagging ( bool &  useBTag_,
edm::Handle< std::vector< pat::Jet > > &  jets,
std::vector< int > &  combi,
std::string &  bTagAlgo_,
double &  minBTagValueBJets_,
double &  maxBTagValueNonBJets_ 
) [private]

Definition at line 155 of file TtSemiLepKinFitProducer.h.

References TtSemiLepEvtPartons::HadB, analyzePatCleaning_cfg::jets, TtSemiLepEvtPartons::LepB, TtSemiLepEvtPartons::LightQ, and TtSemiLepEvtPartons::LightQBar.

                                                                                                                                         {
  
  if( !useBTag_ ) {
    return true;
  }
  if( useBTag_ &&
      (*jets)[combi[TtSemiLepEvtPartons::HadB     ]].bDiscriminator(bTagAlgo_) >= minBTagValueBJet_ &&
      (*jets)[combi[TtSemiLepEvtPartons::LepB     ]].bDiscriminator(bTagAlgo_) >= minBTagValueBJet_ &&
      (*jets)[combi[TtSemiLepEvtPartons::LightQ   ]].bDiscriminator(bTagAlgo_) <  maxBTagValueNonBJet_ &&
      (*jets)[combi[TtSemiLepEvtPartons::LightQBar]].bDiscriminator(bTagAlgo_) <  maxBTagValueNonBJet_ ) {
    return true;
  }
  else{
    return false;
  }
}
template<typename LeptonCollection >
TtSemiLepKinFitter::Param TtSemiLepKinFitProducer< LeptonCollection >::param ( unsigned  val) [private]

Definition at line 391 of file TtSemiLepKinFitProducer.h.

References Exception, TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, and query::result.

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

{
  TtSemiLepKinFitter::Param result;
  switch(val){
  case TtSemiLepKinFitter::kEMom       : result=TtSemiLepKinFitter::kEMom;       break;
  case TtSemiLepKinFitter::kEtEtaPhi   : result=TtSemiLepKinFitter::kEtEtaPhi;   break;
  case TtSemiLepKinFitter::kEtThetaPhi : result=TtSemiLepKinFitter::kEtThetaPhi; break;
  default: 
    throw cms::Exception("Configuration") 
      << "Chosen jet parametrization is not supported: " << val << "\n";
    break;
  }
  return result;
} 
template<typename LeptonCollection >
void TtSemiLepKinFitProducer< LeptonCollection >::produce ( edm::Event evt,
const edm::EventSetup setup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 174 of file TtSemiLepKinFitProducer.h.

References TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::Chi2, edm::Event::getByLabel(), TtSemiLepEvtPartons::HadB, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::HadB, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::HadP, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::HadQ, i, UserOptions_cff::idx, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::JetCombi, analyzePatCleaning_cfg::jets, TtSemiLepEvtPartons::LepB, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::LepB, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::LepL, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::LepN, TtSemiLepEvtPartons::LightQ, TtSemiLepEvtPartons::LightQBar, match(), stdcomb::next_combination(), nPartons, TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::Prob, edm::Event::put(), query::result, ntuplemaker::status, and TtSemiLepKinFitProducer< LeptonCollection >::KinFitResult::Status.

{
  std::auto_ptr< std::vector<pat::Particle> > pPartonsHadP( new std::vector<pat::Particle> );
  std::auto_ptr< std::vector<pat::Particle> > pPartonsHadQ( new std::vector<pat::Particle> );
  std::auto_ptr< std::vector<pat::Particle> > pPartonsHadB( new std::vector<pat::Particle> );
  std::auto_ptr< std::vector<pat::Particle> > pPartonsLepB( new std::vector<pat::Particle> );
  std::auto_ptr< std::vector<pat::Particle> > pLeptons    ( new std::vector<pat::Particle> );
  std::auto_ptr< std::vector<pat::Particle> > pNeutrinos  ( new std::vector<pat::Particle> );

  std::auto_ptr< std::vector<std::vector<int> > > pCombi ( new std::vector<std::vector<int> > );
  std::auto_ptr< std::vector<double>            > pChi2  ( new std::vector<double> );
  std::auto_ptr< std::vector<double>            > pProb  ( new std::vector<double> );
  std::auto_ptr< std::vector<int>               > pStatus( new std::vector<int> );

  std::auto_ptr<int> pJetsConsidered(new int);

  edm::Handle<std::vector<pat::Jet> > jets;
  evt.getByLabel(jets_, jets);

  edm::Handle<std::vector<pat::MET> > mets;
  evt.getByLabel(mets_, mets);

  edm::Handle<LeptonCollection> leps;
  evt.getByLabel(leps_, leps);

  const unsigned int nPartons = 4;

  std::vector<int> match;
  bool invalidMatch = false;
  if(useOnlyMatch_) {
    *pJetsConsidered = nPartons;
    edm::Handle<std::vector<std::vector<int> > > matchHandle;
    evt.getByLabel(match_, matchHandle);
    match = *(matchHandle->begin());
    // check if match is valid
    if(match.size()!=nPartons) invalidMatch=true;
    else {
      for(unsigned int idx=0; idx<match.size(); ++idx) {
        if(match[idx]<0 || match[idx]>=(int)jets->size()) {
          invalidMatch=true;
          break;
        }
      }
    }
  }

  // -----------------------------------------------------
  // skip events with no appropriate lepton candidate in
  // or empty MET or less jets than partons or invalid match
  // -----------------------------------------------------

  if( leps->empty() || mets->empty() || jets->size()<nPartons || invalidMatch ) {
    // the kinFit getters return empty objects here
    pPartonsHadP->push_back( fitter->fittedHadP()     );
    pPartonsHadQ->push_back( fitter->fittedHadQ()     );
    pPartonsHadB->push_back( fitter->fittedHadB()     );
    pPartonsLepB->push_back( fitter->fittedLepB()     );
    pLeptons    ->push_back( fitter->fittedLepton()   );
    pNeutrinos  ->push_back( fitter->fittedNeutrino() );
    // indices referring to the jet combination
    std::vector<int> invalidCombi;
    for(unsigned int i = 0; i < nPartons; ++i) 
      invalidCombi.push_back( -1 );
    pCombi->push_back( invalidCombi );
    // chi2
    pChi2->push_back( -1. );
    // chi2 probability
    pProb->push_back( -1. );
    // status of the fitter
    pStatus->push_back( -1 );
    // number of jets
    *pJetsConsidered = jets->size();
    // feed out all products
    evt.put(pCombi);
    evt.put(pPartonsHadP, "PartonsHadP");
    evt.put(pPartonsHadQ, "PartonsHadQ");
    evt.put(pPartonsHadB, "PartonsHadB");
    evt.put(pPartonsLepB, "PartonsLepB");
    evt.put(pLeptons    , "Leptons"    );
    evt.put(pNeutrinos  , "Neutrinos"  );
    evt.put(pChi2       , "Chi2"       );
    evt.put(pProb       , "Prob"       );
    evt.put(pStatus     , "Status"     );
    evt.put(pJetsConsidered, "NumberOfConsideredJets");
    return;
  }

  // -----------------------------------------------------
  // analyze different jet combinations using the KinFitter
  // (or only a given jet combination if useOnlyMatch=true)
  // -----------------------------------------------------
  
  std::vector<int> jetIndices;
  if(!useOnlyMatch_) {
    for(unsigned int i=0; i<jets->size(); ++i){
      if(maxNJets_ >= (int) nPartons && maxNJets_ == (int) i) {
        *pJetsConsidered = i;
        break;
      }
      jetIndices.push_back(i);
    }
  }
  
  std::vector<int> combi;
  for(unsigned int i=0; i<nPartons; ++i) {
    if(useOnlyMatch_) combi.push_back( match[i] );
    else combi.push_back(i);
  }

  std::list<KinFitResult> FitResultList;

  do{
    for(int cnt = 0; cnt < TMath::Factorial( combi.size() ); ++cnt){
      // take into account indistinguishability of the two jets from the hadr. W decay,
      // reduces combinatorics by a factor of 2
      if( (combi[TtSemiLepEvtPartons::LightQ] < combi[TtSemiLepEvtPartons::LightQBar]
         || useOnlyMatch_ ) && doBTagging(useBTag_, jets, combi, bTagAlgo_, minBTagValueBJet_, maxBTagValueNonBJet_) ){
        
        std::vector<pat::Jet> jetCombi;
        jetCombi.resize(nPartons);
        jetCombi[TtSemiLepEvtPartons::LightQ   ] = (*jets)[combi[TtSemiLepEvtPartons::LightQ   ]];
        jetCombi[TtSemiLepEvtPartons::LightQBar] = (*jets)[combi[TtSemiLepEvtPartons::LightQBar]];
        jetCombi[TtSemiLepEvtPartons::HadB     ] = (*jets)[combi[TtSemiLepEvtPartons::HadB     ]];
        jetCombi[TtSemiLepEvtPartons::LepB     ] = (*jets)[combi[TtSemiLepEvtPartons::LepB     ]];

        // do the kinematic fit
        const int status = fitter->fit(jetCombi, (*leps)[0], (*mets)[0]);

        if( status == 0 ) { // only take into account converged fits
          KinFitResult result;
          result.Status = status;
          result.Chi2 = fitter->fitS();
          result.Prob = fitter->fitProb();
          result.HadB = fitter->fittedHadB();
          result.HadP = fitter->fittedHadP();
          result.HadQ = fitter->fittedHadQ();
          result.LepB = fitter->fittedLepB();
          result.LepL = fitter->fittedLepton();
          result.LepN = fitter->fittedNeutrino();
          result.JetCombi = combi;

          FitResultList.push_back(result);
        }

      }
      if(useOnlyMatch_) break; // don't go through combinatorics if useOnlyMatch was chosen
      next_permutation( combi.begin(), combi.end() );
    }
    if(useOnlyMatch_) break; // don't go through combinatorics if useOnlyMatch was chosen
  }
  while(stdcomb::next_combination( jetIndices.begin(), jetIndices.end(), combi.begin(), combi.end() ));

  // sort results w.r.t. chi2 values
  FitResultList.sort();
  
  // -----------------------------------------------------
  // feed out result
  // starting with the JetComb having the smallest chi2
  // -----------------------------------------------------

  if( (unsigned)FitResultList.size() < 1 ) { // in case no fit results were stored in the list (all fits aborted)
    pPartonsHadP->push_back( fitter->fittedHadP()     );
    pPartonsHadQ->push_back( fitter->fittedHadQ()     );
    pPartonsHadB->push_back( fitter->fittedHadB()     );
    pPartonsLepB->push_back( fitter->fittedLepB()     );
    pLeptons    ->push_back( fitter->fittedLepton()   );
    pNeutrinos  ->push_back( fitter->fittedNeutrino() );
    // indices referring to the jet combination
    std::vector<int> invalidCombi;
    for(unsigned int i = 0; i < nPartons; ++i) 
      invalidCombi.push_back( -1 );
    pCombi->push_back( invalidCombi );
    // chi2
    pChi2->push_back( -1. );
    // chi2 probability
    pProb->push_back( -1. );
    // status of the fitter
    pStatus->push_back( -1 );
  }
  else {
    unsigned int iComb = 0;
    for(typename std::list<KinFitResult>::const_iterator result = FitResultList.begin(); result != FitResultList.end(); ++result) {
      if(maxNComb_ >= 1 && iComb == (unsigned int) maxNComb_) break;
      iComb++;
      // partons
      pPartonsHadP->push_back( result->HadP );
      pPartonsHadQ->push_back( result->HadQ );
      pPartonsHadB->push_back( result->HadB );
      pPartonsLepB->push_back( result->LepB );
      // lepton
      pLeptons->push_back( result->LepL );
      // neutrino
      pNeutrinos->push_back( result->LepN );
      // indices referring to the jet combination
      pCombi->push_back( result->JetCombi );
      // chi2
      pChi2->push_back( result->Chi2 );
      // chi2 probability
      pProb->push_back( result->Prob );
      // status of the fitter
      pStatus->push_back( result->Status );
    }
  }
  evt.put(pCombi);
  evt.put(pPartonsHadP, "PartonsHadP");
  evt.put(pPartonsHadQ, "PartonsHadQ");
  evt.put(pPartonsHadB, "PartonsHadB");
  evt.put(pPartonsLepB, "PartonsLepB");
  evt.put(pLeptons    , "Leptons"    );
  evt.put(pNeutrinos  , "Neutrinos"  );
  evt.put(pChi2       , "Chi2"       );
  evt.put(pProb       , "Prob"       );
  evt.put(pStatus     , "Status"     );
  evt.put(pJetsConsidered, "NumberOfConsideredJets");
}

Member Data Documentation

template<typename LeptonCollection >
std::vector<edm::ParameterSet> TtSemiLepKinFitProducer< LeptonCollection >::bResolutions_ [private]
template<typename LeptonCollection >
std::string TtSemiLepKinFitProducer< LeptonCollection >::bTagAlgo_ [private]

input tag for b-tagging algorithm

Definition at line 42 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
std::vector<unsigned> TtSemiLepKinFitProducer< LeptonCollection >::constraints_ [private]
template<typename LeptonCollection >
TtSemiLepKinFitter* TtSemiLepKinFitProducer< LeptonCollection >::fitter [private]
template<typename LeptonCollection >
std::vector<double> TtSemiLepKinFitProducer< LeptonCollection >::jetEnergyResolutionEtaBinning_ [private]
template<typename LeptonCollection >
std::vector<double> TtSemiLepKinFitProducer< LeptonCollection >::jetEnergyResolutionScaleFactors_ [private]

scale factors for jet energy resolution

Definition at line 68 of file TtSemiLepKinFitProducer.h.

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

template<typename LeptonCollection >
unsigned int TtSemiLepKinFitProducer< LeptonCollection >::jetParam_ [private]
template<typename LeptonCollection >
edm::InputTag TtSemiLepKinFitProducer< LeptonCollection >::jets_ [private]

Definition at line 34 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
unsigned int TtSemiLepKinFitProducer< LeptonCollection >::lepParam_ [private]
template<typename LeptonCollection >
std::vector<edm::ParameterSet> TtSemiLepKinFitProducer< LeptonCollection >::lepResolutions_ [private]
template<typename LeptonCollection >
edm::InputTag TtSemiLepKinFitProducer< LeptonCollection >::leps_ [private]

Definition at line 35 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
edm::InputTag TtSemiLepKinFitProducer< LeptonCollection >::match_ [private]

Definition at line 38 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepKinFitProducer< LeptonCollection >::maxBTagValueNonBJet_ [private]

max value of bTag for a non-b-jet

Definition at line 46 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepKinFitProducer< LeptonCollection >::maxDeltaS_ [private]

maximal chi2 equivalent

Definition at line 57 of file TtSemiLepKinFitProducer.h.

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

template<typename LeptonCollection >
double TtSemiLepKinFitProducer< LeptonCollection >::maxF_ [private]

maximal deviation for contstraints

Definition at line 59 of file TtSemiLepKinFitProducer.h.

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

template<typename LeptonCollection >
int TtSemiLepKinFitProducer< LeptonCollection >::maxNComb_ [private]

maximal number of combinations to be written to the event

Definition at line 52 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
int TtSemiLepKinFitProducer< LeptonCollection >::maxNJets_ [private]

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

Definition at line 50 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
unsigned int TtSemiLepKinFitProducer< LeptonCollection >::maxNrIter_ [private]

maximal number of iterations to be performed for the fit

Definition at line 55 of file TtSemiLepKinFitProducer.h.

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

template<typename LeptonCollection >
unsigned int TtSemiLepKinFitProducer< LeptonCollection >::metParam_ [private]
template<typename LeptonCollection >
std::vector<edm::ParameterSet> TtSemiLepKinFitProducer< LeptonCollection >::metResolutions_ [private]
template<typename LeptonCollection >
edm::InputTag TtSemiLepKinFitProducer< LeptonCollection >::mets_ [private]

Definition at line 36 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepKinFitProducer< LeptonCollection >::minBTagValueBJet_ [private]

min value of bTag for a b-jet

Definition at line 44 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
double TtSemiLepKinFitProducer< LeptonCollection >::mTop_ [private]
template<typename LeptonCollection >
double TtSemiLepKinFitProducer< LeptonCollection >::mW_ [private]
template<typename LeptonCollection >
std::vector<edm::ParameterSet> TtSemiLepKinFitProducer< LeptonCollection >::udscResolutions_ [private]

config-file-based object resolutions

Definition at line 71 of file TtSemiLepKinFitProducer.h.

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

template<typename LeptonCollection >
bool TtSemiLepKinFitProducer< LeptonCollection >::useBTag_ [private]

switch to tell whether to use b-tagging or not

Definition at line 48 of file TtSemiLepKinFitProducer.h.

template<typename LeptonCollection >
bool TtSemiLepKinFitProducer< LeptonCollection >::useOnlyMatch_ [private]

switch to use only a combination given by another hypothesis

Definition at line 40 of file TtSemiLepKinFitProducer.h.