#include <TtSemiLepKinFitter.h>
Public Types | |
enum | Constraint { kWHadMass = 1, kWLepMass, kTopHadMass, kTopLepMass, kNeutrinoMass, kEqualTopMasses, kSumPt } |
supported constraints More... | |
Public Member Functions | |
TtSemiEvtSolution | addKinFitInfo (TtSemiEvtSolution *asol) |
add kin fit information to the old event solution (in for legacy reasons) | |
int | fit (const TLorentzVector &p4HadP, const TLorentzVector &p4HadQ, const TLorentzVector &p4HadB, const TLorentzVector &p4LepB, const TLorentzVector &p4Lepton, const TLorentzVector &p4Neutrino, const int leptonCharge, const CovarianceMatrix::ObjectType leptonType) |
kinematic fit interface for plain 4-vecs | |
template<class LeptonType > | |
int | fit (const std::vector< pat::Jet > &jets, const pat::Lepton< LeptonType > &leps, const pat::MET &met) |
kinematic fit interface for PAT objects | |
int | fit (const TLorentzVector &p4HadP, const TLorentzVector &p4HadQ, const TLorentzVector &p4HadB, const TLorentzVector &p4LepB, const TLorentzVector &p4Lepton, const TLorentzVector &p4Neutrino, const TMatrixD &covHadP, const TMatrixD &covHadQ, const TMatrixD &covHadB, const TMatrixD &covLepB, const TMatrixD &covLepton, const TMatrixD &covNeutrino, const int leptonCharge) |
common core of the fit interface | |
const pat::Particle | fittedHadB () const |
return hadronic b quark candidate | |
const pat::Particle | fittedHadP () const |
return hadronic light quark candidate | |
const pat::Particle | fittedHadQ () const |
return hadronic light quark candidate | |
const pat::Particle | fittedLepB () const |
return leptonic b quark candidate | |
const pat::Particle | fittedLepton () const |
return lepton candidate | |
const pat::Particle | fittedNeutrino () const |
return neutrino candidate | |
TtSemiLepKinFitter () | |
default constructor | |
TtSemiLepKinFitter (Param jetParam, Param lepParam, Param metParam, int maxNrIter, double maxDeltaS, double maxF, const std::vector< Constraint > &constraints, double mW=80.4, double mTop=173., const std::vector< edm::ParameterSet > *udscResolutions=0, const std::vector< edm::ParameterSet > *bResolutions=0, const std::vector< edm::ParameterSet > *lepResolutions=0, const std::vector< edm::ParameterSet > *metResolutions=0, const std::vector< double > *jetEnergyResolutionScaleFactors=0, const std::vector< double > *jetEnergyResolutionEtaBinning=0) | |
constructor initialized with built-in types and class enum's custom parameters | |
~TtSemiLepKinFitter () | |
default destructor | |
Private Member Functions | |
void | printSetup () const |
print fitter setup | |
void | setupConstraints () |
initialize constraints | |
void | setupFitter () |
setup fitter | |
void | setupJets () |
initialize jet inputs | |
void | setupLeptons () |
initialize lepton inputs | |
Private Attributes | |
const std::vector < edm::ParameterSet > * | bResolutions_ |
bool | constrainSumPt_ |
internally use simple boolean for this constraint to reduce the per-event computing time | |
std::vector< Constraint > | constrList_ |
vector of constraints to be used | |
CovarianceMatrix * | covM_ |
object used to construct the covariance matrices for the individual particles | |
pat::Particle | fittedHadB_ |
output particles | |
pat::Particle | fittedHadP_ |
pat::Particle | fittedHadQ_ |
pat::Particle | fittedLepB_ |
pat::Particle | fittedLepton_ |
pat::Particle | fittedNeutrino_ |
TAbsFitParticle * | hadB_ |
input particles | |
TAbsFitParticle * | hadP_ |
TAbsFitParticle * | hadQ_ |
const std::vector< double > * | jetEnergyResolutionEtaBinning_ |
const std::vector< double > * | jetEnergyResolutionScaleFactors_ |
scale factors for the jet energy resolution | |
Param | jetParam_ |
jet parametrization | |
TAbsFitParticle * | lepB_ |
Param | lepParam_ |
lepton parametrization | |
const std::vector < edm::ParameterSet > * | lepResolutions_ |
TAbsFitParticle * | lepton_ |
std::map< Constraint, TFitConstraintM * > | massConstr_ |
supported constraints | |
Param | metParam_ |
met parametrization | |
const std::vector < edm::ParameterSet > * | metResolutions_ |
TAbsFitParticle * | neutrino_ |
TFitConstraintEp * | sumPxConstr_ |
TFitConstraintEp * | sumPyConstr_ |
const std::vector < edm::ParameterSet > * | udscResolutions_ |
resolutions |
Definition at line 30 of file TtSemiLepKinFitter.h.
supported constraints
Definition at line 35 of file TtSemiLepKinFitter.h.
{ kWHadMass = 1, kWLepMass, kTopHadMass, kTopLepMass, kNeutrinoMass, kEqualTopMasses, kSumPt };
TtSemiLepKinFitter::TtSemiLepKinFitter | ( | ) | [explicit] |
default constructor
default configuration is: Parametrization kEMom, Max iterations = 200, deltaS<= 5e-5, maxF<= 1e-4, no constraints
Definition at line 15 of file TtSemiLepKinFitter.cc.
References setupFitter().
: TopKinFitter(), hadB_(0), hadP_(0), hadQ_(0), lepB_(0), lepton_(0), neutrino_(0), udscResolutions_(0), bResolutions_(0), lepResolutions_(0), metResolutions_(0), jetEnergyResolutionScaleFactors_(0), jetEnergyResolutionEtaBinning_(0), jetParam_(kEMom), lepParam_(kEMom), metParam_(kEMom) { setupFitter(); }
TtSemiLepKinFitter::TtSemiLepKinFitter | ( | Param | jetParam, |
Param | lepParam, | ||
Param | metParam, | ||
int | maxNrIter, | ||
double | maxDeltaS, | ||
double | maxF, | ||
const std::vector< Constraint > & | constraints, | ||
double | mW = 80.4 , |
||
double | mTop = 173. , |
||
const std::vector< edm::ParameterSet > * | udscResolutions = 0 , |
||
const std::vector< edm::ParameterSet > * | bResolutions = 0 , |
||
const std::vector< edm::ParameterSet > * | lepResolutions = 0 , |
||
const std::vector< edm::ParameterSet > * | metResolutions = 0 , |
||
const std::vector< double > * | jetEnergyResolutionScaleFactors = 0 , |
||
const std::vector< double > * | jetEnergyResolutionEtaBinning = 0 |
||
) | [explicit] |
constructor initialized with built-in types and class enum's custom parameters
Definition at line 25 of file TtSemiLepKinFitter.cc.
References setupFitter().
: TopKinFitter(maxNrIter, maxDeltaS, maxF, mW, mTop), hadB_(0), hadP_(0), hadQ_(0), lepB_(0), lepton_(0), neutrino_(0), udscResolutions_(udscResolutions), bResolutions_(bResolutions), lepResolutions_(lepResolutions), metResolutions_(metResolutions), jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors), jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning), jetParam_(jetParam), lepParam_(lepParam), metParam_(metParam), constrList_(constraints) { setupFitter(); }
TtSemiLepKinFitter::~TtSemiLepKinFitter | ( | ) |
default destructor
Definition at line 43 of file TtSemiLepKinFitter.cc.
References covM_, hadB_, hadP_, hadQ_, lepB_, lepton_, massConstr_, neutrino_, sumPxConstr_, and sumPyConstr_.
{ delete hadB_; delete hadP_; delete hadQ_; delete lepB_; delete lepton_; delete neutrino_; delete covM_; for(std::map<Constraint, TFitConstraintM*>::iterator it = massConstr_.begin(); it != massConstr_.end(); ++it) delete it->second; delete sumPxConstr_; delete sumPyConstr_; }
TtSemiEvtSolution TtSemiLepKinFitter::addKinFitInfo | ( | TtSemiEvtSolution * | asol | ) |
add kin fit information to the old event solution (in for legacy reasons)
Definition at line 292 of file TtSemiLepKinFitter.cc.
References fit(), TopKinFitter::fitProb(), fittedHadB(), fittedHadP(), fittedHadQ(), fittedLepB(), fittedLepton(), fittedNeutrino(), TopKinFitter::fitter_, TtSemiEvtSolution::getCalHadb(), TtSemiEvtSolution::getCalHadp(), TtSemiEvtSolution::getCalHadq(), TtSemiEvtSolution::getCalLepb(), TtSemiEvtSolution::getCalLepe(), TtSemiEvtSolution::getCalLepm(), TtSemiEvtSolution::getCalLepn(), TtSemiEvtSolution::getDecay(), TKinFitter::getStatus(), TtSemiLepDaughter::HadB, analyzePatCleaning_cfg::jets, TtSemiLepDaughter::LepB, TtFullHadDaughter::LightQ, TtFullHadDaughter::LightQBar, TtSemiEvtSolution::setFitHadb(), TtSemiEvtSolution::setFitHadp(), TtSemiEvtSolution::setFitHadq(), TtSemiEvtSolution::setFitLepb(), TtSemiEvtSolution::setFitLepl(), TtSemiEvtSolution::setFitLepn(), and TtSemiEvtSolution::setProbChi2().
Referenced by TtSemiEvtSolutionMaker::produce().
{ TtSemiEvtSolution fitsol(*asol); std::vector<pat::Jet> jets; jets.resize(4); jets[TtSemiLepEvtPartons::LightQ ] = fitsol.getCalHadp(); jets[TtSemiLepEvtPartons::LightQBar] = fitsol.getCalHadq(); jets[TtSemiLepEvtPartons::HadB ] = fitsol.getCalHadb(); jets[TtSemiLepEvtPartons::LepB ] = fitsol.getCalLepb(); // perform the fit, either using the electron or the muon if(fitsol.getDecay() == "electron") fit( jets, fitsol.getCalLepe(), fitsol.getCalLepn() ); if(fitsol.getDecay() == "muon" ) fit( jets, fitsol.getCalLepm(), fitsol.getCalLepn() ); // add fitted information to the solution if (fitter_->getStatus() == 0) { // fill the fitted particles fitsol.setFitHadb( fittedHadB() ); fitsol.setFitHadp( fittedHadP() ); fitsol.setFitHadq( fittedHadQ() ); fitsol.setFitLepb( fittedLepB() ); fitsol.setFitLepl( fittedLepton() ); fitsol.setFitLepn( fittedNeutrino() ); // store the fit's chi2 probability fitsol.setProbChi2( fitProb() ); } return fitsol; }
int TtSemiLepKinFitter::fit | ( | const std::vector< pat::Jet > & | jets, |
const pat::Lepton< LeptonType > & | leps, | ||
const pat::MET & | met | ||
) |
kinematic fit interface for PAT objects
Definition at line 189 of file TtSemiLepKinFitter.cc.
References edm::errors::Configuration, covM_, reco::LeafCandidate::energy(), reco::LeafCandidate::et(), TtSemiLepDaughter::HadB, jetParam_, TtSemiLepDaughter::LepB, lepParam_, TtFullHadDaughter::LightQ, TtFullHadDaughter::LightQBar, metParam_, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), and CovarianceMatrix::setupMatrix().
Referenced by addKinFitInfo(), and fit().
{ if( jets.size()<4 ) throw edm::Exception( edm::errors::Configuration, "Cannot run the TtSemiLepKinFitter with less than 4 jets" ); // get jets in right order const pat::Jet hadP = jets[TtSemiLepEvtPartons::LightQ ]; const pat::Jet hadQ = jets[TtSemiLepEvtPartons::LightQBar]; const pat::Jet hadB = jets[TtSemiLepEvtPartons::HadB ]; const pat::Jet lepB = jets[TtSemiLepEvtPartons::LepB ]; // initialize particles const TLorentzVector p4HadP( hadP.px(), hadP.py(), hadP.pz(), hadP.energy() ); const TLorentzVector p4HadQ( hadQ.px(), hadQ.py(), hadQ.pz(), hadQ.energy() ); const TLorentzVector p4HadB( hadB.px(), hadB.py(), hadB.pz(), hadB.energy() ); const TLorentzVector p4LepB( lepB.px(), lepB.py(), lepB.pz(), lepB.energy() ); const TLorentzVector p4Lepton ( lepton.px(), lepton.py(), lepton.pz(), lepton.energy() ); const TLorentzVector p4Neutrino( neutrino.px(), neutrino.py(), 0, neutrino.et() ); // initialize covariance matrices TMatrixD covHadP = covM_->setupMatrix(hadP, jetParam_); TMatrixD covHadQ = covM_->setupMatrix(hadQ, jetParam_); TMatrixD covHadB = covM_->setupMatrix(hadB, jetParam_, "bjets"); TMatrixD covLepB = covM_->setupMatrix(lepB, jetParam_, "bjets"); TMatrixD covLepton = covM_->setupMatrix(lepton , lepParam_); TMatrixD covNeutrino = covM_->setupMatrix(neutrino, metParam_); // now do the part that is fully independent of PAT features return fit(p4HadP, p4HadQ, p4HadB, p4LepB, p4Lepton, p4Neutrino, covHadP, covHadQ, covHadB, covLepB, covLepton, covNeutrino, lepton.charge()); }
int TtSemiLepKinFitter::fit | ( | const TLorentzVector & | p4HadP, |
const TLorentzVector & | p4HadQ, | ||
const TLorentzVector & | p4HadB, | ||
const TLorentzVector & | p4LepB, | ||
const TLorentzVector & | p4Lepton, | ||
const TLorentzVector & | p4Neutrino, | ||
const int | leptonCharge, | ||
const CovarianceMatrix::ObjectType | leptonType | ||
) |
kinematic fit interface for plain 4-vecs
Definition at line 222 of file TtSemiLepKinFitter.cc.
References covM_, fit(), jetParam_, CovarianceMatrix::kBJet, CovarianceMatrix::kMet, CovarianceMatrix::kUdscJet, lepParam_, metParam_, and CovarianceMatrix::setupMatrix().
{ // initialize covariance matrices TMatrixD covHadP = covM_->setupMatrix(p4HadP, CovarianceMatrix::kUdscJet, jetParam_); TMatrixD covHadQ = covM_->setupMatrix(p4HadQ, CovarianceMatrix::kUdscJet, jetParam_); TMatrixD covHadB = covM_->setupMatrix(p4HadB, CovarianceMatrix::kBJet, jetParam_); TMatrixD covLepB = covM_->setupMatrix(p4LepB, CovarianceMatrix::kBJet, jetParam_); TMatrixD covLepton = covM_->setupMatrix(p4Lepton , leptonType , lepParam_); TMatrixD covNeutrino = covM_->setupMatrix(p4Neutrino, CovarianceMatrix::kMet , metParam_); // now do the part that is fully independent of PAT features return fit(p4HadP, p4HadQ, p4HadB, p4LepB, p4Lepton, p4Neutrino, covHadP, covHadQ, covHadB, covLepB, covLepton, covNeutrino, leptonCharge); }
int TtSemiLepKinFitter::fit | ( | const TLorentzVector & | p4HadP, |
const TLorentzVector & | p4HadQ, | ||
const TLorentzVector & | p4HadB, | ||
const TLorentzVector & | p4LepB, | ||
const TLorentzVector & | p4Lepton, | ||
const TLorentzVector & | p4Neutrino, | ||
const TMatrixD & | covHadP, | ||
const TMatrixD & | covHadQ, | ||
const TMatrixD & | covHadB, | ||
const TMatrixD & | covLepB, | ||
const TMatrixD & | covLepton, | ||
const TMatrixD & | covNeutrino, | ||
const int | leptonCharge | ||
) |
common core of the fit interface
Definition at line 239 of file TtSemiLepKinFitter.cc.
References constrainSumPt_, TKinFitter::fit(), fittedHadB_, fittedHadP_, fittedHadQ_, fittedLepB_, fittedLepton_, fittedNeutrino_, TopKinFitter::fitter_, TAbsFitParticle::getCurr4Vec(), TKinFitter::getStatus(), hadB_, hadP_, hadQ_, lepB_, lepton_, neutrino_, TFitConstraintEp::setConstraint(), TAbsFitParticle::setCovMatrix(), TAbsFitParticle::setIni4Vec(), sumPxConstr_, and sumPyConstr_.
{ // set the kinematics of the objects to be fitted hadP_->setIni4Vec( &p4HadP ); hadQ_->setIni4Vec( &p4HadQ ); hadB_->setIni4Vec( &p4HadB ); lepB_->setIni4Vec( &p4LepB ); lepton_->setIni4Vec( &p4Lepton ); neutrino_->setIni4Vec( &p4Neutrino ); hadP_->setCovMatrix( &covHadP ); hadQ_->setCovMatrix( &covHadQ ); hadB_->setCovMatrix( &covHadB ); lepB_->setCovMatrix( &covLepB ); lepton_ ->setCovMatrix( &covLepton ); neutrino_->setCovMatrix( &covNeutrino ); if(constrainSumPt_){ // setup Px and Py constraint for curent event configuration so that sum Pt will be conserved sumPxConstr_->setConstraint( p4HadP.Px() + p4HadQ.Px() + p4HadB.Px() + p4LepB.Px() + p4Lepton.Px() + p4Neutrino.Px() ); sumPyConstr_->setConstraint( p4HadP.Py() + p4HadQ.Py() + p4HadB.Py() + p4LepB.Py() + p4Lepton.Py() + p4Neutrino.Py() ); } // now do the fit fitter_->fit(); // read back the resulting particles if the fit converged if(fitter_->getStatus()==0){ // read back jet kinematics fittedHadP_= pat::Particle(reco::LeafCandidate(0, math::XYZTLorentzVector(hadP_->getCurr4Vec()->X(), hadP_->getCurr4Vec()->Y(), hadP_->getCurr4Vec()->Z(), hadP_->getCurr4Vec()->E()), math::XYZPoint())); fittedHadQ_= pat::Particle(reco::LeafCandidate(0, math::XYZTLorentzVector(hadQ_->getCurr4Vec()->X(), hadQ_->getCurr4Vec()->Y(), hadQ_->getCurr4Vec()->Z(), hadQ_->getCurr4Vec()->E()), math::XYZPoint())); fittedHadB_= pat::Particle(reco::LeafCandidate(0, math::XYZTLorentzVector(hadB_->getCurr4Vec()->X(), hadB_->getCurr4Vec()->Y(), hadB_->getCurr4Vec()->Z(), hadB_->getCurr4Vec()->E()), math::XYZPoint())); fittedLepB_= pat::Particle(reco::LeafCandidate(0, math::XYZTLorentzVector(lepB_->getCurr4Vec()->X(), lepB_->getCurr4Vec()->Y(), lepB_->getCurr4Vec()->Z(), lepB_->getCurr4Vec()->E()), math::XYZPoint())); // read back lepton kinematics fittedLepton_= pat::Particle(reco::LeafCandidate(leptonCharge, math::XYZTLorentzVector(lepton_->getCurr4Vec()->X(), lepton_->getCurr4Vec()->Y(), lepton_->getCurr4Vec()->Z(), lepton_->getCurr4Vec()->E()), math::XYZPoint())); // read back the MET kinematics fittedNeutrino_= pat::Particle(reco::LeafCandidate(0, math::XYZTLorentzVector(neutrino_->getCurr4Vec()->X(), neutrino_->getCurr4Vec()->Y(), neutrino_->getCurr4Vec()->Z(), neutrino_->getCurr4Vec()->E()), math::XYZPoint())); } return fitter_->getStatus(); }
const pat::Particle TtSemiLepKinFitter::fittedHadB | ( | ) | const [inline] |
return hadronic b quark candidate
Definition at line 64 of file TtSemiLepKinFitter.h.
References fittedHadB_, TopKinFitter::fitter_, and TKinFitter::getStatus().
Referenced by addKinFitInfo().
{ return (fitter_->getStatus()==0 ? fittedHadB_ : pat::Particle()); };
const pat::Particle TtSemiLepKinFitter::fittedHadP | ( | ) | const [inline] |
return hadronic light quark candidate
Definition at line 66 of file TtSemiLepKinFitter.h.
References fittedHadP_, TopKinFitter::fitter_, and TKinFitter::getStatus().
Referenced by addKinFitInfo().
{ return (fitter_->getStatus()==0 ? fittedHadP_ : pat::Particle()); };
const pat::Particle TtSemiLepKinFitter::fittedHadQ | ( | ) | const [inline] |
return hadronic light quark candidate
Definition at line 68 of file TtSemiLepKinFitter.h.
References fittedHadQ_, TopKinFitter::fitter_, and TKinFitter::getStatus().
Referenced by addKinFitInfo().
{ return (fitter_->getStatus()==0 ? fittedHadQ_ : pat::Particle()); };
const pat::Particle TtSemiLepKinFitter::fittedLepB | ( | ) | const [inline] |
return leptonic b quark candidate
Definition at line 70 of file TtSemiLepKinFitter.h.
References fittedLepB_, TopKinFitter::fitter_, and TKinFitter::getStatus().
Referenced by addKinFitInfo().
{ return (fitter_->getStatus()==0 ? fittedLepB_ : pat::Particle()); };
const pat::Particle TtSemiLepKinFitter::fittedLepton | ( | ) | const [inline] |
return lepton candidate
Definition at line 72 of file TtSemiLepKinFitter.h.
References fittedLepton_, TopKinFitter::fitter_, and TKinFitter::getStatus().
Referenced by addKinFitInfo().
{ return (fitter_->getStatus()==0 ? fittedLepton_ : pat::Particle()); };
const pat::Particle TtSemiLepKinFitter::fittedNeutrino | ( | ) | const [inline] |
return neutrino candidate
Definition at line 74 of file TtSemiLepKinFitter.h.
References fittedNeutrino_, TopKinFitter::fitter_, and TKinFitter::getStatus().
Referenced by addKinFitInfo().
{ return (fitter_->getStatus()==0 ? fittedNeutrino_ : pat::Particle()); };
void TtSemiLepKinFitter::printSetup | ( | ) | const [private] |
print fitter setup
Definition at line 58 of file TtSemiLepKinFitter.cc.
References constrList_, i, jetParam_, kEqualTopMasses, kNeutrinoMass, kSumPt, kTopHadMass, kTopLepMass, kWHadMass, kWLepMass, lepParam_, TopKinFitter::maxDeltaS_, TopKinFitter::maxF_, TopKinFitter::maxNrIter_, metParam_, TopKinFitter::mTop_, TopKinFitter::mW_, and TopKinFitter::param().
Referenced by setupFitter().
{ std::stringstream constr; for(unsigned int i=0; i<constrList_.size(); ++i){ switch(constrList_[i]){ case kWHadMass : constr << " * hadronic W-mass (" << mW_ << " GeV) \n"; break; case kWLepMass : constr << " * leptonic W-mass (" << mW_ << " GeV) \n"; break; case kTopHadMass : constr << " * hadronic t-mass (" << mTop_ << " GeV) \n"; break; case kTopLepMass : constr << " * leptonic t-mass (" << mTop_ << " GeV) \n"; break; case kNeutrinoMass : constr << " * neutrino mass (0 GeV) \n"; break; case kEqualTopMasses : constr << " * equal t-masses \n"; break; case kSumPt : constr << " * summed transverse momentum \n"; break; } } edm::LogVerbatim( "TtSemiLepKinFitter" ) << "\n" << "+++++++++++ TtSemiLepKinFitter Setup ++++++++++++ \n" << " Parametrization: \n" << " * jet : " << param(jetParam_) << "\n" << " * lep : " << param(lepParam_) << "\n" << " * met : " << param(metParam_) << "\n" << " Constraints: \n" << constr.str() << " Max(No iterations): " << maxNrIter_ << "\n" << " Max(deltaS) : " << maxDeltaS_ << "\n" << " Max(F) : " << maxF_ << "\n" << "+++++++++++++++++++++++++++++++++++++++++++++++++ \n"; }
void TtSemiLepKinFitter::setupConstraints | ( | ) | [private] |
initialize constraints
Definition at line 128 of file TtSemiLepKinFitter.cc.
References TFitConstraintEp::addParticles(), constrainSumPt_, constrList_, spr::find(), hadB_, hadP_, hadQ_, kEqualTopMasses, kNeutrinoMass, kSumPt, kTopHadMass, kTopLepMass, kWHadMass, kWLepMass, lepB_, lepton_, massConstr_, TopKinFitter::mTop_, TopKinFitter::mW_, neutrino_, TFitConstraintEp::pX, TFitConstraintEp::pY, sumPxConstr_, and sumPyConstr_.
Referenced by setupFitter().
{ massConstr_[kWHadMass ] = new TFitConstraintM("WMassHad", "WMassHad", 0, 0, mW_ ); massConstr_[kWLepMass ] = new TFitConstraintM("WMassLep", "WMassLep", 0, 0, mW_ ); massConstr_[kTopHadMass ] = new TFitConstraintM("TopMassHad", "TopMassHad", 0, 0, mTop_); massConstr_[kTopLepMass ] = new TFitConstraintM("TopMassLep", "TopMassLep", 0, 0, mTop_); massConstr_[kNeutrinoMass ] = new TFitConstraintM("NeutrinoMass", "NeutrinoMass", 0, 0, 0.); massConstr_[kEqualTopMasses] = new TFitConstraintM("EqualTopMasses","EqualTopMasses",0, 0, 0.); sumPxConstr_ = new TFitConstraintEp("SumPx", "SumPx", 0, TFitConstraintEp::pX, 0.); sumPyConstr_ = new TFitConstraintEp("SumPy", "SumPy", 0, TFitConstraintEp::pY, 0.); massConstr_[kWHadMass ]->addParticles1(hadP_, hadQ_ ); massConstr_[kWLepMass ]->addParticles1(lepton_, neutrino_); massConstr_[kTopHadMass ]->addParticles1(hadP_, hadQ_, hadB_); massConstr_[kTopLepMass ]->addParticles1(lepton_, neutrino_, lepB_); massConstr_[kNeutrinoMass ]->addParticle1 (neutrino_); massConstr_[kEqualTopMasses]->addParticles1(hadP_, hadQ_, hadB_); massConstr_[kEqualTopMasses]->addParticles2(lepton_, neutrino_, lepB_); sumPxConstr_->addParticles(lepton_, neutrino_, hadP_, hadQ_, hadB_, lepB_); sumPyConstr_->addParticles(lepton_, neutrino_, hadP_, hadQ_, hadB_, lepB_); if(std::find(constrList_.begin(), constrList_.end(), kSumPt)!=constrList_.end()) constrainSumPt_ = true; constrainSumPt_ = false; }
void TtSemiLepKinFitter::setupFitter | ( | ) | [private] |
setup fitter
Definition at line 154 of file TtSemiLepKinFitter.cc.
References TKinFitter::addConstraint(), TKinFitter::addMeasParticle(), bResolutions_, constrainSumPt_, constrList_, covM_, TopKinFitter::fitter_, hadB_, hadP_, hadQ_, i, jetEnergyResolutionEtaBinning_, jetEnergyResolutionScaleFactors_, kSumPt, lepB_, lepResolutions_, lepton_, massConstr_, metResolutions_, neutrino_, printSetup(), setupConstraints(), setupJets(), setupLeptons(), sumPxConstr_, sumPyConstr_, and udscResolutions_.
Referenced by TtSemiLepKinFitter().
{ printSetup(); setupJets(); setupLeptons(); setupConstraints(); // add measured particles fitter_->addMeasParticle(hadB_); fitter_->addMeasParticle(hadP_); fitter_->addMeasParticle(hadQ_); fitter_->addMeasParticle(lepB_); fitter_->addMeasParticle(lepton_); fitter_->addMeasParticle(neutrino_); // add constraints for(unsigned int i=0; i<constrList_.size(); i++){ if(constrList_[i]!=kSumPt) fitter_->addConstraint(massConstr_[constrList_[i]]); } if(constrainSumPt_) { fitter_->addConstraint(sumPxConstr_); fitter_->addConstraint(sumPyConstr_); } // initialize helper class used to bring the resolutions into covariance matrices if(udscResolutions_->size() && bResolutions_->size() && lepResolutions_->size() && metResolutions_->size()) covM_ = new CovarianceMatrix(*udscResolutions_, *bResolutions_, *lepResolutions_, *metResolutions_, *jetEnergyResolutionScaleFactors_, *jetEnergyResolutionEtaBinning_); else covM_ = new CovarianceMatrix(); }
void TtSemiLepKinFitter::setupJets | ( | ) | [private] |
initialize jet inputs
Definition at line 87 of file TtSemiLepKinFitter.cc.
References hadB_, hadP_, hadQ_, jetParam_, TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, and lepB_.
Referenced by setupFitter().
{ TMatrixD empty3x3(3,3); TMatrixD empty4x4(4,4); switch(jetParam_){ // setup jets according to parameterization case kEMom : hadB_= new TFitParticleEMomDev ("Jet1", "Jet1", 0, &empty4x4); hadP_= new TFitParticleEMomDev ("Jet2", "Jet2", 0, &empty4x4); hadQ_= new TFitParticleEMomDev ("Jet3", "Jet3", 0, &empty4x4); lepB_= new TFitParticleEMomDev ("Jet4", "Jet4", 0, &empty4x4); break; case kEtEtaPhi : hadB_= new TFitParticleEtEtaPhi ("Jet1", "Jet1", 0, &empty3x3); hadP_= new TFitParticleEtEtaPhi ("Jet2", "Jet2", 0, &empty3x3); hadQ_= new TFitParticleEtEtaPhi ("Jet3", "Jet3", 0, &empty3x3); lepB_= new TFitParticleEtEtaPhi ("Jet4", "Jet4", 0, &empty3x3); break; case kEtThetaPhi : hadB_= new TFitParticleEtThetaPhi("Jet1", "Jet1", 0, &empty3x3); hadP_= new TFitParticleEtThetaPhi("Jet2", "Jet2", 0, &empty3x3); hadQ_= new TFitParticleEtThetaPhi("Jet3", "Jet3", 0, &empty3x3); lepB_= new TFitParticleEtThetaPhi("Jet4", "Jet4", 0, &empty3x3); break; } }
void TtSemiLepKinFitter::setupLeptons | ( | ) | [private] |
initialize lepton inputs
Definition at line 113 of file TtSemiLepKinFitter.cc.
References TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, lepParam_, lepton_, metParam_, and neutrino_.
Referenced by setupFitter().
{ TMatrixD empty3x3(3,3); switch(lepParam_){ // setup lepton according to parameterization case kEMom : lepton_ = new TFitParticleEScaledMomDev("Lepton", "Lepton", 0, &empty3x3); break; case kEtEtaPhi : lepton_ = new TFitParticleEtEtaPhi ("Lepton", "Lepton", 0, &empty3x3); break; case kEtThetaPhi : lepton_ = new TFitParticleEtThetaPhi ("Lepton", "Lepton", 0, &empty3x3); break; } switch(metParam_){ // setup neutrino according to parameterization case kEMom : neutrino_= new TFitParticleEScaledMomDev("Neutrino", "Neutrino", 0, &empty3x3); break; case kEtEtaPhi : neutrino_= new TFitParticleEtEtaPhi ("Neutrino", "Neutrino", 0, &empty3x3); break; case kEtThetaPhi : neutrino_= new TFitParticleEtThetaPhi ("Neutrino", "Neutrino", 0, &empty3x3); break; } }
const std::vector<edm::ParameterSet>* TtSemiLepKinFitter::bResolutions_ [private] |
Definition at line 100 of file TtSemiLepKinFitter.h.
Referenced by setupFitter().
bool TtSemiLepKinFitter::constrainSumPt_ [private] |
internally use simple boolean for this constraint to reduce the per-event computing time
Definition at line 128 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), and setupFitter().
std::vector<Constraint> TtSemiLepKinFitter::constrList_ [private] |
vector of constraints to be used
Definition at line 126 of file TtSemiLepKinFitter.h.
Referenced by printSetup(), setupConstraints(), and setupFitter().
CovarianceMatrix* TtSemiLepKinFitter::covM_ [private] |
object used to construct the covariance matrices for the individual particles
Definition at line 107 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupFitter(), and ~TtSemiLepKinFitter().
pat::Particle TtSemiLepKinFitter::fittedHadB_ [private] |
output particles
Definition at line 113 of file TtSemiLepKinFitter.h.
Referenced by fit(), and fittedHadB().
pat::Particle TtSemiLepKinFitter::fittedHadP_ [private] |
Definition at line 114 of file TtSemiLepKinFitter.h.
Referenced by fit(), and fittedHadP().
pat::Particle TtSemiLepKinFitter::fittedHadQ_ [private] |
Definition at line 115 of file TtSemiLepKinFitter.h.
Referenced by fit(), and fittedHadQ().
pat::Particle TtSemiLepKinFitter::fittedLepB_ [private] |
Definition at line 116 of file TtSemiLepKinFitter.h.
Referenced by fit(), and fittedLepB().
Definition at line 117 of file TtSemiLepKinFitter.h.
Referenced by fit(), and fittedLepton().
Definition at line 118 of file TtSemiLepKinFitter.h.
Referenced by fit(), and fittedNeutrino().
TAbsFitParticle* TtSemiLepKinFitter::hadB_ [private] |
input particles
Definition at line 92 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), setupJets(), and ~TtSemiLepKinFitter().
TAbsFitParticle* TtSemiLepKinFitter::hadP_ [private] |
Definition at line 93 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), setupJets(), and ~TtSemiLepKinFitter().
TAbsFitParticle* TtSemiLepKinFitter::hadQ_ [private] |
Definition at line 94 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), setupJets(), and ~TtSemiLepKinFitter().
const std::vector<double>* TtSemiLepKinFitter::jetEnergyResolutionEtaBinning_ [private] |
Definition at line 105 of file TtSemiLepKinFitter.h.
Referenced by setupFitter().
const std::vector<double>* TtSemiLepKinFitter::jetEnergyResolutionScaleFactors_ [private] |
scale factors for the jet energy resolution
Definition at line 104 of file TtSemiLepKinFitter.h.
Referenced by setupFitter().
Param TtSemiLepKinFitter::jetParam_ [private] |
jet parametrization
Definition at line 120 of file TtSemiLepKinFitter.h.
Referenced by fit(), printSetup(), and setupJets().
TAbsFitParticle* TtSemiLepKinFitter::lepB_ [private] |
Definition at line 95 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), setupJets(), and ~TtSemiLepKinFitter().
Param TtSemiLepKinFitter::lepParam_ [private] |
lepton parametrization
Definition at line 122 of file TtSemiLepKinFitter.h.
Referenced by fit(), printSetup(), and setupLeptons().
const std::vector<edm::ParameterSet>* TtSemiLepKinFitter::lepResolutions_ [private] |
Definition at line 101 of file TtSemiLepKinFitter.h.
Referenced by setupFitter().
TAbsFitParticle* TtSemiLepKinFitter::lepton_ [private] |
Definition at line 96 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), setupLeptons(), and ~TtSemiLepKinFitter().
std::map<Constraint, TFitConstraintM*> TtSemiLepKinFitter::massConstr_ [private] |
supported constraints
Definition at line 109 of file TtSemiLepKinFitter.h.
Referenced by setupConstraints(), setupFitter(), and ~TtSemiLepKinFitter().
Param TtSemiLepKinFitter::metParam_ [private] |
met parametrization
Definition at line 124 of file TtSemiLepKinFitter.h.
Referenced by fit(), printSetup(), and setupLeptons().
const std::vector<edm::ParameterSet>* TtSemiLepKinFitter::metResolutions_ [private] |
Definition at line 102 of file TtSemiLepKinFitter.h.
Referenced by setupFitter().
TAbsFitParticle* TtSemiLepKinFitter::neutrino_ [private] |
Definition at line 97 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), setupLeptons(), and ~TtSemiLepKinFitter().
TFitConstraintEp* TtSemiLepKinFitter::sumPxConstr_ [private] |
Definition at line 110 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), and ~TtSemiLepKinFitter().
TFitConstraintEp* TtSemiLepKinFitter::sumPyConstr_ [private] |
Definition at line 111 of file TtSemiLepKinFitter.h.
Referenced by fit(), setupConstraints(), setupFitter(), and ~TtSemiLepKinFitter().
const std::vector<edm::ParameterSet>* TtSemiLepKinFitter::udscResolutions_ [private] |