File indexing completed on 2023-10-25 09:54:51
0001 #ifndef L1Trigger_L1TCalorimeter_L1TStage2CaloLayer2Comp
0002 #define L1Trigger_L1TCalorimeter_L1TStage2CaloLayer2Comp
0003
0004 #include "FWCore/Framework/interface/MakerMacros.h"
0005
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009
0010 #include "DataFormats/L1Trigger/interface/EGamma.h"
0011 #include "DataFormats/L1Trigger/interface/Jet.h"
0012 #include "DataFormats/L1Trigger/interface/EtSum.h"
0013 #include "DataFormats/L1Trigger/interface/Tau.h"
0014
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/MakerMacros.h"
0017 #include "FWCore/Framework/interface/stream/EDProducer.h"
0018
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020
0021 #include "CondFormats/L1TObjects/interface/CaloParams.h"
0022 #include "CondFormats/DataRecord/interface/L1TCaloParamsRcd.h"
0023
0024 #include "DataFormats/L1TCalorimeter/interface/CaloTower.h"
0025 #include "DataFormats/L1TCalorimeter/interface/CaloCluster.h"
0026 #include "DataFormats/L1Trigger/interface/EGamma.h"
0027 #include "DataFormats/L1Trigger/interface/Tau.h"
0028 #include "DataFormats/L1Trigger/interface/Jet.h"
0029 #include "DataFormats/L1Trigger/interface/EtSum.h"
0030
0031
0032 #include "L1Trigger/L1TCalorimeter/interface/CaloTools.h"
0033
0034
0035 #include "L1Trigger/L1TCalorimeter/interface/Stage2Layer2JetAlgorithmFirmware.h"
0036 #include "L1Trigger/L1TCalorimeter/interface/AccumulatingSort.h"
0037 #include "L1Trigger/L1TCalorimeter/interface/BitonicSort.h"
0038 #include "DataFormats/Math/interface/LorentzVector.h"
0039
0040 namespace l1t {
0041 bool operator>(const l1t::Jet &a, l1t::Jet &b) { return a.hwPt() > b.hwPt(); }
0042 bool operator>(const l1t::EGamma &a, l1t::EGamma &b) { return a.hwPt() > b.hwPt(); }
0043 bool operator>(const l1t::Tau &a, l1t::Tau &b) { return a.hwPt() > b.hwPt(); }
0044 }
0045
0046
0047
0048
0049
0050
0051
0052 class L1TStage2CaloLayer2Comp : public edm::stream::EDProducer<> {
0053 public:
0054
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064
0065
0066
0067 L1TStage2CaloLayer2Comp(const edm::ParameterSet &ps);
0068
0069 protected:
0070
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082 void produce(edm::Event &, const edm::EventSetup &) override;
0083
0084 private:
0085
0086
0087
0088
0089
0090
0091
0092
0093
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105
0106
0107 bool compareJets(const edm::Handle<l1t::JetBxCollection> &dataCol, const edm::Handle<l1t::JetBxCollection> &emulCol);
0108
0109
0110
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122
0123
0124
0125
0126
0127
0128
0129
0130
0131
0132
0133 bool compareEGs(const edm::Handle<l1t::EGammaBxCollection> &dataCol,
0134 const edm::Handle<l1t::EGammaBxCollection> &emulCol);
0135
0136
0137
0138
0139
0140
0141
0142
0143
0144
0145
0146
0147
0148
0149
0150
0151
0152
0153
0154
0155
0156
0157
0158
0159
0160 bool compareTaus(const edm::Handle<l1t::TauBxCollection> &dataCol, const edm::Handle<l1t::TauBxCollection> &emulCol);
0161
0162
0163
0164
0165
0166
0167
0168
0169
0170
0171
0172
0173
0174
0175
0176
0177
0178
0179
0180
0181 bool compareSums(const edm::Handle<l1t::EtSumBxCollection> &dataCol,
0182 const edm::Handle<l1t::EtSumBxCollection> &emulCol);
0183
0184 void accuSort(std::vector<l1t::Jet> *jets);
0185 void accuSort(std::vector<l1t::Tau> *taus);
0186 void accuSort(std::vector<l1t::EGamma> *egs);
0187 void accuSort(std::vector<l1t::L1Candidate> &jets);
0188
0189 void dumpEventToFile();
0190 void dumpEventToEDM(edm::Event &e);
0191
0192
0193 edm::EDGetToken calol2JetCollectionData;
0194 edm::EDGetToken calol2JetCollectionEmul;
0195 edm::EDGetToken calol2EGammaCollectionData;
0196 edm::EDGetToken calol2EGammaCollectionEmul;
0197 edm::EDGetToken calol2TauCollectionData;
0198 edm::EDGetToken calol2TauCollectionEmul;
0199 edm::EDGetToken calol2EtSumCollectionData;
0200 edm::EDGetToken calol2EtSumCollectionEmul;
0201 edm::EDGetToken calol2CaloTowerCollectionData;
0202 edm::EDGetToken calol2CaloTowerCollectionEmul;
0203
0204
0205 edm::Handle<l1t::JetBxCollection> jetDataCol;
0206 edm::Handle<l1t::JetBxCollection> jetEmulCol;
0207 edm::Handle<l1t::EGammaBxCollection> egDataCol;
0208 edm::Handle<l1t::EGammaBxCollection> egEmulCol;
0209 edm::Handle<l1t::TauBxCollection> tauDataCol;
0210 edm::Handle<l1t::TauBxCollection> tauEmulCol;
0211 edm::Handle<l1t::EtSumBxCollection> sumDataCol;
0212 edm::Handle<l1t::EtSumBxCollection> sumEmulCol;
0213 edm::Handle<l1t::CaloTowerBxCollection> caloTowerDataCol;
0214 edm::Handle<l1t::CaloTowerBxCollection> caloTowerEmulCol;
0215
0216 const unsigned int currBx = 0;
0217
0218 bool dumpTowers = false;
0219 bool dumpWholeEvent = false;
0220 };
0221
0222 L1TStage2CaloLayer2Comp::L1TStage2CaloLayer2Comp(const edm::ParameterSet &ps)
0223 : calol2JetCollectionData(
0224 consumes<l1t::JetBxCollection>(ps.getParameter<edm::InputTag>("calol2JetCollectionData"))),
0225 calol2JetCollectionEmul(
0226 consumes<l1t::JetBxCollection>(ps.getParameter<edm::InputTag>("calol2JetCollectionEmul"))),
0227 calol2EGammaCollectionData(
0228 consumes<l1t::EGammaBxCollection>(ps.getParameter<edm::InputTag>("calol2EGammaCollectionData"))),
0229 calol2EGammaCollectionEmul(
0230 consumes<l1t::EGammaBxCollection>(ps.getParameter<edm::InputTag>("calol2EGammaCollectionEmul"))),
0231 calol2TauCollectionData(
0232 consumes<l1t::TauBxCollection>(ps.getParameter<edm::InputTag>("calol2TauCollectionData"))),
0233 calol2TauCollectionEmul(
0234 consumes<l1t::TauBxCollection>(ps.getParameter<edm::InputTag>("calol2TauCollectionEmul"))),
0235 calol2EtSumCollectionData(
0236 consumes<l1t::EtSumBxCollection>(ps.getParameter<edm::InputTag>("calol2EtSumCollectionData"))),
0237 calol2EtSumCollectionEmul(
0238 consumes<l1t::EtSumBxCollection>(ps.getParameter<edm::InputTag>("calol2EtSumCollectionEmul"))),
0239 calol2CaloTowerCollectionData(
0240 consumes<l1t::CaloTowerBxCollection>(ps.getParameter<edm::InputTag>("calol2CaloTowerCollectionData"))),
0241 calol2CaloTowerCollectionEmul(
0242 consumes<l1t::CaloTowerBxCollection>(ps.getParameter<edm::InputTag>("calol2CaloTowerCollectionEmul"))) {
0243 produces<l1t::JetBxCollection>("dataJet").setBranchAlias("dataJets");
0244 produces<l1t::JetBxCollection>("emulJet").setBranchAlias("emulJets");
0245 produces<l1t::EGammaBxCollection>("dataEg").setBranchAlias("dataEgs");
0246 produces<l1t::EGammaBxCollection>("emulEg").setBranchAlias("emulEgs");
0247 produces<l1t::TauBxCollection>("dataTau").setBranchAlias("dataTaus");
0248 produces<l1t::TauBxCollection>("emulTau").setBranchAlias("emulTaus");
0249 produces<l1t::EtSumBxCollection>("dataEtSum").setBranchAlias("dataEtSums");
0250 produces<l1t::EtSumBxCollection>("emulEtSum").setBranchAlias("emulEtSums");
0251 produces<l1t::CaloTowerBxCollection>("dataCaloTower").setBranchAlias("dataCaloTowers");
0252 produces<l1t::CaloTowerBxCollection>("emulCaloTower").setBranchAlias("emulCaloTowers");
0253
0254 dumpTowers = ps.getParameter<bool>("dumpTowers");
0255 dumpWholeEvent = ps.getParameter<bool>("dumpWholeEvent");
0256 }
0257
0258 void L1TStage2CaloLayer2Comp::produce(edm::Event &e, const edm::EventSetup &c) {
0259 bool eventGood = true;
0260
0261
0262 e.getByToken(calol2JetCollectionData, jetDataCol);
0263 e.getByToken(calol2JetCollectionEmul, jetEmulCol);
0264 e.getByToken(calol2EGammaCollectionData, egDataCol);
0265 e.getByToken(calol2EGammaCollectionEmul, egEmulCol);
0266 e.getByToken(calol2TauCollectionData, tauDataCol);
0267 e.getByToken(calol2TauCollectionEmul, tauEmulCol);
0268 e.getByToken(calol2EtSumCollectionData, sumDataCol);
0269 e.getByToken(calol2EtSumCollectionEmul, sumEmulCol);
0270 e.getByToken(calol2CaloTowerCollectionData, caloTowerDataCol);
0271 e.getByToken(calol2CaloTowerCollectionEmul, caloTowerEmulCol);
0272
0273 edm::LogProblem("l1tcalol2ebec") << "Processing event " << e.id() << std::endl;
0274
0275 if (caloTowerDataCol->size() == 0) {
0276 edm::LogProblem("l1tcalol2ebec") << "Empty caloTowers. Skipping event " << std::endl;
0277 return;
0278 }
0279
0280 if (!compareJets(jetDataCol, jetEmulCol)) {
0281 eventGood = false;
0282 }
0283
0284 if (!compareEGs(egDataCol, egEmulCol)) {
0285 eventGood = false;
0286 }
0287
0288 if (!compareTaus(tauDataCol, tauEmulCol)) {
0289 eventGood = false;
0290 }
0291
0292 if (!compareSums(sumDataCol, sumEmulCol)) {
0293 eventGood = false;
0294 }
0295
0296 if (!eventGood) {
0297 if (dumpWholeEvent) {
0298
0299 dumpEventToFile();
0300 }
0301
0302
0303 dumpEventToEDM(e);
0304
0305 edm::LogProblem("l1tcalol2ebec") << "Bad event! " << std::endl;
0306 }
0307 }
0308
0309
0310 bool L1TStage2CaloLayer2Comp::compareJets(const edm::Handle<l1t::JetBxCollection> &dataCol,
0311 const edm::Handle<l1t::JetBxCollection> &emulCol) {
0312 bool eventGood = true;
0313
0314 std::vector<l1t::Jet> *jets = new std::vector<l1t::Jet>;
0315 std::vector<l1t::Jet>::iterator dataIt;
0316 l1t::JetBxCollection::const_iterator dataBxIt = dataCol->begin(currBx);
0317 l1t::JetBxCollection::const_iterator emulBxIt = emulCol->begin(currBx);
0318
0319 if (dataCol->size(currBx) > 0) {
0320
0321 while (true) {
0322 jets->emplace_back(*dataBxIt);
0323
0324 dataBxIt++;
0325 if (dataBxIt == dataCol->end(currBx))
0326 break;
0327 }
0328
0329 accuSort(jets);
0330 dataIt = jets->begin();
0331 }
0332
0333
0334 if (jets->size() != emulCol->size(currBx)) {
0335 edm::LogProblem("l1tcalol2ebec") << "Jet collection size difference: "
0336 << "data = " << jets->size() << " emul = " << emulCol->size(currBx) << std::endl;
0337 return false;
0338 }
0339
0340 if (dataIt != jets->end() || emulBxIt != emulCol->end(currBx)) {
0341 int nNeg = 0;
0342
0343 while (true) {
0344 bool posGood = true;
0345 bool etGood = true;
0346 bool qualGood = true;
0347
0348
0349 if (dataIt->hwPt() != emulBxIt->hwPt()) {
0350 etGood = false;
0351 eventGood = false;
0352 }
0353
0354
0355 if (dataIt->hwPhi() != emulBxIt->hwPhi()) {
0356 posGood = false;
0357 }
0358
0359
0360 if (dataIt->hwEta() != emulBxIt->hwEta()) {
0361 posGood = false;
0362 }
0363
0364
0365 if (dataIt->hwQual() != emulBxIt->hwQual()) {
0366 qualGood = false;
0367 eventGood = false;
0368 }
0369
0370
0371 if (etGood && !posGood) {
0372 l1t::JetBxCollection::const_iterator emulItCheckSort;
0373 l1t::JetBxCollection::const_iterator dataItCheckSort;
0374 for (emulItCheckSort = emulCol->begin(currBx); emulItCheckSort != emulCol->end(currBx); ++emulItCheckSort) {
0375 for (dataItCheckSort = jets->begin(); dataItCheckSort != jets->end(); ++dataItCheckSort) {
0376 if (dataItCheckSort->hwEta() < 0)
0377 ++nNeg;
0378
0379 if (dataIt->hwPt() == emulItCheckSort->hwPt() && dataIt->hwPhi() == emulItCheckSort->hwPhi() &&
0380 dataIt->hwEta() == emulItCheckSort->hwEta())
0381 posGood = true;
0382 }
0383 }
0384 }
0385
0386 if (etGood && dataIt->hwEta() > 0 && ((distance(dataIt, jets->end()) - nNeg) < 5))
0387 posGood = true;
0388 if (etGood && dataIt->hwEta() < 0 && (distance(dataIt, jets->end()) < 5))
0389 posGood = true;
0390
0391 if (!posGood)
0392 eventGood = false;
0393
0394
0395 if (!etGood || !posGood || !qualGood) {
0396 edm::LogProblem("l1tcalol2ebec") << "Jet Problem (data emul): "
0397 << "\tEt = " << dataIt->hwPt() << " " << emulBxIt->hwPt()
0398 << "\teta = " << dataIt->hwEta() << " " << emulBxIt->hwEta()
0399 << "\tphi = " << dataIt->hwPhi() << " " << emulBxIt->hwPhi()
0400 << "\tqual = " << dataIt->hwQual() << " " << emulBxIt->hwQual() << std::endl;
0401 }
0402
0403
0404 ++dataIt;
0405 ++emulBxIt;
0406
0407 if (dataIt == jets->end() || emulBxIt == emulCol->end(currBx))
0408 break;
0409 }
0410 } else {
0411 if (!jets->empty() || emulCol->size(currBx) != 0)
0412 return false;
0413 }
0414
0415
0416
0417 return eventGood;
0418 }
0419
0420
0421 bool L1TStage2CaloLayer2Comp::compareEGs(const edm::Handle<l1t::EGammaBxCollection> &dataCol,
0422 const edm::Handle<l1t::EGammaBxCollection> &emulCol) {
0423 bool eventGood = true;
0424
0425 std::vector<l1t::EGamma> *egs = new std::vector<l1t::EGamma>;
0426 std::vector<l1t::EGamma>::iterator dataIt;
0427 l1t::EGammaBxCollection::const_iterator dataBxIt = dataCol->begin(currBx);
0428 l1t::EGammaBxCollection::const_iterator emulBxIt = emulCol->begin(currBx);
0429
0430 if (dataCol->size(currBx) > 0) {
0431
0432 while (true) {
0433 egs->emplace_back(*dataBxIt);
0434
0435 dataBxIt++;
0436 if (dataBxIt == dataCol->end(currBx))
0437 break;
0438 }
0439 accuSort(egs);
0440 dataIt = egs->begin();
0441 }
0442
0443
0444 if (egs->size() != emulCol->size(currBx)) {
0445 edm::LogProblem("l1tcalol2ebec") << "EG collection size difference: "
0446 << "data = " << egs->size() << " emul = " << emulCol->size(currBx) << std::endl;
0447 return false;
0448 }
0449
0450
0451 if (dataIt != egs->end() || emulBxIt != emulCol->end(currBx)) {
0452 int nNeg = 0;
0453
0454 while (true) {
0455 bool posGood = true;
0456 bool etGood = true;
0457
0458
0459 if (dataIt->hwPt() != emulBxIt->hwPt()) {
0460 etGood = false;
0461 eventGood = false;
0462 }
0463
0464
0465 if (dataIt->hwPhi() != emulBxIt->hwPhi()) {
0466 posGood = false;
0467 }
0468
0469
0470 if (dataIt->hwEta() != emulBxIt->hwEta()) {
0471 posGood = false;
0472 }
0473
0474
0475 if (etGood && !posGood) {
0476 l1t::EGammaBxCollection::const_iterator emulItCheckSort;
0477 l1t::EGammaBxCollection::const_iterator dataItCheckSort;
0478 for (emulItCheckSort = emulCol->begin(currBx); emulItCheckSort != emulCol->end(currBx); ++emulItCheckSort) {
0479 for (dataItCheckSort = egs->begin(); dataItCheckSort != egs->end(); ++dataItCheckSort) {
0480 if (dataItCheckSort->hwEta() < 0)
0481 ++nNeg;
0482
0483 if (dataIt->hwPt() == emulItCheckSort->hwPt() && dataIt->hwPhi() == emulItCheckSort->hwPhi() &&
0484 dataIt->hwEta() == emulItCheckSort->hwEta())
0485 posGood = true;
0486 }
0487 }
0488 }
0489
0490 if (etGood && dataIt->hwEta() > 0 && ((distance(dataIt, egs->end()) - nNeg) < 5))
0491 posGood = true;
0492 if (etGood && dataIt->hwEta() < 0 && (distance(dataIt, egs->end()) < 5))
0493 posGood = true;
0494
0495 if (!posGood)
0496 eventGood = false;
0497
0498
0499 if (!posGood || !etGood) {
0500 edm::LogProblem("l1tcalol2ebec") << "EG Problem (data emul): "
0501 << "\tEt = " << dataIt->hwPt() << " " << emulBxIt->hwPt()
0502 << "\teta = " << dataIt->hwEta() << " " << emulBxIt->hwEta()
0503 << "\tphi = " << dataIt->hwPhi() << " " << emulBxIt->hwPhi() << std::endl;
0504 }
0505
0506
0507 ++dataIt;
0508 ++emulBxIt;
0509
0510 if (dataIt == egs->end() || emulBxIt == emulCol->end(currBx))
0511 break;
0512 }
0513 } else {
0514 if (!egs->empty() || emulCol->size(currBx) != 0)
0515 return false;
0516 }
0517
0518
0519
0520 return eventGood;
0521 }
0522
0523
0524 bool L1TStage2CaloLayer2Comp::compareTaus(const edm::Handle<l1t::TauBxCollection> &dataCol,
0525 const edm::Handle<l1t::TauBxCollection> &emulCol) {
0526 bool eventGood = true;
0527
0528 std::vector<l1t::Tau> *taus = new std::vector<l1t::Tau>;
0529 std::vector<l1t::Tau>::iterator dataIt;
0530 l1t::TauBxCollection::const_iterator dataBxIt = dataCol->begin(currBx);
0531 l1t::TauBxCollection::const_iterator emulBxIt = emulCol->begin(currBx);
0532
0533 if (dataCol->size(currBx) > 0) {
0534
0535 while (true) {
0536 taus->emplace_back(*dataBxIt);
0537
0538 dataBxIt++;
0539 if (dataBxIt == dataCol->end(currBx))
0540 break;
0541 }
0542 accuSort(taus);
0543 dataIt = taus->begin();
0544 }
0545
0546
0547 if (taus->size() != emulCol->size(currBx)) {
0548 edm::LogProblem("l1tcalol2ebec") << "Tau collection size difference: "
0549 << "data = " << taus->size() << " emul = " << emulCol->size(currBx) << std::endl;
0550 return false;
0551 }
0552
0553
0554 if (dataIt != taus->end() || emulBxIt != emulCol->end(currBx)) {
0555 int nNeg = 0;
0556
0557 while (true) {
0558 bool posGood = true;
0559 bool etGood = true;
0560
0561
0562 if (dataIt->hwPt() != emulBxIt->hwPt()) {
0563 etGood = false;
0564 eventGood = false;
0565 }
0566
0567
0568 if (dataIt->hwPhi() != emulBxIt->hwPhi()) {
0569 posGood = false;
0570 }
0571
0572
0573 if (dataIt->hwEta() != emulBxIt->hwEta()) {
0574 posGood = false;
0575 }
0576
0577
0578 if (etGood && !posGood) {
0579 l1t::TauBxCollection::const_iterator emulItCheckSort;
0580 l1t::TauBxCollection::const_iterator dataItCheckSort;
0581 for (emulItCheckSort = emulCol->begin(currBx); emulItCheckSort != emulCol->end(currBx); ++emulItCheckSort) {
0582 for (dataItCheckSort = taus->begin(); dataItCheckSort != taus->end(); ++dataItCheckSort) {
0583 if (dataItCheckSort->hwEta() < 0)
0584 ++nNeg;
0585
0586 if (dataIt->hwPt() == emulItCheckSort->hwPt() && dataIt->hwPhi() == emulItCheckSort->hwPhi() &&
0587 dataIt->hwEta() == emulItCheckSort->hwEta())
0588 posGood = true;
0589 }
0590 }
0591 }
0592
0593 if (etGood && dataIt->hwEta() > 0 && ((distance(dataIt, taus->end()) - nNeg) < 5))
0594 posGood = true;
0595 if (etGood && dataIt->hwEta() < 0 && (distance(dataIt, taus->end()) < 5))
0596 posGood = true;
0597
0598 if (!posGood)
0599 eventGood = false;
0600
0601
0602 if (!posGood || !etGood) {
0603 edm::LogProblem("l1tcalol2ebec") << "Tau Problem (data emul): "
0604 << "\tEt = " << dataIt->hwPt() << " " << emulBxIt->hwPt()
0605 << "\teta = " << dataIt->hwEta() << " " << emulBxIt->hwEta()
0606 << "\tphi = " << dataIt->hwPhi() << " " << emulBxIt->hwPhi() << std::endl;
0607 }
0608
0609
0610 ++dataIt;
0611 ++emulBxIt;
0612
0613 if (dataIt == taus->end() || emulBxIt == emulCol->end(currBx))
0614 break;
0615 }
0616 } else {
0617 if (!taus->empty() || emulCol->size(currBx) != 0)
0618 return false;
0619 }
0620
0621
0622
0623 return eventGood;
0624 }
0625
0626
0627 bool L1TStage2CaloLayer2Comp::compareSums(const edm::Handle<l1t::EtSumBxCollection> &dataCol,
0628 const edm::Handle<l1t::EtSumBxCollection> &emulCol) {
0629 bool eventGood = true;
0630
0631 double dataEt = 0;
0632 double emulEt = 0;
0633
0634
0635
0636 if ((dataCol->isEmpty(currBx)) || (emulCol->isEmpty(currBx)) || (dataCol->size(currBx) != emulCol->size(currBx))) {
0637 edm::LogProblem("l1tcalol2ebec") << "EtSum collection size difference: "
0638 << "data = " << dataCol->size(currBx) << " emul = " << emulCol->size(currBx)
0639 << std::endl;
0640
0641 return false;
0642 }
0643
0644 for (unsigned int i = 0; i < emulCol->size(currBx); i++) {
0645 l1t::EtSum const &emulSum = emulCol->at(currBx, i);
0646 l1t::EtSum const &dataSum = dataCol->at(currBx, l1t::CaloTools::emul_to_data_sum_index_map[i]);
0647
0648 if (emulSum.getType() != dataSum.getType()) {
0649 edm::LogProblem("l1tcalol2ebec") << "EtSum type problem (data emul): " << dataSum.getType() << " "
0650 << emulSum.getType() << std::endl;
0651 }
0652
0653 dataEt = dataSum.hwPt();
0654 emulEt = emulSum.hwPt();
0655
0656 if (dataEt != emulEt) {
0657 eventGood = false;
0658 edm::LogProblem("l1tcalol2ebec") << "EtSum problem (data emul):\tType = " << emulSum.getType()
0659 << "\tEt = " << dataEt << " " << emulEt << std::endl;
0660 }
0661 }
0662
0663
0664
0665 return eventGood;
0666 }
0667
0668
0669 void L1TStage2CaloLayer2Comp::accuSort(std::vector<l1t::Jet> *jets) {
0670 math::PtEtaPhiMLorentzVector emptyP4;
0671 l1t::Jet tempJet(emptyP4, 0, 0, 0, 0);
0672 std::vector<std::vector<l1t::Jet> > jetEtaPos(41, std::vector<l1t::Jet>(18, tempJet));
0673 std::vector<std::vector<l1t::Jet> > jetEtaNeg(41, std::vector<l1t::Jet>(18, tempJet));
0674
0675 for (unsigned int iJet = 0; iJet < jets->size(); iJet++) {
0676 if (jets->at(iJet).hwEta() > 0)
0677 jetEtaPos.at(jets->at(iJet).hwEta() - 1).at((jets->at(iJet).hwPhi() - 1) / 4) = jets->at(iJet);
0678 else
0679 jetEtaNeg.at(-(jets->at(iJet).hwEta() + 1)).at((jets->at(iJet).hwPhi() - 1) / 4) = jets->at(iJet);
0680 }
0681
0682 AccumulatingSort<l1t::Jet> etaPosSorter(7);
0683 AccumulatingSort<l1t::Jet> etaNegSorter(7);
0684 std::vector<l1t::Jet> accumEtaPos;
0685 std::vector<l1t::Jet> accumEtaNeg;
0686
0687 for (int ieta = 0; ieta < 41; ++ieta) {
0688
0689 std::vector<l1t::Jet>::iterator start_, end_;
0690 start_ = jetEtaPos.at(ieta).begin();
0691 end_ = jetEtaPos.at(ieta).end();
0692 BitonicSort<l1t::Jet>(down, start_, end_);
0693 etaPosSorter.Merge(jetEtaPos.at(ieta), accumEtaPos);
0694
0695
0696 start_ = jetEtaNeg.at(ieta).begin();
0697 end_ = jetEtaNeg.at(ieta).end();
0698 BitonicSort<l1t::Jet>(down, start_, end_);
0699 etaNegSorter.Merge(jetEtaNeg.at(ieta), accumEtaNeg);
0700 }
0701
0702
0703 if (accumEtaPos.at(6).hwPt() == accumEtaPos.at(5).hwPt() && accumEtaPos.at(6).hwEta() == accumEtaPos.at(5).hwEta() &&
0704 accumEtaPos.at(6).hwPhi() > accumEtaPos.at(5).hwPhi()) {
0705 accumEtaPos.at(5) = accumEtaPos.at(6);
0706 }
0707 if (accumEtaNeg.at(6).hwPt() == accumEtaNeg.at(5).hwPt() && accumEtaNeg.at(6).hwEta() == accumEtaNeg.at(5).hwEta() &&
0708 accumEtaNeg.at(6).hwPhi() > accumEtaNeg.at(5).hwPhi()) {
0709 accumEtaNeg.at(5) = accumEtaNeg.at(6);
0710 }
0711
0712
0713 accumEtaPos.resize(6);
0714 accumEtaNeg.resize(6);
0715
0716
0717 jets->clear();
0718 for (const l1t::Jet &accjet : accumEtaPos) {
0719 if (accjet.hwPt() > 0)
0720 jets->push_back(accjet);
0721 }
0722 for (const l1t::Jet &accjet : accumEtaNeg) {
0723 if (accjet.hwPt() > 0)
0724 jets->push_back(accjet);
0725 }
0726 }
0727
0728
0729 void L1TStage2CaloLayer2Comp::accuSort(std::vector<l1t::EGamma> *egs) {
0730 math::PtEtaPhiMLorentzVector emptyP4;
0731 l1t::EGamma tempEGamma(emptyP4, 0, 0, 0, 0);
0732 std::vector<std::vector<l1t::EGamma> > jetEtaPos(41, std::vector<l1t::EGamma>(18, tempEGamma));
0733 std::vector<std::vector<l1t::EGamma> > jetEtaNeg(41, std::vector<l1t::EGamma>(18, tempEGamma));
0734
0735 for (unsigned int iEGamma = 0; iEGamma < egs->size(); iEGamma++) {
0736 if (egs->at(iEGamma).hwEta() > 0)
0737 jetEtaPos.at(egs->at(iEGamma).hwEta() - 1).at((egs->at(iEGamma).hwPhi() - 1) / 4) = egs->at(iEGamma);
0738 else
0739 jetEtaNeg.at(-(egs->at(iEGamma).hwEta() + 1)).at((egs->at(iEGamma).hwPhi() - 1) / 4) = egs->at(iEGamma);
0740 }
0741
0742 AccumulatingSort<l1t::EGamma> etaPosSorter(7);
0743 AccumulatingSort<l1t::EGamma> etaNegSorter(7);
0744 std::vector<l1t::EGamma> accumEtaPos;
0745 std::vector<l1t::EGamma> accumEtaNeg;
0746
0747 for (int ieta = 0; ieta < 41; ++ieta) {
0748
0749 std::vector<l1t::EGamma>::iterator start_, end_;
0750 start_ = jetEtaPos.at(ieta).begin();
0751 end_ = jetEtaPos.at(ieta).end();
0752 BitonicSort<l1t::EGamma>(down, start_, end_);
0753 etaPosSorter.Merge(jetEtaPos.at(ieta), accumEtaPos);
0754
0755
0756 start_ = jetEtaNeg.at(ieta).begin();
0757 end_ = jetEtaNeg.at(ieta).end();
0758 BitonicSort<l1t::EGamma>(down, start_, end_);
0759 etaNegSorter.Merge(jetEtaNeg.at(ieta), accumEtaNeg);
0760 }
0761
0762
0763 if (accumEtaPos.at(6).hwPt() == accumEtaPos.at(5).hwPt() && accumEtaPos.at(6).hwEta() == accumEtaPos.at(5).hwEta() &&
0764 accumEtaPos.at(6).hwPhi() > accumEtaPos.at(5).hwPhi()) {
0765 accumEtaPos.at(5) = accumEtaPos.at(6);
0766 }
0767 if (accumEtaNeg.at(6).hwPt() == accumEtaNeg.at(5).hwPt() && accumEtaNeg.at(6).hwEta() == accumEtaNeg.at(5).hwEta() &&
0768 accumEtaNeg.at(6).hwPhi() > accumEtaNeg.at(5).hwPhi()) {
0769 accumEtaNeg.at(5) = accumEtaNeg.at(6);
0770 }
0771
0772
0773 accumEtaPos.resize(6);
0774 accumEtaNeg.resize(6);
0775
0776
0777 egs->clear();
0778 for (const l1t::EGamma &accjet : accumEtaPos) {
0779 if (accjet.hwPt() > 0)
0780 egs->push_back(accjet);
0781 }
0782 for (const l1t::EGamma &accjet : accumEtaNeg) {
0783 if (accjet.hwPt() > 0)
0784 egs->push_back(accjet);
0785 }
0786 }
0787
0788
0789 void L1TStage2CaloLayer2Comp::accuSort(std::vector<l1t::Tau> *taus) {
0790 math::PtEtaPhiMLorentzVector emptyP4;
0791 l1t::Tau tempTau(emptyP4, 0, 0, 0, 0);
0792 std::vector<std::vector<l1t::Tau> > jetEtaPos(41, std::vector<l1t::Tau>(18, tempTau));
0793 std::vector<std::vector<l1t::Tau> > jetEtaNeg(41, std::vector<l1t::Tau>(18, tempTau));
0794
0795 for (unsigned int iTau = 0; iTau < taus->size(); iTau++) {
0796 if (taus->at(iTau).hwEta() > 0)
0797 jetEtaPos.at(taus->at(iTau).hwEta() - 1).at((taus->at(iTau).hwPhi() - 1) / 4) = taus->at(iTau);
0798 else
0799 jetEtaNeg.at(-(taus->at(iTau).hwEta() + 1)).at((taus->at(iTau).hwPhi() - 1) / 4) = taus->at(iTau);
0800 }
0801
0802 AccumulatingSort<l1t::Tau> etaPosSorter(7);
0803 AccumulatingSort<l1t::Tau> etaNegSorter(7);
0804 std::vector<l1t::Tau> accumEtaPos;
0805 std::vector<l1t::Tau> accumEtaNeg;
0806
0807 for (int ieta = 0; ieta < 41; ++ieta) {
0808
0809 std::vector<l1t::Tau>::iterator start_, end_;
0810 start_ = jetEtaPos.at(ieta).begin();
0811 end_ = jetEtaPos.at(ieta).end();
0812 BitonicSort<l1t::Tau>(down, start_, end_);
0813 etaPosSorter.Merge(jetEtaPos.at(ieta), accumEtaPos);
0814
0815
0816 start_ = jetEtaNeg.at(ieta).begin();
0817 end_ = jetEtaNeg.at(ieta).end();
0818 BitonicSort<l1t::Tau>(down, start_, end_);
0819 etaNegSorter.Merge(jetEtaNeg.at(ieta), accumEtaNeg);
0820 }
0821
0822
0823 if (accumEtaPos.at(6).hwPt() == accumEtaPos.at(5).hwPt() && accumEtaPos.at(6).hwEta() == accumEtaPos.at(5).hwEta() &&
0824 accumEtaPos.at(6).hwPhi() > accumEtaPos.at(5).hwPhi()) {
0825 accumEtaPos.at(5) = accumEtaPos.at(6);
0826 }
0827 if (accumEtaNeg.at(6).hwPt() == accumEtaNeg.at(5).hwPt() && accumEtaNeg.at(6).hwEta() == accumEtaNeg.at(5).hwEta() &&
0828 accumEtaNeg.at(6).hwPhi() > accumEtaNeg.at(5).hwPhi()) {
0829 accumEtaNeg.at(5) = accumEtaNeg.at(6);
0830 }
0831
0832
0833 accumEtaPos.resize(6);
0834 accumEtaNeg.resize(6);
0835
0836
0837 taus->clear();
0838 for (const l1t::Tau &accjet : accumEtaPos) {
0839 if (accjet.hwPt() > 0)
0840 taus->push_back(accjet);
0841 }
0842 for (const l1t::Tau &accjet : accumEtaNeg) {
0843 if (accjet.hwPt() > 0)
0844 taus->push_back(accjet);
0845 }
0846 }
0847
0848 void L1TStage2CaloLayer2Comp::dumpEventToFile() {
0849 edm::LogProblem("l1tcalol2ebec") << "==== Problems found, dumping full event contents ====" << std::endl;
0850
0851 edm::LogProblem("l1tcalol2ebec") << "==== Event contents in data: ====" << std::endl;
0852
0853 if (dumpTowers) {
0854 edm::LogProblem("l1tcalol2ebec") << "==== Towers: ====" << std::endl;
0855
0856 for (auto tower = caloTowerDataCol->begin(currBx); tower != caloTowerDataCol->end(currBx); ++tower)
0857 edm::LogProblem("l1tcalol2ebec") << "Tower: Et = " << tower->hwPt() << ", "
0858 << "eta = " << tower->hwEta() << ", "
0859 << "phi = " << tower->hwPhi() << std::endl;
0860 }
0861
0862 edm::LogProblem("l1tcalol2ebec") << "==== Jets: ====" << std::endl;
0863 for (auto jet = jetDataCol->begin(currBx); jet != jetDataCol->end(currBx); ++jet)
0864 edm::LogProblem("l1tcalol2ebec") << "Jet: Et = " << jet->hwPt() << ", "
0865 << "eta = " << jet->hwEta() << ", "
0866 << "phi = " << jet->hwPhi() << std::endl;
0867
0868 edm::LogProblem("l1tcalol2ebec") << "==== EGs: ====" << std::endl;
0869 for (auto eg = egDataCol->begin(currBx); eg != egDataCol->end(currBx); ++eg)
0870 edm::LogProblem("l1tcalol2ebec") << "EG: Et = " << eg->hwPt() << ", "
0871 << "eta = " << eg->hwEta() << ", "
0872 << "phi = " << eg->hwPhi() << std::endl;
0873
0874 edm::LogProblem("l1tcalol2ebec") << "==== Taus: ====" << std::endl;
0875 for (auto tau = tauDataCol->begin(currBx); tau != tauDataCol->end(currBx); ++tau)
0876 edm::LogProblem("l1tcalol2ebec") << "Tau: Et = " << tau->hwPt() << ", "
0877 << "eta = " << tau->hwEta() << ", "
0878 << "phi = " << tau->hwPhi() << std::endl;
0879
0880 edm::LogProblem("l1tcalol2ebec") << "==== Sums: ====" << std::endl;
0881 for (auto sum = sumDataCol->begin(currBx); sum != sumDataCol->end(currBx); ++sum)
0882 edm::LogProblem("l1tcalol2ebec") << "Sum: type = " << sum->getType() << " "
0883 << "Et = " << sum->hwPt() << ", "
0884 << "eta = " << sum->hwEta() << ", "
0885 << "phi = " << sum->hwPhi() << std::endl;
0886
0887 edm::LogProblem("l1tcalol2ebec") << "==== Event contents in emul: ====" << std::endl;
0888
0889 if (dumpTowers) {
0890 edm::LogProblem("l1tcalol2ebec") << "==== Towers: ====" << std::endl;
0891
0892 for (auto tower = caloTowerEmulCol->begin(currBx); tower != caloTowerEmulCol->end(currBx); ++tower)
0893 edm::LogProblem("l1tcalol2ebec") << "Tower: Et = " << tower->hwPt() << ", "
0894 << "eta = " << tower->hwEta() << ", "
0895 << "phi = " << tower->hwPhi() << std::endl;
0896 }
0897
0898 edm::LogProblem("l1tcalol2ebec") << "==== Jets: ====" << std::endl;
0899 for (auto jet = jetEmulCol->begin(currBx); jet != jetEmulCol->end(currBx); ++jet)
0900 edm::LogProblem("l1tcalol2ebec") << "Jet: Et = " << jet->hwPt() << ", "
0901 << "eta = " << jet->hwEta() << ", "
0902 << "phi = " << jet->hwPhi() << std::endl;
0903
0904 edm::LogProblem("l1tcalol2ebec") << "==== EGs: ====" << std::endl;
0905 for (auto eg = egEmulCol->begin(currBx); eg != egEmulCol->end(currBx); ++eg)
0906 edm::LogProblem("l1tcalol2ebec") << "EG: Et = " << eg->hwPt() << ", "
0907 << "eta = " << eg->hwEta() << ", "
0908 << "phi = " << eg->hwPhi() << std::endl;
0909
0910 edm::LogProblem("l1tcalol2ebec") << "==== Taus: ====" << std::endl;
0911 for (auto tau = tauEmulCol->begin(currBx); tau != tauEmulCol->end(currBx); ++tau)
0912 edm::LogProblem("l1tcalol2ebec") << "Tau: Et = " << tau->hwPt() << ", "
0913 << "eta = " << tau->hwEta() << ", "
0914 << "phi = " << tau->hwPhi() << std::endl;
0915
0916 edm::LogProblem("l1tcalol2ebec") << "==== Sums: ====" << std::endl;
0917 for (auto sum = sumEmulCol->begin(currBx); sum != sumEmulCol->end(currBx); ++sum)
0918 edm::LogProblem("l1tcalol2ebec") << "Sum: type = " << sum->getType() << " "
0919 << "Et = " << sum->hwPt() << ", "
0920 << "eta = " << sum->hwEta() << ", "
0921 << "phi = " << sum->hwPhi() << std::endl;
0922 }
0923
0924 void L1TStage2CaloLayer2Comp::dumpEventToEDM(edm::Event &e) {
0925
0926 std::unique_ptr<l1t::JetBxCollection> mpjets_data(new l1t::JetBxCollection(0, currBx, currBx));
0927 std::unique_ptr<l1t::JetBxCollection> mpjets_emul(new l1t::JetBxCollection(0, currBx, currBx));
0928
0929 for (auto jet = jetDataCol->begin(currBx); jet != jetDataCol->end(currBx); ++jet)
0930 mpjets_data->push_back(0, (*jet));
0931 for (auto jet = jetEmulCol->begin(currBx); jet != jetEmulCol->end(currBx); ++jet)
0932 mpjets_emul->push_back(0, (*jet));
0933
0934 e.put(std::move(mpjets_data), "dataJet");
0935 e.put(std::move(mpjets_emul), "emulJet");
0936
0937
0938 std::unique_ptr<l1t::EGammaBxCollection> mpEGammas_data(new l1t::EGammaBxCollection(0, currBx, currBx));
0939 std::unique_ptr<l1t::EGammaBxCollection> mpEGammas_emul(new l1t::EGammaBxCollection(0, currBx, currBx));
0940
0941 for (auto eg = egDataCol->begin(currBx); eg != egDataCol->end(currBx); ++eg)
0942 mpEGammas_data->push_back(0, (*eg));
0943 for (auto eg = egEmulCol->begin(currBx); eg != egEmulCol->end(currBx); ++eg)
0944 mpEGammas_emul->push_back(0, (*eg));
0945
0946 e.put(std::move(mpEGammas_data), "dataEg");
0947 e.put(std::move(mpEGammas_emul), "emulEg");
0948
0949
0950 std::unique_ptr<l1t::TauBxCollection> mptaus_data(new l1t::TauBxCollection(0, currBx, currBx));
0951 std::unique_ptr<l1t::TauBxCollection> mptaus_emul(new l1t::TauBxCollection(0, currBx, currBx));
0952
0953 for (auto tau = tauDataCol->begin(currBx); tau != tauDataCol->end(currBx); ++tau)
0954 mptaus_data->push_back(0, (*tau));
0955 for (auto tau = tauEmulCol->begin(currBx); tau != tauEmulCol->end(currBx); ++tau)
0956 mptaus_emul->push_back(0, (*tau));
0957
0958 e.put(std::move(mptaus_data), "dataTau");
0959 e.put(std::move(mptaus_emul), "emulTau");
0960
0961
0962 std::unique_ptr<l1t::EtSumBxCollection> mpsums_data(new l1t::EtSumBxCollection(0, currBx, currBx));
0963 std::unique_ptr<l1t::EtSumBxCollection> mpsums_emul(new l1t::EtSumBxCollection(0, currBx, currBx));
0964
0965 for (auto sum = sumDataCol->begin(currBx); sum != sumDataCol->end(currBx); ++sum)
0966 mpsums_data->push_back(0, (*sum));
0967 for (auto sum = sumEmulCol->begin(currBx); sum != sumEmulCol->end(currBx); ++sum)
0968 mpsums_emul->push_back(0, (*sum));
0969
0970 e.put(std::move(mpsums_data), "dataEtSum");
0971 e.put(std::move(mpsums_emul), "emulEtSum");
0972
0973
0974 std::unique_ptr<l1t::CaloTowerBxCollection> mptowers_data(new l1t::CaloTowerBxCollection(0, currBx, currBx));
0975 std::unique_ptr<l1t::CaloTowerBxCollection> mptowers_emul(new l1t::CaloTowerBxCollection(0, currBx, currBx));
0976
0977 for (auto tower = caloTowerDataCol->begin(currBx); tower != caloTowerDataCol->end(currBx); ++tower)
0978 mptowers_data->push_back(0, (*tower));
0979 for (auto tower = caloTowerEmulCol->begin(currBx); tower != caloTowerEmulCol->end(currBx); ++tower)
0980 mptowers_emul->push_back(0, (*tower));
0981
0982 e.put(std::move(mptowers_data), "dataCaloTower");
0983 e.put(std::move(mptowers_emul), "emulCaloTower");
0984 }
0985
0986 DEFINE_FWK_MODULE(L1TStage2CaloLayer2Comp);
0987
0988 #endif