CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
TtHadEvtSolutionMaker Class Reference
Inheritance diagram for TtHadEvtSolutionMaker:
edm::stream::EDProducer<>

Public Member Functions

void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 TtHadEvtSolutionMaker (const edm::ParameterSet &iConfig)
 constructor More...
 
 ~TtHadEvtSolutionMaker () override
 destructor More...
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Attributes

bool addLRJetComb_
 
bool addLRSignalSel_
 
std::vector< unsigned int > constraints_
 
bool doKinFit_
 
edm::EDGetTokenT< TtGenEventgenEvtToken_
 
int jetCorrScheme_
 
int jetParam_
 
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
 
std::string lrJetCombFile_
 
std::vector< int > lrJetCombObs_
 
std::string lrSignalSelFile_
 
std::vector< int > lrSignalSelObs_
 
int matchingAlgo_
 
bool matchToGenEvt_
 
double maxDeltaS_
 
double maxDist_
 
double maxF_
 
int maxNrIter_
 
TtFullHadKinFittermyKinFitter
 
TtHadLRJetCombCalcmyLRJetCombCalc
 
TtHadLRJetCombObservablesmyLRJetCombObservables
 
TtHadLRSignalSelCalcmyLRSignalSelCalc
 
TtHadLRSignalSelObservablesmyLRSignalSelObservables
 
TtHadSimpleBestJetCombmySimpleBestJetComb
 
bool useDeltaR_
 
bool useMaxDist_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 26 of file TtHadEvtSolutionMaker.cc.

Constructor & Destructor Documentation

◆ TtHadEvtSolutionMaker()

TtHadEvtSolutionMaker::TtHadEvtSolutionMaker ( const edm::ParameterSet iConfig)
explicit

constructor

Definition at line 59 of file TtHadEvtSolutionMaker.cc.

References addLRJetComb_, addLRSignalSel_, constraints_, doKinFit_, genEvtToken_, edm::ParameterSet::getParameter(), HLT_2022v15_cff::InputTag, jetCorrScheme_, jetParam_, jetSrcToken_, lrJetCombFile_, lrJetCombObs_, lrSignalSelFile_, lrSignalSelObs_, matchingAlgo_, matchToGenEvt_, maxDeltaS_, maxDist_, maxF_, maxNrIter_, myKinFitter, myLRJetCombCalc, myLRJetCombObservables, myLRSignalSelCalc, myLRSignalSelObservables, mySimpleBestJetComb, AlCaHLTBitMon_QueryRunRegistry::string, useDeltaR_, and useMaxDist_.

59  {
60  // configurables
61  jetSrcToken_ = consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("jetSource"));
62  jetCorrScheme_ = iConfig.getParameter<int>("jetCorrectionScheme");
63  doKinFit_ = iConfig.getParameter<bool>("doKinFit");
64  addLRSignalSel_ = iConfig.getParameter<bool>("addLRSignalSel");
65  lrSignalSelObs_ = iConfig.getParameter<std::vector<int> >("lrSignalSelObs");
66  lrSignalSelFile_ = iConfig.getParameter<std::string>("lrSignalSelFile");
67  addLRJetComb_ = iConfig.getParameter<bool>("addLRJetComb");
68  lrJetCombObs_ = iConfig.getParameter<std::vector<int> >("lrJetCombObs");
69  lrJetCombFile_ = iConfig.getParameter<std::string>("lrJetCombFile");
70  maxNrIter_ = iConfig.getParameter<int>("maxNrIter");
71  maxDeltaS_ = iConfig.getParameter<double>("maxDeltaS");
72  maxF_ = iConfig.getParameter<double>("maxF");
73  jetParam_ = iConfig.getParameter<int>("jetParametrisation");
74  constraints_ = iConfig.getParameter<std::vector<unsigned int> >("constraints");
75  matchToGenEvt_ = iConfig.getParameter<bool>("matchToGenEvt");
76  matchingAlgo_ = iConfig.getParameter<bool>("matchingAlgorithm");
77  useMaxDist_ = iConfig.getParameter<bool>("useMaximalDistance");
78  useDeltaR_ = iConfig.getParameter<bool>("useDeltaR");
79  maxDist_ = iConfig.getParameter<double>("maximalDistance");
80  genEvtToken_ = mayConsume<TtGenEvent>(edm::InputTag("genEvt"));
81 
82  // define kinfitter
83  if (doKinFit_) {
85  }
86 
87  // define jet combinations related calculators
91  if (addLRJetComb_)
93 
94  // instantiate signal selection calculator
95  if (addLRSignalSel_)
97 
98  // define what will be produced
99  produces<std::vector<TtHadEvtSolution> >();
100 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
TtHadLRSignalSelObservables * myLRSignalSelObservables
TtHadSimpleBestJetComb * mySimpleBestJetComb
TtFullHadKinFitter * myKinFitter
std::vector< int > lrJetCombObs_
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
Steering class for the overall hadronic top likelihood.
TtHadLRJetCombCalc * myLRJetCombCalc
std::vector< unsigned int > constraints_
Class to calculate the jet combination LR value and purity from a root-file with fit functions...
Based on the TtSemiSimpleBestJetComb.by: Jan Heyninck version: TtSemiSimpleBestJetComb.cc,v 1.2 2007/06/09 01:17:40 lowette Exp.
edm::EDGetTokenT< TtGenEvent > genEvtToken_
TtHadLRJetCombObservables * myLRJetCombObservables
TtHadLRSignalSelCalc * myLRSignalSelCalc
std::vector< int > lrSignalSelObs_

