28 const std::vector<edm::ParameterSet>* udscResolutions,
29 const std::vector<edm::ParameterSet>* bResolutions,
30 const std::vector<edm::ParameterSet>* lepResolutions,
31 const std::vector<edm::ParameterSet>* metResolutions,
32 const std::vector<double>* jetEnergyResolutionScaleFactors,
33 const std::vector<double>* jetEnergyResolutionEtaBinning)
35 udscResolutions_(udscResolutions),
36 bResolutions_(bResolutions),
37 lepResolutions_(lepResolutions),
38 metResolutions_(metResolutions),
39 jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors),
40 jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning),
44 constrList_(constraints) {
51 std::stringstream constr;
55 constr <<
" * hadronic W-mass (" <<
mW_ <<
" GeV) \n";
58 constr <<
" * leptonic W-mass (" << mW_ <<
" GeV) \n";
61 constr <<
" * hadronic t-mass (" <<
mTop_ <<
" GeV) \n";
64 constr <<
" * leptonic t-mass (" << mTop_ <<
" GeV) \n";
67 constr <<
" * neutrino mass (0 GeV) \n";
70 constr <<
" * equal t-masses \n";
73 constr <<
" * summed transverse momentum \n";
78 <<
"+++++++++++ TtSemiLepKinFitter Setup ++++++++++++ \n"
79 <<
" Parametrization: \n"
84 << constr.str() <<
" Max(No iterations): " <<
maxNrIter_ <<
"\n"
86 <<
" Max(F) : " <<
maxF_ <<
"\n"
87 <<
"+++++++++++++++++++++++++++++++++++++++++++++++++ \n";
91 TMatrixD empty3x3(3, 3);
92 TMatrixD empty4x4(4, 4);
95 hadB_ = std::make_unique<TFitParticleEMomDev>(
"Jet1",
"Jet1",
nullptr, &empty4x4);
96 hadP_ = std::make_unique<TFitParticleEMomDev>(
"Jet2",
"Jet2",
nullptr, &empty4x4);
97 hadQ_ = std::make_unique<TFitParticleEMomDev>(
"Jet3",
"Jet3",
nullptr, &empty4x4);
98 lepB_ = std::make_unique<TFitParticleEMomDev>(
"Jet4",
"Jet4",
nullptr, &empty4x4);
101 hadB_ = std::make_unique<TFitParticleEtEtaPhi>(
"Jet1",
"Jet1",
nullptr, &empty3x3);
102 hadP_ = std::make_unique<TFitParticleEtEtaPhi>(
"Jet2",
"Jet2",
nullptr, &empty3x3);
103 hadQ_ = std::make_unique<TFitParticleEtEtaPhi>(
"Jet3",
"Jet3",
nullptr, &empty3x3);
104 lepB_ = std::make_unique<TFitParticleEtEtaPhi>(
"Jet4",
"Jet4",
nullptr, &empty3x3);
107 hadB_ = std::make_unique<TFitParticleEtThetaPhi>(
"Jet1",
"Jet1",
nullptr, &empty3x3);
108 hadP_ = std::make_unique<TFitParticleEtThetaPhi>(
"Jet2",
"Jet2",
nullptr, &empty3x3);
109 hadQ_ = std::make_unique<TFitParticleEtThetaPhi>(
"Jet3",
"Jet3",
nullptr, &empty3x3);
110 lepB_ = std::make_unique<TFitParticleEtThetaPhi>(
"Jet4",
"Jet4",
nullptr, &empty3x3);
116 TMatrixD empty3x3(3, 3);
119 lepton_ = std::make_unique<TFitParticleEScaledMomDev>(
"Lepton",
"Lepton",
nullptr, &empty3x3);
122 lepton_ = std::make_unique<TFitParticleEtEtaPhi>(
"Lepton",
"Lepton",
nullptr, &empty3x3);
125 lepton_ = std::make_unique<TFitParticleEtThetaPhi>(
"Lepton",
"Lepton",
nullptr, &empty3x3);
130 neutrino_ = std::make_unique<TFitParticleEScaledMomDev>(
"Neutrino",
"Neutrino",
nullptr, &empty3x3);
133 neutrino_ = std::make_unique<TFitParticleEtEtaPhi>(
"Neutrino",
"Neutrino",
nullptr, &empty3x3);
136 neutrino_ = std::make_unique<TFitParticleEtThetaPhi>(
"Neutrino",
"Neutrino",
nullptr, &empty3x3);
144 massConstr_[
kTopHadMass] = std::make_unique<TFitConstraintM>(
"TopMassHad",
"TopMassHad",
nullptr,
nullptr,
mTop_);
145 massConstr_[
kTopLepMass] = std::make_unique<TFitConstraintM>(
"TopMassLep",
"TopMassLep",
nullptr,
nullptr,
mTop_);
146 massConstr_[
kNeutrinoMass] = std::make_unique<TFitConstraintM>(
"NeutrinoMass",
"NeutrinoMass",
nullptr,
nullptr, 0.);
148 std::make_unique<TFitConstraintM>(
"EqualTopMasses",
"EqualTopMasses",
nullptr,
nullptr, 0.);
201 covM_ = std::make_unique<CovarianceMatrix>();
205 const TLorentzVector& p4HadQ,
206 const TLorentzVector& p4HadB,
207 const TLorentzVector& p4LepB,
208 const TLorentzVector& p4Lepton,
209 const TLorentzVector& p4Neutrino,
210 const int leptonCharge,
217 TMatrixD covLepton =
covM_->setupMatrix(p4Lepton, leptonType,
lepParam_);
237 const TLorentzVector& p4HadQ,
238 const TLorentzVector& p4HadB,
239 const TLorentzVector& p4LepB,
240 const TLorentzVector& p4Lepton,
241 const TLorentzVector& p4Neutrino,
242 const TMatrixD& covHadP,
243 const TMatrixD& covHadQ,
244 const TMatrixD& covHadB,
245 const TMatrixD& covLepB,
246 const TMatrixD& covLepton,
247 const TMatrixD& covNeutrino,
248 const int leptonCharge) {
250 hadP_->setIni4Vec(&p4HadP);
251 hadQ_->setIni4Vec(&p4HadQ);
252 hadB_->setIni4Vec(&p4HadB);
253 lepB_->setIni4Vec(&p4LepB);
254 lepton_->setIni4Vec(&p4Lepton);
257 hadP_->setCovMatrix(&covHadP);
258 hadQ_->setCovMatrix(&covHadQ);
259 hadB_->setCovMatrix(&covHadB);
260 lepB_->setCovMatrix(&covLepB);
261 lepton_->setCovMatrix(&covLepton);
266 sumPxConstr_->setConstraint(p4HadP.Px() + p4HadQ.Px() + p4HadB.Px() + p4LepB.Px() + p4Lepton.Px() +
268 sumPyConstr_->setConstraint(p4HadP.Py() + p4HadQ.Py() + p4HadB.Py() + p4LepB.Py() + p4Lepton.Py() +
276 if (
fitter_->getStatus() == 0) {
281 hadP_->getCurr4Vec()->X(),
hadP_->getCurr4Vec()->Y(),
hadP_->getCurr4Vec()->Z(),
hadP_->getCurr4Vec()->E()),
286 hadQ_->getCurr4Vec()->X(),
hadQ_->getCurr4Vec()->Y(),
hadQ_->getCurr4Vec()->Z(),
hadQ_->getCurr4Vec()->E()),
291 hadB_->getCurr4Vec()->X(),
hadB_->getCurr4Vec()->Y(),
hadB_->getCurr4Vec()->Z(),
hadB_->getCurr4Vec()->E()),
296 lepB_->getCurr4Vec()->X(),
lepB_->getCurr4Vec()->Y(),
lepB_->getCurr4Vec()->Z(),
lepB_->getCurr4Vec()->E()),
321 std::vector<pat::Jet>
jets;
329 if (fitsol.
getDecay() ==
"electron")
335 if (
fitter_->getStatus() == 0) {
Log< level::Info, true > LogVerbatim
void setFitHadq(const pat::Particle &aFitHadq)
const pat::Particle fittedHadP() const
return hadronic light quark candidate
void setFitHadb(const pat::Particle &aFitHadb)
std::map< Constraint, std::unique_ptr< TFitConstraintM > > massConstr_
supported constraints
std::unique_ptr< TAbsFitParticle > hadQ_
Param
supported parameterizations
const std::vector< edm::ParameterSet > * udscResolutions_
resolutions
const pat::Particle fittedNeutrino() const
return neutrino candidate
const pat::Particle fittedLepB() const
return leptonic b quark candidate
pat::Jet getCalHadq() const
std::unique_ptr< TKinFitter > fitter_
kinematic fitter
std::unique_ptr< TAbsFitParticle > neutrino_
bool constrainSumPt_
internally use simple boolean for this constraint to reduce the per-event computing time ...
std::string param(const Param ¶m) const
convert Param to human readable form
pat::MET getCalLepn() const
std::string getDecay() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
int maxNrIter_
maximal allowed number of iterations to be used for the fit
std::unique_ptr< TAbsFitParticle > lepton_
pat::Particle fittedHadP_
double fitProb() const
return fit probability
pat::Particle fittedLepB_
Param lepParam_
lepton parametrization
pat::Particle fittedHadB_
output particles
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
pat::Electron getCalLepe() const
const std::vector< edm::ParameterSet > * bResolutions_
int fit(const std::vector< pat::Jet > &jets, const pat::Lepton< LeptonType > &leps, const pat::MET &met)
kinematic fit interface for PAT objects
void setupJets()
initialize jet inputs
const std::vector< edm::ParameterSet > * lepResolutions_
void printSetup() const
print fitter setup
void setFitLepb(const pat::Particle &aFitLepb)
pat::Particle fittedLepton_
void setFitHadp(const pat::Particle &aFitHadp)
const std::vector< double > * jetEnergyResolutionScaleFactors_
scale factors for the jet energy resolution
Param jetParam_
jet parametrization
const pat::Particle fittedHadB() const
return hadronic b quark candidate
std::unique_ptr< TFitConstraintEp > sumPyConstr_
double mW_
W mass value used for constraints.
TtSemiLepKinFitter()
default constructor
Param metParam_
met parametrization
pat::Particle fittedHadQ_
pat::Muon getCalLepm() const
double maxDeltaS_
maximal allowed chi2 (not normalized to degrees of freedom)
pat::Jet getCalHadb() const
TtSemiEvtSolution addKinFitInfo(TtSemiEvtSolution *asol)
add kin fit information to the old event solution (in for legacy reasons)
~TtSemiLepKinFitter()
default destructor
void setupConstraints()
initialize constraints
std::vector< Constraint > constrList_
vector of constraints to be used
std::unique_ptr< TAbsFitParticle > lepB_
void setupLeptons()
initialize lepton inputs
void setFitLepn(const pat::Particle &aFitLepn)
XYZPointD XYZPoint
point in space with cartesian internal representation
Analysis-level particle class.
std::unique_ptr< TFitConstraintEp > sumPxConstr_
std::unique_ptr< TAbsFitParticle > hadB_
input particles
const pat::Particle fittedHadQ() const
return hadronic light quark candidate
const pat::Particle fittedLepton() const
return lepton candidate
void setProbChi2(double c)
pat::Jet getCalLepb() const
void setupFitter()
setup fitter
pat::Particle fittedNeutrino_
const std::vector< double > * jetEnergyResolutionEtaBinning_
void setFitLepl(const pat::Particle &aFitLepl)
std::unique_ptr< CovarianceMatrix > covM_
object used to construct the covariance matrices for the individual particles
std::unique_ptr< TAbsFitParticle > hadP_
double mTop_
top mass value used for constraints
double maxF_
maximal allowed distance from constraints
pat::Jet getCalHadp() const
const std::vector< edm::ParameterSet > * metResolutions_