File indexing completed on 2024-04-06 12:25:45
0001 #include "RecoLocalCalo/EcalRecProducers/test/EcalCompactTrigPrimProducerTest.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005
0006 #include <iostream>
0007 #include <iomanip>
0008 #include <sstream>
0009
0010 void EcalCompactTrigPrimProducerTest::analyze(edm::Event const& event, edm::EventSetup const& c) {
0011 edm::Handle<EcalTrigPrimDigiCollection> hTpDigis;
0012 event.getByToken(tpDigiToken_, hTpDigis);
0013
0014 const EcalTrigPrimDigiCollection* trigPrims = hTpDigis.product();
0015
0016 edm::Handle<EcalTrigPrimCompactColl> hTpRecs;
0017 event.getByToken(tpRecToken_, hTpRecs);
0018
0019 const EcalTrigPrimCompactColl* trigPrimRecs = hTpRecs.product();
0020
0021 int nTps = 0;
0022 err_ = false;
0023 for (EcalTrigPrimDigiCollection::const_iterator trigPrim = trigPrims->begin(); trigPrim != trigPrims->end();
0024 ++trigPrim) {
0025 const EcalTrigTowerDetId& ttId = trigPrim->id();
0026
0027 if ((trigPrim->sample(trigPrim->sampleOfInterest()).raw() & 0x1FFF) != (trigPrimRecs->raw(ttId) & 0x1FFF)) {
0028 err("Different TP (0x") << std::hex << trigPrim->sample(trigPrim->sampleOfInterest()).raw() << " -- "
0029 << trigPrimRecs->raw(ttId) << std::dec << ") for " << ttId << "\n";
0030 }
0031 if (trigPrim->compressedEt() != trigPrimRecs->compressedEt(ttId))
0032 err("\tDifferent compressed Et\n");
0033 if (trigPrim->fineGrain() != trigPrimRecs->fineGrain(ttId))
0034 err("\tDifferent FGVB\n") << ttId << "\n";
0035 if (trigPrim->ttFlag() != trigPrimRecs->ttFlag(ttId))
0036 err("\tDifferent compressed TTF\n") << ttId << "\n";
0037 if (trigPrim->l1aSpike() != trigPrimRecs->l1aSpike(ttId))
0038 err("\tDifferent compressed L1Spike flag\n");
0039 if (trigPrim->compressedEt() != 0)
0040 ++nCompressEt_;
0041 if (trigPrim->fineGrain() != 0)
0042 ++nFineGrain_;
0043 if (trigPrim->ttFlag() != 0)
0044 ++nTTF_;
0045 if (trigPrim->l1aSpike() != 0)
0046 ++nL1aSpike_;
0047 ++nTps;
0048 }
0049 if (nTps != 4032)
0050 err("Unexpected number of TPs: ") << nTps << "\n";
0051
0052 if (!err_)
0053 std::cout << "Validation of compact trigger primitive collection " << (err_ ? "failed" : "succeeded") << "\n";
0054
0055 if (err_) {
0056 std::cout << "Cannot check compact to legacy collection convertion because of previous failure\n";
0057 } else {
0058
0059 EcalTrigPrimDigiCollection col2;
0060 trigPrimRecs->toEcalTrigPrimDigiCollection(col2);
0061 if (col2.size() != trigPrims->size()) {
0062 err("Collection size error!\n");
0063 err_ = true;
0064 } else {
0065 EcalTrigPrimDigiCollection::const_iterator trigPrim2 = col2.begin();
0066 for (EcalTrigPrimDigiCollection::const_iterator trigPrim = trigPrims->begin();
0067 trigPrim != trigPrims->end() && !err_;
0068 ++trigPrim, ++trigPrim2) {
0069 if ((trigPrim->sample(trigPrim->sampleOfInterest()).raw() &
0070 0x1FFF)
0071 != (trigPrim2->sample(0).raw() & 0x1FFF))
0072 err("Trig prim differs: ") << *trigPrim << " (" << std::hex
0073 << (trigPrim->sample(trigPrim->sampleOfInterest()).raw()) << std::dec << ") "
0074 << " -- " << *trigPrim2 << " (" << std::hex
0075 << (trigPrim2->sample(trigPrim2->sampleOfInterest()).raw()) << std::dec << ") "
0076 << "\n";
0077 }
0078 }
0079 std::cout << "Validation of compact-to-legacy trigger primitive collection conversion "
0080 << (err_ ? "failed" : "succeeded") << "\n";
0081 }
0082
0083
0084 edm::Handle<EcalTrigPrimDigiCollection> hSkimTpDigis;
0085 event.getByToken(tpSkimToken_, hSkimTpDigis);
0086 const EcalTrigPrimDigiCollection* skimTrigPrims = hSkimTpDigis.product();
0087 err_ = false;
0088 if (skimTrigPrims->size() == 0)
0089 err("Skimmed TP collection is empty!");
0090 std::stringstream tpList;
0091 for (EcalTrigPrimDigiCollection::const_iterator skimTrigPrim = skimTrigPrims->begin();
0092 skimTrigPrim != skimTrigPrims->end();
0093 ++skimTrigPrim) {
0094 tpList << "\t- detid=" << skimTrigPrim->id().rawId() << " ieta=" << skimTrigPrim->id().ieta()
0095 << " iphi=" << skimTrigPrim->id().iphi() << "\n";
0096 EcalTrigPrimDigiCollection::const_iterator origTrigPrim = trigPrims->find(skimTrigPrim->id());
0097 if (origTrigPrim == trigPrims->end())
0098 err("Skimmed TP ") << skimTrigPrim->id() << " not found in original TP collection!\n";
0099 else {
0100 if (skimTrigPrim->size() != origTrigPrim->size()) {
0101 std::cout << "TP from skimmed colletion has " << skimTrigPrim->size() << " sample(s), "
0102 << "while TP from original collection has " << origTrigPrim->size() << " sample(s)!\n";
0103 }
0104 bool oneSample = false;
0105 if (skimTrigPrim->size() == 1 || origTrigPrim->size() == 1) {
0106 std::cout << "\tComparing only the \"sample of interest\"!\n";
0107 oneSample = true;
0108 }
0109 bool eq = true;
0110 if (oneSample) {
0111 const int skimSample = skimTrigPrim->sampleOfInterest();
0112 const int origSample = origTrigPrim->sampleOfInterest();
0113
0114 eq = ((skimTrigPrim->sample(skimSample).raw() & 0x1FFF) == (origTrigPrim->sample(origSample).raw() & 0x1FFF));
0115 } else if (skimTrigPrim->size() == origTrigPrim->size()) {
0116 for (int iS = 0; iS < skimTrigPrim->size(); ++iS) {
0117
0118 if ((skimTrigPrim->sample(iS).raw() & 0x1FFF) != (origTrigPrim->sample(iS).raw() & 0x1FFF))
0119 eq = false;
0120 }
0121 } else {
0122 err_ = true;
0123 }
0124 if (!eq)
0125 err("Skimmed Trig prim differs from original one: ") << *skimTrigPrim << " -- " << *origTrigPrim << "\n";
0126 }
0127 }
0128 std::cout << "List if TT in skimmed TP collection: " << tpList.str() << "\n";
0129 std::cout << "Validation of skimmed trigger primitive collection " << (err_ ? "failed" : "succeeded") << "\n";
0130 }
0131
0132 std::ostream& EcalCompactTrigPrimProducerTest::err(const char* mess) {
0133 err_ = true;
0134 std::cout << mess;
0135 return std::cout;
0136 }
0137
0138 EcalCompactTrigPrimProducerTest::~EcalCompactTrigPrimProducerTest() {
0139 std::cout << "# of non-null compressed Et: " << nCompressEt_ << "\n";
0140 std::cout << "# of non-null FGVB: " << nFineGrain_ << "\n";
0141 std::cout << "# of non-null TTF: " << nTTF_ << "\n";
0142 std::cout << "# of non-null L1ASpike: " << nL1aSpike_ << "\n";
0143 }
0144
0145 #include "FWCore/Framework/interface/MakerMacros.h"
0146 DEFINE_FWK_MODULE(EcalCompactTrigPrimProducerTest);