CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
31  public:
34 
35  public:
39  std::vector<TtFullHadKinFitter::Constraint> intToConstraint(const std::vector<unsigned int>& constraints);
41  TtFullHadKinFitter(int jetParam, int maxNrIter, double maxDeltaS, double maxF, const std::vector<unsigned int>& constraints,
42  double mW=80.4, double mTop=173.,
43  const std::vector<edm::ParameterSet>* udscResolutions=0,
44  const std::vector<edm::ParameterSet>* bResolutions =0,
45  const std::vector<double>* jetEnergyResolutionScaleFactors=0,
46  const std::vector<double>* jetEnergyResolutionEtaBinning =0);
48  TtFullHadKinFitter(Param jetParam, int maxNrIter, double maxDeltaS, double maxF, const std::vector<Constraint>& constraints,
49  double mW=80.4, double mTop=173.,
50  const std::vector<edm::ParameterSet>* udscResolutions=0,
51  const std::vector<edm::ParameterSet>* bResolutions =0,
52  const std::vector<double>* jetEnergyResolutionScaleFactors=0,
53  const std::vector<double>* jetEnergyResolutionEtaBinning =0);
56 
58  int fit(const std::vector<pat::Jet>& jets);
60  const pat::Particle fittedB() const { return (fitter_->getStatus()==0 ? fittedB_ : pat::Particle()); };
62  const pat::Particle fittedBBar() const { return (fitter_->getStatus()==0 ? fittedBBar_ : pat::Particle()); };
64  const pat::Particle fittedLightQ() const { return (fitter_->getStatus()==0 ? fittedLightQ_ : pat::Particle()); };
68  const pat::Particle fittedLightP() const { return (fitter_->getStatus()==0 ? fittedLightP_ : pat::Particle()); };
73 
74  private:
76  void printSetup() const;
78  void setupFitter();
80  void setupJets();
82  void setupConstraints();
83 
84  private:
93  const std::vector<edm::ParameterSet>* udscResolutions_;
94  const std::vector<edm::ParameterSet>* bResolutions_;
96  const std::vector<double>* jetEnergyResolutionScaleFactors_;
97  const std::vector<double>* jetEnergyResolutionEtaBinning_;
99  std::map<Constraint, TFitConstraintM*> massConstr_;
110  std::vector<Constraint> constraints_;
111 
114 
115  public:
116 
118  struct KinFitResult {
119  int Status;
120  double Chi2;
121  double Prob;
128  std::vector<int> JetCombi;
129  bool operator< (const KinFitResult& rhs) { return Chi2 < rhs.Chi2; };
130  };
131 
133  class KinFit {
134 
135  public:
136 
138  KinFit();
140  KinFit(bool useBTagging, unsigned int bTags, std::string bTagAlgo, double minBTagValueBJet, double maxBTagValueNonBJet,
141  const std::vector<edm::ParameterSet>& udscResolutions, const std::vector<edm::ParameterSet>& bResolutions, const std::vector<double>& jetEnergyResolutionScaleFactors,
142  const std::vector<double>& jetEnergyResolutionEtaBinning, std::string jetCorrectionLevel, int maxNJets, int maxNComb,
143  unsigned int maxNrIter, double maxDeltaS, double maxF, unsigned int jetParam, const std::vector<unsigned>& constraints, double mW, double mTop);
145  ~KinFit();
146 
148  void setBTagging(bool useBTagging, unsigned int bTags, std::string bTagAlgo, double minBTagValueBJet, double maxBTagValueNonBJet){
149  useBTagging_ = useBTagging;
150  bTags_ = bTags;
151  bTagAlgo_ = bTagAlgo;
152  minBTagValueBJet_ = minBTagValueBJet;
153  maxBTagValueNonBJet_ = maxBTagValueNonBJet;
154  }
156  void setResolutions(const std::vector<edm::ParameterSet>& udscResolutions, const std::vector<edm::ParameterSet>& bResolutions,
157  const std::vector<double>& jetEnergyResolutionScaleFactors, const std::vector<double>& jetEnergyResolutionEtaBinning){
158  udscResolutions_ = udscResolutions;
159  bResolutions_ = bResolutions;
160  jetEnergyResolutionScaleFactors_ = jetEnergyResolutionScaleFactors;
161  jetEnergyResolutionEtaBinning_ = jetEnergyResolutionEtaBinning;
162  }
164  void setFitter(int maxNJets, unsigned int maxNrIter, double maxDeltaS, double maxF,
165  unsigned int jetParam, const std::vector<unsigned>& constraints, double mW, double mTop){
166  maxNJets_ = maxNJets;
167  maxNrIter_ = maxNrIter;
168  maxDeltaS_ = maxDeltaS;
169  maxF_ = maxF;
170  jetParam_ = jetParam;
172  mW_ = mW;
173  mTop_ = mTop;
174  }
176  void setJEC(std::string jetCorrectionLevel){
177  jetCorrectionLevel_ = jetCorrectionLevel;
178  }
180  void setUseOnlyMatch(bool useOnlyMatch){
181  useOnlyMatch_ = useOnlyMatch;
182  }
184  void setMatch(const std::vector<int>& match){
185  match_ = match;
186  }
188  void setMatchInvalidity(bool invalidMatch){
189  invalidMatch_ = invalidMatch;
190  }
192  void setOutput(int maxNComb){
193  maxNComb_ = maxNComb;
194  }
195 
197  std::list<TtFullHadKinFitter::KinFitResult> fit(const std::vector<pat::Jet>& jets);
198 
199  private:
200 
201  // helper function for b-tagging
202  bool doBTagging(const std::vector<pat::Jet>& jets, const unsigned int& bJetCounter, std::vector<int>& combi);
204  pat::Jet corJet(const pat::Jet& jet, const std::string& quarkType);
205 
206  // convert unsigned to Param
207  TtFullHadKinFitter::Param param(unsigned int configParameter);
208  // convert unsigned int to Constraint
209  TtFullHadKinFitter::Constraint constraint(unsigned int configParameter);
210  // convert vector of unsigned int's to vector of Contraint's
211  std::vector<TtFullHadKinFitter::Constraint> constraints(const std::vector<unsigned int>& configParameters);
212 
218  unsigned int bTags_;
226  std::vector<edm::ParameterSet> udscResolutions_, bResolutions_;
228  std::vector<double> jetEnergyResolutionScaleFactors_;
229  std::vector<double> jetEnergyResolutionEtaBinning_;
237  unsigned int maxNrIter_;
239  double maxDeltaS_;
241  double maxF_;
243  unsigned int jetParam_;
245  std::vector<unsigned> constraints_;
247  double mW_;
249  double mTop_;
253  std::vector<int> match_;
256 
259 
260  };
261 };
262 
263 #endif
TAbsFitParticle * lightQ_
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)
CovarianceMatrix * covM_
get object resolutions and put them into a matrix
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)
TAbsFitParticle * lightP_
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
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
TtFullHadKinFitter * fitter
kinematic fit interface
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
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_
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
TAbsFitParticle * b_
input particles
Int_t getStatus()
Definition: TKinFitter.h:41
double minBTagValueBJet_
min value of bTag for a b-jet
std::map< Constraint, TFitConstraintM * > massConstr_
supported constraints
std::vector< TtFullHadKinFitter::Constraint > constraints(const std::vector< unsigned int > &configParameters)
pat::Particle fittedLightP_
TAbsFitParticle * lightQBar_
const pat::Particle fittedLightQ() const
return fitted light quark candidate
pat::Particle fittedBBar_
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
TAbsFitParticle * bBar_
~KinFit()
default destructor
Analysis-level particle class.
Definition: Particle.h:32
std::string bTagAlgo_
input tag for b-tagging algorithm
const std::vector< edm::ParameterSet > * udscResolutions_
resolutions
Analysis-level calorimeter jet class.
Definition: Jet.h:73
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
const pat::Particle fittedLightP() const
return fitted light quark candidate
void setupJets()
initialize jet inputs
TKinFitter * fitter_
kinematic fitter
Definition: TopKinFitter.h:46
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
bool operator<(const KinFitResult &rhs)
double maxBTagValueNonBJet_
max value of bTag for a non-b-jet
Definition: Chi2.h:17
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
TAbsFitParticle * lightPBar_
void setMatchInvalidity(bool invalidMatch)
set the validity of a match