◆ ~TtHadEvtSolutionMaker()

TtHadEvtSolutionMaker::~TtHadEvtSolutionMaker ( )
override

destructor

Definition at line 103 of file TtHadEvtSolutionMaker.cc.

References addLRJetComb_, addLRSignalSel_, doKinFit_, myKinFitter, myLRJetCombCalc, myLRJetCombObservables, myLRSignalSelCalc, myLRSignalSelObservables, and mySimpleBestJetComb.

103  {
104  if (doKinFit_) {
105  delete myKinFitter;
106  }
107  delete mySimpleBestJetComb;
109  delete myLRJetCombObservables;
110  if (addLRSignalSel_)
111  delete myLRSignalSelCalc;
112  if (addLRJetComb_)
113  delete myLRJetCombCalc;
114 }
TtHadLRSignalSelObservables * myLRSignalSelObservables
TtHadSimpleBestJetComb * mySimpleBestJetComb
TtFullHadKinFitter * myKinFitter
TtHadLRJetCombCalc * myLRJetCombCalc
TtHadLRJetCombObservables * myLRJetCombObservables
TtHadLRSignalSelCalc * myLRSignalSelCalc

Member Function Documentation

◆ produce()

void TtHadEvtSolutionMaker::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 116 of file TtHadEvtSolutionMaker.cc.

References TtFullHadKinFitter::addKinFitInfo(), addLRJetComb_, addLRSignalSel_, gather_cfg::cout, doKinFit_, TtGenEvtProducer_cfi::genEvt, genEvtToken_, EgammaValidation_cff::genp, JetPartonMatching::getDistanceForParton(), JetPartonMatching::getMatchForParton(), JetPartonMatching::getSumDistances(), mps_fire::i, iEvent, dqmiolumiharvest::j, jetCorrScheme_, jetParam_, PDWG_EXODelayedJetMET_cff::jets, jetSrcToken_, dqmdumpme::k, matchingAlgo_, matchToGenEvt_, maxDist_, eostools::move(), myKinFitter, AlCaHLTBitMon_ParallelJobs::p, submitPVResolutionJobs::q, alignCSCRings::s, TtHadEvtSolution::setJetParametrisation(), useDeltaR_, and useMaxDist_.

