#include <TopQuarkAnalysis/TopEventProducers/interface/TtDilepEvtSolutionMaker.h>
Public Member Functions | |
virtual void | produce (edm::Event &iEvent, const edm::EventSetup &iSetup) |
TtDilepEvtSolutionMaker (const edm::ParameterSet &iConfig) | |
~TtDilepEvtSolutionMaker () | |
Private Member Functions | |
bool | HasPositiveCharge (const reco::Candidate *) const |
bool | LepDiffCharge (const reco::Candidate *, const reco::Candidate *) const |
bool | PTComp (const reco::Candidate *, const reco::Candidate *) const |
Private Attributes | |
bool | calcTopMass_ |
bool | eeChannel_ |
edm::InputTag | electronSource_ |
bool | emuChannel_ |
bool | etauChannel_ |
edm::InputTag | evtSource_ |
int | jetCorrScheme_ |
edm::InputTag | jetSource_ |
bool | matchToGenEvt_ |
edm::InputTag | metSource_ |
bool | mumuChannel_ |
edm::InputTag | muonSource_ |
bool | mutauChannel_ |
TtDilepLRSignalSelObservables * | myLRSignalSelObservables |
unsigned int | nrCombJets_ |
edm::InputTag | tauSource_ |
bool | tautauChannel_ |
double | tmassbegin_ |
double | tmassend_ |
double | tmassstep_ |
bool | useMCforBest_ |
Definition at line 13 of file TtDilepEvtSolutionMaker.h.
TtDilepEvtSolutionMaker::TtDilepEvtSolutionMaker | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 18 of file TtDilepEvtSolutionMaker.cc.
References calcTopMass_, eeChannel_, electronSource_, emuChannel_, etauChannel_, evtSource_, edm::ParameterSet::getParameter(), jetCorrScheme_, jetSource_, matchToGenEvt_, metSource_, mumuChannel_, muonSource_, mutauChannel_, myLRSignalSelObservables, nrCombJets_, tauSource_, tautauChannel_, tmassbegin_, tmassend_, tmassstep_, and useMCforBest_.
00019 { 00020 // configurables 00021 electronSource_ = iConfig.getParameter<edm::InputTag>("electronSource"); 00022 muonSource_ = iConfig.getParameter<edm::InputTag>("muonSource"); 00023 tauSource_ = iConfig.getParameter<edm::InputTag>("tauSource"); 00024 metSource_ = iConfig.getParameter<edm::InputTag>("metSource"); 00025 jetSource_ = iConfig.getParameter<edm::InputTag>("jetSource"); 00026 jetCorrScheme_ = iConfig.getParameter<int> ("jetCorrectionScheme"); 00027 evtSource_ = iConfig.getParameter<edm::InputTag>("evtSource"); 00028 nrCombJets_ = iConfig.getParameter<unsigned int> ("nrCombJets"); 00029 matchToGenEvt_ = iConfig.getParameter<bool> ("matchToGenEvt"); 00030 calcTopMass_ = iConfig.getParameter<bool> ("calcTopMass"); 00031 useMCforBest_ = iConfig.getParameter<bool> ("bestSolFromMC"); 00032 eeChannel_ = iConfig.getParameter<bool> ("eeChannel"); 00033 emuChannel_ = iConfig.getParameter<bool> ("emuChannel"); 00034 mumuChannel_ = iConfig.getParameter<bool> ("mumuChannel"); 00035 mutauChannel_ = iConfig.getParameter<bool> ("mutauChannel"); 00036 etauChannel_ = iConfig.getParameter<bool> ("etauChannel"); 00037 tautauChannel_ = iConfig.getParameter<bool> ("tautauChannel"); 00038 tmassbegin_ = iConfig.getParameter<double> ("tmassbegin"); 00039 tmassend_ = iConfig.getParameter<double> ("tmassend"); 00040 tmassstep_ = iConfig.getParameter<double> ("tmassstep"); 00041 00042 // define what will be produced 00043 produces<std::vector<TtDilepEvtSolution> >(); 00044 00045 myLRSignalSelObservables = new TtDilepLRSignalSelObservables(); 00046 myLRSignalSelObservables->jetSource(jetSource_); 00047 }
TtDilepEvtSolutionMaker::~TtDilepEvtSolutionMaker | ( | ) |
bool TtDilepEvtSolutionMaker::HasPositiveCharge | ( | const reco::Candidate * | l | ) | const [inline, private] |
Definition at line 55 of file TtDilepEvtSolutionMaker.h.
References reco::Particle::charge().
Referenced by produce().
00056 { 00057 return (l->charge() > 0); 00058 }
bool TtDilepEvtSolutionMaker::LepDiffCharge | ( | const reco::Candidate * | l1, | |
const reco::Candidate * | l2 | |||
) | const [inline, private] |
Definition at line 50 of file TtDilepEvtSolutionMaker.h.
References reco::Particle::charge().
Referenced by produce().
void TtDilepEvtSolutionMaker::produce | ( | edm::Event & | iEvent, | |
const edm::EventSetup & | iSetup | |||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 53 of file TtDilepEvtSolutionMaker.cc.
References TtFullLepKinSolver::addKinSolInfo(), calcTopMass_, GenMuonPlsPt100GeV_cfg::cout, eeChannel_, electronSource_, emuChannel_, lat::endl(), etauChannel_, evtSource_, find(), TtGenEvtProducer_cfi::genEvt, edm::Event::getByLabel(), TtDilepEvtSolution::getJetResidual(), HasPositiveCharge(), int, metsig::jet, jetCorrScheme_, pfTauBenchmarkGeneric_cfi::jets, jetSource_, LepDiffCharge(), matchToGenEvt_, metSource_, mumuChannel_, muons_cfi::muons, muonSource_, mutauChannel_, nrCombJets_, NULL, PTComp(), edm::Event::put(), s, TtDilepEvtSolution::setB(), TtDilepEvtSolution::setBbar(), TtDilepEvtSolution::setElectronm(), TtDilepEvtSolution::setElectronp(), TtDilepEvtSolution::setGenEvt(), TtDilepEvtSolution::setJetCorrectionScheme(), TtDilepEvtSolution::setMET(), TtDilepEvtSolution::setMuonm(), TtDilepEvtSolution::setMuonp(), TtDilepEvtSolution::setTaum(), TtDilepEvtSolution::setTaup(), metsig::tau, tauSource_, tautauChannel_, tmassbegin_, tmassend_, tmassstep_, useMCforBest_, and TtFullLepKinSolver::useWeightFromMC().
00054 { 00055 using namespace edm; 00056 Handle<std::vector<pat::Tau> > taus; 00057 iEvent.getByLabel(tauSource_, taus); 00058 Handle<std::vector<pat::Muon> > muons; 00059 iEvent.getByLabel(muonSource_, muons); 00060 Handle<std::vector<pat::Electron> > electrons; 00061 iEvent.getByLabel(electronSource_, electrons); 00062 Handle<std::vector<pat::MET> > mets; 00063 iEvent.getByLabel(metSource_, mets); 00064 Handle<std::vector<pat::Jet> > jets; 00065 iEvent.getByLabel(jetSource_, jets); 00066 00067 int selMuonp = -1, selMuonm = -1; 00068 int selElectronp = -1, selElectronm = -1; 00069 int selTaup = -1, selTaum = -1; 00070 bool leptonFound = false; 00071 bool mumu = false; 00072 bool emu = false; 00073 bool ee = false; 00074 bool etau = false; 00075 bool mutau = false; 00076 bool tautau = false; 00077 bool leptonFoundEE = false; 00078 bool leptonFoundMM = false; 00079 bool leptonFoundTT = false; 00080 bool leptonFoundEpMm = false; 00081 bool leptonFoundEmMp = false; 00082 bool leptonFoundEpTm = false; 00083 bool leptonFoundEmTp = false; 00084 bool leptonFoundMpTm = false; 00085 bool leptonFoundMmTp = false; 00086 bool jetsFound = false; 00087 bool METFound = false; 00088 std::vector<int> JetVetoByTaus; 00089 00090 //select MET (TopMET vector is sorted on ET) 00091 if(mets->size()>=1) { METFound = true; } 00092 00093 // If we have electrons and muons available, 00094 // build a solutions with electrons and muons. 00095 if (muons->size() + electrons->size() >=2) { 00096 // select leptons 00097 if (electrons->size() == 0) mumu = true; 00098 else if (muons->size() == 0) ee = true; 00099 else if (electrons->size() == 1) { 00100 if (muons->size() == 1) emu = true; 00101 else if (PTComp(&(*electrons)[0], &(*muons)[1])) emu = true; 00102 else mumu = true; 00103 } 00104 else if (electrons->size() > 1) { 00105 if (PTComp(&(*electrons)[1], &(*muons)[0])) ee = true; 00106 else if (muons->size() == 1) emu = true; 00107 else if (PTComp(&(*electrons)[0], &(*muons)[1])) emu = true; 00108 else mumu = true; 00109 } 00110 if (ee) { 00111 if (LepDiffCharge(&(*electrons)[0], &(*electrons)[1])) { 00112 leptonFound = true; 00113 leptonFoundEE = true; 00114 if (HasPositiveCharge(&(*electrons)[0])) { 00115 selElectronp = 0; 00116 selElectronm = 1; 00117 } else { 00118 selElectronp = 1; 00119 selElectronm = 0; 00120 } 00121 } 00122 } 00123 else if (emu) { 00124 if (LepDiffCharge(&(*electrons)[0], &(*muons)[0])) { 00125 leptonFound = true; 00126 if (HasPositiveCharge(&(*electrons)[0])) { 00127 leptonFoundEpMm = true; 00128 selElectronp = 0; 00129 selMuonm = 0; 00130 } else { 00131 leptonFoundEmMp = true; 00132 selMuonp = 0; 00133 selElectronm = 0; 00134 } 00135 } 00136 } 00137 else if (mumu) { 00138 if (LepDiffCharge(&(*muons)[0], &(*muons)[1])) { 00139 leptonFound = true; 00140 leptonFoundMM = true; 00141 if (HasPositiveCharge(&(*muons)[0])) { 00142 selMuonp = 0; 00143 selMuonm = 1; 00144 } else { 00145 selMuonp = 1; 00146 selMuonm = 0; 00147 } 00148 } 00149 } 00150 //select Jets (TopJet vector is sorted on ET) 00151 if(jets->size()>=2) { jetsFound = true; } 00152 } 00153 // If a tau is needed to have two leptons, then only consider the taus. 00154 // This is the minimal modification of the dilept selection that includes taus, 00155 // since we are considering taus only when no other solution exist. 00156 else if(muons->size() + electrons->size()==1 && taus->size()>0) { 00157 // select leptons 00158 if(muons->size()==1) { 00159 mutau = true; 00160 // depending on the muon charge, set the right muon index and specify channel 00161 int expectedCharge = - muons->begin()->charge(); 00162 int* tauIdx = NULL; 00163 if (expectedCharge<0) { 00164 selMuonp = 0; 00165 tauIdx = &selTaum; 00166 leptonFoundMpTm = true; 00167 } else { 00168 selMuonm = 0; 00169 tauIdx = &selTaup; 00170 leptonFoundMmTp = true; 00171 } 00172 // loop over the vector of taus to find the ones 00173 // that have the charge opposite to the muon one, and do not match in eta-phi 00174 std::vector<std::vector<pat::Tau>::const_iterator> subset1; 00175 for(std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau ) { 00176 if(tau->charge()*expectedCharge>=0 && DeltaR<reco::Particle>()(*tau,*(muons->begin()))>0.1) { 00177 *tauIdx = tau-taus->begin(); 00178 leptonFound = true; 00179 subset1.push_back(tau); 00180 } 00181 } 00182 // if there are more than one tau with ecalIsol==0, take the smallest E/P 00183 float iso = 999.; 00184 for(std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin(); tau < subset1.end(); ++tau) { 00185 if((*tau)->isCaloTau() && (*tau)->isolationTracksPtSum()<iso) { 00186 *tauIdx = *tau - taus->begin(); 00187 iso = (*tau)->isolationTracksPtSum(); 00188 } 00189 if((*tau)->isPFTau() && (*tau)->isolationPFChargedHadrCandsPtSum()<iso) { 00190 *tauIdx = *tau - taus->begin(); 00191 iso = (*tau)->isolationPFChargedHadrCandsPtSum(); 00192 } 00193 } 00194 00195 // check that one combination has been found 00196 if(!leptonFound) { leptonFoundMpTm = false; leptonFoundMmTp = false; } 00197 // discard the jet that matches the tau (if one) 00198 if(leptonFound) { 00199 for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet<jets->end(); ++jet) { 00200 if(DeltaR<reco::Particle>()(*(taus->begin()+*tauIdx),*jet)<0.1) { 00201 JetVetoByTaus.push_back(jet-jets->begin()); 00202 } 00203 } 00204 } 00205 } 00206 else { 00207 etau = true; 00208 // depending on the electron charge, set the right electron index and specify channel 00209 int expectedCharge = - electrons->begin()->charge(); 00210 int* tauIdx = NULL; 00211 if (expectedCharge<0) { 00212 selElectronp = 0; 00213 tauIdx = &selTaum; 00214 leptonFoundEpTm = true; 00215 } else { 00216 selElectronm = 0; 00217 tauIdx = &selTaup; 00218 leptonFoundEmTp = true; 00219 } 00220 // loop over the vector of taus to find the ones 00221 // that have the charge opposite to the muon one, and do not match in eta-phi 00222 std::vector<std::vector<pat::Tau>::const_iterator> subset1; 00223 for(std::vector<pat::Tau>::const_iterator tau = taus->begin(); tau < taus->end(); ++tau ) { 00224 if(tau->charge()*expectedCharge>=0 && DeltaR<reco::Particle>()(*tau,*(electrons->begin()))>0.1) { 00225 *tauIdx = tau-taus->begin(); 00226 leptonFound = true; 00227 subset1.push_back(tau); 00228 } 00229 } 00230 // if there are more than one tau with ecalIsol==0, take the smallest E/P 00231 float iso = 999.; 00232 for(std::vector<std::vector<pat::Tau>::const_iterator>::const_iterator tau = subset1.begin(); tau < subset1.end(); ++tau) { 00233 if((*tau)->isCaloTau() && (*tau)->isolationTracksPtSum()<iso) { 00234 *tauIdx = *tau - taus->begin(); 00235 iso = (*tau)->isolationTracksPtSum(); 00236 } 00237 if((*tau)->isPFTau() && (*tau)->isolationPFChargedHadrCandsPtSum()<iso) { 00238 *tauIdx = *tau - taus->begin(); 00239 iso = (*tau)->isolationPFChargedHadrCandsPtSum(); 00240 } 00241 } 00242 00243 // check that one combination has been found 00244 if(!leptonFound) { leptonFoundEpTm = false; leptonFoundEmTp = false; } 00245 // discard the jet that matches the tau (if one) 00246 if(leptonFound) { 00247 for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet<jets->end(); ++jet) { 00248 if(DeltaR<reco::Particle>()(*(taus->begin()+*tauIdx),*jet)<0.1) { 00249 JetVetoByTaus.push_back(jet-jets->begin()); 00250 } 00251 } 00252 } 00253 } 00254 // select Jets (TopJet vector is sorted on ET) 00255 jetsFound = ((jets->size()-JetVetoByTaus.size())>=2); 00256 } else if(taus->size()>1) { 00257 tautau = true; 00258 if(LepDiffCharge(&(*taus)[0],&(*taus)[1])) { 00259 leptonFound = true; 00260 leptonFoundTT = true; 00261 if(HasPositiveCharge(&(*taus)[0])) { 00262 selTaup = 0; 00263 selTaum = 1; 00264 } 00265 else { 00266 selTaup = 1; 00267 selTaum = 0; 00268 } 00269 } 00270 for(std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet<jets->end(); ++jet) { 00271 if(DeltaR<reco::Particle>()((*taus)[0],*jet)<0.1 || DeltaR<reco::Particle>()((*taus)[1],*jet)<0.1) { 00272 JetVetoByTaus.push_back(jet-jets->begin()); 00273 } 00274 } 00275 // select Jets (TopJet vector is sorted on ET) 00276 jetsFound = ((jets->size()-JetVetoByTaus.size())>=2); 00277 } 00278 00279 // Check that the above work makes sense 00280 if(int(ee)+int(emu)+int(mumu)+int(etau)+int(mutau)+int(tautau)>1) 00281 std::cout << "[TtDilepEvtSolutionMaker]: " 00282 << "Lepton selection criteria uncorrectly defined" << std::endl; 00283 00284 bool correctLepton = (leptonFoundEE && eeChannel_) || 00285 ((leptonFoundEmMp || leptonFoundEpMm) && emuChannel_) || 00286 (leptonFoundMM && mumuChannel_) || 00287 ((leptonFoundMmTp || leptonFoundMpTm) && mutauChannel_)|| 00288 ((leptonFoundEmTp || leptonFoundEpTm) && etauChannel_) || 00289 (leptonFoundTT && tautauChannel_) ; 00290 00291 std::vector<TtDilepEvtSolution> * evtsols = new std::vector<TtDilepEvtSolution>(); 00292 if(correctLepton && METFound && jetsFound) { 00293 // protect against reading beyond array boundaries while discounting vetoed jets 00294 unsigned int nrCombJets = 0; 00295 unsigned int numberOfJets = 0; 00296 for(; nrCombJets<jets->size() && numberOfJets<nrCombJets_; ++nrCombJets) { 00297 if(find(JetVetoByTaus.begin(),JetVetoByTaus.end(),int(nrCombJets))==JetVetoByTaus.end()) ++numberOfJets; 00298 } 00299 // consider all permutations 00300 for (unsigned int ib = 0; ib < nrCombJets; ib++) { 00301 // skipped jet vetoed during components-flagging. 00302 if(find(JetVetoByTaus.begin(),JetVetoByTaus.end(),int(ib))!=JetVetoByTaus.end())continue; 00303 // second loop of the permutations 00304 for (unsigned int ibbar = 0; ibbar < nrCombJets; ibbar++) { 00305 // avoid the diagonal: b and bbar must be distinct jets 00306 if(ib==ibbar) continue; 00307 // skipped jet vetoed during components-flagging. 00308 if(find(JetVetoByTaus.begin(),JetVetoByTaus.end(),int(ibbar))!=JetVetoByTaus.end())continue; 00309 // Build and save a solution 00310 TtDilepEvtSolution asol; 00311 asol.setJetCorrectionScheme(jetCorrScheme_); 00312 double xconstraint = 0, yconstraint = 0; 00313 // Set e+ in the event 00314 if (leptonFoundEE || leptonFoundEpMm || leptonFoundEpTm) { 00315 asol.setElectronp(electrons, selElectronp); 00316 xconstraint += (*electrons)[selElectronp].px(); 00317 yconstraint += (*electrons)[selElectronp].py(); 00318 } 00319 // Set e- in the event 00320 if (leptonFoundEE || leptonFoundEmMp || leptonFoundEmTp) { 00321 asol.setElectronm(electrons, selElectronm); 00322 xconstraint += (*electrons)[selElectronm].px(); 00323 yconstraint += (*electrons)[selElectronm].py(); 00324 } 00325 // Set mu+ in the event 00326 if (leptonFoundMM || leptonFoundEmMp || leptonFoundMpTm) { 00327 asol.setMuonp(muons, selMuonp); 00328 xconstraint += (*muons)[selMuonp].px(); 00329 yconstraint += (*muons)[selMuonp].py(); 00330 } 00331 // Set mu- in the event 00332 if (leptonFoundMM || leptonFoundEpMm || leptonFoundMmTp) { 00333 asol.setMuonm(muons, selMuonm); 00334 xconstraint += (*muons)[selMuonm].px(); 00335 yconstraint += (*muons)[selMuonm].py(); 00336 } 00337 // Set tau- in the event 00338 if (leptonFoundEpTm || leptonFoundMpTm || leptonFoundTT) { 00339 asol.setTaum(taus, selTaum); 00340 xconstraint += (*taus)[selTaum].px(); 00341 yconstraint += (*taus)[selTaum].py(); 00342 } 00343 // Set tau+ in the event 00344 if (leptonFoundEmTp || leptonFoundMmTp || leptonFoundTT) { 00345 asol.setTaup(taus, selTaup); 00346 xconstraint += (*taus)[selTaup].px(); 00347 yconstraint += (*taus)[selTaup].py(); 00348 } 00349 // Set Jets/MET in the event 00350 asol.setB(jets, ib); 00351 asol.setBbar(jets, ibbar); 00352 asol.setMET(mets, 0); 00353 xconstraint += (*jets)[ib].px() + (*jets)[ibbar].px() + (*mets)[0].px(); 00354 yconstraint += (*jets)[ib].py() + (*jets)[ibbar].py() + (*mets)[0].py(); 00355 // if asked for, match the event solutions to the gen Event 00356 if(matchToGenEvt_){ 00357 Handle<TtGenEvent> genEvt; 00358 iEvent.getByLabel (evtSource_,genEvt); 00359 asol.setGenEvt(genEvt); 00360 } 00361 // If asked, use the kin fitter to compute the top mass 00362 if (calcTopMass_) { 00363 TtFullLepKinSolver solver(tmassbegin_, tmassend_, tmassstep_, xconstraint, yconstraint); 00364 solver.useWeightFromMC(useMCforBest_); 00365 asol = solver.addKinSolInfo(&asol); 00366 } 00367 00368 // these lines calculate the observables to be used in the TtDilepSignalSelection LR 00369 (*myLRSignalSelObservables)(asol, iEvent); 00370 00371 evtsols->push_back(asol); 00372 } 00373 } 00374 // flag the best solution (MC matching) 00375 if(matchToGenEvt_){ 00376 double bestSolDR = 9999.; 00377 int bestSol = -1; 00378 double dR = 0.; 00379 for(size_t s=0; s<evtsols->size(); s++) { 00380 dR = (*evtsols)[s].getJetResidual(); 00381 if(dR<bestSolDR) { bestSolDR = dR; bestSol = s; } 00382 } 00383 if(bestSol!=-1) (*evtsols)[bestSol].setBestSol(true); 00384 } 00385 // put the result in the event 00386 std::auto_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols); 00387 iEvent.put(pOut); 00388 } else { 00389 // no solution: put a dummy solution in the event 00390 TtDilepEvtSolution asol; 00391 evtsols->push_back(asol); 00392 std::auto_ptr<std::vector<TtDilepEvtSolution> > pOut(evtsols); 00393 iEvent.put(pOut); 00394 } 00395 }
bool TtDilepEvtSolutionMaker::PTComp | ( | const reco::Candidate * | l1, | |
const reco::Candidate * | l2 | |||
) | const [inline, private] |
Definition at line 45 of file TtDilepEvtSolutionMaker.h.
References reco::Particle::pt().
Referenced by produce().
bool TtDilepEvtSolutionMaker::calcTopMass_ [private] |
Definition at line 39 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::eeChannel_ [private] |
Definition at line 40 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
Definition at line 31 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::emuChannel_ [private] |
Definition at line 40 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::etauChannel_ [private] |
Definition at line 40 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
Definition at line 36 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
int TtDilepEvtSolutionMaker::jetCorrScheme_ [private] |
Definition at line 37 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
Definition at line 35 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::matchToGenEvt_ [private] |
Definition at line 39 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
Definition at line 34 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::mumuChannel_ [private] |
Definition at line 40 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
Definition at line 32 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::mutauChannel_ [private] |
Definition at line 40 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
unsigned int TtDilepEvtSolutionMaker::nrCombJets_ [private] |
Definition at line 38 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
Definition at line 33 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::tautauChannel_ [private] |
Definition at line 40 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
double TtDilepEvtSolutionMaker::tmassbegin_ [private] |
Definition at line 41 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
double TtDilepEvtSolutionMaker::tmassend_ [private] |
Definition at line 41 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
double TtDilepEvtSolutionMaker::tmassstep_ [private] |
Definition at line 41 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().
bool TtDilepEvtSolutionMaker::useMCforBest_ [private] |
Definition at line 39 of file TtDilepEvtSolutionMaker.h.
Referenced by produce(), and TtDilepEvtSolutionMaker().