00001 #ifndef TtSemiLepKinFitter_h
00002 #define TtSemiLepKinFitter_h
00003
00004 #include <vector>
00005
00006 #include "TMath.h"
00007 #include "TMatrixD.h"
00008 #include "TLorentzVector.h"
00009
00010 #include "DataFormats/PatCandidates/interface/Lepton.h"
00011 #include "PhysicsTools/KinFitter/interface/TKinFitter.h"
00012
00013 #include "AnalysisDataFormats/TopObjects/interface/TtSemiEvtSolution.h"
00014
00015 class TAbsFitParticle;
00016 class TFitConstraintM;
00017
00018
00019 class TtSemiLepKinFitter {
00020
00021 public:
00022
00023
00024 enum Constraint
00025 { kWHadMass = 1,
00026 kWLepMass,
00027 kTopHadMass,
00028 kTopLepMass,
00029 kNeutrinoMass
00030 };
00031
00032 enum Param
00033 { kEMom,
00034 kEtEtaPhi,
00035 kEtThetaPhi
00036 };
00037
00038 public:
00039
00040
00041 explicit TtSemiLepKinFitter();
00042
00043 explicit TtSemiLepKinFitter(Param jetParam,
00044 Param lepParam,
00045 Param metParam,
00046 int maxNrIter,
00047 double maxDeltaS,
00048 double maxF,
00049 std::vector<Constraint> constraints);
00050
00051 ~TtSemiLepKinFitter();
00052
00053
00054 template <class LeptonType> int fit(const std::vector<pat::Jet> & jets, const pat::Lepton<LeptonType>& leps, const pat::MET& met);
00055
00056
00057 const pat::Particle fittedHadB() const { return (fitter_->getStatus()==0 ? fittedHadB_ : pat::Particle()); };
00058 const pat::Particle fittedHadP() const { return (fitter_->getStatus()==0 ? fittedHadP_ : pat::Particle()); };
00059 const pat::Particle fittedHadQ() const { return (fitter_->getStatus()==0 ? fittedHadQ_ : pat::Particle()); };
00060 const pat::Particle fittedLepB() const { return (fitter_->getStatus()==0 ? fittedLepB_ : pat::Particle()); };
00061 const pat::Particle fittedLepton() const { return (fitter_->getStatus()==0 ? fittedLepton_ : pat::Particle()); };
00062 const pat::Particle fittedNeutrino() const { return (fitter_->getStatus()==0 ? fittedNeutrino_ : pat::Particle()); };
00063
00064
00065 double fitS() const { return fitter_->getS(); };
00066 int fitNrIter() const { return fitter_->getNbIter(); };
00067 double fitProb() const { return TMath::Prob(fitter_->getS(), fitter_->getNDF()); };
00068
00069
00070 TtSemiEvtSolution addKinFitInfo(TtSemiEvtSolution* asol);
00071
00072 private:
00073
00074 void printSetup();
00075 void setupFitter();
00076 void setupJets();
00077 void setupLeptons();
00078 void setupConstraints();
00079
00080 std::string param(Param& param){
00081 std::string parName;
00082 switch(param){
00083 case kEMom : parName="EMom"; break;
00084 case kEtEtaPhi : parName="EtEtaPhi"; break;
00085 case kEtThetaPhi : parName="EtThetaPhi"; break;
00086 }
00087 return parName;
00088 }
00089 std::vector<float> translateCovM(TMatrixD &);
00090
00091 private:
00092
00093
00094 TKinFitter* fitter_;
00095
00096
00097 TAbsFitParticle* hadB_;
00098 TAbsFitParticle* hadP_;
00099 TAbsFitParticle* hadQ_;
00100 TAbsFitParticle* lepB_;
00101 TAbsFitParticle* lepton_;
00102 TAbsFitParticle* neutrino_;
00103
00104
00105 std::map<Constraint, TFitConstraintM*> massConstr_;
00106
00107
00108 pat::Particle fittedHadB_;
00109 pat::Particle fittedHadP_;
00110 pat::Particle fittedHadQ_;
00111 pat::Particle fittedLepB_;
00112 pat::Particle fittedLepton_;
00113 pat::Particle fittedNeutrino_;
00114
00115
00116 Param jetParam_, lepParam_, metParam_;
00117 int maxNrIter_;
00118 double maxDeltaS_;
00119 double maxF_;
00120 std::vector<Constraint> constrList_;
00121 };
00122
00123 #endif