116  {
117  // TopObject Selection
118  // Select Jets
119 
120  bool jetsFound = false;
122  iEvent.getByToken(jetSrcToken_, jets);
123 
124  if (jets->size() >= 6)
125  jetsFound = true;
126 
127  // Build Event solutions according to the ambiguity in the jet combination
128  // Note, hardcoded to only run through the 6 most energetic jets - could be changed ....
129 
130  std::vector<TtHadEvtSolution>* evtsols = new std::vector<TtHadEvtSolution>();
131  if (jetsFound) {
132  for (unsigned int p = 0; p < 3; p++) { // loop over light jet p
133  for (unsigned int q = p + 1; q < 4; q++) { // loop over light jet q
134  for (unsigned int j = q + 1; j < 5; j++) { // loop over light jet j
135  for (unsigned int k = j + 1; k < 6; k++) { // loop over light jet k
136  for (unsigned int bh = 0; bh != jets->size(); bh++) { //loop over hadronic b-jet1
137  if (!(bh == p || bh == q || bh == j || bh == k)) {
138  for (unsigned int bbarh = 0; bbarh != jets->size(); bbarh++) { //loop over hadronic b-jet2
139  if (!(bbarh == p || bbarh == q || bbarh == j || bbarh == k) && !(bbarh == bh)) {
140  // Make event solutions for all possible combinations of the 4 light
141  // jets and 2 possible b-jets, not including the option of the b's being swapped.
142  // Hadp,Hadq is one pair, Hadj,Hadk the other
143  std::vector<TtHadEvtSolution> asol;
144  asol.resize(3);
145  //[p][q][b] and [j][k][bbar]
146  asol[0].setJetCorrectionScheme(jetCorrScheme_);
147  asol[0].setHadp(jets, p);
148  asol[0].setHadq(jets, q);
149  asol[0].setHadj(jets, j);
150  asol[0].setHadk(jets, k);
151  asol[0].setHadb(jets, bh);
152  asol[0].setHadbbar(jets, bbarh);
153 
154  //[p][j][b] and [q][k][bbar]
155  asol[1].setJetCorrectionScheme(jetCorrScheme_);
156  asol[1].setHadp(jets, p);
157  asol[1].setHadq(jets, j);
158  asol[1].setHadj(jets, q);
159  asol[1].setHadk(jets, k);
160  asol[1].setHadb(jets, bh);
161  asol[1].setHadbbar(jets, bbarh);
162 
163  //[p][k][b] and [j][q][bbar]
164  asol[2].setJetCorrectionScheme(jetCorrScheme_);
165  asol[2].setHadp(jets, p);
166  asol[2].setHadq(jets, k);
167  asol[2].setHadj(jets, j);
168  asol[2].setHadk(jets, q);
169  asol[2].setHadb(jets, bh);
170  asol[2].setHadbbar(jets, bbarh);
171 
172  if (doKinFit_) {
173  for (unsigned int i = 0; i != asol.size(); i++) {
174  asol[i] = myKinFitter->addKinFitInfo(&(asol[i]));
176  }
177 
178  } else {
179  std::cout << "Fitting needed to decide on best solution, enable fitting!" << std::endl;
180  }
181  // these lines calculate the observables to be used in the TtHadSignalSelection LR
182 
183  for (unsigned int i = 0; i != asol.size(); i++) {
184  (*myLRSignalSelObservables)(asol[i]);
185  }
186  // if asked for, calculate with these observable values the LRvalue and
187  // (depending on the configuration) probability this event is signal
188  if (addLRSignalSel_) {
189  for (unsigned int i = 0; i != asol.size(); i++) {
190  (*myLRSignalSelCalc)(asol[i]);
191  }
192  }
193 
194  // these lines calculate the observables to be used in the TtHadJetCombination LR
195  for (unsigned int i = 0; i != asol.size(); i++) {
196  (*myLRJetCombObservables)(asol[i]);
197  }
198  // if asked for, calculate with these observable values the LRvalue and
199  // (depending on the configuration) probability a jet combination is correct
200  if (addLRJetComb_) {
201  for (unsigned int i = 0; i != asol.size(); i++) {
202  (*myLRJetCombCalc)(asol[i]);
203  }
204  }
205  //std::cout<<"SignalSelLRval = "<<asol.getLRSignalEvtLRval()<<" JetCombProb = "<<asol.getLRSignalEvtProb()<<std::endl;
206  //std::cout<<"JetCombLRval = "<<asol.getLRJetCombLRval()<<" JetCombProb = "<<asol.getLRJetCombProb()<<std::endl;
207 
208  // fill solution to vector with all possible solutions
209  for (unsigned int i = 0; i != asol.size(); i++) {
210  evtsols->push_back(asol[i]);
211  }
212  }
213  }
214  }
215  }
216  }
217  }
218  }
219  }
220 
221  // add TtHadSimpleBestJetComb to solutions
222  int simpleBestJetComb = (*mySimpleBestJetComb)(*evtsols);
223 
224  for (size_t s = 0; s < evtsols->size(); s++) {
225  (*evtsols)[s].setSimpleBestJetComb(simpleBestJetComb);
226  // if asked for, match the event solutions to the gen Event
227  if (matchToGenEvt_) {
228  int bestSolution = -999;
229  int bestSolutionChangeW1Q = -999;
230  int bestSolutionChangeW2Q = -999;
232  iEvent.getByToken(genEvtToken_, genEvt);
233  std::vector<const reco::Candidate*> quarks;
234  const reco::Candidate& genp = *(genEvt->daughterQuarkOfWPlus());
235  const reco::Candidate& genq = *(genEvt->daughterQuarkBarOfWPlus());
236  const reco::Candidate& genb = *(genEvt->b());
237  const reco::Candidate& genj = *(genEvt->daughterQuarkOfWMinus());
238  const reco::Candidate& genk = *(genEvt->daughterQuarkBarOfWMinus());
239  const reco::Candidate& genbbar = *(genEvt->bBar());
240  quarks.push_back(&genp);
241  quarks.push_back(&genq);
242  quarks.push_back(&genb);
243  quarks.push_back(&genj);
244  quarks.push_back(&genk);
245  quarks.push_back(&genbbar);
246  std::vector<const reco::Candidate*> jets;
247  for (size_t s = 0; s < evtsols->size(); s++) {
248  jets.clear();
249  const reco::Candidate& jetp = (*evtsols)[s].getRecHadp();
250  const reco::Candidate& jetq = (*evtsols)[s].getRecHadq();
251  const reco::Candidate& jetbh = (*evtsols)[s].getRecHadb();
252  const reco::Candidate& jetj = (*evtsols)[s].getRecHadj();
253  const reco::Candidate& jetk = (*evtsols)[s].getRecHadk();
254  const reco::Candidate& jetbbar = (*evtsols)[s].getRecHadbbar();
255  jets.push_back(&jetp);
256  jets.push_back(&jetq);
257  jets.push_back(&jetbh);
258  jets.push_back(&jetj);
259  jets.push_back(&jetk);
260  jets.push_back(&jetbbar);
262  (*evtsols)[s].setGenEvt(genEvt);
263  (*evtsols)[s].setMCBestSumAngles(aMatch.getSumDistances());
264  (*evtsols)[s].setMCBestAngleHadp(aMatch.getDistanceForParton(0));
265  (*evtsols)[s].setMCBestAngleHadq(aMatch.getDistanceForParton(1));
266  (*evtsols)[s].setMCBestAngleHadb(aMatch.getDistanceForParton(2));
267  (*evtsols)[s].setMCBestAngleHadb(aMatch.getDistanceForParton(2));
268  (*evtsols)[s].setMCBestAngleHadj(aMatch.getDistanceForParton(3));
269  (*evtsols)[s].setMCBestAngleHadk(aMatch.getDistanceForParton(4));
270  (*evtsols)[s].setMCBestAngleHadbbar(aMatch.getDistanceForParton(5));
271 
272  // Check match - checking if two light quarks are swapped wrt matched gen particle
273  if ((aMatch.getMatchForParton(2) == 2 && aMatch.getMatchForParton(5) == 5) ||
274  (aMatch.getMatchForParton(2) == 5 && aMatch.getMatchForParton(5) == 2)) { // check b-jets
275 
276  if (aMatch.getMatchForParton(3) == 3 && aMatch.getMatchForParton(4) == 4) { //check light jets
277  bestSolutionChangeW2Q = 0;
278  if (aMatch.getMatchForParton(0) == 0 && aMatch.getMatchForParton(1) == 1) {
279  bestSolution = s;
280  bestSolutionChangeW1Q = 0;
281  } else {
282  if (aMatch.getMatchForParton(0) == 1 && aMatch.getMatchForParton(1) == 0) {
283  bestSolution = s;
284  bestSolutionChangeW1Q = 1;
285  }
286  }
287  } else {
288  if (aMatch.getMatchForParton(2) == 3 && aMatch.getMatchForParton(3) == 2) { // or check if swapped
289  bestSolutionChangeW2Q = 1;
290  if (aMatch.getMatchForParton(0) == 1 && aMatch.getMatchForParton(1) == 0) {
291  bestSolution = s;
292  bestSolutionChangeW1Q = 1;
293  } else {
294  if (aMatch.getMatchForParton(0) == 0 && aMatch.getMatchForParton(1) == 1) {
295  bestSolution = s;
296  bestSolutionChangeW1Q = 0;
297  }
298  }
299  }
300  if (aMatch.getMatchForParton(2) == 2 && aMatch.getMatchForParton(3) == 3) {
301  bestSolutionChangeW2Q = 0;
302  if (aMatch.getMatchForParton(0) == 0 && aMatch.getMatchForParton(1) == 1) {
303  bestSolution = s;
304  bestSolutionChangeW1Q = 0;
305  } else if (aMatch.getMatchForParton(0) == 1 && aMatch.getMatchForParton(1) == 0) {
306  bestSolution = s;
307  bestSolutionChangeW1Q = 1;
308  }
309  }
310  }
311  }
312  for (size_t s = 0; s < evtsols->size(); s++) {
313  (*evtsols)[s].setMCBestJetComb(bestSolution);
314  (*evtsols)[s].setMCChangeW1Q(bestSolutionChangeW1Q);
315  (*evtsols)[s].setMCChangeW2Q(bestSolutionChangeW2Q);
316  }
317  }
318  } // end matchEvt
319  }
320  //store the vector of solutions to the event
321 
322  std::unique_ptr<std::vector<TtHadEvtSolution> > pOut(evtsols);
323  iEvent.put(std::move(pOut));
324  } else { //end loop jet/MET found
325  std::cout << "No calibrated solutions built, because only " << jets->size() << " were present";
326 
327  std::unique_ptr<std::vector<TtHadEvtSolution> > pOut(evtsols);
328  iEvent.put(std::move(pOut));
329  }
330 }
genp
produce generated paricles in acceptance #
TtFullHadKinFitter * myKinFitter
TtHadEvtSolution addKinFitInfo(TtHadEvtSolution *asol)
add kin fit information to the old event solution (in for legacy reasons)
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
edm::EDGetTokenT< TtGenEvent > genEvtToken_
int iEvent
Definition: GenABIO.cc:224
void setJetParametrisation(int jp)
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ addLRJetComb_

