1 #ifndef TtFullLepKinSolutionProducer_h
2 #define TtFullLepKinSolutionProducer_h
9 #include "TLorentzVector.h"
35 bool operator()(std::pair<double, int>
a, std::pair<double, int>
b) {
return a.first > b.first; }
53 return (l1->
pt() > l2->
pt());
83 produces<std::vector<std::vector<int> > >();
84 produces<std::vector<reco::LeafCandidate> >(
"fullLepNeutrinos");
85 produces<std::vector<reco::LeafCandidate> >(
"fullLepNeutrinoBars");
86 produces<std::vector<double> >(
"solWeight");
87 produces<bool>(
"isWrongCharge");
100 std::vector<std::vector<int> > idcsV;
101 std::vector<reco::LeafCandidate> nusV;
102 std::vector<reco::LeafCandidate> nuBarsV;
103 std::vector<std::pair<double, int> > weightsV;
106 std::unique_ptr<std::vector<std::vector<int> > > pIdcs(
new std::vector<std::vector<int> >);
107 std::unique_ptr<std::vector<reco::LeafCandidate> > pNus(
new std::vector<reco::LeafCandidate>);
108 std::unique_ptr<std::vector<reco::LeafCandidate> > pNuBars(
new std::vector<reco::LeafCandidate>);
109 std::unique_ptr<std::vector<double> > pWeight(
new std::vector<double>);
110 std::unique_ptr<bool> pWrongCharge(
new bool);
121 int selMuon1 = -1, selMuon2 = -1;
122 int selElectron1 = -1, selElectron2 = -1;
126 bool isWrongCharge =
false;
127 bool jetsFound =
false;
128 bool METFound =
false;
129 bool electronsFound =
false;
130 bool electronMuonFound =
false;
131 bool muonsFound =
false;
134 if (jets->size() >= 2) {
139 if (!mets->empty()) {
145 if (muons->size() + electrons->size() >= 2) {
147 if (electrons->empty())
149 else if (muons->empty())
151 else if (electrons->size() == 1) {
152 if (muons->size() == 1)
154 else if (
PTComp(&(*electrons)[0], &(*muons)[1]))
158 }
else if (electrons->size() > 1) {
159 if (
PTComp(&(*electrons)[1], &(*muons)[0]))
161 else if (muons->size() == 1)
163 else if (
PTComp(&(*electrons)[0], &(*muons)[1]))
177 electronsFound =
true;
179 isWrongCharge =
true;
185 electronMuonFound =
true;
187 isWrongCharge =
true;
200 isWrongCharge =
true;
205 *pWrongCharge = isWrongCharge;
208 if (
int(ee) + int(emu) + int(mumu) > 1) {
209 edm::LogWarning(
"TtFullLepKinSolutionProducer") <<
"Lepton selection criteria uncorrectly defined";
213 bool correctLeptons =
220 if (correctLeptons && METFound && jetsFound) {
225 if (jets->size() <
static_cast<unsigned int>(stop) || stop < 0)
232 for (
int ib = 0;
ib < stop;
ib++) {
234 for (
int ibbar = 0; ibbar < stop; ibbar++) {
239 std::vector<int> idcs;
243 idcs.push_back(ibbar);
245 TLorentzVector LV_l1;
246 TLorentzVector LV_l2;
247 TLorentzVector LV_b = TLorentzVector((*jets)[
ib].correctedJet(
jetCorrLevel_,
"bottom").px(),
251 TLorentzVector LV_bbar = TLorentzVector((*jets)[ibbar].correctedJet(
jetCorrLevel_,
"bottom").px(),
256 double xconstraint = 0, yconstraint = 0;
259 idcs.push_back(selElectron1);
260 LV_l1.SetXYZT((*electrons)[selElectron1].px(),
261 (*electrons)[selElectron1].py(),
262 (*electrons)[selElectron1].pz(),
263 (*electrons)[selElectron1].
energy());
264 xconstraint += (*electrons)[selElectron1].px();
265 yconstraint += (*electrons)[selElectron1].py();
267 idcs.push_back(selElectron2);
268 LV_l2.SetXYZT((*electrons)[selElectron2].px(),
269 (*electrons)[selElectron2].py(),
270 (*electrons)[selElectron2].pz(),
271 (*electrons)[selElectron2].
energy());
272 xconstraint += (*electrons)[selElectron2].px();
273 yconstraint += (*electrons)[selElectron2].py();
280 if (!isWrongCharge) {
282 idcs.push_back(selElectron1);
283 LV_l1.SetXYZT((*electrons)[selElectron1].px(),
284 (*electrons)[selElectron1].py(),
285 (*electrons)[selElectron1].pz(),
286 (*electrons)[selElectron1].
energy());
287 xconstraint += (*electrons)[selElectron1].px();
288 yconstraint += (*electrons)[selElectron1].py();
293 idcs.push_back(selMuon1);
294 LV_l2.SetXYZT((*muons)[selMuon1].px(),
295 (*muons)[selMuon1].py(),
296 (*muons)[selMuon1].pz(),
297 (*muons)[selMuon1].
energy());
298 xconstraint += (*muons)[selMuon1].px();
299 yconstraint += (*muons)[selMuon1].py();
303 idcs.push_back(selMuon1);
304 LV_l1.SetXYZT((*muons)[selMuon1].px(),
305 (*muons)[selMuon1].py(),
306 (*muons)[selMuon1].pz(),
307 (*muons)[selMuon1].
energy());
308 xconstraint += (*muons)[selMuon1].px();
309 yconstraint += (*muons)[selMuon1].py();
311 idcs.push_back(selElectron1);
312 LV_l2.SetXYZT((*electrons)[selElectron1].px(),
313 (*electrons)[selElectron1].py(),
314 (*electrons)[selElectron1].pz(),
315 (*electrons)[selElectron1].
energy());
316 xconstraint += (*electrons)[selElectron1].px();
317 yconstraint += (*electrons)[selElectron1].py();
323 idcs.push_back(selElectron1);
324 LV_l1.SetXYZT((*electrons)[selElectron1].px(),
325 (*electrons)[selElectron1].py(),
326 (*electrons)[selElectron1].pz(),
327 (*electrons)[selElectron1].
energy());
328 xconstraint += (*electrons)[selElectron1].px();
329 yconstraint += (*electrons)[selElectron1].py();
333 idcs.push_back(selMuon1);
334 LV_l2.SetXYZT((*muons)[selMuon1].px(),
335 (*muons)[selMuon1].py(),
336 (*muons)[selMuon1].pz(),
337 (*muons)[selMuon1].
energy());
338 xconstraint += (*muons)[selMuon1].px();
339 yconstraint += (*muons)[selMuon1].py();
345 idcs.push_back(selElectron1);
346 LV_l2.SetXYZT((*electrons)[selElectron1].px(),
347 (*electrons)[selElectron1].py(),
348 (*electrons)[selElectron1].pz(),
349 (*electrons)[selElectron1].
energy());
350 xconstraint += (*electrons)[selElectron1].px();
351 yconstraint += (*electrons)[selElectron1].py();
355 idcs.push_back(selMuon1);
356 LV_l1.SetXYZT((*muons)[selMuon1].px(),
357 (*muons)[selMuon1].py(),
358 (*muons)[selMuon1].pz(),
359 (*muons)[selMuon1].
energy());
360 xconstraint += (*muons)[selMuon1].px();
361 yconstraint += (*muons)[selMuon1].py();
370 idcs.push_back(selMuon1);
372 (*muons)[selMuon1].px(), (*muons)[selMuon1].py(), (*muons)[selMuon1].pz(), (*muons)[selMuon1].
energy());
373 xconstraint += (*muons)[selMuon1].px();
374 yconstraint += (*muons)[selMuon1].py();
376 idcs.push_back(selMuon2);
378 (*muons)[selMuon2].px(), (*muons)[selMuon2].py(), (*muons)[selMuon2].pz(), (*muons)[selMuon2].
energy());
379 xconstraint += (*muons)[selMuon2].px();
380 yconstraint += (*muons)[selMuon2].py();
383 xconstraint += (*jets)[
ib].px() + (*jets)[ibbar].px() + (*mets)[0].px();
384 yconstraint += (*jets)[
ib].py() + (*jets)[ibbar].py() + (*mets)[0].py();
393 idcsV.push_back(idcs);
400 weightsV.push_back(std::make_pair(nuSol.
weight, nSol));
408 if (weightsV.empty()) {
410 std::vector<int> idcs;
412 for (
int i = 0;
i < 6; ++
i)
415 idcsV.push_back(idcs);
416 weightsV.push_back(std::make_pair(-1, 0));
420 nuBarsV.push_back(nuBar);
424 int weightL = weightsV.
size();
425 int nuL = nusV.size();
426 int nuBarL = nuBarsV.size();
427 int idxL = idcsV.size();
429 if (weightL != nuL || weightL != nuBarL || weightL != idxL) {
431 <<
"Output vectors are of different length:"
432 <<
"\n weight: " << weightL <<
"\n nu: " << nuL <<
"\n nubar: " << nuBarL <<
"\n idcs: " << idxL;
436 if (weightsV.size() > 1) {
437 sort(weightsV.begin(), weightsV.end(),
Compare());
445 for (
int i = 0;
i < stop; ++
i) {
446 pWeight->push_back(weightsV[
i].
first);
447 pNus->push_back(nusV[weightsV[
i].
second]);
448 pNuBars->push_back(nuBarsV[weightsV[
i].second]);
449 pIdcs->push_back(idcsV[weightsV[
i].second]);
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
virtual double pt() const =0
transverse momentum
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
virtual int charge() const =0
electric charge
edm::EDGetTokenT< std::vector< pat::MET > > metsToken_
edm::EDGetTokenT< std::vector< pat::Jet > > jetsToken_
std::string jetCorrLevel_
T getParameter(std::string const &) const
TtFullLepKinSolutionProducer(const edm::ParameterSet &iConfig)
edm::EDGetTokenT< std::vector< pat::Electron > > electronsToken_
Log< level::Warning, false > LogWarning
reco::LeafCandidate neutrinoBar
TtFullLepKinSolver * solver