#include <TtFullHadKinFitProducer.h>
Public Member Functions | |
TtFullHadKinFitProducer (const edm::ParameterSet &cfg) | |
default constructor | |
~TtFullHadKinFitProducer () | |
default destructor | |
Public Attributes | |
TtFullHadKinFitter::KinFit * | kinFitter |
kinematic fit interface | |
Private Member Functions | |
virtual void | produce (edm::Event &event, const edm::EventSetup &setup) |
produce fitted object collections and meta data describing fit quality | |
Private Attributes | |
std::vector< edm::ParameterSet > | bResolutions_ |
std::string | bTagAlgo_ |
input tag for b-tagging algorithm | |
unsigned int | bTags_ |
minimal number of b-jets | |
std::vector< unsigned > | constraints_ |
numbering of different possible kinematic constraints | |
double | energyResolutionSmearFactor_ |
smearing factor for jet energy resolutions | |
std::string | jetCorrectionLevel_ |
correction level for jets | |
unsigned int | jetParam_ |
numbering of different possible jet parametrizations | |
edm::InputTag | jets_ |
input tag for jets | |
edm::InputTag | match_ |
input tag for matches (in case the fit should be performed on certain matches) | |
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 | |
double | minBTagValueBJet_ |
min value of bTag for a b-jet | |
double | mTop_ |
top mass value used for constraints | |
double | mW_ |
W mass value used for constraints. | |
std::vector< edm::ParameterSet > | udscResolutions_ |
store the resolutions for the jets | |
bool | useBTagging_ |
switch to tell whether to use b-tagging or not | |
bool | useOnlyMatch_ |
Definition at line 22 of file TtFullHadKinFitProducer.h.
TtFullHadKinFitProducer::TtFullHadKinFitProducer | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
default constructor
Definition at line 6 of file TtFullHadKinFitProducer.cc.
References bResolutions_, bTagAlgo_, bTags_, constraints_, energyResolutionSmearFactor_, Exception, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), jetCorrectionLevel_, jetParam_, kinFitter, maxBTagValueNonBJet_, maxDeltaS_, maxF_, maxNComb_, maxNJets_, maxNrIter_, minBTagValueBJet_, mTop_, mW_, udscResolutions_, and useBTagging_.
: jets_ (cfg.getParameter<edm::InputTag>("jets")), match_ (cfg.getParameter<edm::InputTag>("match")), useOnlyMatch_ (cfg.getParameter<bool>("useOnlyMatch")), bTagAlgo_ (cfg.getParameter<std::string>("bTagAlgo")), minBTagValueBJet_ (cfg.getParameter<double>("minBTagValueBJet")), maxBTagValueNonBJet_ (cfg.getParameter<double>("maxBTagValueNonBJet")), useBTagging_ (cfg.getParameter<bool>("useBTagging")), bTags_ (cfg.getParameter<unsigned int>("bTags")), jetCorrectionLevel_ (cfg.getParameter<std::string>("jetCorrectionLevel")), maxNJets_ (cfg.getParameter<int>("maxNJets")), maxNComb_ (cfg.getParameter<int>("maxNComb")), maxNrIter_ (cfg.getParameter<unsigned int>("maxNrIter")), maxDeltaS_ (cfg.getParameter<double>("maxDeltaS")), maxF_ (cfg.getParameter<double>("maxF")), jetParam_ (cfg.getParameter<unsigned>("jetParametrisation")), constraints_ (cfg.getParameter<std::vector<unsigned> >("constraints")), mW_ (cfg.getParameter<double>("mW" )), mTop_ (cfg.getParameter<double>("mTop")), energyResolutionSmearFactor_(cfg.getParameter<double>("energyResolutionSmearFactor")) { if(cfg.exists("udscResolutions") && cfg.exists("bResolutions")){ udscResolutions_ = cfg.getParameter <std::vector<edm::ParameterSet> >("udscResolutions"); bResolutions_ = cfg.getParameter <std::vector<edm::ParameterSet> >("bResolutions"); } else if(cfg.exists("udscResolutions") || cfg.exists("bResolutions")){ if(cfg.exists("udscResolutions")) throw cms::Exception("WrongConfig") << "Parameter 'bResolutions' is needed if parameter 'udscResolutions' is defined!\n"; else throw cms::Exception("WrongConfig") << "Parameter 'udscResolutions' is needed if parameter 'bResolutions' is defined!\n"; } // define kinematic fit interface kinFitter = new TtFullHadKinFitter::KinFit(useBTagging_, bTags_, bTagAlgo_, minBTagValueBJet_, maxBTagValueNonBJet_, udscResolutions_, bResolutions_, energyResolutionSmearFactor_ , jetCorrectionLevel_, maxNJets_, maxNComb_, maxNrIter_, maxDeltaS_, maxF_, jetParam_, constraints_, mW_, mTop_); // produces the following collections produces< std::vector<pat::Particle> >("PartonsB"); produces< std::vector<pat::Particle> >("PartonsBBar"); produces< std::vector<pat::Particle> >("PartonsLightQ"); produces< std::vector<pat::Particle> >("PartonsLightQBar"); produces< std::vector<pat::Particle> >("PartonsLightP"); produces< std::vector<pat::Particle> >("PartonsLightPBar"); produces< std::vector<std::vector<int> > >(); produces< std::vector<double> >("Chi2"); produces< std::vector<double> >("Prob"); produces< std::vector<int> >("Status"); }
TtFullHadKinFitProducer::~TtFullHadKinFitProducer | ( | ) |
default destructor
Definition at line 57 of file TtFullHadKinFitProducer.cc.
References kinFitter.
{ delete kinFitter; }
void TtFullHadKinFitProducer::produce | ( | edm::Event & | event, |
const edm::EventSetup & | setup | ||
) | [private, virtual] |
produce fitted object collections and meta data describing fit quality
set match to be used
set the validity of a match
Implements edm::EDProducer.
Definition at line 64 of file TtFullHadKinFitProducer.cc.
References TtFullHadKinFitter::KinFit::fit(), analyzePatCleaning_cfg::jets, jets_, kinFitter, match(), match_, maxNComb_, nPartons, TtFullHadKinFitter::KinFit::setMatch(), TtFullHadKinFitter::KinFit::setMatchInvalidity(), TtFullHadKinFitter::KinFit::setUseOnlyMatch(), and useOnlyMatch_.
{ // get jet collection edm::Handle<std::vector<pat::Jet> > jets; event.getByLabel(jets_, jets); // get match in case that useOnlyMatch_ is true std::vector<int> match; bool invalidMatch=false; if(useOnlyMatch_) { kinFitter->setUseOnlyMatch(true); // in case that only a ceratin match should be used, get match here edm::Handle<std::vector<std::vector<int> > > matches; event.getByLabel(match_, matches); match = *(matches->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; } } } kinFitter->setMatch(match); } kinFitter->setMatchInvalidity(invalidMatch); std::list<TtFullHadKinFitter::KinFitResult> fitResults = kinFitter->fit(*jets); // pointer for output collections std::auto_ptr< std::vector<pat::Particle> > pPartonsB( new std::vector<pat::Particle> ); std::auto_ptr< std::vector<pat::Particle> > pPartonsBBar( new std::vector<pat::Particle> ); std::auto_ptr< std::vector<pat::Particle> > pPartonsLightQ ( new std::vector<pat::Particle> ); std::auto_ptr< std::vector<pat::Particle> > pPartonsLightQBar( new std::vector<pat::Particle> ); std::auto_ptr< std::vector<pat::Particle> > pPartonsLightP ( new std::vector<pat::Particle> ); std::auto_ptr< std::vector<pat::Particle> > pPartonsLightPBar( new std::vector<pat::Particle> ); // pointer for meta information 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> ); unsigned int iComb = 0; for(std::list<TtFullHadKinFitter::KinFitResult>::const_iterator res = fitResults.begin(); res != fitResults.end(); ++res){ if(maxNComb_>=1 && iComb==(unsigned int)maxNComb_){ break; } ++iComb; pPartonsB ->push_back( res->B ); pPartonsBBar ->push_back( res->BBar ); pPartonsLightQ ->push_back( res->LightQ ); pPartonsLightQBar->push_back( res->LightQBar ); pPartonsLightP ->push_back( res->LightP ); pPartonsLightPBar->push_back( res->LightPBar ); pCombi ->push_back( res->JetCombi ); pChi2 ->push_back( res->Chi2 ); pProb ->push_back( res->Prob ); pStatus->push_back( res->Status ); } event.put(pCombi); event.put(pPartonsB , "PartonsB" ); event.put(pPartonsBBar , "PartonsBBar" ); event.put(pPartonsLightQ , "PartonsLightQ" ); event.put(pPartonsLightQBar, "PartonsLightQBar"); event.put(pPartonsLightP , "PartonsLightP" ); event.put(pPartonsLightPBar, "PartonsLightPBar"); event.put(pChi2 , "Chi2" ); event.put(pProb , "Prob" ); event.put(pStatus , "Status" ); }
std::vector<edm::ParameterSet> TtFullHadKinFitProducer::bResolutions_ [private] |
Definition at line 73 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
std::string TtFullHadKinFitProducer::bTagAlgo_ [private] |
input tag for b-tagging algorithm
Definition at line 43 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
unsigned int TtFullHadKinFitProducer::bTags_ [private] |
minimal number of b-jets
Definition at line 51 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
std::vector<unsigned> TtFullHadKinFitProducer::constraints_ [private] |
numbering of different possible kinematic constraints
Definition at line 67 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
double TtFullHadKinFitProducer::energyResolutionSmearFactor_ [private] |
smearing factor for jet energy resolutions
Definition at line 75 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
std::string TtFullHadKinFitProducer::jetCorrectionLevel_ [private] |
correction level for jets
Definition at line 53 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
unsigned int TtFullHadKinFitProducer::jetParam_ [private] |
numbering of different possible jet parametrizations
Definition at line 65 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
edm::InputTag TtFullHadKinFitProducer::jets_ [private] |
input tag for jets
Definition at line 36 of file TtFullHadKinFitProducer.h.
Referenced by produce().
kinematic fit interface
Definition at line 80 of file TtFullHadKinFitProducer.h.
Referenced by produce(), TtFullHadKinFitProducer(), and ~TtFullHadKinFitProducer().
edm::InputTag TtFullHadKinFitProducer::match_ [private] |
input tag for matches (in case the fit should be performed on certain matches)
Definition at line 38 of file TtFullHadKinFitProducer.h.
Referenced by produce().
double TtFullHadKinFitProducer::maxBTagValueNonBJet_ [private] |
max value of bTag for a non-b-jet
Definition at line 47 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
double TtFullHadKinFitProducer::maxDeltaS_ [private] |
maximal chi2 equivalent
Definition at line 61 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
double TtFullHadKinFitProducer::maxF_ [private] |
maximal deviation for contstraints
Definition at line 63 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
int TtFullHadKinFitProducer::maxNComb_ [private] |
maximal number of combinations to be written to the event
Definition at line 57 of file TtFullHadKinFitProducer.h.
Referenced by produce(), and TtFullHadKinFitProducer().
int TtFullHadKinFitProducer::maxNJets_ [private] |
maximal number of jets (-1 possible to indicate 'all')
Definition at line 55 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
unsigned int TtFullHadKinFitProducer::maxNrIter_ [private] |
maximal number of iterations to be performed for the fit
Definition at line 59 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
double TtFullHadKinFitProducer::minBTagValueBJet_ [private] |
min value of bTag for a b-jet
Definition at line 45 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
double TtFullHadKinFitProducer::mTop_ [private] |
top mass value used for constraints
Definition at line 71 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
double TtFullHadKinFitProducer::mW_ [private] |
W mass value used for constraints.
Definition at line 69 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
std::vector<edm::ParameterSet> TtFullHadKinFitProducer::udscResolutions_ [private] |
store the resolutions for the jets
Definition at line 73 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
bool TtFullHadKinFitProducer::useBTagging_ [private] |
switch to tell whether to use b-tagging or not
Definition at line 49 of file TtFullHadKinFitProducer.h.
Referenced by TtFullHadKinFitProducer().
bool TtFullHadKinFitProducer::useOnlyMatch_ [private] |
switch to tell whether all possible combinations should be used for the fit or only a certain combination
Definition at line 41 of file TtFullHadKinFitProducer.h.
Referenced by produce().