bool TtHadEvtSolutionMaker::addLRJetComb_
private

Definition at line 39 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ addLRSignalSel_

bool TtHadEvtSolutionMaker::addLRSignalSel_
private

Definition at line 39 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ constraints_

std::vector<unsigned int> TtHadEvtSolutionMaker::constraints_
private

Definition at line 47 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ doKinFit_

bool TtHadEvtSolutionMaker::doKinFit_
private

Definition at line 39 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ genEvtToken_

edm::EDGetTokenT<TtGenEvent> TtHadEvtSolutionMaker::genEvtToken_
private

Definition at line 48 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ jetCorrScheme_

int TtHadEvtSolutionMaker::jetCorrScheme_
private

Definition at line 37 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ jetParam_

int TtHadEvtSolutionMaker::jetParam_
private

Definition at line 45 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ jetSrcToken_

edm::EDGetTokenT<std::vector<pat::Jet> > TtHadEvtSolutionMaker::jetSrcToken_
private

Definition at line 36 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ lrJetCombFile_

std::string TtHadEvtSolutionMaker::lrJetCombFile_
private

Definition at line 38 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ lrJetCombObs_

std::vector<int> TtHadEvtSolutionMaker::lrJetCombObs_
private

Definition at line 46 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ lrSignalSelFile_

