CMS 3D CMS Logo

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

Public Member Functions

TtSemiLepKinFitter::Constraint constraint (unsigned)
 
std::vector< TtSemiLepKinFitter::Constraintconstraints (std::vector< unsigned > &)
 
TtSemiLepKinFitter::Param param (unsigned)
 
void produce (edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 TtSemiEvtSolutionMaker (const edm::ParameterSet &iConfig)
 constructor More...
 
 ~TtSemiEvtSolutionMaker () 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 > constraints_
 
bool doKinFit_
 
edm::EDGetTokenT< std::vector< pat::Electron > > electronSrcToken_
 
edm::EDGetTokenT< TtGenEventgenEvtToken_
 
int jetCorrScheme_
 
int jetParam_
 
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
 
int lepParam_
 
std::string leptonFlavour_
 
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_
 
int metParam_
 
edm::EDGetTokenT< std::vector< pat::MET > > metSrcToken_
 
edm::EDGetTokenT< std::vector< pat::Muon > > muonSrcToken_
 
TtSemiLepKinFittermyKinFitter
 
TtSemiLRJetCombCalcmyLRJetCombCalc
 
TtSemiLRJetCombObservablesmyLRJetCombObservables
 
TtSemiLRSignalSelCalcmyLRSignalSelCalc
 
TtSemiLRSignalSelObservablesmyLRSignalSelObservables
 
TtSemiSimpleBestJetCombmySimpleBestJetComb
 
unsigned int nrCombJets_
 
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 22 of file TtSemiEvtSolutionMaker.cc.

Constructor & Destructor Documentation

◆ TtSemiEvtSolutionMaker()

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

constructor

Definition at line 66 of file TtSemiEvtSolutionMaker.cc.

References addLRJetComb_, addLRSignalSel_, constraints(), constraints_, doKinFit_, electronSrcToken_, contentValuesFiles::fullPath, genEvtToken_, edm::ParameterSet::getParameter(), ProducerED_cfi::InputTag, jetCorrScheme_, jetParam_, jetSrcToken_, lepParam_, leptonFlavour_, lrJetCombFile_, lrJetCombObs_, lrSignalSelFile_, lrSignalSelObs_, matchingAlgo_, matchToGenEvt_, maxDeltaS_, maxDist_, maxF_, maxNrIter_, metParam_, metSrcToken_, muonSrcToken_, myKinFitter, myLRJetCombCalc, myLRJetCombObservables, myLRSignalSelCalc, myLRSignalSelObservables, mySimpleBestJetComb, nrCombJets_, param(), AlCaHLTBitMon_QueryRunRegistry::string, useDeltaR_, and useMaxDist_.

66  {
67  // configurables
68  electronSrcToken_ = mayConsume<std::vector<pat::Electron> >(iConfig.getParameter<edm::InputTag>("electronSource"));
69  muonSrcToken_ = mayConsume<std::vector<pat::Muon> >(iConfig.getParameter<edm::InputTag>("muonSource"));
70  metSrcToken_ = consumes<std::vector<pat::MET> >(iConfig.getParameter<edm::InputTag>("metSource"));
71  jetSrcToken_ = consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("jetSource"));
72  leptonFlavour_ = iConfig.getParameter<std::string>("leptonFlavour");
73  jetCorrScheme_ = iConfig.getParameter<int>("jetCorrectionScheme");
74  nrCombJets_ = iConfig.getParameter<unsigned int>("nrCombJets");
75  doKinFit_ = iConfig.getParameter<bool>("doKinFit");
76  addLRSignalSel_ = iConfig.getParameter<bool>("addLRSignalSel");
77  lrSignalSelObs_ = iConfig.getParameter<std::vector<int> >("lrSignalSelObs");
78  lrSignalSelFile_ = iConfig.getParameter<std::string>("lrSignalSelFile");
79  addLRJetComb_ = iConfig.getParameter<bool>("addLRJetComb");
80  lrJetCombObs_ = iConfig.getParameter<std::vector<int> >("lrJetCombObs");
81  lrJetCombFile_ = iConfig.getParameter<std::string>("lrJetCombFile");
82  maxNrIter_ = iConfig.getParameter<int>("maxNrIter");
83  maxDeltaS_ = iConfig.getParameter<double>("maxDeltaS");
84  maxF_ = iConfig.getParameter<double>("maxF");
85  jetParam_ = iConfig.getParameter<int>("jetParametrisation");
86  lepParam_ = iConfig.getParameter<int>("lepParametrisation");
87  metParam_ = iConfig.getParameter<int>("metParametrisation");
88  constraints_ = iConfig.getParameter<std::vector<unsigned> >("constraints");
89  matchToGenEvt_ = iConfig.getParameter<bool>("matchToGenEvt");
90  matchingAlgo_ = iConfig.getParameter<int>("matchingAlgorithm");
91  useMaxDist_ = iConfig.getParameter<bool>("useMaximalDistance");
92  useDeltaR_ = iConfig.getParameter<bool>("useDeltaR");
93  maxDist_ = iConfig.getParameter<double>("maximalDistance");
94  genEvtToken_ = mayConsume<TtGenEvent>(edm::InputTag("genEvt"));
95 
96  // define kinfitter
97  if (doKinFit_) {
100  }
101 
102  // define jet combinations related calculators
106  if (addLRJetComb_)
108 
109  // instantiate signal selection calculator
110  if (addLRSignalSel_)
112 
113  // define what will be produced
114  produces<std::vector<TtSemiEvtSolution> >();
115 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
TtSemiLRJetCombCalc * myLRJetCombCalc
edm::EDGetTokenT< std::vector< pat::Electron > > electronSrcToken_
TtSemiLepKinFitter * myKinFitter
TtSemiLRJetCombObservables * myLRJetCombObservables
TtSemiSimpleBestJetComb * mySimpleBestJetComb
std::vector< unsigned > constraints_
TtSemiLepKinFitter::Param param(unsigned)
edm::EDGetTokenT< std::vector< pat::Muon > > muonSrcToken_
Simple method to get the correct jet combination in semileptonic ttbar events.
Class to calculate the jet combination LR value and purity from a root-file with fit functions...
edm::EDGetTokenT< std::vector< pat::MET > > metSrcToken_
Class to calculate the jet combination LR value and purity from a root-file with fit functions...
Steering class for the overall top-lepton likelihood.
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
edm::EDGetTokenT< TtGenEvent > genEvtToken_
std::vector< int > lrSignalSelObs_
std::vector< TtSemiLepKinFitter::Constraint > constraints(std::vector< unsigned > &)
TtSemiLRSignalSelObservables * myLRSignalSelObservables
std::vector< int > lrJetCombObs_
TtSemiLRSignalSelCalc * myLRSignalSelCalc

◆ ~TtSemiEvtSolutionMaker()

TtSemiEvtSolutionMaker::~TtSemiEvtSolutionMaker ( )
override

destructor

Definition at line 118 of file TtSemiEvtSolutionMaker.cc.

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

118  {
119  if (doKinFit_)
120  delete myKinFitter;
121  delete mySimpleBestJetComb;
123  delete myLRJetCombObservables;
124  if (addLRSignalSel_)
125  delete myLRSignalSelCalc;
126  if (addLRJetComb_)
127  delete myLRJetCombCalc;
128 }
TtSemiLRJetCombCalc * myLRJetCombCalc
TtSemiLepKinFitter * myKinFitter
TtSemiLRJetCombObservables * myLRJetCombObservables
TtSemiSimpleBestJetComb * mySimpleBestJetComb
TtSemiLRSignalSelObservables * myLRSignalSelObservables
TtSemiLRSignalSelCalc * myLRSignalSelCalc

Member Function Documentation

◆ constraint()

TtSemiLepKinFitter::Constraint TtSemiEvtSolutionMaker::constraint ( unsigned  val)

Definition at line 352 of file TtSemiEvtSolutionMaker.cc.

References Exception, TtSemiLepKinFitter::kNeutrinoMass, TtSemiLepKinFitter::kTopHadMass, TtSemiLepKinFitter::kTopLepMass, TtSemiLepKinFitter::kWHadMass, TtSemiLepKinFitter::kWLepMass, mps_fire::result, and heppy_batch::val.

Referenced by constraints().

◆ constraints()

std::vector< TtSemiLepKinFitter::Constraint > TtSemiEvtSolutionMaker::constraints ( std::vector< unsigned > &  val)

Definition at line 377 of file TtSemiEvtSolutionMaker.cc.

References constraint(), mps_fire::i, mps_fire::result, and heppy_batch::val.

Referenced by TtSemiEvtSolutionMaker().

377  {
378  std::vector<TtSemiLepKinFitter::Constraint> result;
379  for (unsigned i = 0; i < val.size(); ++i) {
380  result.push_back(constraint(val[i]));
381  }
382  return result;
383 }
TtSemiLepKinFitter::Constraint constraint(unsigned)

◆ param()

TtSemiLepKinFitter::Param TtSemiEvtSolutionMaker::param ( unsigned  val)

Definition at line 333 of file TtSemiEvtSolutionMaker.cc.

References Exception, TopKinFitter::kEMom, TopKinFitter::kEtEtaPhi, TopKinFitter::kEtThetaPhi, mps_fire::result, and heppy_batch::val.

Referenced by TtSemiEvtSolutionMaker().

333  {
335  switch (val) {
338  break;
341  break;
344  break;
345  default:
346  throw cms::Exception("WrongConfig") << "Chosen jet parametrization is not supported: " << val << "\n";
347  break;
348  }
349  return result;
350 }
Param
supported parameterizations
Definition: TopKinFitter.h:22

◆ produce()

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

Definition at line 130 of file TtSemiEvtSolutionMaker.cc.

References TtSemiLepKinFitter::addKinFitInfo(), addLRJetComb_, addLRSignalSel_, doKinFit_, pwdgSkimBPark_cfi::electrons, electronSrcToken_, TtGenEvtProducer_cfi::genEvt, genEvtToken_, EgammaValidation_cff::genp, JetPartonMatching::getDistanceForParton(), JetPartonMatching::getMatchForParton(), JetPartonMatching::getSumDistances(), iEvent, jetCorrScheme_, jetParam_, PDWG_EXODelayedJetMET_cff::jets, jetSrcToken_, lepParam_, leptonFlavour_, matchingAlgo_, matchToGenEvt_, maxDist_, metParam_, singleTopDQM_cfi::mets, metSrcToken_, eostools::move(), DiMuonV_cfg::muons, muonSrcToken_, myKinFitter, TtDilepEvtSolProducer_cfi::nrCombJets, nrCombJets_, AlCaHLTBitMon_ParallelJobs::p, submitPVResolutionJobs::q, alignCSCRings::s, TtSemiEvtSolution::setElectron(), TtSemiEvtSolution::setGenEvt(), TtSemiEvtSolution::setHadb(), TtSemiEvtSolution::setHadp(), TtSemiEvtSolution::setHadq(), TtSemiEvtSolution::setJetCorrectionScheme(), TtSemiEvtSolution::setJetParametrisation(), TtSemiEvtSolution::setLepb(), TtSemiEvtSolution::setLeptonParametrisation(), TtSemiEvtSolution::setMuon(), TtSemiEvtSolution::setNeutrino(), TtSemiEvtSolution::setNeutrinoParametrisation(), useDeltaR_, and useMaxDist_.

130  {
131  //
132  // TopObject Selection
133  //
134 
135  // select lepton (the TtLepton vectors are, for the moment, sorted on pT)
136  bool leptonFound = false;
138  if (leptonFlavour_ == "muon") {
139  iEvent.getByToken(muonSrcToken_, muons);
140  if (!muons->empty())
141  leptonFound = true;
142  }
144  if (leptonFlavour_ == "electron") {
145  iEvent.getByToken(electronSrcToken_, electrons);
146  if (!electrons->empty())
147  leptonFound = true;
148  }
149 
150  // select MET (TopMET vector is sorted on ET)
151  bool metFound = false;
153  iEvent.getByToken(metSrcToken_, mets);
154  if (!mets->empty())
155  metFound = true;
156 
157  // select Jets
158  bool jetsFound = false;
160  iEvent.getByToken(jetSrcToken_, jets);
161  if (jets->size() >= 4)
162  jetsFound = true;
163 
164  //
165  // Build Event solutions according to the ambiguity in the jet combination
166  //
167  std::vector<TtSemiEvtSolution>* evtsols = new std::vector<TtSemiEvtSolution>();
168  if (leptonFound && metFound && jetsFound) {
169  // protect against reading beyond array boundaries
170  unsigned int nrCombJets = nrCombJets_; // do not overwrite nrCombJets_
171  if (jets->size() < nrCombJets)
172  nrCombJets = jets->size();
173  // loop over all jets
174  for (unsigned int p = 0; p < nrCombJets; p++) {
175  for (unsigned int q = 0; q < nrCombJets; q++) {
176  for (unsigned int bh = 0; bh < nrCombJets; bh++) {
177  if (q > p && !(bh == p || bh == q)) {
178  for (unsigned int bl = 0; bl < nrCombJets; bl++) {
179  if (!(bl == p || bl == q || bl == bh)) {
180  TtSemiEvtSolution asol;
182  if (leptonFlavour_ == "muon")
183  asol.setMuon(muons, 0);
184  if (leptonFlavour_ == "electron")
185  asol.setElectron(electrons, 0);
186  asol.setNeutrino(mets, 0);
187  asol.setHadp(jets, p);
188  asol.setHadq(jets, q);
189  asol.setHadb(jets, bh);
190  asol.setLepb(jets, bl);
191  if (doKinFit_) {
192  asol = myKinFitter->addKinFitInfo(&asol);
193  // just to keep a record in the event (drop? -> present in provenance anyway...)
197  }
198  if (matchToGenEvt_) {
200  iEvent.getByToken(genEvtToken_, genEvt);
201  if (genEvt->numberOfBQuarks() ==
202  2 && // FIXME: in rare cases W->bc decay, resulting in a wrong filled genEvt leading to a segmentation fault
203  genEvt->numberOfLeptons() ==
204  1) { // FIXME: temporary solution to avoid crash in JetPartonMatching for non semi-leptonic events
205  asol.setGenEvt(genEvt);
206  }
207  }
208  // these lines calculate the observables to be used in the TtSemiSignalSelection LR
209  (*myLRSignalSelObservables)(asol, *jets);
210 
211  // if asked for, calculate with these observable values the LRvalue and
212  // (depending on the configuration) probability this event is signal
213  // FIXME: DO WE NEED TO DO THIS FOR EACH SOLUTION??? (S.L.20/8/07)
214  if (addLRSignalSel_)
215  (*myLRSignalSelCalc)(asol);
216 
217  // these lines calculate the observables to be used in the TtSemiJetCombination LR
218  //(*myLRJetCombObservables)(asol);
219 
220  (*myLRJetCombObservables)(asol, iEvent);
221 
222  // if asked for, calculate with these observable values the LRvalue and
223  // (depending on the configuration) probability a jet combination is correct
224  if (addLRJetComb_)
225  (*myLRJetCombCalc)(asol);
226 
227  //std::cout<<"SignalSelLRval = "<<asol.getLRSignalEvtLRval()<<" JetCombProb = "<<asol.getLRSignalEvtProb()<<std::endl;
228  //std::cout<<"JetCombLRval = "<<asol.getLRJetCombLRval()<<" JetCombProb = "<<asol.getLRJetCombProb()<<std::endl;
229 
230  // fill solution to vector
231  asol.setupHyp();
232  evtsols->push_back(asol);
233  }
234  }
235  }
236  }
237  }
238  }
239 
240  // if asked for, match the event solutions to the gen Event
241  if (matchToGenEvt_) {
242  int bestSolution = -999;
243  int bestSolutionChangeWQ = -999;
245  iEvent.getByToken(genEvtToken_, genEvt);
246  if (genEvt->numberOfBQuarks() ==
247  2 && // FIXME: in rare cases W->bc decay, resulting in a wrong filled genEvt leading to a segmentation fault
248  genEvt->numberOfLeptons() ==
249  1) { // FIXME: temporary solution to avoid crash in JetPartonMatching for non semi-leptonic events
250  std::vector<const reco::Candidate*> quarks;
251  const reco::Candidate& genp = *(genEvt->hadronicDecayQuark());
252  const reco::Candidate& genq = *(genEvt->hadronicDecayQuarkBar());
253  const reco::Candidate& genbh = *(genEvt->hadronicDecayB());
254  const reco::Candidate& genbl = *(genEvt->leptonicDecayB());
255  quarks.push_back(&genp);
256  quarks.push_back(&genq);
257  quarks.push_back(&genbh);
258  quarks.push_back(&genbl);
259  std::vector<const reco::Candidate*> recjets;
260  for (size_t s = 0; s < evtsols->size(); s++) {
261  recjets.clear();
262  const reco::Candidate& jetp = (*evtsols)[s].getRecHadp();
263  const reco::Candidate& jetq = (*evtsols)[s].getRecHadq();
264  const reco::Candidate& jetbh = (*evtsols)[s].getRecHadb();
265  const reco::Candidate& jetbl = (*evtsols)[s].getRecLepb();
266  recjets.push_back(&jetp);
267  recjets.push_back(&jetq);
268  recjets.push_back(&jetbh);
269  recjets.push_back(&jetbl);
270  JetPartonMatching aMatch(quarks, recjets, matchingAlgo_, useMaxDist_, useDeltaR_, maxDist_);
271  (*evtsols)[s].setGenEvt(genEvt);
272  (*evtsols)[s].setMCBestSumAngles(aMatch.getSumDistances());
273  (*evtsols)[s].setMCBestAngleHadp(aMatch.getDistanceForParton(0));
274  (*evtsols)[s].setMCBestAngleHadq(aMatch.getDistanceForParton(1));
275  (*evtsols)[s].setMCBestAngleHadb(aMatch.getDistanceForParton(2));
276  (*evtsols)[s].setMCBestAngleLepb(aMatch.getDistanceForParton(3));
277  if (aMatch.getMatchForParton(2) == 2 && aMatch.getMatchForParton(3) == 3) {
278  if (aMatch.getMatchForParton(0) == 0 && aMatch.getMatchForParton(1) == 1) {
279  bestSolution = s;
280  bestSolutionChangeWQ = 0;
281  } else if (aMatch.getMatchForParton(0) == 1 && aMatch.getMatchForParton(1) == 0) {
282  bestSolution = s;
283  bestSolutionChangeWQ = 1;
284  }
285  }
286  }
287  }
288  for (size_t s = 0; s < evtsols->size(); s++) {
289  (*evtsols)[s].setMCBestJetComb(bestSolution);
290  (*evtsols)[s].setMCChangeWQ(bestSolutionChangeWQ);
291  }
292  }
293 
294  // add TtSemiSimpleBestJetComb to solutions
295  int simpleBestJetComb = (*mySimpleBestJetComb)(*evtsols);
296  for (size_t s = 0; s < evtsols->size(); s++)
297  (*evtsols)[s].setSimpleBestJetComb(simpleBestJetComb);
298 
299  // choose the best jet combination according to LR value
300  if (addLRJetComb_ && !evtsols->empty()) {
301  float bestLRVal = -1000000;
302  int bestSol = (*evtsols)[0].getLRBestJetComb(); // duplicate the default
303  for (size_t s = 0; s < evtsols->size(); s++) {
304  if ((*evtsols)[s].getLRJetCombLRval() > bestLRVal) {
305  bestLRVal = (*evtsols)[s].getLRJetCombLRval();
306  bestSol = s;
307  }
308  }
309  for (size_t s = 0; s < evtsols->size(); s++) {
310  (*evtsols)[s].setLRBestJetComb(bestSol);
311  }
312  }
313 
314  //store the vector of solutions to the event
315  std::unique_ptr<std::vector<TtSemiEvtSolution> > pOut(evtsols);
316  iEvent.put(std::move(pOut));
317 
318  } else {
319  /*
320  std::cout<<"No calibrated solutions built, because: ";
321  if(jets->size()<4) std::cout<<"nr sel jets < 4"<<std::endl;
322  if(leptonFlavour_ == "muon" && muons->size() == 0) std::cout<<"no good muon candidate"<<std::endl;
323  if(leptonFlavour_ == "electron" && electrons->size() == 0) std::cout<<"no good electron candidate"<<std::endl;
324  if(mets->size() == 0) std::cout<<"no MET reconstruction"<<std::endl;
325  */
326  // TtSemiEvtSolution asol;
327  // evtsols->push_back(asol);
328  std::unique_ptr<std::vector<TtSemiEvtSolution> > pOut(evtsols);
329  iEvent.put(std::move(pOut));
330  }
331 }
void setLepb(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
edm::EDGetTokenT< std::vector< pat::Electron > > electronSrcToken_
genp
produce generated paricles in acceptance #
TtSemiLepKinFitter * myKinFitter
void setGenEvt(const edm::Handle< TtGenEvent > &aGenEvt)
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
edm::EDGetTokenT< std::vector< pat::Muon > > muonSrcToken_
void setHadp(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
void setLeptonParametrisation(int lp)
void setMuon(const edm::Handle< std::vector< pat::Muon > > &muon, int i)
int iEvent
Definition: GenABIO.cc:224
void setJetCorrectionScheme(int scheme)
edm::EDGetTokenT< std::vector< pat::MET > > metSrcToken_
void setJetParametrisation(int jp)
void setElectron(const edm::Handle< std::vector< pat::Electron > > &elec, int i)
TtSemiEvtSolution addKinFitInfo(TtSemiEvtSolution *asol)
add kin fit information to the old event solution (in for legacy reasons)
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
void setNeutrinoParametrisation(int mp)
edm::EDGetTokenT< TtGenEvent > genEvtToken_
void setNeutrino(const edm::Handle< std::vector< pat::MET > > &met, int i)
void setHadb(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
def move(src, dest)
Definition: eostools.py:511
void setHadq(const edm::Handle< std::vector< pat::Jet > > &jet, int i)

Member Data Documentation

◆ addLRJetComb_

bool TtSemiEvtSolutionMaker::addLRJetComb_
private

◆ addLRSignalSel_

bool TtSemiEvtSolutionMaker::addLRSignalSel_
private

◆ constraints_

std::vector<unsigned> TtSemiEvtSolutionMaker::constraints_
private

Definition at line 54 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ doKinFit_

bool TtSemiEvtSolutionMaker::doKinFit_
private

◆ electronSrcToken_

edm::EDGetTokenT<std::vector<pat::Electron> > TtSemiEvtSolutionMaker::electronSrcToken_
private

Definition at line 38 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ genEvtToken_

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

Definition at line 55 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ jetCorrScheme_

int TtSemiEvtSolutionMaker::jetCorrScheme_
private

Definition at line 43 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ jetParam_

int TtSemiEvtSolutionMaker::jetParam_
private

Definition at line 52 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ jetSrcToken_

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

Definition at line 41 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ lepParam_

int TtSemiEvtSolutionMaker::lepParam_
private

Definition at line 52 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ leptonFlavour_

std::string TtSemiEvtSolutionMaker::leptonFlavour_
private

Definition at line 42 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ lrJetCombFile_

std::string TtSemiEvtSolutionMaker::lrJetCombFile_
private

Definition at line 45 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ lrJetCombObs_

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

Definition at line 53 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ lrSignalSelFile_

std::string TtSemiEvtSolutionMaker::lrSignalSelFile_
private

Definition at line 45 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ lrSignalSelObs_

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

Definition at line 53 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ matchingAlgo_

int TtSemiEvtSolutionMaker::matchingAlgo_
private

Definition at line 47 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ matchToGenEvt_

bool TtSemiEvtSolutionMaker::matchToGenEvt_
private

Definition at line 46 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ maxDeltaS_

double TtSemiEvtSolutionMaker::maxDeltaS_
private

Definition at line 51 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ maxDist_

double TtSemiEvtSolutionMaker::maxDist_
private

Definition at line 49 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ maxF_

double TtSemiEvtSolutionMaker::maxF_
private

Definition at line 51 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ maxNrIter_

int TtSemiEvtSolutionMaker::maxNrIter_
private

Definition at line 50 of file TtSemiEvtSolutionMaker.cc.

Referenced by TtSemiEvtSolutionMaker().

◆ metParam_

int TtSemiEvtSolutionMaker::metParam_
private

Definition at line 52 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ metSrcToken_

edm::EDGetTokenT<std::vector<pat::MET> > TtSemiEvtSolutionMaker::metSrcToken_
private

Definition at line 40 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ muonSrcToken_

edm::EDGetTokenT<std::vector<pat::Muon> > TtSemiEvtSolutionMaker::muonSrcToken_
private

Definition at line 39 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ myKinFitter

TtSemiLepKinFitter* TtSemiEvtSolutionMaker::myKinFitter
private

◆ myLRJetCombCalc

TtSemiLRJetCombCalc* TtSemiEvtSolutionMaker::myLRJetCombCalc
private

Definition at line 60 of file TtSemiEvtSolutionMaker.cc.

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

◆ myLRJetCombObservables

TtSemiLRJetCombObservables* TtSemiEvtSolutionMaker::myLRJetCombObservables
private

Definition at line 59 of file TtSemiEvtSolutionMaker.cc.

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

◆ myLRSignalSelCalc

TtSemiLRSignalSelCalc* TtSemiEvtSolutionMaker::myLRSignalSelCalc
private

Definition at line 62 of file TtSemiEvtSolutionMaker.cc.

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

◆ myLRSignalSelObservables

TtSemiLRSignalSelObservables* TtSemiEvtSolutionMaker::myLRSignalSelObservables
private

Definition at line 61 of file TtSemiEvtSolutionMaker.cc.

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

◆ mySimpleBestJetComb

TtSemiSimpleBestJetComb* TtSemiEvtSolutionMaker::mySimpleBestJetComb
private

Definition at line 58 of file TtSemiEvtSolutionMaker.cc.

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

◆ nrCombJets_

unsigned int TtSemiEvtSolutionMaker::nrCombJets_
private

Definition at line 44 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ useDeltaR_

bool TtSemiEvtSolutionMaker::useDeltaR_
private

Definition at line 48 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().

◆ useMaxDist_

bool TtSemiEvtSolutionMaker::useMaxDist_
private

Definition at line 48 of file TtSemiEvtSolutionMaker.cc.

Referenced by produce(), and TtSemiEvtSolutionMaker().