File indexing completed on 2024-04-06 12:31:22
0001 #ifndef TtFullHadKinFitter_h
0002 #define TtFullHadKinFitter_h
0003
0004 #include <vector>
0005
0006 #include "TLorentzVector.h"
0007
0008 #include "AnalysisDataFormats/TopObjects/interface/TtHadEvtSolution.h"
0009
0010 #include "TopQuarkAnalysis/TopKinFitter/interface/CovarianceMatrix.h"
0011 #include "TopQuarkAnalysis/TopKinFitter/interface/TopKinFitter.h"
0012
0013 #include "PhysicsTools/JetMCUtils/interface/combination.h"
0014
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016
0017 class TAbsFitParticle;
0018 class TFitConstraintM;
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 class TtFullHadKinFitter : public TopKinFitter {
0030 public:
0031
0032 enum Constraint { kWPlusMass = 1, kWMinusMass, kTopMass, kTopBarMass, kEqualTopMasses };
0033
0034 public:
0035
0036 TtFullHadKinFitter();
0037
0038 std::vector<TtFullHadKinFitter::Constraint> intToConstraint(const std::vector<unsigned int>& constraints);
0039
0040 TtFullHadKinFitter(int jetParam,
0041 int maxNrIter,
0042 double maxDeltaS,
0043 double maxF,
0044 const std::vector<unsigned int>& constraints,
0045 double mW = 80.4,
0046 double mTop = 173.,
0047 const std::vector<edm::ParameterSet>* udscResolutions = nullptr,
0048 const std::vector<edm::ParameterSet>* bResolutions = nullptr,
0049 const std::vector<double>* jetEnergyResolutionScaleFactors = nullptr,
0050 const std::vector<double>* jetEnergyResolutionEtaBinning = nullptr);
0051
0052 TtFullHadKinFitter(Param jetParam,
0053 int maxNrIter,
0054 double maxDeltaS,
0055 double maxF,
0056 const std::vector<Constraint>& constraints,
0057 double mW = 80.4,
0058 double mTop = 173.,
0059 const std::vector<edm::ParameterSet>* udscResolutions = nullptr,
0060 const std::vector<edm::ParameterSet>* bResolutions = nullptr,
0061 const std::vector<double>* jetEnergyResolutionScaleFactors = nullptr,
0062 const std::vector<double>* jetEnergyResolutionEtaBinning = nullptr);
0063
0064 ~TtFullHadKinFitter();
0065
0066
0067 int fit(const std::vector<pat::Jet>& jets);
0068
0069 const pat::Particle fittedB() const { return (fitter_->getStatus() == 0 ? fittedB_ : pat::Particle()); };
0070
0071 const pat::Particle fittedBBar() const { return (fitter_->getStatus() == 0 ? fittedBBar_ : pat::Particle()); };
0072
0073 const pat::Particle fittedLightQ() const { return (fitter_->getStatus() == 0 ? fittedLightQ_ : pat::Particle()); };
0074
0075 const pat::Particle fittedLightQBar() const {
0076 return (fitter_->getStatus() == 0 ? fittedLightQBar_ : pat::Particle());
0077 };
0078
0079 const pat::Particle fittedLightP() const { return (fitter_->getStatus() == 0 ? fittedLightP_ : pat::Particle()); };
0080
0081 const pat::Particle fittedLightPBar() const {
0082 return (fitter_->getStatus() == 0 ? fittedLightPBar_ : pat::Particle());
0083 };
0084
0085 TtHadEvtSolution addKinFitInfo(TtHadEvtSolution* asol);
0086
0087 private:
0088
0089 void printSetup() const;
0090
0091 void setupFitter();
0092
0093 void setupJets();
0094
0095 void setupConstraints();
0096
0097 private:
0098
0099 std::unique_ptr<TAbsFitParticle> b_;
0100 std::unique_ptr<TAbsFitParticle> bBar_;
0101 std::unique_ptr<TAbsFitParticle> lightQ_;
0102 std::unique_ptr<TAbsFitParticle> lightQBar_;
0103 std::unique_ptr<TAbsFitParticle> lightP_;
0104 std::unique_ptr<TAbsFitParticle> lightPBar_;
0105
0106 const std::vector<edm::ParameterSet>* udscResolutions_ = nullptr;
0107 const std::vector<edm::ParameterSet>* bResolutions_ = nullptr;
0108
0109 const std::vector<double>* jetEnergyResolutionScaleFactors_ = nullptr;
0110 const std::vector<double>* jetEnergyResolutionEtaBinning_ = nullptr;
0111
0112 std::map<Constraint, std::unique_ptr<TFitConstraintM>> massConstr_;
0113
0114 pat::Particle fittedB_;
0115 pat::Particle fittedBBar_;
0116 pat::Particle fittedLightQ_;
0117 pat::Particle fittedLightQBar_;
0118 pat::Particle fittedLightP_;
0119 pat::Particle fittedLightPBar_;
0120
0121 Param jetParam_;
0122
0123 std::vector<Constraint> constraints_;
0124
0125
0126 std::unique_ptr<CovarianceMatrix> covM_;
0127
0128 public:
0129
0130 struct KinFitResult {
0131 int Status;
0132 double Chi2;
0133 double Prob;
0134 pat::Particle B;
0135 pat::Particle BBar;
0136 pat::Particle LightQ;
0137 pat::Particle LightQBar;
0138 pat::Particle LightP;
0139 pat::Particle LightPBar;
0140 std::vector<int> JetCombi;
0141 bool operator<(const KinFitResult& rhs) { return Chi2 < rhs.Chi2; };
0142 };
0143
0144
0145 class KinFit {
0146 public:
0147
0148 KinFit();
0149
0150 KinFit(bool useBTagging,
0151 unsigned int bTags,
0152 std::string bTagAlgo,
0153 double minBTagValueBJet,
0154 double maxBTagValueNonBJet,
0155 const std::vector<edm::ParameterSet>& udscResolutions,
0156 const std::vector<edm::ParameterSet>& bResolutions,
0157 const std::vector<double>& jetEnergyResolutionScaleFactors,
0158 const std::vector<double>& jetEnergyResolutionEtaBinning,
0159 std::string jetCorrectionLevel,
0160 int maxNJets,
0161 int maxNComb,
0162 unsigned int maxNrIter,
0163 double maxDeltaS,
0164 double maxF,
0165 unsigned int jetParam,
0166 const std::vector<unsigned>& constraints,
0167 double mW,
0168 double mTop);
0169
0170 ~KinFit();
0171
0172
0173 void setBTagging(bool useBTagging,
0174 unsigned int bTags,
0175 std::string bTagAlgo,
0176 double minBTagValueBJet,
0177 double maxBTagValueNonBJet) {
0178 useBTagging_ = useBTagging;
0179 bTags_ = bTags;
0180 bTagAlgo_ = bTagAlgo;
0181 minBTagValueBJet_ = minBTagValueBJet;
0182 maxBTagValueNonBJet_ = maxBTagValueNonBJet;
0183 }
0184
0185 void setResolutions(const std::vector<edm::ParameterSet>& udscResolutions,
0186 const std::vector<edm::ParameterSet>& bResolutions,
0187 const std::vector<double>& jetEnergyResolutionScaleFactors,
0188 const std::vector<double>& jetEnergyResolutionEtaBinning) {
0189 udscResolutions_ = udscResolutions;
0190 bResolutions_ = bResolutions;
0191 jetEnergyResolutionScaleFactors_ = jetEnergyResolutionScaleFactors;
0192 jetEnergyResolutionEtaBinning_ = jetEnergyResolutionEtaBinning;
0193 }
0194
0195 void setFitter(int maxNJets,
0196 unsigned int maxNrIter,
0197 double maxDeltaS,
0198 double maxF,
0199 unsigned int jetParam,
0200 const std::vector<unsigned>& constraints,
0201 double mW,
0202 double mTop) {
0203 maxNJets_ = maxNJets;
0204 maxNrIter_ = maxNrIter;
0205 maxDeltaS_ = maxDeltaS;
0206 maxF_ = maxF;
0207 jetParam_ = jetParam;
0208 constraints_ = constraints;
0209 mW_ = mW;
0210 mTop_ = mTop;
0211 }
0212
0213 void setJEC(std::string jetCorrectionLevel) { jetCorrectionLevel_ = jetCorrectionLevel; }
0214
0215 void setUseOnlyMatch(bool useOnlyMatch) { useOnlyMatch_ = useOnlyMatch; }
0216
0217 void setMatch(const std::vector<int>& match) { match_ = match; }
0218
0219 void setMatchInvalidity(bool invalidMatch) { invalidMatch_ = invalidMatch; }
0220
0221 void setOutput(int maxNComb) { maxNComb_ = maxNComb; }
0222
0223
0224 std::list<TtFullHadKinFitter::KinFitResult> fit(const std::vector<pat::Jet>& jets);
0225
0226 private:
0227
0228 bool doBTagging(const std::vector<pat::Jet>& jets, const unsigned int& bJetCounter, std::vector<int>& combi);
0229
0230 pat::Jet corJet(const pat::Jet& jet, const std::string& quarkType);
0231
0232
0233 TtFullHadKinFitter::Param param(unsigned int configParameter);
0234
0235 TtFullHadKinFitter::Constraint constraint(unsigned int configParameter);
0236
0237 std::vector<TtFullHadKinFitter::Constraint> constraints(const std::vector<unsigned int>& configParameters);
0238
0239
0240
0241
0242 bool useBTagging_;
0243
0244 unsigned int bTags_;
0245
0246 std::string bTagAlgo_;
0247
0248 double minBTagValueBJet_;
0249
0250 double maxBTagValueNonBJet_;
0251
0252 std::vector<edm::ParameterSet> udscResolutions_, bResolutions_;
0253
0254 std::vector<double> jetEnergyResolutionScaleFactors_;
0255 std::vector<double> jetEnergyResolutionEtaBinning_;
0256
0257 std::string jetCorrectionLevel_;
0258
0259 int maxNJets_;
0260
0261 int maxNComb_;
0262
0263 unsigned int maxNrIter_;
0264
0265 double maxDeltaS_;
0266
0267 double maxF_;
0268
0269 unsigned int jetParam_;
0270
0271 std::vector<unsigned> constraints_;
0272
0273 double mW_;
0274
0275 double mTop_;
0276
0277 bool useOnlyMatch_;
0278
0279 std::vector<int> match_;
0280
0281 bool invalidMatch_;
0282
0283
0284 std::unique_ptr<TtFullHadKinFitter> fitter;
0285 };
0286 };
0287
0288 #endif