1 #ifndef TtFullLepKinSolutionProducer_h 2 #define TtFullLepKinSolutionProducer_h 9 #include "TLorentzVector.h" 39 bool operator()(std::pair<double, int>
a, std::pair<double, int>
b){
40 return a.first > b.first;
60 return (l1->
pt() > l2->
pt());
93 produces<std::vector<std::vector<int> > > ();
94 produces<std::vector<reco::LeafCandidate> >(
"fullLepNeutrinos");
95 produces<std::vector<reco::LeafCandidate> >(
"fullLepNeutrinoBars");
96 produces<std::vector<double> >(
"solWeight");
97 produces<bool>(
"isWrongCharge");
117 std::vector<std::vector<int> > idcsV;
118 std::vector<reco::LeafCandidate> nusV;
119 std::vector<reco::LeafCandidate> nuBarsV;
120 std::vector<std::pair<double, int> > weightsV;
123 std::unique_ptr<std::vector<std::vector<int> > > pIdcs(
new std::vector<std::vector<int> >);
124 std::unique_ptr<std::vector<reco::LeafCandidate> > pNus(
new std::vector<reco::LeafCandidate>);
125 std::unique_ptr<std::vector<reco::LeafCandidate> > pNuBars(
new std::vector<reco::LeafCandidate>);
126 std::unique_ptr<std::vector<double> > pWeight(
new std::vector<double>);
127 std::unique_ptr<bool> pWrongCharge(
new bool);
138 int selMuon1 = -1, selMuon2 = -1;
139 int selElectron1 = -1, selElectron2 = -1;
143 bool isWrongCharge =
false;
144 bool jetsFound =
false;
145 bool METFound =
false;
146 bool electronsFound =
false;
147 bool electronMuonFound =
false;
148 bool muonsFound =
false;
151 if(jets->size()>=2) { jetsFound =
true; }
154 if(!mets->empty()) { METFound =
true; }
158 if(muons->size() + electrons->size() >=2) {
160 if(electrons->empty()) mumu =
true;
161 else if(muons->empty()) ee =
true;
162 else if(electrons->size() == 1) {
163 if(muons->size() == 1) emu =
true;
164 else if(
PTComp(&(*electrons)[0], &(*muons)[1])) emu =
true;
167 else if(electrons->size() > 1) {
168 if(
PTComp(&(*electrons)[1], &(*muons)[0])) ee =
true;
169 else if(muons->size() == 1) emu =
true;
170 else if(
PTComp(&(*electrons)[0], &(*muons)[1])) emu =
true;
183 electronsFound =
true;
184 if(!
LepDiffCharge(&(*electrons)[0], &(*electrons)[1])) isWrongCharge =
true;
191 electronMuonFound =
true;
192 if(!
LepDiffCharge(&(*electrons)[0], &(*muons)[0])) isWrongCharge =
true;
206 if(!
LepDiffCharge(&(*muons)[0], &(*muons)[1])) isWrongCharge =
true;
211 *pWrongCharge = isWrongCharge;
214 if(
int(ee)+
int(emu)+
int(mumu)>1) {
215 edm::LogWarning(
"TtFullLepKinSolutionProducer") <<
"Lepton selection criteria uncorrectly defined";
223 if(correctLeptons && METFound && jetsFound) {
229 if(jets->size()<
static_cast<unsigned int>(stop) || stop<0) stop=jets->size();
235 for (
int ib = 0;
ib<stop;
ib++) {
237 for (
int ibbar = 0; ibbar<stop; ibbar++) {
239 if(
ib==ibbar)
continue;
241 std::vector<int> idcs;
245 idcs.push_back(ibbar);
247 TLorentzVector LV_l1;
248 TLorentzVector LV_l2;
249 TLorentzVector LV_b = TLorentzVector((*jets)[
ib].correctedJet(
jetCorrLevel_,
"bottom").px(),
253 TLorentzVector LV_bbar = TLorentzVector((*jets)[ibbar].correctedJet(
jetCorrLevel_,
"bottom").px(),
258 double xconstraint = 0, yconstraint = 0;
261 idcs.push_back(selElectron1);
262 LV_l1.SetXYZT((*electrons)[selElectron1].px(),
263 (*electrons)[selElectron1].py(),
264 (*electrons)[selElectron1].pz(),
265 (*electrons)[selElectron1].
energy());
266 xconstraint += (*electrons)[selElectron1].px();
267 yconstraint += (*electrons)[selElectron1].py();
269 idcs.push_back(selElectron2);
270 LV_l2.SetXYZT((*electrons)[selElectron2].px(),
271 (*electrons)[selElectron2].py(),
272 (*electrons)[selElectron2].pz(),
273 (*electrons)[selElectron2].
energy());
274 xconstraint += (*electrons)[selElectron2].px();
275 yconstraint += (*electrons)[selElectron2].py();
284 idcs.push_back(selElectron1);
285 LV_l1.SetXYZT((*electrons)[selElectron1].px(),
286 (*electrons)[selElectron1].py(),
287 (*electrons)[selElectron1].pz(),
288 (*electrons)[selElectron1].
energy());
289 xconstraint += (*electrons)[selElectron1].px();
290 yconstraint += (*electrons)[selElectron1].py();
295 idcs.push_back(selMuon1);
296 LV_l2.SetXYZT((*muons)[selMuon1].px(),
297 (*muons)[selMuon1].py(),
298 (*muons)[selMuon1].pz(),
299 (*muons)[selMuon1].
energy());
300 xconstraint += (*muons)[selMuon1].px();
301 yconstraint += (*muons)[selMuon1].py();
306 idcs.push_back(selMuon1);
307 LV_l1.SetXYZT((*muons)[selMuon1].px(),
308 (*muons)[selMuon1].py(),
309 (*muons)[selMuon1].pz(),
310 (*muons)[selMuon1].
energy());
311 xconstraint += (*muons)[selMuon1].px();
312 yconstraint += (*muons)[selMuon1].py();
314 idcs.push_back(selElectron1);
315 LV_l2.SetXYZT((*electrons)[selElectron1].px(),
316 (*electrons)[selElectron1].py(),
317 (*electrons)[selElectron1].pz(),
318 (*electrons)[selElectron1].
energy());
319 xconstraint += (*electrons)[selElectron1].px();
320 yconstraint += (*electrons)[selElectron1].py();
327 idcs.push_back(selElectron1);
328 LV_l1.SetXYZT((*electrons)[selElectron1].px(),
329 (*electrons)[selElectron1].py(),
330 (*electrons)[selElectron1].pz(),
331 (*electrons)[selElectron1].
energy());
332 xconstraint += (*electrons)[selElectron1].px();
333 yconstraint += (*electrons)[selElectron1].py();
337 idcs.push_back(selMuon1);
338 LV_l2.SetXYZT((*muons)[selMuon1].px(),
339 (*muons)[selMuon1].py(),
340 (*muons)[selMuon1].pz(),
341 (*muons)[selMuon1].
energy());
342 xconstraint += (*muons)[selMuon1].px();
343 yconstraint += (*muons)[selMuon1].py();
350 idcs.push_back(selElectron1);
351 LV_l2.SetXYZT((*electrons)[selElectron1].px(),
352 (*electrons)[selElectron1].py(),
353 (*electrons)[selElectron1].pz(),
354 (*electrons)[selElectron1].
energy());
355 xconstraint += (*electrons)[selElectron1].px();
356 yconstraint += (*electrons)[selElectron1].py();
360 idcs.push_back(selMuon1);
361 LV_l1.SetXYZT((*muons)[selMuon1].px(),
362 (*muons)[selMuon1].py(),
363 (*muons)[selMuon1].pz(),
364 (*muons)[selMuon1].
energy());
365 xconstraint += (*muons)[selMuon1].px();
366 yconstraint += (*muons)[selMuon1].py();
375 idcs.push_back(selMuon1);
376 LV_l1.SetXYZT((*muons)[selMuon1].px(),
377 (*muons)[selMuon1].py(),
378 (*muons)[selMuon1].pz(),
379 (*muons)[selMuon1].
energy());
380 xconstraint += (*muons)[selMuon1].px();
381 yconstraint += (*muons)[selMuon1].py();
383 idcs.push_back(selMuon2);
384 LV_l2.SetXYZT((*muons)[selMuon2].px(),
385 (*muons)[selMuon2].py(),
386 (*muons)[selMuon2].pz(),
387 (*muons)[selMuon2].
energy());
388 xconstraint += (*muons)[selMuon2].px();
389 yconstraint += (*muons)[selMuon2].py();
392 xconstraint += (*jets)[
ib].px() + (*jets)[ibbar].px() + (*mets)[0].px();
393 yconstraint += (*jets)[
ib].py() + (*jets)[ibbar].py() + (*mets)[0].py();
402 idcsV.push_back(idcs);
409 weightsV.push_back(std::make_pair(nuSol.
weight, nSol));
417 if(weightsV.empty()){
419 std::vector<int> idcs;
420 for(
int i=0;
i<6; ++
i)
423 idcsV.push_back(idcs);
424 weightsV.push_back(std::make_pair(-1,0));
428 nuBarsV.push_back(nuBar);
432 int weightL = weightsV.
size();
433 int nuL = nusV.size();
434 int nuBarL = nuBarsV.size();
435 int idxL = idcsV.size();
437 if(weightL!=nuL || weightL!=nuBarL || weightL!=idxL){
438 edm::LogWarning(
"TtFullLepKinSolutionProducer") <<
"Output vectors are of different length:" 439 <<
"\n weight: " << weightL
441 <<
"\n nubar: " << nuBarL
442 <<
"\n idcs: " << idxL;
446 if(weightsV.size()>1){
447 sort(weightsV.begin(), weightsV.end(),
Compare());
454 for(
int i=0;
i<stop; ++
i){
455 pWeight->push_back(weightsV[
i].
first);
456 pNus ->push_back(nusV[weightsV[
i].
second]);
457 pNuBars->push_back(nuBarsV[weightsV[
i].second]);
458 pIdcs ->push_back(idcsV[weightsV[
i].second]);
T getParameter(std::string const &) const
bool operator()(std::pair< double, int > a, std::pair< double, int > b)
edm::EDGetTokenT< std::vector< pat::Muon > > muonsToken_
bool LepDiffCharge(const reco::Candidate *, const reco::Candidate *) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
NeutrinoSolution getNuSolution(const TLorentzVector &LV_l, const TLorentzVector &LV_l_, const TLorentzVector &LV_b, const TLorentzVector &LV_b_)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
bool HasPositiveCharge(const reco::Candidate *) const
~TtFullLepKinSolutionProducer() override
std::vector< double > nupars_
void SetConstraints(const double xx=0, const double yy=0)
U second(std::pair< T, U > const &p)
void produce(edm::Event &evt, const edm::EventSetup &iSetup) override
bool PTComp(const reco::Candidate *, const reco::Candidate *) const
reco::LeafCandidate neutrino
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
virtual double pt() const =0
transverse momentum
std::string jetCorrLevel_
TtFullLepKinSolutionProducer(const edm::ParameterSet &iConfig)
edm::EDGetTokenT< std::vector< pat::Electron > > electronsToken_
virtual int charge() const =0
electric charge
reco::LeafCandidate neutrinoBar
TtFullLepKinSolver * solver