CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TtFullHadKinFitter.h
Go to the documentation of this file.
1 #ifndef TtFullHadKinFitter_h
2 #define TtFullHadKinFitter_h
3 
4 #include <vector>
5 
6 #include "TLorentzVector.h"
7 
9 
12 
14 
16 
17 class TAbsFitParticle;
18 class TFitConstraintM;
19 
20 /*
21  \class TtFullHadKinFitter TtFullHadKinFitter.h "TopQuarkAnalysis/TopKinFitter/interface/TtFullHadKinFitter.h"
22 
23  \brief one line description to be added here...
24 
25  text to be added here...
26 
27 **/
28 
30 public:
33 
34 public:
38  std::vector<TtFullHadKinFitter::Constraint> intToConstraint(const std::vector<unsigned int>& constraints);
40  TtFullHadKinFitter(int jetParam,
41  int maxNrIter,
42  double maxDeltaS,
43  double maxF,
44  const std::vector<unsigned int>& constraints,
45  double mW = 80.4,
46  double mTop = 173.,
47  const std::vector<edm::ParameterSet>* udscResolutions = nullptr,
48  const std::vector<edm::ParameterSet>* bResolutions = nullptr,
49  const std::vector<double>* jetEnergyResolutionScaleFactors = nullptr,
50  const std::vector<double>* jetEnergyResolutionEtaBinning = nullptr);
52  TtFullHadKinFitter(Param jetParam,
53  int maxNrIter,
54  double maxDeltaS,
55  double maxF,
56  const std::vector<Constraint>& constraints,
57  double mW = 80.4,
58  double mTop = 173.,
59  const std::vector<edm::ParameterSet>* udscResolutions = nullptr,
60  const std::vector<edm::ParameterSet>* bResolutions = nullptr,
61  const std::vector<double>* jetEnergyResolutionScaleFactors = nullptr,
62  const std::vector<double>* jetEnergyResolutionEtaBinning = nullptr);
65 
67  int fit(const std::vector<pat::Jet>& jets);
69  const pat::Particle fittedB() const { return (fitter_->getStatus() == 0 ? fittedB_ : pat::Particle()); };
71  const pat::Particle fittedBBar() const { return (fitter_->getStatus() == 0 ? fittedBBar_ : pat::Particle()); };
73  const pat::Particle fittedLightQ() const { return (fitter_->getStatus() == 0 ? fittedLightQ_ : pat::Particle()); };
76  return (fitter_->getStatus() == 0 ? fittedLightQBar_ : pat::Particle());
77  };
79  const pat::Particle fittedLightP() const { return (fitter_->getStatus() == 0 ? fittedLightP_ : pat::Particle()); };
82  return (fitter_->getStatus() == 0 ? fittedLightPBar_ : pat::Particle());
83  };
86 
87 private:
89  void printSetup() const;
91  void setupFitter();
93  void setupJets();
95  void setupConstraints();
96 
97 private:
99  std::unique_ptr<TAbsFitParticle> b_;
100  std::unique_ptr<TAbsFitParticle> bBar_;
101  std::unique_ptr<TAbsFitParticle> lightQ_;
102  std::unique_ptr<TAbsFitParticle> lightQBar_;
103  std::unique_ptr<TAbsFitParticle> lightP_;
104  std::unique_ptr<TAbsFitParticle> lightPBar_;
106  const std::vector<edm::ParameterSet>* udscResolutions_ = nullptr;
107  const std::vector<edm::ParameterSet>* bResolutions_ = nullptr;
109  const std::vector<double>* jetEnergyResolutionScaleFactors_ = nullptr;
110  const std::vector<double>* jetEnergyResolutionEtaBinning_ = nullptr;
112  std::map<Constraint, std::unique_ptr<TFitConstraintM>> massConstr_;
123  std::vector<Constraint> constraints_;
124 
126  std::unique_ptr<CovarianceMatrix> covM_;
127 
128 public:
130  struct KinFitResult {
131  int Status;
132  double Chi2;
133  double Prob;
140  std::vector<int> JetCombi;
141  bool operator<(const KinFitResult& rhs) { return Chi2 < rhs.Chi2; };
142  };
143 
145  class KinFit {
146  public:
148  KinFit();
150  KinFit(bool useBTagging,
151  unsigned int bTags,
152  std::string bTagAlgo,
153  double minBTagValueBJet,
154  double maxBTagValueNonBJet,
155  const std::vector<edm::ParameterSet>& udscResolutions,
156  const std::vector<edm::ParameterSet>& bResolutions,
157  const std::vector<double>& jetEnergyResolutionScaleFactors,
158  const std::vector<double>& jetEnergyResolutionEtaBinning,
159  std::string jetCorrectionLevel,
160  int maxNJets,
161  int maxNComb,
162  unsigned int maxNrIter,
163  double maxDeltaS,
164  double maxF,
165  unsigned int jetParam,
166  const std::vector<unsigned>& constraints,
167  double mW,
168  double mTop);
170  ~KinFit();
171 
173  void setBTagging(bool useBTagging,
174  unsigned int bTags,
175  std::string bTagAlgo,
176  double minBTagValueBJet,
177  double maxBTagValueNonBJet) {
178  useBTagging_ = useBTagging;
179  bTags_ = bTags;
180  bTagAlgo_ = bTagAlgo;
181  minBTagValueBJet_ = minBTagValueBJet;
182  maxBTagValueNonBJet_ = maxBTagValueNonBJet;
183  }
185  void setResolutions(const std::vector<edm::ParameterSet>& udscResolutions,
186  const std::vector<edm::ParameterSet>& bResolutions,
187  const std::vector<double>& jetEnergyResolutionScaleFactors,
188  const std::vector<double>& jetEnergyResolutionEtaBinning) {
189  udscResolutions_ = udscResolutions;
190  bResolutions_ = bResolutions;
191  jetEnergyResolutionScaleFactors_ = jetEnergyResolutionScaleFactors;
192  jetEnergyResolutionEtaBinning_ = jetEnergyResolutionEtaBinning;
193  }
195  void setFitter(int maxNJets,
196  unsigned int maxNrIter,
197  double maxDeltaS,
198  double maxF,
199  unsigned int jetParam,
200  const std::vector<unsigned>& constraints,
201  double mW,
202  double mTop) {
203  maxNJets_ = maxNJets;
204  maxNrIter_ = maxNrIter;
205  maxDeltaS_ = maxDeltaS;
206  maxF_ = maxF;
207  jetParam_ = jetParam;
209  mW_ = mW;
210  mTop_ = mTop;
211  }
213  void setJEC(std::string jetCorrectionLevel) { jetCorrectionLevel_ = jetCorrectionLevel; }
215  void setUseOnlyMatch(bool useOnlyMatch) { useOnlyMatch_ = useOnlyMatch; }
217  void setMatch(const std::vector<int>& match) { match_ = match; }
219  void setMatchInvalidity(bool invalidMatch) { invalidMatch_ = invalidMatch; }
221  void setOutput(int maxNComb) { maxNComb_ = maxNComb; }
222 
224  std::list<TtFullHadKinFitter::KinFitResult> fit(const std::vector<pat::Jet>& jets);
225 
226  private:
227  // helper function for b-tagging
228  bool doBTagging(const std::vector<pat::Jet>& jets, const unsigned int& bJetCounter, std::vector<int>& combi);
230  pat::Jet corJet(const pat::Jet& jet, const std::string& quarkType);
231 
232  // convert unsigned to Param
233  TtFullHadKinFitter::Param param(unsigned int configParameter);
234  // convert unsigned int to Constraint
235  TtFullHadKinFitter::Constraint constraint(unsigned int configParameter);
236  // convert vector of unsigned int's to vector of Contraint's
237  std::vector<TtFullHadKinFitter::Constraint> constraints(const std::vector<unsigned int>& configParameters);
238 
244  unsigned int bTags_;
252  std::vector<edm::ParameterSet> udscResolutions_, bResolutions_;
254  std::vector<double> jetEnergyResolutionScaleFactors_;
255  std::vector<double> jetEnergyResolutionEtaBinning_;
263  unsigned int maxNrIter_;
265  double maxDeltaS_;
267  double maxF_;
269  unsigned int jetParam_;
271  std::vector<unsigned> constraints_;
273  double mW_;
275  double mTop_;
279  std::vector<int> match_;
282 
284  std::unique_ptr<TtFullHadKinFitter> fitter;
285  };
286 };
287 
288 #endif
void setResolutions(const std::vector< edm::ParameterSet > &udscResolutions, const std::vector< edm::ParameterSet > &bResolutions, const std::vector< double > &jetEnergyResolutionScaleFactors, const std::vector< double > &jetEnergyResolutionEtaBinning)
set resolutions
double maxF_
maximal deviation for contstraints
const std::vector< double > * jetEnergyResolutionEtaBinning_
double maxDeltaS_
maximal chi2 equivalent
TtFullHadKinFitter::Param param(unsigned int configParameter)
Param
supported parameterizations
Definition: TopKinFitter.h:22
void setupConstraints()
initialize constraints
void setFitter(int maxNJets, unsigned int maxNrIter, double maxDeltaS, double maxF, unsigned int jetParam, const std::vector< unsigned > &constraints, double mW, double mTop)
set parameters for fitter
TtHadEvtSolution addKinFitInfo(TtHadEvtSolution *asol)
add kin fit information to the old event solution (in for legacy reasons)
std::unique_ptr< TAbsFitParticle > lightPBar_
void setBTagging(bool useBTagging, unsigned int bTags, std::string bTagAlgo, double minBTagValueBJet, double maxBTagValueNonBJet)
set all parameters for b-tagging
int maxNComb_
maximal number of combinations to be written to the event
TtFullHadKinFitter::Constraint constraint(unsigned int configParameter)
const pat::Particle fittedLightPBar() const
return fitted light quark candidate
std::unique_ptr< TKinFitter > fitter_
kinematic fitter
Definition: TopKinFitter.h:49
pat::Particle fittedLightPBar_
void setUseOnlyMatch(bool useOnlyMatch)
set useOnlyMatch
pat::Particle fittedB_
output particles
int fit(const std::vector< pat::Jet > &jets)
kinematic fit interface
void printSetup() const
print fitter setup
double mW_
W mass value used for constraints.
std::vector< unsigned > constraints_
numbering of different possible kinematic constraints
pat::Particle fittedLightQBar_
const std::vector< double > * jetEnergyResolutionScaleFactors_
scale factors for the jet energy resolution
std::list< TtFullHadKinFitter::KinFitResult > fit(const std::vector< pat::Jet > &jets)
do the fitting and return fit result
std::vector< double > jetEnergyResolutionScaleFactors_
scale factors for the jet energy resolution
const pat::Particle fittedBBar() const
return fitted b quark candidate
std::vector< edm::ParameterSet > udscResolutions_
store the resolutions for the jets
class that does the fitting
std::vector< edm::ParameterSet > bResolutions_
std::string jetCorrectionLevel_
correction level for jets
std::map< Constraint, std::unique_ptr< TFitConstraintM > > massConstr_
supported constraints
void setJEC(std::string jetCorrectionLevel)
set jec level
bool doBTagging(const std::vector< pat::Jet > &jets, const unsigned int &bJetCounter, std::vector< int > &combi)
pat::Jet corJet(const pat::Jet &jet, const std::string &quarkType)
helper function to construct the proper corrected jet for its corresponding quarkType ...
std::vector< double > jetEnergyResolutionEtaBinning_
std::unique_ptr< TAbsFitParticle > lightP_
unsigned int bTags_
minimal number of b-jets
Constraint
supported constraints
void setMatch(const std::vector< int > &match)
set match to be used
unsigned int maxNrIter_
maximal number of iterations to be performed for the fit
std::vector< int > match_
the combination that should be used
std::vector< TtFullHadKinFitter::Constraint > intToConstraint(const std::vector< unsigned int > &constraints)
used to convert vector of int&#39;s to vector of constraints (just used in TtFullHadKinFitter(int, int, double, double, std::vector&lt;unsigned int&gt;))
~TtFullHadKinFitter()
default destructor
const std::vector< edm::ParameterSet > * bResolutions_
vector< PseudoJet > jets
double minBTagValueBJet_
min value of bTag for a b-jet
std::vector< TtFullHadKinFitter::Constraint > constraints(const std::vector< unsigned int > &configParameters)
pat::Particle fittedLightP_
const pat::Particle fittedLightQ() const
return fitted light quark candidate
pat::Particle fittedBBar_
std::unique_ptr< TAbsFitParticle > b_
input particles
double mTop_
top mass value used for constraints
std::vector< Constraint > constraints_
vector of constraints to be used
Param jetParam_
jet parametrization
pat::Particle fittedLightQ_
void setupFitter()
setup fitter
~KinFit()
default destructor
Analysis-level particle class.
Definition: Particle.h:30
std::string bTagAlgo_
input tag for b-tagging algorithm
const std::vector< edm::ParameterSet > * udscResolutions_
resolutions
Analysis-level calorimeter jet class.
Definition: Jet.h:77
void setOutput(int maxNComb)
set number of combinations of output
unsigned int jetParam_
numbering of different possible jet parametrizations
bool useOnlyMatch_
fit or only a certain combination
int maxNJets_
maximal number of jets (-1 possible to indicate &#39;all&#39;)
TtFullHadKinFitter()
default constructor
std::unique_ptr< CovarianceMatrix > covM_
get object resolutions and put them into a matrix
const pat::Particle fittedLightP() const
return fitted light quark candidate
void setupJets()
initialize jet inputs
std::unique_ptr< TAbsFitParticle > lightQBar_
std::unique_ptr< TtFullHadKinFitter > fitter
kinematic fit interface
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::unique_ptr< TAbsFitParticle > lightQ_
bool operator<(const KinFitResult &rhs)
double maxBTagValueNonBJet_
max value of bTag for a non-b-jet
Definition: Chi2.h:15
bool invalidMatch_
match is invalid
const pat::Particle fittedLightQBar() const
return fitted light quark candidate
const pat::Particle fittedB() const
return fitted b quark candidate
KinFit()
default constructor
std::unique_ptr< TAbsFitParticle > bBar_
void setMatchInvalidity(bool invalidMatch)
set the validity of a match