File indexing completed on 2023-03-17 11:25:44
0001
0002 #include <memory>
0003
0004
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/one/EDProducer.h"
0007
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012
0013 #include "DataFormats/Common/interface/DetSetVector.h"
0014 #include "SimDataFormats/TrackerDigiSimLink/interface/StripCompactDigiSimLinks.h"
0015 #include "SimDataFormats/TrackerDigiSimLink/interface/StripDigiSimLink.h"
0016
0017 #ifdef SCDSL_DEBUG
0018 #define DEBUG(X) X
0019 #else
0020 #define DEBUG(X)
0021 #endif
0022
0023 class StripCompactDigiSimLinksProducer : public edm::one::EDProducer<> {
0024 public:
0025 StripCompactDigiSimLinksProducer(const edm::ParameterSet &iConfig);
0026 ~StripCompactDigiSimLinksProducer() override;
0027
0028 void produce(edm::Event &, const edm::EventSetup &) override;
0029
0030 private:
0031 edm::InputTag src_;
0032 uint32_t maxHoleSize_;
0033 };
0034
0035 StripCompactDigiSimLinksProducer::StripCompactDigiSimLinksProducer(const edm::ParameterSet &iConfig)
0036 : src_(iConfig.getParameter<edm::InputTag>("src")), maxHoleSize_(iConfig.getParameter<uint32_t>("maxHoleSize")) {
0037 produces<StripCompactDigiSimLinks>();
0038 }
0039
0040 StripCompactDigiSimLinksProducer::~StripCompactDigiSimLinksProducer() {}
0041
0042 void StripCompactDigiSimLinksProducer::produce(edm::Event &iEvent, const edm::EventSetup &) {
0043 using namespace edm;
0044 Handle<DetSetVector<StripDigiSimLink>> src;
0045 iEvent.getByLabel(src_, src);
0046
0047 StripCompactDigiSimLinks::Filler output;
0048
0049 int previousStrip;
0050
0051
0052 int previousLinkStrip;
0053
0054 std::vector<StripCompactDigiSimLinks::key_type> thisStripSignals;
0055 std::vector<StripCompactDigiSimLinks::key_type> previousStripSignals;
0056
0057 for (auto const &det : *src) {
0058 DEBUG(std::cerr << "\n\nProcessing detset " << det.detId() << ", size = " << det.size() << std::endl;)
0059 previousStrip = -2;
0060
0061 previousLinkStrip = -2;
0062
0063 thisStripSignals.clear();
0064 previousStripSignals.clear();
0065 for (DetSet<StripDigiSimLink>::const_iterator it = det.begin(), ed = det.end(); it != ed; ++it) {
0066 DEBUG(std::cerr << " processing digiSimLink on strip " << it->channel() << " left by particle "
0067 << it->SimTrackId() << ", event " << it->eventId().rawId() << std::endl;)
0068 if (int(it->channel()) != previousLinkStrip) {
0069 previousStrip = previousLinkStrip;
0070 DEBUG(std::cerr << " strip changed!" << std::endl;)
0071 swap(thisStripSignals, previousStripSignals);
0072 thisStripSignals.clear();
0073 }
0074 DEBUG(std::cerr << " previous strip " << previousStrip << ", previous link strip " << previousLinkStrip
0075 << std::endl;)
0076
0077
0078
0079
0080
0081
0082 StripCompactDigiSimLinks::key_type key(it->eventId(), it->SimTrackId());
0083 bool alreadyClusterized = false;
0084 if (int(it->channel()) == previousStrip + 1) {
0085 DEBUG(std::cerr << " on next strip" << std::endl;)
0086 if (std::find(previousStripSignals.begin(), previousStripSignals.end(), key) != previousStripSignals.end()) {
0087 alreadyClusterized = true;
0088 DEBUG(std::cerr << " and part of previous cluster" << std::endl;)
0089 }
0090 }
0091 if (!alreadyClusterized) {
0092 DEBUG(std::cerr << " clusterize!" << std::endl;)
0093 unsigned int size = 1;
0094 int myLastStrip = it->channel();
0095 for (DetSet<StripDigiSimLink>::const_iterator it2 = it + 1; it2 < ed; ++it2) {
0096 DEBUG(std::cerr << " digiSimLink on strip " << it2->channel() << " left by particle " << it2->SimTrackId()
0097 << ", event " << it2->eventId().rawId() << std::endl;)
0098 if ((it2->channel() - myLastStrip) > maxHoleSize_ + 1) {
0099 DEBUG(std::cerr << " found hole of size " << (it2->channel() - myLastStrip) << ", stopping."
0100 << std::endl;)
0101 break;
0102 }
0103 if ((it2->eventId() == key.first) && (it2->SimTrackId() == key.second)) {
0104 size++;
0105 DEBUG(std::cerr << " extending cluster, now size = " << size << std::endl;)
0106 myLastStrip = it2->channel();
0107 }
0108 }
0109 output.insert(key, StripCompactDigiSimLinks::HitRecord(det.detId(), it->channel(), size));
0110 }
0111 if (int(it->channel()) != previousLinkStrip) {
0112 previousLinkStrip = it->channel();
0113 }
0114 thisStripSignals.push_back(key);
0115 DEBUG(std::cerr << " ending state " << previousStrip << ", previous link strip " << previousLinkStrip
0116 << std::endl;)
0117
0118
0119
0120
0121
0122
0123 DEBUG(std::cerr << std::endl;)
0124 }
0125 }
0126
0127 std::unique_ptr<StripCompactDigiSimLinks> ptr(new StripCompactDigiSimLinks(output));
0128 iEvent.put(std::move(ptr));
0129 }
0130
0131 DEFINE_FWK_MODULE(StripCompactDigiSimLinksProducer);