Go to the documentation of this file.00001 #ifndef TtSemiLepKinFitter_h
00002 #define TtSemiLepKinFitter_h
00003
00004 #include <vector>
00005
00006 #include "TLorentzVector.h"
00007
00008 #include "DataFormats/PatCandidates/interface/Lepton.h"
00009
00010 #include "AnalysisDataFormats/TopObjects/interface/TtSemiEvtSolution.h"
00011
00012 #include "TopQuarkAnalysis/TopKinFitter/interface/TopKinFitter.h"
00013 #include "TopQuarkAnalysis/TopKinFitter/interface/CovarianceMatrix.h"
00014
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016
00017 class TAbsFitParticle;
00018 class TFitConstraintM;
00019 class TFitConstraintEp;
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030 class TtSemiLepKinFitter : public TopKinFitter {
00031
00032 public:
00033
00035 enum Constraint { kWHadMass = 1, kWLepMass, kTopHadMass, kTopLepMass, kNeutrinoMass, kEqualTopMasses, kSumPt };
00036
00037 public:
00039 explicit TtSemiLepKinFitter();
00041 explicit TtSemiLepKinFitter(Param jetParam, Param lepParam, Param metParam, int maxNrIter, double maxDeltaS, double maxF,
00042 const std::vector<Constraint>& constraints, double mW=80.4, double mTop=173.,
00043 const std::vector<edm::ParameterSet>* udscResolutions=0,
00044 const std::vector<edm::ParameterSet>* bResolutions =0,
00045 const std::vector<edm::ParameterSet>* lepResolutions =0,
00046 const std::vector<edm::ParameterSet>* metResolutions =0,
00047 const std::vector<double>* jetEnergyResolutionScaleFactors=0,
00048 const std::vector<double>* jetEnergyResolutionEtaBinning =0);
00050 ~TtSemiLepKinFitter();
00051
00053 template <class LeptonType> int fit(const std::vector<pat::Jet>& jets, const pat::Lepton<LeptonType>& leps, const pat::MET& met);
00055 int fit(const TLorentzVector& p4HadP, const TLorentzVector& p4HadQ, const TLorentzVector& p4HadB, const TLorentzVector& p4LepB,
00056 const TLorentzVector& p4Lepton, const TLorentzVector& p4Neutrino, const int leptonCharge, const CovarianceMatrix::ObjectType leptonType);
00058 int fit(const TLorentzVector& p4HadP, const TLorentzVector& p4HadQ, const TLorentzVector& p4HadB, const TLorentzVector& p4LepB,
00059 const TLorentzVector& p4Lepton, const TLorentzVector& p4Neutrino,
00060 const TMatrixD& covHadP, const TMatrixD& covHadQ, const TMatrixD& covHadB, const TMatrixD& covLepB,
00061 const TMatrixD& covLepton, const TMatrixD& covNeutrino,
00062 const int leptonCharge);
00064 const pat::Particle fittedHadB() const { return (fitter_->getStatus()==0 ? fittedHadB_ : pat::Particle()); };
00066 const pat::Particle fittedHadP() const { return (fitter_->getStatus()==0 ? fittedHadP_ : pat::Particle()); };
00068 const pat::Particle fittedHadQ() const { return (fitter_->getStatus()==0 ? fittedHadQ_ : pat::Particle()); };
00070 const pat::Particle fittedLepB() const { return (fitter_->getStatus()==0 ? fittedLepB_ : pat::Particle()); };
00072 const pat::Particle fittedLepton() const { return (fitter_->getStatus()==0 ? fittedLepton_ : pat::Particle()); };
00074 const pat::Particle fittedNeutrino() const { return (fitter_->getStatus()==0 ? fittedNeutrino_ : pat::Particle()); };
00076 TtSemiEvtSolution addKinFitInfo(TtSemiEvtSolution* asol);
00077
00078 private:
00080 void printSetup() const;
00082 void setupFitter();
00084 void setupJets();
00086 void setupLeptons();
00088 void setupConstraints();
00089
00090 private:
00092 TAbsFitParticle* hadB_;
00093 TAbsFitParticle* hadP_;
00094 TAbsFitParticle* hadQ_;
00095 TAbsFitParticle* lepB_;
00096 TAbsFitParticle* lepton_;
00097 TAbsFitParticle* neutrino_;
00099 const std::vector<edm::ParameterSet>* udscResolutions_;
00100 const std::vector<edm::ParameterSet>* bResolutions_;
00101 const std::vector<edm::ParameterSet>* lepResolutions_;
00102 const std::vector<edm::ParameterSet>* metResolutions_;
00104 const std::vector<double>* jetEnergyResolutionScaleFactors_;
00105 const std::vector<double>* jetEnergyResolutionEtaBinning_;
00107 CovarianceMatrix* covM_;
00109 std::map<Constraint, TFitConstraintM*> massConstr_;
00110 TFitConstraintEp* sumPxConstr_;
00111 TFitConstraintEp* sumPyConstr_;
00113 pat::Particle fittedHadB_;
00114 pat::Particle fittedHadP_;
00115 pat::Particle fittedHadQ_;
00116 pat::Particle fittedLepB_;
00117 pat::Particle fittedLepton_;
00118 pat::Particle fittedNeutrino_;
00120 Param jetParam_;
00122 Param lepParam_;
00124 Param metParam_;
00126 std::vector<Constraint> constrList_;
00128 bool constrainSumPt_;
00129 };
00130
00131 #endif