std::string TtHadEvtSolutionMaker::lrSignalSelFile_
private

Definition at line 38 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ lrSignalSelObs_

std::vector<int> TtHadEvtSolutionMaker::lrSignalSelObs_
private

Definition at line 46 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ matchingAlgo_

int TtHadEvtSolutionMaker::matchingAlgo_
private

Definition at line 40 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ matchToGenEvt_

bool TtHadEvtSolutionMaker::matchToGenEvt_
private

Definition at line 39 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ maxDeltaS_

double TtHadEvtSolutionMaker::maxDeltaS_
private

Definition at line 44 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ maxDist_

double TtHadEvtSolutionMaker::maxDist_
private

Definition at line 42 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ maxF_

double TtHadEvtSolutionMaker::maxF_
private

Definition at line 44 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ maxNrIter_

int TtHadEvtSolutionMaker::maxNrIter_
private

Definition at line 43 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker().

◆ myKinFitter

TtFullHadKinFitter* TtHadEvtSolutionMaker::myKinFitter
private

Definition at line 50 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ myLRJetCombCalc

TtHadLRJetCombCalc* TtHadEvtSolutionMaker::myLRJetCombCalc
private

Definition at line 53 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ myLRJetCombObservables

TtHadLRJetCombObservables* TtHadEvtSolutionMaker::myLRJetCombObservables
private

Definition at line 52 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ myLRSignalSelCalc

TtHadLRSignalSelCalc* TtHadEvtSolutionMaker::myLRSignalSelCalc
private

Definition at line 55 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ myLRSignalSelObservables

TtHadLRSignalSelObservables* TtHadEvtSolutionMaker::myLRSignalSelObservables
private

Definition at line 54 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ mySimpleBestJetComb

TtHadSimpleBestJetComb* TtHadEvtSolutionMaker::mySimpleBestJetComb
private

Definition at line 51 of file TtHadEvtSolutionMaker.cc.

Referenced by TtHadEvtSolutionMaker(), and ~TtHadEvtSolutionMaker().

◆ useDeltaR_

bool TtHadEvtSolutionMaker::useDeltaR_
private

Definition at line 41 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().

◆ useMaxDist_

bool TtHadEvtSolutionMaker::useMaxDist_
private

Definition at line 41 of file TtHadEvtSolutionMaker.cc.

Referenced by produce(), and TtHadEvtSolutionMaker().