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 "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00013
00014 #include "TopQuarkAnalysis/TopKinFitter/interface/TopKinFitter.h"
00015 #include "TopQuarkAnalysis/TopKinFitter/interface/CovarianceMatrix.h"
00016
00017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00018
00019 class TAbsFitParticle;
00020 class TFitConstraintM;
00021 class TFitConstraintEp;
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032 class TtSemiLepKinFitter : public TopKinFitter {
00033
00034 public:
00035
00037 enum Constraint { kWHadMass = 1, kWLepMass, kTopHadMass, kTopLepMass, kNeutrinoMass, kEqualTopMasses, kSumPt };
00038
00039 public:
00041 explicit TtSemiLepKinFitter();
00043 explicit TtSemiLepKinFitter(Param jetParam, Param lepParam, Param metParam, int maxNrIter, double maxDeltaS, double maxF,
00044 const std::vector<Constraint>& constraints, double mW=80.4, double mTop=173.,
00045 const std::vector<edm::ParameterSet>* udscResolutions=0,
00046 const std::vector<edm::ParameterSet>* bResolutions =0,
00047 const std::vector<edm::ParameterSet>* lepResolutions =0,
00048 const std::vector<edm::ParameterSet>* metResolutions =0,
00049 const std::vector<double>* jetEnergyResolutionScaleFactors=0,
00050 const std::vector<double>* jetEnergyResolutionEtaBinning =0);
00052 ~TtSemiLepKinFitter();
00053
00055 template <class LeptonType> int fit(const std::vector<pat::Jet>& jets, const pat::Lepton<LeptonType>& leps, const pat::MET& met);
00057 int fit(const TLorentzVector& p4HadP, const TLorentzVector& p4HadQ, const TLorentzVector& p4HadB, const TLorentzVector& p4LepB,
00058 const TLorentzVector& p4Lepton, const TLorentzVector& p4Neutrino, const int leptonCharge, const CovarianceMatrix::ObjectType leptonType);
00060 int fit(const TLorentzVector& p4HadP, const TLorentzVector& p4HadQ, const TLorentzVector& p4HadB, const TLorentzVector& p4LepB,
00061 const TLorentzVector& p4Lepton, const TLorentzVector& p4Neutrino,
00062 const TMatrixD& covHadP, const TMatrixD& covHadQ, const TMatrixD& covHadB, const TMatrixD& covLepB,
00063 const TMatrixD& covLepton, const TMatrixD& covNeutrino,
00064 const int leptonCharge);
00066 const pat::Particle fittedHadB() const { return (fitter_->getStatus()==0 ? fittedHadB_ : pat::Particle()); };
00068 const pat::Particle fittedHadP() const { return (fitter_->getStatus()==0 ? fittedHadP_ : pat::Particle()); };
00070 const pat::Particle fittedHadQ() const { return (fitter_->getStatus()==0 ? fittedHadQ_ : pat::Particle()); };
00072 const pat::Particle fittedLepB() const { return (fitter_->getStatus()==0 ? fittedLepB_ : pat::Particle()); };
00074 const pat::Particle fittedLepton() const { return (fitter_->getStatus()==0 ? fittedLepton_ : pat::Particle()); };
00076 const pat::Particle fittedNeutrino() const { return (fitter_->getStatus()==0 ? fittedNeutrino_ : pat::Particle()); };
00078 TtSemiEvtSolution addKinFitInfo(TtSemiEvtSolution* asol);
00079
00080 private:
00082 void printSetup() const;
00084 void setupFitter();
00086 void setupJets();
00088 void setupLeptons();
00090 void setupConstraints();
00091
00092 private:
00094 TAbsFitParticle* hadB_;
00095 TAbsFitParticle* hadP_;
00096 TAbsFitParticle* hadQ_;
00097 TAbsFitParticle* lepB_;
00098 TAbsFitParticle* lepton_;
00099 TAbsFitParticle* neutrino_;
00101 const std::vector<edm::ParameterSet>* udscResolutions_;
00102 const std::vector<edm::ParameterSet>* bResolutions_;
00103 const std::vector<edm::ParameterSet>* lepResolutions_;
00104 const std::vector<edm::ParameterSet>* metResolutions_;
00106 const std::vector<double>* jetEnergyResolutionScaleFactors_;
00107 const std::vector<double>* jetEnergyResolutionEtaBinning_;
00109 CovarianceMatrix* covM_;
00111 std::map<Constraint, TFitConstraintM*> massConstr_;
00112 TFitConstraintEp* sumPxConstr_;
00113 TFitConstraintEp* sumPyConstr_;
00115 pat::Particle fittedHadB_;
00116 pat::Particle fittedHadP_;
00117 pat::Particle fittedHadQ_;
00118 pat::Particle fittedLepB_;
00119 pat::Particle fittedLepton_;
00120 pat::Particle fittedNeutrino_;
00122 Param jetParam_;
00124 Param lepParam_;
00126 Param metParam_;
00128 std::vector<Constraint> constrList_;
00130 bool constrainSumPt_;
00131 };
00132
00133 template <class LeptonType>
00134 int TtSemiLepKinFitter::fit(const std::vector<pat::Jet>& jets, const pat::Lepton<LeptonType>& lepton, const pat::MET& neutrino)
00135 {
00136 if( jets.size()<4 )
00137 throw edm::Exception( edm::errors::Configuration, "Cannot run the TtSemiLepKinFitter with less than 4 jets" );
00138
00139
00140 const pat::Jet hadP = jets[TtSemiLepEvtPartons::LightQ ];
00141 const pat::Jet hadQ = jets[TtSemiLepEvtPartons::LightQBar];
00142 const pat::Jet hadB = jets[TtSemiLepEvtPartons::HadB ];
00143 const pat::Jet lepB = jets[TtSemiLepEvtPartons::LepB ];
00144
00145
00146 const TLorentzVector p4HadP( hadP.px(), hadP.py(), hadP.pz(), hadP.energy() );
00147 const TLorentzVector p4HadQ( hadQ.px(), hadQ.py(), hadQ.pz(), hadQ.energy() );
00148 const TLorentzVector p4HadB( hadB.px(), hadB.py(), hadB.pz(), hadB.energy() );
00149 const TLorentzVector p4LepB( lepB.px(), lepB.py(), lepB.pz(), lepB.energy() );
00150 const TLorentzVector p4Lepton ( lepton.px(), lepton.py(), lepton.pz(), lepton.energy() );
00151 const TLorentzVector p4Neutrino( neutrino.px(), neutrino.py(), 0, neutrino.et() );
00152
00153
00154 TMatrixD covHadP = covM_->setupMatrix(hadP, jetParam_);
00155 TMatrixD covHadQ = covM_->setupMatrix(hadQ, jetParam_);
00156 TMatrixD covHadB = covM_->setupMatrix(hadB, jetParam_, "bjets");
00157 TMatrixD covLepB = covM_->setupMatrix(lepB, jetParam_, "bjets");
00158 TMatrixD covLepton = covM_->setupMatrix(lepton , lepParam_);
00159 TMatrixD covNeutrino = covM_->setupMatrix(neutrino, metParam_);
00160
00161
00162 return fit(p4HadP, p4HadQ, p4HadB, p4LepB, p4Lepton, p4Neutrino,
00163 covHadP, covHadQ, covHadB, covLepB, covLepton, covNeutrino,
00164 lepton.charge());
00165 }
00166
00167
00168 #endif