CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/TopQuarkAnalysis/TopKinFitter/interface/TtFullHadKinFitter.h

Go to the documentation of this file.
00001 #ifndef TtFullHadKinFitter_h
00002 #define TtFullHadKinFitter_h
00003 
00004 #include <vector>
00005 
00006 #include "TLorentzVector.h"
00007 
00008 #include "AnalysisDataFormats/TopObjects/interface/TtHadEvtSolution.h"
00009 
00010 #include "TopQuarkAnalysis/TopKinFitter/interface/CovarianceMatrix.h"
00011 #include "TopQuarkAnalysis/TopKinFitter/interface/TopKinFitter.h"
00012 
00013 #include "PhysicsTools/JetMCUtils/interface/combination.h"
00014 
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 
00017 class TAbsFitParticle;
00018 class TFitConstraintM;
00019 
00020 /*
00021   \class   TtFullHadKinFitter TtFullHadKinFitter.h "TopQuarkAnalysis/TopKinFitter/interface/TtFullHadKinFitter.h"
00022   
00023   \brief   one line description to be added here...
00024 
00025   text to be added here...
00026   
00027 **/
00028 
00029 class TtFullHadKinFitter : public TopKinFitter {
00030 
00031  public:
00033   enum Constraint{ kWPlusMass=1, kWMinusMass, kTopMass, kTopBarMass, kEqualTopMasses };
00034   
00035  public:
00037   TtFullHadKinFitter();
00039   std::vector<TtFullHadKinFitter::Constraint> intToConstraint(std::vector<unsigned int> constraints);
00041   TtFullHadKinFitter(int jetParam, int maxNrIter, double maxDeltaS, double maxF, std::vector<unsigned int> constraints,
00042                      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<double>* jetEnergyResolutionScaleFactors=0,
00046                      const std::vector<double>* jetEnergyResolutionEtaBinning  =0);
00048   TtFullHadKinFitter(Param jetParam, int maxNrIter, double maxDeltaS, double maxF, std::vector<Constraint> constraints,
00049                      double mW=80.4, double mTop=173.,
00050                      const std::vector<edm::ParameterSet>* udscResolutions=0, 
00051                      const std::vector<edm::ParameterSet>* bResolutions   =0,
00052                      const std::vector<double>* jetEnergyResolutionScaleFactors=0,
00053                      const std::vector<double>* jetEnergyResolutionEtaBinning  =0);
00055   ~TtFullHadKinFitter();
00056 
00058   int fit(const std::vector<pat::Jet>& jets);
00060   const pat::Particle fittedB() const { return (fitter_->getStatus()==0 ? fittedB_ : pat::Particle()); };
00062   const pat::Particle fittedBBar() const { return (fitter_->getStatus()==0 ? fittedBBar_ : pat::Particle()); };
00064   const pat::Particle fittedLightQ() const { return (fitter_->getStatus()==0 ? fittedLightQ_ : pat::Particle()); };
00066   const pat::Particle fittedLightQBar() const { return (fitter_->getStatus()==0 ? fittedLightQBar_ : pat::Particle()); };
00068   const pat::Particle fittedLightP() const { return (fitter_->getStatus()==0 ? fittedLightP_ : pat::Particle()); };
00070   const pat::Particle fittedLightPBar() const { return (fitter_->getStatus()==0 ? fittedLightPBar_ : pat::Particle()); };
00072   TtHadEvtSolution addKinFitInfo(TtHadEvtSolution * asol);
00073   
00074  private:
00076   void printSetup() const;
00078   void setupFitter();
00080   void setupJets();
00082   void setupConstraints();
00083 
00084  private:
00086   TAbsFitParticle* b_;
00087   TAbsFitParticle* bBar_;
00088   TAbsFitParticle* lightQ_;
00089   TAbsFitParticle* lightQBar_;
00090   TAbsFitParticle* lightP_;
00091   TAbsFitParticle* lightPBar_;
00093   const std::vector<edm::ParameterSet>* udscResolutions_;
00094   const std::vector<edm::ParameterSet>* bResolutions_;
00096   const std::vector<double>* jetEnergyResolutionScaleFactors_;
00097   const std::vector<double>* jetEnergyResolutionEtaBinning_;
00099   std::map<Constraint, TFitConstraintM*> massConstr_;
00101   pat::Particle fittedB_;
00102   pat::Particle fittedBBar_;
00103   pat::Particle fittedLightQ_;
00104   pat::Particle fittedLightQBar_;
00105   pat::Particle fittedLightP_;
00106   pat::Particle fittedLightPBar_;
00108   Param jetParam_;
00110   std::vector<Constraint> constraints_;
00111 
00113   CovarianceMatrix * covM_;
00114 
00115  public:
00116 
00118   struct KinFitResult {
00119     int Status;
00120     double Chi2;
00121     double Prob;
00122     pat::Particle B;
00123     pat::Particle BBar;
00124     pat::Particle LightQ;
00125     pat::Particle LightQBar;
00126     pat::Particle LightP;
00127     pat::Particle LightPBar;
00128     std::vector<int> JetCombi;
00129     bool operator< (const KinFitResult& rhs) { return Chi2 < rhs.Chi2; };
00130   };
00131 
00133   class KinFit {
00134 
00135   public:
00136 
00138     KinFit();
00140     KinFit(bool useBTagging, unsigned int bTags, std::string bTagAlgo, double minBTagValueBJet, double maxBTagValueNonBJet,
00141            std::vector<edm::ParameterSet> udscResolutions, std::vector<edm::ParameterSet> bResolutions, std::vector<double> jetEnergyResolutionScaleFactors,
00142            std::vector<double> jetEnergyResolutionEtaBinning, std::string jetCorrectionLevel, int maxNJets, int maxNComb,
00143            unsigned int maxNrIter, double maxDeltaS, double maxF, unsigned int jetParam, std::vector<unsigned> constraints, double mW, double mTop);
00145     ~KinFit();
00146     
00148     void setBTagging(bool useBTagging, unsigned int bTags, std::string bTagAlgo, double minBTagValueBJet, double maxBTagValueNonBJet){
00149       useBTagging_         = useBTagging;
00150       bTags_               = bTags;
00151       bTagAlgo_            = bTagAlgo;
00152       minBTagValueBJet_    = minBTagValueBJet;
00153       maxBTagValueNonBJet_ = maxBTagValueNonBJet;
00154     }
00156     void setResolutions(std::vector<edm::ParameterSet> udscResolutions, std::vector<edm::ParameterSet> bResolutions,
00157                         std::vector<double> jetEnergyResolutionScaleFactors, std::vector<double> jetEnergyResolutionEtaBinning){
00158       udscResolutions_       = udscResolutions;
00159       bResolutions_          = bResolutions;
00160       jetEnergyResolutionScaleFactors_ = jetEnergyResolutionScaleFactors;
00161       jetEnergyResolutionEtaBinning_   = jetEnergyResolutionEtaBinning;
00162     }
00164     void setFitter(int maxNJets, unsigned int maxNrIter, double maxDeltaS, double maxF,
00165                    unsigned int jetParam, std::vector<unsigned> constraints, double mW, double mTop){
00166       maxNJets_    = maxNJets;
00167       maxNrIter_   = maxNrIter;
00168       maxDeltaS_   = maxDeltaS;
00169       maxF_        = maxF;
00170       jetParam_    = jetParam;
00171       constraints_ = constraints;
00172       mW_          = mW;
00173       mTop_        = mTop;
00174     }
00176     void setJEC(std::string jetCorrectionLevel){
00177       jetCorrectionLevel_ = jetCorrectionLevel;
00178     }
00180     void setUseOnlyMatch(bool useOnlyMatch){
00181       useOnlyMatch_ = useOnlyMatch;
00182     }
00184     void setMatch(std::vector<int> match){
00185       match_ = match;
00186     }
00188     void setMatchInvalidity(bool invalidMatch){
00189       invalidMatch_ = invalidMatch;
00190     }
00192     void setOutput(int maxNComb){
00193       maxNComb_ = maxNComb;
00194     }
00195 
00197     std::list<TtFullHadKinFitter::KinFitResult> fit(const std::vector<pat::Jet>& jets);
00198     
00199   private:
00200 
00201     // helper function for b-tagging
00202     bool doBTagging(const std::vector<pat::Jet>& jets, const unsigned int& bJetCounter, std::vector<int>& combi);
00204     pat::Jet corJet(const pat::Jet& jet, const std::string& quarkType);
00205     
00206     // convert unsigned to Param
00207     TtFullHadKinFitter::Param param(unsigned int configParameter);
00208     // convert unsigned int to Constraint
00209     TtFullHadKinFitter::Constraint constraint(unsigned int configParameter);
00210     // convert vector of unsigned int's to vector of Contraint's
00211     std::vector<TtFullHadKinFitter::Constraint> constraints(std::vector<unsigned int>& configParameters);
00212     
00216     bool useBTagging_;
00218     unsigned int bTags_;
00220     std::string bTagAlgo_;
00222     double minBTagValueBJet_;
00224     double maxBTagValueNonBJet_;
00226     std::vector<edm::ParameterSet> udscResolutions_, bResolutions_;
00228     std::vector<double> jetEnergyResolutionScaleFactors_;
00229     std::vector<double> jetEnergyResolutionEtaBinning_;
00231     std::string jetCorrectionLevel_;
00233     int maxNJets_;
00235     int maxNComb_;
00237     unsigned int maxNrIter_;
00239     double maxDeltaS_;
00241     double maxF_;
00243     unsigned int jetParam_;
00245     std::vector<unsigned> constraints_;
00247     double mW_;
00249     double mTop_;
00251     bool useOnlyMatch_;
00253     std::vector<int> match_;
00255     bool invalidMatch_;
00256 
00258     TtFullHadKinFitter* fitter;
00259  
00260   };
00261 };
00262 
00263 #endif