File indexing completed on 2023-03-17 11:19:45
0001
0002 #include <functional>
0003 #include <numeric>
0004 #include <vector>
0005 #include <iostream>
0006 #include <sstream>
0007
0008
0009 #include "DataFormats/SiStripDigi/interface/SiStripDigi.h"
0010 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0011 #include "CalibFormats/SiStripObjects/interface/SiStripGain.h"
0012 #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
0013 #include "CalibFormats/SiStripObjects/interface/SiStripQuality.h"
0014 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
0015 #include "CalibTracker/Records/interface/SiStripGainRcd.h"
0016 #include "CalibTracker/Records/interface/SiStripQualityRcd.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/Frameworkfwd.h"
0019 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0020 #include "DataFormats/Common/interface/DetSetVector.h"
0021 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0022 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0023 #include "FWCore/Framework/interface/ESHandle.h"
0024 #include "FWCore/Utilities/interface/InputTag.h"
0025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0026
0027 class CompareClusters : public edm::one::EDAnalyzer<> {
0028 typedef edmNew::DetSetVector<SiStripCluster> input_t;
0029
0030 public:
0031 CompareClusters(const edm::ParameterSet& conf);
0032
0033 private:
0034 void analyze(const edm::Event&, const edm::EventSetup&);
0035
0036 void show(uint32_t);
0037 std::string printDigis(uint32_t);
0038 static std::string printCluster(const SiStripCluster&);
0039 static bool identicalDSV(const input_t&, const input_t&);
0040 static bool identicalDetSet(const edmNew::DetSet<SiStripCluster>&, const edmNew::DetSet<SiStripCluster>&);
0041 static bool identicalClusters(const SiStripCluster&, const SiStripCluster&);
0042
0043 const edm::ESGetToken<SiStripGain, SiStripGainRcd> gainESToken;
0044 const edm::ESGetToken<SiStripNoises, SiStripNoisesRcd> noiseESToken;
0045 const edm::ESGetToken<SiStripQuality, SiStripQualityRcd> qualESToken;
0046
0047 std::stringstream message;
0048
0049 const edm::InputTag clusters1, clusters2, digis;
0050 const edm::EDGetTokenT<input_t> clustersToken1, clustersToken2;
0051 const edm::EDGetTokenT<edm::DetSetVector<SiStripDigi>> digisToken;
0052 edm::Handle<input_t> clusterHandle1, clusterHandle2;
0053
0054 edm::Handle<edm::DetSetVector<SiStripDigi>> digiHandle;
0055 edm::ESHandle<SiStripGain> gainHandle;
0056 edm::ESHandle<SiStripNoises> noiseHandle;
0057 edm::ESHandle<SiStripQuality> qualityHandle;
0058 };
0059
0060 CompareClusters::CompareClusters(const edm::ParameterSet& conf)
0061 : gainESToken(esConsumes()),
0062 noiseESToken(esConsumes()),
0063 qualESToken(esConsumes()),
0064 clusters1(conf.getParameter<edm::InputTag>("Clusters1")),
0065 clusters2(conf.getParameter<edm::InputTag>("Clusters2")),
0066 digis(conf.getParameter<edm::InputTag>("Digis")),
0067 clustersToken1(consumes<input_t>(clusters1)),
0068 clustersToken2(consumes<input_t>(clusters2)),
0069 digisToken(consumes<edm::DetSetVector<SiStripDigi>>(digis)) {}
0070
0071 void CompareClusters::analyze(const edm::Event& event, const edm::EventSetup& es) {
0072 event.getByToken(clustersToken1, clusterHandle1);
0073 if (!clusterHandle1.isValid())
0074 throw cms::Exception("Input Not found") << clusters1;
0075 event.getByToken(clustersToken2, clusterHandle2);
0076 if (!clusterHandle2.isValid())
0077 throw cms::Exception("Input Not found") << clusters2;
0078 if (identicalDSV(*clusterHandle1, *clusterHandle2))
0079 return;
0080
0081 {
0082 event.getByToken(digisToken, digiHandle);
0083 if (!digiHandle.isValid())
0084 throw cms::Exception("Input Not found") << digis;
0085 gainHandle = es.getHandle(gainESToken);
0086 noiseHandle = es.getHandle(noiseESToken);
0087 qualityHandle = es.getHandle(qualESToken);
0088 }
0089
0090 input_t::const_iterator set1(clusterHandle1->begin()), end1(clusterHandle1->end()), set2(clusterHandle2->begin()),
0091 end2(clusterHandle2->end());
0092
0093 message.str("");
0094 while (set1 != end1 && set2 != end2) {
0095 while (set1 != end1 && set1->id() < set2->id())
0096 show((set1++)->id());
0097 while (set2 != end2 && set2->id() < set1->id())
0098 show((set2++)->id());
0099 if (set1 != end1 && set2 != end2) {
0100 if (!identicalDetSet(*set1, *set2))
0101 show(set1->id());
0102 set1++;
0103 set2++;
0104 }
0105 }
0106 while (set1 != end1)
0107 show((set1++)->id());
0108 while (set2 != end2)
0109 show((set2++)->id());
0110 edm::LogError("Not Identical") << message.str();
0111
0112 return;
0113 }
0114
0115 void CompareClusters::show(uint32_t id) {
0116 message << std::endl << "detId: " << id << std::endl;
0117 message << "Digis:\n" << printDigis(id);
0118 edmNew::DetSet<SiStripCluster>::const_iterator c1(0), c2(0), end1(0), end2(0);
0119 if (clusterHandle1->find(id) != clusterHandle1->end()) {
0120 c1 = clusterHandle1->find(id)->begin();
0121 end1 = clusterHandle1->find(id)->end();
0122 }
0123 if (clusterHandle2->find(id) != clusterHandle2->end()) {
0124 c2 = clusterHandle2->find(id)->begin();
0125 end2 = clusterHandle2->find(id)->end();
0126 }
0127 message << clusters1.label() << std::endl;
0128 while (c1 != end1)
0129 message << printCluster(*c1++);
0130 message << clusters2.label() << std::endl;
0131 while (c2 != end2)
0132 message << printCluster(*c2++);
0133 }
0134
0135 std::string CompareClusters::printCluster(const SiStripCluster& cluster) {
0136 std::stringstream s;
0137 s << "\t" << cluster.firstStrip() << " [ ";
0138 for (unsigned i = 0; i < cluster.amplitudes().size(); i++) {
0139 s << static_cast<int>(cluster.amplitudes()[i]) << " ";
0140 }
0141 s << "]" << std::endl;
0142 return s.str();
0143 }
0144
0145 std::string CompareClusters::printDigis(uint32_t id) {
0146 std::stringstream s;
0147 SiStripApvGain::Range gainRange = gainHandle->getRange(id);
0148 SiStripNoises::Range noiseRange = noiseHandle->getRange(id);
0149 SiStripQuality::Range qualityRange = qualityHandle->getRange(id);
0150 edm::DetSetVector<SiStripDigi>::const_iterator set = digiHandle->find(id);
0151 if (set != digiHandle->end()) {
0152 for (edm::DetSet<SiStripDigi>::const_iterator digi = set->begin(); digi != set->end(); digi++) {
0153 s << "( " << digi->strip() << ", " << digi->adc() << ", " << noiseHandle->getNoise(digi->strip(), noiseRange)
0154 << ", " << gainHandle->getStripGain(digi->strip(), gainRange) << ", "
0155 << (qualityHandle->IsStripBad(qualityRange, digi->strip()) ? "bad" : "good") << ")\n";
0156 }
0157 }
0158 return s.str();
0159 }
0160
0161 bool CompareClusters::identicalDSV(const input_t& L, const input_t& R) {
0162 return L.size() == R.size() &&
0163 inner_product(L.begin(), L.end(), R.begin(), bool(true), std::logical_and<bool>(), identicalDetSet);
0164 }
0165
0166 bool CompareClusters::identicalDetSet(const edmNew::DetSet<SiStripCluster>& L,
0167 const edmNew::DetSet<SiStripCluster>& R) {
0168 return L.id() == R.id() && L.size() == R.size() &&
0169 inner_product(L.begin(), L.end(), R.begin(), bool(true), std::logical_and<bool>(), identicalClusters);
0170 }
0171
0172 bool CompareClusters::identicalClusters(const SiStripCluster& L, const SiStripCluster& R) {
0173 return L.firstStrip() == R.firstStrip() && L.amplitudes().size() == R.amplitudes().size() &&
0174 inner_product(L.amplitudes().begin(),
0175 L.amplitudes().end(),
0176 R.amplitudes().begin(),
0177 bool(true),
0178 std::logical_and<bool>(),
0179 std::equal_to<uint16_t>());
0180 }
0181
0182 #include "FWCore/Framework/interface/MakerMacros.h"
0183 DEFINE_FWK_MODULE(CompareClusters);